41 string NavierStokesCFE::className =
43 "NavierStokesCFE", NavierStokesCFE::create,
44 "NavierStokes equations in conservative variables.");
46 NavierStokesCFE::NavierStokesCFE(
73 m_session->LoadParameter(
"GasConstant", gasConstant, 287.058);
88 if (
m_session->DefinesParameter(
"thermalConductivity"))
91 "Cannot define both Pr and thermalConductivity.");
93 m_session->LoadParameter(
"thermalConductivity",
109 m_session->LoadSolverInfo(
"DiffusionType", diffName,
"LDGNS");
113 if (
"InteriorPenalty" == diffName)
119 if (
"LDGNS" == diffName ||
"LDGNS3DHomogeneous1D" == diffName)
136 &NavierStokesCFE::GetViscousFluxVectorConservVar<false>,
this);
139 &NavierStokesCFE::GetViscousFluxVectorConservVar<true>,
this);
156 size_t nvariables = inarray.size();
162 for (
size_t i = 0; i < nvariables; ++i)
192 for (
size_t i = 0; i < nvariables; ++i)
194 Vmath::Vadd(npoints, outarrayDiff[i], 1, outarray[i], 1,
205 for (
size_t i = 0; i < nvariables - 1; ++i)
214 m_varConv->GetPressure(inarray, inarrayDiff[0]);
217 m_varConv->GetTemperature(inarray, inarrayDiff[nvariables - 2]);
220 m_varConv->GetVelocityVector(inarray, inarrayDiff);
234 m_varConv->GetTemperature(pFwd, inFwd[nvariables - 2]);
235 m_varConv->GetTemperature(pBwd, inBwd[nvariables - 2]);
237 m_varConv->GetVelocityVector(pFwd, inFwd);
238 m_varConv->GetVelocityVector(pBwd, inBwd);
245 for (
size_t i = 0; i < nvariables; ++i)
247 Vmath::Vadd(npoints, outarrayDiff[i], 1, outarray[i], 1,
269 size_t nScalar = physfield.size();
270 size_t nPts = physfield[0].size();
280 thermalConductivity);
285 Vmath::Vadd(nPts, divVel, 1, derivativesO1[j][j], 1, divVel, 1);
303 Vmath::Vadd(nPts, derivativesO1[i][j], 1, derivativesO1[j][i], 1,
304 viscousTensor[i][j + 1], 1);
306 Vmath::Vmul(nPts, mu, 1, viscousTensor[i][j + 1], 1,
307 viscousTensor[i][j + 1], 1);
312 Vmath::Vadd(nPts, viscousTensor[i][j + 1], 1, divVel, 1,
313 viscousTensor[i][j + 1], 1);
319 viscousTensor[j][i + 1], 1);
331 Vmath::Vvtvp(nPts, physfield[j], 1, viscousTensor[i][j + 1], 1,
354 size_t nScalar = physfield.size();
355 size_t nPts =
m_fields[0]->Get1DScaledTotPoints(OneDptscale);
356 size_t nPts_orig = physfield[0].size();
368 thermalConductivity);
378 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
386 m_fields[0]->PhysInterp1DScaled(OneDptscale, derivativesO1[i][j],
401 Vmath::Vadd(nPts, divVel, 1, deriv_interp[j][j], 1, divVel, 1);
419 Vmath::Vadd(nPts, deriv_interp[i][j], 1, deriv_interp[j][i], 1,
420 out_interp[i][j + 1], 1);
423 out_interp[i][j + 1], 1);
428 Vmath::Vadd(nPts, out_interp[i][j + 1], 1, divVel, 1,
429 out_interp[i][j + 1], 1);
434 out_interp[j][i + 1] = out_interp[i][j + 1];
446 Vmath::Vvtvp(nPts, vel_interp[j], 1, out_interp[i][j + 1], 1,
461 m_fields[0]->PhysGalerkinProjection1DScaled(
462 OneDptscale, out_interp[i][j], viscousTensor[i][j]);
479 size_t nConvectiveFields = consvar.size();
481 size_t nengy = nConvectiveFields - 1;
488 size_t nLengthArray = 0;
491 size_t nBndRegions =
m_fields[nengy]->GetBndCondExpansions().size();
492 for (
size_t j = 0; j < nBndRegions; ++j)
495 ->GetBndConditions()[j]
502 m_fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
503 for (
size_t e = 0; e < nBndEdges; ++e)
505 size_t nBndEdgePts =
m_fields[nengy]
506 ->GetBndCondExpansions()[j]
510 int id2 =
m_fields[0]->GetTrace()->GetPhys_Offset(
511 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
515 m_fields[nengy]->GetBndConditions()[j]->GetUserDefined(),
518 if (nBndEdgePts != nLengthArray)
520 for (
size_t i = 0; i < nConvectiveFields; ++i)
528 nLengthArray = nBndEdgePts;
539 for (
size_t k = 0; k < nConvectiveFields; ++k)
545 m_varConv->GetPressure(bndCons, bndPressure);
547 m_varConv->GetRhoFromPT(bndPressure, bndTotEngy, bndRho);
548 m_varConv->GetEFromRhoP(bndRho, bndPressure, bndIntEndy);
549 m_varConv->GetDynamicEnergy(bndCons, bndTotEngy);
551 Vmath::Vvtvp(nBndEdgePts, bndIntEndy, 1, bndCons[ndens], 1,
552 bndTotEngy, 1, bndTotEngy, 1);
555 tmp = consvar[nengy] + id2, 1);
571 size_t nConvectiveFields = inarray.size();
572 size_t nPts = inaverg[nConvectiveFields - 1].size();
574 for (
size_t i = 0; i < nConvectiveFields - 1; ++i)
576 nonZeroIndex[i] = i + 1;
582 m_varConv->GetTemperature(inaverg, temperature);
585 std::vector<NekDouble> inAvgTmp(nConvectiveFields);
586 std::vector<NekDouble> inTmp(nConvectiveFields);
587 std::vector<NekDouble> outTmp(nConvectiveFields);
588 for (
size_t d = 0; d < nDim; ++d)
590 for (
size_t nderiv = 0; nderiv < nDim; ++nderiv)
592 for (
size_t p = 0;
p < nPts; ++
p)
595 for (
size_t f = 0; f < nConvectiveFields; ++f)
597 inAvgTmp[f] = inaverg[f][
p];
598 inTmp[f] = inarray[f][
p];
602 inAvgTmp.data(), inTmp.data(),
603 mu[
p], outTmp.data());
605 for (
size_t f = 0; f < nConvectiveFields; ++f)
607 outarray[d][f][
p] += normal[nderiv][
p] * outTmp[f];
622 size_t nTracePts = uFwd[0].size();
625 size_t nVariables = uFwd.size();
628 uBwd[nVariables - 1], 1, tAve, 1);
637 for (
size_t i = 0; i < nVariables; ++i)
640 Vmath::Vsub(nTracePts, uFwd[i], 1, uBwd[i], 1, penaltyCoeff[i], 1);
642 if (i < nVariables - 1)
644 Vmath::Vmul(nTracePts, muAve, 1, penaltyCoeff[i], 1,
649 Vmath::Vmul(nTracePts, tcAve, 1, penaltyCoeff[i], 1,
663 size_t nPts = temperature.size();
665 for (
size_t p = 0;
p < nPts; ++
p)
675 size_t nTracePts =
m_fields[0]->GetTrace()->GetTotPoints();
676 if (nPts != nTracePts)
711 size_t nDim = fields[0]->GetCoordim(0);
712 size_t nVar = cnsVar.size();
713 size_t nPts = cnsVar[0].size();
714 size_t nPtsTrc = cnsVarFwd[0].size();
718 primVarBwd(nVar - 1);
720 for (
unsigned short d = 0; d < nVar - 2; ++d)
726 size_t ergLoc = nVar - 2;
734 for (
unsigned short d = 0; d < nVar - 2; ++d)
737 for (
size_t p = 0;
p < nPts; ++
p)
739 primVar[d][
p] = cnsVar[d + 1][
p] / cnsVar[0][
p];
742 for (
size_t p = 0;
p < nPtsTrc; ++
p)
744 primVarFwd[d][
p] = cnsVarFwd[d + 1][
p] / cnsVarFwd[0][
p];
745 primVarBwd[d][
p] = cnsVarBwd[d + 1][
p] / cnsVarBwd[0][
p];
751 for (
unsigned short j = 0; j < nDim; ++j)
754 for (
unsigned short i = 0; i < nVar - 1; ++i)
761 m_diffusion->DiffuseCalcDerivative(fields, primVar, primVarDer, primVarFwd,
776 size_t nDim = pVarDer.size();
777 size_t nPts = div.size();
780 for (
size_t p = 0;
p < nPts; ++
p)
783 for (
unsigned short j = 0; j < nDim; ++j)
785 divTmp += pVarDer[j][j][
p];
793 for (
size_t p = 0;
p < nPts; ++
p)
817 curl0sqr += curl1sqr + curl2sqr;
819 curlSquare[
p] = curl0sqr;
824 for (
size_t p = 0;
p < nPts; ++
p)
833 curlSquare[
p] = curl * curl;
844 std::vector<std::string> &variables)
847 m_session->MatchSolverInfo(
"OutputExtraFields",
"True", extraFields,
true);
850 const int nPhys =
m_fields[0]->GetNpoints();
851 const int nCoeffs =
m_fields[0]->GetNcoeffs();
854 for (
size_t i = 0; i <
m_fields.size(); ++i)
872 m_varConv->GetVelocityVector(cnsVar, velocity);
874 m_varConv->GetTemperature(cnsVar, temperature);
876 m_varConv->GetSoundSpeed(cnsVar, soundspeed);
877 m_varConv->GetMach(cnsVar, soundspeed, mach);
880 m_session->LoadParameter(
"SensorOffset", sensorOffset, 1);
889 string velNames[3] = {
"u",
"v",
"w"};
892 m_fields[0]->FwdTransLocalElmt(velocity[i], velFwd[i]);
893 variables.push_back(velNames[i]);
894 fieldcoeffs.push_back(velFwd[i]);
898 m_fields[0]->FwdTransLocalElmt(temperature, TFwd);
899 m_fields[0]->FwdTransLocalElmt(entropy, sFwd);
900 m_fields[0]->FwdTransLocalElmt(soundspeed, aFwd);
901 m_fields[0]->FwdTransLocalElmt(mach, mFwd);
902 m_fields[0]->FwdTransLocalElmt(sensor, sensFwd);
904 variables.push_back(
"p");
905 variables.push_back(
"T");
906 variables.push_back(
"s");
907 variables.push_back(
"a");
908 variables.push_back(
"Mach");
909 variables.push_back(
"Sensor");
910 fieldcoeffs.push_back(pFwd);
911 fieldcoeffs.push_back(TFwd);
912 fieldcoeffs.push_back(sFwd);
913 fieldcoeffs.push_back(aFwd);
914 fieldcoeffs.push_back(mFwd);
915 fieldcoeffs.push_back(sensFwd);
924 variables.push_back(
"ArtificialVisc");
925 fieldcoeffs.push_back(sensorFwd);
934 for (
size_t i = 0; i <
m_fields.size(); ++i)
938 m_fields[i]->GetFwdBwdTracePhys(cnsVar[i], cnsVarFwd[i],
947 m_fields[0]->FwdTransLocalElmt(div, divFwd);
948 variables.push_back(
"div");
949 fieldcoeffs.push_back(divFwd);
952 m_fields[0]->FwdTransLocalElmt(curlSquare, curlFwd);
953 variables.push_back(
"curl^2");
954 fieldcoeffs.push_back(curlFwd);
960 variables.push_back(
"ArtificialVisc");
961 fieldcoeffs.push_back(muavFwd);
968 if (type ==
"NonSmooth" || type ==
"Physical" || type ==
"Off")
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
NekDouble m_bndEvaluateTime
std::string m_shockCaptureType
SolverUtils::DiffusionSharedPtr m_diffusion
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
VariableConverterSharedPtr m_varConv
void SetBoundaryConditionsBwdWeight()
Set up a weight on physical boundaries for boundary condition applications.
ArtificialDiffusionSharedPtr m_artificialDiffusion
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
virtual ~NavierStokesCFE()
virtual void v_GetViscousFluxVector(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
Return the flux vector for the LDG diffusion problem.
void GetDivCurlFromDvelT(const TensorOfArray3D< NekDouble > &pVarDer, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare)
Get divergence and curl from velocity derivative tensor.
virtual void v_InitObject(bool DeclareField=true) override
Initialization object for CompressibleFlowSystem class.
void SpecialBndTreat(Array< OneD, Array< OneD, NekDouble >> &consvar)
For very special treatment. For general boundaries it does nothing But for WallViscous and WallAdiaba...
void GetDivCurlSquared(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &cnsVar, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare, const Array< OneD, Array< OneD, NekDouble >> &cnsVarFwd, const Array< OneD, Array< OneD, NekDouble >> &cnsVarBwd)
Get divergence and curl squared.
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) override
void GetViscousSymmtrFluxConservVar(const size_t nSpaceDim, const Array< OneD, Array< OneD, NekDouble >> &inaverg, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &normals)
Calculate and return the Symmetric flux in IP method.
void InitObject_Explicit()
NekDouble m_thermalConductivityRef
bool m_is_mu_variable
flag to switch between constant viscosity and Sutherland an enum could be added for more options
bool m_is_shockCaptPhys
flag for shock capturing switch on/off an enum could be added for more options
virtual void v_GetViscousFluxVectorDeAlias(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
Return the flux vector for the LDG diffusion problem.
virtual bool v_SupportsShockCaptType(const std::string type) const override
void GetViscosityAndThermalCondFromTemp(const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &thermalCond)
Update viscosity todo: add artificial viscosity here.
bool m_is_diffIP
flag to switch between IP and LDG an enum could be added for more options
virtual void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables) override
void GetViscousFluxBilinearFormKernel(const unsigned short nDim, const unsigned short FluxDirection, const unsigned short DerivDirection, const T *inaverg, const T *injumpp, const T &mu, T *outarray)
Calculate diffusion flux using the Jacobian form.
virtual void v_GetFluxPenalty(const Array< OneD, const Array< OneD, NekDouble >> &uFwd, const Array< OneD, const Array< OneD, NekDouble >> &uBwd, Array< OneD, Array< OneD, NekDouble >> &penaltyCoeff)
Return the penalty vector for the LDGNS diffusion problem.
void GetViscosityAndThermalCondFromTempKernel(const T &temperature, T &mu, T &thermalCond)
std::string m_ViscosityType
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNpoints()
Base class for unsteady solvers.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DiffusionFactory & GetDiffusionFactory()
EquationSystemFactory & GetEquationSystemFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (scalar times vector plus scalar times vector):
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.
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
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.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.