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 
44 
45 namespace Nektar
46 {
47 namespace LibUtilities
48 {
50 {
51 public:
52  ExplicitTimeIntegrationSchemeSDC(std::string variant, unsigned int order,
53  std::vector<NekDouble> freeParams)
54  : TimeIntegrationSchemeSDC(variant, order, freeParams)
55  {
56  m_name = "ExplicitSpectralDeferredCorrection";
58  }
59 
61  std::string variant, unsigned int order,
62  std::vector<NekDouble> freeParams)
63  {
66  variant, order, freeParams);
67 
68  return p;
69  }
70 
71  static std::string className;
72 
73 protected:
74  LUE virtual void v_InitializeScheme(
75  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
76  const TimeIntegrationSchemeOperators &op) override;
77 
78  LUE virtual void v_ResidualEval(
79  const NekDouble &delta_t, const int n,
80  const TimeIntegrationSchemeOperators &op) override;
81 
82  LUE virtual void v_ResidualEval(
83  const NekDouble &delta_t,
84  const TimeIntegrationSchemeOperators &op) override;
85 
86  LUE virtual void v_ComputeInitialGuess(
87  const NekDouble &delta_t,
88  const TimeIntegrationSchemeOperators &op) override;
89 
90  LUE virtual void v_SDCIterationLoop(
91  const NekDouble &delta_t,
92  const TimeIntegrationSchemeOperators &op) override;
93 
94 }; // end class ExplicitTimeIntegrationSchemeSDC
95 
96 /**
97  * @brief Worker method to initialize the integration scheme.
98  */
100  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
102 {
103  TimeIntegrationSchemeSDC::v_InitializeScheme(deltaT, y_0, time, op);
104 
105  op.DoProjection(m_Y[0], m_Y[0], m_time);
106 }
107 
108 /**
109  * @brief Worker method to compute the residual.
110  */
112  const NekDouble &delta_t, const int n,
114 {
115  // Compute residual
116  op.DoProjection(m_Y[n], m_Y[n], m_time + delta_t * m_points[n]);
117  op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_points[n]);
118 }
119 
121  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
122 {
123  // Compute residual
124  for (unsigned int n = 0; n < m_nQuadPts; ++n)
125  {
126  v_ResidualEval(delta_t, n, op);
127  }
128 }
129 
130 /**
131  * @brief Worker method to compute the initial SDC guess.
132  */
134  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
135 {
136  // Compute initial guess
137  for (unsigned int n = 0; n < m_nQuadPts; ++n)
138  {
139  // Compute residual
140  op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_points[n]);
141 
142  // Use explicit Euler as a first guess
143  if (n < m_nQuadPts - 1)
144  {
145  NekDouble dtn = delta_t * (m_points[n + 1] - m_points[n]);
146  for (unsigned int i = 0; i < m_nvars; ++i)
147  {
148  Vmath::Svtvp(m_npoints, dtn, m_F[n][i], 1, m_Y[n][i], 1,
149  m_Y[n + 1][i], 1);
150  }
151  op.DoProjection(m_Y[n + 1], m_Y[n + 1],
152  m_time + delta_t * m_points[n + 1]);
153  }
154  }
155 }
156 
157 /**
158  * @brief Worker method to compute the SDC iteration.
159  */
161  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
162 {
163  unsigned int kstart = 1;
164  for (unsigned int k = kstart; k < m_order; ++k)
165  {
166  // Update integrated residual
167  UpdateIntegratedResidual(delta_t);
168 
169  // Loop over quadrature points
170  for (unsigned int n = 0; n < m_nQuadPts - 1; ++n)
171  {
172  NekDouble dtn = delta_t * (m_points[n + 1] - m_points[n]);
173 
174  // Update residual if n > 0
175  if (n > 0)
176  {
177  op.DoOdeRhs(m_Y[n], m_Fn, m_time + delta_t * m_points[n]);
178  }
179 
180  // Update solution
181  for (unsigned int i = 0; i < m_nvars; ++i)
182  {
183  // Add Fint contribution
184  Vmath::Vadd(m_npoints, m_Y[n][i], 1, m_Fint[n][i], 1,
185  m_Y[n + 1][i], 1);
186 
187  if (n > 0)
188  {
189  // Add explicit contribution
190  Vmath::Svtvp(m_npoints, m_theta * dtn, m_Fn[i], 1,
191  m_Y[n + 1][i], 1, m_Y[n + 1][i], 1);
192  Vmath::Svtvp(m_npoints, -m_theta * dtn, m_F[n][i], 1,
193  m_Y[n + 1][i], 1, m_Y[n + 1][i], 1);
194 
195  // Copy new rhs value to old
196  Vmath::Vcopy(m_npoints, m_Fn[i], 1, m_F[n][i], 1);
197  }
198  }
199  op.DoProjection(m_Y[n + 1], m_Y[n + 1],
200  m_time + delta_t * m_points[n + 1]);
201  }
202  op.DoOdeRhs(m_Y[m_nQuadPts - 1], m_F[m_nQuadPts - 1],
203  m_time + delta_t * m_points[m_nQuadPts - 1]);
204  }
205 }
206 
207 } // end namespace LibUtilities
208 } // end namespace Nektar
209 
210 #endif
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op) override
Worker method to compute the initial SDC guess.
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.
ExplicitTimeIntegrationSchemeSDC(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op) override
Worker method to compute the SDC iteration.
virtual LUE void v_ResidualEval(const NekDouble &delta_t, const int n, const TimeIntegrationSchemeOperators &op) override
Worker method to compute the residual.
static TimeIntegrationSchemeSharedPtr create(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
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 UpdateIntegratedResidual(const NekDouble &delta_t, const int option=0)
Worker method that compute the integrated flux.
TripleArray m_Fint
Array corresponding to the stage Derivatives.
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:622
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:359
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255