KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::ScalarLaplace< Scalar, dim, isotropic > Class Template Reference

Convenience class for handling diffusion terms in elliptic/parabolic equations. More...

#include <laplace.hh>

Detailed Description

template<class Scalar, int dim, bool isotropic = true>
class Kaskade::ScalarLaplace< Scalar, dim, isotropic >

Convenience class for handling diffusion terms in elliptic/parabolic equations.

This class defines the variational formulation of the diffusion operator \( -\mathrm{div}(\kappa \nabla \cdot)\) for scalar-valued functions, i.e. the variational functional \( u \mapsto \nabla \frac{1}{2}u^T \kappa \nabla u \).

Definition at line 33 of file laplace.hh.

Public Types

using DiffusionTensor = std::conditional_t< isotropic, Scalar, Dune::FieldMatrix< Scalar, dim, dim > >
 The type of the diffusion tensor. More...
 

Public Member Functions

 ScalarLaplace ()
 Constructor The diffusion constant \( kappa \) is initialized to 1, the linearization point \(u_0\) to \( \nabla u_0 = 0 \). More...
 
LaplacesetDiffusionTensor (DiffusionTensor const &kappa_)
 Specifies a scalar diffusion constant. More...
 
LaplacesetLinearizationPoint (Dune::FieldVector< Scalar, dim > const &du0_)
 Specifies the gradient of the linearization point of the Laplace variational functional. More...
 
Scalar d0 () const
 Returns the value of the Laplace variational functional. The value returned is. More...
 
Dune::FieldVector< Scalar, 1 > d1 (VariationalArg< Scalar, dim > const &arg) const
 
Dune::FieldMatrix< Scalar, 1, 1 > d2 (VariationalArg< Scalar, dim > const &arg1, VariationalArg< Scalar, dim > const &arg2) const
 

Member Typedef Documentation

◆ DiffusionTensor

template<class Scalar , int dim, bool isotropic = true>
using Kaskade::ScalarLaplace< Scalar, dim, isotropic >::DiffusionTensor = std::conditional_t<isotropic, Scalar, Dune::FieldMatrix<Scalar,dim,dim> >

The type of the diffusion tensor.

Conceptually, this is always a dim x dim matrix. For isotropic problems, it is just a multiple of the unit matrix. In this case, we use a scalar diffusion constant.

Definition at line 42 of file laplace.hh.

Constructor & Destructor Documentation

◆ ScalarLaplace()

template<class Scalar , int dim, bool isotropic = true>
Kaskade::ScalarLaplace< Scalar, dim, isotropic >::ScalarLaplace ( )
inline

Constructor The diffusion constant \( kappa \) is initialized to 1, the linearization point \(u_0\) to \( \nabla u_0 = 0 \).

Definition at line 49 of file laplace.hh.

Member Function Documentation

◆ d0()

template<class Scalar , int dim, bool isotropic = true>
Scalar Kaskade::ScalarLaplace< Scalar, dim, isotropic >::d0 ( ) const
inline

Returns the value of the Laplace variational functional. The value returned is.

\[ \frac{\kappa}{2} |\nabla u_0|^2, \]

where the gradient \( \nabla u_0 \) of the linearization point \( u_0 \) is set via setLinearizationPoint.

Definition at line 81 of file laplace.hh.

◆ d1()

template<class Scalar , int dim, bool isotropic = true>
Dune::FieldVector< Scalar, 1 > Kaskade::ScalarLaplace< Scalar, dim, isotropic >::d1 ( VariationalArg< Scalar, dim > const &  arg) const
inline

Definition at line 87 of file laplace.hh.

◆ d2()

template<class Scalar , int dim, bool isotropic = true>
Dune::FieldMatrix< Scalar, 1, 1 > Kaskade::ScalarLaplace< Scalar, dim, isotropic >::d2 ( VariationalArg< Scalar, dim > const &  arg1,
VariationalArg< Scalar, dim > const &  arg2 
) const
inline

Definition at line 93 of file laplace.hh.

◆ setDiffusionTensor()

template<class Scalar , int dim, bool isotropic = true>
Laplace & Kaskade::ScalarLaplace< Scalar, dim, isotropic >::setDiffusionTensor ( DiffusionTensor const &  kappa_)
inline

Specifies a scalar diffusion constant.

Parameters
kappathe diffusion constant (has to be symmetric positive definite)

Definition at line 61 of file laplace.hh.

◆ setLinearizationPoint()

template<class Scalar , int dim, bool isotropic = true>
Laplace & Kaskade::ScalarLaplace< Scalar, dim, isotropic >::setLinearizationPoint ( Dune::FieldVector< Scalar, dim > const &  du0_)
inline

Specifies the gradient of the linearization point of the Laplace variational functional.

Definition at line 70 of file laplace.hh.


The documentation for this class was generated from the following file: