Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::DiffusionIP Class Reference

#include <DiffusionIP.h>

Inheritance diagram for Nektar::SolverUtils::DiffusionIP:
[legend]

Public Member Functions

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
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. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
SOLVER_UTILS_EXPORT void 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=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void 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=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Similar with Diffusion::Diffuse(): calculate diffusion flux The difference is in the outarray: it is the coefficients of basis for DiffuseCoeffs() it is the physics on quadrature points for Diffuse() More...
 
SOLVER_UTILS_EXPORT void 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, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void 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, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
SOLVER_UTILS_EXPORT void 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, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
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)
 
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. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &solution_Aver=NullNekDoubleArrayOfArray, Array< OneD, Array< OneD, NekDouble > > &solution_jump=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void 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)
 Add symmetric flux to field in coeff space. More...
 
SOLVER_UTILS_EXPORT void 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)
 Add symmetric flux to field in coeff physical space. More...
 
SOLVER_UTILS_EXPORT void AddSymmFluxIntegralToOffDiag (const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &Fwdflux, TensorOfArray3D< NekDouble > &Bwdflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void FluxVec (TensorOfArray3D< NekDouble > &fluxvector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 
SOLVER_UTILS_EXPORT void SetFluxVector (DiffusionFluxVecCB fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVectorNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxVectorNS (DiffusionFluxVecCBNS fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxPenaltyNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxPenaltyNS (DiffusionFluxPenaltyNS flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionFluxCons (FuncPointerT func, ObjectPointerT obj)
 
void SetDiffusionFluxCons (DiffusionFluxCons flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionFluxConsTrace (FuncPointerT func, ObjectPointerT obj)
 
void SetDiffusionFluxConsTrace (DiffusionFluxCons flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetArtificialDiffusionVector (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetSpecialBndTreat (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetCalcViscosity (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionSymmFluxCons (FuncPointerT func, ObjectPointerT obj)
 
void SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & GetFluxTensor ()
 
SOLVER_UTILS_EXPORT void GetAVmu (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &muvar, Array< OneD, NekDouble > &MuVarTrace)
 Get the mu of artifical viscosity(AV) More...
 
SOLVER_UTILS_EXPORT void 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)
 Get the average and jump value of conservative variables on trace. More...
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal ()
 Get trace normal. More...
 

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionIP ()
 
void GetPenaltyFactor (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get IP penalty factor based on order. More...
 
void GetPenaltyFactorConst (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get a constant IP penalty factor. More...
 
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) More...
 
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) More...
 
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. More...
 
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. More...
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
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)
 
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)
 
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, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
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. More...
 
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)
 
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)
 
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_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)
 
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_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)
 
virtual const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal ()
 
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. More...
 
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) More...
 
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 More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT 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)
 
virtual SOLVER_UTILS_EXPORT 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)
 
virtual SOLVER_UTILS_EXPORT 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 > > &vFwd, const Array< OneD, Array< OneD, NekDouble > > &vBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
virtual SOLVER_UTILS_EXPORT 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)
 
virtual SOLVER_UTILS_EXPORT 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. More...
 
virtual SOLVER_UTILS_EXPORT 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. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump)
 
virtual SOLVER_UTILS_EXPORT 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 SOLVER_UTILS_EXPORT 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)
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual SOLVER_UTILS_EXPORT void v_GetPrimVar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
 Compute primary derivatives. More...
 
SOLVER_UTILS_EXPORT void GetDivCurl (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
 Compute divergence and curl squared. More...
 

Protected Attributes

std::string m_shockCaptureType
 
NekDouble m_IPSymmFluxCoeff
 
NekDouble m_IP2ndDervCoeff
 
NekDouble m_IPPenaltyCoeff
 
Array< OneD, NekDoublem_MuVarTrace
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
Array< OneD, Array< OneD, NekDouble > > m_traceAver
 
Array< OneD, Array< OneD, NekDouble > > m_traceJump
 
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
 Workspace for v_Diffusion. More...
 
TensorOfArray3D< NekDoublem_wspNumDerivBwd
 Workspace for CallTraceNumFlux. More...
 
TensorOfArray3D< NekDoublem_wspNumDerivFwd
 
Array< OneD, NekDoublem_tracBwdWeightAver
 
Array< OneD, NekDoublem_tracBwdWeightJump
 
Array< OneD, NekDoublem_traceNormDirctnElmtLength
 
Array< OneD, NekDoublem_traceNormDirctnElmtLengthRecip
 
LibUtilities::SessionReaderSharedPtr m_session
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
CalcViscosity m_CalcViscosity
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
NekDouble m_time =0.0
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 

Detailed Description

Definition at line 45 of file DiffusionIP.h.

Constructor & Destructor Documentation

◆ DiffusionIP()

Nektar::SolverUtils::DiffusionIP::DiffusionIP ( )
protected

Definition at line 49 of file DiffusionIP.cpp.

50 {
51 }

Referenced by create().

Member Function Documentation

◆ AddSecondDerivToTrace()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Add second derivative term to trace jump (for DDG scheme)

Definition at line 856 of file DiffusionIP.cpp.

863 {
864  Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
865  Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
866  Array<OneD, NekDouble> tmp{nTracePts, 0.0};
867 
868  Array<OneD, Array<OneD, NekDouble>> elmt2ndDerv{nDim};
869  for (int nd1 = 0; nd1 < nDim; ++nd1)
870  {
871  elmt2ndDerv[nd1] = Array<OneD, NekDouble>{nPts, 0.0};
872  }
873 
874  Array<OneD, Array<OneD, NekDouble>> qtmp{3};
875  for (int nd = 0; nd < 3; ++nd)
876  {
877  qtmp[nd] = NullNekDouble1DArray;
878  }
879  for (int nd2 = 0; nd2 < nDim; ++nd2)
880  {
881  qtmp[nd2] = elmt2ndDerv[nd2];
882  }
883 
884  Vmath::Smul(nTracePts, PenaltyFactor2, m_traceNormDirctnElmtLength, 1, tmp,
885  1);
886  // the derivatives are assumed to be exchangable
887  for (int nd1 = 0; nd1 < nDim; ++nd1)
888  {
889  for (int i = 0; i < nConvectiveFields; ++i)
890  {
891  fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
892 
893  for (int nd2 = nd1; nd2 < nDim; ++nd2)
894  {
895  Vmath::Zero(nTracePts, Bwd, 1);
896  fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd,
897  true, true, false);
898  Vmath::Vmul(nTracePts, tmp, 1, Bwd, 1, Bwd, 1);
899  Vmath::Vvtvp(nTracePts, m_traceNormals[nd2], 1, Bwd, 1,
900  numDerivBwd[nd1][i], 1, numDerivBwd[nd1][i], 1);
901  Vmath::Vmul(nTracePts, tmp, 1, Fwd, 1, Fwd, 1);
902  Vmath::Vvtvm(nTracePts, m_traceNormals[nd2], 1, Fwd, 1,
903  numDerivFwd[nd1][i], 1, numDerivFwd[nd1][i], 1);
904  Vmath::Neg(nTracePts, numDerivFwd[nd1][i], 1);
905 
906  if (nd2 != nd1)
907  {
908  Vmath::Vvtvp(nTracePts, m_traceNormals[nd1], 1, Bwd, 1,
909  numDerivBwd[nd2][i], 1, numDerivBwd[nd2][i],
910  1);
911  Vmath::Vvtvm(nTracePts, m_traceNormals[nd1], 1, Fwd, 1,
912  numDerivFwd[nd2][i], 1, numDerivFwd[nd2][i],
913  1);
914  Vmath::Neg(nTracePts, numDerivFwd[nd2][i], 1);
915  }
916  }
917  }
918  }
919 }
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:120
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:109
static Array< OneD, NekDouble > NullNekDouble1DArray
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.
Definition: Vmath.cpp:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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
Definition: Vmath.cpp:513
void Vvtvm(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)
vvtvm (vector times vector plus vector): z = w*x - y
Definition: Vmath.cpp:541
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

References m_traceNormals, m_traceNormDirctnElmtLength, Vmath::Neg(), Nektar::NullNekDouble1DArray, Vmath::Smul(), Vmath::Vmul(), Vmath::Vvtvm(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by CalcTraceNumFlux().

◆ AddSymmFluxIntegralToCoeff()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Add symmetric flux integration to field (in coefficient space)

Definition at line 565 of file DiffusionIP.cpp.

572 {
573  boost::ignore_unused(nTracePts);
574 
575  size_t nCoeffs = outarray[nConvectiveFields - 1].size();
576  Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
577  Array<OneD, Array<OneD, NekDouble>> tmpfield(nDim);
578  for (int i = 0; i < nDim; ++i)
579  {
580  tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
581  }
582  int nv = 0;
583  for (int j = 0; j < nonZeroIndex.size(); ++j)
584  {
585  nv = nonZeroIndex[j];
586  for (int nd = 0; nd < nDim; ++nd)
587  {
588  Vmath::Zero(nPts, tmpfield[nd], 1);
589 
590  fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
591  tracflux[nd][nv], tmpfield[nd]);
592  fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
593  }
594  fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
595  Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
596  }
597 }
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.
Definition: Vmath.cpp:322

References Vmath::Vadd(), and Vmath::Zero().

Referenced by v_AddDiffusionSymmFluxToCoeff().

◆ AddSymmFluxIntegralToPhys()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Add symmetric flux integration to field (in physical space)

Definition at line 599 of file DiffusionIP.cpp.

606 {
607  boost::ignore_unused(nTracePts);
608 
609  size_t nCoeffs = outarray[nConvectiveFields - 1].size();
610  Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
611  Array<OneD, NekDouble> tmpPhysi{nPts, 0.0};
612  Array<OneD, Array<OneD, NekDouble>> tmpfield{nDim};
613  for (int i = 0; i < nDim; ++i)
614  {
615  tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
616  }
617  for (int j = 0; j < nonZeroIndex.size(); ++j)
618  {
619  int nv = nonZeroIndex[j];
620  for (int nd = 0; nd < nDim; ++nd)
621  {
622  Vmath::Zero(nPts, tmpfield[nd], 1);
623 
624  fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
625  tracflux[nd][nv], tmpfield[nd]);
626  fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
627  }
628  fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
629  fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
630  Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
631  }
632 }

References Vmath::Vadd(), and Vmath::Zero().

Referenced by v_AddDiffusionSymmFluxToPhys().

◆ ApplyFluxBndConds()

void Nektar::SolverUtils::DiffusionIP::ApplyFluxBndConds ( const int  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protected

aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic

Definition at line 927 of file DiffusionIP.cpp.

931 {
932  int nengy = nConvectiveFields-1;
933  int cnt;
934  // Compute boundary conditions for Energy
935  cnt = 0;
936  size_t nBndRegions = fields[nengy]->
937  GetBndCondExpansions().size();
938  for (size_t j = 0; j < nBndRegions; ++j)
939  {
940  if (fields[nengy]->GetBndConditions()[j]->
941  GetBoundaryConditionType() ==
943  {
944  continue;
945  }
946 
947  size_t nBndEdges = fields[nengy]->
948  GetBndCondExpansions()[j]->GetExpSize();
949  for (size_t e = 0; e < nBndEdges; ++e)
950  {
951  size_t nBndEdgePts = fields[nengy]->
952  GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
953 
954  std::size_t id2 = fields[0]->GetTrace()->
955  GetPhys_Offset(fields[0]->GetTraceMap()->
956  GetBndCondIDToGlobalTraceID(cnt++));
957 
958  if (fields[0]->GetBndConditions()[j]->
959  GetUserDefined() =="WallAdiabatic")
960  {
961  Vmath::Zero(nBndEdgePts, &flux[nengy][id2], 1);
962  }
963  }
964  }
965 }

References Nektar::SpatialDomains::ePeriodic, and Vmath::Zero().

Referenced by v_DiffuseTraceFlux().

◆ CalcTraceNumFlux()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Calculate numerical flux on traces.

Definition at line 744 of file DiffusionIP.cpp.

755 {
756  size_t nDim = fields[0]->GetCoordim(0);
757  size_t nPts = fields[0]->GetTotPoints();
758  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
759  size_t nConvectiveFields = fields.size();
760 
761  boost::ignore_unused(inarray);
762  const MultiRegions::AssemblyMapDGSharedPtr TraceMap =
763  fields[0]->GetTraceMap();
764 
765 
766  LibUtilities::Timer timer;
767  timer.Start();
768  // with further restructuring this iniziatilization could be eliminated
769  for (int nd = 0; nd < nDim; ++nd)
770  {
771  for (int i = 0; i < nConvectiveFields; ++i)
772  {
773  Vmath::Zero(nTracePts, m_wspNumDerivBwd[nd][i], 1);
774  Vmath::Zero(nTracePts, m_wspNumDerivFwd[nd][i], 1);
775  }
776  }
777 
778  // could this be pre-allocated?
779  Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
780  Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
781 
782  timer.Stop();
783  timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux_alloc",1);
784 
785 
786  timer.Start();
787  if (fabs(PenaltyFactor2) > 1.0E-12)
788  {
789  AddSecondDerivToTrace(nConvectiveFields, nDim, nPts, nTracePts,
790  PenaltyFactor2, fields, qfield, m_wspNumDerivFwd,
792  }
793 
794  timer.Stop();
795  timer.AccumulateRegion("DiffIP:_AddSecondDerivToTrace",1);
796 
797  for (int nd = 0; nd < nDim; ++nd)
798  {
799  for (int i = 0; i < nConvectiveFields; ++i)
800  {
801  // this sequence of operations is really exepensive,
802  // it should be done collectively for all fields together instead of
803  // one by one
804  timer.Start();
805  fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd, true, true,
806  false);
807  timer.Stop();
808  timer.AccumulateRegion("DiffIP:_GetFwdBwdTracePhys",1);
809 
810  for (size_t p = 0; p < nTracePts; ++p)
811  {
812  m_wspNumDerivBwd[nd][i][p] += 0.5 * Bwd[p];
813  m_wspNumDerivFwd[nd][i][p] += 0.5 * Fwd[p];
814  }
815 
816  timer.Start();
817  TraceMap->GetAssemblyCommDG()->PerformExchange(m_wspNumDerivFwd[nd][i],
818  m_wspNumDerivBwd[nd][i]);
819  timer.Stop();
820  timer.AccumulateRegion("DiffIP:_PerformExchange",1);
821 
822  Vmath::Vadd(nTracePts, m_wspNumDerivFwd[nd][i], 1, m_wspNumDerivBwd[nd][i], 1,
823  m_wspNumDerivFwd[nd][i], 1);
824  }
825  }
826 
827  timer.Start();
828  ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
829  solution_jump);
830  timer.Stop();
831  timer.AccumulateRegion("DiffIP:_ConsVarAveJump",1);
832 
833  for (size_t p = 0; p < nTracePts; ++p)
834  {
835  NekDouble PenaltyFactor = m_IPPenaltyCoeff * m_traceNormDirctnElmtLengthRecip[p]; // load 1x
836 
837  for (size_t f = 0; f < nConvectiveFields; ++f)
838  {
839  NekDouble jumpTmp = solution_jump[f][p] * PenaltyFactor; // load 1x
840  for (size_t d = 0; d < nDim; ++d)
841  {
842  NekDouble tmp = m_traceNormals[d][p] * jumpTmp + m_wspNumDerivFwd[d][f][p]; // load 2x
843  m_wspNumDerivFwd[d][f][p] = tmp; // store 1x
844  }
845  }
846  }
847 
848  timer.Start();
849  // Calculate normal viscous flux
850  m_FunctorDiffusionfluxConsTrace(nDim, solution_Aver, m_wspNumDerivFwd, traceflux,
851  nonZeroIndexflux, m_traceNormals, MuVarTrace);
852  timer.Stop();
853  timer.AccumulateRegion("DiffIP:_FunctorDiffFluxConsTrace",1);
854 }
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:471
SOLVER_UTILS_EXPORT void 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)
Get the average and jump value of conservative variables on trace.
Definition: Diffusion.h:446
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:121
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)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
Definition: DiffusionIP.h:115
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:116
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:47
double NekDouble

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddSecondDerivToTrace(), Nektar::SolverUtils::Diffusion::ConsVarAveJump(), Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxConsTrace, m_IPPenaltyCoeff, m_traceNormals, m_traceNormDirctnElmtLengthRecip, m_wspNumDerivBwd, m_wspNumDerivFwd, CellMLToNektar.cellml_metadata::p, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Vmath::Vadd(), and Vmath::Zero().

