Nektar++
CompressibleFlowSystem.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: CompressibleFlowSystem.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: Compressible flow system base class with auxiliary functions
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
39
40namespace Nektar
41{
42
46 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
47{
48}
49
50/**
51 * @brief Initialization object for CompressibleFlowSystem class.
52 */
54{
55 AdvectionSystem::v_InitObject(DeclareFields);
56
57 for (size_t i = 0; i < m_fields.size(); i++)
58 {
59 // Use BwdTrans to make sure initial condition is in solution space
60 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
61 m_fields[i]->UpdatePhys());
62 }
63
66
67 ASSERTL0(m_session->DefinesSolverInfo("UPWINDTYPE"),
68 "No UPWINDTYPE defined in session.");
69
70 // Do not forwards transform initial condition
71 m_homoInitialFwd = false;
72
73 // Set up locations of velocity vector.
76 for (int i = 0; i < m_spacedim; ++i)
77 {
78 m_vecLocs[0][i] = 1 + i;
79 }
80
81 // Loading parameters from session file
83
84 // Setting up advection and diffusion operators
86
87 // Create artificial diffusion with laplacian operator
88 if (m_shockCaptureType == "NonSmooth")
89 {
92 }
93
94 // Forcing terms for the sponge region
96 m_fields, m_fields.size());
97
98 // User-defined boundary conditions
99 size_t cnt = 0;
100 for (size_t n = 0; n < (size_t)m_fields[0]->GetBndConditions().size(); ++n)
101 {
102 std::string type = m_fields[0]->GetBndConditions()[n]->GetUserDefined();
103
104 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
106 {
107 continue;
108 }
109
110 if (!type.empty())
111 {
114 m_spacedim, n, cnt));
115 }
116 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
117 }
118
121
123}
124
125/**
126 * @brief Load CFS parameters from the session file
127 */
129{
130 // Get gamma parameter from session file.
131 m_session->LoadParameter("Gamma", m_gamma, 1.4);
132
133 // Shock capture
134 m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
135
136 // Check if the shock capture type is supported
137 std::string err_msg = "Warning, ShockCaptureType = " + m_shockCaptureType +
138 " is not supported by this solver";
140
141 // Load parameters for exponential filtering
142 m_session->MatchSolverInfo("ExponentialFiltering", "True", m_useFiltering,
143 false);
144 if (m_useFiltering)
145 {
146 m_session->LoadParameter("FilterAlpha", m_filterAlpha, 36);
147 m_session->LoadParameter("FilterExponent", m_filterExponent, 16);
148 m_session->LoadParameter("FilterCutoff", m_filterCutoff, 0);
149 }
150
151 // Load CFL for local time-stepping (for steady state)
152 m_session->MatchSolverInfo("LocalTimeStep", "True", m_useLocalTimeStep,
153 false);
155 {
157 "Local time stepping requires CFL parameter.");
158 }
159}
160
161/**
162 * @brief Create advection and diffusion objects for CFS
163 */
165{
166 // Check if projection type is correct
168 "Unsupported projection type.");
169
170 std::string advName, riemName;
171 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
172
175
177 {
178 m_advObject->SetFluxVector(
180 }
181 else
182 {
184 this);
185 }
186
187 // Setting up Riemann solver for advection operator
188 m_session->LoadSolverInfo("UpwindType", riemName, "Average");
189
192 riemName, m_session);
193
194 // Tell Riemann Solver if doing ALE and provide trace grid velocity
195 riemannSolver->SetALEFlag(m_ALESolver);
196 riemannSolver->SetVector("vgt", &ALEHelper::GetGridVelocityTrace, this);
197
198 // Setting up parameters for advection operator Riemann solver
199 riemannSolver->SetParam("gamma", &CompressibleFlowSystem::GetGamma, this);
200 riemannSolver->SetAuxVec("vecLocs", &CompressibleFlowSystem::GetVecLocs,
201 this);
202 riemannSolver->SetVector("N", &CompressibleFlowSystem::GetNormals, this);
203
204 // Concluding initialisation of advection / diffusion operators
205 m_advObject->SetRiemannSolver(riemannSolver);
206 m_advObject->InitObject(m_session, m_fields);
207}
208
209/**
210 * @brief Compute the right-hand side.
211 */
213 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
214 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
215{
217
218 size_t nvariables = inarray.size();
219 size_t nTracePts = GetTraceTotPoints();
220
221 // This converts our Mu in coefficient space to u in physical space for ALE
222 Array<OneD, Array<OneD, NekDouble>> tmpIn(nvariables);
223 if (m_ALESolver)
224 {
226 }
227 else
228 {
229 tmpIn = inarray;
230 }
231
232 m_bndEvaluateTime = time;
233
234 // Store forwards/backwards space along trace space
235 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
236 Array<OneD, Array<OneD, NekDouble>> Bwd(nvariables);
237
239 {
242 }
243 else
244 {
245 for (size_t i = 0; i < nvariables; ++i)
246 {
247 Fwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
248 Bwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
249 m_fields[i]->GetFwdBwdTracePhys(tmpIn[i], Fwd[i], Bwd[i]);
250 }
251 }
252
253 // Calculate advection
254 timer.Start();
255 DoAdvection(tmpIn, outarray, time, Fwd, Bwd);
256 timer.Stop();
257 timer.AccumulateRegion("DoAdvection");
258
259 // Negate results
260 for (size_t i = 0; i < nvariables; ++i)
261 {
262 Vmath::Neg(outarray[i].size(), outarray[i], 1);
263 }
264
265 // Add diffusion terms
266 timer.Start();
267 DoDiffusion(tmpIn, outarray, Fwd, Bwd);
268 timer.Stop();
269 timer.AccumulateRegion("DoDiffusion");
270
271 // Add forcing terms
272 for (auto &x : m_forcing)
273 {
274 x->Apply(m_fields, tmpIn, outarray, time);
275 }
276
278 {
279 size_t nElements = m_fields[0]->GetExpSize();
280 int nq, offset;
281 NekDouble fac;
283
284 Array<OneD, NekDouble> tstep(nElements, 0.0);
285 GetElmtTimeStep(tmpIn, tstep);
286
287 // Loop over elements
288 for (size_t n = 0; n < nElements; ++n)
289 {
290 nq = m_fields[0]->GetExp(n)->GetTotPoints();
291 offset = m_fields[0]->GetPhys_Offset(n);
292 fac = tstep[n] / m_timestep;
293 for (size_t i = 0; i < nvariables; ++i)
294 {
295 Vmath::Smul(nq, fac, outarray[i] + offset, 1,
296 tmp = outarray[i] + offset, 1);
297 }
298 }
299 }
300}
301
302/**
303 * @brief Compute the projection and call the method for imposing the
304 * boundary conditions in case of discontinuous projection.
305 */
307 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
308 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
309{
310 size_t nvariables = inarray.size();
311
312 // Perform ALE movement
313 if (m_ALESolver)
314 {
316 }
317
318 switch (m_projectionType)
319 {
321 {
322 // Just copy over array
323 for (size_t i = 0; i < nvariables; ++i)
324 {
325 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
326
327 if (m_useFiltering)
328 {
329 m_fields[i]->ExponentialFilter(outarray[i], m_filterAlpha,
332 }
333 }
334 SetBoundaryConditions(outarray, time);
335 break;
336 }
339 {
340 NEKERROR(ErrorUtil::efatal, "No Continuous Galerkin for full "
341 "compressible Navier-Stokes equations");
342 break;
343 }
344 default:
345 NEKERROR(ErrorUtil::efatal, "Unknown projection scheme");
346 break;
347 }
348}
349
350/**
351 * @brief Compute the advection terms for the right-hand side
352 */
354 const Array<OneD, Array<OneD, NekDouble>> &inarray,
355 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
356 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
357 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
358{
359 int nvariables = inarray.size();
361
362 if (m_ALESolver)
363 {
364 auto advWeakDGObject =
365 std::dynamic_pointer_cast<SolverUtils::AdvectionWeakDG>(
367 advWeakDGObject->AdvectCoeffs(nvariables, m_fields, advVel, inarray,
368 outarray, time, pFwd, pBwd);
369 }
370 else
371 {
372 m_advObject->Advect(nvariables, m_fields, advVel, inarray, outarray,
373 time, pFwd, pBwd);
374 }
375}
376
377/**
378 * @brief Add the diffusions terms to the right-hand side
379 */
381 const Array<OneD, Array<OneD, NekDouble>> &inarray,
383 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
384 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
385{
386 v_DoDiffusion(inarray, outarray, pFwd, pBwd);
387}
388
390 Array<OneD, Array<OneD, NekDouble>> &physarray, NekDouble time)
391{
392 size_t nTracePts = GetTraceTotPoints();
393 size_t nvariables = physarray.size();
394
395 // This converts our Mu in coefficient space to u in physical space for ALE
396 Array<OneD, Array<OneD, NekDouble>> tmpIn(nvariables);
397
399 {
400 ALEHelper::ALEDoElmtInvMassBwdTrans(physarray, tmpIn);
401 }
402 else
403 {
404 tmpIn = physarray;
405 }
406
407 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
408 for (size_t i = 0; i < nvariables; ++i)
409 {
410 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
411 m_fields[i]->ExtractTracePhys(tmpIn[i], Fwd[i]);
412 }
413
414 if (!m_bndConds.empty())
415 {
416 // Loop over user-defined boundary conditions
417 for (auto &x : m_bndConds)
418 {
419 x->Apply(Fwd, tmpIn, time);
420 }
421 }
422}
423/**
424 * @brief Set up a weight on physical boundaries for boundary condition
425 * applications
426 */
428{
429 if (m_bndConds.size())
430 {
431 // Loop over user-defined boundary conditions
432 for (auto &x : m_bndConds)
433 {
434 x->ApplyBwdWeight();
435 }
436 }
437}
438
439/**
440 * @brief Return the flux vector for the compressible Euler equations.
441 *
442 * @param physfield Fields.
443 * @param flux Resulting flux.
444 */
446 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
448{
449 size_t nVariables = physfield.size();
450 size_t nPts = physfield[0].size();
451
452 constexpr unsigned short maxVel = 3;
453 constexpr unsigned short maxFld = 5;
454
455 // hardcoding done for performance reasons
456 ASSERTL1(nVariables <= maxFld, "GetFluxVector, hard coded max fields");
457
458 for (size_t p = 0; p < nPts; ++p)
459 {
460 // local storage
461 std::array<NekDouble, maxFld> fieldTmp = {0};
462 std::array<NekDouble, maxVel> velocity = {0};
463
464 // rearrenge and load data
465 for (size_t f = 0; f < nVariables; ++f)
466 {
467 fieldTmp[f] = physfield[f][p]; // load
468 }
469
470 // 1 / rho
471 NekDouble oneOrho = 1.0 / fieldTmp[0];
472
473 for (size_t d = 0; d < m_spacedim; ++d)
474 {
475 // Flux vector for the rho equation
476 flux[0][d][p] = fieldTmp[d + 1]; // store
477 // compute velocity
478 velocity[d] = fieldTmp[d + 1] * oneOrho;
479 }
480
481 NekDouble pressure = m_varConv->GetPressure(fieldTmp.data());
482 NekDouble ePlusP = fieldTmp[m_spacedim + 1] + pressure;
483 for (size_t f = 0; f < m_spacedim; ++f)
484 {
485 // Flux vector for the velocity fields
486 for (size_t d = 0; d < m_spacedim; ++d)
487 {
488 flux[f + 1][d][p] = velocity[d] * fieldTmp[f + 1]; // store
489 }
490
491 // Add pressure to appropriate field
492 flux[f + 1][f][p] += pressure;
493
494 // Flux vector for energy
495 flux[m_spacedim + 1][f][p] = ePlusP * velocity[f]; // store
496 }
497 }
498
499 // @TODO : for each row (3 columns) negative grid velocity component (d * d
500 // + 2 (4 rows)) rho, rhou, rhov, rhow, E,
501 // @TODO : top row is flux for rho etc... each row subtract v_g * conserved
502 // variable for that row... For grid velocity subtract v_g * conserved
503 // variable
504 if (m_ALESolver)
505 {
506 for (int i = 0; i < m_spacedim + 2; ++i)
507 {
508 for (int j = 0; j < m_spacedim; ++j)
509 {
510 for (int k = 0; k < nPts; ++k)
511 {
512 flux[i][j][k] -= physfield[i][k] * m_gridVelocity[j][k];
513 }
514 }
515 }
516 }
517}
518
519/**
520 * @brief Return the flux vector for the compressible Euler equations
521 * by using the de-aliasing technique.
522 *
523 * @param physfield Fields.
524 * @param flux Resulting flux.
525 */
527 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
529{
530 int i, j;
531 size_t nq = physfield[0].size();
532 size_t nVariables = m_fields.size();
533
534 // Factor to rescale 1d points in dealiasing
535 NekDouble OneDptscale = 2;
536 nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
537
540
541 Array<OneD, Array<OneD, NekDouble>> physfield_interp(nVariables);
542 TensorOfArray3D<NekDouble> flux_interp(nVariables);
543
544 for (size_t i = 0; i < nVariables; ++i)
545 {
546 physfield_interp[i] = Array<OneD, NekDouble>(nq);
548 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
549 physfield_interp[i]);
550
551 for (j = 0; j < m_spacedim; ++j)
552 {
553 flux_interp[i][j] = Array<OneD, NekDouble>(nq);
554 }
555 }
556
557 // Flux vector for the rho equation
558 for (i = 0; i < m_spacedim; ++i)
559 {
560 velocity[i] = Array<OneD, NekDouble>(nq);
561
562 // Galerkin project solution back to original space
563 m_fields[0]->PhysGalerkinProjection1DScaled(
564 OneDptscale, physfield_interp[i + 1], flux[0][i]);
565 }
566
567 m_varConv->GetVelocityVector(physfield_interp, velocity);
568 m_varConv->GetPressure(physfield_interp, pressure);
569
570 // Evaluation of flux vector for the velocity fields
571 for (i = 0; i < m_spacedim; ++i)
572 {
573 for (j = 0; j < m_spacedim; ++j)
574 {
575 Vmath::Vmul(nq, velocity[j], 1, physfield_interp[i + 1], 1,
576 flux_interp[i + 1][j], 1);
577 }
578
579 // Add pressure to appropriate field
580 Vmath::Vadd(nq, flux_interp[i + 1][i], 1, pressure, 1,
581 flux_interp[i + 1][i], 1);
582 }
583
584 // Galerkin project solution back to original space
585 for (i = 0; i < m_spacedim; ++i)
586 {
587 for (j = 0; j < m_spacedim; ++j)
588 {
589 m_fields[0]->PhysGalerkinProjection1DScaled(
590 OneDptscale, flux_interp[i + 1][j], flux[i + 1][j]);
591 }
592 }
593
594 // Evaluation of flux vector for energy
595 Vmath::Vadd(nq, physfield_interp[m_spacedim + 1], 1, pressure, 1, pressure,
596 1);
597
598 for (j = 0; j < m_spacedim; ++j)
599 {
600 Vmath::Vmul(nq, velocity[j], 1, pressure, 1,
601 flux_interp[m_spacedim + 1][j], 1);
602
603 // Galerkin project solution back to original space
604 m_fields[0]->PhysGalerkinProjection1DScaled(
605 OneDptscale, flux_interp[m_spacedim + 1][j],
606 flux[m_spacedim + 1][j]);
607 }
608}
609
610/**
611 * @brief Calculate the maximum timestep on each element
612 * subject to CFL restrictions.
613 */
615 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray,
617{
618 size_t nElements = m_fields[0]->GetExpSize();
619
620 // Change value of m_timestep (in case it is set to zero)
621 NekDouble tmp = m_timestep;
622 m_timestep = 1.0;
623
624 Array<OneD, NekDouble> cfl(nElements);
625 cfl = GetElmtCFLVals();
626
627 // Factors to compute the time-step limit
629
630 // Loop over elements to compute the time-step limit for each element
631 for (size_t n = 0; n < nElements; ++n)
632 {
633 tstep[n] = m_cflSafetyFactor * alpha / cfl[n];
634 }
635
636 // Restore value of m_timestep
637 m_timestep = tmp;
638}
639
640/**
641 * @brief Calculate the maximum timestep subject to CFL restrictions.
642 */
644 const Array<OneD, const Array<OneD, NekDouble>> &inarray)
645{
646 int nElements = m_fields[0]->GetExpSize();
647 Array<OneD, NekDouble> tstep(nElements, 0.0);
648
649 GetElmtTimeStep(inarray, tstep);
650
651 // Get the minimum time-step limit and return the time-step
652 NekDouble TimeStep = Vmath::Vmin(nElements, tstep, 1);
653 m_comm->GetSpaceComm()->AllReduce(TimeStep, LibUtilities::ReduceMin);
654
655 return TimeStep;
656}
657
659{
661 if (m_session->DefinesSolverInfo("ICTYPE"))
662 {
663 if (boost::iequals(m_session->GetSolverInfo("ICTYPE"),
664 "IsentropicVortex"))
665 {
666 SolverUtils::AddSummaryItem(s, "Problem Type", "IsentropicVortex");
667 }
668 else if (boost::iequals(m_session->GetSolverInfo("ICTYPE"),
669 "RinglebFlow"))
670 {
671 SolverUtils::AddSummaryItem(s, "Problem Type", "RinglebFlow");
672 }
673 else
674 {
675 NEKERROR(ErrorUtil::efatal, "unknow initial condition");
676 }
677 }
678}
679
680/**
681 * @brief Set up logic for residual calculation.
682 */
684 NekDouble initialtime, bool dumpInitialConditions,
685 [[maybe_unused]] const int domain)
686{
687 if (m_session->DefinesSolverInfo("ICTYPE") &&
688 boost::iequals(m_session->GetSolverInfo("ICTYPE"), "IsentropicVortex"))
689 {
690 // Forward transform to fill the coefficient space
691 for (int i = 0; i < m_fields.size(); ++i)
692 {
693 EvaluateIsentropicVortex(i, m_fields[i]->UpdatePhys(), initialtime);
694 m_fields[i]->SetPhysState(true);
695 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
696 m_fields[i]->UpdateCoeffs());
697 }
698 }
699 else
700 {
701 EquationSystem::v_SetInitialConditions(initialtime, false);
702 m_nchk--; // Note: m_nchk has been incremented in EquationSystem.
703
704 // insert white noise in initial condition
705 NekDouble Noise;
706 int phystot = m_fields[0]->GetTotPoints();
707 Array<OneD, NekDouble> noise(phystot);
708
709 m_session->LoadParameter("Noise", Noise, 0.0);
710 int m_nConvectiveFields = m_fields.size();
711
712 if (Noise > 0.0)
713 {
714 int seed = -m_comm->GetSpaceComm()->GetRank() * m_nConvectiveFields;
715 for (int i = 0; i < m_nConvectiveFields; i++)
716 {
717 Vmath::FillWhiteNoise(phystot, Noise, noise, 1, seed);
718 --seed;
719 Vmath::Vadd(phystot, m_fields[i]->GetPhys(), 1, noise, 1,
720 m_fields[i]->UpdatePhys(), 1);
721 m_fields[i]->FwdTransLocalElmt(m_fields[i]->GetPhys(),
722 m_fields[i]->UpdateCoeffs());
723 }
724 }
725 }
726
727 if (dumpInitialConditions && m_nchk == 0 && m_checksteps &&
728 !m_comm->IsParallelInTime())
729 {
731 }
732 else if (dumpInitialConditions && m_nchk == 0 && m_comm->IsParallelInTime())
733 {
734 std::string newdir = m_sessionName + ".pit";
735 if (!fs::is_directory(newdir))
736 {
737 fs::create_directory(newdir);
738 }
739 if (m_comm->GetTimeComm()->GetRank() == 0)
740 {
741 WriteFld(newdir + "/" + m_sessionName + "_0.fld");
742 }
743 }
744 m_nchk++;
745}
746
748 unsigned int field, Array<OneD, NekDouble> &outfield, const NekDouble time)
749{
750
751 if (!m_session->DefinesFunction("ExactSolution") &&
752 m_session->DefinesSolverInfo("ICTYPE"))
753 {
754 if (boost::iequals(m_session->GetSolverInfo("ICTYPE"),
755 "IsentropicVortex"))
756 {
757 EvaluateIsentropicVortex(field, outfield, time);
758 }
759 else if (boost::iequals(m_session->GetSolverInfo("ICTYPE"),
760 "RinglebFlow"))
761 {
762 GetExactRinglebFlow(field, outfield);
763 }
764 }
765 else
766 {
768 }
769}
770
771/**
772 * @brief Compute the advection velocity in the standard space
773 * for each element of the expansion.
774 */
776 const NekDouble SpeedSoundFactor)
777{
778 size_t nTotQuadPoints = GetTotPoints();
779 size_t n_element = m_fields[0]->GetExpSize();
780 size_t expdim = m_fields[0]->GetGraph()->GetMeshDimension();
781 size_t nfields = m_fields.size();
782 int offset;
784
785 Array<OneD, Array<OneD, NekDouble>> physfields(nfields);
786 for (size_t i = 0; i < nfields; ++i)
787 {
788 physfields[i] = m_fields[i]->GetPhys();
789 }
790
791 Array<OneD, NekDouble> stdV(n_element, 0.0);
792
793 // Getting the velocity vector on the 2D normal space
797 Array<OneD, NekDouble> soundspeed(nTotQuadPoints);
799
800 for (int i = 0; i < m_spacedim; ++i)
801 {
802 velocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
803 stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints, 0.0);
804 stdSoundSpeed[i] = Array<OneD, NekDouble>(nTotQuadPoints, 0.0);
805 }
806
807 m_varConv->GetVelocityVector(physfields, velocity);
808 m_varConv->GetSoundSpeed(physfields, soundspeed);
809
810 // Subtract Ug from the velocity for the ALE formulation
811 for (size_t i = 0; i < m_spacedim; ++i)
812 {
813 Vmath::Vsub(nTotQuadPoints, velocity[i], 1, m_gridVelocity[i], 1,
814 velocity[i], 1);
815 }
816
817 for (int el = 0; el < n_element; ++el)
818 {
819 ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
820 offset = m_fields[0]->GetPhys_Offset(el);
821 int nq = m_fields[0]->GetExp(el)->GetTotPoints();
822
823 const SpatialDomains::GeomFactorsSharedPtr metricInfo =
824 m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo();
825 const Array<TwoD, const NekDouble> &gmat =
826 m_fields[0]
827 ->GetExp(el)
828 ->GetGeom()
829 ->GetMetricInfo()
830 ->GetDerivFactors(ptsKeys);
831
832 // Convert to standard element
833 // consider soundspeed in all directions
834 // (this might overestimate the cfl)
835 if (metricInfo->GetGtype() == SpatialDomains::eDeformed)
836 {
837 // d xi/ dx = gmat = 1/J * d x/d xi
838 for (size_t i = 0; i < expdim; ++i)
839 {
840 Vmath::Vmul(nq, gmat[i], 1, velocity[0] + offset, 1,
841 tmp = stdVelocity[i] + offset, 1);
842 Vmath::Vmul(nq, gmat[i], 1, soundspeed + offset, 1,
843 tmp = stdSoundSpeed[i] + offset, 1);
844 for (size_t j = 1; j < expdim; ++j)
845 {
846 Vmath::Vvtvp(nq, gmat[expdim * j + i], 1,
847 velocity[j] + offset, 1,
848 stdVelocity[i] + offset, 1,
849 tmp = stdVelocity[i] + offset, 1);
850 Vmath::Vvtvp(nq, gmat[expdim * j + i], 1,
851 soundspeed + offset, 1,
852 stdSoundSpeed[i] + offset, 1,
853 tmp = stdSoundSpeed[i] + offset, 1);
854 }
855 }
856 }
857 else
858 {
859 for (size_t i = 0; i < expdim; ++i)
860 {
861 Vmath::Smul(nq, gmat[i][0], velocity[0] + offset, 1,
862 tmp = stdVelocity[i] + offset, 1);
863 Vmath::Smul(nq, gmat[i][0], soundspeed + offset, 1,
864 tmp = stdSoundSpeed[i] + offset, 1);
865 for (size_t j = 1; j < expdim; ++j)
866 {
867 Vmath::Svtvp(nq, gmat[expdim * j + i][0],
868 velocity[j] + offset, 1,
869 stdVelocity[i] + offset, 1,
870 tmp = stdVelocity[i] + offset, 1);
871 Vmath::Svtvp(nq, gmat[expdim * j + i][0],
872 soundspeed + offset, 1,
873 stdSoundSpeed[i] + offset, 1,
874 tmp = stdSoundSpeed[i] + offset, 1);
875 }
876 }
877 }
878
879 NekDouble vel;
880 for (size_t i = 0; i < nq; ++i)
881 {
882 NekDouble pntVelocity = 0.0;
883 for (size_t j = 0; j < expdim; ++j)
884 {
885 // Add sound speed
886 vel = std::abs(stdVelocity[j][offset + i]) +
887 SpeedSoundFactor * std::abs(stdSoundSpeed[j][offset + i]);
888 pntVelocity += vel * vel;
889 }
890 pntVelocity = sqrt(pntVelocity);
891 if (pntVelocity > stdV[el])
892 {
893 stdV[el] = pntVelocity;
894 }
895 }
896 }
897
898 return stdV;
899}
900
901/**
902 * @brief Set the denominator to compute the time step when a cfl
903 * control is employed. This function is no longer used but is still
904 * here for being utilised in the future.
905 *
906 * @param n Order of expansion element by element.
907 */
909{
910 ASSERTL0(n <= 20, "Illegal modes dimension for CFL calculation "
911 "(P has to be less then 20)");
912
913 NekDouble CFLDG[21] = {2.0000, 6.0000, 11.8424, 19.1569, 27.8419,
914 37.8247, 49.0518, 61.4815, 75.0797, 89.8181,
915 105.6700, 122.6200, 140.6400, 159.7300, 179.8500,
916 201.0100, 223.1800, 246.3600, 270.5300, 295.6900,
917 321.8300}; // CFLDG 1D [0-20]
918 NekDouble CFL = 0.0;
919
921 {
922 CFL = CFLDG[n];
923 }
924 else
925 {
926 NEKERROR(ErrorUtil::efatal, "Continuous Galerkin stability "
927 "coefficients not introduced yet.");
928 }
929
930 return CFL;
931}
932
933/**
934 * @brief Compute the vector of denominators to compute the time step
935 * when a cfl control is employed. This function is no longer used but
936 * is still here for being utilised in the future.
937 *
938 * @param ExpOrder Order of expansion element by element.
939 */
941 const Array<OneD, int> &ExpOrder)
942{
943 int i;
944 Array<OneD, NekDouble> returnval(m_fields[0]->GetExpSize(), 0.0);
945 for (i = 0; i < m_fields[0]->GetExpSize(); i++)
946 {
947 returnval[i] = GetStabilityLimit(ExpOrder[i]);
948 }
949 return returnval;
950}
951
953 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
954 std::vector<std::string> &variables)
955{
956 bool extraFields;
957 m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
958 if (extraFields)
959 {
960 const int nPhys = m_fields[0]->GetNpoints();
961 const int nCoeffs = m_fields[0]->GetNcoeffs();
963
964 for (int i = 0; i < m_fields.size(); ++i)
965 {
966 tmp[i] = m_fields[i]->GetPhys();
967 }
968
971 for (int i = 0; i < m_spacedim; ++i)
972 {
973 velocity[i] = Array<OneD, NekDouble>(nPhys);
974 velFwd[i] = Array<OneD, NekDouble>(nCoeffs);
975 }
976
977 Array<OneD, NekDouble> pressure(nPhys), temperature(nPhys);
978 Array<OneD, NekDouble> entropy(nPhys);
979 Array<OneD, NekDouble> soundspeed(nPhys), mach(nPhys);
980 Array<OneD, NekDouble> sensor(nPhys), SensorKappa(nPhys);
981
982 m_varConv->GetVelocityVector(tmp, velocity);
983 m_varConv->GetPressure(tmp, pressure);
984 m_varConv->GetTemperature(tmp, temperature);
985 m_varConv->GetEntropy(tmp, entropy);
986 m_varConv->GetSoundSpeed(tmp, soundspeed);
987 m_varConv->GetMach(tmp, soundspeed, mach);
988
989 int sensorOffset;
990 m_session->LoadParameter("SensorOffset", sensorOffset, 1);
991 m_varConv->GetSensor(m_fields[0], tmp, sensor, SensorKappa,
992 sensorOffset);
993
994 Array<OneD, NekDouble> pFwd(nCoeffs), TFwd(nCoeffs);
995 Array<OneD, NekDouble> sFwd(nCoeffs);
996 Array<OneD, NekDouble> aFwd(nCoeffs), mFwd(nCoeffs);
997 Array<OneD, NekDouble> sensFwd(nCoeffs);
998
999 std::string velNames[3] = {"u", "v", "w"};
1000 for (int i = 0; i < m_spacedim; ++i)
1001 {
1002 m_fields[0]->FwdTransLocalElmt(velocity[i], velFwd[i]);
1003 variables.push_back(velNames[i]);
1004 fieldcoeffs.push_back(velFwd[i]);
1005 }
1006
1007 m_fields[0]->FwdTransLocalElmt(pressure, pFwd);
1008 m_fields[0]->FwdTransLocalElmt(temperature, TFwd);
1009 m_fields[0]->FwdTransLocalElmt(entropy, sFwd);
1010 m_fields[0]->FwdTransLocalElmt(soundspeed, aFwd);
1011 m_fields[0]->FwdTransLocalElmt(mach, mFwd);
1012 m_fields[0]->FwdTransLocalElmt(sensor, sensFwd);
1013
1014 variables.push_back("p");
1015 variables.push_back("T");
1016 variables.push_back("s");
1017 variables.push_back("a");
1018 variables.push_back("Mach");
1019 variables.push_back("Sensor");
1020 fieldcoeffs.push_back(pFwd);
1021 fieldcoeffs.push_back(TFwd);
1022 fieldcoeffs.push_back(sFwd);
1023 fieldcoeffs.push_back(aFwd);
1024 fieldcoeffs.push_back(mFwd);
1025 fieldcoeffs.push_back(sensFwd);
1026
1028 {
1029 // reuse pressure
1030 Array<OneD, NekDouble> sensorFwd(nCoeffs);
1031 m_artificialDiffusion->GetArtificialViscosity(tmp, pressure);
1032 m_fields[0]->FwdTransLocalElmt(pressure, sensorFwd);
1033
1034 variables.push_back("ArtificialVisc");
1035 fieldcoeffs.push_back(sensorFwd);
1036 }
1037
1038 if (m_ALESolver)
1039 {
1040 ExtraFldOutputGridVelocity(fieldcoeffs, variables);
1041 }
1042 }
1043}
1044
1045/**
1046 *
1047 */
1049 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
1051{
1052 m_varConv->GetPressure(physfield, pressure);
1053}
1054
1055/**
1056 *
1057 */
1059 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
1060 Array<OneD, NekDouble> &density)
1061{
1062 density = physfield[0];
1063}
1064
1065/**
1066 *
1067 */
1069 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
1070 Array<OneD, Array<OneD, NekDouble>> &velocity)
1071{
1072 m_varConv->GetVelocityVector(physfield, velocity);
1073}
1074
1077{
1078 const size_t nPoints = GetTotPoints();
1079 const size_t nFields = m_fields.size();
1081 Array<OneD, Array<OneD, NekDouble>> inarray(nFields);
1082 for (size_t i = 0; i < nFields; ++i)
1083 {
1084 rhs[i] = Array<OneD, NekDouble>(nPoints, 0.0);
1085 inarray[i] = m_fields[i]->UpdatePhys();
1086 }
1087
1088 DoOdeRhs(inarray, rhs, m_time);
1089
1090 // Holds L2 errors.
1092 Array<OneD, NekDouble> RHSL2(nFields);
1093 Array<OneD, NekDouble> residual(nFields);
1094
1095 for (size_t i = 0; i < nFields; ++i)
1096 {
1097 tmp = rhs[i];
1098
1099 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
1100 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
1101 }
1102
1103 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
1104
1105 NekDouble onPoints = 1.0 / NekDouble(nPoints);
1106 for (size_t i = 0; i < nFields; ++i)
1107 {
1108 L2[i] = sqrt(residual[i] * onPoints);
1109 }
1110}
1111
1112/**
1113 *
1114 */
1116 unsigned int field, Array<OneD, NekDouble> &outfield, NekDouble time,
1117 const int o)
1118{
1119 NekDouble beta, u0, v0, x0, y0;
1120
1121 int nTotQuadPoints = GetTotPoints();
1122 Array<OneD, NekDouble> x(nTotQuadPoints);
1123 Array<OneD, NekDouble> y(nTotQuadPoints);
1124 Array<OneD, NekDouble> z(nTotQuadPoints);
1126
1127 m_fields[0]->GetCoords(x, y, z);
1128
1129 for (int i = 0; i < m_spacedim + 2; ++i)
1130 {
1131 u[i] = Array<OneD, NekDouble>(nTotQuadPoints);
1132 }
1133 m_session->LoadParameter("IsentropicBeta", beta, 5.0);
1134 m_session->LoadParameter("IsentropicU0", u0, 1.0);
1135 m_session->LoadParameter("IsentropicV0", v0, 0.5);
1136 m_session->LoadParameter("IsentropicX0", x0, 5.0);
1137 m_session->LoadParameter("IsentropicY0", y0, 0.0);
1138
1139 // Flow parameters
1140 NekDouble r, xbar, ybar, tmp;
1141 NekDouble fac = 1.0 / (16.0 * m_gamma * M_PI * M_PI);
1142
1143 // In 3D zero rhow field.
1144 if (m_spacedim == 3)
1145 {
1146 Vmath::Zero(nTotQuadPoints, &u[3][o], 1);
1147 }
1148
1149 // Fill storage
1150 for (int i = 0; i < nTotQuadPoints; ++i)
1151 {
1152 xbar = x[i] - u0 * time - x0;
1153 ybar = y[i] - v0 * time - y0;
1154 r = sqrt(xbar * xbar + ybar * ybar);
1155 tmp = beta * exp(1 - r * r);
1156 u[0][i + o] =
1157 pow(1.0 - (m_gamma - 1.0) * tmp * tmp * fac, 1.0 / (m_gamma - 1.0));
1158 u[1][i + o] = u[0][i + o] * (u0 - tmp * ybar / (2 * M_PI));
1159 u[2][i + o] = u[0][i + o] * (v0 + tmp * xbar / (2 * M_PI));
1160 u[m_spacedim + 1][i + o] =
1161 pow(u[0][i + o], m_gamma) / (m_gamma - 1.0) +
1162 0.5 * (u[1][i + o] * u[1][i + o] + u[2][i + o] * u[2][i + o]) /
1163 u[0][i + o];
1164 }
1165 Vmath::Vcopy(nTotQuadPoints, u[field].data(), 1, outfield.data(), 1);
1166}
1167
1168/**
1169 * @brief Compute the exact solution for the Ringleb flow problem.
1170 */
1172 int field, Array<OneD, NekDouble> &outarray)
1173{
1174 int nTotQuadPoints = GetTotPoints();
1175
1176 Array<OneD, NekDouble> rho(nTotQuadPoints, 100.0);
1177 Array<OneD, NekDouble> rhou(nTotQuadPoints);
1178 Array<OneD, NekDouble> rhov(nTotQuadPoints);
1179 Array<OneD, NekDouble> E(nTotQuadPoints);
1180 Array<OneD, NekDouble> x(nTotQuadPoints);
1181 Array<OneD, NekDouble> y(nTotQuadPoints);
1182 Array<OneD, NekDouble> z(nTotQuadPoints);
1183
1184 m_fields[0]->GetCoords(x, y, z);
1185
1186 // Flow parameters
1187 NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
1188 NekDouble J11, J12, J21, J22, det;
1189 NekDouble Fx, Fy;
1190 NekDouble xi, yi;
1191 NekDouble dV;
1192 NekDouble dtheta;
1193 NekDouble par1;
1194 NekDouble theta = M_PI / 4.0;
1195 NekDouble kExt = 0.7;
1196 NekDouble V = kExt * sin(theta);
1197 NekDouble toll = 1.0e-8;
1198 NekDouble errV = 1.0;
1199 NekDouble errTheta = 1.0;
1200 NekDouble gamma = m_gamma;
1201 NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
1202
1203 for (int i = 0; i < nTotQuadPoints; ++i)
1204 {
1205 while ((abs(errV) > toll) || (abs(errTheta) > toll))
1206 {
1207 VV = V * V;
1208 sint = sin(theta);
1209 c = sqrt(1.0 - gamma_1_2 * VV);
1210 k = V / sint;
1211 phi = 1.0 / k;
1212 pp = phi * phi;
1213 J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
1214 1.0 / (5.0 * c * c * c * c * c) -
1215 0.5 * log((1.0 + c) / (1.0 - c));
1216
1217 r = pow(c, 1.0 / gamma_1_2);
1218 xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
1219 yi = phi / (r * V) * sqrt(1.0 - VV * pp);
1220 par1 = 25.0 - 5.0 * VV;
1221 ss = sint * sint;
1222
1223 Fx = xi - x[i];
1224 Fy = yi - y[i];
1225
1226 J11 =
1227 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) * V +
1228 1562.5 / pow(par1, 2.5) *
1229 (-2.0 / (VV * V) + 4.0 / (VV * V) * ss) +
1230 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
1231 7812.5 / pow(par1, 3.5) * V -
1232 0.25 *
1233 (-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 * pow(par1, 0.5)) -
1234 (1.0 + 0.2 * pow(par1, 0.5)) /
1235 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
1236 pow(par1, 0.5) * V) /
1237 (1.0 + 0.2 * pow(par1, 0.5)) * (1.0 - 0.2 * pow(par1, 0.5));
1238
1239 J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
1240 J21 = -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
1241 pow((1.0 - ss), 0.5) +
1242 78125.0 / V * sint / pow(par1, 3.5) * pow((1.0 - ss), 0.5);
1243
1244 // the matrix is singular when theta = pi/2
1245 if (abs(y[i]) < toll && abs(cos(theta)) < toll)
1246 {
1247 J22 = -39062.5 / pow(par1, 3.5) / V +
1248 3125 / pow(par1, 2.5) / (VV * V) +
1249 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
1250 7812.5 / pow(par1, 3.5) * V -
1251 0.25 *
1252 (-1.0 / pow(par1, 0.5) * V /
1253 (1.0 - 0.2 * pow(par1, 0.5)) -
1254 (1.0 + 0.2 * pow(par1, 0.5)) /
1255 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
1256 pow(par1, 0.5) * V) /
1257 (1.0 + 0.2 * pow(par1, 0.5)) *
1258 (1.0 - 0.2 * pow(par1, 0.5));
1259
1260 // dV = -dV/dx * Fx
1261 dV = -1.0 / J22 * Fx;
1262 dtheta = 0.0;
1263 theta = M_PI / 2.0;
1264 }
1265 else
1266 {
1267 J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
1268 pow((1.0 - ss), 0.5) -
1269 3125.0 / VV * ss / pow(par1, 2.5) / pow((1.0 - ss), 0.5) *
1270 cos(theta);
1271
1272 det = -1.0 / (J11 * J22 - J12 * J21);
1273
1274 // [dV dtheta]' = -[invJ]*[Fx Fy]'
1275 dV = det * (J22 * Fx - J12 * Fy);
1276 dtheta = det * (-J21 * Fx + J11 * Fy);
1277 }
1278
1279 V = V + dV;
1280 theta = theta + dtheta;
1281
1282 errV = abs(dV);
1283 errTheta = abs(dtheta);
1284 }
1285
1286 c = sqrt(1.0 - gamma_1_2 * VV);
1287 r = pow(c, 1.0 / gamma_1_2);
1288
1289 rho[i] = r;
1290 rhou[i] = rho[i] * V * cos(theta);
1291 rhov[i] = rho[i] * V * sin(theta);
1292 P = (c * c) * rho[i] / gamma;
1293 E[i] = P / (gamma - 1.0) +
1294 0.5 * (rhou[i] * rhou[i] / rho[i] + rhov[i] * rhov[i] / rho[i]);
1295
1296 // Resetting the guess value
1297 errV = 1.0;
1298 errTheta = 1.0;
1299 theta = M_PI / 4.0;
1300 V = kExt * sin(theta);
1301 }
1302
1303 switch (field)
1304 {
1305 case 0:
1306 outarray = rho;
1307 break;
1308 case 1:
1309 outarray = rhou;
1310 break;
1311 case 2:
1312 outarray = rhov;
1313 break;
1314 case 3:
1315 outarray = E;
1316 break;
1317 default:
1318 ASSERTL0(false, "Error in variable number!");
1319 break;
1320 }
1321}
1322
1324 int spaceDim, Array<OneD, MultiRegions::ExpListSharedPtr> &fields)
1325{
1326 m_spaceDim = spaceDim;
1327 m_fieldsALE = fields;
1328
1329 // Initialise grid velocities as 0s
1332 for (int i = 0; i < spaceDim; ++i)
1333 {
1334 m_gridVelocity[i] =
1335 Array<OneD, NekDouble>(fields[0]->GetTotPoints(), 0.0);
1337 Array<OneD, NekDouble>(fields[0]->GetTrace()->GetTotPoints(), 0.0);
1338 }
1339
1340 ALEHelper::InitObject(spaceDim, fields);
1341}
1342
1343} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations.
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
void GetElmtTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &tstep)
Calculate the maximum timestep on each element subject to CFL restrictions.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
virtual bool v_SupportsShockCaptType(const std::string type) const =0
MultiRegions::ExpListSharedPtr v_GetPressure() override
Array< OneD, NekDouble > GetStabilityLimitVector(const Array< OneD, int > &ExpOrder)
Function to calculate the stability limit for DG/CG (a vector of them).
void InitAdvection()
Create advection and diffusion objects for CFS.
NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray) override
Calculate the maximum timestep subject to CFL restrictions.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
virtual void v_DoDiffusion(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)=0
void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
Set up logic for residual calculation.
void GetExactRinglebFlow(int field, Array< OneD, NekDouble > &outarray)
Ringleb Flow Test Case.
void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
void v_GetDensity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
void EvaluateIsentropicVortex(unsigned int field, Array< OneD, NekDouble > &outfield, NekDouble time, const int o=0)
Isentropic Vortex Test Case.
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
void GetFluxVectorDeAlias(const Array< OneD, const Array< OneD, NekDouble > > &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations by using the de-aliasing technique.
std::vector< CFSBndCondSharedPtr > m_bndConds
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
VariableConverterSharedPtr m_varConv
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
void SetBoundaryConditionsBwdWeight()
Set up a weight on physical boundaries for boundary condition applications.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
ArtificialDiffusionSharedPtr m_artificialDiffusion
Array< OneD, NekDouble > v_GetMaxStdVelocity(const NekDouble SpeedSoundFactor) override
Compute the advection velocity in the standard space for each element of the expansion.
void DoAdvection(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Compute the advection terms for the right-hand side.
void InitialiseParameters()
Load CFS parameters from the session file.
void DoDiffusion(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Add the diffusions terms to the right-hand side.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2) override
void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time=0.0) override
void v_GetVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsALE
Definition: ALEHelper.h:88
SOLVER_UTILS_EXPORT void InitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Definition: ALEHelper.cpp:48
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Definition: ALEHelper.h:90
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ALEHelper.cpp:148
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
Definition: ALEHelper.cpp:389
SOLVER_UTILS_EXPORT void MoveMesh(const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
Definition: ALEHelper.cpp:168
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace()
Definition: ALEHelper.cpp:289
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
Definition: ALEHelper.h:89
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > GetElmtCFLVals(const bool FlagAcousticCFL=true)
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
int m_spacedim
Spatial dimension (>= expansion dim).
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT int GetExpSize()
std::string m_sessionName
Name of the session.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
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.
enum HomogeneousType m_HomogeneousType
SOLVER_UTILS_EXPORT int GetTotPoints()
int m_nchk
Number of checkpoints written so far.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
int m_checksteps
Number of steps between checkpoints.
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:76
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
@ P
Monomial polynomials .
Definition: BasisType.h:62
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:40
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
ArtificialDiffusionFactory & GetArtificialDiffusionFactory()
Declaration of the artificial diffusion factory singleton.
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.hpp:72
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.hpp:396
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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.hpp:725
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.hpp:366
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:154
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:289
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:294
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285