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