Nektar++
CFLtester.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: CFLtester.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: CFL tester solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 
40 #include <LocalRegions/QuadExp.h>
41 #include <LocalRegions/TriExp.h>
42 
43 namespace Nektar
44 {
45 string CFLtester::className =
47  "CFLtester", CFLtester::create, "Testing CFL restriction");
48 
51  : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
52 {
53 }
54 
55 void CFLtester::v_InitObject(bool DeclareFields)
56 {
57  AdvectionSystem::v_InitObject(DeclareFields);
58 
60  std::vector<std::string> vel;
61  vel.push_back("Vx");
62  vel.push_back("Vy");
63  vel.push_back("Vz");
64  vel.resize(m_spacedim);
65 
66  GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
67 
68  // Type of advection class to be used
69  switch (m_projectionType)
70  {
71  // Discontinuous field
73  {
74  // Do not forwards transform initial condition
75  m_homoInitialFwd = false;
76 
77  // Define the normal velocity fields
78  if (m_fields[0]->GetTrace())
79  {
81  }
82 
83  string advName;
84  string riemName;
85  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
87  advName, advName);
88  m_advObject->SetFluxVector(&CFLtester::GetFluxVector, this);
89  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
92  riemName, m_session);
94  this);
95 
96  m_advObject->SetRiemannSolver(m_riemannSolver);
97  m_advObject->InitObject(m_session, m_fields);
98  break;
99  }
100  // Continuous field
103  {
104  string advName;
105  m_session->LoadSolverInfo("AdvectionType", advName,
106  "NonConservative");
108  advName, advName);
109  m_advObject->SetFluxVector(&CFLtester::GetFluxVector, this);
110  break;
111  }
112  default:
113  {
114  ASSERTL0(false, "Unsupported projection type.");
115  break;
116  }
117  }
118 
120  {
123  }
124  else
125  {
126  ASSERTL0(false, "Implicit unsteady Advection not set up.");
127  }
128 }
129 
131 {
132 }
133 
135  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
136  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
137 {
138  int nvariables = inarray.size();
139  int npoints = GetNpoints();
140 
141  m_advObject->Advect(nvariables, m_fields, m_velocity, inarray, outarray,
142  time);
143  for (int j = 0; j < nvariables; ++j)
144  {
145  Vmath::Neg(npoints, outarray[j], 1);
146  }
147 }
148 
150  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
151  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
152 {
153  int nvariables = inarray.size();
154  SetBoundaryConditions(time);
155 
156  switch (m_projectionType)
157  {
159  {
160  // Just copy over array
161  int npoints = GetNpoints();
162 
163  for (int j = 0; j < nvariables; ++j)
164  {
165  Vmath::Vcopy(npoints, inarray[j], 1, outarray[j], 1);
166  }
167  }
168  break;
171  {
173 
174  for (int j = 0; j < nvariables; ++j)
175  {
176  m_fields[j]->FwdTrans(inarray[j], coeffs);
177  m_fields[j]->BwdTrans_IterPerExp(coeffs, outarray[j]);
178  }
179  break;
180  }
181  default:
182  ASSERTL0(false, "Unknown projection scheme");
183  break;
184  }
185 }
186 
187 /**
188  * @brief Get the normal velocity
189  */
191 {
192  // Number of trace (interface) points
193  int nTracePts = GetTraceNpoints();
194 
195  // Auxiliary variable to compute the normal velocity
196  Array<OneD, NekDouble> tmp(nTracePts);
197 
198  // Reset the normal velocity
199  Vmath::Zero(nTracePts, m_traceVn, 1);
200 
201  for (int i = 0; i < m_velocity.size(); ++i)
202  {
203  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
204 
205  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
206  m_traceVn, 1);
207  }
208 
209  return m_traceVn;
210 }
211 
213  const Array<OneD, Array<OneD, NekDouble>> &physfield,
215 {
216  ASSERTL1(flux[0].size() == m_velocity.size(),
217  "Dimension of flux array and velocity array do not match");
218 
219  int nq = physfield[0].size();
220 
221  for (int i = 0; i < flux.size(); ++i)
222  {
223  for (int j = 0; j < flux[0].size(); ++j)
224  {
225  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
226  }
227  }
228 }
229 
231 {
233 }
234 
236  const Array<OneD, NekDouble> CFL,
237  NekDouble timeCFL)
238 {
239 
240  int n_element = m_fields[0]->GetExpSize();
241 
242  // const NekDouble minLengthStdTri = 1.414213;
243  // const NekDouble minLengthStdQuad = 2.0;
244  // const NekDouble cLambda = 0.2; // Spencer book pag. 317
245 
246  Array<OneD, NekDouble> tstep(n_element, 0.0);
247  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
248  stdVelocity = GetStdVelocity(m_velocity);
249 
250  for (int el = 0; el < n_element; ++el)
251  {
252  int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
253  Array<OneD, NekDouble> one2D(npoints, 1.0);
254  // NekDouble Area = m_fields[0]->GetExp(el)->Integral(one2D);
255  if (std::dynamic_pointer_cast<LocalRegions::TriExp>(
256  m_fields[0]->GetExp(el)))
257  {
258  // tstep[el] =
259  // timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
260  // tstep[el] =
261  // timeCFL*minLengthStdTri/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
262  tstep[el] = CFL[el] / (stdVelocity[el]);
263  }
264  else if (std::dynamic_pointer_cast<LocalRegions::QuadExp>(
265  m_fields[0]->GetExp(el)))
266  {
267  // tstep[el] =
268  // timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
269  // tstep[el] =
270  // timeCFL*minLengthStdQuad/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
271  tstep[el] = CFL[el] / (stdVelocity[el]);
272  }
273  }
274 
275  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
276 
277  return TimeStep;
278 }
279 
281  NekDouble TimeStability)
282 {
283  //================================================================
284  // This function has been created just to test specific problems, hence is
285  // not general and it has been implemented in a rude fashion, as the full
286  // CFLtester class. For real CFL calculations refer to the general
287  // implementation above. (A.Bolis)
288  //================================================================
289 
290  NekDouble TimeStep;
291  NekDouble SpatialStability;
292  int n_elements = m_fields[0]->GetExpSize();
293 
294  // solve ambiguity in windows
295  NekDouble n_elem = n_elements;
296  NekDouble DH = sqrt(n_elem);
297 
298  int H = (int)DH;
299  int P = ExpOrder - 1;
300 
301  //================================================================
302  // Regular meshes
303 
304  // SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
305 
306  //================================================================
307  // Anisotropic meshes
308 
309  if (TimeStability == 1.0)
310  {
311  SpatialStability = EigenvaluesAnaMeshesAB2[H / 2][P - 1];
312  }
313  else if (TimeStability == 2.0)
314  {
315  SpatialStability = EigenvaluesAnaMeshesRK2[H / 2][P - 1];
316  }
317  else if (TimeStability == 2.784)
318  {
319  SpatialStability = EigenvaluesAnaMeshesRK4[H / 2][P - 1];
320  }
321  else
322  {
323  ASSERTL0(false, "Dominant eigenvalues database not present for this "
324  "time-stepping scheme")
325  }
326 
327  //================================================================
328 
329  TimeStep = (TimeStability / SpatialStability) * CFL;
330 
331  //================================================================
332 
333  return TimeStep;
334 }
335 
337  const Array<OneD, Array<OneD, NekDouble>> inarray)
338 {
339  // Checking if the problem is 2D
340  ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
341 
342  int nTotQuadPoints = GetTotPoints();
343  int n_element = m_fields[0]->GetExpSize();
344  int nvel = inarray.size();
345 
346  NekDouble pntVelocity;
347 
348  // Getting the standard velocity vector on the 2D normal space
349  Array<OneD, Array<OneD, NekDouble>> stdVelocity(nvel);
350  Array<OneD, NekDouble> stdV(n_element, 0.0);
351 
352  for (int i = 0; i < nvel; ++i)
353  {
354  stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
355  }
356 
357  if (nvel == 2)
358  {
359  for (int el = 0; el < n_element; ++el)
360  {
362  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
363  int n_points = el2D->GetTotPoints();
364 
366  el2D->GetGeom2D()->GetMetricInfo()->GetJac();
368  el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
369 
370  if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype() ==
372  {
373  for (int i = 0; i < n_points; i++)
374  {
375  stdVelocity[0][i] =
376  gmat[0][i] * inarray[0][i] + gmat[2][i] * inarray[1][i];
377 
378  stdVelocity[1][i] =
379  gmat[1][i] * inarray[0][i] + gmat[3][i] * inarray[1][i];
380  }
381  }
382  else
383  {
384  for (int i = 0; i < n_points; i++)
385  {
386  stdVelocity[0][i] =
387  gmat[0][0] * inarray[0][i] + gmat[2][0] * inarray[1][i];
388 
389  stdVelocity[1][i] =
390  gmat[1][0] * inarray[0][i] + gmat[3][0] * inarray[1][i];
391  }
392  }
393 
394  for (int i = 0; i < n_points; i++)
395  {
396  pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
397  stdVelocity[1][i] * stdVelocity[1][i]);
398 
399  if (pntVelocity > stdV[el])
400  {
401  stdV[el] = pntVelocity;
402  }
403  }
404  }
405  }
406  else
407  {
408  for (int el = 0; el < n_element; ++el)
409  {
411  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
412 
413  int n_points = el3D->GetTotPoints();
414 
416  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
418  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
419 
420  if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype() ==
422  {
423  for (int i = 0; i < n_points; i++)
424  {
425  stdVelocity[0][i] = gmat[0][i] * inarray[0][i] +
426  gmat[3][i] * inarray[1][i] +
427  gmat[6][i] * inarray[2][i];
428 
429  stdVelocity[1][i] = gmat[1][i] * inarray[0][i] +
430  gmat[4][i] * inarray[1][i] +
431  gmat[7][i] * inarray[2][i];
432 
433  stdVelocity[2][i] = gmat[2][i] * inarray[0][i] +
434  gmat[5][i] * inarray[1][i] +
435  gmat[8][i] * inarray[2][i];
436  }
437  }
438  else
439  {
441  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
443  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
444 
445  for (int i = 0; i < n_points; i++)
446  {
447  stdVelocity[0][i] = gmat[0][0] * inarray[0][i] +
448  gmat[3][0] * inarray[1][i] +
449  gmat[6][0] * inarray[2][i];
450 
451  stdVelocity[1][i] = gmat[1][0] * inarray[0][i] +
452  gmat[4][0] * inarray[1][i] +
453  gmat[7][0] * inarray[2][i];
454 
455  stdVelocity[2][i] = gmat[2][0] * inarray[0][i] +
456  gmat[5][0] * inarray[1][i] +
457  gmat[8][0] * inarray[2][i];
458  }
459  }
460 
461  for (int i = 0; i < n_points; i++)
462  {
463  pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
464  stdVelocity[1][i] * stdVelocity[1][i] +
465  stdVelocity[2][i] * stdVelocity[2][i]);
466 
467  if (pntVelocity > stdV[el])
468  {
469  stdV[el] = pntVelocity;
470  }
471  }
472  }
473  }
474 
475  return stdV;
476 }
477 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:158
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Definition: CFLtester.h:157
CFLtester(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: CFLtester.cpp:49
virtual ~CFLtester()
Definition: CFLtester.cpp:130
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Definition: CFLtester.cpp:134
virtual NekDouble v_GetTimeStep(const Array< OneD, int > ExpOrder, const Array< OneD, NekDouble > CFL, NekDouble timeCFL) override
Definition: CFLtester.cpp:235
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: CFLtester.h:142
virtual void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
Definition: CFLtester.cpp:230
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Definition: CFLtester.cpp:149
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
Definition: CFLtester.cpp:55
Array< OneD, NekDouble > m_traceVn
Definition: CFLtester.h:159
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Definition: CFLtester.cpp:190
Array< OneD, NekDouble > GetStdVelocity(const Array< OneD, Array< OneD, NekDouble >> inarray)
Definition: CFLtester.cpp:336
static std::string className
Definition: CFLtester.h:152
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
Definition: CFLtester.cpp:212
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
int m_expdim
Expansion dimension.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT int GetTotPoints()
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:50
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
static NekDouble EigenvaluesAnaMeshesRK2[6][14]
Definition: CFLtester.h:97
static NekDouble EigenvaluesAnaMeshesRK4[6][14]
Definition: CFLtester.h:117
static NekDouble EigenvaluesAnaMeshesAB2[6][14]
Definition: CFLtester.h:77
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:1050
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294