DUNE-ACFEM (2.5.1)

Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess > Class Template Reference

Adaptive fem-scheme for "elliptic" problems. More...

#include <dune/acfem/algorithms/ellipticfemscheme.hh>

+ Collaboration diagram for Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >:

Public Types

typedef DiscreteFunction DiscreteFunctionType
 Type of the discrete solution function.
 
typedef DiscreteFunctionType::GridType GridType
 type of hierarchic grid
 
typedef DiscreteFunctionType::GridPartType GridPartType
 type of the grid view
 
typedef DiscreteFunctionType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
 type of the discrete function space
 
typedef Model ModelImplementationType
 type of the mathematical model
 
typedef DiscreteFunctionSpaceType::FunctionSpaceType FunctionSpaceType
 type of function space
 
typedef SolverSelector< DiscreteFunctionType, OperatorPartsType > SolverSelectorType
 choose type of discrete function, Matrix implementation and solver implementation
 
typedef Fem::NewtonInverseOperator< LinearOperatorType, LinearInverseOperatorType > NonLinearInverseOperatorType
 Non-linear solver.
 
typedef Fem::RestrictProlongDefault< DiscreteFunctionTypeRestrictionProlongationType
 type of restriction/prolongation projection for adaptive simulations (use default here, i.e. More...
 
typedef Fem::AdaptationManager< GridType, RestrictionProlongationTypeAdaptationManagerType
 type of adaptation manager handling adapation and DoF compression
 
typedef ModelType::DirichletBoundaryFunctionType DirichletBoundaryFunctionType
 type of Dirichlet boundary values
 
typedef ModelType::DirichletWeightFunctionType DirichletWeightFunctionType
 type of Dirichlet weight function
 
typedef decltype(std::declval< DirichletBoundaryFunctionType >()/std::declval< typename DirichletWeightFunctionType::GridFunctionType >()) EffectiveDirichletFunctionType
 type of effective Dirichlet values
 
typedef ModelType::DirichletIndicatorType DirichletIndicatorType
 types for various data
 
typedef DirichletConstraints< LinearOperatorType, EffectiveDirichletFunctionType, DirichletIndicatorTypeConstraintsOperatorType
 type of Dirichlet constraints
 
typedef DifferentiableEllipticOperator< LinearOperatorType, OperatorPartsType, ConstraintsOperatorTypeOperatorType
 define differential operator
 
typedef EllipticEstimator< DiscreteFunctionSpaceType, ModelImplementationTypeEstimatorType
 type of error estimator
 
typedef MarkingStrategy< EstimatorTypeMarkingStrategyType
 type of marking strategy
 
typedef InitialGuess InitialGuessType
 Initial guess/exact solution.
 
typedef GridFunctionConverter< InitialGuessType, GridPartType >::WrappedGridFunctionType InitialGuessFunctionType
 The wrapped InitialGuessType (no-op if already instrumented with LocalFunction)
 
typedef std::tuple< DiscreteFunctionType *, const InitialGuessFunctionType * > IOTupleType
 type of input/output tuple
 
typedef DataOutput< GridType, IOTupleTypeDataOutputType
 type of data writer
 

Define the RHS functional.

The default for forces and Neumann is a ZeroGridFunction which will result in a ZeroFunctional which simply will do nothing.

For assembling, both will be added efficently together using FunctionalExpressions.

typedef ModelType::BulkForcesFunctionType BulkForcesFunctionType
 
typedef L2InnerProductFunctional< DiscreteFunctionSpaceType, BulkForcesFunctionType > BulkForcesFunctional
 
typedef ModelType::NeumannBoundaryFunctionType BoundaryFluxFunctionType
 
typedef L2BoundaryFunctional< DiscreteFunctionSpaceType, BoundaryFluxFunctionType, typename ModelType::NeumannIndicatorTypeBoundaryFluxFunctional
 

Public Member Functions

virtual void initialize ()
 initialize solution
 
virtual void solve (bool forceMatrixAssembling=true)
 Solve the system. More...
 
virtual bool mark (const double tolerance)
 mark elements for adaptation
 
virtual double estimate ()
 calculate error estimator
 
virtual void adapt ()
 do the adaptation for a given marking
 
virtual int output ()
 data I/O
 
virtual double residual () const
 calculate residual (in small l^2)
 
virtual double error () const
 calculate L2/H1 error
 
virtual size_t size () const
 return some measure about the number of DOFs in use
 

Protected Member Functions

 EllipticFemSchemeBase (DiscreteFunctionType &solution, const ModelType &model, const InitialGuessType &initialGuess, const std::string &name="ellipt")
 Constructor for the elliptic fem-scheme. More...
 
virtual void nonLinearSolve (DiscreteFunctionType &rhs)
 Run the full Newton-scheme ...
 
virtual void linearSolve (DiscreteFunctionType &rhs, bool forceMatrixAssembling)
 Perform only one step of the Newton scheme for the affine-linear case. More...
 

Detailed Description

template<class DiscreteFunction, class Model, class InitialGuess>
class Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >

Adaptive fem-scheme for "elliptic" problems.

The quotes are due to the fact, that neither the EllipticFemScheme nor the EllipticOperator used in the scheme make any assumptions on the ellipticity except for the solvers which are used inside. But even here in principle the EllipticOperator has the possibility to flag that it is indefinite.

This scheme provides the necessary modules

ESTIMATE - MARK - ADAPT - SOLVE - DATA I/O

which then can be used by an adaptive algorithm for stationary problems, see the exeternal Dune module ellipt-dot-c for an example.

Parameters
[in]DiscreteFunctionType of the discrte function for solution and test functions.
[in]ModelA model satisfying the ModelInterface which describes the operator at a local basis by providing "operator germs", i.e. half of the integrants which constitute the discrete bilinear forms.
[in]InitialGuessAn initial guess in order to start an iterative solver. This can be used in order to pass an "exact solution" for experimental convergence tests.

Member Typedef Documentation

◆ RestrictionProlongationType

template<class DiscreteFunction , class Model , class InitialGuess >
typedef Fem::RestrictProlongDefault<DiscreteFunctionType> Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::RestrictionProlongationType

type of restriction/prolongation projection for adaptive simulations (use default here, i.e.

LagrangeInterpolation)

Constructor & Destructor Documentation

◆ EllipticFemSchemeBase()

template<class DiscreteFunction , class Model , class InitialGuess >
Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::EllipticFemSchemeBase ( DiscreteFunctionType solution,
const ModelType model,
const InitialGuessType initialGuess,
const std::string &  name = "ellipt" 
)
inlineprotected

Constructor for the elliptic fem-scheme.

Parameters
[in]solutionMutable reference to an existing discrete function in order to store the solution.
[in]modelA const reference to something satisfying the ModelInterface.
[in]initialGuessAn initial guess for an iterative solver or an exact solution for experimental convergence tests. EllipticFemScheme::error() will compute the H1-distance from thise function.
[in]nameFancy prefix name for parameters and debugging.

Member Function Documentation

◆ linearSolve()

template<class DiscreteFunction , class Model , class InitialGuess >
virtual void Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::linearSolve ( DiscreteFunctionType rhs,
bool  forceMatrixAssembling 
)
inlineprotectedvirtual

Perform only one step of the Newton scheme for the affine-linear case.

This implies that an affine linear case is really allowed.

Referenced by Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::solve().

◆ solve()

template<class DiscreteFunction , class Model , class InitialGuess >
virtual void Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::solve ( bool  forceMatrixAssembling = true)
inlinevirtual

Solve the system.

Parameters
[in]forceMatrixAssemblingIf the solve step requires assembling of matrices, then by default the application may give a hint to the FemScheme that this might not be necessary. The actual implementation can then decide if it caches the assembled matrix from the last solve() invocation and reuses it. If forceAssemble is true then the implementation has to discard and cached assembled matrices and recompute them.

Implements Dune::ACFem::BasicFemScheme.

References Dune::ACFem::DirichletConstraints< JacobianOperator, BoundaryValues, Indicator >::constrain(), Dune::ACFem::isZero(), Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::linearSolve(), Dune::ACFem::EllipticFemSchemeBase< DiscreteFunction, Model, InitialGuess >::nonLinearSolve(), and Dune::ACFem::DirichletConstraints< JacobianOperator, BoundaryValues, Indicator >::zeroConstrain().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)