DUNE-ACFEM (2.5.1)
Dune::ACFem::EllipticFemScheme< DiscreteFunction, Model, ZeroGridFunction< typename Model::FunctionSpaceType, typename Model::GridPartType > > Class Template Reference
Constructor for the elliptic fem-scheme. More...
#include <dune/acfem/algorithms/ellipticfemscheme.hh>

Public Types | |
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< DiscreteFunctionType > | RestrictionProlongationType |
type of restriction/prolongation projection for adaptive simulations (use default here, i.e. More... | |
typedef Fem::AdaptationManager< GridType, RestrictionProlongationType > | AdaptationManagerType |
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, DirichletIndicatorType > | ConstraintsOperatorType |
type of Dirichlet constraints | |
typedef DifferentiableEllipticOperator< LinearOperatorType, OperatorPartsType, ConstraintsOperatorType > | OperatorType |
define differential operator | |
typedef EllipticEstimator< DiscreteFunctionSpaceType, ModelImplementationType > | EstimatorType |
type of error estimator | |
typedef MarkingStrategy< EstimatorType > | MarkingStrategyType |
type of marking strategy | |
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, IOTupleType > | DataOutputType |
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::NeumannIndicatorType > | BoundaryFluxFunctional |
Public Member Functions | |
EllipticFemScheme (DiscreteFunctionType &solution, const ModelType &model, const InitialGuessType &initialGuess, const std::string &name="ellipt") | |
Constructor for the elliptic fem-scheme. More... | |
EllipticFemScheme (DiscreteFunctionType &solution, const ModelType &model, const std::string &name="ellipt") | |
Constructor for the elliptic fem-scheme. More... | |
virtual void | initialize () |
initialize solution | |
virtual void | solve (bool forceMatrixAssembling=true) |
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 | |
Protected Member Functions | |
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 Dune::ACFem::EllipticFemScheme< DiscreteFunction, Model, ZeroGridFunction< typename Model::FunctionSpaceType, typename Model::GridPartType > >
class Dune::ACFem::EllipticFemScheme< DiscreteFunction, Model, ZeroGridFunction< typename Model::FunctionSpaceType, typename Model::GridPartType > >
Constructor for the elliptic fem-scheme.
- Parameters
-
[in] solution Mutable reference to an existing discrete function in order to store the solution. [in] model A const reference to something satisfying the ModelInterface. [in] initialGuess An 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] name Fancy prefix name for parameters and debugging.
Member Typedef Documentation
◆ RestrictionProlongationType
|
inherited |
type of restriction/prolongation projection for adaptive simulations (use default here, i.e.
LagrangeInterpolation)
Constructor & Destructor Documentation
◆ EllipticFemScheme() [1/2]
template<class DiscreteFunction , class Model >
|
inline |
Constructor for the elliptic fem-scheme.
- Parameters
-
[in] solution Mutable reference to an existing discrete function in order to store the solution. [in] model A const reference to something satisfying the ModelInterface. [in] initialGuess An 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] name Fancy prefix name for parameters and debugging.
◆ EllipticFemScheme() [2/2]
template<class DiscreteFunction , class Model >
|
inline |
Constructor for the elliptic fem-scheme.
- Parameters
-
[in] solution Mutable reference to an existing discrete function in order to store the solution. [in] model A const reference to something satisfying the ModelInterface. [in] name Fancy prefix name for parameters and debugging.
- Note
- This constructor automatically constructs the InitialGuessType as ZeroGridFunction.
Member Function Documentation
◆ linearSolve()
|
inlineprotectedvirtualinherited |
Perform only one step of the Newton scheme for the affine-linear case.
This implies that an affine linear case is really allowed.
◆ solve()
|
inlinevirtualinherited |
Implements Dune::ACFem::BasicFemScheme.
The documentation for this class was generated from the following file:
- dune/acfem/algorithms/ellipticfemscheme.hh
