DUNE-ACFEM (unstable)

Abstract FEM algorithms and schemes. More...

Modules

 SchemeGenerators
 

Classes

class  Dune::ACFem::EllipticFemScheme< DiscreteFunction, Model, InitialGuess, RHSFunctional >
 Adaptive fem-scheme for "elliptic" problems. More...
 
class  Dune::ACFem::BasicFemScheme
 Abstract non-adaptive basic FEM scheme. More...
 
class  Dune::ACFem::AdaptiveFemScheme
 Abstract space adaptative FEM scheme. More...
 
struct  Dune::ACFem::OperatorTypeSelector< DiscreteFunctionSpace, Model, flavour >
 
class  Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >
 Basic parabolic fem-scheme class. More...
 

Enumerations

enum  Dune::ACFem::DGFlavour
 

Functions

static std::pair< double, double > Dune::ACFem::adaptiveAlgorithm (AdaptiveFemScheme &scheme, const std::string &prefix_="")
 Stationary adaptive algorithm like adapt_method_stat() from poor old ALBERTA. More...
 
template<class TimeProvider >
static double Dune::ACFem::adaptiveAlgorithm (TransientAdaptiveFemScheme &scheme, TimeProvider &timeProvider, double endTime, const std::string &prefix_="")
 Space-time adaptive algorithm like adapt_method_instat() from poor old ALBERTA. More...
 

Detailed Description

Abstract FEM algorithms and schemes.

Enumeration Type Documentation

◆ DGFlavour

Choice DG methods

Function Documentation

◆ adaptiveAlgorithm() [1/2]

static std::pair< double, double > Dune::ACFem::adaptiveAlgorithm ( AdaptiveFemScheme scheme,
const std::string &  prefix_ = "" 
)
inlinestatic

Stationary adaptive algorithm like adapt_method_stat() from poor old ALBERTA.

This could be moved into a pre-compiled library.

Parameters
[in,out]schemeThe transient fem-scheme used for the solve-estimate-mark cycles.
[in]prefix_A prefix like "poisson" to tag the varous parameters fetched from the parameter file. In particalar the tolerance is fetched from the Fem::Parameter class with "prefix_.adaptation.tolerance" key. Of course the "." is not added if prefix_ is empty.
Returns
A std::pair with pair.first being the final estimator value and pair.second intentionally being the real error for convergence tests. If no exact solution is known (general case) then the distance to some initial guess of the solution. See EllipticFemScheme.
Todo:
Document the relevant parameters fetched from the parameter file.

References Dune::ACFem::BasicFemScheme::name().

◆ adaptiveAlgorithm() [2/2]

template<class TimeProvider >
static double Dune::ACFem::adaptiveAlgorithm ( TransientAdaptiveFemScheme &  scheme,
TimeProvider &  timeProvider,
double  endTime,
const std::string &  prefix_ = "" 
)
inlinestatic

Space-time adaptive algorithm like adapt_method_instat() from poor old ALBERTA.

This could be moved into a pre-compiled library but for the fact that the TimeProvider is tagged by the GridType.

Note that the implementation of the TransientAdaptiveFemScheme has to take care of always using the time-step-size of the Dune::Fem::TimeProvider implementation, because that may be changed by the adaptive algorithm.

Parameters
[in,out]schemeThe transient fem-scheme used for the solve-estimate-mark cycles.
[in,out]timeProviderThe Dune::Fem::TimeProvider instance, also defines the start-time.The TimeProvider must be already initialized.
[in]endTimeEnd-of-simulation point in time.
[in]prefix_A prefix like "heat" to tag the varous parameters fetched from the parameter file.
Returns
The maximum estimator value from all time-steps.
Todo:
Document the parameter fetched from the parameter file.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)