Nektar++
TimeIntegrationSchemeGEM.h
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////polationMethod
2 //
3 // File: TimeIntegrationSchemeGEM.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 GEM base class
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_GEM
36 #define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_GEM
37 
38 #define LUE LIB_UTILITIES_EXPORT
39 
40 #include <string>
41 
44 
45 namespace Nektar
46 {
47 namespace LibUtilities
48 {
49 
50 ///////////////////////////////////////////////////////////////////////////////
51 /// Class for spectral deferred correction integration.
53 {
54 public:
55  TimeIntegrationSchemeGEM(std::string variant, unsigned int order,
56  std::vector<NekDouble> freeParams)
57  : TimeIntegrationScheme(variant, order, freeParams),
58  m_name("ExtrapolationMethod")
59  {
60  ASSERTL0(variant == "" || variant == "ExplicitEuler" ||
61  variant == "ImplicitEuler" || variant == "IMEXEuler" ||
62  variant == "ExplicitMidpoint" ||
63  variant == "ImplicitMidpoint",
64  "Extrapolation Time integration "
65  "scheme bad variant (ExplicitEuler, ImplicitEuler, "
66  "ExplicitMidpoint)")
67 
68  if (variant == "IMEXEuler")
69  {
70  std::cerr << "WARNING: IMEX Euler extrapolation method has been "
71  "implemented but its use is not recommended as the "
72  "approach is affected by order-reduction problems."
73  << std::endl;
74  }
75 
76  if (variant == "" || variant == "ExplicitEuler" ||
77  variant == "ImplicitEuler" || variant == "IMEXEuler")
78  {
79  ASSERTL0(order >= 1, "Extrapolation Time integration "
80  "scheme bad order numbers (>=1): " +
81  std::to_string(order));
82  }
83  else if (variant == "ExplicitMidpoint" || variant == "ImplicitMidpoint")
84  {
85  ASSERTL0(order >= 2, "Extrapolation Time integration "
86  "scheme bad order numbers (>=2): " +
87  std::to_string(order));
88 
89  ASSERTL0(order % 2 == 0,
90  "Extrapolation Time integration "
91  "scheme bad order numbers (even number): " +
92  std::to_string(order));
93  }
94  m_variant = variant;
95  m_order = order;
96  m_freeParams = freeParams;
97  if (variant == "ExplicitEuler" || variant == "ExplicitMidpoint")
98  {
100  }
101  else if (variant == "ImplicitEuler" || variant == "ImplicitMidpoint")
102  {
104  }
105  else if (variant == "IMEXEuler")
106  {
108  }
109  }
110 
111  /// Destructor
113  {
114  }
115 
117  std::string variant, unsigned int order,
118  std::vector<NekDouble> freeParams)
119  {
122  variant, order, freeParams);
123 
124  return p;
125  }
126 
127  static std::string className;
128 
129 protected:
130  LUE virtual std::string v_GetName() const override;
131  LUE virtual std::string v_GetVariant() const override;
132  LUE virtual unsigned int v_GetOrder() const override;
133  LUE virtual std::vector<NekDouble> v_GetFreeParams() const override;
135  const override;
136  LUE virtual NekDouble v_GetTimeStability() const override;
137  LUE virtual unsigned int v_GetNumIntegrationPhases() const override;
138 
139  /**
140  * \brief Gets the solution vector of the ODE
141  */
142  virtual const TripleArray &v_GetSolutionVector() const override;
143 
144  /**
145  * \brief Sets the solution vector of the ODE
146  */
147  virtual void v_SetSolutionVector(const int Offset,
148  const DoubleArray &y) override;
149 
150  // The worker methods from the base class that are virtual
151  LUE virtual void v_InitializeScheme(
152  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
153  const TimeIntegrationSchemeOperators &op) override;
154 
156  const int timestep, const NekDouble delta_t,
157  const TimeIntegrationSchemeOperators &op) override;
158 
159  LUE virtual void v_print(std::ostream &os) const override;
160  LUE virtual void v_printFull(std::ostream &os) const override;
161 
162  // Variables common to all schemes.
163  std::string m_name;
164  std::string m_variant;
165  std::string m_nQuadType;
166  unsigned int m_order{0};
167  bool m_initialized = false;
168  std::vector<NekDouble> m_freeParams;
170 
172 
173  // Storage of previous states and associated timesteps.
174  TripleArray m_Y; /// Array containing the stage values
175  TripleArray m_T; /// Array containing the solution values
176  TripleArray m_T0; /// Array containing the solution values
177  DoubleArray m_F; /// Array corresponding to the stage Derivatives
178  DoubleArray m_F0; /// Array corresponding to the stage Derivatives
179 
180  int m_nvars{0}; // Number of variables in the integration scheme.
181  int m_npoints{0}; // Number of points in the integration scheme.
182 
183 }; // end class TimeIntegrationSchemeGEM
184 
185 LUE std::ostream &operator<<(std::ostream &os,
186  const TimeIntegrationSchemeGEM &rhs);
187 LUE std::ostream &operator<<(std::ostream &os,
189 
190 } // end namespace LibUtilities
191 } // end namespace Nektar
192 
193 #endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define LUE
Class for spectral deferred correction integration.
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
int m_nvars
Array corresponding to the stage Derivatives.
virtual LUE ConstDoubleArray & v_TimeIntegrate(const int timestep, const NekDouble delta_t, const TimeIntegrationSchemeOperators &op) override
Worker method that performs the time integration.
virtual LUE NekDouble v_GetTimeStability() const override
virtual LUE void v_printFull(std::ostream &os) const override
virtual LUE std::string v_GetName() const override
virtual const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
static TimeIntegrationSchemeSharedPtr create(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
virtual LUE unsigned int v_GetNumIntegrationPhases() const override
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
DoubleArray m_F
Array containing the solution values.
TripleArray m_T
Array containing the stage values.
virtual void v_SetSolutionVector(const int Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
DoubleArray m_F0
Array corresponding to the stage Derivatives.
TripleArray m_T0
Array containing the solution values.
virtual LUE unsigned int v_GetOrder() const override
TimeIntegrationSchemeGEM(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
virtual LUE std::vector< NekDouble > v_GetFreeParams() const override
virtual LUE std::string v_GetVariant() const override
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.
Base class for time integration schemes.
Binds a set of functions for use by time integration schemes.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
@ eImplicit
Fully implicit scheme.
@ eExplicit
Formally explicit scheme.
@ eExtrapolationMethod
Extrapolation scheme.
@ eIMEX
Implicit Explicit General Linear Method.
std::shared_ptr< TimeIntegrationSchemeGEM > TimeIntegrationSchemeGEMSharedPtr
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