Nektar++
TimeIntegrationSchemeFIT.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: TimeIntegrationSchemeFIT.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 FIT base class
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 // Note: The file is named TimeIntegrationSchemeFIT to parallel the
36 // TimeIntegrationSchemeGLM file but the class is named
37 // FractionalInTimeIntegrationScheme so keep with the factory naming
38 // convention.
39 
40 #ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_FIT
41 #define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_FIT
42 
43 #define LUE LIB_UTILITIES_EXPORT
44 
45 #include <string>
46 
48 
49 namespace Nektar
50 {
51 namespace LibUtilities
52 {
53 
54 ///////////////////////////////////////////////////////////////////////////////
55 /// Class for fractional-in-time integration.
57 {
58 public:
59  /// Constructor
60  FractionalInTimeIntegrationScheme(std::string variant, unsigned int order,
61  std::vector<NekDouble> freeParams);
62 
63  /// Destructor
65  {
66  }
67 
68  /// Creator
70  std::string variant, unsigned int order,
71  std::vector<NekDouble> freeParams)
72  {
75  variant, order, freeParams);
76 
77  return p;
78  }
79 
80  static std::string className;
81 
82  // Friend classes
83  LUE friend std::ostream &operator<<(
84  std::ostream &os, const FractionalInTimeIntegrationScheme &rhs);
85  LUE friend std::ostream &operator<<(
86  std::ostream &os,
88 
89 protected:
90  // Access methods from the base class that are virtual
91  LUE virtual std::string v_GetName() const override
92  {
93  return m_name;
94  }
95 
96  LUE virtual std::string v_GetVariant() const override
97  {
98  return m_variant;
99  }
100 
101  LUE virtual unsigned int v_GetOrder() const override
102  {
103  return m_order;
104  }
105 
106  LUE virtual std::vector<NekDouble> v_GetFreeParams() const override
107  {
108  return m_freeParams;
109  }
110 
112  const override
113  {
114  return m_schemeType;
115  }
116 
117  LUE virtual NekDouble v_GetTimeStability() const override
118  {
119  return 1.0;
120  }
121 
122  LUE virtual unsigned int v_GetNumIntegrationPhases() const override
123  {
124  return 1;
125  }
126 
127  /**
128  * \brief Gets the solution vector of the ODE
129  */
130  virtual const TripleArray &v_GetSolutionVector() const override
131  {
132  return m_u;
133  }
134 
135  /**
136  * \brief Sets the solution vector of the ODE
137  */
138  virtual void v_SetSolutionVector(const int Offset,
139  const DoubleArray &y) override
140  {
141  m_u[Offset] = y;
142  }
143 
144  // The worker methods from the base class that are virtual
145  LUE virtual void v_InitializeScheme(
146  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
147  const TimeIntegrationSchemeOperators &op) override;
148 
150  const int timestep, const NekDouble delta_t,
151  const TimeIntegrationSchemeOperators &op) override;
152 
153  LUE virtual void v_print(std::ostream &os) const override;
154  LUE virtual void v_printFull(std::ostream &os) const override;
155 
156  struct Instance
157  {
158  int base;
159 
160  int index; // Index of this instance
161  bool active; // Used to determine if active
162  int activecounter; // counter used to flip active bit
164 
165  // Major storage for auxilliary ODE solutions.
166  // Storage for values of y currently used to update u
168  std::pair<int, int> stage_ind; // Time-step counters indicating the
169  // interval ymain is associated with
170 
171  // Staging allocation
174  int stage_cbase; // This base is halved after the first cycle
176  int stage_fbase; // This base is halved after the first cycle
177 
178  // Ceiling stash allocation
179  int cstash_counter; // Counter used to determine
180  // when to stash
181  int cstash_base; // base for counter
183  std::pair<int, int> cstash_ind; // ind(1) is never used:
184  // it always matches main.ind(1)
185 
186  // Ceiling sandbox allocation
187  bool csandbox_active; // Flag to determine when
188  // stash 2 is utilized
191  std::pair<int, int> csandbox_ind;
192 
193  // Floor stash
196  std::pair<int, int> fstash_ind;
197 
198  // Floor sandbox
203  std::pair<int, int> fsandbox_ind;
204 
205  // Talbot quadrature rule
208 
210 
214  };
215 
216  inline unsigned int modIncrement(const unsigned int counter,
217  const unsigned int base) const;
218 
219  inline unsigned int computeL(const unsigned int base,
220  const unsigned int m) const;
221 
222  inline unsigned int computeQML(const unsigned int base,
223  const unsigned int m);
224 
225  inline unsigned int computeTaus(const unsigned int base,
226  const unsigned int m);
227 
228  void talbotQuadrature(const unsigned int nQuadPts, const NekDouble mu,
229  const NekDouble nu, const NekDouble sigma,
230  ComplexSingleArray &lamb,
231  ComplexSingleArray &w) const;
232 
233  void integralClassInitialize(const unsigned int index,
234  Instance &instance) const;
235 
236  void updateStage(const unsigned int timeStep, Instance &instance);
237 
238  void finalIncrement(const unsigned int timeStep,
240 
241  void integralContribution(const unsigned int timeStep,
242  const unsigned int tauml,
243  const Instance &instance);
244 
245  void timeAdvance(const unsigned int timeStep,
247  Instance &instance, ComplexTripleArray &y);
248 
249  void advanceSandbox(const unsigned int timeStep,
251  Instance &instance);
252 
253  // Variables common to all schemes.
254  std::string m_name;
255  std::string m_variant;
256  unsigned int m_order{0};
257  std::vector<NekDouble> m_freeParams;
258 
260 
261  // Varaibles and methods specific to FIT integration schemes.
263  NekDouble m_T{0}; // Finial time
264  unsigned int m_maxTimeSteps; // Number of time steps.
265  NekDouble m_alpha{0.3}; // Value for exp integration.
266  unsigned int m_base{4}; // "Base" of the algorithm.
267  unsigned int m_nQuadPts{20}; // Number of Talbot quadrature rule points
271 
272  int m_nvars{0}; // Number of variables in the integration scheme.
273  int m_npoints{0}; // Number of points in the integration scheme.
274 
275  unsigned int m_Lmax{0}; // Maxium number of integral groups.
277 
278  // Demarcation integers
280  // Demarcation interval markers
282 
283  // Storage of the initial values.
285  // Storage of the next solution from the final increment.
287  // Storage for the integral contribution.
289  // Storage for the exponential factor in the integral contribution.
291 
292  // Storage of previous states and associated timesteps.
294 
295  // Storage for the stage derivative as the data will be re-used to
296  // update the solution.
298 
299  // J
301 
302  // Ahat array one for each order.
304 
305  // Mulitply the last Ahat array, transposed by J
307 
308 }; // end class FractionalInTimeIntegrator
309 
310 LUE std::ostream &operator<<(std::ostream &os,
312 LUE std::ostream &operator<<(
313  std::ostream &os, const FractionalInTimeIntegrationSchemeSharedPtr &rhs);
314 
315 } // end namespace LibUtilities
316 } // end namespace Nektar
317 
318 #endif
#define LUE
unsigned int modIncrement(const unsigned int counter, const unsigned int base) const
Method that increments the counter then performs mod calculation.
void updateStage(const unsigned int timeStep, Instance &instance)
Method to rearrange of staging/stashing for current time.
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 std::vector< NekDouble > v_GetFreeParams() const override
void integralClassInitialize(const unsigned int index, Instance &instance) const
Method to initialize the integral class.
virtual LUE ConstDoubleArray & v_TimeIntegrate(const int timestep, const NekDouble delta_t, const TimeIntegrationSchemeOperators &op) override
Worker method that performs the time integration.
unsigned int computeQML(const unsigned int base, const unsigned int m)
Method to compute the demarcation integers q_{m, ell}.
virtual void v_SetSolutionVector(const int Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
unsigned int computeL(const unsigned int base, const unsigned int m) const
Method to compute the smallest integer L such that base < 2.
void timeAdvance(const unsigned int timeStep, const TimeIntegrationSchemeOperators &op, Instance &instance, ComplexTripleArray &y)
Method to get the solution to y' = z*y + f(u), using an exponential integrator with implicit order (m...
static TimeIntegrationSchemeSharedPtr create(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
Creator.
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
virtual const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
unsigned int computeTaus(const unsigned int base, const unsigned int m)
Method to compute the demarcation interval marker tau_{m, ell}.
void talbotQuadrature(const unsigned int nQuadPts, const NekDouble mu, const NekDouble nu, const NekDouble sigma, ComplexSingleArray &lamb, ComplexSingleArray &w) const
Method to compute the quadrature rule over Tablot contour.
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
void advanceSandbox(const unsigned int timeStep, const TimeIntegrationSchemeOperators &op, Instance &instance)
Method to update sandboxes to the current time.
virtual LUE void v_printFull(std::ostream &os) const override
FractionalInTimeIntegrationScheme(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
Constructor.
LUE friend std::ostream & operator<<(std::ostream &os, const FractionalInTimeIntegrationScheme &rhs)
void finalIncrement(const unsigned int timeStep, const TimeIntegrationSchemeOperators &op)
Method to approximate the integral.
void integralContribution(const unsigned int timeStep, const unsigned int tauml, const Instance &instance)
Method to get the integral contribution over [taus(i+1) taus(i)]. Stored in m_uInt.
virtual LUE unsigned int v_GetNumIntegrationPhases() const override
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.
@ eFractionalInTime
Fractional in Time scheme.
std::shared_ptr< FractionalInTimeIntegrationScheme > FractionalInTimeIntegrationSchemeSharedPtr
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