DUNE-ACFEM (unstable)

Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits > Class Template Reference

Basic parabolic fem-scheme class. More...

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

+ Collaboration diagram for Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >:

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 TimeProvider TimeProviderType
 access to the simulation time and time-step
 
using ImplicitModelType = ImplicitModel
 type of the mathematical model
 
typedef InitialValue InitialValueType
 Initial value.
 
typedef DiscreteModelType::FunctionSpaceType FunctionSpaceType
 type of function space (scalar functions, \(f: \Omega -> R\))
 
typedef Fem::NewtonInverseOperator< LinearOperatorType, LinearInverseOperatorType > NonLinearInverseOperatorType
 Non-linear solver.
 
typedef Dune::Fem::RestrictProlongDefault< DiscreteFunctionTypeRestrictionProlongationType
 type of restriction/prolongation projection for adaptive simulations (use default here, i.e. More...
 
typedef Dune::Fem::AdaptationManager< GridType, RestrictionProlongationTypeAdaptationManagerType
 type of adaptation manager handling adapation and DoF compression
 
typedef EmptyBlockConstraints< DiscreteFunctionSpaceTypeEmptyConstraintsType
 empty constraints for the explicit operator (old solution is already constrained)
 
using ConstraintsOperatorType = DirichletConstraints< DiscreteFunctionSpaceType, DiscreteModelType >
 type of Dirichlet constraints
 
typedef DifferentiableEllipticOperator< LinearOperatorType, DiscreteModelType, const ConstraintsOperatorType &, QuadratureTraitsImplicitOperatorType
 define differential operator, implicit part
 
typedef EllipticOperator< ExplicitModelType, DiscreteFunctionType, DiscreteFunctionType, EmptyConstraintsType, QuadratureTraitsExplicitOperatorType
 explicit part of differential operator. More...
 
typedef ParabolicEulerEstimator< DiscreteFunctionType, TimeProviderType, DiscreteModelType, ExplicitModelType > EstimatorType
 type of error estimator
 
typedef MarkingStrategy< GridPartTypeMarkingStrategyType
 type of marking strategy for space adaptivity
 
typedef TrueErrorEstimator< InitialValueType, Dune::Fem::L2Norm< GridPartType > > InitialEstimatorType
 intial estimator
 
typedef MarkingStrategy< GridPartTypeInitialMarkingStrategyType
 Initial marking.
 
typedef InitialValueType ExactSolutionFunctionType
 adapter to turn exact solution into a grid function (for visualization)
 
typedef std::tuple< DiscreteFunctionType * > IOTupleType
 type of input/output tuple
 
typedef DataOutput< GridType, IOTupleTypeDataOutputType
 type of data writer (produces VTK data)
 
typedef Dune::Fem::CheckPointer< GridTypeCheckPointerType
 type of check-pointer (dumps unaltered simulation data)
 

Public Member Functions

 ParabolicFemScheme (DiscreteFunctionType &solution, const TimeProviderType &timeProvider, const ImplicitModelType &implicitModel, const ExplicitModelType &explicitModel, const InitialValueType &initialValue, const RHSFunctionalType &rhsFunctional, const std::string name="heat")
 Constructor. More...
 
virtual void initialize ()
 initialize the solution
 
virtual std::string name () const
 name of the Fem Scheme
 
virtual void next ()
 Close the current time-step. More...
 
virtual void solve (bool forceMatrixAssembling=true)
 solve the system
 
virtual bool mark (const double tolerance)
 mark elements for adaptation
 
virtual double estimate ()
 calculate error estimator
 
virtual double timeEstimate ()
 return the most recent time estimate
 
virtual void adapt ()
 do the adaptation for a given marking
 
virtual int output ()
 data I/O
 
virtual bool converged () const
 check whether solver has converged
 
virtual double residual () const
 calculate residual (in small l^2)
 
virtual double error () const
 Calculate L2/H1 error. More...
 
virtual size_t size () const
 return some measure about the number of DOFs in use
 

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 TimeProvider, class ImplicitModel, class ExplicitModel, class InitialValue, class RHSFunctional, template< class > class QuadratureTraits = DefaultQuadratureTraits>
class Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >

Basic parabolic fem-scheme class.

Parameters
[in]DiscreteFunctionThe type of the discrete solution.
[in]TimeProviderThe type of the Dune::Fem::TimeProvider.
[in]ImplicitModelThe type of the implicit model.
[in]ExplicitModelThe type of the explicit model.
[in]InitialValueThe type of the initial value. Additionally, this can also be the exact solution in order to perform experimental convergence tests. ParabolicFemScheme::error() will compute the distance to this function; it will also be dumped to the VTK-files for visualization.
[in]QuadratureTraitsSpecial quadrature rules, for instance to implement mass-lumping. Defaults to DefaultQuadratureTraits. The template argument of the template template-parameter is the GridPartType.

Member Typedef Documentation

◆ ExplicitOperatorType

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
typedef EllipticOperator<ExplicitModelType, DiscreteFunctionType, DiscreteFunctionType, EmptyConstraintsType, QuadratureTraits> Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::ExplicitOperatorType

explicit part of differential operator.

Need not be differentiable as it is only needed for the RHS.

◆ RestrictionProlongationType

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
typedef Dune::Fem::RestrictProlongDefault<DiscreteFunctionType> Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::RestrictionProlongationType

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

LagrangeInterpolation)

Constructor & Destructor Documentation

◆ ParabolicFemScheme()

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::ParabolicFemScheme ( DiscreteFunctionType solution,
const TimeProviderType timeProvider,
const ImplicitModelType implicitModel,
const ExplicitModelType &  explicitModel,
const InitialValueType initialValue,
const RHSFunctionalType &  rhsFunctional,
const std::string  name = "heat" 
)
inline

Constructor.

Parameters
[in,out]solutionStorage for the discrete solution.
[in]timeProviderA Dune::Fem time-provider to communicate time and step-size all over the place.
[in]implicitModelThe part of the operator to be applied at the end of the time interval.
[in]explicitModelThe part of the operator to be applied to the old values.
[in]initialValueGuess what. This can either be a non-discrete function or something which already has a LocalFunction-property. This can also used to pass an "exact" "solution" for experimental convergence tests. If so the function still has to evaluate at the start of the current time-step, in order to be suitable as an initial value.
[in]nameName used for parameters and output.

Member Function Documentation

◆ error()

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
virtual double Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::error ( ) const
inlinevirtual

Calculate L2/H1 error.

This actually computes the distance to the initial value which optionally may be used for experimental convergence tests in order to pass an exact "solution". In this case initialValue may be time-dependent and has to give the correct value at the start of the time-step.

Implements Dune::ACFem::BasicFemScheme.

◆ linearSolve()

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
virtual void Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::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::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::solve().

◆ next()

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits = DefaultQuadratureTraits>
virtual void Dune::ACFem::ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits >::next ( )
inlinevirtual

Close the current time-step.

Here: copy new solution to oldSolution. Note that this does not move the timeProvider to the new time-step.


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 (Dec 27, 23:30, 2024)