DUNE-ACFEM (2.5.1)

Dune::ACFem::OperatorPartsInterface< PartsImpl > Class Template Reference

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

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

Public Types

enum  StructureFlags { isLinear = TraitsType::isLinear , isSymmetric = TraitsType::isSymmetric , isSemiDefinite = TraitsType::isSemiDefinite }
 Static flags for the overall structure of the operator. More...
 
enum  ConstituentFlags
 Provide information about the constituents of the model. More...
 

Public Member Functions

template<class Entity >
void setEntity (const Entity &entity) const
 Per entity initialization, if that is needed. More...
 
template<class Intersection >
bool setIntersection (const Intersection &intersection) const
 Per-intersection initialization for the boundary contributions. More...
 
std::string name () const
 Print a descriptive name for debugging and output.
 
template<class Entity , class Point >
void flux (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More...
 
template<class Entity , class Point >
void linearizedFlux (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 Evaluate the linearized flux in local coordinates. More...
 
template<class Entity , class Point >
void source (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 The zero-order term as function of local coordinates. More...
 
template<class Entity , class Point >
void linearizedSource (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 The linearized source term as function of local coordinates. More...
 
template<class Intersection , class Point >
void robinFlux (const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 The non-linearized Robin-type flux term. More...
 
template<class Intersection , class Point >
void linearizedRobinFlux (const RangeType &uBar, const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 The linearized Robin-type flux term. More...
 
template<class Entity , class Point >
void fluxDivergence (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, const HessianRangeType &hessian, RangeType &result) const
 Compute the point-wise value of the flux-part of the operator, meaning the part of the differential operator which is multiplied by the derivative of the test function. More...
 

Detailed Description

template<class PartsImpl>
class Dune::ACFem::OperatorPartsInterface< PartsImpl >

Interface class for second order elliptic models.

\[ \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]PartsImplThe 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 Enumeration Documentation

◆ ConstituentFlags

Provide information about the constituents of the model.

All other components (like forces, boundary values) can be identified by looking at the data type.

◆ StructureFlags

template<class PartsImpl >
enum Dune::ACFem::OperatorPartsInterface::StructureFlags

Static flags for the overall structure of the operator.

Enumerator
isLinear 

Define to true for the affine-linear case.

isSymmetric 

Define to true for the symmetric case.

isSemiDefinite 

Define to true for the non-indefinite case (non trivial kernel is allowed).

Member Function Documentation

◆ flux()

template<class PartsImpl >
template<class Entity , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::flux ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
JacobianRangeType &  flux 
) const
inline

Evaluate \(A(x, u)\nabla u(x)\) in local coordinates.

This can be interpreted as a diffusive "flux", at least when restricted to some surface. \(\bar u\) denotes the point of linearization for non-linear problems.

Parameters
[in]entityThe mesh element.
[in]xThe point of evaluation, local coordinates.
[in]valueTo allow integration by parts for first order terms and in preparation for non-linear models.
[in]jacobianThe value of \(\nabla\bar\psi\).
[out]fluxThe result of the computation.
Note
"flux" in this sense simply means something which is multiplied by the jacobians of the test-function and covers "diffusive fluxes" as well as advective fluxes.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::flux().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::flux(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedFlux().

◆ fluxDivergence()

template<class PartsImpl >
template<class Entity , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::fluxDivergence ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
const HessianRangeType &  hessian,
RangeType &  result 
) const
inline

Compute the point-wise value of the flux-part of the operator, meaning the part of the differential operator which is multiplied by the derivative of the test function.

Primarily useful for residual error estimators.

Parameters
[in]entityThe mesh element.
[in]xThe point of evaluation, local coordinates.
[in]valueIn preparation for non-linear problems.
[in]jacobianThe value of \(\nabla u\).
[in]hessianThe value of \(D^2u\).
[out]resultThe result of the computation.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::fluxDivergence().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::fluxDivergence().

◆ linearizedFlux()

template<class PartsImpl >
template<class Entity , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedFlux ( const RangeType &  uBar,
const JacobianRangeType &  DuBar,
const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
JacobianRangeType &  flux 
) const
inline

Evaluate the linearized flux in local coordinates.

This can be interpreted as a diffusive "flux", at least when restricted to some surface. \(\bar u\) denotes the point of linearization for non-linear problems.

Parameters
[in]uBarPoint of linearization.
[in]DuBarPoint of linearization.
[in]entityThe mesh element.
[in]xThe point of evaluation, local coordinates.
[in]valueTo allow integration by parts for first order terms and in preparation for non-linear models.
[in]jacobianThe value of \(\nabla\bar\psi\).
[out]fluxThe result of the computation.
Note
"flux" in this sense simply means something which is multiplied by the jacobians of the test-function and covers "diffusive fluxes" as well as advective fluxes.

References Dune::ACFem::asImp(), Dune::ACFem::OperatorPartsInterface< PartsImpl >::flux(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedFlux().

Referenced by Dune::ACFem::DefaultOperatorParts< PartsImpl >::flux(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedFlux().

◆ linearizedRobinFlux()

template<class PartsImpl >
template<class Intersection , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedRobinFlux ( const RangeType &  uBar,
const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inline

The linearized Robin-type flux term.

Parameters
[in]uBarThe point of linearization.
[in]intersectionThe current intersection.
[in]xThe point of evaluation, local coordinates.
[in]unitOuterNormalThe outer normal (outer with respect to the current entity).
[in]valueThe value of u.
[out]resultThe result of the computation.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedRobinFlux().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedRobinFlux(), and Dune::ACFem::DefaultOperatorParts< PartsImpl >::robinFlux().

◆ linearizedSource()

template<class PartsImpl >
template<class Entity , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedSource ( const RangeType &  uBar,
const JacobianRangeType &  DuBar,
const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inline

The linearized source term as function of local coordinates.

Note the "source" in this context includes all terms which are multiplied by the value (not the jacobian) of the test-functions und thus may also include first order terms.

Parameters
[in]uBarThe point of linearization.
[in]DuBarThe jacobian at the point of linearization.
[in]entityThe mesh element.
[in]xThe point of evaluation, local coordinates.
[in]valueThe value of u.
[in]jacobianThe jacobian of u.
[out]resultThe result of the computation.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedSource().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::linearizedSource(), and Dune::ACFem::DefaultOperatorParts< PartsImpl >::source().

◆ robinFlux()

template<class PartsImpl >
template<class Intersection , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::robinFlux ( const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inline

The non-linearized Robin-type flux term.

This is "@f$\alpha@f$" from the Robin boundary condition.

Parameters
[in]intersectionThe current intersection.
[in]xThe point of evaluation, local coordinates.
[in]unitOuterNormalThe outer normal (outer with respect to the current entity).
[in]valueThe value of u.
[out]resultThe result of the computation.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::robinFlux().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::robinFlux().

◆ setEntity()

template<class PartsImpl >
template<class Entity >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::setEntity ( const Entity &  entity) const
inline

Per entity initialization, if that is needed.

Useful for constructing LocalFunction instances if the coefficients depend on discrete data, for example. A cleaner solution would be to provide a "LocalOperatorGerm" class which would be constructed from en EntityType in the same way as the Fem::LocalFunction is constructed from an entity (and a grid-function). OHTO, this would introduce yet another data-structure ...

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::setEntity().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::setEntity().

◆ setIntersection()

template<class PartsImpl >
template<class Intersection >
bool Dune::ACFem::OperatorPartsInterface< PartsImpl >::setIntersection ( const Intersection &  intersection) const
inline

Per-intersection initialization for the boundary contributions.

Several parts of the boundary may carry different Robin-kind boundary conditions. If necessary, per-intersection initializations can be performed here.

Parameters
[in]intersectionThe Intersection implementation for the underlying grid.
Returns
The function must return true if Robin-boundary conditions apply for intersection.
Note
Potential Dirichlet boundary conditions can be ignored here, the higher level code will only try to assemble boundary contributions if Dirichlet-conditions do not apply. It is therefore safe to constantly return true if there is only a split between a Dirichlet boundary part and a Robin boundary part. However, for more complicated cases it is vital that robinFlux() and linearizedRobinFlux() evaluate to zero if setIntersection() return false in order for ModelExpressions to work correctly,

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::setIntersection().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::setIntersection().

◆ source()

template<class PartsImpl >
template<class Entity , class Point >
void Dune::ACFem::OperatorPartsInterface< PartsImpl >::source ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inline

The zero-order term as function of local coordinates.

Has to evaluate \(\nabla\cdot(b(x,u)\,u) + c(x, u)\,u\).

Parameters
[in]entityThe mesh element.
[in]xThe point of evaluation, local coordinates.
[in]valueThe value of u.
[in]jacobianThe jacobian of u.
[out]resultThe result of the computation.

References Dune::ACFem::asImp(), and Dune::ACFem::OperatorPartsInterface< PartsImpl >::source().

Referenced by Dune::ACFem::OperatorPartsInterface< PartsImpl >::source().


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.111.3 (Nov 12, 23:30, 2024)