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

Class for spectral deferred correction integration. More...

#include <TimeIntegrationSchemeSDC.h>

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

Public Member Functions

 TimeIntegrationSchemeSDC (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~TimeIntegrationSchemeSDC () override
 Destructor. More...
 
LUE void SetPFASST (bool pfasst)
 
LUE void SetTime (double time)
 
LUE size_t GetMaxOrder () const
 
LUE bool HasFirstQuadrature () const
 
LUE bool HasLastQuadrature () const
 
LUE size_t GetQuadPtsNumber () const
 
LUE size_t GetNpoints () const
 
LUE size_t GetNvars () const
 
LUE const PointsKeyGetPointsKey () const
 
LUE const DoubleArrayGetFirstQuadratureSolutionVector () const
 
LUE DoubleArrayUpdateFirstQuadratureSolutionVector ()
 
LUE const DoubleArrayGetLastQuadratureSolutionVector () const
 
LUE DoubleArrayUpdateLastQuadratureSolutionVector ()
 
LUE const TripleArrayGetResidualVector () const
 
LUE TripleArrayUpdateResidualVector ()
 
LUE const TripleArrayGetIntegratedResidualVector () const
 
LUE TripleArrayUpdateIntegratedResidualVector ()
 
LUE const TripleArrayGetFAScorrectionVector () const
 
LUE TripleArrayUpdateFAScorrectionVector ()
 
LUE void ResidualEval (const NekDouble &delta_t, const size_t n)
 
LUE void ResidualEval (const NekDouble &delta_t)
 
LUE void ComputeInitialGuess (const NekDouble &delta_t)
 
LUE void SDCIterationLoop (const NekDouble &delta_t)
 
LUE void UpdateFirstQuadrature (void)
 Worker method that update the first quadrature. More...
 
LUE void UpdateLastQuadrature (void)
 Worker method that update the last quadrature. More...
 
LUE void AddFASCorrectionToSFint (void)
 Worker method that add the FASCorrection. More...
 
LUE void UpdateIntegratedResidualSFint (const NekDouble &delta_t)
 Worker method that compute residual integral. More...
 
LUE void UpdateIntegratedResidualQFint (const NekDouble &delta_t)
 
- 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 Public Attributes

static std::string className
 

Protected Member Functions

LUE std::string v_GetName () const override
 
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 NekDouble v_GetTimeStability () const override
 
LUE size_t v_GetNumIntegrationPhases () const override
 
LUE const TripleArrayv_GetSolutionVector () const override
 Gets the solution vector of the ODE. More...
 
LUE TripleArrayv_UpdateSolutionVector () override
 
LUE void v_SetSolutionVector (const size_t Offset, const DoubleArray &y) override
 Sets the solution vector of the ODE. More...
 
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 performs the time integration. More...
 
virtual LUE void v_ResidualEval (const NekDouble &delta_t, const size_t n)
 
virtual LUE void v_ResidualEval (const NekDouble &delta_t)
 
virtual LUE void v_ComputeInitialGuess (const NekDouble &delta_t)
 
virtual LUE void v_SDCIterationLoop (const NekDouble &delta_t)
 
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
 
- 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
 

Protected Attributes

NekDouble m_time
 
std::string m_name
 
std::string m_variant
 
std::vector< NekDoublem_freeParams
 
TimeIntegrationSchemeType m_schemeType {eNoTimeIntegrationSchemeType}
 
TimeIntegrationSchemeOperators m_op
 
PointsKey m_pointsKey
 
SingleArray m_tau
 Object containing quadrature data. More...
 
DoubleArray m_Y_f
 Array containing the quadrature points. More...
 
TripleArray m_Y
 Array containing the last stage values. More...
 
TripleArray m_F
 Array containing the stage values. More...
 
TripleArray m_FAScorr
 Array containing the stage derivatives. More...
 
TripleArray m_SFint
 Array containing the FAS correction term. More...
 
TripleArray m_QFint
 Array containing the integrated residual term. More...
 
SingleArray m_QMat
 Array containing the integrated residual term. More...
 
SingleArray m_interp
 Array containing the integration matrix. More...
 
NekDouble m_theta {1.0}
 Array containing the interpolation coefficients. More...
 
size_t m_ordermin {0}
 SDC parameter. More...
 
size_t m_ordermax {0}
 Minimum order of the integration scheme. More...
 
size_t m_order {0}
 Maximum order of the integration scheme. More...
 
size_t m_nQuadPts {0}
 Order of the integration scheme. More...
 
size_t m_nvars {0}
 Number of quadrature points. More...
 
size_t m_npoints {0}
 Number of variables in the integration scheme. More...
 
bool m_first_quadrature {true}
 Number of points in the integration scheme. More...
 
bool m_last_quadrature {true}
 
bool m_initialized {false}
 
bool m_PFASST {false}
 

Detailed Description

Class for spectral deferred correction integration.

Definition at line 54 of file TimeIntegrationSchemeSDC.h.

Constructor & Destructor Documentation

◆ TimeIntegrationSchemeSDC()

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

Definition at line 57 of file TimeIntegrationSchemeSDC.h.

59 : TimeIntegrationScheme(variant, order, freeParams),
60 m_name("SpectralDeferredCorrection")
61 {
62 ASSERTL0(freeParams.size() == 2,
63 "SDC Time integration scheme invalid number "
64 "of free parameters, expected two "
65 "<theta, number of quadrature>, received " +
66 std::to_string(freeParams.size()));
67
68 ASSERTL0(0.0 <= freeParams[0] && freeParams[0] <= 1.0,
69 "Spectral Deferred Correction Time integration "
70 "scheme bad parameter numbers (0.0 - 1.0): " +
71 std::to_string(freeParams[0]));
72
73 // Set scheme properties
74 m_variant = variant;
75 m_order = order;
76 m_freeParams = freeParams;
79
80 if (m_variant == "Equidistant")
81 {
83 m_variant +
84 " quadrature require quadrature "
85 "numbers (>=1" +
86 "): " + std::to_string(freeParams[1]));
87
88 m_first_quadrature = (m_nQuadPts == 1) ? false : true;
89 m_last_quadrature = (m_nQuadPts == 1) ? false : true;
90 m_ordermin = 1;
92 m_pointsKey = LibUtilities::PointsKey(
94 }
95 else if (m_variant == "GaussLobattoLegendre")
96 {
98 m_variant +
99 " quadrature require quadrature "
100 "numbers (>=2" +
101 "): " + std::to_string(freeParams[1]));
102
103 m_first_quadrature = true;
104 m_last_quadrature = true;
105 m_ordermin = 1;
106 m_ordermax = 2 * m_nQuadPts - 2;
107 m_pointsKey = LibUtilities::PointsKey(
109 }
110 else if (m_variant == "GaussRadauLegendre")
111 {
112 ASSERTL0(2 <= m_nQuadPts,
113 m_variant +
114 " quadrature require quadrature "
115 "numbers (>=1" +
116 "): " + std::to_string(freeParams[1]));
117
118 m_first_quadrature = false;
119 m_last_quadrature = true;
120 m_ordermin = 1;
121 m_ordermax = 2 * m_nQuadPts - 1;
122 m_pointsKey = LibUtilities::PointsKey(
124 }
125 else if (m_variant == "GaussGaussLegendre")
126 {
127 ASSERTL0(1 <= m_nQuadPts,
128 m_variant +
129 " quadrature require quadrature "
130 "numbers (>=1" +
131 "): " + std::to_string(freeParams[1]));
132
133 m_first_quadrature = false;
134 m_last_quadrature = false;
135 m_ordermin = 1;
137 m_pointsKey = LibUtilities::PointsKey(
139 }
140 else
141 {
142 ASSERTL0(false, "unknow variant (quadrature) type");
143 }
144
146 "Spectral Deferred Correction Time integration "
147 "scheme bad order numbers (>=" +
148 std::to_string(m_ordermin) +
149 "): " + std::to_string(m_order));
150
152 "Spectral Deferred Correction Time integration "
153 "scheme bad order numbers (<=" +
154 std::to_string(m_ordermax) +
155 "): " + std::to_string(m_order));
156
157 // Add one extra quadrature points for i.c., if necessary
159 {
161 }
162
163 // Get quadrature points and rescale to [0, 1]
165 size_t offset = m_first_quadrature ? 0 : 1;
166 for (size_t i = offset; i < m_nQuadPts; i++)
167 {
168 NekDouble tau = PointsManager()[m_pointsKey]->GetZ()[i - offset];
169 m_tau[i] = (tau + 1.0) / 2.0;
170 }
171 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LUE TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
size_t m_ordermax
Minimum order of the integration scheme.
size_t m_nQuadPts
Order of the integration scheme.
bool m_first_quadrature
Number of points in the integration scheme.
NekDouble m_theta
Array containing the interpolation coefficients.
size_t m_order
Maximum order of the integration scheme.
SingleArray m_tau
Object containing quadrature data.
PointsManagerT & PointsManager(void)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:49
double NekDouble

References ASSERTL0, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauPLegendre, Nektar::LibUtilities::ePolyEvenlySpaced, m_first_quadrature, m_freeParams, m_last_quadrature, m_nQuadPts, m_order, m_ordermax, m_ordermin, m_pointsKey, m_tau, m_theta, m_variant, and Nektar::LibUtilities::PointsManager().

◆ ~TimeIntegrationSchemeSDC()

Nektar::LibUtilities::TimeIntegrationSchemeSDC::~TimeIntegrationSchemeSDC ( )
inlineoverride

Destructor.

Definition at line 174 of file TimeIntegrationSchemeSDC.h.

175 {
176 }

Member Function Documentation

◆ AddFASCorrectionToSFint()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::AddFASCorrectionToSFint ( void  )

Worker method that add the FASCorrection.

Definition at line 248 of file TimeIntegrationSchemeSDC.cpp.

249{
250 // Add PFASST correction term
251 if (m_PFASST)
252 {
253 for (size_t n = 1; n < m_nQuadPts; ++n)
254 {
255 for (size_t i = 0; i < m_nvars; ++i)
256 {
257 Vmath::Vadd(m_npoints, m_SFint[n][i], 1, m_FAScorr[n][i], 1,
258 m_SFint[n][i], 1);
259 Vmath::Vsub(m_npoints, m_SFint[n][i], 1, m_FAScorr[n - 1][i], 1,
260 m_SFint[n][i], 1);
261 }
262 }
263 }
264}
size_t m_npoints
Number of variables in the integration scheme.
TripleArray m_FAScorr
Array containing the stage derivatives.
TripleArray m_SFint
Array containing the FAS correction term.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References m_FAScorr, m_npoints, m_nQuadPts, m_nvars, m_PFASST, m_SFint, Vmath::Vadd(), and Vmath::Vsub().

Referenced by Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC::v_SDCIterationLoop(), Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC::v_SDCIterationLoop(), and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC::v_SDCIterationLoop().

◆ ComputeInitialGuess()

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::ComputeInitialGuess ( const NekDouble delta_t)
inline

Definition at line 295 of file TimeIntegrationSchemeSDC.h.

296 {
297 v_ComputeInitialGuess(delta_t);
298 }
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t)

References v_ComputeInitialGuess().

Referenced by v_TimeIntegrate().

◆ create()

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

Definition at line 178 of file TimeIntegrationSchemeSDC.h.

180 {
183 variant, order, freeParams);
184
185 return p;
186 }
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.

◆ GetFAScorrectionVector()

LUE const TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetFAScorrectionVector ( ) const
inline

Definition at line 275 of file TimeIntegrationSchemeSDC.h.

276 {
277 return m_FAScorr;
278 }

References m_FAScorr.

◆ GetFirstQuadratureSolutionVector()

LUE const DoubleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetFirstQuadratureSolutionVector ( ) const
inline

Definition at line 235 of file TimeIntegrationSchemeSDC.h.

236 {
237 return m_Y[0];
238 }
TripleArray m_Y
Array containing the last stage values.

References m_Y.

◆ GetIntegratedResidualVector()

LUE const TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetIntegratedResidualVector ( ) const
inline

Definition at line 265 of file TimeIntegrationSchemeSDC.h.

266 {
267 return m_QFint;
268 }
TripleArray m_QFint
Array containing the integrated residual term.

References m_QFint.

◆ GetLastQuadratureSolutionVector()

LUE const DoubleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetLastQuadratureSolutionVector ( ) const
inline

Definition at line 245 of file TimeIntegrationSchemeSDC.h.

246 {
247 return m_last_quadrature ? m_Y[m_nQuadPts - 1] : m_Y_f;
248 }
DoubleArray m_Y_f
Array containing the quadrature points.

References m_last_quadrature, m_nQuadPts, m_Y, and m_Y_f.

◆ GetMaxOrder()

LUE size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetMaxOrder ( ) const
inline

Definition at line 200 of file TimeIntegrationSchemeSDC.h.

201 {
202 return m_ordermax;
203 }

References m_ordermax.

◆ GetNpoints()

LUE size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetNpoints ( ) const
inline

Definition at line 220 of file TimeIntegrationSchemeSDC.h.

221 {
222 return m_npoints;
223 }

References m_npoints.

◆ GetNvars()

LUE size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetNvars ( ) const
inline

Definition at line 225 of file TimeIntegrationSchemeSDC.h.

226 {
227 return m_nvars;
228 }

References m_nvars.

◆ GetPointsKey()

LUE const PointsKey & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetPointsKey ( ) const
inline

Definition at line 230 of file TimeIntegrationSchemeSDC.h.

231 {
232 return m_pointsKey;
233 }

References m_pointsKey.

◆ GetQuadPtsNumber()

LUE size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetQuadPtsNumber ( ) const
inline

Definition at line 215 of file TimeIntegrationSchemeSDC.h.

216 {
217 return m_nQuadPts;
218 }

References m_nQuadPts.

◆ GetResidualVector()

LUE const TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::GetResidualVector ( ) const
inline

Definition at line 255 of file TimeIntegrationSchemeSDC.h.

256 {
257 return m_F;
258 }
TripleArray m_F
Array containing the stage values.

References m_F.

◆ HasFirstQuadrature()

LUE bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::HasFirstQuadrature ( ) const
inline

Definition at line 205 of file TimeIntegrationSchemeSDC.h.

206 {
207 return m_first_quadrature;
208 }

References m_first_quadrature.

◆ HasLastQuadrature()

LUE bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::HasLastQuadrature ( ) const
inline

Definition at line 210 of file TimeIntegrationSchemeSDC.h.

211 {
212 return m_last_quadrature;
213 }

References m_last_quadrature.

◆ ResidualEval() [1/2]

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::ResidualEval ( const NekDouble delta_t)
inline

Definition at line 290 of file TimeIntegrationSchemeSDC.h.

291 {
292 v_ResidualEval(delta_t);
293 }
virtual LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n)

References v_ResidualEval().

◆ ResidualEval() [2/2]

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::ResidualEval ( const NekDouble delta_t,
const size_t  n 
)
inline

Definition at line 285 of file TimeIntegrationSchemeSDC.h.

286 {
287 v_ResidualEval(delta_t, n);
288 }

References v_ResidualEval().

◆ SDCIterationLoop()

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::SDCIterationLoop ( const NekDouble delta_t)
inline

Definition at line 300 of file TimeIntegrationSchemeSDC.h.

301 {
302 v_SDCIterationLoop(delta_t);
303 }
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t)

References v_SDCIterationLoop().

Referenced by v_TimeIntegrate().

◆ SetPFASST()

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::SetPFASST ( bool  pfasst)
inline

Definition at line 190 of file TimeIntegrationSchemeSDC.h.

191 {
192 m_PFASST = pfasst;
193 }

References m_PFASST.

◆ SetTime()

LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::SetTime ( double  time)
inline

Definition at line 195 of file TimeIntegrationSchemeSDC.h.

References m_time.

◆ UpdateFAScorrectionVector()

LUE TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateFAScorrectionVector ( )
inline

Definition at line 280 of file TimeIntegrationSchemeSDC.h.

281 {
282 return m_FAScorr;
283 }

References m_FAScorr.

◆ UpdateFirstQuadrature()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateFirstQuadrature ( void  )

Worker method that update the first quadrature.

Definition at line 217 of file TimeIntegrationSchemeSDC.cpp.

218{
220 for (size_t i = 0; i < m_nvars; ++i)
221 {
222 Vmath::Vcopy(m_npoints, Y_f[i], 1, m_Y[0][i], 1);
223 }
224}
AT< AT< NekDouble > > DoubleArray
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References m_last_quadrature, m_npoints, m_nQuadPts, m_nvars, m_Y, m_Y_f, and Vmath::Vcopy().

Referenced by v_TimeIntegrate().

◆ UpdateFirstQuadratureSolutionVector()

LUE DoubleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateFirstQuadratureSolutionVector ( )
inline

Definition at line 240 of file TimeIntegrationSchemeSDC.h.

241 {
242 return m_Y[0];
243 }

References m_Y.

◆ UpdateIntegratedResidualQFint()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateIntegratedResidualQFint ( const NekDouble delta_t)

Definition at line 299 of file TimeIntegrationSchemeSDC.cpp.

