36 #ifndef NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
37 #define NEKTAR_SOLVERUTILS_DIFFUSIONWEAKDG
50 boost::ignore_unused(diffType);
57 template <
class T,
typename =
typename std::enable_if
59 std::is_floating_point<T>::value ||
64 const size_t nConvectiveFields,
66 const std::vector<T>& vFwd,
67 const std::vector<T>& vBwd,
70 constexpr
int nvelst = 1;
71 int nEngy = nConvectiveFields - 1;
74 T Fweight = 1.0 - Bweight;
75 for (
size_t v = 0; v < nEngy; ++v)
77 aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
83 for (
size_t j = nvelst; j < nveled; ++j)
85 LinternalEngy += vFwd[j] * vFwd[j];
86 RinternalEngy += vBwd[j] * vBwd[j];
87 AinternalEngy += aver[j] * aver[j];
89 LinternalEngy *= -0.5 / vFwd[0];
90 RinternalEngy *= -0.5 / vBwd[0];
91 LinternalEngy += vFwd[nEngy];
92 RinternalEngy += vBwd[nEngy];
94 aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
95 aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
135 const std::size_t nConvectiveFields,
const size_t nDim,
136 const size_t nPts,
const size_t nTracePts,
144 const std::size_t nConvectiveFields,
const size_t nDim,
145 const size_t nPts,
const size_t nTracePts,
153 const std::size_t nConvectiveFields,
const size_t nDim,
162 const std::size_t nConvectiveFields,
177 const std::size_t nConvectiveFields,
185 const std::size_t nConvectiveFields,
193 const std::size_t nConvectiveFields,
220 const int nConvectiveFields,
235 const std::size_t nConvectiveFields,
245 const std::size_t nConvectiveFields,
262 const std::size_t nConvectiveFields,
const size_t npnts,
307 const std::size_t nConvectiveFields,
const size_t nDim,
308 const size_t nPts,
const size_t nTracePts,
316 const int nConvectiveFields,
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
static DiffusionSharedPtr create(std::string diffType)
virtual void v_Diffuse(const std::size_t nConvectiveFields, 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)
Array< OneD, NekDouble > m_tracBwdWeightAver
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
void v_DiffuseTraceFlux(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > &qFwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble > > > &qBwd, const Array< OneD, NekDouble > &MuAVTrace, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &Aver, const Array< OneD, Array< OneD, NekDouble >> &Jump)
void ApplyFluxBndConds(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, Array< OneD, NekDouble > > &flux)
aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic
void GetPenaltyFactor(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get IP penalty factor based on order.
void ConsVarAve(const size_t nConvectiveFields, const T &Bweight, const std::vector< T > &vFwd, const std::vector< T > &vBwd, std::vector< T > &aver)
Calculate the average of conservative variables on traces.
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.
void GetPenaltyFactorConst(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get a constant IP penalty factor.
void AddSecondDerivToTrace(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &numDerivFwd, TensorOfArray3D< NekDouble > &numDerivBwd)
Add second derivative term to trace jump (for DDG scheme)
Array< OneD, NekDouble > m_MuVarTrace
NekDouble m_IPSymmFluxCoeff
virtual const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
NekDouble m_IP2ndDervCoeff
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
NekDouble m_IPPenaltyCoeff
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Add symmetric flux integration to field (in coefficient space)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
virtual void v_ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
virtual void v_AddDiffusionSymmFluxToPhys(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
void DiffuseTraceSymmFlux(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
Calculate symmetric flux on traces interface.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
LibUtilities::SessionReaderSharedPtr m_session
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump, Array< OneD, int > &nonZeroIndexsymm, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &traceSymflux)
Calculate symmetric flux on traces.
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Array< OneD, NekDouble > m_tracBwdWeightJump
virtual void v_AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
virtual void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
virtual void v_DiffuseCoeffs(const std::size_t nConvectiveFields, 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)
void CalcTraceNumFlux(const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble >> &vFwd, const Array< OneD, Array< OneD, NekDouble >> &vBwd, const Array< OneD, NekDouble > &MuVarTrace, Array< OneD, int > &nonZeroIndexflux, TensorOfArray3D< NekDouble > &traceflux, Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump)
Calculate numerical flux on traces.
void AddSymmFluxIntegralToPhys(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Add symmetric flux integration to field (in physical space)
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)
std::string m_shockCaptureType
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
The above copyright notice and this permission notice shall be included.