DUNE-ACFEM (unstable)

Dune::ACFem::PDEModel::ModelFacade< ModelImpl > Class Template Reference

A class defining the "closure" type of all supported model-method and method-call-signatures. More...

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

Public Member Functions

 ModelFacade (const ModelType &model)
 Constructor from Model. More...
 
 ModelFacade (const ModelFacade &other)
 Copy constructor. More...
 
 ModelFacade (ModelFacade &&other)
 Move constructor. More...
 
 operator const ModelType & () const
 Type-cast to underlying model, const version.
 
 operator ModelType & ()
 Type-cast to underlying model.
 
ModelType & model ()
 Getter for wrapped model, non-const version.
 
const ModelType & model () const
 Getter for wrapped model, const version.
 
template<class Entity >
void bind (const Entity &entity)
 Bind to the given entity. More...
 
void unbind ()
 Unbind from the previously bound entity. More...
 
template<class Intersection >
BoundaryConditionsType classifyBoundary (const Intersection &intersection)
 Bind to the given intersection and classify the components w.r.t. More...
 
std::string name () const
 Print a descriptive name for debugging and output.
 
template<class Quadrature >
auto flux (const QuadraturePoint< Quadrature > &x, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More...
 
template<class Quadrature >
auto linearizedFlux (const DomainRangeType &uBar, const DomainJacobianRangeType &DuBar, const QuadraturePoint< Quadrature > &x, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 Evaluate the linearized flux in local coordinates. More...
 
template<class Quadrature >
auto source (const QuadraturePoint< Quadrature > &x, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The zero-order term as function of local coordinates. More...
 
template<class Quadrature >
auto linearizedSource (const DomainRangeType &uBar, const DomainJacobianRangeType &DuBar, const QuadraturePoint< Quadrature > &x, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The linearized source term as function of local coordinates. More...
 
template<class Quadrature >
auto robinFlux (const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The non-linearized Robin-type flux term. More...
 
template<class Quadrature >
auto linearizedRobinFlux (const DomainRangeType &uBar, const DomainJacobianRangeType &DuBar, const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The linearized Robin-type flux term. More...
 
template<class Quadrature >
auto singularFlux (const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The non-linearized singular boundary "flux" term. More...
 
template<class Quadrature >
auto linearizedSingularFlux (const DomainRangeType &uBar, const DomainJacobianRangeType &DuBar, const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal, const DomainRangeType &value, const DomainJacobianRangeType &jacobian) const
 The linearized "singular flux" term. More...
 
template<class Quadrature >
auto fluxDivergence (const QuadraturePoint< Quadrature > &x, const DomainRangeType &value, const DomainJacobianRangeType &jacobian, const DomainHessianRangeType &hessian) 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...
 
template<class Quadrature >
auto dirichlet (const QuadraturePoint< Quadrature > &x, const DomainRangeType &value) const
 Dirichlet values. More...
 
template<class Quadrature >
auto linearizedDirichlet (const DomainRangeType &uBar, const QuadraturePoint< Quadrature > &x, const DomainRangeType &value) const
 Linearized Dirichlet values. More...
 

Detailed Description

template<class ModelImpl>
class Dune::ACFem::PDEModel::ModelFacade< ModelImpl >

A class defining the "closure" type of all supported model-method and method-call-signatures.

ACFem model-implementations are sparse: a particular implementation only needs to implement the very methods and call-signatures it needs, of course according to certain rules, see TO-BE-DOCED-IN-ANOTHER-CONCEPT-DOC-PAGE.

This class fills all missing methods with dummy-implementations and fills the sparse call-signature patterns with dummy parameters and finally forwards everything to the wrapped model.

Note
The class-definition is only a stub in order to inject the documentation of the "model-interface" into Doxygen. The "real" ModelFacade data-type is the return type of the closure() function.

Constructor & Destructor Documentation

◆ ModelFacade() [1/3]

template<class ModelImpl >
Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::ModelFacade ( const ModelType &  model)
inline

Constructor from Model.

This will make a copy of the model.

◆ ModelFacade() [2/3]

template<class ModelImpl >
Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::ModelFacade ( const ModelFacade< ModelImpl > &  other)
inline

Copy constructor.

This will copy the model stored in the other facade.

◆ ModelFacade() [3/3]

template<class ModelImpl >
Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::ModelFacade ( ModelFacade< ModelImpl > &&  other)
inline

Move constructor.

This will move the underlying data from other to this facade.

Member Function Documentation

◆ bind()

◆ classifyBoundary()

template<class ModelImpl >
template<class Intersection >
BoundaryConditionsType Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::classifyBoundary ( const Intersection &  intersection)
inline

Bind to the given intersection and classify the components w.r.t.

to the kind of applicable boundary conditions.

Warning
Note that prior to calling this function the model has to be bound to the inside entity of the given intersection. Failing to do so generates undefined behaviour.
The result of calling the other boundary related methods without binding to an intersection is undefined.
If RESULT.first is false, then the result of calling any of the other boundary related functions is undefined. Philosophically, they should return 0 in this case, but in order to have decent performance they give a damn and just don't care.
If RESULT.first is true, then still you cannot rely on user-friendly behaviour:
  • only if the respective bit of RESULT.second is set to 1, then the Dirichlet value in this compoment is well-defined.
  • only if the respective bit of RESULT.second is set to 0, then the Robin value in this component is well defined.
Parameters
[in]intersectionThe intersection to bind to.
Returns
A tuple. First component is a bool which is true iff any of the boundary related data functions would result in non trivial results. Second component is a bitset of size dimRange which is true if the given component of the system is subject to Dirichlet boundary conditions and false if it is subject to Robin or Neumann boundary conditions. If first is false then the contents of the bitset is undefined.

References Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::model().

Referenced by Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateLocal().

◆ dirichlet()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::dirichlet ( const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value 
) const
inline

Dirichlet values.

Maybe slightly more complicated than necessary, but in order to have a consistent variational formulation we simply model Dirichlet values just the same way as all the other contributions. For ordinary Dirchlet constraints the method should just compute

\[ \text{result} = u(x) - g(x) \]

just the same way as all other non-linearized method also contain the "right hand side" contribution. Things become more interesting when doing arithmetic with models.

Parameters
[in]xThe point of evaluation.
[in]valueThe value of the Ansatz-function at x.
Returns
The result of the computation.

◆ flux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::flux ( const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) 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]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\ps\).
Returns
The 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.

Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateBoundary(), Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateBoundary(), Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateIntersection(), Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateIntersection(), Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateProcessBoundary(), and Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateProcessorBoundary().

◆ fluxDivergence()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::fluxDivergence ( const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian,
const DomainHessianRangeType &  hessian 
) 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]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\).
Returns
The result of the computation.

Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateLocal(), and Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateLocal().

◆ linearizedDirichlet()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::linearizedDirichlet ( const DomainRangeType &  uBar,
const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value 
) const
inline

Linearized Dirichlet values.

For ordinary Dirichlet constraints, this method should simply copy value to result.

Parameters
[in]uBarThe point of linearization, evaluated at x.
[in]xThe point (in space) of evaluation.
[in]valueThe value of the Ansatz function for the linearized problem at x.
Returns
The result of the computation.

◆ linearizedFlux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::linearizedFlux ( const DomainRangeType &  uBar,
const DomainJacobianRangeType &  DuBar,
const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) 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]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\).
Returns
The 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.

◆ linearizedRobinFlux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::linearizedRobinFlux ( const DomainRangeType &  uBar,
const DomainJacobianRangeType &  DuBar,
const QuadraturePoint< Quadrature > &  x,
const DomainType &  unitOuterNormal,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) const
inline

The linearized Robin-type flux term.

Parameters
[in]uBarThe point of linearization.
[in]DuBarThe point of linearization.
[in]xThe point of evaluation, local coordinates.
[in]unitOuterNormalThe outer normal (outer with respect to the current entity).
[in]valueThe value of u.
[in]jacobianThe jacobian of u.
Returns
The result of the computation.

◆ linearizedSingularFlux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::linearizedSingularFlux ( const DomainRangeType &  uBar,
const DomainJacobianRangeType &  DuBar,
const QuadraturePoint< Quadrature > &  x,
const DomainType &  unitOuterNormal,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) const
inline

The linearized "singular flux" term.

Something tested with the gradient of the test functions on the boundary (think about weak capillary boundary conditions, or weak Dirichlet conditions in the spirit of Nitsche)

\[ \int_{\Gamma_s} linearizedSingularFlux(\bar u, \nabla \bar u, x, u, \nabla u)\cdot \nabla \phi do \]

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

◆ linearizedSource()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::linearizedSource ( const DomainRangeType &  uBar,
const DomainJacobianRangeType &  DuBar,
const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) 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]xThe point of evaluation, local coordinates.
[in]valueThe value of u.
[in]jacobianThe jacobian of u.
Returns
The result of the computation.

