DUNE-ACFEM (unstable)
Some utility function in order to conveniently define some standard models without having to go through the "typedef" trouble. More...
Functions | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::BulkLoadFunctionModel< GridFunction >::bind (const Entity &entity) |
void | Dune::ACFem::PDEModel::BulkLoadFunctionModel< GridFunction >::unbind () |
template<class Quadrature > | |
RangeType | Dune::ACFem::PDEModel::BulkLoadFunctionModel< GridFunction >::source (const QuadraturePoint< Quadrature > &x) const |
template<class Fct , std::enable_if_t<(IsWrappableByConstLocalFunction< Fct >::value &&!GridFunction::HasRegularity< Fct, 1 >::value), int > = 0> | |
constexpr auto | Dune::ACFem::PDEModel::bulkLoadFunctionModel (Fct &&f, const std::string &name="") |
Generate a BulkLoadFunctionModel for the "right hand
side". More... | |
JacobianRangeType | Dune::ACFem::PDEModel::DeformationTensorModel< FunctionSpace >::linearizedFlux (const JacobianRangeType &jacobian) const |
Evaluate the linearized flux in local coordinates. More... | |
template<class Point > | |
RangeType | Dune::ACFem::PDEModel::DeformationTensorModel< FunctionSpace >::fluxDivergence (const HessianRangeType &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 Object > | |
auto | Dune::ACFem::PDEModel::deformationTensorModel (const Object &object, const std::string &name="") |
Generate a deformation tensor model fitting the specified object. | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<!IndicatorTraits< Indicator >::emptySupport &&!ExpressionTraits< GridFunction >::isZero > >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<!IndicatorTraits< Indicator >::emptySupport &&!ExpressionTraits< GridFunction >::isZero > >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<!IndicatorTraits< Indicator >::emptySupport &&!ExpressionTraits< GridFunction >::isZero > >::classifyBoundary (const Intersection &intersection) |
Bind to the given intersection and classify the components w.r.t. More... | |
template<class Quadrature > | |
RangeType | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<!IndicatorTraits< Indicator >::emptySupport &&!ExpressionTraits< GridFunction >::isZero > >::dirichlet (const QuadraturePoint< Quadrature > &x, const RangeType &value) const |
Dirichlet values. More... | |
auto | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<!IndicatorTraits< Indicator >::emptySupport &&!ExpressionTraits< GridFunction >::isZero > >::linearizedDirichlet (const RangeType &value) const |
Linearized Dirichlet values. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<(!IndicatorTraits< Indicator >::emptySupport &&ExpressionTraits< GridFunction >::isZero)> >::classifyBoundary (const Intersection &intersection) |
Bind to the given intersection and classify the components w.r.t. More... | |
auto | Dune::ACFem::PDEModel::DirichletBoundaryModel< GridFunction, Indicator, std::enable_if_t<(!IndicatorTraits< Indicator >::emptySupport &&ExpressionTraits< GridFunction >::isZero)> >::linearizedDirichlet (const RangeType &value) const |
Linearized Dirichlet values. More... | |
template<class T , class Indicator , std::enable_if_t< ExpressionTraits< Indicator >::isZero, int > = 0> | |
constexpr auto | Dune::ACFem::PDEModel::dirichletBoundaryModel (T &&values, Indicator &&where, const std::string &name="") |
Generate the zero model for the empty indicator. More... | |
template<class Object , class Indicator = EntireBoundaryIndicator, std::enable_if_t<!IsWrappableByConstLocalFunction< Object >::value, int > = 0> | |
auto | Dune::ACFem::PDEModel::dirichletZeroModel (const Object &object, Indicator &&where=std::decay_t< Indicator >{}, const std::string &name="") |
Generate homogeneous Dirichlet boundary conditions fitting the specified object. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::DivergenceLoadModel< GridFunction >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::DivergenceLoadModel< GridFunction >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Quadrature > | |
RangeType | Dune::ACFem::PDEModel::DivergenceLoadModel< GridFunction >::source (const QuadraturePoint< Quadrature > &x) const |
The zero-order term as function of local coordinates. More... | |
template<class GridFunction > | |
constexpr auto | Dune::ACFem::PDEModel::divergenceLoadModel (GridFunction &&f, const std::string &name="") |
Generate a divergence model which only contributes to the load-vector. More... | |
RangeRangeType | Dune::ACFem::PDEModel::DivergenceModel< FunctionSpace >::linearizedSource (const DomainJacobianRangeType &jacobian) const |
The linearized source term as function of local coordinates. More... | |
template<class Object , std::enable_if_t< Object::FunctionSpaceType::ScalarFunctionSpaceType::dimRange==1, int > = 0> | |
static auto | Dune::ACFem::PDEModel::divergenceModel (const Object &object, const std::string &name="") |
Generate a divergence model from some object which has a function-space. More... | |
template<class GridPartTraits > | |
static auto | Dune::ACFem::PDEModel::divergenceModel (const Fem::GridPartInterface< GridPartTraits > &gridPart, const std::string &name="") |
Generate a DivergenceModel from a GridPart, using ctype as field and dimensionWorld as dimension. | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::classifyBoundary (const Intersection &intersection) |
Bind to the given intersection and classify the components w.r.t. More... | |
JacobianRangeType | Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::flux (const RangeType &value) const |
Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More... | |
template<class Point > | |
JacobianRangeType | Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::linearizedFlux (const RangeType &uBar, const RangeType &value) const |
Evaluate the linearized flux in local coordinates. More... | |
RangeType | Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::fluxDivergence (const RangeType &value, const JacobianRangeType &jacobian) 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 Object > | |
static auto | Dune::ACFem::PDEModel::fluidSelfTransportModel (const Object &object, const std::string &name="") |
Generate a Navier-Stokes non-linearity fitting the given object. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::flux (const QuadraturePoint< Quadrature > &x) const |
Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::fluxDivergence (const QuadraturePoint< Quadrature > &x) 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 Intersection > | |
auto | Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::classifyBoundary (const Intersection &intersection) const |
Bind to the given intersection and classify the components w.r.t. More... | |
template<class GridFunction > | |
constexpr auto | Dune::ACFem::PDEModel::gradientLoadModel (GridFunction &&f, const std::string &name="") |
Generate a Gradient-model as contribution to the load vector from the given grid-function. More... | |
auto | Dune::ACFem::PDEModel::GradientModel< FunctionSpace >::linearizedFlux (const DomainRangeType &value) const |
The linearized source term as function of local coordinates. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::GradientModel< FunctionSpace >::classifyBoundary (const Intersection &intersection) const |
Bind to the given intersection and classify the components w.r.t. More... | |
template<class Object , std::enable_if_t< Object::FunctionSpaceType::ScalarFunctionSpaceType::dimRange==1, int > = 0> | |
static auto | Dune::ACFem::PDEModel::gradientModel (const Object &object, const std::string &name="") |
Generate a gradient model from some object which has a function-space. More... | |
template<class GridPartTraits > | |
static auto | Dune::ACFem::PDEModel::gradientModel (const Fem::GridPartInterface< GridPartTraits > &gridPart, const std::string &name="") |
Generate a GradientModel from a GridPart, using ctype as field and dimensionWorld as dimension. | |
auto | Dune::ACFem::PDEModel::IncompressibleSelfTransportModel< FunctionSpace >::source (const RangeType &value, const JacobianRangeType &jacobian) const |
The zero-order term as function of local coordinates. More... | |
auto | Dune::ACFem::PDEModel::IncompressibleSelfTransportModel< FunctionSpace >::linearizedSource (const RangeType &uBar, const JacobianRangeType &DuBar, const RangeType &value, const JacobianRangeType &jacobian) const |
The linearized source term as function of local coordinates. More... | |
template<class Object > | |
static auto | Dune::ACFem::PDEModel::incompressibleSelfTransportModel (const Object &object, const std::string &name="") |
Generate a Navier-Stokes non-linearity fitting the given object. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::IncompressibleTransportModel< FunctionSpace, GridFunction >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::IncompressibleTransportModel< FunctionSpace, GridFunction >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::IncompressibleTransportModel< FunctionSpace, GridFunction >::linearizedSource (const QuadraturePoint< Quadrature > &x, const JacobianRangeType &jacobian) const |
The linearized source term as function of local coordinates. More... | |
template<class Object , class Velocity > | |
constexpr auto | Dune::ACFem::PDEModel::incompressibleTransportModel (Object &&object, Velocity &&velocity, const std::string &name="") |
Generate an advection-model object. More... | |
auto | Dune::ACFem::PDEModel::LaplacianModel< FunctionSpace >::linearizedFlux (const JacobianRangeType &jacobian) const |
Evaluate the linearized flux in local coordinates. More... | |
auto | Dune::ACFem::PDEModel::LaplacianModel< FunctionSpace >::fluxDivergence (const HessianRangeType &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 Object > | |
static auto | Dune::ACFem::PDEModel::laplacianModel (const Object &object, const std::string &name="") |
Generate Model for a (weak, of course) Laplacian. More... | |
RangeType | Dune::ACFem::PDEModel::MassModel< FunctionSpace >::linearizedSource (const RangeType &value) const |
The linearized source term as function of local coordinates. More... | |
template<class Object > | |
auto | Dune::ACFem::PDEModel::massModel (const Object &object, const std::string &name="") |
Generate a mass model fitting the specified object. More... | |
JacobianRangeType | Dune::ACFem::PDEModel::MeanCurvatureModel< FunctionSpace, Regularization >::flux (const JacobianRangeType &jacobian) const |
Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More... | |
JacobianRangeType | Dune::ACFem::PDEModel::MeanCurvatureModel< FunctionSpace, Regularization >::linearizedFlux (const JacobianRangeType &DuBar, const JacobianRangeType &jacobian) const |
Evaluate the linearized flux in local coordinates. More... | |
RangeType | Dune::ACFem::PDEModel::MeanCurvatureModel< FunctionSpace, Regularization >::fluxDivergence (const JacobianRangeType &jacobian, const HessianRangeType &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 Object , class Regularization = Tensor::FieldVectorTensor<typename Object::FunctionSpaceType::RangeFieldType>> | |
auto | Dune::ACFem::PDEModel::meanCurvatureModel (Regularization &®ularization, const Object &object, const std::string &name="") |
Generate a MeanCurvature-model fitting the specified object. More... | |
std::string | Dune::ACFem::PDEModel::NeumannBoundaryModel< GridFunction, Indicator >::name () const |
template<class Entity > | |
void | Dune::ACFem::PDEModel::NeumannBoundaryModel< GridFunction, Indicator >::bind (const Entity &entity) |
void | Dune::ACFem::PDEModel::NeumannBoundaryModel< GridFunction, Indicator >::unbind () |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::NeumannBoundaryModel< GridFunction, Indicator >::classifyBoundary (const Intersection &intersection) |
template<class Quadrature > | |
RangeType | Dune::ACFem::PDEModel::NeumannBoundaryModel< GridFunction, Indicator >::robinFlux (const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal) const |
The non-linearized Robin-type flux term. More... | |
template<class Fct , class Indicator = EntireBoundaryIndicator, std::enable_if_t<(IsWrappableByConstLocalFunction< Fct >::value &&!GridFunction::HasRegularity< Fct, 1UL >::value), int > = 0> | |
constexpr auto | Dune::ACFem::PDEModel::neumannBoundaryModel (Fct &&values, Indicator &&where=std::decay_t< Indicator >{}, const std::string &name="") |
Generate NeumannBoundaryModel from given grid-function and boundary indicator. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::NitscheDirichletBoundaryModel< Model, PenaltyFunction, Symmetrize >::bind (const Entity &entity) |
Unbind from the previously bound entity. More... | |
void | Dune::ACFem::PDEModel::NitscheDirichletBoundaryModel< Model, PenaltyFunction, Symmetrize >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::NitscheDirichletBoundaryModel< Model, PenaltyFunction, Symmetrize >::classifyBoundary (const Intersection &intersection) |
Bind to the given intersection and classify the components w.r.t. More... | |
template<class Model , class PenaltyFunction , class Symmetrize = SkeletonSymmetrizeDefault<Model>, std::enable_if_t<(IsProperPDEModel< Model >::value &&IsWrappableByConstLocalFunction< PenaltyFunction >::value &&!Expressions::IsPromotedTopLevel< PenaltyFunction >::value &&ModelMethodExists< Model, ModelIntrospection::dirichlet >::value &&IsSign< Symmetrize >::value), int > = 0> | |
auto | Dune::ACFem::PDEModel::nitscheDirichletModel (const Model &m, const PenaltyFunction &p, Symmetrize=Symmetrize{}) |
Wrap an existing model converting any Dirichlet boundary conditions into weak Dirichlet boundary conditions with a general penalty function. | |
template<class Model , class GridPart , class Param = decltype(ModelParameters::nitscheDirichletPenalty()), class Symmetrize = SkeletonSymmetrizeDefault<Model>, std::enable_if_t<(IsPDEModel< Model >::value &&ModelMethodExists< Model, ModelIntrospection::dirichlet >::value &&IsTensor< Param >::value &&IsSign< Symmetrize >::value), int > = 0> | |
auto | Dune::ACFem::PDEModel::nitscheDirichletModel (Model &&m, const GridPart &gridPart, Param &&p=ModelParameters::nitscheDirichletPenalty(), Symmetrize=Symmetrize{}) |
Wrap existing model converting any Dirichlet boundary conditions into weak Dirichlet boundary conditions with a standard penalty function. | |
template<class Model , class... T, std::enable_if_t<(IsProperPDEModel< Model >::value &&!ModelMethodExists< Model, ModelIntrospection::dirichlet >::value), int > = 0> | |
constexpr auto | Dune::ACFem::PDEModel::nitscheDirichletModel (Model &&m, T &&...) |
Just return the original model if it does not have Dirichlet boundary conditions. More... | |
template<class GridFunction , class Param , class Indicator = EntireBoundaryIndicator, class Symmetrize = Sign<1>, std::enable_if_t<(IsWrappableByConstLocalFunction< GridFunction >::value &&IsBoundaryIndicator< Indicator >::value &&IsTensor< Param >::value &&IsSign< Symmetrize >::value), int > = 0> | |
auto | Dune::ACFem::PDEModel::nitscheDirichletModel (GridFunction &&values, Param &&p=ModelParameters::nitscheDirichletPenalty(), Indicator &&where=Indicator(), Symmetrize=Symmetrize{}) |
Create a model supplying weak Dirichet conditions from a given grid-function. | |
auto | Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::flux (const JacobianRangeType &jacobian) const |
Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More... | |
auto | Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::linearizedFlux (const JacobianRangeType &DuBar, const JacobianRangeType &jacobian) const |
Evaluate the linearized flux in local coordinates. More... | |
auto | Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::fluxDivergence (const JacobianRangeType &jacobian, const HessianRangeType &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 Object , class PField > | |
constexpr auto | Dune::ACFem::PDEModel::p_LaplacianModel (PField &&p, const Object &object, const std::string &name="") |
Generate Model for a (weak, of course) p-Laplacian. More... | |
RangeType | Dune::ACFem::PDEModel::P_MassModel< FunctionSpace, PField >::source (const RangeType &value) const |
The zero-order term as function of local coordinates. More... | |
RangeType | Dune::ACFem::PDEModel::P_MassModel< FunctionSpace, PField >::linearizedSource (const RangeType &uBar, const RangeType &value) const |
The linearized source term as function of local coordinates. More... | |
template<class Object , class PField > | |
constexpr auto | Dune::ACFem::PDEModel::p_MassModel (PField &&p, const Object &object, const std::string &name="") |
Generate Model for a (weak, of course) Mass. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::RobinBoundaryModel< FunctionSpace, Indicator >::classifyBoundary (const Intersection &intersection) |
auto | Dune::ACFem::PDEModel::RobinBoundaryModel< FunctionSpace, Indicator >::linearizedRobinFlux (const DomainType &unitOuterNormal, const RangeType &value) const |
template<class Object , class Indicator = EntireBoundaryIndicator> | |
constexpr auto | Dune::ACFem::PDEModel::robinZeroModel (const Object &object, Indicator &&where=Indicator{}, const std::string &name="") |
Generate homogeneous Robin boundary conditions fitting the specified object. More... | |
template<class GridFunction , class Indicator = EntireBoundaryIndicator> | |
constexpr auto | Dune::ACFem::PDEModel::robinBoundaryModel (GridFunction &&values, Indicator &&where=Indicator(), const std::string &name="") |
Generate a RobinBoundaryModel from given grid-function and boundary indicator. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Intersection > | |
auto | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::classifyBoundary (const Intersection &intersection) |
Bind to the given intersection and classify the components w.r.t. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::linearizedFlux (const QuadraturePoint< Quadrature > &x, const RangeType &value) const |
Evaluate the linearized flux in local coordinates. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::fluxDivergence (const QuadraturePoint< Quadrature > &x, const RangeType &value, const JacobianRangeType &jacobian) const |
! More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::TransportModel< FunctionSpace, GridFunction >::linearizedRobinFlux (const QuadraturePoint< Quadrature > &x, const DomainType &unitOuterNormal, const RangeType &value) const |
The linearized Robin-type flux term. More... | |
template<class Object , class GridFunction > | |
constexpr auto | Dune::ACFem::PDEModel::transportModel (const Object &object, GridFunction &&velocity, const std::string &name="") |
Generate an advection-model object. More... | |
template<class GridFunction , class Indicator = EntireBoundaryIndicator, std::enable_if_t< HasTag< GridFunction, Fem::HasLocalFunction >::value, int > = 0> | |
auto | Dune::ACFem::PDEModel::weakDirichletBoundaryModel (GridFunction &&values, Indicator &&where=Indicator{}, const double penalty=Dune::Fem::Parameter::getValue< double >("acfem.dgPenalty")) |
Generate homogeneous WeakDirichlet boundary conditions fitting the specified object. More... | |
template<class Entity > | |
void | Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::bind (const Entity &entity) |
Bind to the given entity. More... | |
void | Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::unbind () |
Unbind from the previously bound entity. More... | |
template<class Quadrature > | |
JacobianRangeType | Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::flux (const QuadraturePoint< Quadrature > &x) const |
Evaluate \(A(x, u)\nabla u(x)\) in local coordinates. More... | |
template<class Quadrature > | |
auto | Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::fluxDivergence (const QuadraturePoint< Quadrature > &x) 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 T , class F = Expressions::Closure, std::enable_if_t< IsPDEModel< T >::value, int > = 0> | |
auto | Dune::ACFem::PDEModel::zeroModel (const T &t, const std::string &name, F closure=F{}) |
Generate a zero model fitting the specified object. More... | |
template<class T , class F = Expressions::Closure, std::enable_if_t< IsPDEModel< T >::value, int > = 0> | |
auto | Dune::ACFem::PDEModel::zeroModel (const T &, F closure=F{}) |
Generate a zero model fitting the specified object. More... | |
template<class T , class F = Expressions::Closure, std::enable_if_t< IsPDEModel< T >::value, int > = 0> | |
auto | Dune::ACFem::PDEModel::zeroModel (F closure=F{}) |
Generate a zero model fitting the specified object. More... | |
Detailed Description
Some utility function in order to conveniently define some standard models without having to go through the "typedef" trouble.
Function Documentation
◆ bind() [1/9]
|
inline |
◆ bind() [2/9]
|
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.
◆ bind() [3/9]
|
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.
◆ bind() [4/9]
|
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.
◆ bind() [5/9]
|
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.
◆ bind() [6/9]
|
inline |
◆ bind() [7/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ bind() [8/9]
|
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.
◆ bind() [9/9]
|
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.
◆ bulkLoadFunctionModel()
|
constexpr |
Generate a BulkLoadFunctionModel for the "right hand side".
If the supplied grid-function is a BindableTensorFunction then first force the differentiability order to 0.
If Fct is a BindableTensorFunction then we first construct a non-differentiable variant and then feed it into the BulkLoadModel.
- Parameters
-
[in] f The L2-function for the RHS. [in] name Optional. If left empty some sensible name is built from f.name() for debugging purposes.
This spares stack-space and complexity.
References Dune::ACFem::Expressions::expressionClosure().
Referenced by main().
◆ classifyBoundary() [1/9]
|
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.
◆ classifyBoundary() [2/9]
|
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.
◆ classifyBoundary() [3/9]
|
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.
◆ classifyBoundary() [4/9]
|
inline |
◆ classifyBoundary() [5/9]
|
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.
◆ classifyBoundary() [6/9]
|
inline |
◆ classifyBoundary() [7/9]
|
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.
◆ classifyBoundary() [8/9]
|
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.
◆ classifyBoundary() [9/9]
|
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.
◆ 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.
◆ dirichletBoundaryModel()
|
constexpr |
Generate the zero model for the empty indicator.
Generate DirichletBoundaryModel from given grid-function and boundary indicator.
- Parameters
-
[in] values The Dirichlet boundary values. [in] where Indicator which decides which part of the boundary is affected
References Dune::ACFem::PDEModel::zeroModel().
Referenced by main().
◆ dirichletZeroModel()
auto Dune::ACFem::PDEModel::dirichletZeroModel | ( | const Object & | object, |
Indicator && | where = std::decay_t<Indicator>{} , |
||
const std::string & | name = "" |
||
) |
Generate homogeneous Dirichlet boundary conditions fitting the specified object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType and a method object.gridPart().
- Parameters
-
[in] object Super-object, the DirichletBoundaryModel generated will be compatible to this object. [in] where Indicator which decides where the Dirichlet b.c. apply.
Examples of suitable objects are
- another model
- Fem::DiscreteFunctionSpace
- Fem::BindableGridFunction
- See also
- CompatibleModel.
Referenced by main().
◆ divergenceLoadModel()
|
constexpr |
Generate a divergence model which only contributes to the load-vector.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
◆ divergenceModel()
|
inlinestatic |
Generate a divergence model from some object which has a function-space.
The generated DivergenceModel uses the dimDomain as its dimension.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
References Dune::ACFem::Expressions::expressionClosure().
◆ fluidSelfTransportModel()
|
inlinestatic |
Generate a Navier-Stokes non-linearity fitting the given object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
References Dune::ACFem::Expressions::expressionClosure().
◆ flux() [1/5]
|
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.
◆ flux() [2/5]
|
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.
References Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::flux(), and Dune::ACFem::Tensor::pow().
Referenced by Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::flux(), and Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::linearizedFlux().
◆ flux() [3/5]
|
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.
References Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::flux().
Referenced by Dune::ACFem::PDEModel::GradientLoadModel< GridFunction >::flux().
◆ flux() [4/5]
|
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.
References Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::flux().
Referenced by Dune::ACFem::PDEModel::WeakDivergenceLoadModel< GridFunction >::flux().
◆ flux() [5/5]
|
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.
References Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::flux().
Referenced by Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::flux(), and Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::linearizedFlux().
◆ fluxDivergence() [1/8]
|
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.
◆ fluxDivergence() [2/8]
|
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.
References Dune::ACFem::ModelBase< FunctionSpace >::dimDomain, and Dune::ACFem::ModelBase< FunctionSpace >::dimRange.
◆ fluxDivergence() [3/8]
|
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.
◆ fluxDivergence() [4/8]
|
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.
References Dune::ACFem::ModelBase< FunctionSpace >::dimDomain, Dune::ACFem::ModelBase< FunctionSpace >::dimRange, and Dune::ACFem::Tensor::pow().
◆ fluxDivergence() [5/8]
|
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.
◆ fluxDivergence() [6/8]
|
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.
This is the strong form, i.e. simply the divergence.
◆ fluxDivergence() [7/8]
|
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.
The implementation also works for transport-velocities which are not divergence free.
References Dune::ACFem::ModelBase< FunctionSpace >::dimDomain.
◆ fluxDivergence() [8/8]
|
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.
◆ gradientLoadModel()
|
constexpr |
Generate a Gradient-model as contribution to the load vector from the given grid-function.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
◆ gradientModel()
|
inlinestatic |
Generate a gradient model from some object which has a function-space.
The generated GradientModel uses the dimDomain as its dimension.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
References Dune::ACFem::Expressions::expressionClosure().
◆ incompressibleSelfTransportModel()
|
inlinestatic |
Generate a Navier-Stokes non-linearity fitting the given object.
This variant moves the derivative to the test function at the cost of introducing a boundary integral.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
References Dune::ACFem::Expressions::expressionClosure().
◆ incompressibleTransportModel()
|
constexpr |
Generate an advection-model object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
- Parameters
-
[in] velocity The advection-velocity. This must already be something with a localFunction() method. It is implicitly assumed that the divergence of this object vanishes.
◆ laplacianModel()
|
inlinestatic |
Generate Model for a (weak, of course) Laplacian.
- Parameters
-
[in] object Something with a public FunctionSpaceType typedef. [in] name An optional name for debugging and pretty-printing.
References Dune::ACFem::Expressions::expressionClosure().
Referenced by main().
◆ linearizedDirichlet() [1/2]
|
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.
◆ linearizedDirichlet() [2/2]
|
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() [1/7]
|
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.
References Dune::ACFem::PDEModel::GradientModel< FunctionSpace >::rangeDimRange.
◆ linearizedFlux() [2/7]
|
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.
◆ linearizedFlux() [3/7]
|
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.
References Dune::ACFem::PDEModel::P_LaplacianModel< FunctionSpace, PField >::flux(), and Dune::ACFem::Tensor::pow().
◆ linearizedFlux() [4/7]
|
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.
◆ linearizedFlux() [5/7]
|
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.
◆ linearizedFlux() [6/7]
|
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.
References Dune::ACFem::ModelBase< FunctionSpace >::dimRange.
◆ linearizedFlux() [7/7]
|
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.
References Dune::ACFem::PDEModel::FluidSelfTransportModel< FunctionSpace >::flux().
◆ linearizedRobinFlux() [1/2]
|
inline |
References Dune::ACFem::zero().
◆ linearizedRobinFlux() [2/2]
|
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.
◆ linearizedSource() [1/5]
|
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.
References Dune::ACFem::PDEModel::DivergenceModel< FunctionSpace >::domainDimRange.
◆ linearizedSource() [2/5]
|
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.
◆ linearizedSource() [3/5]
|
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.
◆ linearizedSource() [4/5]
|
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.
References Dune::ACFem::Tensor::pow().
◆ linearizedSource() [5/5]
|
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.
◆ massModel()
|
inline |
Generate a mass model fitting the specified object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
References Dune::ACFem::Expressions::expressionClosure().
Referenced by Dune::ACFem::l2Projection(), and main().
◆ meanCurvatureModel()
auto Dune::ACFem::PDEModel::meanCurvatureModel | ( | Regularization && | regularization, |
const Object & | object, | ||
const std::string & | name = "" |
||
) |
Generate a MeanCurvature-model fitting the specified object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
- Parameters
-
[in] regularization \(\eta\), see MeanCurvatureModel. For the graph-case this should be 1.0, for the level-set approach this is a regularization parameter in order to be able to cope with fattening.
References Dune::ACFem::Expressions::expressionClosure().
◆ name()
|
inline |
◆ neumannBoundaryModel()
|
constexpr |
Generate NeumannBoundaryModel from given grid-function and boundary indicator.
- Parameters
-
[in] values The Neumann boundary values. [in] where Indicator which decides which part of ]the boundary is affected
◆ nitscheDirichletModel()
|
constexpr |
Just return the original model if it does not have Dirichlet boundary conditions.
This implies that Model is not an expression closure.
References Dune::ACFem::Expressions::expressionClosure().
◆ p_LaplacianModel()
|
constexpr |
Generate Model for a (weak, of course) p-Laplacian.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
- Parameters
-
[in] p The exponent, see P_LaplacianModel. [in] object Something with a public FunctionSpaceType typedef. [in] name An optional name for debugging and pretty-printing.
References Dune::ACFem::Expressions::expressionClosure().
◆ p_MassModel()
|
constexpr |
Generate Model for a (weak, of course) Mass.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
- Parameters
-
[in] p The exponent, see P_MassModel. [in] object Something with a public FunctionSpaceType typedef. [in] name An optional name for debugging and pretty-printing.
References Dune::ACFem::Expressions::expressionClosure().
◆ robinBoundaryModel()
|
constexpr |
Generate a RobinBoundaryModel from given grid-function and boundary indicator.
- Parameters
-
[in] values The Robin boundary values. [in] where Indicator which decides which part of the boundary is affected
◆ robinFlux()
|
inline |
The non-linearized Robin-type flux term.
This is "@f$\alpha@f$" from the Robin boundary condition.
- Parameters
-
[in] intersection The current intersection. [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. [out] result The result of the computation.
◆ robinZeroModel()
|
constexpr |
Generate homogeneous Robin boundary conditions fitting the specified object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType and a method object.gridPart().
- Parameters
-
[in] object Super-object, the RobinBoundaryModel generated will be compatible to this object. [in] where Indicator which decides where the Robin b.c. apply.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
◆ source() [1/4]
|
inline |
◆ source() [2/4]
|
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.
References Dune::ACFem::Tensor::trace().
◆ source() [3/4]
|
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.
References Dune::ACFem::Tensor::pow().
◆ source() [4/4]
|
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.
◆ transportModel()
|
constexpr |
Generate an advection-model object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] object Ignored. We use only the type to get hold of the FunctionSpaceType and the GridPartType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
- Parameters
-
[in] velocity The advection-velocity. This must already be something with a localFunction() method.
◆ unbind() [1/9]
|
inline |
◆ unbind() [2/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [3/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [4/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [5/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [6/9]
|
inline |
◆ unbind() [7/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [8/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ unbind() [9/9]
|
inline |
Unbind from the previously bound entity.
- Warning
- Calling this method on an unbound model may cause undefined behaviour.
◆ weakDirichletBoundaryModel()
auto Dune::ACFem::PDEModel::weakDirichletBoundaryModel | ( | GridFunction && | values, |
Indicator && | where = Indicator{} , |
||
const double | penalty = Dune::Fem::Parameter::getValue<double>("acfem.dgPenalty") |
||
) |
Generate homogeneous WeakDirichlet boundary conditions fitting the specified object.
In order for this to work the data-type of the given object must have the two sub-types Object::FunctionSpaceType and Object::GridPartType and a method object.gridPart().
- Parameters
-
[in] object Super-object, the WeakDirichletBoundaryModel generated will be compatible to this object. [in] where Indicator which decides where the WeakDirichlet b.c. apply.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
◆ zeroModel() [1/3]
auto Dune::ACFem::PDEModel::zeroModel | ( | const T & | , |
F | closure = F{} |
||
) |
Generate a zero model fitting the specified object.
In order for this to work the data-type of the given object must have the sub-type Object::FunctionSpaceType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] t Ignored. We use only the type to get hold of the FunctionSpaceType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
Variant ommitting the name parameter.
◆ zeroModel() [2/3]
auto Dune::ACFem::PDEModel::zeroModel | ( | const T & | t, |
const std::string & | name, | ||
F | closure = F{} |
||
) |
Generate a zero model fitting the specified object.
In order for this to work the data-type of the given object must have the sub-type Object::FunctionSpaceType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] t Ignored. We use only the type to get hold of the FunctionSpaceType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
Referenced by Dune::ACFem::PDEModel::dirichletBoundaryModel(), and main().
◆ zeroModel() [3/3]
auto Dune::ACFem::PDEModel::zeroModel | ( | F | closure = F{} | ) |
Generate a zero model fitting the specified object.
In order for this to work the data-type of the given object must have the sub-type Object::FunctionSpaceType. The actual instance of the object is ignore and simply has to be passed in order that the compiler can deduce the correct types.
- Parameters
-
[in] t Ignored. We use only the type to get hold of the FunctionSpaceType. [in] name Optional. A descriptive name for the generated model. A suitable default will be chosen if omitted.
Examples of suitable objects are something that satisfies the
- ModelInterface (another model)
- Fem::DiscreteFunctionSpaceInterface
- Fem::DiscreteFunctionInterface (including any wrapped or adapted non-discrete function)
- See also
- CompatibleModel.
Variant for use without an instance of T, and without a name.