301{
302 // Zeroing integrated flux
303 for (size_t n = 0; n < m_nQuadPts; ++n)
304 {
305 for (size_t i = 0; i < m_nvars; ++i)
306 {
307 Vmath::Zero(m_npoints, m_QFint[n][i], 1);
308 }
309 }
310
311 // Update integrated flux
312 size_t offset = m_first_quadrature ? 0 : 1;
313 for (size_t n = 0; n < m_nQuadPts; ++n)
314 {
315 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
316 {
317 for (size_t i = 0; i < m_nvars; ++i)
318 {
320 m_npoints, delta_t * m_QMat[n * (m_nQuadPts - offset) + p],
321 m_F[p + offset][i], 1, m_QFint[n][i], 1, m_QFint[n][i], 1);
322 }
323 }
324 }
325}
SingleArray m_QMat
Array containing the integrated residual term.
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References m_F, m_first_quadrature, m_npoints, m_nQuadPts, m_nvars, m_QFint, m_QMat, CellMLToNektar.cellml_metadata::p, Vmath::Svtvp(), and Vmath::Zero().

◆ UpdateIntegratedResidualSFint()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateIntegratedResidualSFint ( const NekDouble delta_t)

Worker method that compute residual integral.

Definition at line 269 of file TimeIntegrationSchemeSDC.cpp.

271{
272 // Zeroing integrated flux
273 for (size_t n = 0; n < m_nQuadPts; ++n)
274 {
275 for (size_t i = 0; i < m_nvars; ++i)
276 {
277 Vmath::Zero(m_npoints, m_SFint[n][i], 1);
278 }
279 }
280
281 // Update integrated flux
282 size_t offset = m_first_quadrature ? 0 : 1;
283 for (size_t n = 1; n < m_nQuadPts; ++n)
284 {
285 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
286 {
287 for (size_t i = 0; i < m_nvars; ++i)
288 {
290 m_npoints,
291 delta_t * (m_QMat[n * (m_nQuadPts - offset) + p] -
292 m_QMat[(n - 1) * (m_nQuadPts - offset) + p]),
293 m_F[p + offset][i], 1, m_SFint[n][i], 1, m_SFint[n][i], 1);
294 }
295 }
296 }
297}

References m_F, m_first_quadrature, m_npoints, m_nQuadPts, m_nvars, m_QMat, m_SFint, CellMLToNektar.cellml_metadata::p, Vmath::Svtvp(), and Vmath::Zero().

Referenced by Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC::v_SDCIterationLoop(), Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC::v_SDCIterationLoop(), and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC::v_SDCIterationLoop().

◆ UpdateIntegratedResidualVector()

LUE TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateIntegratedResidualVector ( )
inline

Definition at line 270 of file TimeIntegrationSchemeSDC.h.

271 {
272 return m_QFint;
273 }

References m_QFint.

◆ UpdateLastQuadrature()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateLastQuadrature ( void  )

Worker method that update the last quadrature.

Definition at line 229 of file TimeIntegrationSchemeSDC.cpp.

230{
232 {
233 for (size_t i = 0; i < m_nvars; ++i)
234 {
236 for (size_t n = 0; n < m_nQuadPts; ++n)
237 {
238 Vmath::Svtvp(m_npoints, m_interp[n], m_Y[n][i], 1, m_Y_f[i], 1,
239 m_Y_f[i], 1);
240 }
241 }
242 }
243}
SingleArray m_interp
Array containing the integration matrix.

References m_interp, m_last_quadrature, m_npoints, m_nQuadPts, m_nvars, m_Y, m_Y_f, Vmath::Svtvp(), and Vmath::Zero().

Referenced by v_TimeIntegrate().

◆ UpdateLastQuadratureSolutionVector()

LUE DoubleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateLastQuadratureSolutionVector ( )
inline

Definition at line 250 of file TimeIntegrationSchemeSDC.h.

251 {
252 return m_last_quadrature ? m_Y[m_nQuadPts - 1] : m_Y_f;
253 }

References m_last_quadrature, m_nQuadPts, m_Y, and m_Y_f.

◆ UpdateResidualVector()

LUE TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::UpdateResidualVector ( )
inline

Definition at line 260 of file TimeIntegrationSchemeSDC.h.

261 {
262 return m_F;
263 }

References m_F.

◆ v_ComputeInitialGuess()

virtual LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_ComputeInitialGuess ( const NekDouble delta_t)
inlineprotectedvirtual

Reimplemented in Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC, Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC, and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC.