◆ robinFlux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::robinFlux ( const QuadraturePoint< Quadrature > &  x,
const DomainType &  unitOuterNormal,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) const
inline

The non-linearized Robin-type flux term.

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

Parameters
[in]xThe point of evaluation, local coordinates.
[in]unitOuterNormalThe outer normal (outer with respect to the current entity).
[in]valueThe value of u.
[in]valueThe jacobian of u.
Returns
The result of the computation.
Note
Normally Robin boundary conditions do not depend on the jacobian. However, in order to support things like consistent weak Dirichlet boundary conditions in the style of Nitsche [2] we need to be a little bit more general. Of course, if the respective implemented Robin boundary conditions really need the jacobian then the model leaves the \(H^1\) context.

Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateBoundary(), and Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateBoundary().

◆ singularFlux()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::singularFlux ( const QuadraturePoint< Quadrature > &  x,
const DomainType &  unitOuterNormal,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) const
inline

The non-linearized singular boundary "flux" term.

Something tested with the gradient of the test functions on the boundary (think about weak capillary boundary conditions, or weak Dirichlet conditions in the spirit of Nitsche)

\[ \int_{\Gamma_s} singularFlux(x, u, \nabla u)\cdot \nabla \phi do \]

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

◆ source()

template<class ModelImpl >
template<class Quadrature >
auto Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::source ( const QuadraturePoint< Quadrature > &  x,
const DomainRangeType &  value,
const DomainJacobianRangeType &  jacobian 
) 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]xThe point of evaluation, local coordinates.
[in]valueThe value of u.
[in]jacobianThe jacobian of u.
Returns
The result of the computation.

Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateLocal(), and Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateLocal().

◆ unbind()

template<class ModelImpl >
void Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::unbind ( )
inline

Unbind from the previously bound entity.

Warning
Calling this method on an unbound model may cause undefined behaviour.

References Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::model().

Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateIntersection().


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 (May 4, 22:30, 2024)