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

#include <DiffusionIP.h>

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

Public Member Functions

template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void ConsVarAve (const size_t nConvectiveFields, const T &Bweight, const std::vector< T > &vFwd, const std::vector< T > &vBwd, std::vector< T > &aver)
 Calculate the average of conservative variables on traces. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Similar with Diffusion::Diffuse(): calculate diffusion flux The difference is in the outarray: it is the coefficients of basis for DiffuseCoeffs() it is the physics on quadrature points for Diffuse() More...
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, NekDouble time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, NekDouble time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
 Diffusion Volume FLux. More...
 
SOLVER_UTILS_EXPORT void DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
 Diffusion term Trace Flux. More...
 
SOLVER_UTILS_EXPORT void DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble >> &solution_Aver=NullNekDoubleArrayOfArray, Array< OneD, Array< OneD, NekDouble >> &solution_jump=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void AddDiffusionSymmFluxToCoeff (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 Add symmetric flux to field in coeff space. More...
 
SOLVER_UTILS_EXPORT void AddDiffusionSymmFluxToPhys (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 Add symmetric flux to field in coeff physical space. More...
 
SOLVER_UTILS_EXPORT void AddSymmFluxIntegralToOffDiag (const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &Fwdflux, TensorOfArray3D< NekDouble > &Bwdflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
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)
 
void SetHomoDerivs (Array< OneD, Array< OneD, NekDouble >> &deriv)
 
virtual TensorOfArray3D< NekDouble > & GetFluxTensor ()
 
SOLVER_UTILS_EXPORT void ConsVarAveJump (const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump)
 Get the average and jump value of conservative variables on trace. More...
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal ()
 Get trace normal. More...
 

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionIP ()
 
void GetPenaltyFactor (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get IP penalty factor based on order. More...
 
void GetPenaltyFactorConst (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
 Get a constant IP penalty factor. More...
 
void AddSymmFluxIntegralToCoeff (const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Add symmetric flux integration to field (in coefficient space) More...
 
void AddSymmFluxIntegralToPhys (const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
 Add symmetric flux integration to field (in physical space) More...
 
void CalcTraceSymFlux (const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump, Array< OneD, int > &nonZeroIndexsymm, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &traceSymflux)
 Calculate symmetric flux on traces. More...
 
void DiffuseTraceSymmFlux (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
 Calculate symmetric flux on traces interface. More...
 
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_DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex) 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...
 
void v_DiffuseTraceFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qfield, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qFwd, const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &qBwd, const Array< OneD, NekDouble > &MuAVTrace, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &Aver, const Array< OneD, Array< OneD, NekDouble >> &Jump)
 
virtual void v_AddDiffusionSymmFluxToCoeff (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
 
virtual void v_AddDiffusionSymmFluxToPhys (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
 
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 void v_ConsVarAveJump (const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump) override
 
const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal () override
 
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...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump)
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble >> &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual SOLVER_UTILS_EXPORT void v_GetPrimVar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
 Compute primary derivatives. More...
 
SOLVER_UTILS_EXPORT void GetDivCurl (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
 Compute divergence and curl squared. More...
 

Protected Attributes

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
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
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
 

Additional Inherited Members

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

Detailed Description

Definition at line 45 of file DiffusionIP.h.

Constructor & Destructor Documentation

◆ DiffusionIP()

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

Definition at line 50 of file DiffusionIP.cpp.

51 {
52 }

Referenced by create().

Member Function Documentation

◆ AddSecondDerivToTrace()

void Nektar::SolverUtils::DiffusionIP::AddSecondDerivToTrace ( const std::size_t  nConvectiveFields,
const size_t  nDim,
const size_t  nPts,
const size_t  nTracePts,
const NekDouble  PenaltyFactor2,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const TensorOfArray3D< NekDouble > &  qfield,
TensorOfArray3D< NekDouble > &  numDerivFwd,
TensorOfArray3D< NekDouble > &  numDerivBwd 
)
protected

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

Definition at line 871 of file DiffusionIP.cpp.

878 {
879  Array<OneD, NekDouble> Fwd(nTracePts, 0.0);
880  Array<OneD, NekDouble> Bwd(nTracePts, 0.0);
881  std::vector<NekDouble> tmp(nTracePts);
882 
883  Array<OneD, Array<OneD, NekDouble>> elmt2ndDerv{nDim};
884  for (int nd1 = 0; nd1 < nDim; ++nd1)
885  {
886  elmt2ndDerv[nd1] = Array<OneD, NekDouble>{nPts, 0.0};
887  }
888 
889  Array<OneD, Array<OneD, NekDouble>> qtmp{3};
890  for (int nd = 0; nd < 3; ++nd)
891  {
892  qtmp[nd] = NullNekDouble1DArray;
893  }
894  for (int nd2 = 0; nd2 < nDim; ++nd2)
895  {
896  qtmp[nd2] = elmt2ndDerv[nd2];
897  }
898 
899  for (int i = 0; i < nTracePts; ++i)
900  {
901  tmp[i] = PenaltyFactor2 * m_traceNormDirctnElmtLength[i];
902  }
903  // the derivatives are assumed to be exchangable
904  for (int nd1 = 0; nd1 < nDim; ++nd1)
905  {
906  for (int i = 0; i < nConvectiveFields; ++i)
907  {
908  fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
909 
910  for (int nd2 = nd1; nd2 < nDim; ++nd2)
911  {
912  Vmath::Zero(nTracePts, Bwd, 1);
913  fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd, true,
914  true, false);
915  for (int p = 0; p < nTracePts; ++p)
916  {
917  Bwd[p] *= tmp[p];
918  numDerivBwd[nd1][i][p] += m_traceNormals[nd2][p] * Bwd[p];
919  Fwd[p] *= tmp[p];
920  numDerivFwd[nd1][i][p] = -(m_traceNormals[nd2][p] * Fwd[p] -
921  numDerivFwd[nd1][i][p]);
922  }
923  if (nd2 != nd1)
924  {
925  for (int p = 0; p < nTracePts; ++p)
926  {
927  numDerivBwd[nd2][i][p] +=
928  m_traceNormals[nd1][p] * Bwd[p];
929  numDerivFwd[nd2][i][p] =
930  -(m_traceNormals[nd1][p] * Fwd[p] -
931  numDerivFwd[nd2][i][p]);
932  }
933  }
934  }
935  }
936  }
937 }
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:110
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:99
static Array< OneD, NekDouble > NullNekDouble1DArray
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492

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

Add symmetric flux integration to field (in coefficient space)

Definition at line 567 of file DiffusionIP.cpp.

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

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

Referenced by v_AddDiffusionSymmFluxToCoeff().

◆ AddSymmFluxIntegralToPhys()

void Nektar::SolverUtils::DiffusionIP::AddSymmFluxIntegralToPhys ( const std::size_t  nConvectiveFields,
const size_t  nDim,
const size_t  nPts,
const size_t  nTracePts,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, const int > &  nonZeroIndex,
TensorOfArray3D< NekDouble > &  tracflux,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
protected

Add symmetric flux integration to field (in physical space)

Definition at line 605 of file DiffusionIP.cpp.

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

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

Referenced by v_AddDiffusionSymmFluxToPhys().

◆ ApplyFluxBndConds()

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

aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic

Definition at line 944 of file DiffusionIP.cpp.

948 {
949  int nengy = nConvectiveFields - 1;
950  int cnt;
951  // Compute boundary conditions for Energy
952  cnt = 0;
953  size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
954  for (size_t j = 0; j < nBndRegions; ++j)
955  {
956  if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
958  {
959  continue;
960  }
961 
962  size_t nBndEdges =
963  fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
964  for (size_t e = 0; e < nBndEdges; ++e)
965  {
966  size_t nBndEdgePts = fields[nengy]
967  ->GetBndCondExpansions()[j]
968  ->GetExp(e)
969  ->GetTotPoints();
970 
971  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
972  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
973 
974  if (fields[0]->GetBndConditions()[j]->GetUserDefined() ==
975  "WallAdiabatic")
976  {
977  Vmath::Zero(nBndEdgePts, &flux[nengy][id2], 1);
978  }
979  }
980  }
981 }

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

Calculate numerical flux on traces.

Definition at line 757 of file DiffusionIP.cpp.

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

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddSecondDerivToTrace(), Nektar::SolverUtils::Diffusion::ConsVarAveJump(), 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 
)
protected

Calculate symmetric flux on traces.

Definition at line 544 of file DiffusionIP.cpp.

550 {
551  size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
552 
553  m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
554  nonZeroIndexsymm, m_traceNormals);
555 
556  for (int nd = 0; nd < nDim; ++nd)
557  {
558  for (int j = 0; j < nonZeroIndexsymm.size(); ++j)
559  {
560  int i = nonZeroIndexsymm[j];
561  Vmath::Smul(nTracePts, -0.5 * m_IPSymmFluxCoeff,
562  traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
563  }
564  }
565 }
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:408
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:248

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

Referenced by DiffuseTraceSymmFlux().

◆ ConsVarAve()

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

Calculate the average of conservative variables on traces.

Definition at line 60 of file DiffusionIP.h.

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

Referenced by v_ConsVarAveJump().

◆ create()

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

Definition at line 48 of file DiffusionIP.h.

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

References DiffusionIP().

◆ DiffuseTraceSymmFlux()

void Nektar::SolverUtils::DiffusionIP::DiffuseTraceSymmFlux ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  qfield,
TensorOfArray3D< NekDouble > &  VolumeFlux,
TensorOfArray3D< NekDouble > &  SymmFlux,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd,
Array< OneD, int > &  nonZeroIndex 
)
protected

Calculate symmetric flux on traces interface.

Definition at line 527 of file DiffusionIP.cpp.

536 {
537  boost::ignore_unused(inarray, qfield, VolumeFlux, pFwd, pBwd);
538  size_t nDim = fields[0]->GetCoordim(0);
539 
540  CalcTraceSymFlux(nConvectiveFields, nDim, m_traceAver, m_traceJump,
541  nonZeroIndex, SymmFlux);
542 }
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:101
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.
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:100

References CalcTraceSymFlux(), m_traceAver, and m_traceJump.

Referenced by v_AddDiffusionSymmFluxToCoeff(), and v_AddDiffusionSymmFluxToPhys().

◆ GetPenaltyFactor()

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

Get IP penalty factor based on order.

Definition at line 640 of file DiffusionIP.cpp.

643 {
644  MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
645  std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
646  tracelist->GetExp();
647  size_t ntotTrac = (*traceExp).size();
648  int nTracPnt, noffset;
649 
650  const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
651  fields[0]->GetLocTraceToTraceMap();
652 
653  const Array<OneD, const Array<OneD, int>> LRAdjExpid =
654  locTraceToTraceMap->GetLeftRightAdjacentExpId();
655  const Array<OneD, const Array<OneD, bool>> LRAdjflag =
656  locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
657 
658  std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
659  fields[0]->GetExp();
660 
661  Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
662 
663  NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
664 
665  int ntmp, numModes;
666 
667  for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
668  {
669  noffset = tracelist->GetPhys_Offset(ntrace);
670  nTracPnt = tracelist->GetTotPoints(ntrace);
671 
672  factorFwdBwd[0] = 0.0;
673  factorFwdBwd[1] = 0.0;
674 
675  for (int nlr = 0; nlr < 2; ++nlr)
676  {
677  if (LRAdjflag[nlr][ntrace])
678  {
679  numModes = 0;
680  for (int nd = 0; nd < spaceDim; nd++)
681  {
682  ntmp = fields[0]
683  ->GetExp(LRAdjExpid[nlr][ntrace])
684  ->GetBasisNumModes(nd);
685  numModes = std::max(ntmp, numModes);
686  }
687  factorFwdBwd[nlr] = (numModes) * (numModes);
688  }
689  }
690 
691  for (int np = 0; np < nTracPnt; ++np)
692  {
693  factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
694  }
695  }
696 }
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr

Referenced by CalcTraceNumFlux().

◆ GetPenaltyFactorConst()

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

Get a constant IP penalty factor.

Definition at line 698 of file DiffusionIP.cpp.

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

References Vmath::Fill(), and m_IPPenaltyCoeff.

◆ v_AddDiffusionSymmFluxToCoeff()

void Nektar::SolverUtils::DiffusionIP::v_AddDiffusionSymmFluxToCoeff ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  qfield,
TensorOfArray3D< NekDouble > &  VolumeFlux,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 458 of file DiffusionIP.cpp.

466 {
467  if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
468  {
469  size_t nDim = fields[0]->GetCoordim(0);
470  size_t nPts = fields[0]->GetTotPoints();
471  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
472  TensorOfArray3D<NekDouble> traceSymflux{nDim};
473  for (int nd = 0; nd < nDim; ++nd)
474  {
475  traceSymflux[nd] =
476  Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
477  for (int j = 0; j < nConvectiveFields; ++j)
478  {
479  traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
480  }
481  }
482  Array<OneD, int> nonZeroIndex;
483  DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
484  VolumeFlux, traceSymflux, pFwd, pBwd,
485  nonZeroIndex);
486 
487  AddSymmFluxIntegralToCoeff(nConvectiveFields, nDim, nPts, nTracePts,
488  fields, nonZeroIndex, traceSymflux,
489  outarray);
490  }
491 }
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Add symmetric flux integration to field (in coefficient space)
void DiffuseTraceSymmFlux(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
Calculate symmetric flux on traces interface.

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

◆ v_AddDiffusionSymmFluxToPhys()

void Nektar::SolverUtils::DiffusionIP::v_AddDiffusionSymmFluxToPhys ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  qfield,
TensorOfArray3D< NekDouble > &  VolumeFlux,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 493 of file DiffusionIP.cpp.

501 {
502  if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
503  {
504  size_t nDim = fields[0]->GetCoordim(0);
505  size_t nPts = fields[0]->GetTotPoints();
506  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
507  TensorOfArray3D<NekDouble> traceSymflux{nDim};
508  for (int nd = 0; nd < nDim; ++nd)
509  {
510  traceSymflux[nd] =
511  Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
512  for (int j = 0; j < nConvectiveFields; ++j)
513  {
514  traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
515  }
516  }
517  Array<OneD, int> nonZeroIndex;
518  DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
519  VolumeFlux, traceSymflux, pFwd, pBwd,
520  nonZeroIndex);
521 
522  AddSymmFluxIntegralToPhys(nConvectiveFields, nDim, nPts, nTracePts,
523  fields, nonZeroIndex, traceSymflux, outarray);
524  }
525 }
void AddSymmFluxIntegralToPhys(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble >> &outarray)
Add symmetric flux integration to field (in physical space)

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

◆ v_ConsVarAveJump()

void Nektar::SolverUtils::DiffusionIP::v_ConsVarAveJump ( const std::size_t  nConvectiveFields,
const size_t  npnts,
const Array< OneD, const Array< OneD, NekDouble >> &  vFwd,
const Array< OneD, const Array< OneD, NekDouble >> &  vBwd,
Array< OneD, Array< OneD, NekDouble >> &  aver,
Array< OneD, Array< OneD, NekDouble >> &  jump 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 706 of file DiffusionIP.cpp.

712 {
713  // ConsVarAve(nConvectiveFields, nPts, vFwd, vBwd, aver);
714 
715  std::vector<NekDouble> vFwdTmp(nConvectiveFields),
716  vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
717  for (size_t p = 0; p < nPts; ++p)
718  {
719  // re-arrange data
720  for (size_t f = 0; f < nConvectiveFields; ++f)
721  {
722  vFwdTmp[f] = vFwd[f][p];
723  vBwdTmp[f] = vBwd[f][p];
724  }
725 
726  ConsVarAve(nConvectiveFields, m_tracBwdWeightAver[p], vFwdTmp, vBwdTmp,
727  averTmp);
728 
729  // store
730  for (size_t f = 0; f < nConvectiveFields; ++f)
731  {
732  aver[f][p] = averTmp[f];
733  }
734  }
735 
736  // if this can be make function of points, the nPts loop can be lifted more
737  m_SpecialBndTreat(aver);
738 
739  // note: here the jump is 2.0*(aver-vFwd)
740  // because Viscous wall use a symmetry value as the Bwd,
741  // not the target one
742 
743  for (size_t f = 0; f < nConvectiveFields; ++f)
744  {
745  for (size_t p = 0; p < nPts; ++p)
746  {
747  NekDouble Bweight = m_tracBwdWeightJump[p];
748  NekDouble Fweight = 2.0 - Bweight;
749 
750  NekDouble tmpF = aver[f][p] - vFwd[f][p];
751  NekDouble tmpB = vBwd[f][p] - aver[f][p];
752  jump[f][p] = tmpF * Fweight + tmpB * Bweight;
753  }
754  }
755 }
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:407
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:108
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:60
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:109

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

◆ v_Diffuse()

void Nektar::SolverUtils::DiffusionIP::v_Diffuse ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
overrideprotectedvirtual

Reimplemented from 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:103

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

394 {
395  boost::ignore_unused(pFwd, pBwd);
396 
397  Array<OneD, Array<OneD, NekDouble>> qtmp{3};
398  size_t nDim = fields[0]->GetCoordim(0);
399  for (int nd = 0; nd < 3; ++nd)
400  {
401  qtmp[nd] = NullNekDouble1DArray;
402  }
403 
404  size_t nConvectiveFields = inarray.size();
405  for (int i = 0; i < nConvectiveFields; ++i)
406  {
407  for (int nd = 0; nd < nDim; ++nd)
408  {
409  qtmp[nd] = qfield[nd][i];
410  }
411  fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
412  }
413 }

References Nektar::NullNekDouble1DArray.

◆ v_DiffuseCoeffs() [1/2]

virtual void Nektar::SolverUtils::DiffusionIP::v_DiffuseCoeffs ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Referenced by v_Diffuse().

◆ v_DiffuseCoeffs() [2/2]

void Nektar::SolverUtils::DiffusionIP::v_DiffuseCoeffs ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, int > &  nonZeroIndex 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 321 of file DiffusionIP.cpp.

329 {
330  int i, j;
331  int nDim = fields[0]->GetCoordim(0);
332  int nPts = fields[0]->GetTotPoints();
333  int nCoeffs = fields[0]->GetNcoeffs();
334  int nTracePts = fields[0]->GetTrace()->GetTotPoints();
335 
336  TensorOfArray3D<NekDouble> elmtFlux(nDim);
337  for (j = 0; j < nDim; ++j)
338  {
339  elmtFlux[j] = Array<OneD, Array<OneD, NekDouble>>(nConvectiveFields);
340  for (i = 0; i < nConvectiveFields; ++i)
341  {
342  elmtFlux[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
343  }
344  }
345 
346  DiffuseVolumeFlux(fields, inarray, qfield, elmtFlux, nonZeroIndex);
347 
348  Array<OneD, Array<OneD, NekDouble>> tmpFluxIprdct(nDim);
349  // volume intergration: the nonZeroIndex indicates which flux is nonzero
350  for (i = 0; i < nonZeroIndex.size(); ++i)
351  {
352  int j = nonZeroIndex[i];
353  for (int k = 0; k < nDim; ++k)
354  {
355  tmpFluxIprdct[k] = elmtFlux[k][j];
356  }
357  fields[j]->IProductWRTDerivBase(tmpFluxIprdct, outarray[j]);
358  Vmath::Neg(nCoeffs, outarray[j], 1);
359  }
360 
361  Array<OneD, Array<OneD, NekDouble>> Traceflux(nConvectiveFields);
362  for (int j = 0; j < nConvectiveFields; ++j)
363  {
364  Traceflux[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
365  }
366 
367  DiffuseTraceFlux(fields, inarray, qfield, elmtFlux, Traceflux, vFwd, vBwd,
368  nonZeroIndex);
369 
370  // release qfield, elmtFlux and muvar;
371  for (j = 0; j < nDim; ++j)
372  {
373  elmtFlux[j] = NullNekDoubleArrayOfArray;
374  }
375 
376  for (i = 0; i < nonZeroIndex.size(); ++i)
377  {
378  int j = nonZeroIndex[i];
379 
380  fields[j]->AddTraceIntegral(Traceflux[j], outarray[j]);
381  fields[j]->SetPhysState(false);
382  }
383 
384  AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray, qfield,
385  elmtFlux, outarray, vFwd, vBwd);
386 }
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:206
SOLVER_UTILS_EXPORT void AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Add symmetric flux to field in coeff space.
Definition: Diffusion.h:254
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:217
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518

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

◆ v_DiffuseTraceFlux() [1/2]

void Nektar::SolverUtils::DiffusionIP::v_DiffuseTraceFlux ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  Qfields,
TensorOfArray3D< NekDouble > &  VolumeFlux,
Array< OneD, Array< OneD, NekDouble >> &  TraceFlux,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd,
Array< OneD, int > &  nonZeroIndex 
)
overrideprotectedvirtual

Diffusion term Trace Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 433 of file DiffusionIP.cpp.

441 {
442  boost::ignore_unused(VolumeFlux);
443 
444  TensorOfArray3D<NekDouble> traceflux3D(1);
445  traceflux3D[0] = TraceFlux;
446 
447  size_t nConvectiveFields = fields.size();
448  LibUtilities::Timer timer;
449  timer.Start();
450  CalcTraceNumFlux(m_IP2ndDervCoeff, fields, inarray, qfield, pFwd, pBwd,
451  nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
452  timer.Stop();
453  timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux", 10);
454 
455  ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
456 }
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 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

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

◆ v_DiffuseTraceFlux() [2/2]

void Nektar::SolverUtils::DiffusionIP::v_DiffuseTraceFlux ( const int  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &  qfield,
Array< OneD, Array< OneD, NekDouble >> &  TraceFlux,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd,
const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &  qFwd,
const Array< OneD, const Array< OneD, Array< OneD, NekDouble >>> &  qBwd,
const Array< OneD, NekDouble > &  MuAVTrace,
Array< OneD, int > &  nonZeroIndex,
const Array< OneD, Array< OneD, NekDouble >> &  Aver,
const Array< OneD, Array< OneD, NekDouble >> &  Jump 
)
protected

◆ v_DiffuseVolumeFlux()

void Nektar::SolverUtils::DiffusionIP::v_DiffuseVolumeFlux ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  qfields,
TensorOfArray3D< NekDouble > &  VolumeFlux,
Array< OneD, int > &  nonZeroIndex 
)
overrideprotectedvirtual

Diffusion Volume Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 415 of file DiffusionIP.cpp.

420 {
421  size_t nDim = fields[0]->GetCoordim(0);
422 
423  Array<OneD, Array<OneD, NekDouble>> tmparray2D = NullNekDoubleArrayOfArray;
424 
425  LibUtilities::Timer timer;
426  timer.Start();
427  m_FunctorDiffusionfluxCons(nDim, inarray, qfield, VolumeFlux, nonZeroIndex,
428  tmparray2D);
429  timer.Stop();
430  timer.AccumulateRegion("DiffIP:_FunctorDiffFluxCons", 10);
431 }
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:405

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

257  {
258  return m_traceNormals;
259  }

References m_traceNormals.

◆ v_InitObject()

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

Reimplemented from 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;
129  Vmath::Smul(nTracePts, 0.25, m_traceNormDirctnElmtLengthRecip, 1,
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:112
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:209
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:553
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:284
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:419
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
protected

Definition at line 96 of file DiffusionIP.h.

Referenced by v_DiffuseTraceFlux(), and v_InitObject().

◆ m_IPPenaltyCoeff

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

Definition at line 97 of file DiffusionIP.h.

Referenced by GetPenaltyFactorConst(), and v_InitObject().

◆ m_IPSymmFluxCoeff

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

◆ m_session

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

Definition at line 112 of file DiffusionIP.h.

Referenced by v_InitObject().

◆ m_tracBwdWeightAver

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

Definition at line 108 of file DiffusionIP.h.

Referenced by v_ConsVarAveJump(), and v_InitObject().

◆ m_tracBwdWeightJump

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

Definition at line 109 of file DiffusionIP.h.

Referenced by v_ConsVarAveJump(), and v_InitObject().

◆ m_traceAver

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

Definition at line 100 of file DiffusionIP.h.

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

◆ m_traceJump

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

Definition at line 101 of file DiffusionIP.h.

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

◆ m_traceNormals

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

◆ m_traceNormDirctnElmtLength

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

Definition at line 110 of file DiffusionIP.h.

Referenced by AddSecondDerivToTrace(), and v_InitObject().

◆ m_traceNormDirctnElmtLengthRecip

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

Definition at line 111 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ m_wspDiff

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

Workspace for v_Diffusion.

Definition at line 103 of file DiffusionIP.h.

Referenced by v_Diffuse(), and v_InitObject().

◆ m_wspNumDerivBwd

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

Workspace for CallTraceNumFlux.

Definition at line 105 of file DiffusionIP.h.

Referenced by CalcTraceNumFlux(), and v_InitObject().

◆ m_wspNumDerivFwd

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

Definition at line 106 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.