Definition at line 353 of file TimeIntegrationSchemeSDC.h.

355 {
356 ASSERTL0(false, "Specific version of spectral deferred correction "
357 "not implemented");
358 }

References ASSERTL0.

Referenced by ComputeInitialGuess().

◆ v_GetFreeParams()

std::vector< NekDouble > Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetFreeParams ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 55 of file TimeIntegrationSchemeSDC.cpp.

56{
57 return m_freeParams;
58}

References m_freeParams.

◆ v_GetIntegrationSchemeType()

TimeIntegrationSchemeType Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetIntegrationSchemeType ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 60 of file TimeIntegrationSchemeSDC.cpp.

62{
63 return m_schemeType;
64}

References m_schemeType.

◆ v_GetName()

std::string Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetName ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 40 of file TimeIntegrationSchemeSDC.cpp.

41{
42 return m_name;
43}

References m_name.

◆ v_GetNumIntegrationPhases()

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetNumIntegrationPhases ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 71 of file TimeIntegrationSchemeSDC.cpp.

72{
73 return 1;
74}

◆ v_GetOrder()

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetOrder ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 50 of file TimeIntegrationSchemeSDC.cpp.

51{
52 return m_order;
53}

References m_order.

◆ v_GetSolutionVector()

const TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetSolutionVector ( ) const
overrideprotectedvirtual

Gets the solution vector of the ODE.

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 76 of file TimeIntegrationSchemeSDC.cpp.

77{
78 return m_Y;
79}

References m_Y.

◆ v_GetTimeStability()

NekDouble Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetTimeStability ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 66 of file TimeIntegrationSchemeSDC.cpp.

67{
68 return 1.0;
69}

◆ v_GetVariant()

std::string Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_GetVariant ( ) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 45 of file TimeIntegrationSchemeSDC.cpp.

46{
47 return m_variant;
48}

References m_variant.

◆ v_InitializeScheme()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_InitializeScheme ( const NekDouble  deltaT,
ConstDoubleArray y_0,
const NekDouble  time,
const TimeIntegrationSchemeOperators op 
)
overrideprotectedvirtual

Worker method to initialize the integration scheme.

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 95 of file TimeIntegrationSchemeSDC.cpp.

98{
99 if (m_initialized)
100 {
101 m_time = time;
102 for (size_t i = 0; i < m_nvars; ++i)
103 {
104 // Store the initial values as the first previous state.
105 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
106 }
107 }
108 else
109 {
110 m_op = op;
111 m_time = time;
112 m_nvars = y_0.size();
113 m_npoints = y_0[0].size();
114
115 // Compute integration matrix.
116 size_t colOffset = m_first_quadrature ? 0 : 1;
117 size_t rowOffset = m_first_quadrature ? 0 : m_nQuadPts - 1;
118 size_t nCols = m_nQuadPts - colOffset;
119 size_t nRows = m_nQuadPts;
120 m_QMat = SingleArray(nRows * nCols, 0.0);
121 Polylib::Qg(&m_QMat[rowOffset], &m_tau[colOffset], nCols);
122
123 // Compute intepolation coefficient.
125 for (size_t i = 0; i < m_nQuadPts; i++)
126 {
127 m_interp[i] = Polylib::hgj(i, 1.0, &m_tau[0], m_nQuadPts, 0.0, 0.0);
128 }
129
130 // Storage of previous states and associated timesteps.
135 for (size_t m = 0; m < m_nQuadPts; ++m)
136 {
137 m_Y[m] = DoubleArray(m_nvars);
138 m_F[m] = DoubleArray(m_nvars);
141 for (size_t i = 0; i < m_nvars; ++i)
142 {
143 m_Y[m][i] = SingleArray(m_npoints, 0.0);
144 m_F[m][i] = SingleArray(m_npoints, 0.0);
145 m_SFint[m][i] = SingleArray(m_npoints, 0.0);
146 m_QFint[m][i] = SingleArray(m_npoints, 0.0);
147 // Store the initial values as the first previous state.
148 if (m == 0)
149 {
150 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[m][i], 1);
151 }
152 }
153 }
154
156 {
158 for (size_t i = 0; i < m_nvars; ++i)
159 {
160 m_Y_f[i] = SingleArray(m_npoints, 0.0);
161 }
162 }
163
164 if (m_PFASST)
165 {
167 for (size_t m = 0; m < m_nQuadPts; ++m)
168 {
170 for (size_t i = 0; i < m_nvars; ++i)
171 {
172 m_FAScorr[m][i] = SingleArray(m_npoints, 0.0);
173 }
174 }
175 }
176
177 m_initialized = true;
178 }
179}
AT< AT< AT< NekDouble > > > TripleArray
void Qg(double *Q, const double *z, const int np)
Compute the Integration Matrix.
Definition: Polylib.cpp:611
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
Definition: Polylib.cpp:914

