Nektar++
TimeIntegrationScheme.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: TimeIntegrationScheme.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 base class
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_LIB_UTILITIES_FOUNDATIONS_TIMEINTEGRATIONSCHEME_H
36 #define NEKTAR_LIB_UTILITIES_FOUNDATIONS_TIMEINTEGRATIONSCHEME_H
37 
38 #include <boost/core/ignore_unused.hpp>
39 
45 
46 namespace Nektar
47 {
48  namespace LibUtilities
49  {
50  // Forward declaration of all classes in this file
51  class TimeIntegrationScheme;
53 
54  // typedefs
55  typedef std::shared_ptr<TimeIntegrationScheme> TimeIntegrationSchemeSharedPtr;
56  typedef std::vector<TimeIntegrationSchemeSharedPtr> TimeIntegrationSchemeVector;
57  typedef std::vector<TimeIntegrationSchemeSharedPtr>::iterator TimeIntegrationSchemeVectorIter;
58  typedef std::shared_ptr<TimeIntegrationSolution> TimeIntegrationSolutionSharedPtr;
59  typedef std::vector<TimeIntegrationSolutionSharedPtr> TimeIntegrationSolutionVector;
60  typedef std::vector<TimeIntegrationSolutionSharedPtr>::iterator TimeIntegrationSolutionVectorIter;
61 
62  // =========================================================================
63  // ==== ENUM LIST OF ALL SUPPORTED INTEGRATION SCHEMES
64  // =========================================================================
66  {
68  eAdamsBashforthOrder1, //!< Adams-Bashforth Forward multi-step scheme of order 1
69  eAdamsBashforthOrder2, //!< Adams-Bashforth Forward multi-step scheme of order 2
70  eAdamsBashforthOrder3, //!< Adams-Bashforth Forward multi-step scheme of order 3
71  eAdamsBashforthOrder4, //!< Adams-Bashforth Forward multi-step scheme of order 4
72  eAdamsMoultonOrder1, //!< Adams-Moulton Forward multi-step scheme of order 1
73  eAdamsMoultonOrder2, //!< Adams-Moulton Forward multi-step scheme of order 2
74  eBDFImplicitOrder1, //!< BDF multi-step scheme of order 1 (implicit)
75  eBDFImplicitOrder2, //!< BDF multi-step scheme of order 2 (implicit)
76  eClassicalRungeKutta4, //!< Runge-Kutta multi-stage scheme 4th order explicit (old name)
77  eRungeKutta4, //!< Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
78  eRungeKutta5, //!< RungeKutta5 method
79  eRungeKutta3_SSP, //!< Nonlinear SSP RungeKutta3 explicit
80  eRungeKutta2_ImprovedEuler, //!< Improved RungeKutta2 explicit (old name meaning Heun's method)
81  eRungeKutta2_SSP, //!< Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
82  eForwardEuler, //!< Forward Euler scheme
83  eBackwardEuler, //!< Backward Euler scheme
84  eIMEXOrder1, //!< IMEX 1st order scheme using Euler Backwards/Euler Forwards
85  eIMEXOrder2, //!< IMEX 2nd order scheme using Backward Different Formula & Extrapolation
86  eIMEXOrder3, //!< IMEX 3rd order scheme using Backward Different Formula & Extrapolation
87  eIMEXOrder4, //!< IMEX 4th order scheme using Backward Different Formula & Extrapolation
88  eMidpoint, //!< midpoint method (old name)
89  eRungeKutta2, //!< Classical RungeKutta2 method (new name for eMidpoint)
90  eDIRKOrder2, //!< Diagonally Implicit Runge Kutta scheme of order 2
91  eDIRKOrder3, //!< Diagonally Implicit Runge Kutta scheme of order 3
92  eCNAB, //!< Crank-Nicolson/Adams-Bashforth Order 2 (CNAB)
93  eIMEXGear, //!< IMEX Gear Order 2
94  eMCNAB, //!< Modified Crank-Nicolson/Adams-Bashforth Order 2 (MCNAB)
95  eIMEXdirk_1_1_1, //!< Forward-Backward Euler IMEX DIRK(1,1,1)
96  eIMEXdirk_1_2_1, //!< Forward-Backward Euler IMEX DIRK(1,2,1)
97  eIMEXdirk_1_2_2, //!< Implicit-Explicit Midpoint IMEX DIRK(1,2,2)
98  eIMEXdirk_2_2_2, //!< L-stable, two stage, second order IMEX DIRK(2,2,2)
99  eIMEXdirk_2_3_2, //!< L-stable, three stage, third order IMEX DIRK(3,4,3)
100  eIMEXdirk_2_3_3, //!< L-stable, two stage, third order IMEX DIRK(2,3,3)
101  eIMEXdirk_3_4_3, //!< L-stable, three stage, third order IMEX DIRK(3,4,3)
102  eIMEXdirk_4_4_3, //!< L-stable, four stage, third order IMEX DIRK(4,4,3)
103  SIZE_TimeIntegrationMethod //!< Length of enum list
104  };
105 
106  const char* const TimeIntegrationMethodMap[] =
107  {
108  "NoTimeIntegrationMethod",
109  "AdamsBashforthOrder1",
110  "AdamsBashforthOrder2",
111  "AdamsBashforthOrder3",
112  "AdamsBashforthOrder4",
113  "AdamsMoultonOrder1",
114  "AdamsMoultonOrder2",
115  "BDFImplicitOrder1",
116  "BDFImplicitOrder2",
117  "ClassicalRungeKutta4",
118  "RungeKutta4",
119  "RungeKutta5",
120  "RungeKutta3_SSP",
121  "RungeKutta2_ImprovedEuler",
122  "RungeKutta2_SSP",
123  "ForwardEuler",
124  "BackwardEuler",
125  "IMEXOrder1",
126  "IMEXOrder2",
127  "IMEXOrder3",
128  "IMEXOrder4",
129  "Midpoint",
130  "RungeKutta2",
131  "DIRKOrder2",
132  "DIRKOrder3",
133  "CNAB",
134  "IMEXGear",
135  "MCNAB",
136  "IMEXdirk_1_1_1",
137  "IMEXdirk_1_2_1",
138  "IMEXdirk_1_2_2",
139  "IMEXdirk_2_2_2",
140  "IMEXdirk_2_3_2",
141  "IMEXdirk_2_3_3",
142  "IMEXdirk_3_4_3",
143  "IMEXdirk_4_4_3",
144  };
145 
147  {
149  eExplicit, //!< Formally explicit scheme
150  eDiagonallyImplicit, //!< Diagonally implicit scheme (e.g. the DIRK schemes)
151  eIMEX, //!< Implicit Explicit General Linear Method
152  eImplicit, //!< Fully implicit scheme
153  };
154 
155  const char* const TimeIntegrationSchemeTypeMap[] =
156  {
157  "NoTimeIntegrationSchemeType",
158  "Explicit",
159  "DiagonallyImplicit",
160  "IMEX",
161  "Implicit"
162  };
163 
164  // =====================================================================
165 
166  // =====================================================================
167  // ==== DEFINITION OF THE CLASS TimeIntegrationSchemeOperators
168  // =====================================================================
170  {
171  public:
174 
175  typedef std::function< void (InArrayType&, OutArrayType&, const NekDouble) > FunctorType1;
176  typedef std::function< void (InArrayType&, OutArrayType&, const NekDouble, const NekDouble) > FunctorType2;
177 
178  typedef const FunctorType1& ConstFunctorType1Ref;
179  typedef const FunctorType2& ConstFunctorType2Ref;
180 
183 
185  m_functors1(4),
186  m_functors2(1)
187  {
188  }
189 
190  template<typename FuncPointerT, typename ObjectPointerT>
191  void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
192  {
193  m_functors1[0] = std::bind(
194  func, obj, std::placeholders::_1, std::placeholders::_2,
195  std::placeholders::_3);
196  }
197 
198  template<typename FuncPointerT, typename ObjectPointerT>
199  void DefineOdeExplicitRhs(FuncPointerT func, ObjectPointerT obj)
200  {
201  m_functors1[1] = std::bind(
202  func, obj, std::placeholders::_1, std::placeholders::_2,
203  std::placeholders::_3);
204  }
205 
206  template<typename FuncPointerT, typename ObjectPointerT>
207  void DefineOdeImplicitRhs(FuncPointerT func, ObjectPointerT obj)
208  {
209  m_functors1[2] = std::bind(
210  func, obj, std::placeholders::_1, std::placeholders::_2,
211  std::placeholders::_3);
212  }
213 
214  template<typename FuncPointerT, typename ObjectPointerT>
215  void DefineProjection(FuncPointerT func, ObjectPointerT obj)
216  {
217  m_functors1[3] = std::bind(
218  func, obj, std::placeholders::_1, std::placeholders::_2,
219  std::placeholders::_3);
220  }
221 
222  template<typename FuncPointerT, typename ObjectPointerT>
223  void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
224  {
225  m_functors2[0] = std::bind(
226  func, obj, std::placeholders::_1, std::placeholders::_2,
227  std::placeholders::_3, std::placeholders::_4);
228  }
229 
230 
231  inline void DoOdeRhs(InArrayType &inarray,
232  OutArrayType &outarray,
233  const NekDouble time) const
234  {
235  ASSERTL1(m_functors1[0],"OdeRhs should be defined for this time integration scheme");
236  m_functors1[0](inarray,outarray,time);
237  }
238 
239  inline void DoOdeExplicitRhs(InArrayType &inarray,
240  OutArrayType &outarray,
241  const NekDouble time) const
242  {
243  ASSERTL1(m_functors1[1],"OdeExplicitRhs should be defined for this time integration scheme");
244  m_functors1[1](inarray,outarray,time);
245  }
246 
247  inline void DoOdeImplicitRhs(InArrayType &inarray,
248  OutArrayType &outarray,
249  const NekDouble time) const
250  {
251  ASSERTL1(m_functors1[2],"OdeImplictRhs should be defined for this time integration scheme");
252  m_functors1[2](inarray,outarray,time);
253  }
254 
255  inline void DoProjection(InArrayType &inarray,
256  OutArrayType &outarray,
257  const NekDouble time) const
258  {
259  ASSERTL1(m_functors1[3],"Projection operation should be defined for this time integration scheme");
260  m_functors1[3](inarray,outarray,time);
261  }
262 
263  inline void DoImplicitSolve(InArrayType &inarray,
264  OutArrayType &outarray,
265  const NekDouble time,
266  const NekDouble lambda) const
267  {
268  ASSERTL1(m_functors2[0],"ImplicitSolve should be defined for this time integration scheme");
269  m_functors2[0](inarray,outarray,time,lambda);
270  }
271 
272 
273  protected:
274  FunctorType1Array m_functors1;
275  FunctorType2Array m_functors2;
276 
277  private:
278 
279  };
280 
281  // =========================================================================
282  // ==== DEFINITION OF THE CLASS TimeIntegrationSchemeKey
283  // =========================================================================
285  {
286  public:
287 
288  struct opLess
289  {
290  LIB_UTILITIES_EXPORT bool operator()(const TimeIntegrationSchemeKey &lhs, const TimeIntegrationSchemeKey &rhs) const;
291  };
292 
294  m_method(method)
295  {
296  }
297 
299  {
300  }
301 
303  {
304  *this = key; // defer to assignment operator
305  }
306 
308  {
309  m_method = key.m_method;
310 
311  return *this;
312  }
313 
315  {
316  return m_method;
317  }
318 
319  inline bool operator==(const TimeIntegrationSchemeKey &key)
320  {
321  return (m_method == key.m_method);
322  }
323 
324  inline bool operator== (const TimeIntegrationSchemeKey *y)
325  {
326  return (*this == *y);
327  }
328 
330  {
331  return (!(*this == y));
332  }
333 
335  {
336  return (!(*this == *y));
337  }
338 
341  LIB_UTILITIES_EXPORT friend bool opLess::operator()(const TimeIntegrationSchemeKey &lhs, const TimeIntegrationSchemeKey &rhs) const;
342 
343  protected:
344  TimeIntegrationMethod m_method; //!< integration method
345 
346  private:
347  // This should never be called
349  {
350  NEKERROR(ErrorUtil::efatal,"Default Constructor for the TimeIntegrationSchemeKey class should not be called");
351  }
352  };
353 
357  LIB_UTILITIES_EXPORT std::ostream& operator<<(std::ostream& os, const TimeIntegrationSchemeKey& rhs);
358 
359  // =====================================================================
360 
361 
362  // =====================================================================
363  // ==== DEFINITION OF THE FUNCTION TimeIntegrationSolution
364  // =====================================================================
365  // This is the interface the user can use to get hold of the different
366  // time integration schemes. It returns you the NekManager
367  // singleton which manages all the time integration schemes...
368  //
372  LIB_UTILITIES_EXPORT TimeIntegrationSchemeManagerT &TimeIntegrationSchemeManager(void);
373  // =====================================================================
374 
375  // =====================================================================
376  // ==== DEFINITION OF THE CLASS TimeIntegrationScheme
377  // =====================================================================
378  class TimeIntegrationScheme
379  {
380  public:
381 
388  typedef std::function< void (ConstDoubleArray&, DoubleArray&, const NekDouble) > FunctorType1;
389  typedef std::function< void (ConstDoubleArray&, DoubleArray&, const NekDouble, const NekDouble) > FunctorType2;
390 
391  public:
392 
394  {
395  }
396 
397  inline const TimeIntegrationSchemeKey& GetIntegrationSchemeKey() const
398  {
399  return m_schemeKey;
400  }
401 
403  {
404  return m_schemeKey.GetIntegrationMethod();
405  }
406 
408  {
409  return m_schemeType;
410  }
411 
412  inline NekDouble A(const unsigned int i, const unsigned int j) const
413  {
414  return m_A[0][i][j];
415  }
416 
417  inline NekDouble B(const unsigned int i, const unsigned int j) const
418  {
419  return m_B[0][i][j];
420  }
421 
422  inline NekDouble U(const unsigned int i, const unsigned int j) const
423  {
424  return m_U[i][j];
425  }
426 
427  inline NekDouble V(const unsigned int i, const unsigned int j) const
428  {
429  return m_V[i][j];
430  }
431 
432  inline NekDouble A_IMEX(const unsigned int i, const unsigned int j) const
433  {
434  return m_A[1][i][j];
435  }
436 
437  inline NekDouble B_IMEX(const unsigned int i, const unsigned int j) const
438  {
439  return m_B[1][i][j];
440  }
441 
442  inline unsigned int GetNsteps(void) const
443  {
444  return m_numsteps;
445  }
446 
447  inline unsigned int GetNstages(void) const
448  {
449  return m_numstages;
450  }
451 
452  inline unsigned int GetNmultiStepValues(void) const
453  {
454  return m_numMultiStepValues;
455  }
456 
457  inline unsigned int GetNmultiStepDerivs(void) const
458  {
459  return m_numMultiStepDerivs;
460  }
461 
462  /**
463  * \brief This function initialises the time integration
464  * scheme
465  *
466  * Given the solution at the initial time level
467  * \f$\boldsymbol{y}(t_0)\f$, this function generates the
468  * vectors \f$\boldsymbol{y}^{[n]}\f$ and \f$t^{[n]}\f$
469  * needed to evaluate the time integration scheme
470  * formulated as a General Linear Method. These vectors
471  * are embedded in an object of the class
472  * TimeIntegrationSolution. This class is the abstraction
473  * of the input and output vectors of the General Linear
474  * Method.
475  *
476  * For single-step methods, this function is trivial as it
477  * just wraps a TimeIntegrationSolution object around the
478  * given initial values and initial time. However, for
479  * multistep methods, actual time stepping is being done
480  * to evaluate the necessary parameters at multiple time
481  * levels needed to start the actual integration.
482  *
483  * \param timestep The size of the timestep, i.e. \f$\Delta t\f$.
484  * \param time on input: the initial time, i.e. \f$t_0\f$.
485  * \param time on output: the new time-level after initialisation, in general this yields
486  * \f$t = t_0 + (r-1)\Delta t\f$ where \f$r\f$ is the number of steps of the multi-step method.
487  * \param nsteps on output: he number of initialisation steps required. In general this corresponds to \f$r-1\f$
488  * where \f$r\f$ is the number of steps of the multi-step method.
489  * \param f an object of the class FuncType, where FuncType should have a method FuncType::ODEforcing
490  * to evaluate the right hand side \f$f(t,\boldsymbol{y})\f$ of the ODE.
491  * \param y_0 the initial value \f$\boldsymbol{y}(t_0)\f$
492  * \return An object of the class TimeIntegrationSolution which represents the vectors \f$\boldsymbol{y}^{[n]}\f$ and \f$t^{[n]}\f$
493  * that can be used to start the actual integration.
494  */
495  LIB_UTILITIES_EXPORT TimeIntegrationSolutionSharedPtr InitializeScheme(const NekDouble timestep,
496  ConstDoubleArray &y_0 ,
497  const NekDouble time ,
498  const TimeIntegrationSchemeOperators &op );
499 
500  /**
501  * \brief Explicit integration of an ODE.
502  *
503  * This function explicitely perfroms a signle integration step of the ODE system:
504  * \f[
505  * \frac{d\boldsymbol{y}}{dt}=\boldsymbol{f}(t,\boldsymbol{y})
506  * \f]
507  *
508  * \param timestep The size of the timestep, i.e. \f$\Delta t\f$.
509  * \param f an object of the class FuncType, where FuncType should have a method FuncType::ODEforcing
510  * to evaluate the right hand side \f$f(t,\boldsymbol{y})\f$ of the ODE.
511  * \param y on input: the vectors \f$\boldsymbol{y}^{[n-1]}\f$ and \f$t^{[n-1]}\f$ (which corresponds to the
512  * solution at the old time level)
513  * \param y on output: the vectors \f$\boldsymbol{y}^{[n]}\f$ and \f$t^{[n]}\f$ (which corresponds to the
514  * solution at the old new level)
515  * \return The actual solution \f$\boldsymbol{y}^{n}\f$ at the new time level
516  * (which in fact is also embedded in the argument y).
517  */
518  LIB_UTILITIES_EXPORT ConstDoubleArray& TimeIntegrate(const NekDouble timestep,
519  TimeIntegrationSolutionSharedPtr &y ,
520  const TimeIntegrationSchemeOperators &op );
521 
523  {
524  return m_timeLevelOffset;
525  }
526 
527 
528  protected:
529  TimeIntegrationSchemeKey m_schemeKey;
531  unsigned int m_numsteps; //< Number of steps in multi-step component.
532  unsigned int m_numstages; //< Number of stages in multi-stage component.
533 
534  bool m_firstStageEqualsOldSolution; //< Optimisation-flag
535  bool m_lastStageEqualsNewSolution; //< Optimisation-flag
536 
537  unsigned int m_numMultiStepValues; // number of entries in input and output vector that correspond
538  // to VALUES at previous time levels
539  unsigned int m_numMultiStepDerivs; // number of entries in input and output vector that correspond
540  // to DERIVATIVES at previous time levels
541  Array<OneD,unsigned int> m_timeLevelOffset; // denotes to which time-level the entries in both
542  // input and output vector correspond, e.g.
543  // INPUT VECTOR --------> m_inputTimeLevelOffset
544  // _ _ _ _
545  // | u^n | | 0 |
546  // | u^{n-1} | | 1 |
547  // | u^{n-2} | -----> | 2 |
548  // | dt f(u^{n-1})| | 1 |
549  // | dt f(u^{n-2})| | 2 |
550  // - - - -
551 
556 
557  private:
558  bool m_initialised; /// bool to identify if array has been initialised
559  int m_nvar; /// The number of variables in integration scheme.
560  int m_npoints; /// The size of inner data which is stored for reuse.
561  DoubleArray m_Y; /// Array containing the stage values
562  DoubleArray m_tmp; /// explicit right hand side of each stage equation
563 
564  TripleArray m_F; /// Array corresponding to the stage Derivatives
565  TripleArray m_F_IMEX; /// Used to store the Explicit stage derivative of IMEX schemes
566 
567  NekDouble m_T; /// Time at the different stages
568 
569  template <typename> friend class Nektar::MemoryManager;
570  LIB_UTILITIES_EXPORT friend TimeIntegrationSchemeManagerT &TimeIntegrationSchemeManager(void);
571 
572  LIB_UTILITIES_EXPORT static std::shared_ptr<TimeIntegrationScheme> Create(const TimeIntegrationSchemeKey &key);
573 
574  TimeIntegrationScheme(const TimeIntegrationSchemeKey &key);
575 
576  // These should never be called
577  TimeIntegrationScheme():m_schemeKey(NullTimeIntegrationSchemeKey)
578  {
579  NEKERROR(ErrorUtil::efatal,"Default Constructor for the TimeIntegrationScheme class should not be called");
580  }
581 
582  TimeIntegrationScheme(const TimeIntegrationScheme &in):m_schemeKey(NullTimeIntegrationSchemeKey)
583  {
584  boost::ignore_unused(in);
585  NEKERROR(ErrorUtil::efatal,"Copy Constructor for the TimeIntegrationScheme class should not be called");
586  }
587 
588  LIB_UTILITIES_EXPORT bool VerifyIntegrationSchemeType(TimeIntegrationSchemeType type,
589  const Array<OneD, const Array<TwoD, NekDouble> >& A,
590  const Array<OneD, const Array<TwoD, NekDouble> >& B,
592  const Array<TwoD, const NekDouble>& V) const;
593 
594  LIB_UTILITIES_EXPORT void TimeIntegrate(const NekDouble timestep,
595  ConstTripleArray &y_old ,
596  ConstSingleArray &t_old ,
597  TripleArray &y_new ,
598  SingleArray &t_new ,
599  const TimeIntegrationSchemeOperators &op );
600 
601 
602  inline int GetFirstDim(ConstTripleArray &y) const
603  {
604  return y[0].num_elements();
605  }
606 
607  inline int GetSecondDim(ConstTripleArray &y) const
608  {
609  return y[0][0].num_elements();
610  }
611 
612  LIB_UTILITIES_EXPORT bool CheckTimeIntegrateArguments(const NekDouble timestep,
613  ConstTripleArray &y_old ,
614  ConstSingleArray &t_old ,
615  TripleArray &y_new ,
616  SingleArray &t_new ,
617  const TimeIntegrationSchemeOperators &op) const;
618 
619  LIB_UTILITIES_EXPORT bool CheckIfFirstStageEqualsOldSolution(const Array<OneD, const Array<TwoD, NekDouble> >& A,
620  const Array<OneD, const Array<TwoD, NekDouble> >& B,
622  const Array<TwoD, const NekDouble>& V) const;
623 
624  LIB_UTILITIES_EXPORT bool CheckIfLastStageEqualsNewSolution(const Array<OneD, const Array<TwoD, NekDouble> >& A,
625  const Array<OneD, const Array<TwoD, NekDouble> >& B,
627  const Array<TwoD, const NekDouble>& V) const;
628 
629 
630 
631  };
632 
633 
634  // =========================================================================
635  // ==== DEFINITION OF THE CLASS TimeIntegrationSolution
636  // =========================================================================
638  {
639  public:
642 
643  // Constructor for single step methods
644  LIB_UTILITIES_EXPORT TimeIntegrationSolution(const TimeIntegrationSchemeKey &key,
645  const DoubleArray& y,
646  const NekDouble time,
647  const NekDouble timestep);
648  // Constructor for multi-step methods
649  LIB_UTILITIES_EXPORT TimeIntegrationSolution(const TimeIntegrationSchemeKey &key,
650  const TripleArray& y,
651  const Array<OneD, NekDouble>& t);
652  LIB_UTILITIES_EXPORT TimeIntegrationSolution(const TimeIntegrationSchemeKey &key,
653  unsigned int nvar,
654  unsigned int npoints);
655  LIB_UTILITIES_EXPORT TimeIntegrationSolution(const TimeIntegrationSchemeKey &key);
656 
657  inline const TimeIntegrationSchemeSharedPtr& GetIntegrationScheme() const
658  {
659  return m_scheme;
660  }
661 
662  inline const TimeIntegrationSchemeKey& GetIntegrationSchemeKey() const
663  {
664  return m_scheme->GetIntegrationSchemeKey();
665  }
666 
668  {
669  return m_scheme->GetIntegrationMethod();
670  }
671 
672  inline const TripleArray& GetSolutionVector() const
673  {
674  return m_solVector;
675  }
676 
677  inline TripleArray& UpdateSolutionVector()
678  {
679  return m_solVector;
680  }
681 
682  inline const DoubleArray& GetSolution() const
683  {
684  return m_solVector[0];
685  }
686 
687  inline DoubleArray& UpdateSolution()
688  {
689  return m_solVector[0];
690  }
691 
693  {
694  return m_t;
695  }
696 
698  {
699  return m_t;
700  }
701 
702  inline NekDouble GetTime() const
703  {
704  return m_t[0];
705  }
706 
707  inline int GetNsteps()
708  {
709  return m_scheme->GetNsteps();
710  }
711 
712  inline int GetFirstDim() const
713  {
714  return m_solVector[0].num_elements();
715  }
716 
717  inline int GetSecondDim() const
718  {
719  return m_solVector[0][0].num_elements();
720  }
721 
722  // Return the number of entries in the solution vector that correspond
723  // to (multi-step) values
724  inline unsigned int GetNvalues(void) const
725  {
726  return m_scheme->GetNmultiStepValues();
727  }
728 
729  // Return the number of entries in the solution vector that correspond
730  // to (multi-step) derivatives
731  inline unsigned int GetNderivs(void) const
732  {
733  return m_scheme->GetNmultiStepDerivs();
734  }
735 
736  // Returns an array which indicates to which time-level the entries in the
737  // solution vector correspond
739  {
740  return m_scheme->GetTimeLevelOffset();
741  }
742 
743  // returns the entry in the solution vector which corresponds to
744  // the (multi-step) value at the time-level with specified offset
745  inline DoubleArray& GetValue(const unsigned int timeLevelOffset)
746  {
747  int nMultiStepVals = m_scheme->GetNmultiStepValues();
748  const Array<OneD, const unsigned int>& offsetvec = GetTimeLevelOffset();
749 
750  for(int i = 0; i < nMultiStepVals; i++)
751  {
752  if( timeLevelOffset == offsetvec[i] )
753  {
754  return m_solVector[i];
755  }
756  }
757  ASSERTL1(false,"The solution vector of this scheme does not contain a value at the requested time-level");
758  return m_solVector[0];
759  }
760 
761  // returns the entry in the solution vector which corresponds to
762  // the (multi-step) derivative at the time-level with specified offset
763  inline DoubleArray& GetDerivative(const unsigned int timeLevelOffset)
764  {
765  int nMultiStepVals = m_scheme->GetNmultiStepValues();
766  int size = m_scheme->GetNsteps();
767  const Array<OneD, const unsigned int>& offsetvec = GetTimeLevelOffset();
768 
769  for(int i = nMultiStepVals; i < size; i++)
770  {
771  if( timeLevelOffset == offsetvec[i] )
772  {
773  return m_solVector[i];
774  }
775  }
776  ASSERTL1(false,"The solution vector of this scheme does not contain a derivative at the requested time-level");
777  return m_solVector[0];
778  }
779 
780  // returns the time associated with the (multi-step) value at the time-level with the
781  // given offset
782  inline NekDouble GetValueTime(const unsigned int timeLevelOffset)
783  {
784  int nMultiStepVals = m_scheme->GetNmultiStepValues();
785  const Array<OneD, const unsigned int>& offsetvec = GetTimeLevelOffset();
786 
787  for(int i = 0; i < nMultiStepVals; i++)
788  {
789  if( timeLevelOffset == offsetvec[i] )
790  {
791  return m_t[i];
792  }
793  }
794  ASSERTL1(false,"The solution vector of this scheme does not contain a value at the requested time-level");
795  return m_t[0];
796  }
797 
798  // sets the (multi-step) value and time in the solution
799  // vector which corresponds to
800  // the value at the time-level with specified offset
801  inline void SetValue(const unsigned int timeLevelOffset, const DoubleArray& y, const NekDouble t)
802  {
803  int nMultiStepVals = m_scheme->GetNmultiStepValues();
804  const Array<OneD, const unsigned int>& offsetvec = GetTimeLevelOffset();
805 
806  for(int i = 0; i < nMultiStepVals; i++)
807  {
808  if( timeLevelOffset == offsetvec[i] )
809  {
810  m_solVector[i] = y;
811  m_t[i] = t;
812  return;
813  }
814  }
815  }
816 
817  // sets the (multi-step) derivative and time in the
818  // solution vector which corresponds to
819  // the derivative at the time-level with specified offset
820  inline void SetDerivative(const unsigned int timeLevelOffset, const DoubleArray& y, const NekDouble timestep)
821  {
822  int nMultiStepVals = m_scheme->GetNmultiStepValues();
823  int size = m_scheme->GetNsteps();
824  const Array<OneD, const unsigned int>& offsetvec = GetTimeLevelOffset();
825 
826  for(int i = nMultiStepVals; i < size; i++)
827  {
828  if( timeLevelOffset == offsetvec[i] )
829  {
830  m_solVector[i] = y;
831  m_t[i] = timestep;
832  return;
833  }
834  }
835  }
836 
837  // sets the soln Vector
838  inline void SetSolVector(const int Offset, const DoubleArray& y)
839  {
840  m_solVector[Offset] = y;
841  }
842 
843  // Rotate the solution vector
844  // (i.e. updating without calculating/inserting new values)
845  inline void RotateSolutionVector()
846  {
847  int i;
848  int nMultiStepVals = m_scheme->GetNmultiStepValues();
849  int size = m_scheme->GetNsteps();
850  for(i = (nMultiStepVals-1); i > 0; i--)
851  {
852  m_solVector[i] = m_solVector[i-1];
853  }
854 
855  for(i = (size-1); i > nMultiStepVals; i--)
856  {
857  m_solVector[i] = m_solVector[i-1];
858  }
859  }
860 
861 
862  private:
863  TimeIntegrationSchemeSharedPtr m_scheme;
864  TripleArray m_solVector;
866  };
867  // =========================================================================
868 
869 
870  LIB_UTILITIES_EXPORT std::ostream& operator<<(std::ostream& os, const TimeIntegrationScheme& rhs);
871  LIB_UTILITIES_EXPORT std::ostream& operator<<(std::ostream& os, const TimeIntegrationSchemeSharedPtr& rhs);
872 
873  // =========================================================================
874 
875  } // end of namespace
876 } // end of namespace
877 
878 #endif //NEKTAR_LIB_UTILITIES_FOUNDATIONS_TIMEINTEGRATIONSCHEME_H
NekDouble A(const unsigned int i, const unsigned int j) const
NekDouble GetValueTime(const unsigned int timeLevelOffset)
TimeIntegrationSchemeKey(const TimeIntegrationMethod &method)
std::vector< TimeIntegrationSchemeSharedPtr > TimeIntegrationSchemeVector
NekDouble B(const unsigned int i, const unsigned int j) const
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
TimeIntegrationSchemeKey & operator=(const TimeIntegrationSchemeKey &key)
BDF multi-step scheme of order 1 (implicit)
Adams-Bashforth Forward multi-step scheme of order 2.
Array< OneD, Array< TwoD, NekDouble > > m_B
static const TimeIntegrationSchemeKey NullTimeIntegrationSchemeKey(eNoTimeIntegrationMethod)
void DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) const
int m_nvar
bool to identify if array has been initialised
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
Runge-Kutta multi-stage scheme 4th order explicit (old name)
Implicit-Explicit Midpoint IMEX DIRK(1,2,2)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
std::function< void(InArrayType &, OutArrayType &, const NekDouble) > FunctorType1
const char *const TimeIntegrationSchemeTypeMap[]
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
DoubleArray m_tmp
Array containing the stage values.
NekDouble V(const unsigned int i, const unsigned int j) const
Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
void DefineOdeExplicitRhs(FuncPointerT func, ObjectPointerT obj)
const TimeIntegrationSchemeKey & GetIntegrationSchemeKey() const
L-stable, four stage, third order IMEX DIRK(4,4,3)
Forward-Backward Euler IMEX DIRK(1,2,1)
void DefineOdeImplicitRhs(FuncPointerT func, ObjectPointerT obj)
bool operator==(const TimeIntegrationSchemeKey &key)
Implicit Explicit General Linear Method.
const char *const TimeIntegrationMethodMap[]
NekManager< TimeIntegrationSchemeKey, TimeIntegrationScheme, TimeIntegrationSchemeKey::opLess > TimeIntegrationSchemeManagerT
Nonlinear SSP RungeKutta3 explicit.
TimeIntegrationMethod GetIntegrationMethod() const
std::vector< TimeIntegrationSolutionSharedPtr >::iterator TimeIntegrationSolutionVectorIter
DoubleArray m_Y
The size of inner data which is stored for reuse.
bool operator!=(const BasisKey &x, const BasisKey &y)
TimeIntegrationSchemeManagerT & TimeIntegrationSchemeManager(void)
Adams-Moulton Forward multi-step scheme of order 2.
std::vector< TimeIntegrationSchemeSharedPtr >::iterator TimeIntegrationSchemeVectorIter
Adams-Bashforth Forward multi-step scheme of order 3.
StandardMatrixTag & lhs
NekDouble A_IMEX(const unsigned int i, const unsigned int j) const
bool operator==(const BasisKey &x, const BasisKey &y)
void SetValue(const unsigned int timeLevelOffset, const DoubleArray &y, const NekDouble t)
std::function< void(ConstDoubleArray &, DoubleArray &, const NekDouble, const NekDouble) > FunctorType2
Crank-Nicolson/Adams-Bashforth Order 2 (CNAB)
const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > ConstTripleArray
DoubleArray & GetDerivative(const unsigned int timeLevelOffset)
IMEX 2nd order scheme using Backward Different Formula & Extrapolation.
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
TripleArray m_F
explicit right hand side of each stage equation
const Array< OneD, const NekDouble > & GetTimeVector() const
Diagonally implicit scheme (e.g. the DIRK schemes)
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
Classical RungeKutta2 method (new name for eMidpoint)
const Array< OneD, const Array< OneD, NekDouble > > InArrayType
Adams-Moulton Forward multi-step scheme of order 1.
void SetDerivative(const unsigned int timeLevelOffset, const DoubleArray &y, const NekDouble timestep)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
Adams-Bashforth Forward multi-step scheme of order 1.
NekDouble U(const unsigned int i, const unsigned int j) const
#define LIB_UTILITIES_EXPORT
std::function< void(ConstDoubleArray &, DoubleArray &, const NekDouble) > FunctorType1
Array< OneD, Array< OneD, NekDouble > > DoubleArray
L-stable, three stage, third order IMEX DIRK(3,4,3)
TimeIntegrationSchemeKey(const TimeIntegrationSchemeKey &key)
void DoProjection(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
int m_npoints
The number of variables in integration scheme.
Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
Improved RungeKutta2 explicit (old name meaning Heun&#39;s method)
Forward-Backward Euler IMEX DIRK(1,1,1)
Array< OneD, Array< OneD, NekDouble > > DoubleArray
double NekDouble
DoubleArray & GetValue(const unsigned int timeLevelOffset)
void SetSolVector(const int Offset, const DoubleArray &y)
NekDouble B_IMEX(const unsigned int i, const unsigned int j) const
bool operator<(const BasisKey &lhs, const BasisKey &rhs)
std::function< void(InArrayType &, OutArrayType &, const NekDouble, const NekDouble) > FunctorType2
TimeIntegrationSchemeType GetIntegrationSchemeType() const
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > TripleArray
TimeIntegrationScheme(const TimeIntegrationScheme &in)
BDF multi-step scheme of order 2 (implicit)
L-stable, two stage, second order IMEX DIRK(2,2,2)
TimeIntegrationMethod m_method
integration method
void DoOdeExplicitRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DoOdeImplicitRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
Array< OneD, Array< OneD, NekDouble > > OutArrayType
Diagonally Implicit Runge Kutta scheme of order 2.
IMEX 3rd order scheme using Backward Different Formula & Extrapolation.
const Array< OneD, const NekDouble > ConstSingleArray
L-stable, three stage, third order IMEX DIRK(3,4,3)
Diagonally Implicit Runge Kutta scheme of order 3.
std::vector< TimeIntegrationSolutionSharedPtr > TimeIntegrationSolutionVector
const TimeIntegrationSchemeKey & GetIntegrationSchemeKey() const
const Array< OneD, const unsigned int > & GetTimeLevelOffset()
L-stable, two stage, third order IMEX DIRK(2,3,3)
IMEX 4th order scheme using Backward Different Formula & Extrapolation.
const TimeIntegrationSchemeSharedPtr & GetIntegrationScheme() const
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
Array< OneD, Array< TwoD, NekDouble > > m_A
Modified Crank-Nicolson/Adams-Bashforth Order 2 (MCNAB)
IMEX 1st order scheme using Euler Backwards/Euler Forwards.
midpoint method (old name)
const Array< OneD, const unsigned int > & GetTimeLevelOffset()
NekDouble m_T
Used to store the Explicit stage derivative of IMEX schemes.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > TripleArray
std::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Adams-Bashforth Forward multi-step scheme of order 4.
TripleArray m_F_IMEX
Array corresponding to the stage Derivatives.
const Array< OneD, const Array< OneD, NekDouble > > ConstDoubleArray