Nektar++
Loading...
Searching...
No Matches
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.
 
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.
 
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.
 
const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal () override
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()=default
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 

Private Member Functions

template<class T , typename = typename std::enable_if< std::is_floating_point_v<T> || tinysimd::is_vector_floating_point_v<T>>::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.
 
void GetPenaltyFactor (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get IP penalty factor based on order.
 
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 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 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.
 
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.
 
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.
 
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)
 
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
 

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.
 
TensorOfArray3D< NekDoublem_wspNumDerivBwd
 Workspace for CallTraceNumFlux.
 
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
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.
 
SOLVER_UTILS_EXPORT void DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
 Diffusion term Trace Flux.
 
SOLVER_UTILS_EXPORT void 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.
 
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)
 
SOLVER_UTILS_EXPORT void SetGridVelocityTrace (Array< OneD, Array< OneD, NekDouble > > &gridVelocityTrace)
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor.
 
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
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
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 49 of file DiffusionIP.cpp.

50{
51}

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

413{
414 if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
415 {
416 size_t nDim = fields[0]->GetCoordim(0);
417 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
418 TensorOfArray3D<NekDouble> traceSymflux{nDim};
419 for (int nd = 0; nd < nDim; ++nd)
420 {
421 traceSymflux[nd] =
422 Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
423 for (int j = 0; j < nConvectiveFields; ++j)
424 {
425 traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
426 }
427 }
428 Array<OneD, int> nonZeroIndex;
429 DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
430 VolumeFlux, traceSymflux, pFwd, pBwd,
431 nonZeroIndex);
432
433 AddSymmFluxIntegralToCoeff(nConvectiveFields, nDim, fields,
434 nonZeroIndex, traceSymflux, outarray);
435 }
436}
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, 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

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

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

Referenced by CalcTraceNumFlux().

◆ AddSymmFluxIntegralToCoeff()

