Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::LibUtilities::IMEXTimeIntegrationScheme Class Reference

#include <IMEXTimeIntegrationSchemes.h>

Inheritance diagram for Nektar::LibUtilities::IMEXTimeIntegrationScheme:
[legend]

Public Member Functions

 IMEXTimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~IMEXTimeIntegrationScheme () override
 
- Public Member Functions inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
LUE void InitializeSecondaryData (TimeIntegrationAlgorithmGLM *phase, NekDouble deltaT) const
 
- Public Member Functions inherited from Nektar::LibUtilities::TimeIntegrationScheme
LUE std::string GetFullName () const
 
LUE std::string GetName () const
 
LUE std::string GetVariant () const
 
LUE size_t GetOrder () const
 
LUE std::vector< NekDoubleGetFreeParams ()
 
LUE TimeIntegrationSchemeType GetIntegrationSchemeType ()
 
LUE NekDouble GetTimeStability () const
 
LUE size_t GetNumIntegrationPhases ()
 
LUE const TripleArrayGetSolutionVector () const
 Gets the solution vector of the ODE. More...
 
LUE TripleArrayUpdateSolutionVector ()
 
LUE void SetSolutionVector (const size_t Offset, const DoubleArray &y)
 Sets the solution vector of the ODE. More...
 
LUE void InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)
 Explicit integration of an ODE. More...
 
LUE ConstDoubleArrayTimeIntegrate (const size_t timestep, const NekDouble delta_t)
 
LUE void print (std::ostream &os) const
 
LUE void printFull (std::ostream &os) const
 

Static Public Member Functions

static TimeIntegrationSchemeSharedPtr create (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
static LUE void SetupSchemeData (TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

LUE std::string v_GetFullName () const override
 
LUE std::string v_GetName () const override
 
LUE NekDouble v_GetTimeStability () const override
 
- Protected Member Functions inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
LUE std::string v_GetVariant () const override
 
LUE size_t v_GetOrder () const override
 
LUE std::vector< NekDoublev_GetFreeParams () const override
 
LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType () const override
 
LUE size_t v_GetNumIntegrationPhases () const override
 
LUE const TripleArrayv_GetSolutionVector () const override
 
LUE TripleArrayv_UpdateSolutionVector () override
 
LUE void v_SetSolutionVector (const size_t Offset, const DoubleArray &y) override
 
LUE void v_InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
 Worker method to initialize the integration scheme. More...
 
LUE ConstDoubleArrayv_TimeIntegrate (const size_t timestep, const NekDouble delta_t) override
 Worker method that actually does the time integration. More...
 
virtual LUE void v_InitializeSecondaryData (TimeIntegrationAlgorithmGLM *phase, NekDouble deltaT) const
 
LUE void v_print (std::ostream &os) const override
 Worker method to print details on the integration scheme. More...
 
LUE void v_printFull (std::ostream &os) const override
 
LUE TimeIntegrationSchemeGLM (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~TimeIntegrationSchemeGLM () override
 
- Protected Member Functions inherited from Nektar::LibUtilities::TimeIntegrationScheme
virtual LUE std::string v_GetFullName () const
 
virtual LUE std::string v_GetName () const =0
 
virtual LUE std::string v_GetVariant () const =0
 
virtual LUE size_t v_GetOrder () const =0
 
virtual LUE std::vector< NekDoublev_GetFreeParams () const =0
 
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType () const =0
 
virtual LUE NekDouble v_GetTimeStability () const =0
 
virtual LUE size_t v_GetNumIntegrationPhases () const =0
 
virtual LUE const TripleArrayv_GetSolutionVector () const =0
 
virtual LUE TripleArrayv_UpdateSolutionVector ()=0
 
virtual LUE void v_SetSolutionVector (const size_t Offset, const DoubleArray &y)=0
 
virtual LUE void v_InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)=0
 
virtual LUE ConstDoubleArrayv_TimeIntegrate (const size_t timestep, const NekDouble delta_t)=0
 
virtual LUE void v_print (std::ostream &os) const =0
 
virtual LUE void v_printFull (std::ostream &os) const =0
 
LUE TimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
LUE TimeIntegrationScheme (const TimeIntegrationScheme &in)=delete
 
virtual ~TimeIntegrationScheme ()=default
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
TimeIntegrationAlgorithmGLMVector m_integration_phases
 
TimeIntegrationSolutionGLMSharedPtr m_solVector
 

Detailed Description

Definition at line 58 of file IMEXTimeIntegrationSchemes.h.

Constructor & Destructor Documentation

◆ IMEXTimeIntegrationScheme()

Nektar::LibUtilities::IMEXTimeIntegrationScheme::IMEXTimeIntegrationScheme ( std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inline

Definition at line 61 of file IMEXTimeIntegrationSchemes.h.

63 : TimeIntegrationSchemeGLM(variant, order, freeParams)
64 {
65 if (variant == "dirk" || variant == "DIRK")
66 {
67 ASSERTL1(freeParams.size() == 2,
68 "IMEX DIRK Time integration scheme invalid number "
69 "of free parameters, expected two "
70 "<implicit stages, explicit stages>, received " +
71 std::to_string(freeParams.size()));
72
73 size_t s = freeParams[0];
74 size_t sigma = freeParams[1];
75
78 new TimeIntegrationAlgorithmGLM(this));
79
80 if (order == 1 && s == 1 && sigma == 1)
81 {
82 // This phase is Forward Backward Euler which has two steps.
85 }
86 else
87 {
89 m_integration_phases[0], order, freeParams);
90 }
91 }
92 else if (variant == "Gear")
93 {
96 new TimeIntegrationAlgorithmGLM(this));
98 new TimeIntegrationAlgorithmGLM(this));
99
101 m_integration_phases[0], 2, std::vector<NekDouble>{2, 2});
104 }
105 else if (variant == "")
106 {
107 // Currently up to 4th order is implemented.
108 ASSERTL1(1 <= order && order <= 4,
109 "IMEX Time integration scheme bad order (1-4): " +
110 std::to_string(order));
111
113
114 for (size_t n = 0; n < order; ++n)
115 {
117 new TimeIntegrationAlgorithmGLM(this));
118 }
119
120 // Last phase
122 m_integration_phases[order - 1], order);
123
124 // Initial phases
125 switch (order)
126 {
127 case 1:
128 // No intial phase.
129 break;
130
131 case 2:
134 break;
135
136 case 3:
139 std::vector<NekDouble>{3, 4});
142 std::vector<NekDouble>{3, 4});
143 break;
144
145 case 4:
148 std::vector<NekDouble>{3, 4});
151 std::vector<NekDouble>{3, 4});
154 std::vector<NekDouble>{3, 4});
155 break;
156
157 default:
158 ASSERTL1(false, "IMEX Time integration scheme bad order: " +
159 std::to_string(order));
160 }
161 }
162 else
163 {
164 ASSERTL1(false, "IMEX Time integration scheme bad variant: " +
165 variant + ". Must be blank, 'dirk' or 'Gear'");
166 }
167 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order)
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order, std::vector< NekDouble > freeParams)
static LUE void SetupSchemeData_1_1_1(TimeIntegrationAlgorithmGLMSharedPtr &phase)
TimeIntegrationAlgorithmGLMVector m_integration_phases
LUE TimeIntegrationSchemeGLM(std::string variant, size_t order, std::vector< NekDouble > freeParams)
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
std::vector< TimeIntegrationAlgorithmGLMSharedPtr > TimeIntegrationAlgorithmGLMVector

References ASSERTL1, Nektar::LibUtilities::TimeIntegrationSchemeGLM::m_integration_phases, Nektar::LibUtilities::IMEXGearTimeIntegrationScheme::SetupSchemeData(), SetupSchemeData(), Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData(), and Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_1_1_1().

◆ ~IMEXTimeIntegrationScheme()

Nektar::LibUtilities::IMEXTimeIntegrationScheme::~IMEXTimeIntegrationScheme ( )
inlineoverride

Definition at line 169 of file IMEXTimeIntegrationSchemes.h.

170 {
171 }

Member Function Documentation

◆ create()

