DUNE-ACFEM (2.5.1)

In the style of the dune-fem-school-generator we define here dynamic problem descriptions and an EllipticModel such that an application can switch between example problems at run-time. More...

Classes

struct  Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >
 A model for a second order elliptic operator. More...
 
class  Dune::ACFem::ProblemInterface< FunctionSpace >
 Problem interface which describes a second order elliptic boundary problem: More...
 
class  Dune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >
 problem interface class for time dependent problem descriptions, i.e. More...
 

Enumerations

enum  Dune::ACFem::ModelTraits< EllipticModel< FunctionSpace, GridPart, Problem > >::ConstituentFlags { Dune::ACFem::ModelTraits< EllipticModel< FunctionSpace, GridPart, Problem > >::hasFlux = true , Dune::ACFem::ModelTraits< EllipticModel< FunctionSpace, GridPart, Problem > >::hasSources = true }
 

Functions

 Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::EllipticModel (const ProblemType &problem)
 Constructor taking problem reference. More...
 
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::flux (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedFlux (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const
 
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::source (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedSource (const RangeType &uBar, const JacobianRangeType &DuBar, const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const
 
template<class Intersection , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::robinFlux (const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 
template<class Intersection , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedRobinFlux (const RangeType &uBar, const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const
 
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::fluxDivergence (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, const HessianRangeType &hessian, RangeType &result) const
 
BulkForcesFunctionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::bulkForcesFunction (const GridPartType &gridPart) const
 Return a grid functin adapter for the external forces, see class Dune::Fem::GridFunctionAdapter.
 
DirichletBoundaryFunctionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::dirichletBoundaryFunction (const GridPartType &gridPart) const
 Return a grid functin adapter for the Dirichlet boundary values, see class Dune::Fem::GridFunctionAdapter.
 
DirichletWeightFunctionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::dirichletWeightFunction (const GridPartType &gridPart) const
 Return the trivial constant one function as Dirichlet weight.
 
NeumannBoundaryFunctionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::neumannBoundaryFunction (const GridPartType &gridPart) const
 Return a grid functin adapter for the Neumann boundary values, see class Dune::Fem::GridFunctionAdapter.
 
ExactSolutionFunctionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::exactSolutionFunction (const GridPartType &gridPart) const
 Return a grid function for a "exact solution" for experimental convergence tests.
 
const ProblemType & Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::problem () const
 return reference to problem
 
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::has (OperatorPartsType what) const
 May be used for optimizations during assembly. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::f (const DomainType &x, RangeType &value) const
 the right hand side data (default = 0)
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::u (const DomainType &x, RangeType &value) const
 the exact solution (default = 0)
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::uJacobian (const DomainType &x, JacobianRangeType &value) const
 the jacobian of the exact solution (default = 0)
 
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isDirichletSegment (const int bndId, const DomainType &center) const
 Classification of the kind of boundary conditions which applies to a boundary segment with the given Id. More...
 
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isNeumannSegment (const int bndId, const DomainType &center) const
 Classification of the kind of boundary conditions which applies to a boundary segment with the given Id. More...
 
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isRobinSegment (const int bndId, const DomainType &center) const
 Classification of the kind of boundary conditions which applies to a boundary segment with the given Id. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::dirichletData (const DomainType &x, RangeType &value) const
 The Dirichlet boundary data. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::neumannData (const DomainType &x, RangeType &value) const
 The Neumann boundary data. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::robinData (const DomainType &x, const RangeType &u, RangeType &value) const
 The Robin boundary data. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::secondOrderCoefficient (const DomainType &x, const JacobianRangeType &gradient, JacobianRangeType &result) const
 This method has to implement the second order term for the weak formulation, it needs to compute. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::secondOrderCoefficient (const DomainType &x, const JacobianRangeType &Du, const HessianRangeType &D2u, RangeType &result) const
 This method has to implement the second order term for the point-wise operator. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::firstOrderCoefficient (const DomainType &x, const RangeType &u, const JacobianRangeType &Du, RangeType &result) const
 First order term with derivative on u. More...
 
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::zeroOrderCoefficient (const DomainType &x, const RangeType &u, RangeType &result) const
 Zero order coefficient. More...
 
void Dune::ACFem::ProblemInterface< FunctionSpace >::FunctionWrapper< id >::evaluate (const DomainType &x, RangeType &ret) const
 evaluate function
 
void Dune::ACFem::ProblemInterface< FunctionSpace >::FunctionWrapper< id >::jacobian (const DomainType &x, JacobianRangeType &jac) const
 jacobian of the function
 
 Dune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >::TransientProblemInterface (const TimeProviderType &timeProvider, double theta=0.0)
 constructor taking time provider
 
double Dune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >::time () const
 return current simulation time
 
double Dune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >::deltaT () const
 return current time step size ( \(\Delta_t\))
 
const TimeViewTypeDune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >::timeView () const
 return reference to Problem's time provider
 
TimeViewTypeDune::ACFem::TransientProblemInterface< FunctionSpace, TimeProvider >::timeView ()
 return reference to Problem's time provider
 

Variables

const ProblemType & Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::problem_
 Pointer to the underlying problem.
 

Problem data-instances

In principle it would be possible to use the GridFunction Wrapper to avoid references to temporaries, but simply constructing the data-functions and keep them here is much safer.

@[

DirichletBoundaryType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::dirichletFunction_
 
NeumannBoundaryType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::neumannFunction_
 
BulkForcesType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::bulkForcesFunction_
 
ExactSolutionType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::solutionFunction_
 
DirichletIndicatorType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::dirichletIndicator_
 
NeumannIndicatorType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::neumannIndicator_
 
RobinIndicatorType Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::robinIndicator_
 

Detailed Description

In the style of the dune-fem-school-generator we define here dynamic problem descriptions and an EllipticModel such that an application can switch between example problems at run-time.

In principle, however, this is not the way to go when writing efficient code for a single, well defined problem. In that case one would rather write the operator "germs" at the ModelInterface level directly.

As the EllipticModel takes a ProblemClass as template parameter, it is of course possible to write a new non-virtual Problem-class which provides the same methods as the ProblemInterface class but does not use virtual functions.

Nevertheless, the classes defined here a convenient way to get a start.

Enumeration Type Documentation

◆ ConstituentFlags

template<class FunctionSpace , class GridPart , class Problem >
enum Dune::ACFem::ModelTraits< EllipticModel< FunctionSpace, GridPart, Problem > >::ConstituentFlags
Enumerator
hasFlux 

non-zero flux()

hasSources 

non-zero sources()

Function Documentation

◆ dirichletData()

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::dirichletData ( const DomainType &  x,
RangeType &  value 
) const
inlinevirtual

The Dirichlet boundary data.

Calls u() by default.

Parameters
[in]xThe point of evaluation
[out]valueThe result.

References Dune::ACFem::ProblemInterface< FunctionSpace >::u().

◆ EllipticModel()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::EllipticModel ( const ProblemType &  problem)
inline

Constructor taking problem reference.

Parameters
[in]problemAn instance of a class derived from class ProblemInterface.

◆ firstOrderCoefficient()

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::firstOrderCoefficient ( const DomainType &  x,
const RangeType &  u,
const JacobianRangeType &  Du,
RangeType &  result 
) const
inlinevirtual

First order term with derivative on u.

For

\[\int\phi\,\nabla\cdot(b\,u) dx\]

this function has to implement

\[\nabla\cdot(b(x)\,u(x)\]

.

Parameters
[in]xPoint of evaluation.
[in]uu(x)
[in]DuFirst derivatives
[out]resultGuess what.

◆ flux()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::flux ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
JacobianRangeType &  flux 
) const
inline

◆ fluxDivergence()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::fluxDivergence ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
const HessianRangeType &  hessian,
RangeType &  result 
) const
inline

◆ has()

template<class FunctionSpace >
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::has ( OperatorPartsType  what) const
inlinevirtual

May be used for optimizations during assembly.

A derived class has to be careful to reimplement this method properly.

Parameters
[in]whatAn Id for the repective part of the operator.
Returns
true if the what is indeed present, for example what(DirichletBoundary) returns true if part of the boundary is subject to Dirichlet boundary conditions and false otherwise.

◆ isDirichletSegment()

template<class FunctionSpace >
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isDirichletSegment ( const int  bndId,
const DomainType &  center 
) const
inlinevirtual

Classification of the kind of boundary conditions which applies to a boundary segment with the given Id.

Parameters
[in]bndIdThe boundary classification from the underlying mesh.
[in]centerThe barycenter of the boundary facet.
Returns
true iff bndId or center belong to a boundary segment subject to Dirichlet conditions.

◆ isNeumannSegment()

template<class FunctionSpace >
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isNeumannSegment ( const int  bndId,
const DomainType &  center 
) const
inlinevirtual

Classification of the kind of boundary conditions which applies to a boundary segment with the given Id.

Parameters
[in]bndIdThe boundary classification from the underlying mesh.
[in]centerThe barycenter of the boundary facet.
Returns
true iff bndId belong to a boundary segment subject to Neumann conditions and neumannData() should be called to compute contributions to the right hand side.

◆ isRobinSegment()

template<class FunctionSpace >
virtual bool Dune::ACFem::ProblemInterface< FunctionSpace >::isRobinSegment ( const int  bndId,
const DomainType &  center 
) const
inlinevirtual

Classification of the kind of boundary conditions which applies to a boundary segment with the given Id.

Parameters
[in]bndIdThe boundary classification from the underlying mesh.
[in]centerThe barycenter of the boundary facet.
Returns
true iff bndId belong to a boundary segment subject to Robin conditions and robinData() should be called to compute contributions to the system matrix.

◆ linearizedFlux()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedFlux ( const RangeType &  uBar,
const JacobianRangeType &  DuBar,
const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
JacobianRangeType &  flux 
) const
inline

◆ linearizedRobinFlux()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Intersection , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedRobinFlux ( const RangeType &  uBar,
const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inline

◆ linearizedSource()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::linearizedSource ( const RangeType &  uBar,
const JacobianRangeType &  DuBar,
const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inline

◆ neumannData()

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::neumannData ( const DomainType &  x,
RangeType &  value 
) const
inlinevirtual

The Neumann boundary data.

This is used to compute contributions to the right hand side.

Parameters
[in]xThe point of evaluation
[out]valueThe result.

◆ robinData()

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::robinData ( const DomainType &  x,
const RangeType &  u,
RangeType &  value 
) const
inlinevirtual

The Robin boundary data.

This is used to compute contributions to the system matrix. This function has to compute

\[ \alpha(x)\,u(x) \]

We pass the boundary mormal as additional argument. The default implementation simply returns 0.

Parameters
[in]xThe point of evaluation
[in]uThe value of u at x.
[out]valueThe result.

◆ robinFlux()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Intersection , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::robinFlux ( const Intersection &  intersection,
const Point &  x,
const DomainType &  unitOuterNormal,
const RangeType &  value,
RangeType &  result 
) const
inline

◆ secondOrderCoefficient() [1/2]

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::secondOrderCoefficient ( const DomainType &  x,
const JacobianRangeType &  Du,
const HessianRangeType &  D2u,
RangeType &  result 
) const
inlinevirtual

This method has to implement the second order term for the point-wise operator.

This is needed by the error estimator.

\[ -\nabla\cdot(A(x)\nabla u) = (-\nabla\cdot A(x))\cdot\nabla u - A(x) D^2 u \]

Parameters
[in]xPoint of evaluation.
[in]DuFirst derivatives u.
[in]D2uSecond derivatives of u.
[out]resultGuess what.

◆ secondOrderCoefficient() [2/2]

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::secondOrderCoefficient ( const DomainType &  x,
const JacobianRangeType &  gradient,
JacobianRangeType &  result 
) const
inlinevirtual

This method has to implement the second order term for the weak formulation, it needs to compute.

\[ A(x)\nabla u(x) \]

Parameters
[in]xPoint of evaluation.
[in]gradientFirst derivatives.
[out]resultGuess what.

References Dune::ACFem::gradient().

◆ source()

template<class FunctionSpace , class GridPart , class Problem = ProblemInterface<FunctionSpace>>
template<class Entity , class Point >
void Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >::source ( const Entity &  entity,
const Point &  x,
const RangeType &  value,
const JacobianRangeType &  jacobian,
RangeType &  result 
) const
inline

◆ zeroOrderCoefficient()

template<class FunctionSpace >
virtual void Dune::ACFem::ProblemInterface< FunctionSpace >::zeroOrderCoefficient ( const DomainType &  x,
const RangeType &  u,
RangeType &  result 
) const
inlinevirtual

Zero order coefficient.

This method has to compute \(c(x)\,u(x)\).

Parameters
[in]xPoint of evaluation.
[in]uu(x)
[out]resultGuess what.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)