DUNE-ACFEM (unstable)
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
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]
|
inline |
Constructor from Model.
This will make a copy of the model.
◆ ModelFacade() [2/3]
|
inline |
Copy constructor.
This will copy the model stored in the other facade.
◆ ModelFacade() [3/3]
|
inline |
Move constructor.
This will move the underlying data from other to this facade.
Member Function Documentation
◆ bind()
|
inline |
Bind to the given entity.
- Parameters
-
[in] entity The entity to bind to.
- Warning
- Calling any other method without first binding the model results in undefined behaviour.
References Dune::ACFem::PDEModel::ModelFacade< ModelImpl >::model().
Referenced by Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateIntersection(), Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateIntersection(), Dune::ACFem::EllipticEstimator< DiscreteFunctionSpace, Model, Norm >::estimateLocal(), and Dune::ACFem::ParabolicEulerEstimator< OldSolutionFunction, TimeProvider, ImplicitModel, ExplicitModel, Norm >::estimateLocal().
◆ classifyBoundary()
|
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] intersection The intersection to bind to.
- Returns
- A tuple. First component is a
bool
which istrue
iff any of the boundary related data functions would result in non trivial results. Second component is a bitset of sizedimRange
which istrue
if the given component of the system is subject to Dirichlet boundary conditions andfalse
if it is subject to Robin or Neumann boundary conditions. Iffirst
isfalse
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()
|
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] x The point of evaluation. [in] value The value of the Ansatz-function at x.
- Returns
- The result of the computation.
◆ flux()
|
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] x The point of evaluation, local coordinates. [in] value To allow integration by parts for first order terms and in preparation for non-linear models. [in] jacobian The 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()
|
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] x The point of evaluation, local coordinates. [in] value In preparation for non-linear problems. [in] jacobian The value of \(\nabla u\). [in] hessian The 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()
|
inline |
Linearized Dirichlet values.
For ordinary Dirichlet constraints, this method should simply copy value to result.
- Parameters
-
[in] uBar The point of linearization, evaluated at x. [in] x The point (in space) of evaluation. [in] value The value of the Ansatz function for the linearized problem at x.
- Returns
- The result of the computation.
◆ linearizedFlux()
|
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] uBar Point of linearization. [in] DuBar Point of linearization. [in] x The point of evaluation, local coordinates. [in] value To allow integration by parts for first order terms and in preparation for non-linear models. [in] jacobian The 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()
|
inline |
The linearized Robin-type flux term.
- Parameters
-
[in] uBar The point of linearization. [in] DuBar The point of linearization. [in] x The point of evaluation, local coordinates. [in] unitOuterNormal The outer normal (outer with respect to the current entity). [in] value The value of u. [in] jacobian The jacobian of u.
- Returns
- The result of the computation.
◆ linearizedSingularFlux()
|
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] uBar The point of linearization. [in] DuBar The jacobian of uBar. [in] x The point of evaluation, local coordinates. [in] value The value of u. [in] jacobian The jacobian of u. [out] result The result of the computation.
◆ linearizedSource()
|
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] uBar The point of linearization. [in] DuBar The jacobian at the point of linearization. [in] x The point of evaluation, local coordinates. [in] value The value of u. [in] jacobian The jacobian of u.
- Returns
- The result of the computation.
◆ robinFlux()
|
inline |
The non-linearized Robin-type flux term.
This is "@f$\alpha@f$" from the Robin boundary condition.
- Parameters
-
[in] x The point of evaluation, local coordinates. [in] unitOuterNormal The outer normal (outer with respect to the current entity). [in] value The value of u. [in] value The 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()
|
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] x The point of evaluation, local coordinates. [in] value The value of u. [in] jacobian The jacobian of u. [out] result The result of the computation.
◆ source()
|
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] x The point of evaluation, local coordinates. [in] value The value of u. [in] jacobian The 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()
|
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:
- dune/acfem/models/modelfacade.hh