References Polylib::hgj(), m_F, m_FAScorr, m_first_quadrature, m_initialized, m_interp, m_last_quadrature, m_npoints, m_nQuadPts, m_nvars, m_op, m_PFASST, m_QFint, m_QMat, m_SFint, m_tau, m_time, m_Y, m_Y_f, Polylib::Qg(), and Vmath::Vcopy().

Referenced by Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC::v_InitializeScheme(), Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC::v_InitializeScheme(), and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC::v_InitializeScheme().

◆ v_print()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_print ( std::ostream &  os) const
overrideprotectedvirtual

Worker method to print details on the integration scheme.

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 330 of file TimeIntegrationSchemeSDC.cpp.

331{
332 os << "Time Integration Scheme: " << GetFullName() << std::endl;
333}

References Nektar::LibUtilities::TimeIntegrationScheme::GetFullName().

◆ v_printFull()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_printFull ( std::ostream &  os) const
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 335 of file TimeIntegrationSchemeSDC.cpp.

336{
337 os << "Time Integration Scheme: " << GetFullName() << std::endl;
338}

References Nektar::LibUtilities::TimeIntegrationScheme::GetFullName().

◆ v_ResidualEval() [1/2]

virtual LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_ResidualEval ( const NekDouble delta_t)
inlineprotectedvirtual

Reimplemented in Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC, Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC, and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC.

Definition at line 347 of file TimeIntegrationSchemeSDC.h.

348 {
349 ASSERTL0(false, "Specific version of spectral deferred correction "
350 "not implemented");
351 }

References ASSERTL0.

◆ v_ResidualEval() [2/2]

virtual LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_ResidualEval ( const NekDouble delta_t,
const size_t  n 
)
inlineprotectedvirtual

Reimplemented in Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC, Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC, and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC.

Definition at line 340 of file TimeIntegrationSchemeSDC.h.

342 {
343 ASSERTL0(false, "Specific version of spectral deferred correction "
344 "not implemented");
345 }

References ASSERTL0.

Referenced by ResidualEval().

◆ v_SDCIterationLoop()

virtual LUE void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_SDCIterationLoop ( const NekDouble delta_t)
inlineprotectedvirtual

Reimplemented in Nektar::LibUtilities::ExplicitTimeIntegrationSchemeSDC, Nektar::LibUtilities::IMEXTimeIntegrationSchemeSDC, and Nektar::LibUtilities::ImplicitTimeIntegrationSchemeSDC.

Definition at line 360 of file TimeIntegrationSchemeSDC.h.

362 {
363 ASSERTL0(false, "Specific version of spectral deferred correction "
364 "not implemented");
365 }

References ASSERTL0.

Referenced by SDCIterationLoop().

◆ v_SetSolutionVector()

void Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_SetSolutionVector ( const size_t  Offset,
const DoubleArray y 
)
overrideprotectedvirtual

Sets the solution vector of the ODE.

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 86 of file TimeIntegrationSchemeSDC.cpp.

88{
89 m_Y[Offset] = y;
90}

References m_Y.

◆ v_TimeIntegrate()

ConstDoubleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_TimeIntegrate ( const size_t  timestep,
const NekDouble  delta_t 
)
overrideprotectedvirtual

Worker method that performs the time integration.

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 184 of file TimeIntegrationSchemeSDC.cpp.

186{
187 for (size_t k = 0; k < m_order; ++k)
188 {
189 // Compute initial guess
190 if (k == 0)
191 {
192 ComputeInitialGuess(delta_t);
193 }
194 // Apply SDC correction loop
195 else
196 {
197 SDCIterationLoop(delta_t);
198 }
199 }
200
201 // Update last quadrature
203
204 // Update first quadrature
206
207 // Update time step
208 m_time += delta_t;
209
210 // Return solution
211 return m_Y[0];
212}
LUE void ComputeInitialGuess(const NekDouble &delta_t)
LUE void UpdateFirstQuadrature(void)
Worker method that update the first quadrature.
LUE void UpdateLastQuadrature(void)
Worker method that update the last quadrature.
LUE void SDCIterationLoop(const NekDouble &delta_t)

