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

#include <BDFImplicitTimeIntegrationSchemes.h>

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

Public Member Functions

 BDFImplicitTimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~BDFImplicitTimeIntegrationScheme () 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_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 57 of file BDFImplicitTimeIntegrationSchemes.h.

Constructor & Destructor Documentation

◆ BDFImplicitTimeIntegrationScheme()

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

Definition at line 60 of file BDFImplicitTimeIntegrationSchemes.h.

62 : TimeIntegrationSchemeGLM(variant, order, freeParams)
63 {
64 // Currently up to 4th order is implemented.
65 ASSERTL1(1 <= order && order <= 4,
66 "BDFImplicit Time integration scheme bad order (1-4): " +
67 std::to_string(order));
68
70
71 for (size_t n = 0; n < order; ++n)
72 {
74 new TimeIntegrationAlgorithmGLM(this));
75 }
76
77 // Next to last phase
78 if (order > 1)
79 {
81 m_integration_phases[order - 2], order - 1);
82 }
83
84 // Last phase
86 m_integration_phases[order - 1], order);
87
88 // Initial phases
89 switch (order)
90 {
91 case 1:
92 // No intial phase.
93 break;
94
95 case 2:
96 // Intial phase set above
97 break;
98
99 case 3:
100 // Order 2
103 break;
104
105 case 4:
106 // Order 3
109 // Order 3
112 break;
113
114 default:
115 ASSERTL1(false,
116 "BDFImplicit Time integration scheme bad order: " +
117 std::to_string(order));
118 }
119 }
#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, size_t order)
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order)
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, SetupSchemeData(), and Nektar::LibUtilities::DIRKTimeIntegrationScheme::SetupSchemeData().

◆ ~BDFImplicitTimeIntegrationScheme()

Nektar::LibUtilities::BDFImplicitTimeIntegrationScheme::~BDFImplicitTimeIntegrationScheme ( )
inlineoverride

Definition at line 121 of file BDFImplicitTimeIntegrationSchemes.h.

122 {
123 }

Member Function Documentation

◆ create()

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

Definition at line 125 of file BDFImplicitTimeIntegrationSchemes.h.

127 {
130 variant, order, freeParams);
131
132 return p;
133 }
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::BDFImplicitTimeIntegrationScheme::SetupSchemeData ( TimeIntegrationAlgorithmGLMSharedPtr phase,
size_t  order 
)
inlinestatic

Definition at line 137 of file BDFImplicitTimeIntegrationSchemes.h.

139 {
140 constexpr NekDouble ABcoefficients[5] = {0.,
141 1., // 1st Order
142 2. / 3., // 2nd Order
143 6. / 11., // 3rd Order
144 12. / 25.}; // 4th Order
145
146 // clang-format off
147 constexpr NekDouble UVcoefficients[5][4] =
148 { { 0., 0., 0., 0. },
149 // 1st Order
150 { 1., 0., 0., 0. },
151 // 2nd Order
152 { 4./ 3., -1./ 3., 0., 0. },
153 // 3rd Order
154 { 18./11., -9./11., 2./11., 0. },
155 // 4th Order
156 { 48./25., -36./25., 16./25., -3./25. } };
157 // clang-format on
158
159 phase->m_schemeType = eDiagonallyImplicit;
160 phase->m_order = order;
161 phase->m_name =
162 std::string("BDFImplicitOrder" + std::to_string(phase->m_order));
163
164 phase->m_numsteps = phase->m_order;
165 phase->m_numstages = 1;
166
167 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
168 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
169
170 phase->m_A[0] =
171 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
172 phase->m_B[0] =
173 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
174 phase->m_U =
175 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 0.0);
176 phase->m_V =
177 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 0.0);
178
179 // Coefficients
180
181 // A/B Coefficient for first row first column
182 phase->m_A[0][0][0] = ABcoefficients[phase->m_order];
183 phase->m_B[0][0][0] = ABcoefficients[phase->m_order];
184
185 // U/V Coefficients for first row additional columns
186 for (size_t n = 0; n < phase->m_order; ++n)
187 {
188 phase->m_U[0][n] = UVcoefficients[phase->m_order][n];
189 phase->m_V[0][n] = UVcoefficients[phase->m_order][n];
190 }
191
192 // V evaluation value shuffling row n column n-1
193 for (size_t n = 1; n < phase->m_order; ++n)
194 {
195 phase->m_V[n][n - 1] = 1.0;
196 }
197
198 phase->m_numMultiStepValues = phase->m_order;
199 phase->m_numMultiStepImplicitDerivs = 0;
200 phase->m_numMultiStepExplicitDerivs = 0;
201 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
202 for (size_t n = 0; n < phase->m_order; ++n)
203 {
204 phase->m_timeLevelOffset[n] = n;
205 }
206
207 phase->CheckAndVerify();
208 }
@ eDiagonallyImplicit
Diagonally implicit scheme (e.g. the DIRK schemes)
double NekDouble

References Nektar::LibUtilities::eDiagonallyImplicit.

Referenced by BDFImplicitTimeIntegrationScheme().

◆ v_GetName()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 211 of file BDFImplicitTimeIntegrationSchemes.h.

212 {
213 return std::string("BDFImplicit");
214 }

◆ v_GetTimeStability()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 216 of file BDFImplicitTimeIntegrationSchemes.h.

217 {
218 return 1.0;
219 }

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::BDFImplicitTimeIntegrationScheme::className
static

Definition at line 135 of file BDFImplicitTimeIntegrationSchemes.h.