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
45{
47{
48public:
49 ImplicitTimeIntegrationSchemeSDC(std::string variant, size_t order,
50 std::vector<NekDouble> freeParams)
51 : TimeIntegrationSchemeSDC(variant, order, freeParams)
52 {
54 "Quadrature type that include the left end point (e.g. "
55 "GaussLobattoLegendre) should not be used for ImplicitSDC");
56
57 m_name = "ImplicitSpectralDeferredCorrection";
59 }
60
62 std::string variant, size_t order, std::vector<NekDouble> freeParams)
63 {
66 variant, order, freeParams);
67
68 return p;
69 }
70
71 static std::string className;
72
73 DoubleArray m_tmp; /// Array for temporary storage
74
75protected:
77 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
78 const TimeIntegrationSchemeOperators &op) override;
79
80 LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override;
81
82 LUE void v_ResidualEval(const NekDouble &delta_t) override;
83
84 LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override;
85
86 LUE void v_SDCIterationLoop(const NekDouble &delta_t) override;
87}; // end class ImplicitTimeIntegrationSchemeSDC
88
89/**
90 * @brief Worker method to initialize the integration scheme.
91 */
93 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
95{
96 if (m_initialized)
97 {
98 m_time = time;
99 for (size_t i = 0; i < m_nvars; ++i)
100 {
101 // Store the initial values as the first previous state.
102 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
103 }
104 }
105 else
106 {
108
110 for (size_t i = 0; i < m_nvars; ++i)
111 {
112 m_tmp[i] = SingleArray(m_npoints, 0.0);
113 }
114 }
115}
116
117/**
118 * @brief Worker method to compute the residual.
119 */
121 [[maybe_unused]] const NekDouble &delta_t, [[maybe_unused]] const size_t n)
122{
123 // Apply time-dependent boundary condition
124 m_op.DoProjection(m_Y[n], m_tmp, m_time + delta_t * m_tau[n]);
125
126 // Compute residual
127 m_op.DoOdeRhs(m_tmp, m_F[n], m_time + delta_t * m_tau[n]);
128}
129
131{
132 for (size_t n = 0; n < m_nQuadPts; ++n)
133 {
134 v_ResidualEval(delta_t, n);
135 }
136}
137
138/**
139 * @brief Worker method to compute the initial SDC guess.
140 */
142 const NekDouble &delta_t)
143{
144 for (size_t n = 0; n < m_nQuadPts; ++n)
145 {
146 if (n != 0)
147 {
148 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
149
150 // Update solution
151 m_op.DoImplicitSolve(m_Y[n - 1], m_Y[n],
152 m_time + delta_t * m_tau[n], dtn);
153
154 // Compute residual from updated solution
155 for (size_t i = 0; i < m_nvars; ++i)
156 {
157 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_Y[n - 1][i], 1,
158 m_F[n][i], 1);
159 Vmath::Smul(m_npoints, 1.0 / dtn, m_F[n][i], 1, m_F[n][i], 1);
160 }
161 }
162 }
163}
164
165/**
166 * @brief Worker method to compute the SDC iteration.
167 */
169 const NekDouble &delta_t)
170{
171 // Update integrated residual
173
174 // Add FAS correction to integrated residual
176
177 // Loop over quadrature points
178 for (size_t n = 1; n < m_nQuadPts; ++n)
179 {
180 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
181
182 // Add rhs terms
183 for (size_t i = 0; i < m_nvars; ++i)
184 {
185 // Add SFint contribution
186 Vmath::Vadd(m_npoints, m_Y[n - 1][i], 1, m_SFint[n][i], 1, m_tmp[i],
187 1);
188
189 // Add implicit contribution
190 Vmath::Svtvp(m_npoints, -m_theta * dtn, m_F[n][i], 1, m_tmp[i], 1,
191 m_tmp[i], 1);
192 }
193
194 // Solve implicit system
195 m_op.DoImplicitSolve(m_tmp, m_Y[n], m_time + delta_t * m_tau[n],
196 m_theta * dtn);
197
198 // Compute residual from updated solution
199 for (size_t i = 0; i < m_nvars; ++i)
200 {
201 Vmath::Vsub(m_npoints, m_Y[n][i], 1, m_tmp[i], 1, m_F[n][i], 1);
202 Vmath::Smul(m_npoints, 1.0 / (m_theta * dtn), m_F[n][i], 1,
203 m_F[n][i], 1);
204 }
205 }
206}
207
208} // namespace Nektar::LibUtilities
209
210#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Array for temporary storage.
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)
LUE void v_SDCIterationLoop(const NekDouble &delta_t) override
Worker method to compute the SDC iteration.
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
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
@ eImplicit
Fully implicit scheme.
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