Nektar++
VariableConverter.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: VariableConverter.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: Auxiliary functions to convert variables in
32// the compressible flow system
33//
34///////////////////////////////////////////////////////////////////////////////
35#include <iomanip>
36#include <iostream>
37
42
43using namespace std;
44
45namespace Nektar
46{
48 const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim,
50 : m_session(pSession), m_spacedim(spaceDim)
51{
52 // Create equation of state object
53 std::string eosType;
54 m_session->LoadSolverInfo("EquationOfState", eosType, "IdealGas");
56
57 // Parameters for dynamic viscosity
58 m_session->LoadParameter("pInf", m_pInf, 101325);
59 m_session->LoadParameter("rhoInf", m_rhoInf, 1.225);
60 m_session->LoadParameter("GasConstant", m_gasConstant, 287.058);
61 m_session->LoadParameter("mu", m_mu, 1.78e-05);
63
64 // Parameters for sensor
65 m_session->LoadParameter("Skappa", m_Skappa, -1.0);
66 m_session->LoadParameter("Kappa", m_Kappa, 0.25);
67
69
70 // Shock sensor
71 m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
72 if (m_shockCaptureType == "Physical")
73 {
74 // Artificial viscosity scaling constant
75 m_session->LoadParameter("mu0", m_mu0, 1.0);
76
79
80 // Check for Modal/Dilatation sensor
81 m_session->LoadSolverInfo("ShockSensorType", m_shockSensorType,
82 "Dilatation");
83
84 // Check for Ducros sensor
85 m_session->LoadSolverInfo("DucrosSensor", m_ducrosSensor, "Off");
86
87 if (m_ducrosSensor != "Off" || m_shockSensorType == "Dilatation")
88 {
89 m_flagCalcDivCurl = true;
90 }
91 }
92 // Load smoothing type.
93 // We only allocate the smoother if a mesh-graph was passed to the
94 // constructor. This is done to prevent the smoother from being allocated
95 // in cases when it won't be used. Otherwise, any class in the code that
96 // allocates a VariableConverter object try to will allocate a smoother,
97 // even if no mesh graph was passed to the constructor.
98 // TODO: the smoother should be separated from the VariableConverter class.
99 m_session->LoadSolverInfo("Smoothing", m_smoothing, "Off");
100 if (m_smoothing == "C0" && pGraph.get() != nullptr)
101 {
104 m_session, pGraph, m_session->GetVariable(0));
105 }
106
107 std::string viscosityType;
108 m_session->LoadSolverInfo("ViscosityType", viscosityType, "Constant");
109 if ("Variable" == viscosityType)
110 {
111 WARNINGL0(
112 m_session->DefinesParameter("Tref"),
113 "The Tref should be given in Kelvin for using the Sutherland's law "
114 "of dynamic viscosity. The default is 288.15. Note the mu or "
115 "Reynolds number should coorespond to this temperature.");
116 m_session->LoadParameter("Tref", m_Tref, 288.15);
117 m_TRatioSutherland = 110.0 / m_Tref;
118 }
119}
120
121/**
122 * @brief Destructor for VariableConverter class.
123 */
125{
126}
127
128/**
129 * @brief Compute the dynamic energy
130 * \f$ e = rho*V^2/2 \f$.
131 */
133 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
135{
136 size_t nPts = physfield[m_spacedim + 1].size();
137 Vmath::Zero(nPts, energy, 1);
138
139 // tmp = (rho * u_i)^2
140 for (size_t i = 0; i < m_spacedim; ++i)
141 {
142 Vmath::Vvtvp(nPts, physfield[i + 1], 1, physfield[i + 1], 1, energy, 1,
143 energy, 1);
144 }
145 // Divide by rho and multiply by 0.5 --> tmp = 0.5 * rho * u^2
146 Vmath::Vdiv(nPts, energy, 1, physfield[0], 1, energy, 1);
147 Vmath::Smul(nPts, 0.5, energy, 1, energy, 1);
148}
149
150/**
151 * @brief Compute the specific internal energy
152 * \f$ e = (E - rho*V^2/2)/rho \f$.
153 */
155 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
157{
158 size_t nPts = physfield[0].size();
159 Array<OneD, NekDouble> tmp(nPts);
160
161 GetDynamicEnergy(physfield, tmp);
162
163 // Calculate rhoe = E - rho*V^2/2
164 Vmath::Vsub(nPts, physfield[m_spacedim + 1], 1, tmp, 1, energy, 1);
165 // Divide by rho
166 Vmath::Vdiv(nPts, energy, 1, physfield[0], 1, energy, 1);
167}
168
169/**
170 * @brief Compute the specific enthalpy \f$ h = e + p/rho \f$.
171 */
173 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
174 Array<OneD, NekDouble> &enthalpy)
175{
176 size_t nPts = physfield[0].size();
177 Array<OneD, NekDouble> energy(nPts, 0.0);
179
180 GetInternalEnergy(physfield, energy);
181 GetPressure(physfield, pressure);
182
183 // Calculate p/rho
184 Vmath::Vdiv(nPts, pressure, 1, physfield[0], 1, enthalpy, 1);
185 // Calculate h = e + p/rho
186 Vmath::Vadd(nPts, energy, 1, enthalpy, 1, enthalpy, 1);
187}
188
189/**
190 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
191 * \f$ \rho\mathbf{v} \f$.
192 *
193 * @param physfield Momentum field.
194 * @param velocity Velocity field.
195 */
197 const Array<OneD, Array<OneD, NekDouble>> &physfield,
199{
200 const size_t nPts = physfield[0].size();
201
202 for (size_t i = 0; i < m_spacedim; ++i)
203 {
204 Vmath::Vdiv(nPts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
205 }
206}
207
208/**
209 * @brief Compute the mach number \f$ M = \| \mathbf{v} \|^2 / c \f$.
210 *
211 * @param physfield Input physical field.
212 * @param soundfield The speed of sound corresponding to physfield.
213 * @param mach The resulting mach number \f$ M \f$.
214 */
216 Array<OneD, NekDouble> &soundspeed,
218{
219 const size_t nPts = physfield[0].size();
220
221 Vmath::Vmul(nPts, physfield[1], 1, physfield[1], 1, mach, 1);
222
223 for (size_t i = 1; i < m_spacedim; ++i)
224 {
225 Vmath::Vvtvp(nPts, physfield[1 + i], 1, physfield[1 + i], 1, mach, 1,
226 mach, 1);
227 }
228
229 Vmath::Vdiv(nPts, mach, 1, physfield[0], 1, mach, 1);
230 Vmath::Vdiv(nPts, mach, 1, physfield[0], 1, mach, 1);
231 Vmath::Vsqrt(nPts, mach, 1, mach, 1);
232
233 Vmath::Vdiv(nPts, mach, 1, soundspeed, 1, mach, 1);
234}
235
236/**
237 * @brief Compute the dynamic viscosity using the Sutherland's law
238 * \f$ \mu = \mu_star * (T / T_star)^3/2 * (1 + C) / (T/T_star + C) \f$,
239 * C : 110. /Tref
240 * Tref : the reference temperature, Tref, should always given in Kelvin,
241 * if non-dimensional should be the reference for non-dimensionalizing
242 * muref : the dynamic viscosity or the 1/Re corresponding to Tref
243 * T_star : m_pInf / (m_rhoInf * m_gasConstant),non-dimensional or dimensional
244 *
245 * WARNING, if this routine is modified the same must be done in the
246 * FieldConvert utility ProcessWSS.cpp (this class should be restructured).
247 *
248 * @param temperature Input temperature.
249 * @param mu The resulting dynamic viscosity.
250 */
253{
254 const size_t nPts = temperature.size();
255
256 for (size_t i = 0; i < nPts; ++i)
257 {
258 mu[i] = GetDynamicViscosity(temperature[i]);
259 }
260}
261
262/**
263 * @brief Compute the dynamic viscosity using the Sutherland's law
264 * \f$ \mu = \mu_star * (T / T_star)^3/2 * (T_star + 110) / (T + 110) \f$,
265 */
267 const Array<OneD, const NekDouble> &temperature,
269{
270 const size_t nPts = temperature.size();
271 NekDouble tmp = 0.0;
272 NekDouble ratio;
273
274 for (size_t i = 0; i < nPts; ++i)
275 {
276 ratio = temperature[i] * m_oneOverT_star;
277 tmp = 0.5 * (ratio + 3.0 * m_TRatioSutherland) /
278 (ratio * (ratio + m_TRatioSutherland));
279 DmuDT[i] = mu[i] * tmp * m_oneOverT_star;
280 }
281}
282
284 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
286{
287 const size_t nPts = physfield[0].size();
288
289 // Getting the velocity vector on the 2D normal space
291
292 Vmath::Zero(Vtot.size(), Vtot, 1);
293
294 for (size_t i = 0; i < m_spacedim; ++i)
295 {
296 velocity[i] = Array<OneD, NekDouble>(nPts);
297 }
298
299 GetVelocityVector(physfield, velocity);
300
301 for (size_t i = 0; i < m_spacedim; ++i)
302 {
303 Vmath::Vvtvp(nPts, velocity[i], 1, velocity[i], 1, Vtot, 1, Vtot, 1);
304 }
305
306 Vmath::Vsqrt(nPts, Vtot, 1, Vtot, 1);
307}
308
309/**
310 * @brief Calculate the pressure using the equation of state.
311 *
312 * @param physfield Input momentum.
313 * @param pressure Computed pressure field.
314 */
316 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
318{
319 size_t nPts = physfield[0].size();
320
321 Array<OneD, NekDouble> energy(nPts);
322 GetInternalEnergy(physfield, energy);
323
324 for (size_t i = 0; i < nPts; ++i)
325 {
326 pressure[i] = m_eos->GetPressure(physfield[0][i], energy[i]);
327 }
328}
329
330/**
331 * @brief Compute the temperature using the equation of state.
332 *
333 * @param physfield Input physical field.
334 * @param temperature The resulting temperature \f$ T \f$.
335 */
337 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
338 Array<OneD, NekDouble> &temperature)
339{
340 size_t nPts = physfield[0].size();
341
342 Array<OneD, NekDouble> energy(nPts);
343 GetInternalEnergy(physfield, energy);
344
345 for (size_t i = 0; i < nPts; ++i)
346 {
347 temperature[i] = m_eos->GetTemperature(physfield[0][i], energy[i]);
348 }
349}
350
351/**
352 * @brief Compute the sound speed using the equation of state.
353 *
354 * @param physfield Input physical field
355 * @param soundspeed The resulting sound speed \f$ c \f$.
356 */
358 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
359 Array<OneD, NekDouble> &soundspeed)
360{
361 size_t nPts = physfield[0].size();
362
363 Array<OneD, NekDouble> energy(nPts);
364 GetInternalEnergy(physfield, energy);
365
366 for (size_t i = 0; i < nPts; ++i)
367 {
368 soundspeed[i] = m_eos->GetSoundSpeed(physfield[0][i], energy[i]);
369 }
370}
371
372/**
373 * @brief Compute the entropy using the equation of state.
374 *
375 * @param physfield Input physical field
376 * @param soundspeed The resulting sound speed \f$ c \f$.
377 */
379 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
380 Array<OneD, NekDouble> &entropy)
381{
382 size_t nPts = physfield[0].size();
383
384 Array<OneD, NekDouble> energy(nPts);
385 GetInternalEnergy(physfield, energy);
386
387 for (size_t i = 0; i < nPts; ++i)
388 {
389 entropy[i] = m_eos->GetEntropy(physfield[0][i], energy[i]);
390 }
391}
392
393/**
394 * @brief Compute \f$ e(rho,p) \f$ using the equation of state.
395 *
396 * @param rho Input density
397 * @param pressure Input pressure
398 * @param energy The resulting internal energy.
399 */
403{
404 size_t nPts = rho.size();
405
406 for (size_t i = 0; i < nPts; ++i)
407 {
408 energy[i] = m_eos->GetEFromRhoP(rho[i], pressure[i]);
409 }
410}
411
412/**
413 * @brief Compute \f$ rho(p,T) \f$ using the equation of state.
414 *
415 * @param pressure Input pressure
416 * @param temperature Input temperature
417 * @param rho The resulting density
418 */
420 const Array<OneD, NekDouble> &temperature,
422{
423 size_t nPts = pressure.size();
424
425 for (size_t i = 0; i < nPts; ++i)
426 {
427 rho[i] = m_eos->GetRhoFromPT(pressure[i], temperature[i]);
428 }
429}
430
433 const Array<OneD, const Array<OneD, NekDouble>> &consVar,
434 const Array<OneD, NekDouble> &div,
435 const Array<OneD, NekDouble> &curlSquared)
436{
437 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
439 {
440 size_t nPts = fields[0]->GetTotPoints();
441 m_muAv = Array<OneD, NekDouble>(nPts, 0.0);
442 m_muAvTrace = Array<OneD, NekDouble>(nTracePts, 0.0);
443 SetElmtMinHP(fields);
444 }
445
446 if (m_shockSensorType == "Modal")
447 {
448 // Get viscosity based on modal sensor
449 GetMuAv(fields, consVar, m_muAv);
450 }
451 else
452 {
453 // Get viscosity based on dilatation sensor
454 GetMuAv(fields, consVar, div, m_muAv);
455 }
456
457 // Apply Ducros sensor
458 if (m_ducrosSensor != "Off")
459 {
460 ApplyDucros(div, curlSquared, m_muAv);
461 }
462
463 // Apply approximate C0 smoothing
464 if (m_smoothing == "C0")
465 {
467 }
468
469 // Set trace AV
470 Array<OneD, NekDouble> muFwd(nTracePts, 0.0), muBwd(nTracePts, 0.0);
471 fields[0]->GetFwdBwdTracePhys(m_muAv, muFwd, muBwd, false, false, false);
472 for (size_t p = 0; p < nTracePts; ++p)
473 {
474 m_muAvTrace[p] = 0.5 * (muFwd[p] + muBwd[p]);
475 }
476}
477
479{
480 ASSERTL1(m_muAv != NullNekDouble1DArray, "m_muAv not set");
481 return m_muAv;
482}
483
485{
486 ASSERTL1(m_muAvTrace != NullNekDouble1DArray, "m_muAvTrace not set");
487 return m_muAvTrace;
488}
489
490/**
491 * @brief Compute an estimate of minimum h/p
492 * for each element of the expansion.
493 */
496{
497 size_t nElements = fields[0]->GetExpSize();
499 {
500 m_hOverP = Array<OneD, NekDouble>(nElements, 1.0);
501 }
502
503 // Determine h/p scaling
504 Array<OneD, int> pOrderElmt = fields[0]->EvalBasisNumModesMaxPerExp();
505 int expdim = fields[0]->GetGraph()->GetMeshDimension();
506 for (size_t e = 0; e < nElements; e++)
507 {
508 NekDouble h = 1.0e+10;
509 switch (expdim)
510 {
511 case 3:
512 {
514 exp3D = fields[0]->GetExp(e)->as<LocalRegions::Expansion3D>();
515 for (int i = 0; i < exp3D->GetNtraces(); ++i)
516 {
517 h = min(
518 h, exp3D->GetGeom3D()->GetEdge(i)->GetVertex(0)->dist(*(
519 exp3D->GetGeom3D()->GetEdge(i)->GetVertex(1))));
520 }
521 break;
522 }
523
524 case 2:
525 {
527 exp2D = fields[0]->GetExp(e)->as<LocalRegions::Expansion2D>();
528 for (int i = 0; i < exp2D->GetNtraces(); ++i)
529 {
530 h = min(
531 h, exp2D->GetGeom2D()->GetEdge(i)->GetVertex(0)->dist(*(
532 exp2D->GetGeom2D()->GetEdge(i)->GetVertex(1))));
533 }
534 break;
535 }
536 case 1:
537 {
539 exp1D = fields[0]->GetExp(e)->as<LocalRegions::Expansion1D>();
540
541 h = min(h, exp1D->GetGeom1D()->GetVertex(0)->dist(
542 *(exp1D->GetGeom1D()->GetVertex(1))));
543
544 break;
545 }
546 default:
547 {
548 ASSERTL0(false, "Dimension out of bound.")
549 }
550 }
551
552 // Store h/p scaling
553 m_hOverP[e] = h / max(pOrderElmt[e] - 1, 1);
554 }
555}
556
558{
559 ASSERTL1(m_hOverP != NullNekDouble1DArray, "m_hOverP not set");
560 return m_hOverP;
561}
562
565 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
567 int offset)
568{
569 NekDouble Skappa;
570 NekDouble order;
572 Array<OneD, int> expOrderElement = field->EvalBasisNumModesMaxPerExp();
573
574 for (int e = 0; e < field->GetExpSize(); e++)
575 {
576 int numModesElement = expOrderElement[e];
577 int nElmtPoints = field->GetExp(e)->GetTotPoints();
578 int physOffset = field->GetPhys_Offset(e);
579 int nElmtCoeffs = field->GetExp(e)->GetNcoeffs();
580 int numCutOff = numModesElement - offset;
581
582 if (numModesElement <= offset)
583 {
584 Vmath::Fill(nElmtPoints, 0.0, tmp = Sensor + physOffset, 1);
585 Vmath::Fill(nElmtPoints, 0.0, tmp = SensorKappa + physOffset, 1);
586 continue;
587 }
588
589 // create vector to save the solution points per element at P = p;
590 Array<OneD, NekDouble> elmtPhys(nElmtPoints,
591 tmp = physarray[0] + physOffset);
592 // Compute coefficients
593 Array<OneD, NekDouble> elmtCoeffs(nElmtCoeffs, 0.0);
594 field->GetExp(e)->FwdTrans(elmtPhys, elmtCoeffs);
595
596 // ReduceOrderCoeffs reduces the polynomial order of the solution
597 // that is represented by the coeffs given as an inarray. This is
598 // done by projecting the higher order solution onto the orthogonal
599 // basis and padding the higher order coefficients with zeros.
600 Array<OneD, NekDouble> reducedElmtCoeffs(nElmtCoeffs, 0.0);
601 field->GetExp(e)->ReduceOrderCoeffs(numCutOff, elmtCoeffs,
602 reducedElmtCoeffs);
603
604 Array<OneD, NekDouble> reducedElmtPhys(nElmtPoints, 0.0);
605 field->GetExp(e)->BwdTrans(reducedElmtCoeffs, reducedElmtPhys);
606
607 NekDouble numerator = 0.0;
608 NekDouble denominator = 0.0;
609
610 // Determining the norm of the numerator of the Sensor
611 Array<OneD, NekDouble> difference(nElmtPoints, 0.0);
612 Vmath::Vsub(nElmtPoints, elmtPhys, 1, reducedElmtPhys, 1, difference,
613 1);
614
615 numerator = Vmath::Dot(nElmtPoints, difference, difference);
616 denominator = Vmath::Dot(nElmtPoints, elmtPhys, elmtPhys);
617
618 NekDouble elmtSensor = sqrt(numerator / denominator);
619 elmtSensor = log10(max(elmtSensor, NekConstants::kNekMachineEpsilon));
620
621 Vmath::Fill(nElmtPoints, elmtSensor, tmp = Sensor + physOffset, 1);
622
623 // Compute reference value for sensor
624 order = max(numModesElement - 1, 1);
625 if (order > 0)
626 {
627 Skappa = m_Skappa - 4.25 * log10(static_cast<NekDouble>(order));
628 }
629 else
630 {
631 Skappa = 0.0;
632 }
633
634 // Compute artificial viscosity
635 NekDouble elmtSensorKappa;
636 if (elmtSensor < (Skappa - m_Kappa))
637 {
638 elmtSensorKappa = 0;
639 }
640 else if (elmtSensor > (Skappa + m_Kappa))
641 {
642 elmtSensorKappa = 1.0;
643 }
644 else
645 {
646 elmtSensorKappa =
647 0.5 * (1 + sin(M_PI * (elmtSensor - Skappa) / (2 * m_Kappa)));
648 }
649 Vmath::Fill(nElmtPoints, elmtSensorKappa,
650 tmp = SensorKappa + physOffset, 1);
651 }
652}
653
654/**
655 * @brief Calculate the physical artificial viscosity based on modal sensor.
656 *
657 * @param consVar Input field.
658 */
661 const Array<OneD, const Array<OneD, NekDouble>> &consVar,
663{
664 size_t nPts = consVar[0].size();
665 // Determine the maximum wavespeed
666 Array<OneD, NekDouble> Lambdas(nPts, 0.0);
667 Array<OneD, NekDouble> soundspeed(nPts, 0.0);
668 Array<OneD, NekDouble> absVelocity(nPts, 0.0);
669 GetSoundSpeed(consVar, soundspeed);
670 GetAbsoluteVelocity(consVar, absVelocity);
671 Vmath::Vadd(nPts, absVelocity, 1, soundspeed, 1, Lambdas, 1);
672
673 // Compute sensor based on rho
674 Array<OneD, NekDouble> Sensor(nPts, 0.0);
675 GetSensor(fields[0], consVar, Sensor, muAv, 1);
676
678 size_t nElmt = fields[0]->GetExpSize();
679 for (size_t e = 0; e < nElmt; ++e)
680 {
681 int physOffset = fields[0]->GetPhys_Offset(e);
682 int nElmtPoints = fields[0]->GetExp(e)->GetTotPoints();
683
684 // Compute the maximum wave speed
685 NekDouble LambdaElmt = 0.0;
686 LambdaElmt = Vmath::Vmax(nElmtPoints, tmp = Lambdas + physOffset, 1);
687
688 // Compute average bounded density
689 NekDouble rhoAve =
690 Vmath::Vsum(nElmtPoints, tmp = consVar[0] + physOffset, 1);
691 rhoAve = rhoAve / nElmtPoints;
692 rhoAve = Smath::Smax(rhoAve, 1.0e-4, 1.0e+4);
693
694 // Scale sensor by coeff, h/p, and density
695 LambdaElmt *= m_mu0 * m_hOverP[e] * rhoAve;
696 Vmath::Smul(nElmtPoints, LambdaElmt, tmp = muAv + physOffset, 1,
697 tmp = muAv + physOffset, 1);
698 }
699}
700
701/**
702 * @brief Calculate the physical artificial viscosity based on dilatation of
703 * velocity vector.
704 *
705 * @param
706 */
709 const Array<OneD, const Array<OneD, NekDouble>> &consVar,
711{
712 size_t nPts = consVar[0].size();
713
714 // Get sound speed
715 // theoretically it should be used the critical sound speed, this
716 // matters for large Mach numbers (above 3.0)
717 Array<OneD, NekDouble> soundSpeed(nPts, 0.0);
718 GetSoundSpeed(consVar, soundSpeed);
719
720 // Get abosolute velocity to compute lambda
721 Array<OneD, NekDouble> absVelocity(nPts, 0.0);
722 GetAbsoluteVelocity(consVar, absVelocity);
723
724 // Loop over elements
725 size_t nElmt = fields[0]->GetExpSize();
726 for (size_t e = 0; e < nElmt; ++e)
727 {
728 int nElmtPoints = fields[0]->GetExp(e)->GetTotPoints();
729 int physOffset = fields[0]->GetPhys_Offset(e);
730 int physEnd = physOffset + nElmtPoints;
731
732 NekDouble hOpTmp = m_hOverP[e];
733
734 // Loop over the points
735 for (int p = physOffset; p < physEnd; ++p)
736 {
737 // Get non-dimensional sensor based on dilatation
738 NekDouble sSpeedTmp = soundSpeed[p];
739 // (only compression waves)
740 NekDouble divTmp = -div[p];
741 divTmp = std::max(divTmp, 0.0);
742 NekDouble sensor = m_mu0 * hOpTmp * divTmp / sSpeedTmp;
743 // Scale to viscosity scale
744 NekDouble rho = consVar[0][p];
745 NekDouble lambda = sSpeedTmp + absVelocity[p];
746 muAv[p] = sensor * rho * lambda * hOpTmp;
747 }
748 }
749}
750
751/**
752 * @brief Apply Ducros (anti-vorticity) sensor averaged over the element.
753 *
754 * @param field Input Field
755 */
757 const Array<OneD, NekDouble> &curlSquare,
759{
760 // machine eps**2
761 NekDouble eps = std::numeric_limits<NekDouble>::epsilon();
762 eps *= eps;
763
764 // loop over points
765 size_t nPts = div.size();
766 for (size_t p = 0; p < nPts; ++p)
767 {
768 NekDouble tmpDiv2 = div[p];
769 tmpDiv2 *= tmpDiv2;
770 NekDouble denDuc = tmpDiv2 + curlSquare[p] + eps;
771 NekDouble Duc = tmpDiv2 / denDuc;
772 // apply
773 muAv[p] *= Duc;
774 }
775}
776
777/**
778 * @brief Make field C0.
779 *
780 * @param field Input Field
781 */
783{
784 // Make sure that the C0 projection operator has been allocated. Note that
785 // the VariableConverter object can be allocated without the C0 smoother.
786 // This is why this check is needed. Ideally, the C0 smoother is separated
787 // from the VariableConverter class.
788 ASSERTL0(m_C0ProjectExp.get() != nullptr,
789 "C0 projection operator not initialized in "
790 "VariableConverter::ApplyC0Smooth()");
791
792 int nCoeffs = m_C0ProjectExp->GetNcoeffs();
793 Array<OneD, NekDouble> muFwd(nCoeffs);
794 Array<OneD, NekDouble> weights(nCoeffs, 1.0);
795 // Assemble global expansion coefficients for viscosity
796 m_C0ProjectExp->FwdTransLocalElmt(field, m_C0ProjectExp->UpdateCoeffs());
797 m_C0ProjectExp->Assemble();
798 Vmath::Vcopy(nCoeffs, m_C0ProjectExp->GetCoeffs(), 1, muFwd, 1);
799 // Global coefficients
800 Vmath::Vcopy(nCoeffs, weights, 1, m_C0ProjectExp->UpdateCoeffs(), 1);
801 // This is the sign vector
802 m_C0ProjectExp->GlobalToLocal();
803 // Get weights
804 m_C0ProjectExp->Assemble();
805 // Divide
806 Vmath::Vdiv(nCoeffs, muFwd, 1, m_C0ProjectExp->GetCoeffs(), 1,
807 m_C0ProjectExp->UpdateCoeffs(), 1);
808 // Get local coefficients
809 m_C0ProjectExp->GlobalToLocal();
810 // Get C0 field
811 m_C0ProjectExp->BwdTrans(m_C0ProjectExp->GetCoeffs(), field);
812}
813
814} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define WARNINGL0(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:242
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
SpatialDomains::Geometry1DSharedPtr GetGeom1D() const
Definition: Expansion1D.h:104
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void ApplyC0Smooth(Array< OneD, NekDouble > &field)
Make field C0.
Array< OneD, NekDouble > & GetAvTrace()
void GetInternalEnergy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &energy)
Compute the specific internal energy .
EquationOfStateSharedPtr m_eos
void GetSensor(const MultiRegions::ExpListSharedPtr &field, const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, NekDouble > &Sensor, Array< OneD, NekDouble > &SensorKappa, int offset=1)
Array< OneD, NekDouble > m_hOverP
h/p scaling
void GetRhoFromPT(const Array< OneD, NekDouble > &pressure, const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &rho)
Compute using the equation of state.
void GetAbsoluteVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &Vtot)
void GetMuAv(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &muAv)
Calculate the physical artificial viscosity based on modal sensor.
void SetAv(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &consVar, const Array< OneD, NekDouble > &div=NullNekDouble1DArray, const Array< OneD, NekDouble > &curlSquared=NullNekDouble1DArray)
Array< OneD, NekDouble > m_muAvTrace
void GetDmuDT(const Array< OneD, const NekDouble > &temperature, const Array< OneD, const NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
Compute the dynamic viscosity using the Sutherland's law ,.
void GetEFromRhoP(const Array< OneD, NekDouble > &rho, const Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &energy)
Compute using the equation of state.
void SetElmtMinHP(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Compute an estimate of minimum h/p for each element of the expansion.
MultiRegions::ContFieldSharedPtr m_C0ProjectExp
void GetPressure(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
Calculate the pressure using the equation of state.
LibUtilities::SessionReaderSharedPtr m_session
NekDouble m_mu0
Shock sensor.
void GetTemperature(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &temperature)
Compute the temperature using the equation of state.
void ApplyDucros(const Array< OneD, NekDouble > &div, const Array< OneD, NekDouble > &curlSquare, Array< OneD, NekDouble > &muAv)
Apply Ducros (anti-vorticity) sensor averaged over the element.
void GetDynamicEnergy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &energy)
Compute the dynamic energy .
void GetSoundSpeed(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &soundspeed)
Compute the sound speed using the equation of state.
VariableConverter(const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim, const SpatialDomains::MeshGraphSharedPtr &pGraph=nullptr)
void GetEntropy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &entropy)
Compute the entropy using the equation of state.
void GetMach(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &soundspeed, Array< OneD, NekDouble > &mach)
Compute the mach number .
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
~VariableConverter()
Destructor for VariableConverter class.
void GetEnthalpy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &enthalpy)
Compute the specific enthalpy .
Array< OneD, NekDouble > & GetElmtMinHP()
Array< OneD, NekDouble > m_muAv
storage
Array< OneD, NekDouble > & GetAv()
void GetDynamicViscosity(const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &mu)
Compute the dynamic viscosity using the Sutherland's law , C : 110. /Tref Tref : the reference temper...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:46
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:50
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:47
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekMachineEpsilon
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
T Smax(const T a, const T b, const T k)
Return the soft max of between two scalars.
Definition: Smath.hpp:53
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
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 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
T Dot(int n, const T *w, const T *x)
dot product
Definition: Vmath.hpp:761
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 Vdiv(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:126
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294