Nektar++
ImplicitTimeIntegrationSchemeSDC.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ImplicitTimeIntegrationSchemeSDC.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_IMPLICIT_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_IMPLICIT_TIME_INTEGRATION_SCHEME_SDC
37
38#define LUE LIB_UTILITIES_EXPORT
39
40#include <string>
41
43
44namespace Nektar
45{
46namespace LibUtilities
47{
49{
50public:
51 ImplicitTimeIntegrationSchemeSDC(std::string variant, size_t order,
52 std::vector<NekDouble> freeParams)
53 : TimeIntegrationSchemeSDC(variant, order, freeParams)
54 {
56 "Quadrature type that include the left end point (e.g. "
57 "GaussLobattoLegendre) should not be used for ImplicitSDC");
58
59 m_name = "ImplicitSpectralDeferredCorrection";
61 }
62
64 std::string variant, size_t order, std::vector<NekDouble> freeParams)
65 {
68 variant, order, freeParams);
69
70 return p;
71 }
72
73 static std::string className;
74
75 DoubleArray m_tmp; /// Array for temporary storage
76
77protected:
78 LUE virtual void v_InitializeScheme(
79 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
80 const TimeIntegrationSchemeOperators &op) override;
81
82 LUE virtual void v_ResidualEval(const NekDouble &delta_t,
83 const size_t n) override;
84
85 LUE virtual void v_ResidualEval(const NekDouble &delta_t) override;
86
87 LUE virtual void v_ComputeInitialGuess(const NekDouble &delta_t) override;
88
89 LUE virtual void v_SDCIterationLoop(const NekDouble &delta_t) override;
90
91}; // end class ImplicitTimeIntegrationSchemeSDC
92
93/**
94 * @brief Worker method to initialize the integration scheme.
95 */
97 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
99{
100 if (m_initialized)
101 {
102 m_time = time;
103 for (size_t i = 0; i < m_nvars; ++i)
104 {
105 // Store the initial values as the first previous state.
106 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
107 }
108 }
109 else
110 {
112
114 for (size_t i = 0; i < m_nvars; ++i)
115 {
116 m_tmp[i] = SingleArray(m_npoints, 0.0);
117 }
118 }
119}
120
121/**
122 * @brief Worker method to compute the residual.
123 */
125 const size_t n)
126{
127 if (n == 0)
128 {
129 // Not implemented, require implicit evaluation for m_F[0].
130 // Quadrature type that include the left end point (e.g.
131 // GaussLobattoLegendre) should not be used.
132 }
133 else
134 {
135 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
136
137 // Update solution
138 m_op.DoImplicitSolve(m_Y[n - 1], m_tmp, m_time + delta_t * m_tau[n],
139 m_theta * dtn);
140
141 // Compute residual from updated solution
142 for (size_t i = 0; i < m_nvars; ++i)
143 {
144 Vmath::Vsub(m_npoints, m_tmp[i], 1, m_Y[n - 1][i], 1, m_F[n][i], 1);
145 Vmath::Smul(m_npoints, 1.0 / (m_theta * dtn), m_F[n][i], 1,
146 m_F[n][i], 1);
147 }
148 }
149}
150
152{
153 for (size_t n = 0; n < m_nQuadPts; ++n)
154 {
155 v_ResidualEval(delta_t, n);
156 }
157}
158
159/**
160 * @brief Worker method to compute the initial SDC guess.
161 */
163 const NekDouble &delta_t)
164{
165 for (size_t n = 0; n < m_nQuadPts; ++n)
166 {
167 if (n == 0)
168 {
169 // Not implemented, require implicit evaluation for m_F[0].
170 // Quadrature type that include the left end point (e.g.
171 // GaussLobattoLegendre) should not be used.
172 }
173 else
174 {
175 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
176
177 // Update solution
178 m_op.DoImplicitSolve(m_Y[n - 1], m_Y[n],
179 m_time + delta_t * m_tau[n], dtn);
180
181 // Compute residual from updated solution
182 for (size_t i = 0; i < m_nvars; ++i)
183 {
184 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_Y[n - 1][i], 1,
185 m_F[n][i], 1);
186 Vmath::Smul(m_npoints, 1.0 / dtn, m_F[n][i], 1, m_F[n][i], 1);
187 }
188 }
189 }
190}
191
192/**
193 * @brief Worker method to compute the SDC iteration.
194 */
196 const NekDouble &delta_t)
197{
198 // Update integrated residual
200
201 // Add FAS correction to integrated residual
203
204 // Loop over quadrature points
205 for (size_t n = 1; n < m_nQuadPts; ++n)
206 {
207 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
208
209 // Add rhs terms
210 for (size_t i = 0; i < m_nvars; ++i)
211 {
212 // Add SFint contribution
213 Vmath::Vadd(m_npoints, m_Y[n - 1][i], 1, m_SFint[n][i], 1, m_tmp[i],
214 1);
215
216 // Add implicit contribution
217 Vmath::Svtvp(m_npoints, -m_theta * dtn, m_F[n][i], 1, m_tmp[i], 1,
218 m_tmp[i], 1);
219 }
220
221 // Solve implicit system
222 m_op.DoImplicitSolve(m_tmp, m_Y[n], m_time + delta_t * m_tau[n],
223 m_theta * dtn);
224
225 // Compute residual from updated solution
226 for (size_t i = 0; i < m_nvars; ++i)
227 {
228 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_tmp[i], 1, m_F[n][i], 1);
229 Vmath::Smul(m_npoints, 1.0 / (m_theta * dtn), m_F[n][i], 1,
230 m_F[n][i], 1);
231 }
232 }
233}
234
235} // end namespace LibUtilities
236} // end namespace Nektar
237
238#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Array for temporary storage.
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override
Worker method to compute the initial SDC guess.
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
ImplicitTimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t) override
Worker method to compute the SDC iteration.
virtual LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override
Worker method to compute the 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
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.
virtual 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
@ eImplicit
Fully implicit scheme.
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:617
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.cpp:354
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.cpp:245
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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.cpp:414