35 #ifndef NEKTAR_LIB_UTILITIES_FOUNDATIONS_TIMEINTEGRATIONSCHEME_H
36 #define NEKTAR_LIB_UTILITIES_FOUNDATIONS_TIMEINTEGRATIONSCHEME_H
38 #include <boost/core/ignore_unused.hpp>
48 namespace LibUtilities
51 class TimeIntegrationScheme;
52 class TimeIntegrationSolution;
108 "NoTimeIntegrationMethod",
109 "AdamsBashforthOrder1",
110 "AdamsBashforthOrder2",
111 "AdamsBashforthOrder3",
112 "AdamsBashforthOrder4",
113 "AdamsMoultonOrder1",
114 "AdamsMoultonOrder2",
117 "ClassicalRungeKutta4",
121 "RungeKutta2_ImprovedEuler",
157 "NoTimeIntegrationSchemeType",
159 "DiagonallyImplicit",
190 template<
typename FuncPo
interT,
typename ObjectPo
interT>
194 func, obj, std::placeholders::_1, std::placeholders::_2,
195 std::placeholders::_3);
198 template<
typename FuncPo
interT,
typename ObjectPo
interT>
202 func, obj, std::placeholders::_1, std::placeholders::_2,
203 std::placeholders::_3);
206 template<
typename FuncPo
interT,
typename ObjectPo
interT>
210 func, obj, std::placeholders::_1, std::placeholders::_2,
211 std::placeholders::_3);
214 template<
typename FuncPo
interT,
typename ObjectPo
interT>
218 func, obj, std::placeholders::_1, std::placeholders::_2,
219 std::placeholders::_3);
222 template<
typename FuncPo
interT,
typename ObjectPo
interT>
226 func, obj, std::placeholders::_1, std::placeholders::_2,
227 std::placeholders::_3, std::placeholders::_4);
243 ASSERTL1(
m_functors1[1],
"OdeExplicitRhs should be defined for this time integration scheme");
251 ASSERTL1(
m_functors1[2],
"OdeImplictRhs should be defined for this time integration scheme");
259 ASSERTL1(
m_functors1[3],
"Projection operation should be defined for this time integration scheme");
268 ASSERTL1(
m_functors2[0],
"ImplicitSolve should be defined for this time integration scheme");
326 return (*
this == *y);
331 return (!(*
this == y));
336 return (!(*
this == *y));
412 inline NekDouble A(
const unsigned int i,
const unsigned int j)
const
417 inline NekDouble B(
const unsigned int i,
const unsigned int j)
const
422 inline NekDouble U(
const unsigned int i,
const unsigned int j)
const
427 inline NekDouble V(
const unsigned int i,
const unsigned int j)
const
584 boost::ignore_unused(in);
604 return y[0].num_elements();
609 return y[0][0].num_elements();
654 unsigned int npoints);
664 return m_scheme->GetIntegrationSchemeKey();
669 return m_scheme->GetIntegrationMethod();
726 return m_scheme->GetNmultiStepValues();
733 return m_scheme->GetNmultiStepDerivs();
740 return m_scheme->GetTimeLevelOffset();
747 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
750 for(
int i = 0; i < nMultiStepVals; i++)
752 if( timeLevelOffset == offsetvec[i] )
757 ASSERTL1(
false,
"The solution vector of this scheme does not contain a value at the requested time-level");
765 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
769 for(
int i = nMultiStepVals; i < size; i++)
771 if( timeLevelOffset == offsetvec[i] )
776 ASSERTL1(
false,
"The solution vector of this scheme does not contain a derivative at the requested time-level");
784 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
787 for(
int i = 0; i < nMultiStepVals; i++)
789 if( timeLevelOffset == offsetvec[i] )
794 ASSERTL1(
false,
"The solution vector of this scheme does not contain a value at the requested time-level");
803 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
806 for(
int i = 0; i < nMultiStepVals; i++)
808 if( timeLevelOffset == offsetvec[i] )
822 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
826 for(
int i = nMultiStepVals; i < size; i++)
828 if( timeLevelOffset == offsetvec[i] )
848 int nMultiStepVals =
m_scheme->GetNmultiStepValues();
850 for(i = (nMultiStepVals-1); i > 0; i--)
855 for(i = (size-1); i > nMultiStepVals; i--)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define LIB_UTILITIES_EXPORT
int m_npoints
The number of variables in integration scheme.
const Array< OneD, const NekDouble > ConstSingleArray
TripleArray m_F
explicit right hand side of each stage equation
std::function< void(ConstDoubleArray &, DoubleArray &, const NekDouble, const NekDouble) > FunctorType2
TimeIntegrationSchemeType m_schemeType
unsigned int GetNmultiStepDerivs(void) const
int GetFirstDim(ConstTripleArray &y) const
bool m_lastStageEqualsNewSolution
NekDouble A_IMEX(const unsigned int i, const unsigned int j) const
static std::shared_ptr< TimeIntegrationScheme > Create(const TimeIntegrationSchemeKey &key)
virtual ~TimeIntegrationScheme()
TimeIntegrationMethod GetIntegrationMethod() const
Array< OneD, unsigned int > m_timeLevelOffset
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > TripleArray
Array< TwoD, NekDouble > m_U
NekDouble B(const unsigned int i, const unsigned int j) const
bool CheckIfFirstStageEqualsOldSolution(const Array< OneD, const Array< TwoD, NekDouble > > &A, const Array< OneD, const Array< TwoD, NekDouble > > &B, const Array< TwoD, const NekDouble > &U, const Array< TwoD, const NekDouble > &V) const
TimeIntegrationSchemeType GetIntegrationSchemeType() const
bool CheckIfLastStageEqualsNewSolution(const Array< OneD, const Array< TwoD, NekDouble > > &A, const Array< OneD, const Array< TwoD, NekDouble > > &B, const Array< TwoD, const NekDouble > &U, const Array< TwoD, const NekDouble > &V) const
int GetSecondDim(ConstTripleArray &y) const
unsigned int m_numMultiStepValues
TripleArray m_F_IMEX
Array corresponding to the stage Derivatives.
TimeIntegrationScheme(const TimeIntegrationScheme &in)
unsigned int GetNstages(void) const
const TimeIntegrationSchemeKey & GetIntegrationSchemeKey() const
bool m_firstStageEqualsOldSolution
Array< TwoD, NekDouble > m_V
TimeIntegrationSchemeKey m_schemeKey
bool VerifyIntegrationSchemeType(TimeIntegrationSchemeType type, const Array< OneD, const Array< TwoD, NekDouble > > &A, const Array< OneD, const Array< TwoD, NekDouble > > &B, const Array< TwoD, const NekDouble > &U, const Array< TwoD, const NekDouble > &V) const
std::function< void(ConstDoubleArray &, DoubleArray &, const NekDouble) > FunctorType1
unsigned int GetNsteps(void) const
Array< OneD, Array< TwoD, NekDouble > > m_B
bool CheckTimeIntegrateArguments(const NekDouble timestep, ConstTripleArray &y_old, ConstSingleArray &t_old, TripleArray &y_new, SingleArray &t_new, const TimeIntegrationSchemeOperators &op) const
Array< OneD, Array< OneD, NekDouble > > DoubleArray
NekDouble B_IMEX(const unsigned int i, const unsigned int j) const
const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > ConstTripleArray
Array< OneD, NekDouble > SingleArray
const Array< OneD, const Array< OneD, NekDouble > > ConstDoubleArray
NekDouble m_T
Used to store the Explicit stage derivative of IMEX schemes.
NekDouble U(const unsigned int i, const unsigned int j) const
ConstDoubleArray & TimeIntegrate(const NekDouble timestep, TimeIntegrationSolutionSharedPtr &y, const TimeIntegrationSchemeOperators &op)
Explicit integration of an ODE.
Array< OneD, Array< TwoD, NekDouble > > m_A
friend TimeIntegrationSchemeManagerT & TimeIntegrationSchemeManager(void)
int m_nvar
bool to identify if array has been initialised
NekDouble A(const unsigned int i, const unsigned int j) const
DoubleArray m_Y
The size of inner data which is stored for reuse.
const Array< OneD, const unsigned int > & GetTimeLevelOffset()
DoubleArray m_tmp
Array containing the stage values.
unsigned int m_numMultiStepDerivs
unsigned int GetNmultiStepValues(void) const
NekDouble V(const unsigned int i, const unsigned int j) const
TimeIntegrationSolutionSharedPtr InitializeScheme(const NekDouble timestep, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)
This function initialises the time integration scheme.
bool operator==(const TimeIntegrationSchemeKey &key)
bool operator!=(const TimeIntegrationSchemeKey &y)
friend bool operator<(const TimeIntegrationSchemeKey &lhs, const TimeIntegrationSchemeKey &rhs)
TimeIntegrationMethod GetIntegrationMethod() const
TimeIntegrationSchemeKey()
TimeIntegrationMethod m_method
integration method
TimeIntegrationSchemeKey & operator=(const TimeIntegrationSchemeKey &key)
virtual ~TimeIntegrationSchemeKey()
TimeIntegrationSchemeKey(const TimeIntegrationSchemeKey &key)
TimeIntegrationSchemeKey(const TimeIntegrationMethod &method)
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineOdeExplicitRhs(FuncPointerT func, ObjectPointerT obj)
const FunctorType2 & ConstFunctorType2Ref
Array< OneD, FunctorType2 > FunctorType2Array
void DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) const
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void DoProjection(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DoOdeExplicitRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
const Array< OneD, const Array< OneD, NekDouble > > InArrayType
FunctorType1Array m_functors1
Array< OneD, FunctorType1 > FunctorType1Array
Array< OneD, Array< OneD, NekDouble > > OutArrayType
FunctorType2Array m_functors2
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DefineOdeImplicitRhs(FuncPointerT func, ObjectPointerT obj)
std::function< void(InArrayType &, OutArrayType &, const NekDouble) > FunctorType1
TimeIntegrationSchemeOperators(void)
const FunctorType1 & ConstFunctorType1Ref
void DoOdeImplicitRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
std::function< void(InArrayType &, OutArrayType &, const NekDouble, const NekDouble) > FunctorType2
DoubleArray & GetDerivative(const unsigned int timeLevelOffset)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > TripleArray
const DoubleArray & GetSolution() const
Array< OneD, NekDouble > m_t
Array< OneD, Array< OneD, NekDouble > > DoubleArray
TimeIntegrationSchemeSharedPtr m_scheme
void RotateSolutionVector()
const Array< OneD, const NekDouble > & GetTimeVector() const
TripleArray & UpdateSolutionVector()
const TimeIntegrationSchemeSharedPtr & GetIntegrationScheme() const
NekDouble GetTime() const
TimeIntegrationMethod GetIntegrationMethod() const
void SetSolVector(const int Offset, const DoubleArray &y)
NekDouble GetValueTime(const unsigned int timeLevelOffset)
unsigned int GetNderivs(void) const
DoubleArray & UpdateSolution()
Array< OneD, NekDouble > & UpdateTimeVector()
const TripleArray & GetSolutionVector() const
DoubleArray & GetValue(const unsigned int timeLevelOffset)
const Array< OneD, const unsigned int > & GetTimeLevelOffset()
unsigned int GetNvalues(void) const
TimeIntegrationSolution(const TimeIntegrationSchemeKey &key, const DoubleArray &y, const NekDouble time, const NekDouble timestep)
void SetDerivative(const unsigned int timeLevelOffset, const DoubleArray &y, const NekDouble timestep)
void SetValue(const unsigned int timeLevelOffset, const DoubleArray &y, const NekDouble t)
const TimeIntegrationSchemeKey & GetIntegrationSchemeKey() const
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
@ eRungeKutta4
Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
@ eIMEXdirk_1_2_1
Forward-Backward Euler IMEX DIRK(1,2,1)
@ eMCNAB
Modified Crank-Nicolson/Adams-Bashforth Order 2 (MCNAB)
@ eAdamsMoultonOrder2
Adams-Moulton Forward multi-step scheme of order 2.
@ eRungeKutta3_SSP
Nonlinear SSP RungeKutta3 explicit.
@ eCNAB
Crank-Nicolson/Adams-Bashforth Order 2 (CNAB)
@ eAdamsMoultonOrder1
Adams-Moulton Forward multi-step scheme of order 1.
@ SIZE_TimeIntegrationMethod
Length of enum list.
@ eIMEXGear
IMEX Gear Order 2.
@ eIMEXdirk_1_2_2
Implicit-Explicit Midpoint IMEX DIRK(1,2,2)
@ eIMEXdirk_4_4_3
L-stable, four stage, third order IMEX DIRK(4,4,3)
@ eNoTimeIntegrationMethod
@ eAdamsBashforthOrder2
Adams-Bashforth Forward multi-step scheme of order 2.
@ eMidpoint
midpoint method (old name)
@ eIMEXdirk_2_3_2
L-stable, three stage, third order IMEX DIRK(3,4,3)
@ eForwardEuler
Forward Euler scheme.
@ eBDFImplicitOrder2
BDF multi-step scheme of order 2 (implicit)
@ eAdamsBashforthOrder4
Adams-Bashforth Forward multi-step scheme of order 4.
@ eIMEXOrder4
IMEX 4th order scheme using Backward Different Formula & Extrapolation.
@ eIMEXOrder1
IMEX 1st order scheme using Euler Backwards/Euler Forwards.
@ eDIRKOrder3
Diagonally Implicit Runge Kutta scheme of order 3.
@ eBackwardEuler
Backward Euler scheme.
@ eRungeKutta5
RungeKutta5 method.
@ eDIRKOrder2
Diagonally Implicit Runge Kutta scheme of order 2.
@ eRungeKutta2_SSP
Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
@ eAdamsBashforthOrder3
Adams-Bashforth Forward multi-step scheme of order 3.
@ eRungeKutta2
Classical RungeKutta2 method (new name for eMidpoint)
@ eBDFImplicitOrder1
BDF multi-step scheme of order 1 (implicit)
@ eIMEXdirk_3_4_3
L-stable, three stage, third order IMEX DIRK(3,4,3)
@ eRungeKutta2_ImprovedEuler
Improved RungeKutta2 explicit (old name meaning Heun's method)
@ eIMEXOrder3
IMEX 3rd order scheme using Backward Different Formula & Extrapolation.
@ eClassicalRungeKutta4
Runge-Kutta multi-stage scheme 4th order explicit (old name)
@ eIMEXdirk_2_3_3
L-stable, two stage, third order IMEX DIRK(2,3,3)
@ eIMEXOrder2
IMEX 2nd order scheme using Backward Different Formula & Extrapolation.
@ eIMEXdirk_1_1_1
Forward-Backward Euler IMEX DIRK(1,1,1)
@ eAdamsBashforthOrder1
Adams-Bashforth Forward multi-step scheme of order 1.
@ eIMEXdirk_2_2_2
L-stable, two stage, second order IMEX DIRK(2,2,2)
std::vector< TimeIntegrationSolutionSharedPtr > TimeIntegrationSolutionVector
bool operator==(const BasisKey &x, const BasisKey &y)
const char *const TimeIntegrationSchemeTypeMap[]
const char *const TimeIntegrationMethodMap[]
std::vector< TimeIntegrationSchemeSharedPtr > TimeIntegrationSchemeVector
TimeIntegrationSchemeType
@ eImplicit
Fully implicit scheme.
@ eExplicit
Formally explicit scheme.
@ eDiagonallyImplicit
Diagonally implicit scheme (e.g. the DIRK schemes)
@ eIMEX
Implicit Explicit General Linear Method.
@ eNoTimeIntegrationSchemeType
bool operator<(const BasisKey &lhs, const BasisKey &rhs)
TimeIntegrationSchemeManagerT & TimeIntegrationSchemeManager(void)
static const TimeIntegrationSchemeKey NullTimeIntegrationSchemeKey(eNoTimeIntegrationMethod)
std::vector< TimeIntegrationSchemeSharedPtr >::iterator TimeIntegrationSchemeVectorIter
NekManager< TimeIntegrationSchemeKey, TimeIntegrationScheme, TimeIntegrationSchemeKey::opLess > TimeIntegrationSchemeManagerT
std::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
std::vector< TimeIntegrationSolutionSharedPtr >::iterator TimeIntegrationSolutionVectorIter
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
bool operator()(const TimeIntegrationSchemeKey &lhs, const TimeIntegrationSchemeKey &rhs) const