DUNE-ACFEM (2.5.1)
A model for a second order elliptic operator. More...
#include <dune/acfem/models/ellipticmodel.hh>
Collaboration diagram for Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >:Public Member Functions | |
| EllipticModel (const ProblemType &problem) | |
| Constructor taking problem reference. More... | |
| template<class Entity , class Point > | |
| void | flux (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, JacobianRangeType &flux) const |
| template<class Entity , class Point > | |
| void | 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 | source (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, RangeType &result) const |
| template<class Entity , class Point > | |
| void | 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 | robinFlux (const Intersection &intersection, const Point &x, const DomainType &unitOuterNormal, const RangeType &value, RangeType &result) const |
| template<class Intersection , class Point > | |
| void | 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 | fluxDivergence (const Entity &entity, const Point &x, const RangeType &value, const JacobianRangeType &jacobian, const HessianRangeType &hessian, RangeType &result) const |
| BulkForcesFunctionType | bulkForcesFunction (const GridPartType &gridPart) const |
| Return a grid functin adapter for the external forces, see class Dune::Fem::GridFunctionAdapter. | |
| DirichletBoundaryFunctionType | dirichletBoundaryFunction (const GridPartType &gridPart) const |
| Return a grid functin adapter for the Dirichlet boundary values, see class Dune::Fem::GridFunctionAdapter. | |
| DirichletWeightFunctionType | dirichletWeightFunction (const GridPartType &gridPart) const |
| Return the trivial constant one function as Dirichlet weight. | |
| NeumannBoundaryFunctionType | neumannBoundaryFunction (const GridPartType &gridPart) const |
| Return a grid functin adapter for the Neumann boundary values, see class Dune::Fem::GridFunctionAdapter. | |
| ExactSolutionFunctionType | exactSolutionFunction (const GridPartType &gridPart) const |
| Return a grid function for a "exact solution" for experimental convergence tests. | |
| const ProblemType & | problem () const |
| return reference to problem | |
| std::string | name () const |
| Print a descriptive name for debugging and output. | |
| TraitsType::template ForcesFunctionalTraits< DiscreteFunctionSpace >::FunctionalType | forcesFunctional (const DiscreteFunctionSpace &space) const |
Detailed Description
struct Dune::ACFem::EllipticModel< FunctionSpace, GridPart, Problem >
A model for a second order elliptic operator.
The model defines the different parts needed to compute the weak form of a second order elliptic operator. The model is constructed in a way that it could be used to implement the linearization of a potentially non-linear operator of the form
\[ \begin{split} -\nabla\cdot (A(x, u)\,\nabla u) + \nabla\cdot(b(x, u)\,u) + c(x, u)\,u &= f(x)\quad \text{ in }\Omega,\\ u &= g_D \text{ on }\Gamma_D,\\ (A(x, u)\nabla u)\cdot\nu + \alpha(x, u)\,u &= g_N \text{ on }\Gamma_R,\\ (A(x, u)\nabla u)\cdot\nu &= g_N \text{ on }\Gamma_N,\\ \end{split} \]
The underlying default problem interface, however, is not suitable for non-linear problems, i.e. it ignores any extra "non-linear" parameters passed to the local functions of the model.
The methods in this class provide all factors of the integrands not involving test functions. The multiplication by the test-functions is added in the class EllipticOperator. The weak formulation then reads:
\[ \int_\Omega (A(x,u)\,\nabla u)\cdot\nabla\phi\,dx + \int_\Omega (\nabla\cdot(b(x, u)\,u) + c(x, u)\,u)\,\phi\,dx + \int_{\Gamma_R} \alpha(x, u)\,u\,\phi\,do - \int_{\Gamma_R\cup\Gamma_N} g_N\,\phi\,do - \int_\Omega f(x)\,\phi\,dx = 0. \]
- Parameters
-
[in] FunctionSpace The function space type where everything lives; in particular defining the rang and domain of each function (Note: this is not the discrete function space). [in] Problem A class describing a linear elliptic problem in global coordinates. Defaults to class ProblemInterface.
Member Function Documentation
◆ forcesFunctional()
|
inlineinherited |
Generate an instance of a class defining a functional which forms part of the force-terms for the model.
The documentation for this struct was generated from the following file:
- dune/acfem/models/ellipticmodel.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden & Uni Heidelberg |
generated with Hugo v0.111.3
(Mar 3, 23:33, 2026)