39#include <boost/algorithm/string/predicate.hpp>
61 std::size_t nDim = pFields[0]->GetCoordim(0);
62 std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
65 if (pSession->DefinesSolverInfo(
"HOMOGENEOUS"))
86 m_session->LoadSolverInfo(
"EquationOfState", eosType,
"IdealGas");
95 std::size_t nElements = pFields[0]->GetExpSize();
97 for (std::size_t e = 0; e < nElements; e++)
100 std::size_t expDim = pFields[0]->GetShapeDimension();
107 for (std::size_t i = 0; i < exp3D->GetNtraces(); ++i)
110 h, exp3D->GetGeom3D()->GetEdge(i)->GetVertex(0)->dist(*(
111 exp3D->GetGeom3D()->GetEdge(i)->GetVertex(1))));
120 for (std::size_t i = 0; i < exp2D->GetNtraces(); ++i)
123 h, exp2D->GetGeom2D()->GetEdge(i)->GetVertex(0)->dist(*(
124 exp2D->GetGeom2D()->GetEdge(i)->GetVertex(1))));
133 h = std::min(h, exp1D->
GetGeom1D()->GetVertex(0)->dist(
134 *(exp1D->GetGeom1D()->GetVertex(1))));
139 ASSERTL0(
false,
"Dimension out of bound.")
147 std::size_t nPts = pFields[0]->GetTotPoints();
150 for (std::size_t e = 0; e < pFields[0]->GetExpSize(); e++)
152 std::size_t nElmtPoints = pFields[0]->GetExp(e)->GetTotPoints();
153 std::size_t physOffset = pFields[0]->GetPhys_Offset(e);
154 Vmath::Fill(nElmtPoints, hEle[e], tmp = hElePts + physOffset, 1);
159 pFields[0]->GetFwdBwdTracePhys(hElePts, Fwd, Bwd);
162 std::size_t nBndRegions = pFields[0]->GetBndCondExpansions().size();
164 for (std::size_t i = 0; i < nBndRegions; ++i)
167 std::size_t nBndEdges =
168 pFields[0]->GetBndCondExpansions()[i]->GetExpSize();
170 if (pFields[0]->GetBndConditions()[i]->GetBoundaryConditionType() ==
177 for (std::size_t e = 0; e < nBndEdges; ++e)
179 std::size_t nBndEdgePts = pFields[0]
180 ->GetBndCondExpansions()[i]
184 std::size_t id2 = pFields[0]->GetTrace()->GetPhys_Offset(
185 pFields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
213 const std::size_t nConvectiveFields,
220 std::size_t nCoeffs = fields[0]->GetNcoeffs();
222 for (std::size_t i = 0; i < nConvectiveFields; ++i)
227 for (std::size_t i = 0; i < nConvectiveFields; ++i)
229 fields[i]->BwdTrans(tmp2[i], outarray[i]);
234 const std::size_t nConvectiveFields,
242 if (fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
248 std::size_t nDim = fields[0]->GetCoordim(0);
249 std::size_t nPts = fields[0]->GetTotPoints();
250 std::size_t nCoeffs = fields[0]->GetNcoeffs();
251 std::size_t nScalars = inarray.size();
252 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
260 for (std::size_t i = 0; i < nScalars; ++i)
275 for (std::size_t i = 0; i < nScalars + 1; ++i)
281 for (std::size_t i = 0; i < nConvectiveFields; ++i)
296 for (std::size_t i = 0; i < nConvectiveFields; ++i)
300 for (std::size_t j = 0; j < nDim; ++j)
304 fields[i]->IProductWRTDerivBase(tmpIn, tmpOut);
308 fields[i]->AddTraceIntegral(viscousFlux[i], tmpOut);
309 fields[i]->SetPhysState(
false);
310 if (!fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
312 fields[i]->MultiplyByElmtInvMass(tmpOut, outarray[i]);
324 std::size_t nDim = fields[0]->GetCoordim(0);
325 std::size_t nCoeffs = fields[0]->GetNcoeffs();
326 std::size_t nScalars = inarray.size();
327 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
328 std::size_t nConvectiveFields = fields.size();
338 for (std::size_t i = 0; i < nScalars; ++i)
346 for (std::size_t j = 0; j < nDim; ++j)
348 for (std::size_t i = 0; i < nScalars; ++i)
350 fields[i]->IProductWRTDerivBase(j, inarray[i], tmp1);
352 fields[i]->AddTraceIntegral(numericalFluxO1[j][i], tmp1);
353 fields[i]->SetPhysState(
false);
354 fields[i]->MultiplyByElmtInvMass(tmp1, tmp1);
355 fields[i]->BwdTrans(tmp1, qfields[j][i]);
361 for (std::size_t i = 0; i < nScalars; ++i)
401 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
402 std::size_t nScalars = inarray.size();
406 for (std::size_t i = 0; i < nScalars; ++i)
408 numflux[i] = {pFwd[i]};
412 if (fields[0]->GetBndCondExpansions().size())
414 ApplyBCsO1(fields, inarray, pFwd, pBwd, numflux);
420 for (std::size_t i = 0; i < nScalars; ++i)
423 numericalFluxO1[j][i], 1);
439 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
440 std::size_t nScalars = inarray.size();
449 for (std::size_t i = 0; i < nScalars; ++i)
455 for (std::size_t i = 0; i < nScalars - 1; ++i)
460 std::size_t nBndRegions = fields[i + 1]->GetBndCondExpansions().size();
461 for (std::size_t j = 0; j < nBndRegions; ++j)
464 ->GetBndConditions()[j]
470 std::size_t nBndEdges =
471 fields[i + 1]->GetBndCondExpansions()[j]->GetExpSize();
472 for (std::size_t e = 0; e < nBndEdges; ++e)
474 std::size_t nBndEdgePts = fields[i + 1]
475 ->GetBndCondExpansions()[j]
480 fields[i + 1]->GetBndCondExpansions()[j]->GetPhys_Offset(e);
482 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
483 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
487 fields[i]->GetBndConditions()[j]->GetUserDefined(),
490 fields[i]->GetBndConditions()[j]->GetUserDefined(),
493 fields[i]->GetBndConditions()[j]->GetUserDefined(),
497 for (
int pt = 0; pt < nBndEdgePts; ++pt)
504 else if (boost::iequals(
505 fields[i]->GetBndConditions()[j]->GetUserDefined(),
508 fields[i]->GetBndConditions()[j]->GetUserDefined(),
517 for (std::size_t k = 0; k < nScalars - 1; ++k)
521 ->GetBndCondExpansions()[j]
522 ->UpdatePhys())[id1],
525 ->GetBndCondExpansions()[j]
526 ->UpdatePhys())[id1],
527 1, &scalarVariables[k][id2], 1);
529 1, &scalarVariables[k][id2], 1,
530 &tmp1[0], 1, &tmp1[0], 1);
532 Vmath::Smul(nBndEdgePts, -1.0, &tmp1[0], 1, &tmp1[0],
536 for (std::size_t k = 0; k < nScalars - 1; ++k)
540 &scalarVariables[k][id2], 1,
541 &scalarVariables[k][id2], 1);
546 ->GetBndConditions()[j]
547 ->GetBoundaryConditionType() ==
553 ->GetBndCondExpansions()[j]
554 ->UpdatePhys())[id1],
557 ->GetBndCondExpansions()[j]
558 ->UpdatePhys())[id1],
559 1, &scalarVariables[i][id2], 1);
564 ->GetBndConditions()[j]
565 ->GetBoundaryConditionType() ==
573 else if ((fields[i]->GetBndConditions()[j])
574 ->GetBoundaryConditionType() ==
577 Vmath::Vcopy(nBndEdgePts, &pFwd[i][id2], 1, &fluxO1[i][id2],
582 Vmath::Vmul(nBndEdgePts, &scalarVariables[i][id2], 1,
583 &scalarVariables[i][id2], 1, &tmp1[id2], 1);
585 Vmath::Smul(nBndEdgePts, 0.5, &tmp1[id2], 1, &tmp1[id2], 1);
587 Vmath::Vadd(nBndEdgePts, &tmp2[id2], 1, &tmp1[id2], 1,
595 std::size_t nBndRegions = fields[nScalars]->GetBndCondExpansions().size();
596 for (std::size_t j = 0; j < nBndRegions; ++j)
599 ->GetBndConditions()[j]
605 std::size_t nBndEdges =
606 fields[nScalars]->GetBndCondExpansions()[j]->GetExpSize();
607 for (std::size_t e = 0; e < nBndEdges; ++e)
609 std::size_t nBndEdgePts = fields[nScalars]
610 ->GetBndCondExpansions()[j]
615 fields[nScalars]->GetBndCondExpansions()[j]->GetPhys_Offset(e);
617 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
618 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
622 fields[nScalars]->GetBndConditions()[j]->GetUserDefined(),
626 &scalarVariables[nScalars - 1][id2], 1);
630 else if (fields[nScalars]
631 ->GetBndConditions()[j]
632 ->GetBoundaryConditionType() ==
637 for (std::size_t n = 0; n < nBndEdgePts; ++n)
640 ->GetBndCondExpansions()[j]
641 ->GetPhys()[id1 + n];
642 ene = fields[nScalars]
643 ->GetBndCondExpansions()[j]
644 ->GetPhys()[id1 + n] /
647 scalarVariables[nScalars - 1][id2 + n] =
648 m_eos->GetTemperature(rho, ene);
654 ->GetBndConditions()[j]
655 ->GetBoundaryConditionType() ==
658 fields[nScalars]->GetBndConditions()[j]->GetUserDefined(),
661 fields[nScalars]->GetBndConditions()[j]->GetUserDefined(),
664 Vmath::Vcopy(nBndEdgePts, &scalarVariables[nScalars - 1][id2],
665 1, &fluxO1[nScalars - 1][id2], 1);
669 else if (((fields[nScalars]->GetBndConditions()[j])
670 ->GetBoundaryConditionType() ==
672 boost::iequals(fields[nScalars]
673 ->GetBndConditions()[j]
676 boost::iequals(fields[nScalars]
677 ->GetBndConditions()[j]
682 &fluxO1[nScalars - 1][id2], 1);
699 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
700 std::size_t nVariables = fields.size();
704 for (std::size_t i = 0; i < nVariables - 1; ++i)
720 std::size_t nDim = fields[0]->GetCoordim(0);
721 for (std::size_t i = 1; i < nVariables; ++i)
723 for (std::size_t j = 0; j < nDim; ++j)
726 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
737 qfluxtemp, 1, qfluxtemp, 1);
740 if (fields[0]->GetBndCondExpansions().size())
742 ApplyBCsO2(fields, i, j, qfield[j][i], qFwd, qBwd, qfluxtemp);
746 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
757 const std::size_t var,
const std::size_t dir,
764 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
766 for (std::size_t i = 0; i < nBndRegions; ++i)
769 std::size_t nBndEdges =
770 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
772 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
779 for (std::size_t e = 0; e < nBndEdges; ++e)
781 std::size_t nBndEdgePts = fields[var]
782 ->GetBndCondExpansions()[i]
786 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
787 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
792 ->GetBndConditions()[i]
793 ->GetBoundaryConditionType() ==
796 fields[var]->GetBndConditions()[i]->GetUserDefined(),
799 fields[var]->GetBndConditions()[i]->GetUserDefined(),
803 &qFwd[id2], 1, &penaltyflux[id2], 1);
807 else if ((fields[var]->GetBndConditions()[i])
808 ->GetBoundaryConditionType() ==
811 ASSERTL0(
false,
"Neumann bcs not implemented for LDGNS");
813 else if (boost::iequals(
814 fields[var]->GetBndConditions()[i]->GetUserDefined(),
817 fields[var]->GetBndConditions()[i]->GetUserDefined(),
828 &qFwd[id2], 1, &penaltyflux[id2], 1);
#define ASSERTL0(condition, msg)
void v_DiffuseCoeffs(const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, 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
Array< OneD, Array< OneD, NekDouble > > m_homoDerivs
void v_Diffuse(const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, 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
Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations:
EquationOfStateSharedPtr m_eos
Equation of system for computing temperature.
void NumericalFluxO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &numericalFluxO1, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Builds the numerical flux for the 1st order derivatives.
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, Array< OneD, NekDouble > > m_traceVel
void NumericalFluxO2(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Build the numerical flux for the 2nd order derivatives.
NekDouble m_C11
Penalty coefficient for LDGNS.
void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Diffusion Flux, calculate the physical derivatives.
TensorOfArray3D< NekDouble > m_viscTensor
void ApplyBCsO2(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
Imposes appropriate bcs for the 2nd order derivatives.
void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
Array< OneD, NekDouble > m_traceOneOverH
h scaling for penalty term
void ApplyBCsO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, Array< OneD, NekDouble > > &flux01)
Imposes appropriate bcs for the 1st order derivatives.
static DiffusionSharedPtr create(std::string diffType)
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
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.
SpatialDomains::Geometry1DSharedPtr GetGeom1D() const
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
DiffusionFluxVecCBNS m_fluxVectorNS
SOLVER_UTILS_EXPORT void DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
DiffusionFluxPenaltyNS m_fluxPenaltyNS
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
DiffusionFactory & GetDiffusionFactory()
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
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 Neg(int n, T *x, const int incx)
Negate x = -x.
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 Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
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.
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)