Referenced by v_DiffuseTraceFlux().

◆ CalcTraceSymFlux()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Calculate symmetric flux on traces.

Definition at line 535 of file DiffusionIP.cpp.

542 {
543  size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
544 
545  for (int i = 0; i < nConvectiveFields; ++i)
546  {
547  Vmath::Smul(nTracePts, m_IPSymmFluxCoeff, solution_jump[i], 1,
548  solution_jump[i], 1);
549  }
550 
551  m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
552  nonZeroIndexsymm, m_traceNormals);
553 
554  for (int i = 0; i < nConvectiveFields; ++i)
555  {
556  MultiRegions::ExpListSharedPtr tracelist = fields[i]->GetTrace();
557  for (int nd = 0; nd < nDim; ++nd)
558  {
559  tracelist->MultiplyByQuadratureMetric(traceSymflux[nd][i],
560  traceSymflux[nd][i]);
561  }
562  }
563 }
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:474
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.

References Nektar::SolverUtils::Diffusion::m_FunctorSymmetricfluxCons, m_IPSymmFluxCoeff, m_traceNormals, and Vmath::Smul().

Referenced by DiffuseTraceSymmFlux().

◆ ConsVarAve()

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
void Nektar::SolverUtils::DiffusionIP::ConsVarAve ( const size_t  nConvectiveFields,
const T &  Bweight,
const std::vector< T > &  vFwd,
const std::vector< T > &  vBwd,
std::vector< T > &  aver 
)
inline

Calculate the average of conservative variables on traces.

Definition at line 63 of file DiffusionIP.h.

69  {
70  constexpr int nvelst = 1;
71  int nEngy = nConvectiveFields - 1;
72  int nveled = nEngy;
73 
74  T Fweight = 1.0 - Bweight;
75  for (size_t v = 0; v < nEngy; ++v)
76  {
77  aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
78  }
79 
80  T LinternalEngy{};
81  T RinternalEngy{};
82  T AinternalEngy{};
83  for (size_t j = nvelst; j < nveled; ++j)
84  {
85  LinternalEngy += vFwd[j] * vFwd[j];
86  RinternalEngy += vBwd[j] * vBwd[j];
87  AinternalEngy += aver[j] * aver[j];
88  }
89  LinternalEngy *= -0.5 / vFwd[0];
90  RinternalEngy *= -0.5 / vBwd[0];
91  LinternalEngy += vFwd[nEngy];
92  RinternalEngy += vBwd[nEngy];
93 
94  aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
95  aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
96 
97  }

Referenced by v_ConsVarAveJump().

◆ create()

static DiffusionSharedPtr Nektar::SolverUtils::DiffusionIP::create ( std::string  diffType)
inlinestatic

Definition at line 48 of file DiffusionIP.h.

49  {
50  boost::ignore_unused(diffType);
51  return DiffusionSharedPtr(new DiffusionIP());
52  }
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627

References DiffusionIP().

◆ DiffuseTraceSymmFlux()

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

Calculate symmetric flux on traces interface.

Definition at line 518 of file DiffusionIP.cpp.

527 {
528  boost::ignore_unused(inarray, qfield, VolumeFlux, pFwd, pBwd);
529  size_t nDim = fields[0]->GetCoordim(0);
530 
531  CalcTraceSymFlux(nConvectiveFields, nDim, fields, m_traceAver, m_traceJump,
532  nonZeroIndex, SymmFlux);
533 }
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:111
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
Definition: DiffusionIP.h:110

References CalcTraceSymFlux(), m_traceAver, and m_traceJump.

Referenced by v_AddDiffusionSymmFluxToCoeff(), and v_AddDiffusionSymmFluxToPhys().

◆ GetPenaltyFactor()

void Nektar::SolverUtils::DiffusionIP::GetPenaltyFactor ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
Array< OneD, NekDouble > &  factor 
)
protected

Get IP penalty factor based on order.

Definition at line 634 of file DiffusionIP.cpp.

637 {
638  MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
639  std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
640  tracelist->GetExp();
641  size_t ntotTrac = (*traceExp).size();
642  int nTracPnt, noffset;
643 
644  const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
645  fields[0]->GetLocTraceToTraceMap();
646 
647  const Array<OneD, const Array<OneD, int>> LRAdjExpid =
648  locTraceToTraceMap->GetLeftRightAdjacentExpId();
649  const Array<OneD, const Array<OneD, bool>> LRAdjflag =
650  locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
651 
652  std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
653  fields[0]->GetExp();
654 
655  Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
656 
657  NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
658 
659  for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
660  {
661  noffset = tracelist->GetPhys_Offset(ntrace);
662  nTracPnt = tracelist->GetTotPoints(ntrace);
663 
664  factorFwdBwd[0] = 0.0;
665  factorFwdBwd[1] = 0.0;
666 
667  for (int nlr = 0; nlr < 2; ++nlr)
668  {
669  if (LRAdjflag[nlr][ntrace])
670  {
671  int numModes = fields[0]->GetNcoeffs(LRAdjExpid[nlr][ntrace]);
672  NekDouble numModesdir =
673  pow(NekDouble(numModes), (1.0 / spaceDim));
674  factorFwdBwd[nlr] = 1.0 * numModesdir * (numModesdir + 1.0);
675  }
676  }
677 
678  for (int np = 0; np < nTracPnt; ++np)
679  {
680  factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
681  }
682  }
683 }
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr

