Nektar++
CFLtester.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyAdvection.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 
38 #include <LocalRegions/TriExp.h>
39 #include <LocalRegions/QuadExp.h>
42 
43 namespace Nektar
44 {
45  string CFLtester::className = GetEquationSystemFactory().RegisterCreatorFunction("CFLtester", CFLtester::create, "Testing CFL restriction");
46 
50  : UnsteadySystem(pSession, pGraph),
51  AdvectionSystem(pSession, pGraph)
52  {
53  }
54 
56  {
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(
86  "AdvectionType", advName, "WeakDG");
88  GetAdvectionFactory().CreateInstance(advName, advName);
89  m_advObject->SetFluxVector(
91  m_session->LoadSolverInfo(
92  "UpwindType", riemName, "Upwind");
95  riemName, m_session);
96  m_riemannSolver->SetScalar(
97  "Vn", &CFLtester::GetNormalVelocity, this);
98 
99  m_advObject->SetRiemannSolver(m_riemannSolver);
100  m_advObject->InitObject(m_session, m_fields);
101  break;
102  }
103  // Continuous field
106  {
107  string advName;
108  m_session->LoadSolverInfo(
109  "AdvectionType", advName, "NonConservative");
111  GetAdvectionFactory().CreateInstance(advName, advName);
112  m_advObject->SetFluxVector(
113  &CFLtester::GetFluxVector, this);
114  break;
115  }
116  default:
117  {
118  ASSERTL0(false, "Unsupported projection type.");
119  break;
120  }
121  }
122 
124  {
127  }
128  else
129  {
130  ASSERTL0(false, "Implicit unsteady Advection not set up.");
131  }
132  }
133 
135  {
136  }
137 
139  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
140  Array<OneD, Array<OneD, NekDouble> >&outarray,
141  const NekDouble time)
142  {
143  int nvariables = inarray.size();
144  int npoints = GetNpoints();
145 
146  m_advObject->Advect(nvariables, m_fields, m_velocity, inarray,
147  outarray, time);
148  for(int j = 0; j < nvariables; ++j)
149  {
150  Vmath::Neg(npoints,outarray[j],1);
151  }
152  }
153 
155  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
156  Array<OneD, Array<OneD, NekDouble> >&outarray,
157  const NekDouble time)
158  {
159  int nvariables = inarray.size();
160  SetBoundaryConditions(time);
161 
162  switch(m_projectionType)
163  {
165  {
166  // Just copy over array
167  int npoints = GetNpoints();
168 
169  for(int j = 0; j < nvariables; ++j)
170  {
171  Vmath::Vcopy(npoints,inarray[j],1,outarray[j],1);
172  }
173  }
174  break;
177  {
179 
180  for(int j = 0; j < nvariables; ++j)
181  {
182  m_fields[j]->FwdTrans(inarray[j],coeffs);
183  m_fields[j]->BwdTrans_IterPerExp(coeffs,outarray[j]);
184  }
185  break;
186  }
187  default:
188  ASSERTL0(false,"Unknown projection scheme");
189  break;
190  }
191  }
192 
193  /**
194  * @brief Get the normal velocity
195  */
197  {
198  // Number of trace (interface) points
199  int nTracePts = GetTraceNpoints();
200 
201  // Auxiliary variable to compute the normal velocity
202  Array<OneD, NekDouble> tmp(nTracePts);
203 
204  // Reset the normal velocity
205  Vmath::Zero(nTracePts, m_traceVn, 1);
206 
207  for (int i = 0; i < m_velocity.size(); ++i)
208  {
209  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
210 
211  Vmath::Vvtvp(nTracePts,
212  m_traceNormals[i], 1,
213  tmp, 1,
214  m_traceVn, 1,
215  m_traceVn, 1);
216  }
217 
218  return m_traceVn;
219  }
220 
222  const Array<OneD, Array<OneD, NekDouble> > &physfield,
224  {
225  ASSERTL1(flux[0].size() == m_velocity.size(),
226  "Dimension of flux array and velocity array do not match");
227 
228  int nq = physfield[0].size();
229 
230  for (int i = 0; i < flux.size(); ++i)
231  {
232  for (int j = 0; j < flux[0].size(); ++j)
233  {
234  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
235  flux[i][j], 1);
236  }
237  }
238  }
239 
241  {
243  }
244 
245 
246 
248  const Array<OneD,int> ExpOrder,
249  const Array<OneD,NekDouble> CFL,
250  NekDouble timeCFL)
251  {
252 
253  int n_element = m_fields[0]->GetExpSize();
254 
255  //const NekDouble minLengthStdTri = 1.414213;
256  //const NekDouble minLengthStdQuad = 2.0;
257  //const NekDouble cLambda = 0.2; // Spencer book pag. 317
258 
259  Array<OneD, NekDouble> tstep (n_element, 0.0);
260  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
261  stdVelocity = GetStdVelocity(m_velocity);
262 
263  for (int el = 0; el < n_element; ++el)
264  {
265  int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
266  Array<OneD, NekDouble> one2D(npoints, 1.0);
267  //NekDouble Area = m_fields[0]->GetExp(el)->Integral(one2D);
268  if(std::dynamic_pointer_cast<LocalRegions::TriExp>(m_fields[0]->GetExp(el)))
269  {
270  //tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
271  //tstep[el] = timeCFL*minLengthStdTri/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
272  tstep[el] = CFL[el]/(stdVelocity[el]);
273  }
274  else if(std::dynamic_pointer_cast<LocalRegions::QuadExp>(m_fields[0]->GetExp(el)))
275  {
276  //tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
277  //tstep[el] = timeCFL*minLengthStdQuad/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
278  tstep[el] = CFL[el]/(stdVelocity[el]);
279  }
280  }
281 
282  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
283 
284  return TimeStep;
285  }
286 
287 
288 
290  int ExpOrder,
291  NekDouble CFL,
292  NekDouble TimeStability)
293  {
294  //================================================================
295  // This function has been created just to test specific problems, hence is not general
296  // and it has been implemented in a rude fashion, as the full CFLtester class.
297  // For real CFL calculations refer to the general implementation above. (A.Bolis)
298  //================================================================
299 
300  NekDouble TimeStep;
301  NekDouble SpatialStability;
302  int n_elements = m_fields[0]->GetExpSize();
303 
304  //solve ambiguity in windows
305  NekDouble n_elem = n_elements;
306  NekDouble DH = sqrt(n_elem);
307 
308  int H = (int)DH;
309  int P = ExpOrder-1;
310 
311  //================================================================
312  // Regular meshes
313 
314  //SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
315 
316  //================================================================
317  // Anisotropic meshes
318 
319  if (TimeStability == 1.0)
320  {
321  SpatialStability = EigenvaluesAnaMeshesAB2[H/2][P-1];
322  }
323  else if (TimeStability == 2.0)
324  {
325  SpatialStability = EigenvaluesAnaMeshesRK2[H/2][P-1];
326  }
327  else if (TimeStability == 2.784)
328  {
329  SpatialStability = EigenvaluesAnaMeshesRK4[H/2][P-1];
330  }
331  else
332  {
333  ASSERTL0(false,"Dominant eigenvalues database not present for this time-stepping scheme")
334  }
335 
336  //================================================================
337 
338  TimeStep = (TimeStability/SpatialStability)*CFL;
339 
340  //================================================================
341 
342  return TimeStep;
343  }
344 
345 
346 
348  const Array<OneD, Array<OneD,NekDouble> > inarray)
349  {
350  // Checking if the problem is 2D
351  ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
352 
353  int nTotQuadPoints = GetTotPoints();
354  int n_element = m_fields[0]->GetExpSize();
355  int nvel = inarray.size();
356 
357  NekDouble pntVelocity;
358 
359  // Getting the standard velocity vector on the 2D normal space
360  Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
361  Array<OneD, NekDouble> stdV(n_element, 0.0);
362 
363  for (int i = 0; i < nvel; ++i)
364  {
365  stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
366  }
367 
368  if (nvel == 2)
369  {
370  for (int el = 0; el < n_element; ++el)
371  {
373  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
374  int n_points = el2D->GetTotPoints();
375 
377  el2D->GetGeom2D()->GetMetricInfo()->GetJac();
379  el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
380 
381  if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype()
383  {
384  for (int i = 0; i < n_points; i++)
385  {
386  stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
387  + gmat[2][i]*inarray[1][i];
388 
389  stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
390  + gmat[3][i]*inarray[1][i];
391  }
392  }
393  else
394  {
395  for (int i = 0; i < n_points; i++)
396  {
397  stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
398  + gmat[2][0]*inarray[1][i];
399 
400  stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
401  + gmat[3][0]*inarray[1][i];
402  }
403  }
404 
405 
406  for (int i = 0; i < n_points; i++)
407  {
408  pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
409  + stdVelocity[1][i]*stdVelocity[1][i]);
410 
411  if (pntVelocity>stdV[el])
412  {
413  stdV[el] = pntVelocity;
414  }
415  }
416  }
417  }
418  else
419  {
420  for (int el = 0; el < n_element; ++el)
421  {
423  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
424 
425  int n_points = el3D->GetTotPoints();
426 
428  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
430  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
431 
432  if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype()
434  {
435  for (int i = 0; i < n_points; i++)
436  {
437  stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
438  + gmat[3][i]*inarray[1][i]
439  + gmat[6][i]*inarray[2][i];
440 
441  stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
442  + gmat[4][i]*inarray[1][i]
443  + gmat[7][i]*inarray[2][i];
444 
445  stdVelocity[2][i] = gmat[2][i]*inarray[0][i]
446  + gmat[5][i]*inarray[1][i]
447  + gmat[8][i]*inarray[2][i];
448  }
449  }
450  else
451  {
453  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
455  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
456 
457  for (int i = 0; i < n_points; i++)
458  {
459  stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
460  + gmat[3][0]*inarray[1][i]
461  + gmat[6][0]*inarray[2][i];
462 
463  stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
464  + gmat[4][0]*inarray[1][i]
465  + gmat[7][0]*inarray[2][i];
466 
467  stdVelocity[2][i] = gmat[2][0]*inarray[0][i]
468  + gmat[5][0]*inarray[1][i]
469  + gmat[8][0]*inarray[2][i];
470  }
471  }
472 
473  for (int i = 0; i < n_points; i++)
474  {
475  pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
476  + stdVelocity[1][i]*stdVelocity[1][i]
477  + stdVelocity[2][i]*stdVelocity[2][i]);
478 
479  if (pntVelocity > stdV[el])
480  {
481  stdV[el] = pntVelocity;
482  }
483  }
484  }
485  }
486 
487  return stdV;
488  }
489 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102
virtual void v_InitObject()
Init object for UnsteadySystem class.
Definition: CFLtester.cpp:55
Array< OneD, NekDouble > GetStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Definition: CFLtester.cpp:347
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Definition: CFLtester.h:101
CFLtester(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: CFLtester.cpp:47
virtual ~CFLtester()
Definition: CFLtester.cpp:134
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: CFLtester.h:86
virtual void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
Definition: CFLtester.cpp:240
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:154
Array< OneD, NekDouble > m_traceVn
Definition: CFLtester.h:103
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Definition: CFLtester.cpp:196
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: CFLtester.cpp:221
static std::string className
Definition: CFLtester.h:96
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:138
virtual NekDouble v_GetTimeStep(const Array< OneD, int > ExpOrder, const Array< OneD, NekDouble > CFL, NekDouble timeCFL)
Definition: CFLtester.cpp:247
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
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()
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)
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:134
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:49
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:46
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
static NekDouble EigenvaluesAnaMeshesRK2[6][14]
Definition: CFLtester.h:65
static NekDouble EigenvaluesAnaMeshesRK4[6][14]
Definition: CFLtester.h:73
static NekDouble EigenvaluesAnaMeshesAB2[6][14]
Definition: CFLtester.h:57
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:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:992
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:513
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
P
Definition: main.py:133
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267