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 ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 
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) override
 
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) override
 
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) override
 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) override
 Diffusion term Trace Flux. More...
 
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) override
 Diffusion Flux, calculate the physical derivatives. More...
 
virtual 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 45 of file DiffusionIP.h.

Constructor & Destructor Documentation

◆ DiffusionIP()

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

Definition at line 50 of file DiffusionIP.cpp.

51{
52}

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 391 of file DiffusionIP.cpp.

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

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

801{
802 Array<OneD, NekDouble> Fwd(nTracePts, 0.0);
803 Array<OneD, NekDouble> Bwd(nTracePts, 0.0);
804 std::vector<NekDouble> tmp(nTracePts);
805
806 Array<OneD, Array<OneD, NekDouble>> elmt2ndDerv{nDim};
807 for (int nd1 = 0; nd1 < nDim; ++nd1)
808 {
809 elmt2ndDerv[nd1] = Array<OneD, NekDouble>{nPts, 0.0};
810 }
811
812 Array<OneD, Array<OneD, NekDouble>> qtmp{3};
813 for (int nd = 0; nd < 3; ++nd)
814 {
815 qtmp[nd] = NullNekDouble1DArray;
816 }
817 for (int nd2 = 0; nd2 < nDim; ++nd2)
818 {
819 qtmp[nd2] = elmt2ndDerv[nd2];
820 }
821
822 for (int i = 0; i < nTracePts; ++i)
823 {
824 tmp[i] = PenaltyFactor2 * m_traceNormDirctnElmtLength[i];
825 }
826 // the derivatives are assumed to be exchangable
827 for (int nd1 = 0; nd1 < nDim; ++nd1)
828 {
829 for (int i = 0; i < nConvectiveFields; ++i)
830 {
831 fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
832
833 for (int nd2 = nd1; nd2 < nDim; ++nd2)
834 {
835 Vmath::Zero(nTracePts, Bwd, 1);
836 fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd, true,
837 true, false);
838 for (int p = 0; p < nTracePts; ++p)
839 {
840 Bwd[p] *= tmp[p];
841 numDerivBwd[nd1][i][p] += m_traceNormals[nd2][p] * Bwd[p];
842 Fwd[p] *= tmp[p];
843 numDerivFwd[nd1][i][p] = -(m_traceNormals[nd2][p] * Fwd[p] -
844 numDerivFwd[nd1][i][p]);
845 }
846 if (nd2 != nd1)
847 {
848 for (int p = 0; p < nTracePts; ++p)
849 {
850 numDerivBwd[nd2][i][p] +=
851 m_traceNormals[nd1][p] * Bwd[p];
852 numDerivFwd[nd2][i][p] =
853 -(m_traceNormals[nd1][p] * Fwd[p] -
854 numDerivFwd[nd2][i][p]);
855 }
856 }
857 }
858 }
859 }
860}
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:125
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:114
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

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 500 of file DiffusionIP.cpp.

507{
508 boost::ignore_unused(nTracePts);
509
510 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
511 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
512 Array<OneD, Array<OneD, NekDouble>> tmpfield(nDim);
513 for (int i = 0; i < nDim; ++i)
514 {
515 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
516 }
517 int nv = 0;
518 for (int j = 0; j < nonZeroIndex.size(); ++j)
519 {
520 nv = nonZeroIndex[j];
521 MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
522 for (int nd = 0; nd < nDim; ++nd)
523 {
524 Vmath::Zero(nPts, tmpfield[nd], 1);
525
526 tracelist->MultiplyByQuadratureMetric(tracflux[nd][nv],
527 tracflux[nd][nv]);
528
529 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
530 tracflux[nd][nv], tmpfield[nd]);
531 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
532 }
533 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
534 Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
535 }
536}
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.cpp:354

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 538 of file DiffusionIP.cpp.

545{
546 boost::ignore_unused(nTracePts);
547
548 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
549 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
550 Array<OneD, NekDouble> tmpPhysi{nPts, 0.0};
551 Array<OneD, Array<OneD, NekDouble>> tmpfield{nDim};
552 for (int i = 0; i < nDim; ++i)
553 {
554 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
555 }
556 for (int j = 0; j < nonZeroIndex.size(); ++j)
557 {
558 int nv = nonZeroIndex[j];
559 for (int nd = 0; nd < nDim; ++nd)
560 {
561 Vmath::Zero(nPts, tmpfield[nd], 1);
562
563 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
564 tracflux[nd][nv], tmpfield[nd]);
565 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
566 }
567 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
568 fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
569 Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
570 }
571}

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 867 of file DiffusionIP.cpp.

871{
872 int nengy = nConvectiveFields - 1;
873 int cnt;
874 // Compute boundary conditions for Energy
875 cnt = 0;
876 size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
877 for (size_t j = 0; j < nBndRegions; ++j)
878 {
879 if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
881 {
882 continue;
883 }
884
885 size_t nBndEdges =
886 fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
887 for (size_t e = 0; e < nBndEdges; ++e)
888 {
889 size_t nBndEdgePts = fields[nengy]
890 ->GetBndCondExpansions()[j]
891 ->GetExp(e)
892 ->GetTotPoints();
893
894 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
895 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
896
897 if (fields[0]->GetBndConditions()[j]->GetUserDefined() ==
898 "WallAdiabatic")
899 {
900 Vmath::Zero(nBndEdgePts, &flux[nengy][id2], 1);
901 }
902 }
903 }
904}

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 680 of file DiffusionIP.cpp.

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

483{
484 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
485
486 m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
487 nonZeroIndexsymm, m_traceNormals);
488
489 for (int nd = 0; nd < nDim; ++nd)
490 {
491 for (int j = 0; j < nonZeroIndexsymm.size(); ++j)
492 {
493 int i = nonZeroIndexsymm[j];
494 Vmath::Smul(nTracePts, -0.5 * m_IPSymmFluxCoeff,
495 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
496 }
497 }
498}
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:356
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:245

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 133 of file DiffusionIP.h.

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

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 631 of file DiffusionIP.cpp.

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

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 48 of file DiffusionIP.h.

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

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 460 of file DiffusionIP.cpp.

469{
470 boost::ignore_unused(inarray, qfield, VolumeFlux, pFwd, pBwd);
471 size_t nDim = fields[0]->GetCoordim(0);
472
473 CalcTraceSymFlux(nConvectiveFields, nDim, m_traceAver, m_traceJump,
474 nonZeroIndex, SymmFlux);
475}
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:116
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:115
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 573 of file DiffusionIP.cpp.

576{
577 MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
578 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
579 tracelist->GetExp();
580 size_t ntotTrac = (*traceExp).size();
581 int nTracPnt, noffset;
582
583 const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
584 fields[0]->GetLocTraceToTraceMap();
585
586 const Array<OneD, const Array<OneD, int>> LRAdjExpid =
587 locTraceToTraceMap->GetLeftRightAdjacentExpId();
588 const Array<OneD, const Array<OneD, bool>> LRAdjflag =
589 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
590
591 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
592 fields[0]->GetExp();
593
594 Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
595
596 NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
597
598 int ntmp, numModes;
599
600 for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
601 {
602 noffset = tracelist->GetPhys_Offset(ntrace);
603 nTracPnt = tracelist->GetTotPoints(ntrace);
604
605 factorFwdBwd[0] = 0.0;
606 factorFwdBwd[1] = 0.0;
607
608 for (int nlr = 0; nlr < 2; ++nlr)
609 {
610 if (LRAdjflag[nlr][ntrace])
611 {
612 numModes = 0;
613 for (int nd = 0; nd < spaceDim; nd++)
614 {
615 ntmp = fields[0]
616 ->GetExp(LRAdjExpid[nlr][ntrace])
617 ->GetBasisNumModes(nd);
618 numModes = std::max(ntmp, numModes);
619 }
620 factorFwdBwd[nlr] = (numModes) * (numModes);
621 }
622 }
623
624 for (int np = 0; np < nTracPnt; ++np)
625 {
626 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
627 }
628 }
629}
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 175 of file DiffusionIP.cpp.

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

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 321 of file DiffusionIP.cpp.

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

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 196 of file DiffusionIP.cpp.

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

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 366 of file DiffusionIP.cpp.

374{
375 boost::ignore_unused(VolumeFlux);
376
377 TensorOfArray3D<NekDouble> traceflux3D(1);
378 traceflux3D[0] = TraceFlux;
379
380 size_t nConvectiveFields = fields.size();
381 LibUtilities::Timer timer;
382 timer.Start();
383 CalcTraceNumFlux(m_IP2ndDervCoeff, fields, inarray, qfield, pFwd, pBwd,
384 nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
385 timer.Stop();
386 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux", 10);
387
388 ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
389}
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 348 of file DiffusionIP.cpp.

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

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxCons, 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 ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 103 of file DiffusionIP.h.

105 {
106 return m_traceNormals;
107 }

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 54 of file DiffusionIP.cpp.

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

Referenced by v_DiffuseTraceFlux(), and v_InitObject().

◆ m_IPPenaltyCoeff

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

Definition at line 112 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 127 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_tracBwdWeightAver

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

Definition at line 123 of file DiffusionIP.h.

Referenced by ConsVarAveJump(), and v_InitObject().

◆ m_tracBwdWeightJump

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

Definition at line 124 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 115 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 116 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 125 of file DiffusionIP.h.

Referenced by AddSecondDerivToTrace(), and v_InitObject().

◆ m_traceNormDirctnElmtLengthRecip

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

Definition at line 126 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 118 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 120 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ m_wspNumDerivFwd

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

Definition at line 121 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:198
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:48
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

Definition at line 54 of file DiffusionIP.h.