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
41#include <LocalRegions/TriExp.h>
42
43namespace Nektar
44{
47 "CFLtester", CFLtester::create, "Testing CFL restriction");
48
51 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
52{
53}
54
55void 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();
155
156 switch (m_projectionType)
157 {
159 {
160 // Just copy over array
161 if (inarray != outarray)
162 {
163 int npoints = GetNpoints();
164
165 for (int j = 0; j < nvariables; ++j)
166 {
167 Vmath::Vcopy(npoints, inarray[j], 1, outarray[j], 1);
168 }
169 }
170 }
171 break;
174 {
176
177 for (int j = 0; j < nvariables; ++j)
178 {
179 m_fields[j]->FwdTrans(inarray[j], coeffs);
180 m_fields[j]->BwdTrans_IterPerExp(coeffs, outarray[j]);
181 }
182 break;
183 }
184 default:
185 ASSERTL0(false, "Unknown projection scheme");
186 break;
187 }
188}
189
190/**
191 * @brief Get the normal velocity
192 */
194{
195 // Number of trace (interface) points
196 int nTracePts = GetTraceNpoints();
197
198 // Auxiliary variable to compute the normal velocity
199 Array<OneD, NekDouble> tmp(nTracePts);
200
201 // Reset the normal velocity
202 Vmath::Zero(nTracePts, m_traceVn, 1);
203
204 for (int i = 0; i < m_velocity.size(); ++i)
205 {
206 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
207
208 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
209 m_traceVn, 1);
210 }
211
212 return m_traceVn;
213}
214
216 const Array<OneD, Array<OneD, NekDouble>> &physfield,
218{
219 ASSERTL1(flux[0].size() == m_velocity.size(),
220 "Dimension of flux array and velocity array do not match");
221
222 int nq = physfield[0].size();
223
224 for (int i = 0; i < flux.size(); ++i)
225 {
226 for (int j = 0; j < flux[0].size(); ++j)
227 {
228 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
229 }
230 }
231}
232
234{
236}
237
239 const Array<OneD, NekDouble> CFL,
240 NekDouble timeCFL)
241{
242
243 int n_element = m_fields[0]->GetExpSize();
244
245 Array<OneD, NekDouble> tstep(n_element, 0.0);
246 Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
247 stdVelocity = GetStdVelocity(m_velocity);
248
249 for (int el = 0; el < n_element; ++el)
250 {
251 int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
252 Array<OneD, NekDouble> one2D(npoints, 1.0);
253 if (std::dynamic_pointer_cast<LocalRegions::TriExp>(
254 m_fields[0]->GetExp(el)))
255 {
256 tstep[el] = CFL[el] / (stdVelocity[el]);
257 }
258 else if (std::dynamic_pointer_cast<LocalRegions::QuadExp>(
259 m_fields[0]->GetExp(el)))
260 {
261 tstep[el] = CFL[el] / (stdVelocity[el]);
262 }
263 }
264
265 NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
266
267 return TimeStep;
268}
269
271 NekDouble TimeStability)
272{
273 //================================================================
274 // This function has been created just to test specific problems, hence is
275 // not general and it has been implemented in a rude fashion, as the full
276 // CFLtester class. For real CFL calculations refer to the general
277 // implementation above. (A.Bolis)
278 //================================================================
279
280 NekDouble TimeStep;
281 NekDouble SpatialStability;
282 int n_elements = m_fields[0]->GetExpSize();
283
284 // solve ambiguity in windows
285 NekDouble n_elem = n_elements;
286 NekDouble DH = sqrt(n_elem);
287
288 int H = (int)DH;
289 int P = ExpOrder - 1;
290
291 //================================================================
292 // Regular meshes
293
294 // SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
295
296 //================================================================
297 // Anisotropic meshes
298
299 if (TimeStability == 1.0)
300 {
301 SpatialStability = EigenvaluesAnaMeshesAB2[H / 2][P - 1];
302 }
303 else if (TimeStability == 2.0)
304 {
305 SpatialStability = EigenvaluesAnaMeshesRK2[H / 2][P - 1];
306 }
307 else if (TimeStability == 2.784)
308 {
309 SpatialStability = EigenvaluesAnaMeshesRK4[H / 2][P - 1];
310 }
311 else
312 {
313 ASSERTL0(false, "Dominant eigenvalues database not present for this "
314 "time-stepping scheme")
315 }
316
317 //================================================================
318
319 TimeStep = (TimeStability / SpatialStability) * CFL;
320
321 //================================================================
322
323 return TimeStep;
324}
325
327 const Array<OneD, Array<OneD, NekDouble>> inarray)
328{
329 // Checking if the problem is 2D
330 ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
331
332 int nTotQuadPoints = GetTotPoints();
333 int n_element = m_fields[0]->GetExpSize();
334 int nvel = inarray.size();
335
336 NekDouble pntVelocity;
337
338 // Getting the standard velocity vector on the 2D normal space
339 Array<OneD, Array<OneD, NekDouble>> stdVelocity(nvel);
340 Array<OneD, NekDouble> stdV(n_element, 0.0);
341
342 for (int i = 0; i < nvel; ++i)
343 {
344 stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
345 }
346
347 if (nvel == 2)
348 {
349 for (int el = 0; el < n_element; ++el)
350 {
352 m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
353 int n_points = el2D->GetTotPoints();
354
356 el2D->GetGeom2D()->GetMetricInfo()->GetJac();
358 el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
359
360 if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype() ==
362 {
363 for (int i = 0; i < n_points; i++)
364 {
365 stdVelocity[0][i] =
366 gmat[0][i] * inarray[0][i] + gmat[2][i] * inarray[1][i];
367
368 stdVelocity[1][i] =
369 gmat[1][i] * inarray[0][i] + gmat[3][i] * inarray[1][i];
370 }
371 }
372 else
373 {
374 for (int i = 0; i < n_points; i++)
375 {
376 stdVelocity[0][i] =
377 gmat[0][0] * inarray[0][i] + gmat[2][0] * inarray[1][i];
378
379 stdVelocity[1][i] =
380 gmat[1][0] * inarray[0][i] + gmat[3][0] * inarray[1][i];
381 }
382 }
383
384 for (int i = 0; i < n_points; i++)
385 {
386 pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
387 stdVelocity[1][i] * stdVelocity[1][i]);
388
389 if (pntVelocity > stdV[el])
390 {
391 stdV[el] = pntVelocity;
392 }
393 }
394 }
395 }
396 else
397 {
398 for (int el = 0; el < n_element; ++el)
399 {
401 m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
402
403 int n_points = el3D->GetTotPoints();
404
406 el3D->GetGeom3D()->GetMetricInfo()->GetJac();
408 el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
409
410 if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype() ==
412 {
413 for (int i = 0; i < n_points; i++)
414 {
415 stdVelocity[0][i] = gmat[0][i] * inarray[0][i] +
416 gmat[3][i] * inarray[1][i] +
417 gmat[6][i] * inarray[2][i];
418
419 stdVelocity[1][i] = gmat[1][i] * inarray[0][i] +
420 gmat[4][i] * inarray[1][i] +
421 gmat[7][i] * inarray[2][i];
422
423 stdVelocity[2][i] = gmat[2][i] * inarray[0][i] +
424 gmat[5][i] * inarray[1][i] +
425 gmat[8][i] * inarray[2][i];
426 }
427 }
428 else
429 {
431 el3D->GetGeom3D()->GetMetricInfo()->GetJac();
433 el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
434
435 for (int i = 0; i < n_points; i++)
436 {
437 stdVelocity[0][i] = gmat[0][0] * inarray[0][i] +
438 gmat[3][0] * inarray[1][i] +
439 gmat[6][0] * inarray[2][i];
440
441 stdVelocity[1][i] = gmat[1][0] * inarray[0][i] +
442 gmat[4][0] * inarray[1][i] +
443 gmat[7][0] * inarray[2][i];
444
445 stdVelocity[2][i] = gmat[2][0] * inarray[0][i] +
446 gmat[5][0] * inarray[1][i] +
447 gmat[8][0] * inarray[2][i];
448 }
449 }
450
451 for (int i = 0; i < n_points; i++)
452 {
453 pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
454 stdVelocity[1][i] * stdVelocity[1][i] +
455 stdVelocity[2][i] * stdVelocity[2][i]);
456
457 if (pntVelocity > stdV[el])
458 {
459 stdV[el] = pntVelocity;
460 }
461 }
462 }
463 }
464
465 return stdV;
466}
467} // 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
Array< OneD, NekDouble > GetStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Definition: CFLtester.cpp:326
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
virtual NekDouble v_GetTimeStep(const Array< OneD, int > ExpOrder, const Array< OneD, NekDouble > CFL, NekDouble timeCFL) override
Definition: CFLtester.cpp:238
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:233
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:193
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: CFLtester.cpp:215
static std::string className
Definition: CFLtester.h:152
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:134
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
@ P
Monomial polynomials .
Definition: BasisType.h:64
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:48
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:51
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:176
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:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:1045
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:569
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294