39 #include <boost/core/ignore_unused.hpp>
40 #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"))
83 m_session->LoadSolverInfo(
"EquationOfState",
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)
109 h = std::min(h, exp3D->GetGeom3D()->GetEdge(i)->GetVertex(0)->
110 dist(*(exp3D->GetGeom3D()->GetEdge(i)->GetVertex(1))));
119 for(std::size_t i = 0; i < exp2D->GetNtraces(); ++i)
121 h = std::min(h, exp2D->GetGeom2D()->GetEdge(i)->GetVertex(0)->
122 dist(*(exp2D->GetGeom2D()->GetEdge(i)->GetVertex(1))));
131 h = std::min(h, exp1D->
GetGeom1D()->GetVertex(0)->
132 dist(*(exp1D->GetGeom1D()->GetVertex(1))));
137 ASSERTL0(
false,
"Dimension out of bound.")
145 std::size_t nPts = pFields[0]->GetTotPoints();
148 for (std::size_t e = 0; e < pFields[0]->GetExpSize(); e++)
150 std::size_t nElmtPoints = pFields[0]->GetExp(e)->GetTotPoints();
151 std::size_t physOffset = pFields[0]->GetPhys_Offset(e);
153 tmp = hElePts + physOffset, 1);
158 pFields[0]->GetFwdBwdTracePhys(hElePts, Fwd, Bwd);
161 std::size_t nBndRegions = pFields[0]->GetBndCondExpansions().size();
163 for (std::size_t i = 0; i < nBndRegions; ++i)
166 std::size_t nBndEdges = pFields[0]->
167 GetBndCondExpansions()[i]->GetExpSize();
169 if (pFields[0]->GetBndConditions()[i]->GetBoundaryConditionType()
176 for (std::size_t e = 0; e < nBndEdges; ++e)
178 std::size_t nBndEdgePts = pFields[0]->
179 GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
181 std::size_t id2 = pFields[0]->GetTrace()->
182 GetPhys_Offset(pFields[0]->GetTraceMap()->
183 GetBndCondIDToGlobalTraceID(cnt++));
211 const std::size_t nConvectiveFields,
218 std::size_t nCoeffs = fields[0]->GetNcoeffs();
220 for (std::size_t i = 0; i < nConvectiveFields; ++i)
225 for (std::size_t i = 0; i < nConvectiveFields; ++i)
227 fields[i]->BwdTrans (tmp2[i], outarray[i]);
232 const std::size_t nConvectiveFields,
239 std::size_t nDim = fields[0]->GetCoordim(0);
240 std::size_t nPts = fields[0]->GetTotPoints();
241 std::size_t nCoeffs = fields[0]->GetNcoeffs();
242 std::size_t nScalars = inarray.size();
243 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
252 for (std::size_t i = 0; i < nScalars; ++i)
267 for (std::size_t i = 0; i < nScalars+1; ++i)
273 for (std::size_t i = 0; i < nConvectiveFields; ++i)
288 for (std::size_t i = 0; i < nConvectiveFields; ++i)
292 for (std::size_t j = 0; j < nDim; ++j)
296 fields[i]->IProductWRTDerivBase(tmpIn, tmpOut);
300 fields[i]->AddTraceIntegral (viscousFlux[i], tmpOut);
301 fields[i]->SetPhysState (
false);
302 fields[i]->MultiplyByElmtInvMass(tmpOut, outarray[i]);
313 std::size_t nDim = fields[0]->GetCoordim(0);
314 std::size_t nCoeffs = fields[0]->GetNcoeffs();
315 std::size_t nScalars = inarray.size();
316 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
317 std::size_t nConvectiveFields = fields.size();
327 for (std::size_t i = 0; i < nScalars; ++i)
335 for (std::size_t j = 0; j < nDim; ++j)
337 for (std::size_t i = 0; i < nScalars; ++i)
339 fields[i]->IProductWRTDerivBase (j, inarray[i], tmp1);
341 fields[i]->AddTraceIntegral (numericalFluxO1[j][i], tmp1);
342 fields[i]->SetPhysState (
false);
343 fields[i]->MultiplyByElmtInvMass(tmp1, tmp1);
344 fields[i]->BwdTrans (tmp1, qfields[j][i]);
350 for (std::size_t i = 0; i < nScalars; ++i)
365 boost::ignore_unused(fields, nonZeroIndex);
380 boost::ignore_unused(inarray, qfields, nonZeroIndex);
395 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
396 std::size_t nScalars = inarray.size();
400 for (std::size_t i = 0; i < nScalars; ++i)
402 numflux[i] = {pFwd[i]};
406 if (fields[0]->GetBndCondExpansions().size())
408 ApplyBCsO1(fields, inarray, pFwd, pBwd, numflux);
414 for (std::size_t i = 0; i < nScalars; ++i)
417 numericalFluxO1[j][i], 1);
433 boost::ignore_unused(pBwd);
435 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
436 std::size_t nScalars = inarray.size();
445 for (std::size_t i = 0; i < nScalars; ++i)
451 for (std::size_t i = 0; i < nScalars-1; ++i)
456 std::size_t nBndRegions = fields[i+1]->
457 GetBndCondExpansions().size();
458 for (std::size_t j = 0; j < nBndRegions; ++j)
460 if (fields[i+1]->GetBndConditions()[j]->
466 std::size_t nBndEdges = fields[i+1]->
467 GetBndCondExpansions()[j]->GetExpSize();
468 for (std::size_t e = 0; e < nBndEdges; ++e)
470 std::size_t nBndEdgePts = fields[i+1]->
471 GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
473 std::size_t id1 = fields[i+1]->
474 GetBndCondExpansions()[j]->GetPhys_Offset(e);
476 std::size_t id2 = fields[0]->GetTrace()->
477 GetPhys_Offset(fields[0]->GetTraceMap()->
478 GetBndCondIDToGlobalTraceID(cnt++));
480 if (boost::iequals(fields[i]->GetBndConditions()[j]->
481 GetUserDefined(),
"WallViscous") ||
482 boost::iequals(fields[i]->GetBndConditions()[j]->
483 GetUserDefined(),
"WallAdiabatic"))
486 Vmath::Zero(nBndEdgePts, &scalarVariables[i][id2], 1);
490 boost::iequals(fields[i]->GetBndConditions()[j]->
491 GetUserDefined(),
"Wall") ||
492 boost::iequals(fields[i]->GetBndConditions()[j]->
493 GetUserDefined(),
"Symmetry"))
501 for (std::size_t k = 0; k < nScalars-1; ++k)
504 &(fields[k+1]->GetBndCondExpansions()[j]->
505 UpdatePhys())[id1], 1,
506 &(fields[0]->GetBndCondExpansions()[j]->
507 UpdatePhys())[id1], 1,
508 &scalarVariables[k][id2], 1);
511 &scalarVariables[k][id2], 1,
520 for (std::size_t k = 0; k < nScalars-1; ++k)
525 &scalarVariables[k][id2], 1,
526 &scalarVariables[k][id2], 1);
530 else if (fields[i]->GetBndConditions()[j]->
531 GetBoundaryConditionType() ==
536 &(fields[i+1]->GetBndCondExpansions()[j]->
537 UpdatePhys())[id1], 1,
538 &(fields[0]->GetBndCondExpansions()[j]->
539 UpdatePhys())[id1], 1,
540 &scalarVariables[i][id2], 1);
544 if (fields[i]->GetBndConditions()[j]->
545 GetBoundaryConditionType() ==
549 &scalarVariables[i][id2], 1,
554 else if ((fields[i]->GetBndConditions()[j])->
555 GetBoundaryConditionType() ==
565 &scalarVariables[i][id2], 1,
566 &scalarVariables[i][id2], 1,
583 std::size_t nBndRegions = fields[nScalars]->
584 GetBndCondExpansions().size();
585 for (std::size_t j = 0; j < nBndRegions; ++j)
587 if (fields[nScalars]->GetBndConditions()[j]->
588 GetBoundaryConditionType() ==
594 std::size_t nBndEdges = fields[nScalars]->
595 GetBndCondExpansions()[j]->GetExpSize();
596 for (std::size_t e = 0; e < nBndEdges; ++e)
598 std::size_t nBndEdgePts = fields[nScalars]->
599 GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
601 std::size_t id1 = fields[nScalars]->
602 GetBndCondExpansions()[j]->GetPhys_Offset(e);
604 std::size_t id2 = fields[0]->GetTrace()->
605 GetPhys_Offset(fields[0]->GetTraceMap()->
606 GetBndCondIDToGlobalTraceID(cnt++));
609 if (boost::iequals(fields[nScalars]->GetBndConditions()[j]->
610 GetUserDefined(),
"WallViscous"))
614 &scalarVariables[nScalars-1][id2], 1);
618 else if (fields[nScalars]->GetBndConditions()[j]->
619 GetBoundaryConditionType() ==
624 for (std::size_t n = 0; n < nBndEdgePts; ++n)
627 GetBndCondExpansions()[j]->
629 ene = fields[nScalars]->
630 GetBndCondExpansions()[j]->
631 GetPhys()[id1 +n] / rho - tmp2[id2+n];
632 scalarVariables[nScalars-1][id2+n] =
633 m_eos->GetTemperature(rho, ene);
638 if (fields[nScalars]->GetBndConditions()[j]->
640 !boost::iequals(fields[nScalars]->GetBndConditions()[j]
641 ->GetUserDefined(),
"WallAdiabatic"))
644 &scalarVariables[nScalars-1][id2], 1,
645 &fluxO1[nScalars-1][id2], 1);
650 else if (((fields[nScalars]->GetBndConditions()[j])->
652 boost::iequals(fields[nScalars]->GetBndConditions()[j]->
653 GetUserDefined(),
"WallAdiabatic"))
656 &pFwd[nScalars-1][id2], 1,
657 &fluxO1[nScalars-1][id2], 1);
675 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
676 std::size_t nVariables = fields.size();
680 for (std::size_t i = 0; i < nVariables-1; ++i)
696 std::size_t nDim = fields[0]->GetCoordim(0);
697 for (std::size_t i = 1; i < nVariables; ++i)
700 for (std::size_t j = 0; j < nDim; ++j)
703 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
714 qfluxtemp, 1, qfluxtemp, 1);
717 if (fields[0]->GetBndCondExpansions().size())
719 ApplyBCsO2(fields, i, j, qfield[j][i], qFwd, qBwd, qfluxtemp);
723 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
735 const std::size_t var,
736 const std::size_t dir,
742 boost::ignore_unused(qfield, qBwd);
745 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
747 for (std::size_t i = 0; i < nBndRegions; ++i)
750 std::size_t nBndEdges = fields[var]->
751 GetBndCondExpansions()[i]->GetExpSize();
753 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType()
760 for (std::size_t e = 0; e < nBndEdges; ++e)
762 std::size_t nBndEdgePts = fields[var]->
763 GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
765 std::size_t id2 = fields[0]->GetTrace()->
766 GetPhys_Offset(fields[0]->GetTraceMap()->
767 GetBndCondIDToGlobalTraceID(cnt++));
771 if (fields[var]->GetBndConditions()[i]->
773 && !boost::iequals(fields[var]->GetBndConditions()[i]->
774 GetUserDefined(),
"WallAdiabatic"))
779 &penaltyflux[id2], 1);
783 else if ((fields[var]->GetBndConditions()[i])->
787 "Neumann bcs not implemented for LDGNS");
798 else if (boost::iequals(fields[var]->GetBndConditions()[i]->
799 GetUserDefined(),
"WallAdiabatic"))
811 &penaltyflux[id2], 1);
828 int nDim = fields[0]->GetCoordim(0);
829 int nPts = fields[0]->GetTotPoints();
830 for(
int i = 0; i < nDim; ++i)
833 primVar[i] = inarray[i];
#define ASSERTL0(condition, msg)
Array< OneD, Array< OneD, NekDouble > > m_homoDerivs
EquationOfStateSharedPtr m_eos
Equation of system for computing temperature.
virtual void v_GetPrimVar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
Compute primary variables.
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
virtual 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)
Diffusion Volume Flux.
virtual 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)
Diffusion term Trace Flux.
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.
virtual 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)
TensorOfArray3D< NekDouble > m_viscTensor
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
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.
virtual 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)
Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations:
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
virtual 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)
Diffusion Flux, calculate the physical derivatives.
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 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.
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.
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)
DiffusionFluxVecCBNS m_fluxVectorNS
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()
The above copyright notice and this permission notice shall be included.
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)
vvtvvtp (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/y.
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)