Nektar++
TimeIntegrationAlgorithmGLM.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: TimeIntegrationAlgorithmGLM.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2019 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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Header file of time integration scheme data class
33 //
34 // The TimeIntegrationAlgorithmGLM class should only be used by the
35 // TimeIntegrationSchemeGLM class.
36 //
37 ///////////////////////////////////////////////////////////////////////////////
38 
39 #ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_ALGORITHM_GLM
40 #define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_ALGORITHM_GLM
41 
42 #define LUE LIB_UTILITIES_EXPORT
43 
46 
47 namespace Nektar
48 {
49 namespace LibUtilities
50 {
51 
53 {
54 public:
56  : m_parent(parent)
57  {
58  }
59 
61  {
62  }
63 
64  /**
65  * \brief This function initialises the time integration scheme.
66  *
67  * Given the solution at the initial time level \f$\boldsymbol{y}(t_0)\f$,
68  * this function generates the vectors \f$\boldsymbol{y}^{[n]}\f$ and
69  * \f$t^{[n]}\f$ needed to evaluate the time integration scheme formulated
70  * as a General Linear Method. These vectors are embedded in an object of
71  * the class TimeIntegrationSolutionGLM. This class is the abstraction of
72  * the input and output vectors of the General Linear Method.
73  *
74  * For single-step methods, this function is trivial as it just wraps a
75  * TimeIntegrationSolutionGLM object around the given initial values and
76  * initial time. However, for multistep methods, actual time stepping is
77  * being done to evaluate the necessary parameters at multiple time levels
78  * needed to start the actual integration.
79  *
80  * \param timestep The size of the timestep, i.e. \f$\Delta t\f$.
81  * \param time on input: the initial time, i.e. \f$t_0\f$.
82  * \param time on output: the new time-level after initialisation, in
83  * general this yields
84  * \f$t = t_0 + (r-1)\Delta t\f$ where \f$r\f$ is the number of steps
85  * of the multi-step method.
86  * \param nsteps on output: he number of initialisation steps required. In
87  * general this corresponds to \f$r-1\f$
88  * where \f$r\f$ is the number of steps of the multi-step method.
89  * \param f an object of the class FuncType, where FuncType should have a
90  * method FuncType::ODEforcing
91  * to evaluate the right hand side \f$f(t,\boldsymbol{y})\f$ of the
92  * ODE.
93  * \param y_0 the initial value \f$\boldsymbol{y}(t_0)\f$
94  * \return An object of the class TimeIntegrationSolutionGLM which
95  * represents the vectors \f$\boldsymbol{y}^{[n]}\f$ and \f$t^{[n]}\f$ that
96  * can be used to start the actual integration.
97  */
99  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
101 
102  /**
103  * \brief Explicit integration of an ODE.
104  *
105  * This function explicitely perfroms a single integration step of the ODE
106  * system:
107  * \f[
108  * \frac{d\boldsymbol{y}}{dt}=\boldsymbol{f}(t,\boldsymbol{y})
109  * \f]
110  *
111  * \param deltaT The size of the timestep, i.e. \f$\Delta t\f$.
112  * \param f an object of the class FuncType, where FuncType should have a
113  * method FuncType::ODEforcing
114  * to evaluate the right hand side \f$f(t,\boldsymbol{y})\f$ of the
115  * ODE.
116  * \param y on input: the vectors \f$\boldsymbol{y}^{[n-1]}\f$ and
117  * \f$t^{[n-1]}\f$ (which corresponds to the
118  * solution at the old time level)
119  * \param y on output: the vectors \f$\boldsymbol{y}^{[n]}\f$ and
120  * \f$t^{[n]}\f$ (which corresponds to the
121  * solution at the old new level)
122  * \return The actual solution \f$\boldsymbol{y}^{n}\f$ at the new time
123  * level
124  * (which in fact is also embedded in the argument y).
125  */
129 
130  // Does the actual multi-stage multi-step integration.
131  LUE void TimeIntegrate(const NekDouble deltaT, ConstTripleArray &y_old,
132  ConstSingleArray &t_old, TripleArray &y_new,
133  SingleArray &t_new,
135 
137  {
138  return m_schemeType;
139  }
140 
142  {
146  };
147 
148  inline unsigned int GetNmultiStepValues() const
149  {
150  return m_numMultiStepValues;
151  }
152 
153  inline unsigned int GetNmultiStepImplicitDerivs() const
154  {
156  }
157 
158  inline unsigned int GetNmultiStepDerivs() const
159  {
160  return m_numMultiStepDerivs;
161  }
162 
164  {
165  return m_timeLevelOffset;
166  }
167 
168  // Variables - all public for easy access when setting up the phase.
169  /// Parent scheme object
171 
172  std::string m_name;
173  std::string m_variant;
174  unsigned int m_order{0};
175  std::vector<NekDouble> m_freeParams;
176 
177  // Type of time integration scheme (Explicit, Implicit, IMEX, etc)
179 
180  unsigned int m_numstages{0}; //< Number of stages in multi-stage component.
181  unsigned int m_numsteps{0}; //< Number of steps in this integration phase
182 
183  unsigned int m_numMultiStepValues{0}; // number of entries in input and
184  // output vector that correspond
185  // to VALUES at previous time levels
187  0}; // number of entries in input and
188  // output vector that correspond
189  // to implicit DERIVATIVES at previous levels
190  unsigned int m_numMultiStepDerivs{0}; // number of entries in input and
191  // output vector that correspond
192  // to DERIVATIVES at previous levels
194  m_timeLevelOffset; // denotes to which time-level the entries in both
195  // input and output vector correspond, e.g.
196  // INPUT VECTOR --------> m_inputTimeLevelOffset
197  // _ _ _ _
198  // | u^n | | 0 |
199  // | u^{n-1} | | 1 |
200  // | u^{n-2} | -----> | 2 |
201  // | dt f(u^{n-1})| | 1 |
202  // | dt f(u^{n-2})| | 2 |
203  // - - - -
204 
209 
210  // Arrays used for the exponential integrators.
212 
215 
218 
219  bool m_initialised{false}; /// bool to identify array initialization
220 
221  // Values uses for exponential integration
222  NekDouble m_lastDeltaT{0}; /// Last delta T
223  NekDouble m_lastNVars{0}; /// Last number of vars
224 
225  int m_nvars; ///< The number of variables in integration scheme.
226  int m_npoints; ///< The size of inner data which is stored for reuse.
227 
228  // Friend classes
229  LUE friend std::ostream &operator<<(std::ostream &os,
230  const TimeIntegrationScheme &rhs);
231  LUE friend std::ostream &operator<<(
232  std::ostream &os, const TimeIntegrationSchemeSharedPtr &rhs);
233 
234  LUE friend std::ostream &operator<<(std::ostream &os,
235  const TimeIntegrationAlgorithmGLM &rhs);
236  LUE friend std::ostream &operator<<(
237  std::ostream &os, const TimeIntegrationAlgorithmGLMSharedPtr &rhs);
238 
239 private:
240  DoubleArray m_Y; /// Array containing the stage values
241  DoubleArray m_tmp; /// Explicit RHS of each stage equation
242 
243  TripleArray m_F; /// Array corresponding to the stage Derivatives
244  TripleArray m_F_IMEX; /// Used to store the Explicit stage derivative of
245  /// IMEX schemes
246 
247  NekDouble m_T{0}; /// Time at the different stages
248 
249  bool m_firstStageEqualsOldSolution{false}; //< Optimisation-flag
250  bool m_lastStageEqualsNewSolution{false}; //< Optimisation-flag
251 
255 
256  bool CheckTimeIntegrateArguments( // const NekDouble timestep,
257  ConstTripleArray &y_old, ConstSingleArray &t_old, TripleArray &y_new,
258  SingleArray &t_new, const TimeIntegrationSchemeOperators &op) const;
259 
260  // Helpers
261  inline NekDouble A(const unsigned int i, const unsigned int j) const
262  {
263  return m_A[0][i][j];
264  }
265  // Overloaded accessor for GLM parameter A; both exponential and
266  // conventional schemes are currently supported
267  inline NekDouble A(const unsigned int k, const unsigned int i,
268  const unsigned int j) const
269  {
270  if (m_schemeType == eExponential)
271  {
272  // For exponential schemes, the parameter A is specific for
273  // each variable k
274  return m_A_phi[k][i][j];
275  }
276  else
277  {
278  return m_A[0][i][j];
279  }
280  }
281  inline NekDouble B(const unsigned int i, const unsigned int j) const
282  {
283  return m_B[0][i][j];
284  }
285  // Overloaded accessor for GLM parameter B; both exponential and
286  // conventional schemes are currently supported
287  inline NekDouble B(const unsigned int k, const unsigned int i,
288  const unsigned int j) const
289  {
290  if (m_schemeType == eExponential)
291  {
292  // For exponential schemes, the parameter B is specific for
293  // each variable k
294  return m_B_phi[k][i][j];
295  }
296  else
297  {
298  return m_B[0][i][j];
299  }
300  }
301  inline NekDouble U(const unsigned int i, const unsigned int j) const
302  {
303  return m_U[i][j];
304  }
305  // Overloaded accessor for GLM parameter U; both exponential and
306  // conventional schemes are currently supported
307  inline NekDouble U(const unsigned int k, const unsigned int i,
308  const unsigned int j) const
309  {
310  if (m_schemeType == eExponential)
311  {
312  // For exponential schemes, the parameter U is specific for
313  // each variable k
314  return m_U_phi[k][i][j];
315  }
316  else
317  {
318  return m_U[i][j];
319  }
320  }
321  inline NekDouble V(const unsigned int i, const unsigned int j) const
322  {
323  return m_V[i][j];
324  }
325  // Overloaded accessor for GLM parameter V; both exponential and
326  // conventional schemes are currently supported
327  inline NekDouble V(const unsigned int k, const unsigned int i,
328  const unsigned int j) const
329  {
330  if (m_schemeType == eExponential)
331  {
332  // For exponential schemes, the parameter V is specific for
333  // each variable k
334  return m_V_phi[k][i][j];
335  }
336  else
337  {
338  return m_V[i][j];
339  }
340  }
341 
342  inline NekDouble A_IMEX(const unsigned int i, const unsigned int j) const
343  {
344  return m_A[1][i][j];
345  }
346  inline NekDouble B_IMEX(const unsigned int i, const unsigned int j) const
347  {
348  return m_B[1][i][j];
349  }
350 
351  inline unsigned int GetNsteps(void) const
352  {
353  return m_numsteps;
354  }
355 
356  inline unsigned int GetNstages(void) const
357  {
358  return m_numstages;
359  }
360 
361  inline int GetFirstDim(ConstTripleArray &y) const
362  {
363  return y[0].size();
364  }
365  inline int GetSecondDim(ConstTripleArray &y) const
366  {
367  return y[0][0].size();
368  }
369 
370 }; // end class TimeIntegrationAlgorithmGLM
371 
372 LUE std::ostream &operator<<(std::ostream &os,
373  const TimeIntegrationAlgorithmGLM &rhs);
374 
375 LUE std::ostream &operator<<(std::ostream &os,
377 
378 // =========================================================================
379 
380 } // end of namespace LibUtilities
381 } // end of namespace Nektar
382 
383 #endif
NekDouble A(const unsigned int i, const unsigned int j) const
NekDouble U(const unsigned int k, const unsigned int i, const unsigned int j) const
NekDouble A(const unsigned int k, const unsigned int i, const unsigned int j) const
NekDouble B(const unsigned int k, const unsigned int i, const unsigned int j) const
TripleArray m_F
Explicit RHS of each stage equation.
bool CheckTimeIntegrateArguments(ConstTripleArray &y_old, ConstSingleArray &t_old, TripleArray &y_new, SingleArray &t_new, const TimeIntegrationSchemeOperators &op) const
int m_npoints
The size of inner data which is stored for reuse.
LUE friend std::ostream & operator<<(std::ostream &os, const TimeIntegrationScheme &rhs)
DoubleArray m_tmp
Array containing the stage values.
NekDouble B_IMEX(const unsigned int i, const unsigned int j) const
NekDouble V(const unsigned int i, const unsigned int j) const
NekDouble V(const unsigned int k, const unsigned int i, const unsigned int j) const
const TimeIntegrationScheme * m_parent
Parent scheme object.
NekDouble m_T
Used to store the Explicit stage derivative of IMEX schemes.
const Array< OneD, const unsigned int > & GetTimeLevelOffset() const
NekDouble A_IMEX(const unsigned int i, const unsigned int j) const
NekDouble B(const unsigned int i, const unsigned int j) const
TimeIntegrationAlgorithmGLM(const TimeIntegrationScheme *parent)
LUE ConstDoubleArray & TimeIntegrate(const NekDouble deltaT, TimeIntegrationSolutionGLMSharedPtr &y, const TimeIntegrationSchemeOperators &op)
Explicit integration of an ODE.
LUE TimeIntegrationSolutionGLMSharedPtr InitializeData(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)
This function initialises the time integration scheme.
NekDouble m_lastDeltaT
bool to identify array initialization
NekDouble U(const unsigned int i, const unsigned int j) const
TripleArray m_F_IMEX
Array corresponding to the stage Derivatives.
Base class for time integration schemes.
Binds a set of functions for use by time integration schemes.
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
@ eExponential
Exponential scheme.
std::shared_ptr< TimeIntegrationSolutionGLM > TimeIntegrationSolutionGLMSharedPtr
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble