KASKADE 7 development version
|
Extrapolated linearly implicit Euler method. More...
#include <limex.hh>
Extrapolated linearly implicit Euler method.
This class implements the extrapolated linearly implicit Euler method for integrating time-dependent evolution problems. The implementation follows Deuflhard/Bornemann Chapter 6.4.3.
This class implements the extrapolated linearly implicit Euler method for integrating time-dependent evolution problems. The implementation follows Deuflhard/Bornemann Chapter 6.4.3.
Eq | a parabolic equation, a model of the ParabolicEquation interface. |
This class implements the extrapolated linearly implicit Euler method for integrating time-dependent evolution problems. The implementation follows Deuflhard/Bornemann Chapter 6.4.3.
Public Types | |
typedef Eq | EvolutionEquation |
typedef EvolutionEquation::AnsatzVars::VariableSet | State |
typedef Eq | EvolutionEquation |
typedef EvolutionEquation::AnsatzVars::VariableSet | State |
typedef Eq | EvolutionEquation |
typedef EvolutionEquation::AnsatzVars::VariableSet | State |
Public Member Functions | |
Limex (GridManager< typename EvolutionEquation::AnsatzVars::Grid > &gridManager, EvolutionEquation &eq_, typename EvolutionEquation::AnsatzVars const &ansatzVars_, std::vector< std::pair< double, double > > const &tolX) | |
State const & | step (State const &x, double dt, int order) |
double | estimateError (State const &x, int i, int j) const |
template<class OutStream > | |
void | reportTime (OutStream &out) const |
void | advanceTime (double dt) |
Limex (GridManager< typename EvolutionEquation::AnsatzVars::Grid > &gridManager_, EvolutionEquation &eq_, typename EvolutionEquation::AnsatzVars const &ansatzVars_, DirectType st=DirectType::MUMPS, PrecondType precondType_=PrecondType::ILUK, int verbosity_=1) | |
State const & | step (State const &x, double dt, int order, std::vector< std::pair< double, double > > const &tolX) |
Computes a state increment that advances the given state in time. The time in the given evolution equation is increased by dt. More... | |
std::vector< std::pair< double, double > > | estimateError (State const &x, int i, int j) const |
template<class OutStream > | |
void | reportTime (OutStream &out) const |
void | advanceTime (double dt) |
Limex (GridManager< typename EvolutionEquation::AnsatzVars::Grid > &gridManager_, EvolutionEquation &eq_, typename EvolutionEquation::AnsatzVars const &ansatzVars_, DirectType st=DirectType::UMFPACK) | |
State const & | step (State const &x, double dt, int order, std::vector< std::pair< double, double > > const &tolX) |
State const & | step (State const &x, double dt, int order, std::vector< std::pair< double, double > > const &tolX, std::vector< std::vector< bool > > &refinements) |
std::vector< std::pair< double, double > > | estimateError (State const &x, int i, int j) const |
template<class OutStream > | |
void | reportTime (OutStream &out) const |
void | advanceTime (double dt) |
Public Attributes | |
ExtrapolationTableau< State > | extrap |
double | rhsAssemblyTime |
double | matrixAssemblyTime |
double | factorizationTime |
double | solutionTime |
double | precAssemblyTime |
double | elapsedTimeSinceReset |
double | initSolutionTime |
double | estimateTime |
DirectType | directType |
PrecondType | precondType |
int | verbosity |
double | adaptivityTime |
DirectType | solverType |
typedef Eq Kaskade::Limex< Eq >::EvolutionEquation |
typedef Eq Kaskade::Limex< Eq >::EvolutionEquation |
Definition at line 58 of file limexWithoutJens.hh.
typedef Eq Kaskade::Limex< Eq >::EvolutionEquation |
Definition at line 73 of file limexWithoutJensHierarchicEst.hh.
typedef EvolutionEquation::AnsatzVars::VariableSet Kaskade::Limex< Eq >::State |
typedef EvolutionEquation::AnsatzVars::VariableSet Kaskade::Limex< Eq >::State |
Definition at line 59 of file limexWithoutJens.hh.
typedef EvolutionEquation::AnsatzVars::VariableSet Kaskade::Limex< Eq >::State |
Definition at line 74 of file limexWithoutJensHierarchicEst.hh.
|
inline |
|
inline |
Constructs an ODE integrator. The arguments eq and ansatzVars have to exist during the lifetime of the integrator.
Definition at line 70 of file limexWithoutJens.hh.
|
inline |
Constructs an ODE integrator. The arguments eq and ansatzVars have to exist during the lifetime of the integrator.
Definition at line 106 of file limexWithoutJensHierarchicEst.hh.
|
inline |
|
inline |
Definition at line 424 of file limexWithoutJens.hh.
|
inline |
Definition at line 705 of file limexWithoutJensHierarchicEst.hh.
|
inline |
|
inline |
Estimates the time discretization error of the previously computed step by taking the difference between the diagonal and subdiagonal extrapolation values of maximal order. This requires that order>1 has been given for the last step.
Definition at line 399 of file limexWithoutJens.hh.
|
inline |
Estimates the time discretization error of the previously computed step by taking the difference between the diagonal and subdiagonal extrapolation values of maximal order. This requires that order>1 has been given for the last step.
Definition at line 683 of file limexWithoutJensHierarchicEst.hh.
|
inline |
|
inline |
Definition at line 413 of file limexWithoutJens.hh.
|
inline |
Definition at line 697 of file limexWithoutJensHierarchicEst.hh.
|
inline |
Computes a state increment that advances the given state in time. The time in the given evolution equation is increased by dt.
x | the initial state to be evolved |
dt | the time step |
order | the extrapolation order |
Definition at line 66 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::step().
|
inline |
Computes a state increment that advances the given state in time. The time in the given evolution equation is increased by dt.
x | the initial state to be evolved |
dt | the time step |
order | the extrapolation order >= 0 (0 corresponds to the linearly implicit Euler) |
tolX |
Definition at line 98 of file limexWithoutJens.hh.
|
inline |
In order to maintain compatibility with existing code, this overload is needed, as non const reference parameters (refinements) can not have default values. All ist does is call the original method with a temporary parameter.
Definition at line 122 of file limexWithoutJensHierarchicEst.hh.
|
inline |
Computes a state increment that advances the given state in time. The time in the given evolution equation is increased by dt.
x | the initial state to be evolved |
dt | the time step |
order | the extrapolation order >= 0 (0 corresponds to the linearly implicit Euler) |
tolX | |
refinements | keep track of the cells marked for refinement |
Definition at line 146 of file limexWithoutJensHierarchicEst.hh.
double Kaskade::Limex< Eq >::adaptivityTime |
Definition at line 718 of file limexWithoutJensHierarchicEst.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
DirectType Kaskade::Limex< Eq >::directType |
Definition at line 440 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::elapsedTimeSinceReset |
Definition at line 439 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::estimateTime |
Definition at line 439 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
ExtrapolationTableau< State > Kaskade::Limex< Eq >::extrap |
Definition at line 171 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::estimateError(), and Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::factorizationTime |
Definition at line 172 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::initSolutionTime |
Definition at line 439 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::matrixAssemblyTime |
Definition at line 172 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::precAssemblyTime |
Definition at line 438 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
PrecondType Kaskade::Limex< Eq >::precondType |
Definition at line 441 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::rhsAssemblyTime |
Definition at line 172 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
double Kaskade::Limex< Eq >::solutionTime |
Definition at line 172 of file limex.hh.
Referenced by Kaskade::Limex< Eq >::reportTime(), and Kaskade::Limex< Eq >::step().
DirectType Kaskade::Limex< Eq >::solverType |
Definition at line 719 of file limexWithoutJensHierarchicEst.hh.
Referenced by Kaskade::Limex< Eq >::step().
int Kaskade::Limex< Eq >::verbosity |
Definition at line 442 of file limexWithoutJens.hh.
Referenced by Kaskade::Limex< Eq >::step().