Nektar++
IMEXTimeIntegrationSchemeSDC.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: IMEXTimeIntegrationSchemeSDC.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Header file of time integration scheme SDC base class
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_IMEX_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_IMEX_TIME_INTEGRATION_SCHEME_SDC
37
38#define LUE LIB_UTILITIES_EXPORT
39
40#include <string>
41
43
45{
47{
48public:
49 IMEXTimeIntegrationSchemeSDC(std::string variant, size_t order,
50 std::vector<NekDouble> freeParams)
51 : TimeIntegrationSchemeSDC(variant, order, freeParams)
52 {
53
55 "Quadrature type that include the left end point (e.g. "
56 "GaussLobattoLegendre) should not be used for IMEXSDC");
57
58 std::cerr
59 << "WARNING: IMEX Spectral Deferred Correction method has been "
60 "implemented but its use is not recommended as the "
61 "approach is affected by order-reduction problems."
62 << std::endl;
63
64 m_name = "IMEXSpectralDeferredCorrection";
66 }
67
69 std::string variant, size_t order, std::vector<NekDouble> freeParams)
70 {
73 variant, order, freeParams);
74
75 return p;
76 }
77
78 static std::string className;
79
80 TripleArray m_Fexp; /// Array corresponding to the stage derivatives
81 TripleArray m_Fimp; /// Array corresponding to the stage derivatives
82 DoubleArray m_tmp; /// Array for temporary storage
83
84protected:
86 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
87 const TimeIntegrationSchemeOperators &op) override;
88
89 LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override;
90
91 LUE void v_ResidualEval(const NekDouble &delta_t) override;
92
93 LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override;
94
95 LUE void v_SDCIterationLoop(const NekDouble &delta_t) override;
96
97private:
98 void ComputeTotalResidual(const size_t n);
99
100}; // end class IMEXTimeIntegrationSchemeSDC
101
102/**
103 * @brief Worker method to initialize the integration scheme.
104 */
106 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
108{
109 if (m_initialized)
110 {
111 m_time = time;
112 for (size_t i = 0; i < m_nvars; ++i)
113 {
114 // Store the initial values as the first previous state.
115 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
116 }
117 }
118 else
119 {
121
124 for (size_t m = 0; m < m_nQuadPts; ++m)
125 {
128 for (size_t i = 0; i < m_nvars; ++i)
129 {
130 m_Fexp[m][i] = SingleArray(m_npoints, 0.0);
131 m_Fimp[m][i] = SingleArray(m_npoints, 0.0);
132 }
133 }
134
136 for (size_t i = 0; i < m_nvars; ++i)
137 {
138 m_tmp[i] = SingleArray(m_npoints, 0.0);
139 }
140 }
141}
142
143/**
144 * @brief Worker method to compute the residual.
145 */
147 [[maybe_unused]] const NekDouble &delta_t, [[maybe_unused]] const size_t n)
148{
149 ASSERTL0(false, "v_ResidualEval not implemented for IMEX SDC");
150}
151
153{
154 for (size_t n = 0; n < m_nQuadPts; ++n)
155 {
156 v_ResidualEval(delta_t, n);
157 }
158}
159
160/**
161 * @brief Worker method to compute the initial SDC guess.
162 */
164 const NekDouble &delta_t)
165{
166 for (size_t n = 0; n < m_nQuadPts; ++n)
167 {
168 if (n == 0)
169 {
170 // Apply time-dependent boundary condition
171 m_op.DoProjection(m_Y[0], m_Y[0], m_time);
172 }
173 else
174 {
175 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
176
177 // Add explicit contribution to rhs
178 for (size_t i = 0; i < m_nvars; ++i)
179 {
180 Vmath::Svtvp(m_npoints, dtn, m_Fexp[n - 1][i], 1, m_Y[n - 1][i],
181 1, m_tmp[i], 1);
182 }
183
184 // Solve implicit system from rhs
185 m_op.DoImplicitSolve(m_tmp, m_Y[n], m_time + delta_t * m_tau[n],
186 dtn);
187
188 // Compute implicit flux from updated solution
189 for (size_t i = 0; i < m_nvars; ++i)
190 {
191 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_tmp[i], 1, m_Fimp[n][i],
192 1);
193 Vmath::Smul(m_npoints, 1.0 / dtn, m_Fimp[n][i], 1, m_Fimp[n][i],
194 1);
195 }
196 }
197
198 // Compute explicit residual
199 m_op.DoOdeRhs(m_Y[n], m_Fexp[n], m_time + delta_t * m_tau[n]);
200
201 // Compute total residual
203 }
204}
205
206/**
207 * @brief Worker method to compute the SDC iteration.
208 */
210{
211 // Update integrated residual
213
214 // Add FAS correction to integrated residual
216
217 // Loop over quadrature points
218 for (size_t n = 1; n < m_nQuadPts; ++n)
219 {
220 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
221
222 // Add rhs terms
223 for (size_t i = 0; i < m_nvars; ++i)
224 {
225 // Add SFint contribution to rhs
226 Vmath::Vadd(m_npoints, m_Y[n - 1][i], 1, m_SFint[n][i], 1, m_tmp[i],
227 1);
228
229 // Add explicit contribution to rhs
230 if (n > 1)
231 {
232 Vmath::Svtvp(m_npoints, m_theta * dtn, m_Fexp[n - 1][i], 1,
233 m_tmp[i], 1, m_tmp[i], 1);
234 }
235
236 // Add explicit contribution to SFint
237 if (n < m_nQuadPts - 1)
238 {
239 NekDouble dtnp = delta_t * (m_tau[n + 1] - m_tau[n]);
240 Vmath::Svtvp(m_npoints, -m_theta * dtnp, m_Fexp[n][i], 1,
241 m_SFint[n + 1][i], 1, m_SFint[n + 1][i], 1);
242 }
243
244 // Add implicit contribution to rhs
245 Vmath::Svtvp(m_npoints, -m_theta * dtn, m_Fimp[n][i], 1, m_tmp[i],
246 1, m_tmp[i], 1);
247 }
248
249 // Solve implicit system from rhs
250 m_op.DoImplicitSolve(m_tmp, m_Y[n], m_time + delta_t * m_tau[n],
251 m_theta * dtn);
252
253 // Compute implicit residual from updated solution
254 for (size_t i = 0; i < m_nvars; ++i)
255 {
256 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_tmp[i], 1, m_Fimp[n][i], 1);
257 Vmath::Smul(m_npoints, 1.0 / (m_theta * dtn), m_Fimp[n][i], 1,
258 m_Fimp[n][i], 1);
259 }
260
261 // Compute explicit residual
262 m_op.DoOdeRhs(m_Y[n], m_Fexp[n], m_time + delta_t * m_tau[n]);
263
264 // Compute total residual
266 }
267}
268
269/**
270 * @brief Worker method to compute the total residual.
271 */
273{
274 for (size_t i = 0; i < m_nvars; ++i)
275 {
276 Vmath::Vadd(m_npoints, m_Fimp[n][i], 1, m_Fexp[n][i], 1, m_F[n][i], 1);
277 }
278}
279
280} // namespace Nektar::LibUtilities
281
282#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override
Worker method to compute the initial SDC guess.
LUE void v_SDCIterationLoop(const NekDouble &delta_t) override
Worker method to compute the SDC iteration.
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Array for temporary storage.
DoubleArray m_tmp
Array corresponding to the stage derivatives.
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
IMEXTimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override
Worker method to compute the residual.
TripleArray m_Fimp
Array corresponding to the stage derivatives.
void ComputeTotalResidual(const size_t n)
Worker method to compute the total residual.
Binds a set of functions for use by time integration schemes.
void DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) const
void DoProjection(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
Class for spectral deferred correction integration.
LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t)
Worker method that compute residual integral.
TripleArray m_Y
Array containing the last stage values.
size_t m_nQuadPts
Order of the integration scheme.
LUE void AddFASCorrectionToSFint(void)
Worker method that add the FASCorrection.
bool m_first_quadrature
Number of points in the integration scheme.
size_t m_npoints
Number of variables in the integration scheme.
TripleArray m_SFint
Array containing the FAS correction term.
NekDouble m_theta
Array containing the interpolation coefficients.
SingleArray m_tau
Object containing quadrature data.
TripleArray m_F
Array containing the stage values.
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
AT< AT< NekDouble > > DoubleArray
@ eIMEX
Implicit Explicit General Linear Method.
AT< AT< AT< NekDouble > > > TripleArray
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
double NekDouble
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 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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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