43 "NavierStokes equations in conservative variables.");
68 m_session->LoadParameter(
"GasConstant", gasConstant, 287.058);
83 if (
m_session->DefinesParameter(
"thermalConductivity"))
86 "Cannot define both Pr and thermalConductivity.");
88 m_session->LoadParameter(
"thermalConductivity",
103 std::string diffName;
104 m_session->LoadSolverInfo(
"DiffusionType", diffName,
"LDGNS");
108 if (
"InteriorPenalty" == diffName)
114 if (
"LDGNS" == diffName ||
"LDGNS3DHomogeneous1D" == diffName)
131 &NavierStokesCFE::GetViscousFluxVectorConservVar<false>,
this);
134 &NavierStokesCFE::GetViscousFluxVectorConservVar<true>,
this);
153 size_t nvariables = inarray.size();
163 for (
size_t i = 0; i < nvariables; ++i)
204 for (
size_t i = 0; i < nvariables; ++i)
206 Vmath::Vadd(npointsOut, outarrayDiff[i], 1, outarray[i], 1,
217 for (
size_t i = 0; i < nvariables - 1; ++i)
225 m_varConv->GetTemperature(inarray, inarrayDiff[nvariables - 2]);
228 m_varConv->GetVelocityVector(inarray, inarrayDiff);
239 m_varConv->GetTemperature(pFwd, inFwd[nvariables - 2]);
240 m_varConv->GetTemperature(pBwd, inBwd[nvariables - 2]);
242 m_varConv->GetVelocityVector(pFwd, inFwd);
243 m_varConv->GetVelocityVector(pBwd, inBwd);
250 outarrayDiff, inFwd, inBwd);
255 outarrayDiff, inFwd, inBwd);
258 for (
size_t i = 0; i < nvariables; ++i)
260 Vmath::Vadd(npointsOut, outarrayDiff[i], 1, outarray[i], 1,
282 size_t nScalar = physfield.size();
283 size_t nPts = physfield[0].size();
293 thermalConductivity);
298 Vmath::Vadd(nPts, divVel, 1, derivativesO1[j][j], 1, divVel, 1);
316 Vmath::Vadd(nPts, derivativesO1[i][j], 1, derivativesO1[j][i], 1,
317 viscousTensor[i][j + 1], 1);
319 Vmath::Vmul(nPts, mu, 1, viscousTensor[i][j + 1], 1,
320 viscousTensor[i][j + 1], 1);
325 Vmath::Vadd(nPts, viscousTensor[i][j + 1], 1, divVel, 1,
326 viscousTensor[i][j + 1], 1);
332 viscousTensor[j][i + 1], 1);
344 Vmath::Vvtvp(nPts, physfield[j], 1, viscousTensor[i][j + 1], 1,
367 size_t nScalar = physfield.size();
368 size_t nPts =
m_fields[0]->Get1DScaledTotPoints(OneDptscale);
369 size_t nPts_orig = physfield[0].size();
381 thermalConductivity);
391 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
399 m_fields[0]->PhysInterp1DScaled(OneDptscale, derivativesO1[i][j],
414 Vmath::Vadd(nPts, divVel, 1, deriv_interp[j][j], 1, divVel, 1);
432 Vmath::Vadd(nPts, deriv_interp[i][j], 1, deriv_interp[j][i], 1,
433 out_interp[i][j + 1], 1);
436 out_interp[i][j + 1], 1);
441 Vmath::Vadd(nPts, out_interp[i][j + 1], 1, divVel, 1,
442 out_interp[i][j + 1], 1);
447 out_interp[j][i + 1] = out_interp[i][j + 1];
459 Vmath::Vvtvp(nPts, vel_interp[j], 1, out_interp[i][j + 1], 1,
474 m_fields[0]->PhysGalerkinProjection1DScaled(
475 OneDptscale, out_interp[i][j], viscousTensor[i][j]);
492 size_t nConvectiveFields = consvar.size();
494 size_t nengy = nConvectiveFields - 1;
501 size_t nLengthArray = 0;
504 size_t nBndRegions =
m_fields[nengy]->GetBndCondExpansions().size();
505 for (
size_t j = 0; j < nBndRegions; ++j)
508 ->GetBndConditions()[j]
515 m_fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
516 for (
size_t e = 0; e < nBndEdges; ++e)
518 size_t nBndEdgePts =
m_fields[nengy]
519 ->GetBndCondExpansions()[j]
523 int id2 =
m_fields[0]->GetTrace()->GetPhys_Offset(
524 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
528 m_fields[nengy]->GetBndConditions()[j]->GetUserDefined(),
531 if (nBndEdgePts != nLengthArray)
533 for (
size_t i = 0; i < nConvectiveFields; ++i)
541 nLengthArray = nBndEdgePts;
552 for (
size_t k = 0; k < nConvectiveFields; ++k)
558 m_varConv->GetPressure(bndCons, bndPressure);
560 m_varConv->GetRhoFromPT(bndPressure, bndTotEngy, bndRho);
561 m_varConv->GetEFromRhoP(bndRho, bndPressure, bndIntEndy);
562 m_varConv->GetDynamicEnergy(bndCons, bndTotEngy);
564 Vmath::Vvtvp(nBndEdgePts, bndIntEndy, 1, bndCons[ndens], 1,
565 bndTotEngy, 1, bndTotEngy, 1);
568 tmp = consvar[nengy] + id2, 1);
584 size_t nConvectiveFields = inarray.size();
585 size_t nPts = inaverg[nConvectiveFields - 1].size();
587 for (
size_t i = 0; i < nConvectiveFields - 1; ++i)
589 nonZeroIndex[i] = i + 1;
595 m_varConv->GetTemperature(inaverg, temperature);
598 std::vector<NekDouble> inAvgTmp(nConvectiveFields);
599 std::vector<NekDouble> inTmp(nConvectiveFields);
600 std::vector<NekDouble> outTmp(nConvectiveFields);
601 for (
size_t d = 0; d < nDim; ++d)
603 for (
size_t nderiv = 0; nderiv < nDim; ++nderiv)
605 for (
size_t p = 0; p < nPts; ++p)
608 for (
size_t f = 0; f < nConvectiveFields; ++f)
610 inAvgTmp[f] = inaverg[f][p];
611 inTmp[f] = inarray[f][p];
615 inAvgTmp.data(), inTmp.data(),
616 mu[p], outTmp.data());
618 for (
size_t f = 0; f < nConvectiveFields; ++f)
620 outarray[d][f][p] += normal[nderiv][p] * outTmp[f];
635 size_t nTracePts = uFwd[0].size();
638 size_t nVariables = uFwd.size();
641 uBwd[nVariables - 1], 1, tAve, 1);
650 for (
size_t i = 0; i < nVariables; ++i)
653 Vmath::Vsub(nTracePts, uFwd[i], 1, uBwd[i], 1, penaltyCoeff[i], 1);
655 if (i < nVariables - 1)
657 Vmath::Vmul(nTracePts, muAve, 1, penaltyCoeff[i], 1,
662 Vmath::Vmul(nTracePts, tcAve, 1, penaltyCoeff[i], 1,
676 size_t nPts = temperature.size();
678 for (
size_t p = 0; p < nPts; ++p)
688 size_t nTracePts =
m_fields[0]->GetTrace()->GetTotPoints();
689 if (nPts != nTracePts)
724 size_t nDim = fields[0]->GetCoordim(0);
725 size_t nVar = cnsVar.size();
726 size_t nPts = cnsVar[0].size();
727 size_t nPtsTrc = cnsVarFwd[0].size();
731 primVarBwd(nVar - 1);
733 for (
unsigned short d = 0; d < nVar - 2; ++d)
739 size_t ergLoc = nVar - 2;
747 for (
unsigned short d = 0; d < nVar - 2; ++d)
750 for (
size_t p = 0; p < nPts; ++p)
752 primVar[d][p] = cnsVar[d + 1][p] / cnsVar[0][p];
755 for (
size_t p = 0; p < nPtsTrc; ++p)
757 primVarFwd[d][p] = cnsVarFwd[d + 1][p] / cnsVarFwd[0][p];
758 primVarBwd[d][p] = cnsVarBwd[d + 1][p] / cnsVarBwd[0][p];
764 for (
unsigned short j = 0; j < nDim; ++j)
767 for (
unsigned short i = 0; i < nVar - 1; ++i)
774 m_diffusion->DiffuseCalcDerivative(fields, primVar, primVarDer, primVarFwd,
789 size_t nDim = pVarDer.size();
790 size_t nPts = div.size();
793 for (
size_t p = 0; p < nPts; ++p)
796 for (
unsigned short j = 0; j < nDim; ++j)
798 divTmp += pVarDer[j][j][p];
806 for (
size_t p = 0; p < nPts; ++p)
830 curl0sqr += curl1sqr + curl2sqr;
832 curlSquare[p] = curl0sqr;
837 for (
size_t p = 0; p < nPts; ++p)
846 curlSquare[p] = curl * curl;
857 std::vector<std::string> &variables)
860 m_session->MatchSolverInfo(
"OutputExtraFields",
"True", extraFields,
true);
863 const int nPhys =
m_fields[0]->GetNpoints();
864 const int nCoeffs =
m_fields[0]->GetNcoeffs();
867 for (
size_t i = 0; i <
m_fields.size(); ++i)
885 m_varConv->GetVelocityVector(cnsVar, velocity);
886 m_varConv->GetPressure(cnsVar, pressure);
887 m_varConv->GetTemperature(cnsVar, temperature);
889 m_varConv->GetSoundSpeed(cnsVar, soundspeed);
890 m_varConv->GetMach(cnsVar, soundspeed, mach);
893 m_session->LoadParameter(
"SensorOffset", sensorOffset, 1);
902 std::string velNames[3] = {
"u",
"v",
"w"};
905 m_fields[0]->FwdTransLocalElmt(velocity[i], velFwd[i]);
906 variables.push_back(velNames[i]);
907 fieldcoeffs.push_back(velFwd[i]);
910 m_fields[0]->FwdTransLocalElmt(pressure, pFwd);
911 m_fields[0]->FwdTransLocalElmt(temperature, TFwd);
912 m_fields[0]->FwdTransLocalElmt(entropy, sFwd);
913 m_fields[0]->FwdTransLocalElmt(soundspeed, aFwd);
914 m_fields[0]->FwdTransLocalElmt(mach, mFwd);
915 m_fields[0]->FwdTransLocalElmt(sensor, sensFwd);
917 variables.push_back(
"p");
918 variables.push_back(
"T");
919 variables.push_back(
"s");
920 variables.push_back(
"a");
921 variables.push_back(
"Mach");
922 variables.push_back(
"Sensor");
923 fieldcoeffs.push_back(pFwd);
924 fieldcoeffs.push_back(TFwd);
925 fieldcoeffs.push_back(sFwd);
926 fieldcoeffs.push_back(aFwd);
927 fieldcoeffs.push_back(mFwd);
928 fieldcoeffs.push_back(sensFwd);
935 m_fields[0]->FwdTransLocalElmt(pressure, sensorFwd);
937 variables.push_back(
"ArtificialVisc");
938 fieldcoeffs.push_back(sensorFwd);
947 for (
size_t i = 0; i <
m_fields.size(); ++i)
951 m_fields[i]->GetFwdBwdTracePhys(cnsVar[i], cnsVarFwd[i],
960 m_fields[0]->FwdTransLocalElmt(div, divFwd);
961 variables.push_back(
"div");
962 fieldcoeffs.push_back(divFwd);
965 m_fields[0]->FwdTransLocalElmt(curlSquare, curlFwd);
966 variables.push_back(
"curl^2");
967 fieldcoeffs.push_back(curlFwd);
973 variables.push_back(
"ArtificialVisc");
974 fieldcoeffs.push_back(muavFwd);
987 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
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.
NavierStokesCFE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void SpecialBndTreat(Array< OneD, Array< OneD, NekDouble > > &consvar)
For very special treatment. For general boundaries it does nothing But for WallViscous and WallAdiaba...
void GetDivCurlFromDvelT(const TensorOfArray3D< NekDouble > &pVarDer, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare)
Get divergence and curl from velocity derivative tensor.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void v_InitObject(bool DeclareField=true) override
Initialization object for CompressibleFlowSystem class.
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.
void InitObject_Explicit()
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
NekDouble m_thermalConductivityRef
bool m_is_mu_variable
flag to switch between constant viscosity and Sutherland an enum could be added for more options
void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
bool m_is_shockCaptPhys
flag for shock capturing switch on/off an enum could be added for more options
bool v_SupportsShockCaptType(const std::string type) const 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 GetViscosityAndThermalCondFromTemp(const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &thermalCond)
Update viscosity todo: add artificial viscosity here.
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.
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.
static std::string className
bool m_is_diffIP
flag to switch between IP and LDG an enum could be added for more options
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.
void GetViscosityAndThermalCondFromTempKernel(const T &temperature, T &mu, T &thermalCond)
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.
std::string m_ViscosityType
SOLVER_UTILS_EXPORT void ExtraFldOutputGrid(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
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 GetTraceTotPoints()
Base class for unsteady solvers.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
DiffusionFactory & GetDiffusionFactory()
EquationSystemFactory & GetEquationSystemFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
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)
Svtsvtp (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.