References ComputeInitialGuess(), m_order, m_time, m_Y, SDCIterationLoop(), UpdateFirstQuadrature(), and UpdateLastQuadrature().

◆ v_UpdateSolutionVector()

TripleArray & Nektar::LibUtilities::TimeIntegrationSchemeSDC::v_UpdateSolutionVector ( )
overrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 81 of file TimeIntegrationSchemeSDC.cpp.

82{
83 return m_Y;
84}

References m_Y.

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::TimeIntegrationSchemeSDC::className
static

Definition at line 188 of file TimeIntegrationSchemeSDC.h.

◆ m_F

TripleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_F
protected

◆ m_FAScorr

TripleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_FAScorr
protected

Array containing the stage derivatives.

Definition at line 384 of file TimeIntegrationSchemeSDC.h.

Referenced by AddFASCorrectionToSFint(), GetFAScorrectionVector(), UpdateFAScorrectionVector(), and v_InitializeScheme().

◆ m_first_quadrature

bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_first_quadrature {true}
protected

◆ m_freeParams

std::vector<NekDouble> Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_freeParams
protected

Definition at line 374 of file TimeIntegrationSchemeSDC.h.

Referenced by TimeIntegrationSchemeSDC(), and v_GetFreeParams().

◆ m_initialized

bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_initialized {false}
protected

◆ m_interp

SingleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_interp
protected

Array containing the integration matrix.

Definition at line 388 of file TimeIntegrationSchemeSDC.h.

Referenced by UpdateLastQuadrature(), and v_InitializeScheme().

◆ m_last_quadrature

bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_last_quadrature {true}
protected

◆ m_name

std::string Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_name
protected

◆ m_npoints

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_npoints {0}
protected

◆ m_nQuadPts

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_nQuadPts {0}
protected

◆ m_nvars

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_nvars {0}
protected

◆ m_op

TimeIntegrationSchemeOperators Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_op
protected

◆ m_order

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_order {0}
protected

Maximum order of the integration scheme.

Definition at line 394 of file TimeIntegrationSchemeSDC.h.

Referenced by TimeIntegrationSchemeSDC(), v_GetOrder(), and v_TimeIntegrate().

◆ m_ordermax

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_ordermax {0}
protected

Minimum order of the integration scheme.

Definition at line 393 of file TimeIntegrationSchemeSDC.h.

Referenced by GetMaxOrder(), and TimeIntegrationSchemeSDC().

◆ m_ordermin

size_t Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_ordermin {0}
protected

SDC parameter.

Definition at line 392 of file TimeIntegrationSchemeSDC.h.

Referenced by TimeIntegrationSchemeSDC().

◆ m_PFASST

bool Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_PFASST {false}
protected

◆ m_pointsKey

PointsKey Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_pointsKey
protected

Definition at line 379 of file TimeIntegrationSchemeSDC.h.

Referenced by GetPointsKey(), and TimeIntegrationSchemeSDC().

◆ m_QFint

TripleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_QFint
protected

Array containing the integrated residual term.

Definition at line 386 of file TimeIntegrationSchemeSDC.h.

Referenced by GetIntegratedResidualVector(), UpdateIntegratedResidualQFint(), UpdateIntegratedResidualVector(), and v_InitializeScheme().

◆ m_QMat

SingleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_QMat
protected

Array containing the integrated residual term.

Definition at line 387 of file TimeIntegrationSchemeSDC.h.

Referenced by UpdateIntegratedResidualQFint(), UpdateIntegratedResidualSFint(), and v_InitializeScheme().

◆ m_schemeType

TimeIntegrationSchemeType Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_schemeType {eNoTimeIntegrationSchemeType}
protected

◆ m_SFint

TripleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_SFint
protected

◆ m_tau

SingleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_tau
protected

◆ m_theta

NekDouble Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_theta {1.0}
protected

◆ m_time

NekDouble Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_time
protected

◆ m_variant

std::string Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_variant
protected

Definition at line 373 of file TimeIntegrationSchemeSDC.h.

Referenced by TimeIntegrationSchemeSDC(), and v_GetVariant().

◆ m_Y

TripleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_Y
protected

◆ m_Y_f

DoubleArray Nektar::LibUtilities::TimeIntegrationSchemeSDC::m_Y_f
protected