◆ GetPenaltyFactorConst()

void Nektar::SolverUtils::DiffusionIP::GetPenaltyFactorConst ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
Array< OneD, NekDouble > &  factor 
)
protected

Get a constant IP penalty factor.

Definition at line 685 of file DiffusionIP.cpp.

688 {
689  boost::ignore_unused(fields);
690  Vmath::Fill(factor.size(), m_IPPenaltyCoeff, factor, 1);
691 }
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References Vmath::Fill(), and m_IPPenaltyCoeff.

◆ v_AddDiffusionSymmFluxToCoeff()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 448 of file DiffusionIP.cpp.

456 {
457  if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
458  {
459  size_t nDim = fields[0]->GetCoordim(0);
460  size_t nPts = fields[0]->GetTotPoints();
461  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
462  TensorOfArray3D<NekDouble> traceSymflux{nDim};
463  for (int nd = 0; nd < nDim; ++nd)
464  {
465  traceSymflux[nd] =
466  Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
467  for (int j = 0; j < nConvectiveFields; ++j)
468  {
469  traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
470  }
471  }
472  Array<OneD, int> nonZeroIndex;
473  DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
474  VolumeFlux, traceSymflux, pFwd, pBwd,
475  nonZeroIndex);
476 
477  AddSymmFluxIntegralToCoeff(nConvectiveFields, nDim, nPts, nTracePts,
478  fields, nonZeroIndex, traceSymflux,
479  outarray);
480  }
481 }
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)
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.

References AddSymmFluxIntegralToCoeff(), DiffuseTraceSymmFlux(), and m_IPSymmFluxCoeff.

◆ v_AddDiffusionSymmFluxToPhys()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 483 of file DiffusionIP.cpp.

492 {
493  if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
494  {
495  size_t nDim = fields[0]->GetCoordim(0);
496  size_t nPts = fields[0]->GetTotPoints();
497  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
498  TensorOfArray3D<NekDouble> traceSymflux{nDim};
499  for (int nd = 0; nd < nDim; ++nd)
500  {
501  traceSymflux[nd] =
502  Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
503  for (int j = 0; j < nConvectiveFields; ++j)
504  {
505  traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
506  }
507  }
508  Array<OneD, int> nonZeroIndex;
509  DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
510  VolumeFlux, traceSymflux, pFwd, pBwd,
511  nonZeroIndex);
512 
513  AddSymmFluxIntegralToPhys(nConvectiveFields, nDim, nPts, nTracePts,
514  fields, nonZeroIndex, traceSymflux, outarray);
515  }
516 }
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)

References AddSymmFluxIntegralToPhys(), DiffuseTraceSymmFlux(), and m_IPSymmFluxCoeff.

◆ v_ConsVarAveJump()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 693 of file DiffusionIP.cpp.

699 {
700  // ConsVarAve(nConvectiveFields, nPts, vFwd, vBwd, aver);
701 
702  std::vector<NekDouble> vFwdTmp(nConvectiveFields),
703  vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields) ;
704  for (size_t p = 0; p < nPts; ++p)
705  {
706  // re-arrange data
707  for (size_t f = 0; f < nConvectiveFields; ++f)
708  {
709  vFwdTmp[f] = vFwd[f][p];
710  vBwdTmp[f] = vBwd[f][p];
711  }
712 
713  ConsVarAve(nConvectiveFields, m_tracBwdWeightAver[p], vFwdTmp, vBwdTmp,
714  averTmp);
715 
716  // store
717  for (size_t f = 0; f < nConvectiveFields; ++f)
718  {
719  aver[f][p] = averTmp[f];
720  }
721  }
722 
723  // if this can be make function of points, the nPts loop can be lifted more
724  m_SpecialBndTreat(aver);
725 
726  // note: here the jump is 2.0*(aver-vFwd)
727  // because Viscous wall use a symmetry value as the Bwd,
728  // not the target one
729 
730  for (size_t f = 0; f < nConvectiveFields; ++f)
731  {
732  for (size_t p = 0; p < nPts; ++p)
733  {
734  NekDouble Bweight = m_tracBwdWeightJump[p];
735  NekDouble Fweight = 2.0 - Bweight;
736 
737  NekDouble tmpF = aver[f][p] - vFwd[f][p];
738  NekDouble tmpB = vBwd[f][p] - aver[f][p];
739  jump[f][p] = tmpF * Fweight + tmpB * Bweight;
740  }
741  }
742 }
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:472
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:118
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.
Definition: DiffusionIP.h:63
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:119

References ConsVarAve(), Nektar::SolverUtils::Diffusion::m_SpecialBndTreat, m_tracBwdWeightAver, m_tracBwdWeightJump, and CellMLToNektar.cellml_metadata::p.

◆ v_Diffuse()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 154 of file DiffusionIP.cpp.

161 {
162 
163  LibUtilities::Timer timer;
164  timer.Start();
165  DiffusionIP::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, m_wspDiff,
166  pFwd, pBwd);
167  for (int i = 0; i < nConvectiveFields; ++i)
168  {
169  fields[i]->BwdTrans(m_wspDiff[i], outarray[i]);
170  }
171  timer.Stop();
172  timer.AccumulateRegion("Diffusion IP");
173 }
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)
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Definition: DiffusionIP.h:113

References Nektar::LibUtilities::Timer::AccumulateRegion(), m_wspDiff, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and v_DiffuseCoeffs().

◆ v_DiffuseCalcDerivative()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 368 of file DiffusionIP.cpp.

374 {
375  size_t nConvectiveFields = fields.size();
376  boost::ignore_unused(pFwd, pBwd);
377 
378  size_t nDim = fields[0]->GetCoordim(0);
379 
380  Array<OneD, Array<OneD, NekDouble>> qtmp{3};
381  for (int nd = 0; nd < 3; ++nd)
382  {
383  qtmp[nd] = NullNekDouble1DArray;
384  }
385  for (int i = 0; i < nConvectiveFields; ++i)
386  {
387  for (int nd = 0; nd < nDim; ++nd)
388  {
389  qtmp[nd] = qfield[nd][i];
390  }
391  fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
392  }
393 }

References Nektar::NullNekDouble1DArray.

◆ v_DiffuseCoeffs() [1/2]

virtual void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Referenced by v_Diffuse().

◆ v_DiffuseCoeffs() [2/2]

void Nektar::SolverUtils::DiffusionIP::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,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, int > &  nonZeroIndex 
)
protectedvirtual

Definition at line 300 of file DiffusionIP.cpp.

309 {
310  int i, j;
311  int nDim = fields[0]->GetCoordim(0);
312  int nPts = fields[0]->GetTotPoints();
313  int nCoeffs = fields[0]->GetNcoeffs();
314  int nTracePts = fields[0]->GetTrace()->GetTotPoints();
315 
316  TensorOfArray3D<NekDouble> elmtFlux(nDim);
317  for (j = 0; j < nDim; ++j)
318  {
319  elmtFlux[j] = Array<OneD, Array<OneD, NekDouble>>(nConvectiveFields);
320  for (i = 0; i < nConvectiveFields; ++i)
321  {
322  elmtFlux[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
323  }
324  }
325 
326  DiffuseVolumeFlux(fields,inarray,qfield,elmtFlux,nonZeroIndex);
327 
328  Array<OneD, Array<OneD, NekDouble>> tmpFluxIprdct(nDim);
329  // volume intergration: the nonZeroIndex indicates which flux is nonzero
330  for(i = 0; i < nonZeroIndex.size(); ++i)
331  {
332  int j = nonZeroIndex[i];
333  for (int k = 0; k < nDim; ++k)
334  {
335  tmpFluxIprdct[k] = elmtFlux[k][j];
336  }
337  fields[j]->IProductWRTDerivBase(tmpFluxIprdct,outarray[j]);
338  Vmath::Neg (nCoeffs, outarray[j], 1);
339  }
340 
341  Array<OneD, Array<OneD, NekDouble > > Traceflux(nConvectiveFields);
342  for (int j = 0; j < nConvectiveFields; ++j)
343  {
344  Traceflux[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
345  }
346 
347  DiffuseTraceFlux(fields,inarray,qfield,elmtFlux,
348  Traceflux,vFwd,vBwd,nonZeroIndex);
349 
350  // release qfield, elmtFlux and muvar;
351  for (j = 0; j < nDim; ++j)
352  {
353  elmtFlux[j] = NullNekDoubleArrayOfArray;
354  }
355 
356  for(i = 0; i < nonZeroIndex.size(); ++i)
357  {
358  int j = nonZeroIndex[i];
359 
360  fields[j]->AddTraceIntegral (Traceflux[j], outarray[j]);
361  fields[j]->SetPhysState (false);
362  }
363 
364  AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray,qfield,
365  elmtFlux, outarray, vFwd, vBwd);
366 }
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.
Definition: Diffusion.h:249
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.
Definition: Diffusion.h:236
SOLVER_UTILS_EXPORT void 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)
Add symmetric flux to field in coeff space.
Definition: Diffusion.h:285
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References Nektar::SolverUtils::Diffusion::AddDiffusionSymmFluxToCoeff(), Nektar::SolverUtils::Diffusion::DiffuseTraceFlux(), Nektar::SolverUtils::Diffusion::DiffuseVolumeFlux(), Vmath::Neg(), and Nektar::NullNekDoubleArrayOfArray.

◆ v_DiffuseTraceFlux() [1/2]

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Definition at line 421 of file DiffusionIP.cpp.

429 {
430  boost::ignore_unused(VolumeFlux);
431 
432  TensorOfArray3D<NekDouble> traceflux3D(1);
433  traceflux3D[0] = TraceFlux;
434 
435  size_t nConvectiveFields = fields.size();
436 
437  LibUtilities::Timer timer;
438  timer.Start();
440  fields, inarray, qfield, pFwd, pBwd, m_MuVarTrace,
441  nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
442  timer.Stop();
443  timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux",1);
444 
445  ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
446 }
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
Array< OneD, NekDouble > m_MuVarTrace
Definition: DiffusionIP.h:108
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.

References Nektar::LibUtilities::Timer::AccumulateRegion(), ApplyFluxBndConds(), CalcTraceNumFlux(), m_IP2ndDervCoeff, m_MuVarTrace, m_traceAver, m_traceJump, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

◆ v_DiffuseTraceFlux() [2/2]

void Nektar::SolverUtils::DiffusionIP::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 
)
protected

◆ v_DiffuseVolumeFlux()

void Nektar::SolverUtils::DiffusionIP::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 
)
protectedvirtual

Diffusion Volume Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 395 of file DiffusionIP.cpp.

400 {
401  size_t nDim = fields[0]->GetCoordim(0);
402  size_t nPts = fields[0]->GetTotPoints();
403 
404  Array<OneD, NekDouble> muvar = NullNekDouble1DArray;
406  {
407  muvar = Array<OneD, NekDouble>{nPts, 0.0};
408  GetAVmu(fields, inarray, muvar, m_MuVarTrace);
409  }
410 
411  Array<OneD, Array<OneD, NekDouble>> tmparray2D = NullNekDoubleArrayOfArray;
412 
413  LibUtilities::Timer timer;
414  timer.Start();
415  m_FunctorDiffusionfluxCons(nDim, inarray, qfield, VolumeFlux, nonZeroIndex,
416  tmparray2D, muvar);
417  timer.Stop();
418  timer.AccumulateRegion("DiffIP:_FunctorDiffFluxCons",1);
419 }
SOLVER_UTILS_EXPORT void GetAVmu(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &muvar, Array< OneD, NekDouble > &MuVarTrace)
Get the mu of artifical viscosity(AV)
Definition: Diffusion.cpp:115
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
Definition: Diffusion.h:469
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:470

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::SolverUtils::Diffusion::GetAVmu(), Nektar::SolverUtils::Diffusion::m_ArtificialDiffusionVector, Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxCons, m_MuVarTrace, Nektar::NullNekDouble1DArray, Nektar::NullNekDoubleArrayOfArray, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

◆ v_GetTraceNormal()

virtual const Array<OneD, const Array<OneD, NekDouble> >& Nektar::SolverUtils::DiffusionIP::v_GetTraceNormal ( )
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 268 of file DiffusionIP.h.

269  {
270  return m_traceNormals;
271  }

References m_traceNormals.

◆ v_InitObject()

void Nektar::SolverUtils::DiffusionIP::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 53 of file DiffusionIP.cpp.

56 {
57  m_session = pSession;
58 
59  m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
60 
61  m_session->LoadParameter("IPSymmFluxCoeff", m_IPSymmFluxCoeff, 0.0); // 0.5
62 
63  m_session->LoadParameter("IP2ndDervCoeff", m_IP2ndDervCoeff,
64  0.0); // 1.0/12.0
65 
66  m_session->LoadParameter("IPPenaltyCoeff", m_IPPenaltyCoeff,
67  4.0); // 1.0/12.0
68 
69  // Setting up the normals
70  size_t nDim = pFields[0]->GetCoordim(0);
71  size_t nVariable = pFields.size();
72  size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
73 
74  m_traceNormals = Array<OneD, Array<OneD, NekDouble>>{nDim};
75  for (int i = 0; i < nDim; ++i)
76  {
77  m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
78  }
79  m_traceAver = Array<OneD, Array<OneD, NekDouble>>{nVariable};
80  m_traceJump = Array<OneD, Array<OneD, NekDouble>>{nVariable};
81  for (int i = 0; i < nVariable; ++i)
82  {
83  m_traceAver[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
84  m_traceJump[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
85  }
86 
87  pFields[0]->GetTrace()->GetNormals(m_traceNormals);
88  Array<OneD, NekDouble> lengthFwd{nTracePts, 0.0};
89  Array<OneD, NekDouble> lengthBwd{nTracePts, 0.0};
90  pFields[0]->GetTrace()->GetElmtNormalLength(lengthFwd, lengthBwd);
91  pFields[0]->PeriodicBwdCopy(lengthFwd, lengthBwd);
92 
94  pFields[0]->GetTraceMap();
95  TraceMap->GetAssemblyCommDG()->PerformExchange(lengthFwd, lengthBwd);
96 
97  Vmath::Vadd(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthFwd, 1);
98  m_traceNormDirctnElmtLength = lengthFwd;
100  Vmath::Sdiv(nTracePts, 1.0, m_traceNormDirctnElmtLength, 1,
102 
103  m_tracBwdWeightAver = Array<OneD, NekDouble>{nTracePts, 0.0};
104  m_tracBwdWeightJump = Array<OneD, NekDouble>{nTracePts, 0.0};
105  pFields[0]->GetBwdWeight(m_tracBwdWeightAver, m_tracBwdWeightJump);
106  Array<OneD, NekDouble> tmpBwdWeight{nTracePts, 0.0};
107  Array<OneD, NekDouble> tmpBwdWeightJump{nTracePts, 0.0};
108  for (int i = 1; i < nVariable; ++i)
109  {
110  pFields[i]->GetBwdWeight(tmpBwdWeight, tmpBwdWeightJump);
111  Vmath::Vsub(nTracePts, tmpBwdWeight, 1, m_tracBwdWeightAver, 1,
112  tmpBwdWeight, 1);
113  Vmath::Vabs(nTracePts, tmpBwdWeight, 1, tmpBwdWeight, 1);
114  NekDouble norm = 0.0;
115  for (int j = 0; j < nTracePts; ++j)
116  {
117  norm += tmpBwdWeight[j];
118  }
119  ASSERTL0(norm < 1.0E-11,
120  "different BWD for different variable not coded yet");
121  }
122 
125  {
126  m_MuVarTrace = Array<OneD, NekDouble>{nTracePts, 0.0};
127  }
128 
129  // workspace for v_diffuse
130  size_t nCoeffs = pFields[0]->GetNcoeffs();
131  m_wspDiff = Array<OneD, Array<OneD, NekDouble>>{nVariable};
132  for (int i = 0; i < nVariable; ++i)
133  {
134  m_wspDiff[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
135  }
136 
137  // workspace for callnumtraceflux
138  m_wspNumDerivBwd = TensorOfArray3D<NekDouble>{nDim};
139  m_wspNumDerivFwd = TensorOfArray3D<NekDouble>{nDim};
140  for (int nd = 0; nd < nDim; ++nd)
141  {
142  m_wspNumDerivBwd[nd] =
143  Array<OneD, Array<OneD, NekDouble>>{nVariable};
144  m_wspNumDerivFwd[nd] =
145  Array<OneD, Array<OneD, NekDouble>>{nVariable};
146  for (int i = 0; i < nVariable; ++i)
147  {
148  m_wspNumDerivBwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
149  m_wspNumDerivFwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
150  }
151  }
152 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:122
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:493
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:291
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.
Definition: Vmath.cpp:372

References ASSERTL0, Nektar::SolverUtils::Diffusion::m_ArtificialDiffusionVector, m_IP2ndDervCoeff, m_IPPenaltyCoeff, m_IPSymmFluxCoeff, m_MuVarTrace, m_session, m_shockCaptureType, m_tracBwdWeightAver, m_tracBwdWeightJump, m_traceAver, m_traceJump, m_traceNormals, m_traceNormDirctnElmtLength, m_traceNormDirctnElmtLengthRecip, m_wspDiff, m_wspNumDerivBwd, m_wspNumDerivFwd, Nektar::NullNekDouble1DArray, Vmath::Sdiv(), Vmath::Vabs(), Vmath::Vadd(), and Vmath::Vsub().

Member Data Documentation

◆ m_IP2ndDervCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IP2ndDervCoeff
protected

Definition at line 105 of file DiffusionIP.h.

Referenced by v_DiffuseTraceFlux(), and v_InitObject().

◆ m_IPPenaltyCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IPPenaltyCoeff
protected

Definition at line 106 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), GetPenaltyFactorConst(), and v_InitObject().

◆ m_IPSymmFluxCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IPSymmFluxCoeff
protected

◆ m_MuVarTrace

Array<OneD, NekDouble> Nektar::SolverUtils::DiffusionIP::m_MuVarTrace
protected

Definition at line 108 of file DiffusionIP.h.

Referenced by v_DiffuseTraceFlux(), v_DiffuseVolumeFlux(), and v_InitObject().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionIP::m_session
protected

Definition at line 122 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_shockCaptureType

std::string Nektar::SolverUtils::DiffusionIP::m_shockCaptureType
protected

Definition at line 103 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_tracBwdWeightAver

Array<OneD, NekDouble> Nektar::SolverUtils::DiffusionIP::m_tracBwdWeightAver
protected

Definition at line 118 of file DiffusionIP.h.

Referenced by v_ConsVarAveJump(), and v_InitObject().

◆ m_tracBwdWeightJump

Array<OneD, NekDouble> Nektar::SolverUtils::DiffusionIP::m_tracBwdWeightJump
protected

Definition at line 119 of file DiffusionIP.h.

Referenced by v_ConsVarAveJump(), and v_InitObject().

◆ m_traceAver

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionIP::m_traceAver
protected

Definition at line 110 of file DiffusionIP.h.

Referenced by DiffuseTraceSymmFlux(), v_DiffuseTraceFlux(), and v_InitObject().

◆ m_traceJump

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionIP::m_traceJump
protected

Definition at line 111 of file DiffusionIP.h.

Referenced by DiffuseTraceSymmFlux(), v_DiffuseTraceFlux(), and v_InitObject().

◆ m_traceNormals

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionIP::m_traceNormals
protected

◆ m_traceNormDirctnElmtLength

Array<OneD, NekDouble> Nektar::SolverUtils::DiffusionIP::m_traceNormDirctnElmtLength
protected

Definition at line 120 of file DiffusionIP.h.

Referenced by AddSecondDerivToTrace(), and v_InitObject().

◆ m_traceNormDirctnElmtLengthRecip

Array<OneD, NekDouble> Nektar::SolverUtils::DiffusionIP::m_traceNormDirctnElmtLengthRecip
protected

Definition at line 121 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ m_wspDiff

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionIP::m_wspDiff
protected

Workspace for v_Diffusion.

Definition at line 113 of file DiffusionIP.h.

Referenced by v_Diffuse(), and v_InitObject().

◆ m_wspNumDerivBwd

TensorOfArray3D<NekDouble> Nektar::SolverUtils::DiffusionIP::m_wspNumDerivBwd
protected

Workspace for CallTraceNumFlux.

Definition at line 115 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ m_wspNumDerivFwd

TensorOfArray3D<NekDouble> Nektar::SolverUtils::DiffusionIP::m_wspNumDerivFwd
protected

Definition at line 116 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ type

std::string Nektar::SolverUtils::DiffusionIP::type
static
Initial value:
"InteriorPenalty", DiffusionIP::create)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:48
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

Definition at line 54 of file DiffusionIP.h.