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
49namespace Nektar
50{
51namespace LibUtilities
52{
53
54///////////////////////////////////////////////////////////////////////////////
55/// Class for fractional-in-time integration.
57{
58public:
59 /// Constructor
60 FractionalInTimeIntegrationScheme(std::string variant, size_t order,
61 std::vector<NekDouble> freeParams);
62
63 /// Destructor
65 {
66 }
67
68 /// Creator
70 std::string variant, size_t order, std::vector<NekDouble> freeParams)
71 {
74 variant, order, freeParams);
75
76 return p;
77 }
78
79 static std::string className;
80
81 // Friend classes
82 LUE friend std::ostream &operator<<(
83 std::ostream &os, const FractionalInTimeIntegrationScheme &rhs);
84 LUE friend std::ostream &operator<<(
85 std::ostream &os,
87
88protected:
89 // Access methods from the base class that are virtual
90 LUE virtual std::string v_GetName() const override
91 {
92 return m_name;
93 }
94
95 LUE virtual std::string v_GetVariant() const override
96 {
97 return m_variant;
98 }
99
100 LUE virtual size_t v_GetOrder() const override
101 {
102 return m_order;
103 }
104
105 LUE virtual std::vector<NekDouble> v_GetFreeParams() const override
106 {
107 return m_freeParams;
108 }
109
111 const override
112 {
113 return m_schemeType;
114 }
115
116 LUE virtual NekDouble v_GetTimeStability() const override
117 {
118 return 1.0;
119 }
120
121 LUE virtual size_t v_GetNumIntegrationPhases() const override
122 {
123 return 1;
124 }
125
126 /**
127 * \brief Gets the solution vector of the ODE
128 */
129 virtual const TripleArray &v_GetSolutionVector() const override
130 {
131 return m_u;
132 }
134 {
135 return m_u;
136 }
137
138 /**
139 * \brief Sets the solution vector of the ODE
140 */
141 virtual void v_SetSolutionVector(const size_t Offset,
142 const DoubleArray &y) override
143 {
144 m_u[Offset] = y;
145 }
146
147 // The worker methods from the base class that are virtual
148 LUE virtual void v_InitializeScheme(
149 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
150 const TimeIntegrationSchemeOperators &op) override;
151
153 const size_t timestep, const NekDouble delta_t) override;
154
155 LUE virtual void v_print(std::ostream &os) const override;
156 LUE virtual void v_printFull(std::ostream &os) const override;
157
158 struct Instance
159 {
160 size_t base;
161
162 size_t index; // Index of this instance
163 bool active; // Used to determine if active
164 size_t activecounter; // counter used to flip active bit
166
167 // Major storage for auxilliary ODE solutions.
168 // Storage for values of y currently used to update u
170 std::pair<size_t, size_t>
171 stage_ind; // Time-step counters indicating the
172 // interval ymain is associated with
173
174 // Staging allocation
177 size_t stage_cbase; // This base is halved after the first cycle
179 size_t stage_fbase; // This base is halved after the first cycle
180
181 // Ceiling stash allocation
182 size_t cstash_counter; // Counter used to determine
183 // when to stash
184 size_t cstash_base; // base for counter
186 std::pair<size_t, size_t> cstash_ind; // ind(1) is never used:
187 // it always matches main.ind(1)
188
189 // Ceiling sandbox allocation
190 bool csandbox_active; // Flag to determine when
191 // stash 2 is utilized
194 std::pair<size_t, size_t> csandbox_ind;
195
196 // Floor stash
199 std::pair<size_t, size_t> fstash_ind;
200
201 // Floor sandbox
206 std::pair<size_t, size_t> fsandbox_ind;
207
208 // Talbot quadrature rule
211
213
217 };
218
219 inline size_t modIncrement(const size_t counter, const size_t base) const;
220
221 inline size_t computeL(const size_t base, const size_t m) const;
222
223 inline size_t computeQML(const size_t base, const size_t m);
224
225 inline size_t computeTaus(const size_t base, const size_t m);
226
227 void talbotQuadrature(const size_t nQuadPts, const NekDouble mu,
228 const NekDouble nu, const NekDouble sigma,
229 ComplexSingleArray &lamb,
230 ComplexSingleArray &w) const;
231
232 void integralClassInitialize(const size_t index, Instance &instance) const;
233
234 void updateStage(const size_t timeStep, Instance &instance);
235
236 void finalIncrement(const size_t timeStep);
237
238 void integralContribution(const size_t timeStep, const size_t tauml,
239 const Instance &instance);
240
241 void timeAdvance(const size_t timeStep, Instance &instance,
243
244 void advanceSandbox(const size_t timeStep, Instance &instance);
245
246 // Variables common to all schemes.
247 std::string m_name;
248 std::string m_variant;
249 size_t m_order{0};
250 std::vector<NekDouble> m_freeParams;
251
254
255 // Varaibles and methods specific to FIT integration schemes.
257 NekDouble m_T{0}; // Finial time
258 size_t m_maxTimeSteps; // Number of time steps.
259 NekDouble m_alpha{0.3}; // Value for exp integration.
260 size_t m_base{4}; // "Base" of the algorithm.
261 size_t m_nQuadPts{20}; // Number of Talbot quadrature rule points
265
266 size_t m_nvars{0}; // Number of variables in the integration scheme.
267 size_t m_npoints{0}; // Number of points in the integration scheme.
268
269 size_t m_Lmax{0}; // Maxium number of integral groups.
271
272 // Demarcation integers
274 // Demarcation interval markers
276
277 // Storage of the initial values.
279 // Storage of the next solution from the final increment.
281 // Storage for the integral contribution.
283 // Storage for the exponential factor in the integral contribution.
285
286 // Storage of previous states and associated timesteps.
288
289 // Storage for the stage derivative as the data will be re-used to
290 // update the solution.
292
293 // J
295
296 // Ahat array one for each order.
298
299 // Multiply the last Ahat array, transposed by J
301
302}; // end class FractionalInTimeIntegrator
303
304LUE std::ostream &operator<<(std::ostream &os,
306LUE std::ostream &operator<<(
307 std::ostream &os, const FractionalInTimeIntegrationSchemeSharedPtr &rhs);
308
309} // end namespace LibUtilities
310} // end namespace Nektar
311
312#endif
#define LUE
size_t computeTaus(const size_t base, const size_t m)
Method to compute the demarcation interval marker tau_{m, ell}.
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 const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
void timeAdvance(const size_t timeStep, 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, size_t order, std::vector< NekDouble > freeParams)
Creator.
void updateStage(const size_t timeStep, Instance &instance)
Method to rearrange of staging/stashing for current time.
void advanceSandbox(const size_t timeStep, Instance &instance)
Method to update sandboxes to the current time.
void integralClassInitialize(const size_t index, Instance &instance) const
Method to initialize the integral class.
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
size_t modIncrement(const size_t counter, const size_t base) const
Method that increments the counter then performs mod calculation.
void finalIncrement(const size_t timeStep)
Method to approximate the integral.
virtual LUE std::vector< NekDouble > v_GetFreeParams() const override
void integralContribution(const size_t timeStep, const size_t tauml, const Instance &instance)
Method to get the integral contribution over [taus(i+1) taus(i)]. Stored in m_uInt.
void talbotQuadrature(const size_t nQuadPts, const NekDouble mu, const NekDouble nu, const NekDouble sigma, ComplexSingleArray &lamb, ComplexSingleArray &w) const
Method to compute the quadrature rule over Tablot contour.
size_t computeL(const size_t base, const size_t m) const
Method to compute the smallest integer L such that base < 2.
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
virtual LUE void v_printFull(std::ostream &os) const override
virtual LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
FractionalInTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
Constructor.
size_t computeQML(const size_t base, const size_t m)
Method to compute the demarcation integers q_{m, ell}.
virtual void v_SetSolutionVector(const size_t Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
LUE friend std::ostream & operator<<(std::ostream &os, const FractionalInTimeIntegrationScheme &rhs)
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
std::vector< double > w(NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble