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

#include <DiffusionIP.h>

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

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionIP ()
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 
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) override
 
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) override
 
void v_DiffuseVolumeFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
 Diffusion Volume Flux. More...
 
void v_DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex) override
 Diffusion term Trace Flux. More...
 
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) override
 Diffusion Flux, calculate the physical derivatives. More...
 
const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal () override
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)=0
 
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)=0
 
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_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_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 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 void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal ()
 

Private 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...
 
void GetPenaltyFactor (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get IP penalty factor based on order. More...
 
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)
 
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, 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 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)
 
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)
 
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...
 
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, 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...
 

Private Attributes

NekDouble m_IPSymmFluxCoeff
 
NekDouble m_IP2ndDervCoeff
 
NekDouble m_IPPenaltyCoeff
 
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
 

Additional Inherited Members

- 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 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=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, 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 SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
SOLVER_UTILS_EXPORT TensorOfArray3D< NekDouble > & GetFluxTensor ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal ()
 Get trace normal. More...
 
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 SetSpecialBndTreat (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionSymmFluxCons (FuncPointerT func, ObjectPointerT obj)
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 
DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
NekDouble m_time = 0.0
 

Detailed Description

Definition at line 43 of file DiffusionIP.h.

Constructor & Destructor Documentation

◆ DiffusionIP()

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

Definition at line 48 of file DiffusionIP.cpp.

49{
50}

Referenced by create().

Member Function Documentation

◆ AddDiffusionSymmFluxToCoeff()

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

Definition at line 386 of file DiffusionIP.cpp.

394{
395 if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
396 {
397 size_t nDim = fields[0]->GetCoordim(0);
398 size_t nPts = fields[0]->GetTotPoints();
399 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
400 TensorOfArray3D<NekDouble> traceSymflux{nDim};
401 for (int nd = 0; nd < nDim; ++nd)
402 {
403 traceSymflux[nd] =
404 Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
405 for (int j = 0; j < nConvectiveFields; ++j)
406 {
407 traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
408 }
409 }
410 Array<OneD, int> nonZeroIndex;
411 DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
412 VolumeFlux, traceSymflux, pFwd, pBwd,
413 nonZeroIndex);
414
415 AddSymmFluxIntegralToCoeff(nConvectiveFields, nDim, nPts, nTracePts,
416 fields, nonZeroIndex, traceSymflux,
417 outarray);
418 }
419}
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.

Referenced by v_DiffuseCoeffs().

◆ AddDiffusionSymmFluxToPhys()

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

Definition at line 421 of file DiffusionIP.cpp.

429{
430 if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
431 {
432 size_t nDim = fields[0]->GetCoordim(0);
433 size_t nPts = fields[0]->GetTotPoints();
434 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
435 TensorOfArray3D<NekDouble> traceSymflux{nDim};
436 for (int nd = 0; nd < nDim; ++nd)
437 {
438 traceSymflux[nd] =
439 Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
440 for (int j = 0; j < nConvectiveFields; ++j)
441 {
442 traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
443 }
444 }
445 Array<OneD, int> nonZeroIndex;
446 DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
447 VolumeFlux, traceSymflux, pFwd, pBwd,
448 nonZeroIndex);
449
450 AddSymmFluxIntegralToPhys(nConvectiveFields, nDim, nPts, nTracePts,
451 fields, nonZeroIndex, traceSymflux, outarray);
452 }
453}
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.

◆ 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 
)
private

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

Definition at line 784 of file DiffusionIP.cpp.

791{
792 Array<OneD, NekDouble> Fwd(nTracePts, 0.0);
793 Array<OneD, NekDouble> Bwd(nTracePts, 0.0);
794 std::vector<NekDouble> tmp(nTracePts);
795
796 Array<OneD, Array<OneD, NekDouble>> elmt2ndDerv{nDim};
797 for (int nd1 = 0; nd1 < nDim; ++nd1)
798 {
799 elmt2ndDerv[nd1] = Array<OneD, NekDouble>{nPts, 0.0};
800 }
801
802 Array<OneD, Array<OneD, NekDouble>> qtmp{3};
803 for (int nd = 0; nd < 3; ++nd)
804 {
805 qtmp[nd] = NullNekDouble1DArray;
806 }
807 for (int nd2 = 0; nd2 < nDim; ++nd2)
808 {
809 qtmp[nd2] = elmt2ndDerv[nd2];
810 }
811
812 for (int i = 0; i < nTracePts; ++i)
813 {
814 tmp[i] = PenaltyFactor2 * m_traceNormDirctnElmtLength[i];
815 }
816 // the derivatives are assumed to be exchangable
817 for (int nd1 = 0; nd1 < nDim; ++nd1)
818 {
819 for (int i = 0; i < nConvectiveFields; ++i)
820 {
821 fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
822
823 for (int nd2 = nd1; nd2 < nDim; ++nd2)
824 {
825 Vmath::Zero(nTracePts, Bwd, 1);
826 fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd, true,
827 true, false);
828 for (int p = 0; p < nTracePts; ++p)
829 {
830 Bwd[p] *= tmp[p];
831 numDerivBwd[nd1][i][p] += m_traceNormals[nd2][p] * Bwd[p];
832 Fwd[p] *= tmp[p];
833 numDerivFwd[nd1][i][p] = -(m_traceNormals[nd2][p] * Fwd[p] -
834 numDerivFwd[nd1][i][p]);
835 }
836 if (nd2 != nd1)
837 {
838 for (int p = 0; p < nTracePts; ++p)
839 {
840 numDerivBwd[nd2][i][p] +=
841 m_traceNormals[nd1][p] * Bwd[p];
842 numDerivFwd[nd2][i][p] =
843 -(m_traceNormals[nd1][p] * Fwd[p] -
844 numDerivFwd[nd2][i][p]);
845 }
846 }
847 }
848 }
849 }
850}
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References m_traceNormals, m_traceNormDirctnElmtLength, Nektar::NullNekDouble1DArray, CellMLToNektar.cellml_metadata::p, 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 
)
private

Add symmetric flux integration to field (in coefficient space)

Definition at line 495 of file DiffusionIP.cpp.

502{
503 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
504 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
505 Array<OneD, Array<OneD, NekDouble>> tmpfield(nDim);
506 for (int i = 0; i < nDim; ++i)
507 {
508 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
509 }
510 int nv = 0;
511 for (int j = 0; j < nonZeroIndex.size(); ++j)
512 {
513 nv = nonZeroIndex[j];
514 MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
515 for (int nd = 0; nd < nDim; ++nd)
516 {
517 Vmath::Zero(nPts, tmpfield[nd], 1);
518
519 tracelist->MultiplyByQuadratureMetric(tracflux[nd][nv],
520 tracflux[nd][nv]);
521
522 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
523 tracflux[nd][nv], tmpfield[nd]);
524 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
525 }
526 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
527 Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
528 }
529}
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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.hpp:180

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

Referenced by 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 
)
private

Add symmetric flux integration to field (in physical space)

Definition at line 531 of file DiffusionIP.cpp.

538{
539 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
540 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
541 Array<OneD, NekDouble> tmpPhysi{nPts, 0.0};
542 Array<OneD, Array<OneD, NekDouble>> tmpfield{nDim};
543 for (int i = 0; i < nDim; ++i)
544 {
545 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
546 }
547 for (int j = 0; j < nonZeroIndex.size(); ++j)
548 {
549 int nv = nonZeroIndex[j];
550 for (int nd = 0; nd < nDim; ++nd)
551 {
552 Vmath::Zero(nPts, tmpfield[nd], 1);
553
554 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
555 tracflux[nd][nv], tmpfield[nd]);
556 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
557 }
558 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
559 fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
560 Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
561 }
562}

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

Referenced by AddDiffusionSymmFluxToPhys().

◆ ApplyFluxBndConds()

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

aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic

Definition at line 857 of file DiffusionIP.cpp.

861{
862 int nengy = nConvectiveFields - 1;
863 int cnt;
864 // Compute boundary conditions for Energy
865 cnt = 0;
866 size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
867 for (size_t j = 0; j < nBndRegions; ++j)
868 {
869 if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
871 {
872 continue;
873 }
874
875 size_t nBndEdges =
876 fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
877 for (size_t e = 0; e < nBndEdges; ++e)
878 {
879 size_t nBndEdgePts = fields[nengy]
880 ->GetBndCondExpansions()[j]
881 ->GetExp(e)
882 ->GetTotPoints();
883
884 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
885 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
886
887 if (fields[0]->GetBndConditions()[j]->GetUserDefined() ==
888 "WallAdiabatic")
889 {
890 Vmath::Zero(nBndEdgePts, &flux[nengy][id2], 1);
891 }
892 }
893 }
894}

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,
Array< OneD, int > &  nonZeroIndexflux,
TensorOfArray3D< NekDouble > &  traceflux,
Array< OneD, Array< OneD, NekDouble > > &  solution_Aver,
Array< OneD, Array< OneD, NekDouble > > &  solution_jump 
)
private

