35 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_NAVIERSTOKESCFE_H
36 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_NAVIERSTOKESCFE_H
108 const size_t nSpaceDim,
125 virtual void v_InitObject(
bool DeclareField =
true)
override;
171 std::vector<std::string> &variables)
override;
173 template <
class T,
typename =
typename std::enable_if<
174 std::is_floating_point<T>::value ||
177 T &mu, T &thermalCond)
181 thermalCond = tRa * mu;
184 template <
class T,
typename =
typename std::enable_if<
185 std::is_floating_point<T>::value ||
192 mu =
m_varConv->GetDynamicViscosity(temperature);
208 template <
class T,
typename =
typename std::enable_if<
209 std::is_floating_point<T>::value ||
212 const unsigned short nDim,
const unsigned short FluxDirection,
213 const unsigned short DerivDirection,
215 const T *inaverg,
const T *injumpp,
const T &mu, T *outarray)
218 unsigned short nDim_plus_one = nDim + 1;
219 unsigned short FluxDirection_plus_one = FluxDirection + 1;
220 unsigned short DerivDirection_plus_one = DerivDirection + 1;
223 NekDouble one_minus_gammaoPr = 1.0 - gammaoPr;
226 constexpr
NekDouble TwoThird = 2. * OneThird;
227 constexpr
NekDouble FourThird = 4. * OneThird;
229 if (DerivDirection == FluxDirection)
235 T oneOrho = 1.0 / inaverg[0];
237 std::array<T, 3> u = {{0.0, 0.0, 0.0}};
238 std::array<T, 3> u2 = {{0.0, 0.0, 0.0}};
240 for (
unsigned short d = 0; d < nDim; ++d)
242 u[d] = inaverg[d + 1] * oneOrho;
248 T E_minus_u2sum = inaverg[nDim_plus_one];
249 E_minus_u2sum *= oneOrho;
250 E_minus_u2sum -= u2sum;
257 T tmp1 = OneThird * u2[FluxDirection] + u2sum;
258 tmp1 += gammaoPr * E_minus_u2sum;
261 T tmp2 = gammaoPr * injumpp[nDim_plus_one] - tmp1;
265 for (
unsigned short d = 0; d < nDim; ++d)
267 unsigned short d_plus_one = d + 1;
269 T outTmpD = injumpp[d_plus_one] - u[d] * injumpp[0];
272 T tmp3 = one_minus_gammaoPr * u[d];
273 outTmpE += tmp3 * injumpp[d_plus_one];
275 if (d == FluxDirection)
277 outTmpD *= FourThird;
278 T tmp4 = OneThird * u[FluxDirection];
279 outTmpE += tmp4 * injumpp[FluxDirection_plus_one];
282 outarray[d_plus_one] = outTmpD;
287 outarray[nDim_plus_one] = outTmpE;
295 T oneOrho = 1.0 / inaverg[0];
297 std::array<T, 3> u, u2;
299 for (
unsigned short d = 0; d < nDim; ++d)
301 unsigned short d_plus_one = d + 1;
302 u[d] = inaverg[d_plus_one] * oneOrho;
307 outarray[d_plus_one] = 0.0;
311 T E_minus_u2sum = inaverg[nDim_plus_one];
312 E_minus_u2sum *= oneOrho;
313 E_minus_u2sum -= u2sum;
320 T tmp1 = u[DerivDirection] * injumpp[0] -
321 injumpp[DerivDirection_plus_one];
323 outarray[FluxDirection_plus_one] = nu * tmp1;
326 injumpp[FluxDirection_plus_one] - u[FluxDirection] * injumpp[0];
327 outarray[DerivDirection_plus_one] = nu * tmp1;
329 tmp1 = OneThird * u[FluxDirection] * u[DerivDirection];
333 TwoThird * u[FluxDirection] * injumpp[DerivDirection_plus_one];
337 tmp1 = u[DerivDirection] * injumpp[FluxDirection_plus_one] - tmp1;
338 outarray[nDim_plus_one] = nu * tmp1;
346 template <
bool IS_TRACE>
353 size_t nConvectiveFields = inarray.size();
354 size_t nPts = inarray[0].size();
355 size_t n_nonZero = nConvectiveFields - 1;
358 constexpr
unsigned short nOutMax = 3 - 2 * IS_TRACE;
359 constexpr
unsigned short nVarMax = 5;
360 constexpr
unsigned short nDimMax = 3;
368 m_varConv->GetTemperature(inarray, temperature);
370 thermalConductivity);
375 size_t sizeVec = (nPts / vec_t::width) * vec_t::width;
378 for (;
p < sizeVec;
p += vec_t::width)
381 alignas(vec_t::alignment) std::array<vec_t, nVarMax> inTmp,
383 alignas(vec_t::alignment) std::array<vec_t, nDimMax> normalTmp;
384 alignas(vec_t::alignment) std::array<vec_t, nVarMax * nOutMax>
388 for (
size_t f = 0; f < nConvectiveFields; ++f)
398 for (
size_t d = 0; d < nDim; ++d)
400 outArrTmp[f + nConvectiveFields * d] = 0.0;
406 for (
size_t d = 0; d < nDim; ++d)
416 for (
size_t nderiv = 0; nderiv < nDim; ++nderiv)
419 for (
size_t f = 0; f < nConvectiveFields; ++f)
421 qfieldsTmp[f].load(&(qfields[nderiv][f][
p]),
425 for (
size_t d = 0; d < nDim; ++d)
428 nDim, d, nderiv, inTmp.data(), qfieldsTmp.data(), muV,
433 for (
size_t f = 0; f < nConvectiveFields; ++f)
435 outArrTmp[f] += normalTmp[d] * outTmp[f];
440 for (
size_t f = 0; f < nConvectiveFields; ++f)
442 outArrTmp[f + nConvectiveFields * d] += outTmp[f];
451 for (
size_t f = 0; f < nConvectiveFields; ++f)
458 for (
size_t d = 0; d < nDim; ++d)
460 for (
size_t f = 0; f < nConvectiveFields; ++f)
462 outArrTmp[f + nConvectiveFields * d].store(
470 for (;
p < nPts; ++
p)
472 std::array<NekDouble, nVarMax> inTmp, qfieldsTmp, outTmp;
473 std::array<NekDouble, nDimMax> normalTmp;
474 std::array<NekDouble, nVarMax * nOutMax> outArrTmp{{}};
476 for (
size_t f = 0; f < nConvectiveFields; ++f)
478 inTmp[f] = inarray[f][
p];
486 for (
size_t d = 0; d < nDim; ++d)
488 outArrTmp[f + nConvectiveFields * d] = 0.0;
495 for (
size_t d = 0; d < nDim; ++d)
497 normalTmp[d] = normal[d][
p];
504 for (
size_t nderiv = 0; nderiv < nDim; ++nderiv)
507 for (
size_t f = 0; f < nConvectiveFields; ++f)
509 qfieldsTmp[f] = qfields[nderiv][f][
p];
512 for (
size_t d = 0; d < nDim; ++d)
515 nDim, d, nderiv, inTmp.data(), qfieldsTmp.data(), muS,
520 for (
size_t f = 0; f < nConvectiveFields; ++f)
522 outArrTmp[f] += normalTmp[d] * outTmp[f];
527 for (
size_t f = 0; f < nConvectiveFields; ++f)
529 outArrTmp[f + nConvectiveFields * d] += outTmp[f];
538 for (
size_t f = 0; f < nConvectiveFields; ++f)
540 outarray[0][f][
p] = outArrTmp[f];
545 for (
size_t d = 0; d < nDim; ++d)
547 for (
size_t f = 0; f < nConvectiveFields; ++f)
550 outArrTmp[f + nConvectiveFields * d];
558 for (
int i = 1; i < n_nonZero + 1; ++i)
560 nonZeroIndex[n_nonZero - i] = nConvectiveFields - i;
VariableConverterSharedPtr m_varConv
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void CalcViscosity(const Array< OneD, const Array< OneD, NekDouble >> &inaverg, Array< OneD, NekDouble > &mu)
virtual ~NavierStokesCFE()
EquationOfStateSharedPtr m_eos
Equation of system for computing temperature.
NavierStokesCFE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
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.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
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 GetPhysicalAV(const Array< OneD, const Array< OneD, NekDouble >> &physfield)
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 GetArtificialViscosity(const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &muav)
static std::string className
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()
std::string m_physicalSensorType
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.
MultiRegions::ContFieldSharedPtr m_C0ProjectExp
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.
void C0Smooth(Array< OneD, NekDouble > &field)
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 Ducros(Array< OneD, NekDouble > &field)
void GetViscousFluxVectorConservVar(const size_t nDim, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &normal)
Return the flux vector for the IP diffusion problem, based on conservative variables.
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 GetViscosityFromTempKernel(const T &temperature, T &mu)
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)
void GetViscousFluxVectorConservVar(const size_t nDim, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex=NullInt1DArray, const Array< OneD, Array< OneD, NekDouble >> &normal=NullNekDoubleArrayOfArray)
std::string m_ViscosityType
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
static Array< OneD, int > NullInt1DArray
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
tinysimd::simd< NekDouble > vec_t
static constexpr struct tinysimd::is_not_aligned_t is_not_aligned
typename abi< ScalarType, width >::type simd