DUNE-ACFEM (unstable)

SchemeGenerators

Functions

template<class DiscreteFunction , class Model , class InitialGuess , class RHSFunctional , std::enable_if_t<(IsWrappableByConstLocalFunction< InitialGuess >::value &&IsLinearFunctional< RHSFunctional >::value), int > = 0>
auto Dune::ACFem::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. More...
 
template<class DiscreteFunction , class Model , class InitialGuess , std::enable_if_t< IsWrappableByConstLocalFunction< InitialGuess >::value, int > = 0>
auto Dune::ACFem::ellipticFemScheme (DiscreteFunction &solution, const Model &model, const InitialGuess &initialGuess, const std::string name="acfem.schemes.elliptic")
 Adaptive fem-scheme for "elliptic" problems. More...
 
template<class DiscreteFunction , class Model , class RHSFunctional , std::enable_if_t< IsLinearFunctional< RHSFunctional >::value, int > = 0>
auto Dune::ACFem::ellipticFemScheme (DiscreteFunction &solution, const Model &model, const RHSFunctional &rhsFunctional, const std::string name="acfem.schemes.elliptic")
 Adaptive fem-scheme for "elliptic" problems. More...
 
template<class DiscreteFunction , class Model >
auto Dune::ACFem::ellipticFemScheme (DiscreteFunction &solution, const Model &model, const std::string name="acfem.schemes.elliptic")
 Adaptive fem-scheme for "elliptic" problems. More...
 
template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits>
std::enable_if_t< IsLinearFunctional< RHSFunctional >::value, ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits > > Dune::ACFem::parabolicFemScheme (DiscreteFunction &solution, const TimeProvider &timeProvider, const ImplicitModel &implicitModel, const ExplicitModel &explicitModel, const InitialValue &initialValue, const RHSFunctional &rhsFunctional, const QuadratureTraits< typename DiscreteFunction::DiscreteFunctionSpaceType::GridPartType > &quadTraits, const std::string &name="acfem.schemes.parabolic")
 Basic parabolic fem-scheme class. More...
 
template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , template< class > class QuadratureTraits>
auto Dune::ACFem::parabolicFemScheme (DiscreteFunction &solution, const TimeProvider &timeProvider, const ImplicitModel &implicitModel, const ExplicitModel &explicitModel, const InitialValue &initialValue, const QuadratureTraits< typename DiscreteFunction::DiscreteFunctionSpaceType::GridPartType > &quadTraits, const std::string &name="acfem.schemes.parabolic")
 Basic parabolic fem-scheme class. More...
 
template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional >
std::enable_if_t< IsLinearFunctional< RHSFunctional >::value, ParabolicFemScheme< DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, DefaultQuadratureTraits > > Dune::ACFem::parabolicFemScheme (DiscreteFunction &solution, const TimeProvider &timeProvider, const ImplicitModel &implicitModel, const ExplicitModel &explicitModel, const InitialValue &initialValue, const RHSFunctional &rhsFunctional, const std::string &name="acfem.schemes.parabolic")
 Basic parabolic fem-scheme class. More...
 
template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue >
auto Dune::ACFem::parabolicFemScheme (DiscreteFunction &solution, const TimeProvider &timeProvider, const ImplicitModel &implicitModel, const ExplicitModel &explicitModel, const InitialValue &initialValue, const std::string &name="acfem.schemes.parabolic")
 Basic parabolic fem-scheme class. More...
 

Detailed Description

Function Documentation

◆ ellipticFemScheme() [1/4]

template<class DiscreteFunction , class Model , class InitialGuess , class RHSFunctional , std::enable_if_t<(IsWrappableByConstLocalFunction< InitialGuess >::value &&IsLinearFunctional< RHSFunctional >::value), int > = 0>
auto Dune::ACFem::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.

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.

Referenced by Dune::ACFem::l2Projection().

◆ ellipticFemScheme() [2/4]

template<class DiscreteFunction , class Model , class InitialGuess , std::enable_if_t< IsWrappableByConstLocalFunction< InitialGuess >::value, int > = 0>
auto Dune::ACFem::ellipticFemScheme ( DiscreteFunction &  solution,
const Model &  model,
const InitialGuess &  initialGuess,
const std::string  name = "acfem.schemes.elliptic" 
)

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.

◆ ellipticFemScheme() [3/4]

template<class DiscreteFunction , class Model , class RHSFunctional , std::enable_if_t< IsLinearFunctional< RHSFunctional >::value, int > = 0>
auto Dune::ACFem::ellipticFemScheme ( DiscreteFunction &  solution,
const Model &  model,
const RHSFunctional &  rhsFunctional,
const std::string  name = "acfem.schemes.elliptic" 
)

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.

◆ ellipticFemScheme() [4/4]

template<class DiscreteFunction , class Model >
auto Dune::ACFem::ellipticFemScheme ( DiscreteFunction &  solution,
const Model &  model,
const std::string  name = "acfem.schemes.elliptic" 
)

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.

◆ parabolicFemScheme() [1/4]

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , template< class > class QuadratureTraits>
auto Dune::ACFem::parabolicFemScheme ( DiscreteFunction &  solution,
const TimeProvider &  timeProvider,
const ImplicitModel &  implicitModel,
const ExplicitModel &  explicitModel,
const InitialValue &  initialValue,
const QuadratureTraits< typename DiscreteFunction::DiscreteFunctionSpaceType::GridPartType > &  quadTraits,
const std::string &  name = "acfem.schemes.parabolic" 
)

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.

◆ parabolicFemScheme() [2/4]

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional , template< class > class QuadratureTraits>
std::enable_if_t<IsLinearFunctional<RHSFunctional>::value, ParabolicFemScheme<DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, QuadratureTraits> > Dune::ACFem::parabolicFemScheme ( DiscreteFunction &  solution,
const TimeProvider &  timeProvider,
const ImplicitModel &  implicitModel,
const ExplicitModel &  explicitModel,
const InitialValue &  initialValue,
const RHSFunctional &  rhsFunctional,
const QuadratureTraits< typename DiscreteFunction::DiscreteFunctionSpaceType::GridPartType > &  quadTraits,
const std::string &  name = "acfem.schemes.parabolic" 
)

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.

◆ parabolicFemScheme() [3/4]

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue , class RHSFunctional >
std::enable_if_t<IsLinearFunctional<RHSFunctional>::value, ParabolicFemScheme<DiscreteFunction, TimeProvider, ImplicitModel, ExplicitModel, InitialValue, RHSFunctional, DefaultQuadratureTraits> > Dune::ACFem::parabolicFemScheme ( DiscreteFunction &  solution,
const TimeProvider &  timeProvider,
const ImplicitModel &  implicitModel,
const ExplicitModel &  explicitModel,
const InitialValue &  initialValue,
const RHSFunctional &  rhsFunctional,
const std::string &  name = "acfem.schemes.parabolic" 
)

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.

◆ parabolicFemScheme() [4/4]

template<class DiscreteFunction , class TimeProvider , class ImplicitModel , class ExplicitModel , class InitialValue >
auto Dune::ACFem::parabolicFemScheme ( DiscreteFunction &  solution,
const TimeProvider &  timeProvider,
const ImplicitModel &  implicitModel,
const ExplicitModel &  explicitModel,
const InitialValue &  initialValue,
const std::string &  name = "acfem.schemes.parabolic" 
)

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.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 25, 22:37, 2024)