DUNE-ACFEM (unstable)

Generate Model-Building-Blocks Conveniently

Some utility function in order to conveniently define some standard models without having to go through the "typedef" trouble. More...

Functions

template<class Fct , std::enable_if_t<(true &&IsWrappableByConstLocalFunction< Fct >::value &&!GridFunction::HasRegularity< Fct, 1 >::value), int > = 0>
constexpr auto Dune::ACFem::PDEModel::bulkLoadFunctionModel (Fct &&f, const std::string &name="") noexcept
 Generate a BulkLoadFunctionModel for the "right hand side". 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. 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="") noexcept
 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 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...
 
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 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 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...
 
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.
 
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 Object , class Velocity >
constexpr auto Dune::ACFem::PDEModel::incompressibleTransportModel (Object &&object, Velocity &&velocity, const std::string &name="")
 Generate an advection-model object. 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...
 
template<class Object >
auto Dune::ACFem::PDEModel::massModel (const Object &object, const std::string &name="")
 Generate a mass model fitting the specified object. More...
 
template<class Object , class Regularization = Tensor::FieldVectorTensor<typename Object::FunctionSpaceType::RangeFieldType>>
auto Dune::ACFem::PDEModel::meanCurvatureModel (Regularization &&regularization, const Object &object, const std::string &name="")
 Generate a MeanCurvature-model fitting the specified object. More...
 
template<class Fct , class Indicator = EntireBoundaryIndicator, std::enable_if_t<(true &&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="") noexcept
 Generate NeumannBoundaryModel from given grid-function and boundary indicator. More...
 
template<class Model , class PenaltyFunction , class Symmetrize = SkeletonSymmetrizeDefault<Model>, std::enable_if_t<(true &&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<(true &&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<(true &&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<(true &&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.
 
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...
 
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 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 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 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

◆ bulkLoadFunctionModel()

template<class Fct , std::enable_if_t<(true &&IsWrappableByConstLocalFunction< Fct >::value &&!GridFunction::HasRegularity< Fct, 1 >::value), int > = 0>
constexpr auto Dune::ACFem::PDEModel::bulkLoadFunctionModel ( Fct &&  f,
const std::string &  name = "" 
)
noexcept

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]fThe L2-function for the RHS.
[in]nameOptional. 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().

◆ deformationTensorModel()

template<class Object >
auto Dune::ACFem::PDEModel::deformationTensorModel ( const Object &  object,
const std::string &  name = "" 
)
inline

Generate a deformation tensor model fitting the specified object.

References Dune::ACFem::Expressions::expressionClosure().

◆ dirichletBoundaryModel()

template<class T , class Indicator , std::enable_if_t< ExpressionTraits< Indicator >::isZero, int > = 0>
constexpr auto Dune::ACFem::PDEModel::dirichletBoundaryModel ( T &&  values,
Indicator &&  where = std::decay_t< Indicator >{},
const std::string &  name = "" 
)
noexcept

Generate the zero model for the empty indicator.

Generate DirichletBoundaryModel from given grid-function and boundary indicator.

Parameters
[in]valuesThe Dirichlet boundary values.
[in]whereIndicator which decides which part of the boundary is affected

References Dune::ACFem::PDEModel::zeroModel().

Referenced by main().

◆ dirichletZeroModel()

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.

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]objectSuper-object, the DirichletBoundaryModel generated will be compatible to this object.
[in]whereIndicator 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()

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.

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

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 = "" 
)
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

template<class Object >
static auto Dune::ACFem::PDEModel::fluidSelfTransportModel ( const Object &  object,
const std::string &  name = "" 
)
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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().

◆ gradientLoadModel()

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.

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

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 = "" 
)
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

template<class Object >
static auto Dune::ACFem::PDEModel::incompressibleSelfTransportModel ( const Object &  object,
const std::string &  name = "" 
)
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

template<class Object , class Velocity >
constexpr auto Dune::ACFem::PDEModel::incompressibleTransportModel ( Object &&  object,
Velocity &&  velocity,
const std::string &  name = "" 
)

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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]velocityThe advection-velocity. This must already be something with a localFunction() method. It is implicitly assumed that the divergence of this object vanishes.

◆ laplacianModel()

template<class Object >
static auto Dune::ACFem::PDEModel::laplacianModel ( const Object &  object,
const std::string &  name = "" 
)
inlinestatic

Generate Model for a (weak, of course) Laplacian.

Parameters
[in]objectSomething with a public FunctionSpaceType typedef.
[in]nameAn optional name for debugging and pretty-printing.

References Dune::ACFem::Expressions::expressionClosure().

Referenced by main().

◆ massModel()

template<class Object >
auto Dune::ACFem::PDEModel::massModel ( const Object &  object,
const std::string &  name = "" 
)
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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()

template<class Object , class Regularization = Tensor::FieldVectorTensor<typename Object::FunctionSpaceType::RangeFieldType>>
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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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().

◆ neumannBoundaryModel()

template<class Fct , class Indicator = EntireBoundaryIndicator, std::enable_if_t<(true &&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 = "" 
)
noexcept

Generate NeumannBoundaryModel from given grid-function and boundary indicator.

Parameters
[in]valuesThe Neumann boundary values.
[in]whereIndicator which decides which part of ]the boundary is affected

◆ nitscheDirichletModel()

template<class Model , class... T, std::enable_if_t<(true &&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.

This implies that Model is not an expression closure.

References Dune::ACFem::Expressions::expressionClosure().

◆ p_LaplacianModel()

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.

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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]pThe exponent, see P_LaplacianModel.
[in]objectSomething with a public FunctionSpaceType typedef.
[in]nameAn optional name for debugging and pretty-printing.

References Dune::ACFem::Expressions::expressionClosure().

◆ p_MassModel()

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.

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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]pThe exponent, see P_MassModel.
[in]objectSomething with a public FunctionSpaceType typedef.
[in]nameAn optional name for debugging and pretty-printing.

References Dune::ACFem::Expressions::expressionClosure().

◆ robinBoundaryModel()

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.

Parameters
[in]valuesThe Robin boundary values.
[in]whereIndicator which decides which part of the boundary is affected

◆ robinZeroModel()

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.

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]objectSuper-object, the RobinBoundaryModel generated will be compatible to this object.
[in]whereIndicator 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.

◆ transportModel()

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.

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]objectIgnored. We use only the type to get hold of the FunctionSpaceType and the GridPartType.
[in]nameOptional. 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]velocityThe advection-velocity. This must already be something with a localFunction() method.

◆ weakDirichletBoundaryModel()

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.

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]objectSuper-object, the WeakDirichletBoundaryModel generated will be compatible to this object.
[in]whereIndicator 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]

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,
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]tIgnored. We use only the type to get hold of the FunctionSpaceType.
[in]nameOptional. 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() [2/3]

template<class T , class F = Expressions::Closure, std::enable_if_t< IsPDEModel< T >::value, int > = 0>
auto Dune::ACFem::PDEModel::zeroModel ( const T &  ,
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]tIgnored. We use only the type to get hold of the FunctionSpaceType.
[in]nameOptional. 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() [3/3]

template<class T , class F = Expressions::Closure, std::enable_if_t< IsPDEModel< T >::value, int > = 0>
auto Dune::ACFem::PDEModel::zeroModel ( 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]tIgnored. We use only the type to get hold of the FunctionSpaceType.
[in]nameOptional. 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.

Creative Commons License   |  Legal Statements / Impressum  |  generated with Hugo v0.55.6 (Jul 16, 22:22, 2019)