void Nektar::SolverUtils::DiffusionIP::AddSymmFluxIntegralToCoeff ( const std::size_t  nConvectiveFields,
const size_t  nDim,
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 478 of file DiffusionIP.cpp.

484{
485 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
486 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
487 Array<OneD, Array<OneD, NekDouble>> locTracesDotNorm(nDim);
488 int nv = 0;
489 for (int j = 0; j < nonZeroIndex.size(); ++j)
490 {
491 nv = nonZeroIndex[j];
492 auto tracelist = fields[nv]->GetTrace();
493 auto locTraceToTraceMap = fields[nv]->GetLocTraceToTraceMap();
494 auto DisContField =
495 std::dynamic_pointer_cast<MultiRegions::DisContField>(fields[nv]);
496
497 // Use p2p interp
498 int nLocTracePts = DisContField->GetLocElmtTrace()->GetTotPoints();
499 int nTotPhys = DisContField->GetTotPoints();
500
501 // For locTrace, Fwd and Bwd are clustered, so we need to get start
502 // of Bwd
503 Array<OneD, NekDouble> locTraceBwd;
504 Array<OneD, NekDouble> wsp(nLocTracePts);
505 Array<OneD, NekDouble> wsp2(nLocTracePts);
506 Array<OneD, Array<OneD, NekDouble>> tmpfields(nDim);
507
508 for (int nd = 0; nd < nDim; ++nd)
509 {
510 locTracesDotNorm[nd] = Array<OneD, NekDouble>(nLocTracePts, 0.0);
511 // interp trace to locTrace and unshuffle to loc elmt trace
512 locTraceToTraceMap->InterpTraceToLocTrace(0, tracflux[nd][nv],
513 wsp2);
514 locTraceToTraceMap->UnshuffleLocTraces(0, wsp2, wsp);
515 // interp trace to locTrace and unshuffle to loc elmt trace
516 locTraceBwd = wsp2 + locTraceToTraceMap->GetNFwdLocTracePts();
517 locTraceToTraceMap->InterpTraceToLocTrace(1, tracflux[nd][nv],
518 locTraceBwd);
519 locTraceToTraceMap->UnshuffleLocTraces(1, locTraceBwd, wsp);
520 // multiply by loc wJ
521 DisContField->GetLocElmtTrace()->MultiplyByQuadratureMetric(wsp,
522 wsp);
523 // Reshuffle LocElmtTrace to LocTraces
524 locTraceToTraceMap->ReshuffleLocTracesForInterp(
525 0, wsp, locTracesDotNorm[nd]);
526 locTraceBwd =
527 locTracesDotNorm[nd] + locTraceToTraceMap->GetNFwdLocTracePts();
528 locTraceToTraceMap->ReshuffleLocTracesForInterp(1, wsp,
529 locTraceBwd);
530 // add locTraceFlux back to tmpfields
531 tmpfields[nd] = Array<OneD, NekDouble>(nTotPhys, 0.0);
532 locTraceToTraceMap->AddLocTracesToField(locTracesDotNorm[nd],
533 tmpfields[nd]);
534 // To use GLwithEnds here, we must handle division by zero
535 DisContField->DivideByQuadratureMetric(tmpfields[nd],
536 tmpfields[nd]);
537 }
538
539 DisContField->IProductWRTDerivBase(tmpfields, tmpCoeff);
540 Vmath::Vadd(nCoeffs, outarray[nv], 1, tmpCoeff, 1, outarray[nv], 1);
541 }
542}
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 Nektar::MultiRegions::ExpList::DivideByQuadratureMetric(), Nektar::MultiRegions::DisContField::GetLocElmtTrace(), Nektar::MultiRegions::ExpList::GetTotPoints(), Nektar::MultiRegions::ExpList::IProductWRTDerivBase(), and Vmath::Vadd().

Referenced by AddDiffusionSymmFluxToCoeff().

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

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

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

661{
662 size_t nDim = fields[0]->GetCoordim(0);
663 size_t nPts = fields[0]->GetTotPoints();
664 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
665 size_t nConvectiveFields = fields.size();
666
668 fields[0]->GetTraceMap();
669 const MultiRegions::InterfaceMapDGSharedPtr InterfaceMap =
670 fields[0]->GetInterfaceMap();
671
672 LibUtilities::Timer timer;
673 timer.Start();
674 // with further restructuring this iniziatilization could be eliminated
675 for (int nd = 0; nd < nDim; ++nd)
676 {
677 for (int i = 0; i < nConvectiveFields; ++i)
678 {
679 Vmath::Zero(nTracePts, m_wspNumDerivBwd[nd][i], 1);
680 Vmath::Zero(nTracePts, m_wspNumDerivFwd[nd][i], 1);
681 }
682 }
683
684 // could this be pre-allocated?
685 Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
686 Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
687
688 timer.Stop();
689 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux_alloc", 10);
690
691 timer.Start();
692 if (fabs(PenaltyFactor2) > 1.0E-12)
693 {
694 AddSecondDerivToTrace(nConvectiveFields, nDim, nPts, nTracePts,
695 PenaltyFactor2, fields, qfield, m_wspNumDerivFwd,
697 }
698
699 timer.Stop();
700 timer.AccumulateRegion("DiffIP:_AddSecondDerivToTrace", 10);
701 for (int nd = 0; nd < nDim; ++nd)
702 {
703 for (int i = 0; i < nConvectiveFields; ++i)
704 {
705 // this sequence of operations is really exepensive,
706 // it should be done collectively for all fields together instead of
707 // one by one
708 timer.Start();
709 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd, true, true,
710 false);
711 timer.Stop();
712 timer.AccumulateRegion("DiffIP:_GetFwdBwdTracePhys", 10);
713
714 for (size_t p = 0; p < nTracePts; ++p)
715 {
716 m_wspNumDerivBwd[nd][i][p] += 0.5 * Bwd[p];
717 m_wspNumDerivFwd[nd][i][p] += 0.5 * Fwd[p];
718 }
719
720 timer.Start();
721 TraceMap->GetAssemblyCommDG()->PerformExchange(
722 m_wspNumDerivFwd[nd][i], m_wspNumDerivBwd[nd][i]);
723 InterfaceMap->ExchangeTrace(m_wspNumDerivFwd[nd][i],
724 m_wspNumDerivBwd[nd][i]);
725 timer.Stop();
726 timer.AccumulateRegion("DiffIP:_PerformExchange", 10);
727 }
728 }
729
730 // Rotate for rotated periodic boundary conditions
731 fields[0]->PeriodicDeriveBwdRot(m_wspNumDerivBwd);
732 // Rotate for sector rotation along interface
733 if (fields[0]->GetGraph()->GetMovement()->GetSectorRotateFlag())
734 {
735 fields[0]->RotLocalBwdDeriveTrace(m_wspNumDerivBwd);
736 }
737
738 for (int nd = 0; nd < nDim; ++nd)
739 {
740
741 for (int i = 0; i < nConvectiveFields; ++i)
742 {
743 Vmath::Vadd(nTracePts, m_wspNumDerivFwd[nd][i], 1,
744 m_wspNumDerivBwd[nd][i], 1, m_wspNumDerivFwd[nd][i], 1);
745 }
746 }
747
748 timer.Start();
749 ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
750 solution_jump);
751 timer.Stop();
752 timer.AccumulateRegion("DiffIP:_ConsVarAveJump", 10);
753
754 Array<OneD, NekDouble> penaltyCoeff(nTracePts, 0.0);
755 GetPenaltyFactor(fields, penaltyCoeff);
756 for (size_t p = 0; p < nTracePts; ++p)
757 {
758 NekDouble PenaltyFactor =
759 penaltyCoeff[p] * m_traceNormDirctnElmtLengthRecip[p]; // load 1x
760
761 for (size_t f = 0; f < nConvectiveFields; ++f)
762 {
763 NekDouble jumpTmp = solution_jump[f][p] * PenaltyFactor; // load 1x
764 for (size_t d = 0; d < nDim; ++d)
765 {
766 NekDouble tmp = m_traceNormals[d][p] * jumpTmp +
767 m_wspNumDerivFwd[d][f][p]; // load 2x
768 m_wspNumDerivFwd[d][f][p] = tmp; // store 1x
769 }
770 }
771 }
772
773 timer.Start();
774 // Calculate normal viscous flux
776 traceflux, nonZeroIndexflux,
778 timer.Stop();
779 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxConsTrace", 10);
780}
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition Diffusion.h:364
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
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.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
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
std::shared_ptr< InterfaceMapDG > InterfaceMapDGSharedPtr
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddSecondDerivToTrace(), ConsVarAveJump(), GetPenaltyFactor(), Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxConsTrace, m_traceNormals, m_traceNormDirctnElmtLengthRecip, m_wspNumDerivBwd, m_wspNumDerivFwd, 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 455 of file DiffusionIP.cpp.

461{
462 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
463
464 m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
465 nonZeroIndexsymm, m_traceNormals);
466
467 for (int nd = 0; nd < nDim; ++nd)
468 {
469 for (int j = 0; j < nonZeroIndexsymm.size(); ++j)
470 {
471 int i = nonZeroIndexsymm[j];
472 Vmath::Smul(nTracePts, -0.5 * m_IPSymmFluxCoeff,
473 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
474 }
475 }
476}
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition Diffusion.h:366
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_v<T> || tinysimd::is_vector_floating_point_v<T>>::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 602 of file DiffusionIP.cpp.

608{
609 std::vector<NekDouble> vFwdTmp(nConvectiveFields),
610 vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
611 for (size_t p = 0; p < nPts; ++p)
612 {
613 // re-arrange data
614 for (size_t f = 0; f < nConvectiveFields; ++f)
615 {
616 vFwdTmp[f] = vFwd[f][p];
617 vBwdTmp[f] = vBwd[f][p];
618 }
619
620 ConsVarAve(nConvectiveFields, m_tracBwdWeightAver[p], vFwdTmp, vBwdTmp,
621 averTmp);
622
623 // store
624 for (size_t f = 0; f < nConvectiveFields; ++f)
625 {
626 aver[f][p] = averTmp[f];
627 }
628 }
629
630 // if this can be make function of points, the nPts loop can be lifted more
631 m_SpecialBndTreat(aver);
632
633 // note: here the jump is 2.0*(aver-vFwd)
634 // because Viscous wall use a symmetry value as the Bwd,
635 // not the target one
636
637 for (size_t f = 0; f < nConvectiveFields; ++f)
638 {
639 for (size_t p = 0; p < nPts; ++p)
640 {
642 NekDouble Fweight = 2.0 - Bweight;
643
644 NekDouble tmpF = aver[f][p] - vFwd[f][p];
645 NekDouble tmpB = vBwd[f][p] - aver[f][p];
646 jump[f][p] = tmpF * Fweight + tmpB * Bweight;
647 }
648 }
649}
SpecialBndTreat m_SpecialBndTreat
Definition Diffusion.h:365
Array< OneD, NekDouble > m_tracBwdWeightAver
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.
Array< OneD, NekDouble > m_tracBwdWeightJump

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

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:55

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

