DUNE-ACFEM (2.5.1)

Dune::ACFem::ModelInterface< ModelType > Class Template Reference

Interface class for second order elliptic models. More...

#include <dune/acfem/models/modelinterface.hh>

Classes

struct  ForcesFunctionalTraits
 Helper structure to define an optional functional as part of the "right hand side". More...
 

Public Types

Right-hand-side functional definitions

Bulk and boundary contributions to the "Right-Hand-Side".

The ModelInterface allows for "Right-Hand-Side" functionals of the form

\[ \langle\rho\chi\sigma, \phi\rangle = \int_\Omega f\,\phi + \int_{\partial\Omega} g\,\phi + \langle\Psi,\phi\rangle \]

\(f\) is modelled by the ModelInterface::BulkForcesFunctionType, \(g\) is modelled by the ModelInterface::NeumannBoundaryFunctionType and the \(H^1\)-functional \(\Psi\) is modelled by ForcesFunctionalTraits::FunctionalType.

The default is to keep everything at zero.

Note
Functionals cannot not be operated on in a natural way by \(L^\infty\)-functions, only a vector space structure is "natural". Therefore ModelExpressions templates which model multiplication by functions will yield undefined results when applied to PDE-Models with a LinearFunctional as "right hand side" (maybe rejected by the compiler, maybe simply forming something which was most likely not intended).
typedef TraitsType::BulkForcesFunctionType BulkForcesFunctionType
 A function modelling "force" terms in the bulk-phase. More...
 
typedef TraitsType::NeumannBoundaryFunctionType NeumannBoundaryFunctionType
 A function modelling "force" terms in the bulk-phase. More...
 
typedef TraitsType::NeumannIndicatorType NeumannIndicatorType
 A function modelling "force" terms in the bulk-phase. More...
 
Dirichlet boundary definitions
typedef TraitsType::DirichletBoundaryFunctionType DirichletBoundaryFunctionType
 A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType. More...
 
typedef TraitsType::DirichletWeightFunctionType DirichletWeightFunctionType
 A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType. More...
 
typedef TraitsType::DirichletIndicatorType DirichletIndicatorType
 A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType. More...
 

Public Member Functions

std::string name () const
 Print a descriptive name for debugging and output.
 
OperatorPartsType operatorParts () const
 Return the integral kernels for the bilinear form.
 
BulkForcesFunctionType bulkForcesFunction (const GridPartType &gridPart) const
 Generate an instance of a class defining the bulk-forces the model is subject to. More...
 
template<class DiscreteFunctionSpace >
TraitsType::template ForcesFunctionalTraits< DiscreteFunctionSpace >::FunctionalType forcesFunctional (const DiscreteFunctionSpace &space) const
 Generate an instance of a class defining a functional which forms part of the force-terms for the model.
 
DirichletBoundaryFunctionType dirichletBoundaryFunction (const GridPartType &gridPart) const
 Generate an instance of a class defining Dirichlet boundary values as a Fem grid-function.
 
DirichletWeightFunctionType dirichletWeightFunction (const GridPartType &gridPart) const
 Generate an instance of a class defining a "left hand side" weight on the Dirichlet boundary conditions, DefaultModelTraits::DirichletWeightFunctionType.
 
NeumannBoundaryFunctionType neumannBoundaryFunction (const GridPartType &gridPart) const
 Generate an instance of a class defining Neumann boundary values as a Fem grid-function.
 
DirichletIndicatorType dirichletIndicator () const
 Generate an object to identify parts of the boundary subject to Dirichlet boundary conditions. More...
 
NeumannIndicatorType neumannIndicator () const
 Generate an object to identify parts of the boundary subject to Neumann boundary conditions. More...
 

Detailed Description

template<class ModelType>
class Dune::ACFem::ModelInterface< ModelType >

Interface class for second order elliptic models.

The ModelInterface defines all the different parts needed to compute the weak form of a second order elliptic operator of the form

\[ \begin{split} -\nabla\cdot (A(x, u, \nabla u)\,\nabla u) + \nabla\cdot(b(x, u)\,u) + c(x, u)\,u &= f(x) + \Psi\quad \text{ in }\Omega,\\ u &= g_D \text{ on }\Gamma_D,\\ (A(x, u, \nabla u)\nabla u)\cdot\nu + \alpha(x, u)\,u &= g_N \text{ on }\Gamma_R,\\ (A(x, u, \nabla u)\nabla u)\cdot\nu &= g_N \text{ on }\Gamma_N,\\ \end{split} \]

In particular, \(\Psi\in H^{-1}\) models a functional, \(f\in L^2\) is the usual "right hand side".

The methods in this class provide all factors of the integrands not involving test functions. The multiplication by the test-functions is added in the class EllipticOperator. The weak formulation then reads:

\[ \int_\Omega (A(x,u,\nabla u)\,\nabla u)\cdot\nabla\phi\,dx + \int_\Omega (\nabla\cdot(b(x, u)\,u) + c(x, u)\,u)\,\phi\,dx + \int_{\Gamma_R} \alpha(x, u)\,u\,\phi\,do - \int_{\Gamma_N} g_N\,\phi\,do - \int_\Omega f(x)\,\phi\,dx - \langle\Psi,\,\phi\rangle = 0. \]

As the ModelInterface is potentially non-linear one has the option to leave the ModelInterface::BulkForcesFunctionType and ModelInterface::NeumannBoundaryFunctionType at zero and instead stuff the "right-hand-side" functions into the ModelInterface::robinFlux() and ModelInterface::source() methods. In principle even the functional \(\Psi\) could be expressed by a \(L^2\) scalar-product (don't shout at me: as the FEM-space is finite dimensional this is of course possible and one way to implement such a functional).

Nevertheless the different forces and other right-hand-side components are there and will be used by EllipticFemScheme and ParabolicFemScheme. Generally, the implementation checks for zero-objects (see ZeroExpression, ZeroGridFunction, ZeroFunctional, ZeroModel); the decision about which component has to be taken into account is taken at compile-time.

Parameters
[in]ModelTypeThe implementation.
Bug:
The setEntity() method is evil, because it in principle prevents models from being tagged as "const", or is the reason for the evil use of the "mutable" keyword in several places. Way out would be an "OperatorGerm"-class which is constructed from a given entity and/or intersection.
Todo:
In order to support systems there should be a DomainFunctionSpace (ansatz-functions) and a RangeFunctionSpace (test-functions). Example would be velocity-pressure for mixed discretizations of Stokes-like problems. Although it would be possible to stuff such a problem into the existing model-problem this does not feel like being desirable.

Member Typedef Documentation

◆ BulkForcesFunctionType

template<class ModelType >
typedef TraitsType::BulkForcesFunctionType Dune::ACFem::ModelInterface< ModelType >::BulkForcesFunctionType

A function modelling "force" terms in the bulk-phase.

◆ DirichletBoundaryFunctionType

template<class ModelType >
typedef TraitsType::DirichletBoundaryFunctionType Dune::ACFem::ModelInterface< ModelType >::DirichletBoundaryFunctionType

A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType.

See also
BoundarySupportedFunction, in particular the "factory" function boundarySupportedFunction() and BoundaryFunctionExpressions for the set of supported algebraic operations for such functions.

◆ DirichletIndicatorType

template<class ModelType >
typedef TraitsType::DirichletIndicatorType Dune::ACFem::ModelInterface< ModelType >::DirichletIndicatorType

A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType.

See also
BoundarySupportedFunction, in particular the "factory" function boundarySupportedFunction() and BoundaryFunctionExpressions for the set of supported algebraic operations for such functions.

◆ DirichletWeightFunctionType

template<class ModelType >
typedef TraitsType::DirichletWeightFunctionType Dune::ACFem::ModelInterface< ModelType >::DirichletWeightFunctionType

A BoundarySupportedFunction which must be sub-ordinate to the DirichletIndicatorType.

See also
BoundarySupportedFunction, in particular the "factory" function boundarySupportedFunction() and BoundaryFunctionExpressions for the set of supported algebraic operations for such functions.

◆ NeumannBoundaryFunctionType

template<class ModelType >
typedef TraitsType::NeumannBoundaryFunctionType Dune::ACFem::ModelInterface< ModelType >::NeumannBoundaryFunctionType

A function modelling "force" terms in the bulk-phase.

◆ NeumannIndicatorType

template<class ModelType >
typedef TraitsType::NeumannIndicatorType Dune::ACFem::ModelInterface< ModelType >::NeumannIndicatorType

A function modelling "force" terms in the bulk-phase.

Member Function Documentation

◆ bulkForcesFunction()

template<class ModelType >
BulkForcesFunctionType Dune::ACFem::ModelInterface< ModelType >::bulkForcesFunction ( const GridPartType &  gridPart) const
inline

Generate an instance of a class defining the bulk-forces the model is subject to.

The return object is a Fem grid-function. This can be used for the case where the "right hand side" is an L2-function.

References Dune::ACFem::asImp().

Referenced by gridWalkTest().

◆ dirichletIndicator()

template<class ModelType >
DirichletIndicatorType Dune::ACFem::ModelInterface< ModelType >::dirichletIndicator ( ) const
inline

Generate an object to identify parts of the boundary subject to Dirichlet boundary conditions.

The return value has to obey the BoundaryIndicatorInterface.

References Dune::ACFem::asImp().

Referenced by gridWalkTest().

◆ neumannIndicator()

template<class ModelType >
NeumannIndicatorType Dune::ACFem::ModelInterface< ModelType >::neumannIndicator ( ) const
inline

Generate an object to identify parts of the boundary subject to Neumann boundary conditions.

The return value has to obey the BoundaryIndicatorInterface.

References Dune::ACFem::asImp().

Referenced by gridWalkTest(), and Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, QuadratureTraits >::solve().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)