static TimeIntegrationSchemeSharedPtr Nektar::LibUtilities::IMEXTimeIntegrationScheme::create ( std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inlinestatic

Definition at line 173 of file IMEXTimeIntegrationSchemes.h.

175 {
178 variant, order, freeParams);
179
180 return p;
181 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ SetupSchemeData()

static LUE void Nektar::LibUtilities::IMEXTimeIntegrationScheme::SetupSchemeData ( TimeIntegrationAlgorithmGLMSharedPtr phase,
size_t  order 
)
inlinestatic

Definition at line 185 of file IMEXTimeIntegrationSchemes.h.

187 {
188 constexpr NekDouble ABcoefficients[5] = {0.,
189 1., // 1st Order
190 2. / 3., // 2nd Order
191 6. / 11., // 3rd Order
192 12. / 25.}; // 4th Order
193
194 // Nsteps = 2 * order
195
196 // clang-format off
197 constexpr NekDouble UVcoefficients[5][8] =
198 { { 0., 0., 0., 0.,
199 0., 0., 0., 0. },
200 // 1st Order
201 { 1., 1., 0., 0.,
202 0., 0., 0., 0. },
203 // 2nd Order
204 { 4./ 3., -1./ 3., 4./3., -2./ 3.,
205 0., 0., 0., 0. },
206 // 3rd Order
207 { 18./11., -9./11., 2./11., 18./11.,
208 -18./11., 6./11., 0., 0. },
209 // 4th Order
210 { 48./25., -36./25., 16./25., -3./25.,
211 48./25., -72./25., 48./25., -12./25. } };
212 // clang-format on
213
214 phase->m_schemeType = eIMEX;
215 phase->m_order = order;
216 phase->m_name =
217 std::string("IMEXOrder" + std::to_string(phase->m_order));
218
219 phase->m_numsteps = 2 * phase->m_order;
220 phase->m_numstages = 1;
221
222 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(2);
223 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(2);
224
225 phase->m_A[0] =
226 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
227 phase->m_B[0] =
228 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
229 phase->m_A[1] =
230 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
231 phase->m_B[1] =
232 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
233 phase->m_U =
234 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 0.0);
235 phase->m_V =
236 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 0.0);
237
238 // Coefficients
239 phase->m_B[1][phase->m_order][0] = 1.0;
240
241 phase->m_A[0][0][0] = ABcoefficients[phase->m_order];
242 phase->m_B[0][0][0] = ABcoefficients[phase->m_order];
243
244 for (size_t n = 0; n < 2 * phase->m_order; ++n)
245 {
246 phase->m_U[0][n] = UVcoefficients[phase->m_order][n];
247 phase->m_V[0][n] = UVcoefficients[phase->m_order][n];
248 }
249
250 // V evaluation value shuffling row n column n-1
251 for (size_t n = 1; n < 2 * phase->m_order; ++n)
252 {
253 if (n != phase->m_order)
254 {
255 phase->m_V[n][n - 1] = 1.0; // constant 1
256 }
257 }
258
259 phase->m_numMultiStepValues = phase->m_order;
260 phase->m_numMultiStepImplicitDerivs = 0;
261 phase->m_numMultiStepExplicitDerivs = phase->m_order;
262
263 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
264
265 // Values and derivatives needed.
266 for (size_t n = 0; n < phase->m_order; ++n)
267 {
268 phase->m_timeLevelOffset[n] = n;
269 phase->m_timeLevelOffset[phase->m_order + n] = n;
270 }
271
272 phase->CheckAndVerify();
273 }
@ eIMEX
Implicit Explicit General Linear Method.
double NekDouble

References Nektar::LibUtilities::eIMEX.

Referenced by IMEXTimeIntegrationScheme().

◆ v_GetFullName()

LUE std::string Nektar::LibUtilities::IMEXTimeIntegrationScheme::v_GetFullName ( ) const
inlineoverrideprotectedvirtual

◆ v_GetName()

LUE std::string Nektar::LibUtilities::IMEXTimeIntegrationScheme::v_GetName ( ) const
inlineoverrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 281 of file IMEXTimeIntegrationSchemes.h.

282 {
283 return std::string("IMEX");
284 }

◆ v_GetTimeStability()

LUE NekDouble Nektar::LibUtilities::IMEXTimeIntegrationScheme::v_GetTimeStability ( ) const
inlineoverrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 286 of file IMEXTimeIntegrationSchemes.h.

287 {
288 return 1.0;
289 }

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::IMEXTimeIntegrationScheme::className
static

Definition at line 183 of file IMEXTimeIntegrationSchemes.h.