Nektar++
ExplicitTimeIntegrationSchemeSDC.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExplicitTimeIntegrationSchemeSDC.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_EXPLICIT_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_EXPLICIT_TIME_INTEGRATION_SCHEME_SDC
37
38#define LUE LIB_UTILITIES_EXPORT
39
40#include <string>
41
43
45{
47{
48public:
49 ExplicitTimeIntegrationSchemeSDC(std::string variant, size_t order,
50 std::vector<NekDouble> freeParams)
51 : TimeIntegrationSchemeSDC(variant, order, freeParams)
52 {
53 m_name = "ExplicitSpectralDeferredCorrection";
55 }
56
58 std::string variant, size_t order, std::vector<NekDouble> freeParams)
59 {
62 variant, order, freeParams);
63
64 return p;
65 }
66
67 static std::string className;
68
69protected:
71 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
72 const TimeIntegrationSchemeOperators &op) override;
73
74 LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override;
75
76 LUE void v_ResidualEval(const NekDouble &delta_t) override;
77
78 LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override;
79
80 LUE void v_SDCIterationLoop(const NekDouble &delta_t) override;
81}; // end class ExplicitTimeIntegrationSchemeSDC
82
83/**
84 * @brief Worker method to initialize the integration scheme.
85 */
87 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
89{
91}
92
93/**
94 * @brief Worker method to compute the residual.
95 */
97 const size_t n)
98{
99 m_op.DoProjection(m_Y[n], m_Y[n], m_time + delta_t * m_tau[n]);
100 m_op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_tau[n]);
101}
102
104{
105 for (size_t n = 0; n < m_nQuadPts; ++n)
106 {
107 v_ResidualEval(delta_t, n);
108 }
109}
110
111/**
112 * @brief Worker method to compute the initial SDC guess.
113 */
115 const NekDouble &delta_t)
116{
117 for (size_t n = 0; n < m_nQuadPts; ++n)
118 {
119 // Use explicit Euler as a first guess
120 if (n > 0)
121 {
122 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
123 for (size_t i = 0; i < m_nvars; ++i)
124 {
125 Vmath::Svtvp(m_npoints, dtn, m_F[n - 1][i], 1, m_Y[n - 1][i], 1,
126 m_Y[n][i], 1);
127 }
128 }
129
130 // Compute residual
131 m_op.DoProjection(m_Y[n], m_Y[n], m_time + delta_t * m_tau[n]);
132 m_op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_tau[n]);
133 }
134}
135
136/**
137 * @brief Worker method to compute the SDC iteration.
138 */
140 const NekDouble &delta_t)
141{
142 // Update integrated residual
144
145 // Add FAS correction to integrated residual
147
148 // Loop over quadrature points
149 for (size_t n = 1; n < m_nQuadPts; ++n)
150 {
151 // Update solution
152 for (size_t i = 0; i < m_nvars; ++i)
153 {
154 // Add SFint contribution to solution
155 Vmath::Vadd(m_npoints, m_Y[n - 1][i], 1, m_SFint[n][i], 1,
156 m_Y[n][i], 1);
157
158 // Add explicit contribution to solution
159 if (n > 1)
160 {
161 NekDouble dtn = delta_t * (m_tau[n] - m_tau[n - 1]);
162 Vmath::Svtvp(m_npoints, m_theta * dtn, m_F[n - 1][i], 1,
163 m_Y[n][i], 1, m_Y[n][i], 1);
164 }
165
166 // Add explicit contribution to SFint
167 if (n < m_nQuadPts - 1)
168 {
169 NekDouble dtnp = delta_t * (m_tau[n + 1] - m_tau[n]);
170 Vmath::Svtvp(m_npoints, -m_theta * dtnp, m_F[n][i], 1,
171 m_SFint[n + 1][i], 1, m_SFint[n + 1][i], 1);
172 }
173 }
174
175 // Compute residual
176 m_op.DoProjection(m_Y[n], m_Y[n], m_time + delta_t * m_tau[n]);
177 m_op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_tau[n]);
178 }
179}
180
181} // namespace Nektar::LibUtilities
182
183#endif
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override
Worker method to compute the initial SDC guess.
ExplicitTimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(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 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.
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.
@ eExplicit
Formally explicit 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