Calculate numerical flux on traces.

Definition at line 671 of file DiffusionIP.cpp.

681{
682 size_t nDim = fields[0]->GetCoordim(0);
683 size_t nPts = fields[0]->GetTotPoints();
684 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
685 size_t nConvectiveFields = fields.size();
686
688 fields[0]->GetTraceMap();
689
690 LibUtilities::Timer timer;
691 timer.Start();
692 // with further restructuring this iniziatilization could be eliminated
693 for (int nd = 0; nd < nDim; ++nd)
694 {
695 for (int i = 0; i < nConvectiveFields; ++i)
696 {
697 Vmath::Zero(nTracePts, m_wspNumDerivBwd[nd][i], 1);
698 Vmath::Zero(nTracePts, m_wspNumDerivFwd[nd][i], 1);
699 }
700 }
701
702 // could this be pre-allocated?
703 Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
704 Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
705
706 timer.Stop();
707 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux_alloc", 10);
708
709 timer.Start();
710 if (fabs(PenaltyFactor2) > 1.0E-12)
711 {
712 AddSecondDerivToTrace(nConvectiveFields, nDim, nPts, nTracePts,
713 PenaltyFactor2, fields, qfield, m_wspNumDerivFwd,
715 }
716
717 timer.Stop();
718 timer.AccumulateRegion("DiffIP:_AddSecondDerivToTrace", 10);
719
720 for (int nd = 0; nd < nDim; ++nd)
721 {
722 for (int i = 0; i < nConvectiveFields; ++i)
723 {
724 // this sequence of operations is really exepensive,
725 // it should be done collectively for all fields together instead of
726 // one by one
727 timer.Start();
728 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd, true, true,
729 false);
730 timer.Stop();
731 timer.AccumulateRegion("DiffIP:_GetFwdBwdTracePhys", 10);
732
733 for (size_t p = 0; p < nTracePts; ++p)
734 {
735 m_wspNumDerivBwd[nd][i][p] += 0.5 * Bwd[p];
736 m_wspNumDerivFwd[nd][i][p] += 0.5 * Fwd[p];
737 }
738
739 timer.Start();
740 TraceMap->GetAssemblyCommDG()->PerformExchange(
741 m_wspNumDerivFwd[nd][i], m_wspNumDerivBwd[nd][i]);
742 timer.Stop();
743 timer.AccumulateRegion("DiffIP:_PerformExchange", 10);
744
745 Vmath::Vadd(nTracePts, m_wspNumDerivFwd[nd][i], 1,
746 m_wspNumDerivBwd[nd][i], 1, m_wspNumDerivFwd[nd][i], 1);
747 }
748 }
749
750 timer.Start();
751 ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
752 solution_jump);
753 timer.Stop();
754 timer.AccumulateRegion("DiffIP:_ConsVarAveJump", 10);
755
756 Array<OneD, NekDouble> penaltyCoeff(nTracePts, 0.0);
757 GetPenaltyFactor(fields, penaltyCoeff);
758 for (size_t p = 0; p < nTracePts; ++p)
759 {
760 NekDouble PenaltyFactor =
761 penaltyCoeff[p] * m_traceNormDirctnElmtLengthRecip[p]; // load 1x
762
763 for (size_t f = 0; f < nConvectiveFields; ++f)
764 {
765 NekDouble jumpTmp = solution_jump[f][p] * PenaltyFactor; // load 1x
766 for (size_t d = 0; d < nDim; ++d)
767 {
768 NekDouble tmp = m_traceNormals[d][p] * jumpTmp +
769 m_wspNumDerivFwd[d][f][p]; // load 2x
770 m_wspNumDerivFwd[d][f][p] = tmp; // store 1x
771 }
772 }
773 }
774
775 timer.Start();
776 // Calculate normal viscous flux
778 traceflux, nonZeroIndexflux,
780 timer.Stop();
781 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxConsTrace", 10);
782}
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:350
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:121
void GetPenaltyFactor(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get IP penalty factor based on order.
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
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)
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:46
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddSecondDerivToTrace(), ConsVarAveJump(), Nektar::UnitTests::d(), GetPenaltyFactor(), Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxConsTrace, 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, Array< OneD, NekDouble > > &  solution_Aver,
Array< OneD, Array< OneD, NekDouble > > &  solution_jump,
Array< OneD, int > &  nonZeroIndexsymm,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  traceSymflux 
)
private

Calculate symmetric flux on traces.

Definition at line 472 of file DiffusionIP.cpp.

478{
479 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
480
481 m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
482 nonZeroIndexsymm, m_traceNormals);
483
484 for (int nd = 0; nd < nDim; ++nd)
485 {
486 for (int j = 0; j < nonZeroIndexsymm.size(); ++j)
487 {
488 int i = nonZeroIndexsymm[j];
489 Vmath::Smul(nTracePts, -0.5 * m_IPSymmFluxCoeff,
490 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
491 }
492 }
493}
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:352
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.hpp:100

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 
)
inlineprivate

Calculate the average of conservative variables on traces.

Definition at line 128 of file DiffusionIP.h.

131 {
132 constexpr int nvelst = 1;
133 int nEngy = nConvectiveFields - 1;
134 int nveled = nEngy;
135
136 T Fweight = 1.0 - Bweight;
137 for (size_t v = 0; v < nEngy; ++v)
138 {
139 aver[v] = Fweight * vFwd[v] + Bweight * vBwd[v];
140 }
141
142 T LinternalEngy{};
143 T RinternalEngy{};
144 T AinternalEngy{};
145 for (size_t j = nvelst; j < nveled; ++j)
146 {
147 LinternalEngy += vFwd[j] * vFwd[j];
148 RinternalEngy += vBwd[j] * vBwd[j];
149 AinternalEngy += aver[j] * aver[j];
150 }
151 LinternalEngy *= -0.5 / vFwd[0];
152 RinternalEngy *= -0.5 / vBwd[0];
153 LinternalEngy += vFwd[nEngy];
154 RinternalEngy += vBwd[nEngy];
155
156 aver[nEngy] = Fweight * LinternalEngy + Bweight * RinternalEngy;
157 aver[nEngy] += AinternalEngy * (0.5 / aver[0]);
158 }

Referenced by ConsVarAveJump().

◆ ConsVarAveJump()

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

Definition at line 622 of file DiffusionIP.cpp.

628{
629 std::vector<NekDouble> vFwdTmp(nConvectiveFields),
630 vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
631 for (size_t p = 0; p < nPts; ++p)
632 {
633 // re-arrange data
634 for (size_t f = 0; f < nConvectiveFields; ++f)
635 {
636 vFwdTmp[f] = vFwd[f][p];
637 vBwdTmp[f] = vBwd[f][p];
638 }
639
640 ConsVarAve(nConvectiveFields, m_tracBwdWeightAver[p], vFwdTmp, vBwdTmp,
641 averTmp);
642
643 // store
644 for (size_t f = 0; f < nConvectiveFields; ++f)
645 {
646 aver[f][p] = averTmp[f];
647 }
648 }
649
650 // if this can be make function of points, the nPts loop can be lifted more
651 m_SpecialBndTreat(aver);
652
653 // note: here the jump is 2.0*(aver-vFwd)
654 // because Viscous wall use a symmetry value as the Bwd,
655 // not the target one
656
657 for (size_t f = 0; f < nConvectiveFields; ++f)
658 {
659 for (size_t p = 0; p < nPts; ++p)
660 {
662 NekDouble Fweight = 2.0 - Bweight;
663
664 NekDouble tmpF = aver[f][p] - vFwd[f][p];
665 NekDouble tmpB = vBwd[f][p] - aver[f][p];
666 jump[f][p] = tmpF * Fweight + tmpB * Bweight;
667 }
668 }
669}
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:351
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:128
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.

Referenced by CalcTraceNumFlux().

◆ create()

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

Definition at line 46 of file DiffusionIP.h.

47 {
48 return DiffusionSharedPtr(new DiffusionIP());
49 }
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:54

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 
)
private

Calculate symmetric flux on traces interface.

Definition at line 455 of file DiffusionIP.cpp.

465{
466 size_t nDim = fields[0]->GetCoordim(0);
467
468 CalcTraceSymFlux(nConvectiveFields, nDim, m_traceAver, m_traceJump,
469 nonZeroIndex, SymmFlux);
470}
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:111
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:110
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, 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.

References CalcTraceSymFlux(), m_traceAver, and m_traceJump.

Referenced by AddDiffusionSymmFluxToCoeff(), and AddDiffusionSymmFluxToPhys().

◆ GetPenaltyFactor()

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

Get IP penalty factor based on order.

Definition at line 564 of file DiffusionIP.cpp.

567{
568 MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
569 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
570 tracelist->GetExp();
571 size_t ntotTrac = (*traceExp).size();
572 int nTracPnt, noffset;
573
574 const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
575 fields[0]->GetLocTraceToTraceMap();
576
577 const Array<OneD, const Array<OneD, int>> LRAdjExpid =
578 locTraceToTraceMap->GetLeftRightAdjacentExpId();
579 const Array<OneD, const Array<OneD, bool>> LRAdjflag =
580 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
581
582 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
583 fields[0]->GetExp();
584
585 Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
586
587 NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
588
589 int ntmp, numModes;
590
591 for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
592 {
593 noffset = tracelist->GetPhys_Offset(ntrace);
594 nTracPnt = tracelist->GetTotPoints(ntrace);
595
596 factorFwdBwd[0] = 0.0;
597 factorFwdBwd[1] = 0.0;
598
599 for (int nlr = 0; nlr < 2; ++nlr)
600 {
601 if (LRAdjflag[nlr][ntrace])
602 {
603 numModes = 0;
604 for (int nd = 0; nd < spaceDim; nd++)
605 {
606 ntmp = fields[0]
607 ->GetExp(LRAdjExpid[nlr][ntrace])
608 ->GetBasisNumModes(nd);
609 numModes = std::max(ntmp, numModes);
610 }
611 factorFwdBwd[nlr] = (numModes) * (numModes);
612 }
613 }
614
615 for (int np = 0; np < nTracPnt; ++np)
616 {
617 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
618 }
619 }
620}
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr

Referenced by CalcTraceNumFlux().

◆ 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 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Diffusion.

Definition at line 173 of file DiffusionIP.cpp.

180{
181
182 LibUtilities::Timer timer;
183 timer.Start();
184 DiffusionIP::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, m_wspDiff,
185 pFwd, pBwd);
186 for (int i = 0; i < nConvectiveFields; ++i)
187 {
188 fields[i]->BwdTrans(m_wspDiff[i], outarray[i]);
189 }
190 timer.Stop();
191 timer.AccumulateRegion("Diffusion IP");
192}
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) override
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 > &  qfields,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
overrideprotectedvirtual

Diffusion Flux, calculate the physical derivatives.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 319 of file DiffusionIP.cpp.

325{
326 Array<OneD, Array<OneD, NekDouble>> qtmp{3};
327 size_t nDim = fields[0]->GetCoordim(0);
328 for (int nd = 0; nd < 3; ++nd)
329 {
330 qtmp[nd] = NullNekDouble1DArray;
331 }
332
333 size_t nConvectiveFields = inarray.size();
334 for (int i = 0; i < nConvectiveFields; ++i)
335 {
336 for (int nd = 0; nd < nDim; ++nd)
337 {
338 qtmp[nd] = qfield[nd][i];
339 }
340 fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
341 }
342}

References Nektar::NullNekDouble1DArray.

◆ v_DiffuseCoeffs()

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 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 194 of file DiffusionIP.cpp.

201{
202 LibUtilities::Timer timer;
203 timer.Start();
204
205 size_t nDim = fields[0]->GetCoordim(0);
206 size_t nPts = fields[0]->GetTotPoints();
207 size_t nCoeffs = fields[0]->GetNcoeffs();
208 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
209
210 // pre-allocate this?
211 Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
212 Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
213 TensorOfArray3D<NekDouble> elmtFlux{nDim};
214 TensorOfArray3D<NekDouble> qfield{nDim};
215
216 for (int j = 0; j < nDim; ++j)
217 {
218 qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
219 elmtFlux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
220 for (int i = 0; i < nConvectiveFields; ++i)
221 {
222 qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
223 }
224 for (int i = 0; i < nConvectiveFields; ++i)
225 {
226 elmtFlux[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
227 }
228 }
229
230 // pre-allocate this?
231 Array<OneD, Array<OneD, NekDouble>> vFwd{nConvectiveFields};
232 Array<OneD, Array<OneD, NekDouble>> vBwd{nConvectiveFields};
233 // when does this happen?
235 {
236 for (int i = 0; i < nConvectiveFields; ++i)
237 {
238 vFwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
239 vBwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
240 }
241 for (int i = 0; i < nConvectiveFields; ++i)
242 {
243 fields[i]->GetFwdBwdTracePhys(inarray[i], vFwd[i], vBwd[i]);
244 }
245 }
246 else
247 {
248 for (int i = 0; i < nConvectiveFields; ++i)
249 {
250 vFwd[i] = pFwd[i];
251 vBwd[i] = pBwd[i];
252 }
253 }
254
255 DiffuseCalcDerivative(fields, inarray, qfield, vFwd, vBwd);
256 Array<OneD, int> nonZeroIndex;
257 DiffuseVolumeFlux(fields, inarray, qfield, elmtFlux, nonZeroIndex);
258 timer.Stop();
259 timer.AccumulateRegion("Diffusion:Volumeflux", 10);
260 timer.Start();
261
262 // pre-allocate this?
263 Array<OneD, Array<OneD, NekDouble>> tmpFluxIprdct{nDim};
264 // volume intergration: the nonZeroIndex indicates which flux is nonzero
265 for (int i = 0; i < nonZeroIndex.size(); ++i)
266 {
267 int j = nonZeroIndex[i];
268 for (int k = 0; k < nDim; ++k)
269 {
270 tmpFluxIprdct[k] = elmtFlux[k][j];
271 }
272 fields[j]->IProductWRTDerivBase(tmpFluxIprdct, outarray[j]);
273 Vmath::Neg(nCoeffs, outarray[j], 1);
274 }
275 timer.Stop();
276 timer.AccumulateRegion("Diffusion:IPWRTDB", 10);
277
278 // release qfield, elmtFlux and muvar;
279 timer.Start();
280 for (int j = 0; j < nDim; ++j)
281 {
282 elmtFlux[j] = NullNekDoubleArrayOfArray;
283 }
284
285 // pre-allocate this?
286 Array<OneD, Array<OneD, NekDouble>> Traceflux{nConvectiveFields};
287 for (int j = 0; j < nConvectiveFields; ++j)
288 {
289 Traceflux[j] = Array<OneD, NekDouble>{nTracePts, 0.0};
290 }
291
292 DiffuseTraceFlux(fields, inarray, qfield, elmtFlux, Traceflux, vFwd, vBwd,
293 nonZeroIndex);
294 timer.Stop();
295 timer.AccumulateRegion("Diffusion:TraceFlux", 10);
296
297 for (int i = 0; i < nonZeroIndex.size(); ++i)
298 {
299 int j = nonZeroIndex[i];
300
301 fields[j]->AddTraceIntegral(Traceflux[j], outarray[j]);
302 fields[j]->SetPhysState(false);
303 }
304
305 AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray, qfield,
306 elmtFlux, outarray, vFwd, vBwd);
307
308 timer.Start();
309 for (int i = 0; i < nonZeroIndex.size(); ++i)
310 {
311 int j = nonZeroIndex[i];
312
313 fields[j]->MultiplyByElmtInvMass(outarray[j], outarray[j]);
314 }
315 timer.Stop();
316 timer.AccumulateRegion("DiffIP:Diffusion Coeff", 10);
317}
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)
Definition: Diffusion.h:201
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:214
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:225
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)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddDiffusionSymmFluxToCoeff(), Nektar::SolverUtils::Diffusion::DiffuseCalcDerivative(), Nektar::SolverUtils::Diffusion::DiffuseTraceFlux(), Nektar::SolverUtils::Diffusion::DiffuseVolumeFlux(), Vmath::Neg(), Nektar::NullNekDoubleArrayOfArray, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by v_Diffuse().

◆ v_DiffuseTraceFlux()

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 
)
overrideprotectedvirtual

Diffusion term Trace Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 362 of file DiffusionIP.cpp.

371{
372 TensorOfArray3D<NekDouble> traceflux3D(1);
373 traceflux3D[0] = TraceFlux;
374
375 size_t nConvectiveFields = fields.size();
376 LibUtilities::Timer timer;
377 timer.Start();
378 CalcTraceNumFlux(m_IP2ndDervCoeff, fields, inarray, qfield, pFwd, pBwd,
379 nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
380 timer.Stop();
381 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux", 10);
382
383 ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
384}
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 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, 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_traceAver, m_traceJump, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

◆ 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 
)
overrideprotectedvirtual

Diffusion Volume Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 344 of file DiffusionIP.cpp.

349{
350 size_t nDim = fields[0]->GetCoordim(0);
351
352 Array<OneD, Array<OneD, NekDouble>> tmparray2D = NullNekDoubleArrayOfArray;
353
354 LibUtilities::Timer timer;
355 timer.Start();
356 m_FunctorDiffusionfluxCons(nDim, inarray, qfield, VolumeFlux, nonZeroIndex,
357 tmparray2D);
358 timer.Stop();
359 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxCons", 10);
360}
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:349

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

◆ v_GetTraceNormal()

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

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 99 of file DiffusionIP.h.

100 {
101 return m_traceNormals;
102 }

References m_traceNormals.

◆ v_InitObject()

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

Implements Nektar::SolverUtils::Diffusion.

Definition at line 52 of file DiffusionIP.cpp.

55{
56 m_session = pSession;
57
58 m_session->LoadParameter("IPSymmFluxCoeff", m_IPSymmFluxCoeff,
59 0.0); // SIP=1.0; NIP=-1.0; IIP=0.0
60
61 m_session->LoadParameter("IP2ndDervCoeff", m_IP2ndDervCoeff,
62 0.0); // 1.0/12.0
63
64 m_session->LoadParameter("IPPenaltyCoeff", m_IPPenaltyCoeff,
65 4.0); // 1.0/12.0
66
67 // Setting up the normals
68 size_t nDim = pFields[0]->GetCoordim(0);
69 size_t nVariable = pFields.size();
70 size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
71
72 m_traceNormals = Array<OneD, Array<OneD, NekDouble>>{nDim};
73 for (int i = 0; i < nDim; ++i)
74 {
75 m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
76 }
77 m_traceAver = Array<OneD, Array<OneD, NekDouble>>{nVariable};
78 m_traceJump = Array<OneD, Array<OneD, NekDouble>>{nVariable};
79 for (int i = 0; i < nVariable; ++i)
80 {
81 m_traceAver[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
82 m_traceJump[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
83 }
84
85 pFields[0]->GetTrace()->GetNormals(m_traceNormals);
86
87 // Compute the length of the elements in the boundary-normal direction
88 // The function "GetElmtNormalLength" returns half the length of the left
89 // and right adjacent element in the lengthFwd and lengthBwd arrays. If
90 // the element belongs to a boundary (including periodic boundaries) or
91 // a parallel interface, the Bwd array will contain zeros.
92 Array<OneD, NekDouble> lengthFwd{nTracePts, 0.0};
93 Array<OneD, NekDouble> lengthBwd{nTracePts, 0.0};
94 pFields[0]->GetTrace()->GetElmtNormalLength(lengthFwd, lengthBwd);
95
96 // Copy Fwd to Bwd on parallel interfaces
97 // TODO: Move this into GetElmtNormalLength()
98 pFields[0]->GetTraceMap()->GetAssemblyCommDG()->PerformExchange(lengthFwd,
99 lengthBwd);
100 // Copy Fwd to Bwd on periodic interfaces
101 // TODO: Move this into GetElmtNormalLength()
102 pFields[0]->PeriodicBwdCopy(lengthFwd, lengthBwd);
103 // Scale the length by 0.5 on boundaries, and copy Fwd into Bwd
104 // Notes:
105 // - It is not quite clear why we need to scale by 0.5 on the boundaries
106 // - The current implementation is not perfect, it would be nicer to call
107 // a function similar to DiscontField::v_FillBwdWithBoundCond() with
108 // PutFwdInBwdOnBCs = true. If we wouldn't do the scaling by 0.5, this
109 // function could have been used.
110 for (int i = 0; i < nTracePts; ++i)
111 {
112 if (std::abs(lengthBwd[i]) < NekConstants::kNekMachineEpsilon)
113 {
114 lengthFwd[i] *= 0.5;
115 lengthBwd[i] = lengthFwd[i];
116 }
117 }
118
119 // Compute the average element normal length along the edge
120 Array<OneD, NekDouble> lengthsum(nTracePts, 0.0);
121 Array<OneD, NekDouble> lengthmul(nTracePts, 0.0);
122 Vmath::Vadd(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthsum, 1);
123 Vmath::Vmul(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthmul, 1);
124 Vmath::Vdiv(nTracePts, lengthsum, 1, lengthmul, 1, lengthFwd, 1);
125 m_traceNormDirctnElmtLength = lengthsum;
129
130 m_tracBwdWeightAver = Array<OneD, NekDouble>{nTracePts, 0.0};
131 m_tracBwdWeightJump = Array<OneD, NekDouble>{nTracePts, 0.0};
132 pFields[0]->GetBwdWeight(m_tracBwdWeightAver, m_tracBwdWeightJump);
133 Array<OneD, NekDouble> tmpBwdWeight{nTracePts, 0.0};
134 Array<OneD, NekDouble> tmpBwdWeightJump{nTracePts, 0.0};
135 for (int i = 1; i < nVariable; ++i)
136 {
137 pFields[i]->GetBwdWeight(tmpBwdWeight, tmpBwdWeightJump);
138 Vmath::Vsub(nTracePts, tmpBwdWeight, 1, m_tracBwdWeightAver, 1,
139 tmpBwdWeight, 1);
140 Vmath::Vabs(nTracePts, tmpBwdWeight, 1, tmpBwdWeight, 1);
141 NekDouble norm = 0.0;
142 for (int j = 0; j < nTracePts; ++j)
143 {
144 norm += tmpBwdWeight[j];
145 }
146 ASSERTL0(norm < 1.0E-11,
147 "different BWD for different variable not coded yet");
148 }
149
150 // workspace for v_diffuse
151 size_t nCoeffs = pFields[0]->GetNcoeffs();
152 m_wspDiff = Array<OneD, Array<OneD, NekDouble>>{nVariable};
153 for (int i = 0; i < nVariable; ++i)
154 {
155 m_wspDiff[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
156 }
157
158 // workspace for callnumtraceflux
159 m_wspNumDerivBwd = TensorOfArray3D<NekDouble>{nDim};
160 m_wspNumDerivFwd = TensorOfArray3D<NekDouble>{nDim};
161 for (int nd = 0; nd < nDim; ++nd)
162 {
163 m_wspNumDerivBwd[nd] = Array<OneD, Array<OneD, NekDouble>>{nVariable};
164 m_wspNumDerivFwd[nd] = Array<OneD, Array<OneD, NekDouble>>{nVariable};
165 for (int i = 0; i < nVariable; ++i)
166 {
167 m_wspNumDerivBwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
168 m_wspNumDerivFwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
169 }
170 }
171}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:122
static const NekDouble kNekMachineEpsilon
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.hpp:72
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
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.
Definition: Vmath.hpp:126
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.hpp:220
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298

References tinysimd::abs(), ASSERTL0, Nektar::NekConstants::kNekMachineEpsilon, m_IP2ndDervCoeff, m_IPPenaltyCoeff, m_IPSymmFluxCoeff, m_session, m_tracBwdWeightAver, m_tracBwdWeightJump, m_traceAver, m_traceJump, m_traceNormals, m_traceNormDirctnElmtLength, m_traceNormDirctnElmtLengthRecip, m_wspDiff, m_wspNumDerivBwd, m_wspNumDerivFwd, Vmath::Smul(), Vmath::Vabs(), Vmath::Vadd(), Vmath::Vdiv(), Vmath::Vmul(), and Vmath::Vsub().

Member Data Documentation

◆ m_IP2ndDervCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IP2ndDervCoeff
private

Definition at line 106 of file DiffusionIP.h.

Referenced by v_DiffuseTraceFlux(), and v_InitObject().

◆ m_IPPenaltyCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IPPenaltyCoeff
private

Definition at line 107 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_IPSymmFluxCoeff

NekDouble Nektar::SolverUtils::DiffusionIP::m_IPSymmFluxCoeff
private

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionIP::m_session
private

Definition at line 122 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_tracBwdWeightAver

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

Definition at line 118 of file DiffusionIP.h.

Referenced by ConsVarAveJump(), and v_InitObject().

◆ m_tracBwdWeightJump

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

Definition at line 119 of file DiffusionIP.h.

Referenced by ConsVarAveJump(), and v_InitObject().

◆ m_traceAver

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

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
private

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
private

◆ m_traceNormDirctnElmtLength

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

Definition at line 120 of file DiffusionIP.h.

Referenced by AddSecondDerivToTrace(), and v_InitObject().

◆ m_traceNormDirctnElmtLengthRecip

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

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
private

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
private

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
private

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, "Interior Penalty")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:46
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39

Definition at line 51 of file DiffusionIP.h.