Nektar++
TimeIntegrationSchemeSDC.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: TimeIntegrationSchemeSDC.cpp
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: implementation of time integration scheme SDC class
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 
38 namespace Nektar
39 {
40 namespace LibUtilities
41 {
42 
44 {
45  return m_name;
46 }
47 
49 {
50  return m_variant;
51 }
52 
54 {
55  return m_order;
56 }
57 
58 std::vector<NekDouble> TimeIntegrationSchemeSDC::v_GetFreeParams() const
59 {
60  return m_freeParams;
61 }
62 
64  const
65 {
66  return m_schemeType;
67 }
68 
70 {
71  return 1.0;
72 }
73 
75 {
76  return 1;
77 }
78 
79 /**
80  * \brief Gets the solution vector of the ODE
81  */
83 {
84  return m_Y;
85 }
86 
87 /**
88  * \brief Sets the solution vector of the ODE
89  */
91  const DoubleArray &y)
92 {
93  m_Y[Offset] = y;
94 }
95 
96 /**
97  * @brief Worker method to initialize the integration scheme.
98  */
100  const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
102 {
103  boost::ignore_unused(op);
104  boost::ignore_unused(deltaT);
105 
106  if (m_initialized)
107  {
108  for (unsigned int i = 0; i < m_nvars; ++i)
109  {
110  // Store the initial values as the first previous state.
111  Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
112  }
113  }
114  else
115  {
116  m_time = time;
117  m_nvars = y_0.size();
118  m_npoints = y_0[0].size();
119 
120  // Compute integration matrix
122 
123  // Use first quadrature point in the integral
124  if (m_variant == "Equidistant" || m_variant == "GaussLobattoLegendre" ||
125  m_variant == "GaussLobattoChebyshev")
126  {
127  Polylib::Qg(&m_QMat[0], &m_points[0], m_nQuadPts, 0);
128  }
129  else if (m_variant == "GaussRadauLegendre" ||
130  m_variant == "GaussRadauChebyshev")
131  {
133  }
134 
135  // Storage of previous states and associated timesteps.
138 
139  for (unsigned int m = 0; m < m_nQuadPts; ++m)
140  {
141  m_Y[m] = DoubleArray(m_nvars);
142  m_F[m] = DoubleArray(m_nvars);
143 
144  for (unsigned int i = 0; i < m_nvars; ++i)
145  {
146  m_Y[m][i] = SingleArray(m_npoints, 0.0);
147  m_F[m][i] = SingleArray(m_npoints, 0.0);
148 
149  // Store the initial values as the first previous state.
150  if (m == 0)
151  {
152  Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[m][i], 1);
153  }
154  }
155  }
156 
158  for (unsigned int m = 0; m < m_nQuadPts - 1; ++m)
159  {
160  m_Fint[m] = DoubleArray(m_nvars);
161 
162  for (unsigned int i = 0; i < m_nvars; ++i)
163  {
164  m_Fint[m][i] = SingleArray(m_npoints, 0.0);
165  }
166  }
167 
170  for (unsigned int i = 0; i < m_nvars; ++i)
171  {
172  m_Fn[i] = SingleArray(m_npoints, 0.0);
173  m_tmp[i] = SingleArray(m_npoints, 0.0);
174  }
175 
176  m_initialized = true;
177  }
178 }
179 
180 /**
181  * @brief Worker method that performs the time integration.
182  */
184  const int timestep, const NekDouble delta_t,
186 {
187 
188  boost::ignore_unused(timestep);
189 
190  // 1. Compute initial guess
191  ComputeInitialGuess(delta_t, op);
192 
193  // 2. Return initial guess if m_order = 1
194  if (m_order == 1)
195  {
196  // 2.1. Copy final solution
197  for (unsigned int i = 0; i < m_nvars; ++i)
198  {
199  Vmath::Vcopy(m_npoints, m_Y[m_nQuadPts - 1][i], 1, m_Y[0][i], 1);
200  }
201 
202  // 2.2. Update time step
203  m_time += delta_t;
204 
205  // 2.3. Return solution
206  return m_Y[m_nQuadPts - 1];
207  }
208 
209  // 3. Apply SDC correction loop
210  v_SDCIterationLoop(delta_t, op);
211 
212  // 4. Get solution
213  // 4.1. Copy final solution
214  for (unsigned int i = 0; i < m_nvars; ++i)
215  {
216  Vmath::Vcopy(m_npoints, m_Y[m_nQuadPts - 1][i], 1, m_Y[0][i], 1);
217  }
218 
219  // 4.2. Update time step
220  m_time += delta_t;
221 
222  // 4.3. Return solution
223  return m_Y[m_nQuadPts - 1];
224 }
225 
226 /**
227  * @brief Worker method that compute the integrated flux.
228  */
230  const NekDouble &delta_t, const int option)
231 {
232  // Zeroing integrated flux
233  for (unsigned int n = 0; n < m_nQuadPts - 1; ++n)
234  {
235  for (unsigned int i = 0; i < m_nvars; ++i)
236  {
237  Vmath::Zero(m_npoints, m_Fint[n][i], 1);
238  }
239  }
240 
241  // Update integrated flux
242  for (unsigned int n = 0; n < m_nQuadPts - 1; ++n)
243  {
244  for (unsigned int p = 0; p < m_nQuadPts; ++p)
245  {
246  for (unsigned int i = 0; i < m_nvars; ++i)
247  {
248  if (option == 0)
249  {
251  delta_t * (m_QMat[(n + 1) * m_nQuadPts + p] -
252  m_QMat[n * m_nQuadPts + p]),
253  m_F[p][i], 1, m_Fint[n][i], 1, m_Fint[n][i],
254  1);
255  }
256  else if (option == 1)
257  {
258  Vmath::Svtvp(
259  m_npoints, delta_t * m_QMat[(n + 1) * m_nQuadPts + p],
260  m_F[p][i], 1, m_Fint[n][i], 1, m_Fint[n][i], 1);
261  }
262  }
263  }
264  }
265 }
266 
267 /**
268  * @brief Worker method to print details on the integration scheme
269  */
270 void TimeIntegrationSchemeSDC::v_print(std::ostream &os) const
271 {
272  os << "Time Integration Scheme: " << GetFullName() << std::endl;
273 }
274 
275 void TimeIntegrationSchemeSDC::v_printFull(std::ostream &os) const
276 {
277  os << "Time Integration Scheme: " << GetFullName() << std::endl;
278 }
279 
280 // Friend Operators
281 std::ostream &operator<<(std::ostream &os, const TimeIntegrationSchemeSDC &rhs)
282 {
283  rhs.print(os);
284 
285  return os;
286 }
287 
288 std::ostream &operator<<(std::ostream &os,
290 {
291  os << *rhs.get();
292 
293  return os;
294 }
295 
296 } // end namespace LibUtilities
297 } // namespace Nektar
Binds a set of functions for use by time integration schemes.
Class for spectral deferred correction integration.
virtual LUE unsigned int v_GetOrder() const override
LUE void UpdateIntegratedResidual(const NekDouble &delta_t, const int option=0)
Worker method that compute the integrated flux.
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 TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
virtual LUE std::string v_GetName() const override
LUE void ComputeInitialGuess(const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
virtual LUE void v_SetSolutionVector(const int Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
virtual LUE void v_printFull(std::ostream &os) const override
TripleArray m_Fint
Array corresponding to the stage Derivatives.
virtual LUE std::vector< NekDouble > v_GetFreeParams() const override
TripleArray m_F
Array containing the stage values.
virtual LUE unsigned int v_GetNumIntegrationPhases() 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.
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t, const TimeIntegrationSchemeOperators &op)
virtual LUE std::string v_GetVariant() const override
virtual LUE NekDouble v_GetTimeStability() const override
virtual LUE const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
AT< AT< NekDouble > > DoubleArray
std::shared_ptr< TimeIntegrationSchemeSDC > TimeIntegrationSchemeSDCSharedPtr
AT< AT< AT< NekDouble > > > TripleArray
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Qg(double *Q, const double *z, const int np, const int offset)
Compute the Integration Matrix.
Definition: Polylib.cpp:584
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255