DUNE-ACFEM (unstable)

Auto-generated dune-acfem Documentation

Goals of Dune::ACFem

Dune::ACFem started by extracting the examples fem-poisson and fem-heat from the fem-school-generator. While the school code works great, it has some potential to factor out common code to form a general support library. As we are primarily interested in solving diffusion dominated problems in the context of conforming finite element methods, we called this "factored out" module ACFem for Adaptive Continuous FEM. Example applications can be found a couple of small Dune-modules. In the style of the dune-fem-howto there are two stand-alone modules ellipt-dot-c and heat-dot-c which primarily consist of a single source file and a header and do the same as the howto- or school-code does with a PDE-model based on dynamic polymorphism.

More efficient are the static-interface examples acfem-reactiondiffusion, acfem-quasilinear and acfem-heat which use expression-template techniques in order to form quite complicated models from a set of basic predefined "atom-models" in a convenient way (see the example below).

ACFem's relation to Dune::Fem and the dune-fem-howto

Dune::ACFem is kind of a add-on module to Dune::Fem. Dune::ACFem is not meant or even suitable as replacement for the dune-fem-school or some other fem-school. ACFem hides lots of stuf (like: how do I iterate over the mesh, how do I compute boundary integrals, how do I perform numerical integration at all?) behind somewhat cooked interface and implementation classes. So it is not meant as "school code". However, once having worked myself through the school-code, I found it more convenient to factor out common stuff to a central place, a template library.

Introductory Model-Expression Example

The following code is the principal part of a quasi-linear reaction-diffusion example which finally reproduces an approximation to \(u(x) = e^{-10\,|x|^2}\). The zero-order term exhibits a quadratic non-linearity. Besides the usual setup-procedure (construct the mesh, choose discrete function type, define data-types) this is all what has to be provided by the application (aka "the user"). This is taken from a Dune-module called acfem-reactiondiffusion. Further example programs are acfem-quasilinear and acfem-heat which implement respective model-problems in the same manner.

Quasi-linear Reaction-Diffusion Source-Code

using namespace ACFem;
template<class HGridType>
std::pair<double, double> algorithm(HGridType &grid)
{
const unsigned dimDomain = FunctionSpaceType::dimDomain;
GridPartType gridPart(grid);
DiscreteFunctionSpaceType discreteSpace(gridPart);
DiscreteFunctionType solution("solution", discreteSpace);
// Define some basic ingredients
auto zeroFct = zeroFunction(solution);
auto X = identityGridFunction(gridPart);
auto X0 = coordinateGridFunction(gridPart, Dune::ACFem::X<0>());
auto X1 = coordinateGridFunction(gridPart, Dune::ACFem::X<1>());
const double C = 10.0;
auto exactSolution = exp(-C*sqr(X));
// boundary data
auto gN = BoundaryIdIndicator(1) * zeroFct;
auto gR = BoundaryIdIndicator(2) * (-exp(-C*(sqr(X0)+1.)));
auto gD = !(BoundaryIdIndicator(1) || BoundaryIdIndicator(2)) * exactSolution;
// boundary models
auto Dbc = dirichletBoundaryModel(gD);
auto Nbc = neumannBoundaryModel(gN);
auto Rbc = robinBoundaryModel(gR);
// bulk contributions
auto F = 2.0 * C * dimDomain * exactSolution;
auto DU_DPhi = laplacianModel(discreteSpace);
auto U_U_Phi = p_MassModel(3., discreteSpace);
// now define the discrete model ...
auto pdeModel = (DU_DPhi + 4.0*(C*C)*sqr(X)*exp(C*sqr(X))*U_U_Phi - F
+
Dbc + Nbc + C*Rbc);
// create adaptive scheme, with exact solution for testing
// purposes. The exact solution defaults to the ZeroGridFunction if
// not specified.
auto scheme = ellipticFemScheme(solution, pdeModel, exactSolution);
return adaptiveAlgorithm(scheme);
}
static std::pair< double, double > adaptiveAlgorithm(AdaptiveFemScheme &scheme, const std::string &prefix_="")
Stationary adaptive algorithm like adapt_method_stat() from poor old ALBERTA.
Definition: stationaryadaptivealgorithm.hh:42
constexpr auto dirichletBoundaryModel(T &&values, Indicator &&where, const std::string &name="")
Generate the zero model for the empty indicator.
Definition: dirichletmodel.hh:237
constexpr auto p_MassModel(PField &&p, const Object &object, const std::string &name="")
Generate Model for a (weak, of course) Mass.
Definition: pmassmodel.hh:128
static auto laplacianModel(const Object &object, const std::string &name="")
Generate Model for a (weak, of course) Laplacian.
Definition: laplacianmodel.hh:95
constexpr auto robinBoundaryModel(GridFunction &&values, Indicator &&where=Indicator(), const std::string &name="")
Generate a RobinBoundaryModel from given grid-function and boundary indicator.
Definition: robinmodel.hh:199
constexpr auto neumannBoundaryModel(Fct &&values, Indicator &&where=std::decay_t< Indicator >{}, const std::string &name="")
Generate NeumannBoundaryModel from given grid-function and boundary indicator.
Definition: neumannmodel.hh:189
auto ellipticFemScheme(DiscreteFunction &solution, const Model &model, const InitialGuess &initialGuess, const RHSFunctional &rhsFunctional, const std::string name="acfem.schemes.elliptic")
Adaptive fem-scheme for "elliptic" problems.
Definition: ellipticfemscheme.hh:542

Quasi-linear Reaction-Diffusion Explanations

The above example implements an elliptic toy problem with somewhat complicated boundary conditions in order to test the modular discrete-model frame-work of Dune::ACFem.

The "exact solution" \(u(x) = e^{-C\,|x|^2}\) is used in order to construct all data. The discretization will reproduce the "solution" on the unit-square with boundary ids like defined by the following DGF-file:

DGF
Interval
0 0
1 1
1 1
#
GridParameter
% longest or arbitrary (see DGF docu)
refinementedge longest
overlap 0
#
BoundaryDomain
default 3
1 0 0 1 0 % lower boundary
2 0 1 1 1 % upper boundary
#
#

The code imposes Dirichlet boundary condition on the left and right boundaries (boundary id 3), homogeneous boundary conditions on the lower boundary (boundary id 1) and Robin-boundary conditions

\[ \frac{\partial u}{\partial \nu} = C\,(-e^{-C\,|x|^2} - u) \]

on the upper boundary (boundary id 2). The bulk-equation reads:

\[ -\Delta\,u + 4\,C^2\,e^{C\,|x|^2}\,|x|^2\,|u|^2 = 20\,d\,e^{-C\,|x|^2} \]

After defining auxiliary variables the discrete model is finally constructed in symbolic notation

auto pdeModel = (DU_DPhi + 4.0*(C*C)*sqr(X)*exp(C*sqr(X))*U_U_Phi - F
+
Dbc + Nbc + C*Rbc);

and then passed via an EllipticFemScheme to an adaptiveAlgorithm():

SchemeType scheme(solution, pdeModel, exactSolution);
return adaptiveAlgorithm(scheme);
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 19, 22:31, 2024)