448{
449 size_t nDim = fields[0]->GetCoordim(0);
450
451 CalcTraceSymFlux(nConvectiveFields, nDim, m_traceAver, m_traceJump,
452 nonZeroIndex, SymmFlux);
453}
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Array< OneD, Array< OneD, NekDouble > > m_traceAver
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().

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

547{
548 MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
549 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
550 tracelist->GetExp();
551 size_t ntotTrac = (*traceExp).size();
552 int nTracPnt, noffset;
553
554 const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
555 fields[0]->GetLocTraceToTraceMap();
556
557 const Array<OneD, const Array<OneD, int>> LRAdjExpid =
558 locTraceToTraceMap->GetLeftRightAdjacentExpId();
559 const Array<OneD, const Array<OneD, bool>> LRAdjflag =
560 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
561
562 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
563 fields[0]->GetExp();
564
565 Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
566
567 NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
568
569 int ntmp, numModes;
570
571 for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
572 {
573 noffset = tracelist->GetPhys_Offset(ntrace);
574 nTracPnt = tracelist->GetTotPoints(ntrace);
575
576 factorFwdBwd[0] = 0.0;
577 factorFwdBwd[1] = 0.0;
578
579 for (int nlr = 0; nlr < 2; ++nlr)
580 {
581 if (LRAdjflag[nlr][ntrace])
582 {
583 numModes = 0;
584 for (int nd = 0; nd < spaceDim; nd++)
585 {
586 ntmp = fields[0]
587 ->GetExp(LRAdjExpid[nlr][ntrace])
588 ->GetBasisNumModes(nd);
589 numModes = std::max(ntmp, numModes);
590 }
591 factorFwdBwd[nlr] = (numModes) * (numModes);
592 }
593 }
594
595 for (int np = 0; np < nTracPnt; ++np)
596 {
597 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
598 }
599 }
600}
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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 174 of file DiffusionIP.cpp.

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

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

344{
345 Array<OneD, Array<OneD, NekDouble>> qtmp{3};
346 size_t nDim = fields[0]->GetCoordim(0);
347 for (int nd = 0; nd < 3; ++nd)
348 {
349 qtmp[nd] = NullNekDouble1DArray;
350 }
351
352 size_t nConvectiveFields = inarray.size();
353 for (int i = 0; i < nConvectiveFields; ++i)
354 {
355 for (int nd = 0; nd < nDim; ++nd)
356 {
357 qtmp[nd] = qfield[nd][i];
358 }
359 fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
360 }
361}

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

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

390{
391 TensorOfArray3D<NekDouble> traceflux3D(1);
392 traceflux3D[0] = TraceFlux;
393
394 size_t nConvectiveFields = fields.size();
395 LibUtilities::Timer timer;
396 timer.Start();
397 CalcTraceNumFlux(m_IP2ndDervCoeff, fields, inarray, qfield, pFwd, pBwd,
398 nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
399 timer.Stop();
400 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux", 10);
401
402 ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
403}
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 363 of file DiffusionIP.cpp.

368{
369 size_t nDim = fields[0]->GetCoordim(0);
370
371 Array<OneD, Array<OneD, NekDouble>> tmparray2D = NullNekDoubleArrayOfArray;
372
373 LibUtilities::Timer timer;
374 timer.Start();
375 m_FunctorDiffusionfluxCons(nDim, inarray, qfield, VolumeFlux, nonZeroIndex,
376 tmparray2D);
377 timer.Stop();
378 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxCons", 10);
379}
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition Diffusion.h:363

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

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

References 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

Definition at line 105 of file DiffusionIP.h.

Referenced by AddDiffusionSymmFluxToCoeff(), CalcTraceSymFlux(), and v_InitObject().

◆ 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.
static DiffusionSharedPtr create(std::string diffType)
Definition DiffusionIP.h:46
DiffusionFactory & GetDiffusionFactory()
Definition Diffusion.cpp:39

Definition at line 51 of file DiffusionIP.h.