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 
44 
45 namespace Nektar
46 {
47 namespace LibUtilities
48 {
50 {
51 public:
52  ImplicitTimeIntegrationSchemeSDC(std::string variant, unsigned int order,
53  std::vector<NekDouble> freeParams)
54  : TimeIntegrationSchemeSDC(variant, order, freeParams)
55  {
56  ASSERTL0(variant == "GaussRadauLegendre",
57  "only GaussRadauLegendre variant "
58  "(quadrature) type available for ImplicitSDC");
59  ASSERTL0(0.0 < freeParams[0],
60  "Spectral Deferred Correction Time integration "
61  "scheme bad parameter numbers (> 0.0): " +
62  std::to_string(freeParams[0]));
63 
64  m_name = "ImplicitSpectralDeferredCorrection";
66  }
67 
69  std::string variant, unsigned int order,
70  std::vector<NekDouble> freeParams)
71  {
74  variant, order, freeParams);
75 
76  return p;
77  }
78 
79  static std::string className;
80 
81 protected:
82  LUE virtual void v_InitializeScheme(
83  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
84  const TimeIntegrationSchemeOperators &op) override;
85 
86  LUE virtual void v_ResidualEval(
87  const NekDouble &delta_t, const int n,
88  const TimeIntegrationSchemeOperators &op) override;
89 
90  LUE virtual void v_ResidualEval(
91  const NekDouble &delta_t,
92  const TimeIntegrationSchemeOperators &op) override;
93 
94  LUE virtual void v_ComputeInitialGuess(
95  const NekDouble &delta_t,
96  const TimeIntegrationSchemeOperators &op) override;
97 
98  LUE virtual void v_SDCIterationLoop(
99  const NekDouble &delta_t,
100  const TimeIntegrationSchemeOperators &op) override;
101 
102 }; // end class ImplicitTimeIntegrationSchemeSDC
103 
104 /**
105  * @brief Worker method to initialize the integration scheme.
106  */
108  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
110 {
111  TimeIntegrationSchemeSDC::v_InitializeScheme(deltaT, y_0, time, op);
112 }
113 
114 /**
115  * @brief Worker method to compute the residual.
116  */
118  const NekDouble &delta_t, const int n,
120 {
121  // Compute residual
122  op.DoOdeRhs(m_Y[n], m_F[n], m_time + delta_t * m_points[n]);
123 }
124 
126  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
127 {
128  // Compute residual
129  for (unsigned int n = 0; n < m_nQuadPts; ++n)
130  {
131  v_ResidualEval(delta_t, n, op);
132  }
133 }
134 
135 /**
136  * @brief Worker method to compute the initial SDC guess.
137  */
139  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
140 {
141  // Compute initial guess
142  for (unsigned int n = 0; n < m_nQuadPts - 1; ++n)
143  {
144  // Use implicit Euler as a first guess
145  NekDouble dtn = delta_t * (m_points[n + 1] - m_points[n]);
146  op.DoImplicitSolve(m_Y[n], m_Y[n + 1],
147  m_time + delta_t * m_points[n + 1], dtn);
148 
149  // Compute residual from updated solution
150  for (unsigned int i = 0; i < m_nvars; ++i)
151  {
152  Vmath::Vsub(m_npoints, m_Y[n + 1][i], 1, m_Y[n][i], 1,
153  m_F[n + 1][i], 1);
154  Vmath::Smul(m_npoints, 1.0 / dtn, m_F[n + 1][i], 1, m_F[n + 1][i],
155  1);
156  }
157  }
158 }
159 
160 /**
161  * @brief Worker method to compute the SDC iteration.
162  */
164  const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
165 {
166  unsigned int kstart = 1;
167  for (unsigned int k = kstart; k < m_order; ++k)
168  {
169  // Update integrated residual
170  UpdateIntegratedResidual(delta_t);
171 
172  // Loop over quadrature points
173  for (unsigned int n = 0; n < m_nQuadPts - 1; ++n)
174  {
175  NekDouble dtn = delta_t * (m_points[n + 1] - m_points[n]);
176 
177  // Update solution
178  for (unsigned int i = 0; i < m_nvars; ++i)
179  {
180  // Add Fint contribution
181  Vmath::Vadd(m_npoints, m_Y[n][i], 1, m_Fint[n][i], 1, m_tmp[i],
182  1);
183 
184  // Add implicit contribution
185  Vmath::Svtvp(m_npoints, -m_theta * dtn, m_F[n + 1][i], 1,
186  m_tmp[i], 1, m_tmp[i], 1);
187  }
188 
189  // Solve implicit system
190  op.DoImplicitSolve(m_tmp, m_Y[n + 1],
191  m_time + delta_t * m_points[n + 1],
192  m_theta * dtn);
193 
194  // Compute residual from updated solution
195  for (unsigned int i = 0; i < m_nvars; ++i)
196  {
197  Vmath::Vsub(m_npoints, m_Y[n + 1][i], 1, m_tmp[i], 1,
198  m_F[n + 1][i], 1);
199  Vmath::Smul(m_npoints, 1.0 / (m_theta * dtn), m_F[n + 1][i], 1,
200  m_F[n + 1][i], 1);
201  }
202  }
203  }
204 }
205 
206 } // end namespace LibUtilities
207 } // end namespace Nektar
208 
209 #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
Worker method to initialize the integration scheme.
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.
ImplicitTimeIntegrationSchemeSDC(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op) override
Worker method to compute the initial SDC guess.
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 DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) 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.
@ 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: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 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:248
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:419