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

#include <Diffusion.h>

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

Public Member Functions

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...
 

Public Attributes

Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 

Protected Member Functions

virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
virtual SOLVER_UTILS_EXPORT void v_Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &vFwd, const Array< OneD, Array< OneD, NekDouble >> &vBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
virtual SOLVER_UTILS_EXPORT void v_ConsVarAveJump (const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump)
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseCalcDerivative (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 Diffusion Flux, calculate the physical derivatives. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseVolumeFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex)
 Diffusion Volume Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &Qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
 Diffusion term Trace Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble >> &solution_Aver, Array< OneD, Array< OneD, NekDouble >> &solution_jump)
 
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToCoeff (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToPhys (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble >> &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal ()
 
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

DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
NekDouble m_time = 0.0
 

Detailed Description

Definition at line 127 of file Diffusion.h.

Constructor & Destructor Documentation

◆ ~Diffusion()

Diffusion::~Diffusion ( )
inlinevirtual

Definition at line 135 of file Diffusion.h.

135 {};

Member Function Documentation

◆ AddDiffusionSymmFluxToCoeff()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Add symmetric flux to field in coeff space.

Definition at line 254 of file Diffusion.h.

263  {
264  v_AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray,
265  qfield, VolumeFlux, outarray, pFwd, pBwd);
266  }
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Definition: Diffusion.h:496

References v_AddDiffusionSymmFluxToCoeff().

Referenced by Nektar::SolverUtils::DiffusionIP::v_DiffuseCoeffs().

◆ AddDiffusionSymmFluxToPhys()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Add symmetric flux to field in coeff physical space.

Definition at line 269 of file Diffusion.h.

278  {
279  v_AddDiffusionSymmFluxToPhys(nConvectiveFields, fields, inarray, qfield,
280  VolumeFlux, outarray, pFwd, pBwd);
281  }
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToPhys(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Definition: Diffusion.h:509

References v_AddDiffusionSymmFluxToPhys().

◆ AddSymmFluxIntegralToOffDiag()

void Diffusion::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 
)

Definition at line 111 of file Diffusion.cpp.

118 {
119  boost::ignore_unused(nTracePts);
120 
121  int nCoeffs = outarray[nConvectiveFields - 1].size();
122  Array<OneD, NekDouble> tmpCoeff(nCoeffs, 0.0);
123  Array<OneD, Array<OneD, NekDouble>> tmpfield(nDim);
124 
125  for (int i = 0; i < nDim; i++)
126  {
127  tmpfield[i] = Array<OneD, NekDouble>(nPts, 0.0);
128  }
129  int nv = 0;
130 
131  for (int j = 0; j < nonZeroIndex.size(); j++)
132  {
133  nv = nonZeroIndex[j];
134  MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
135  for (int nd = 0; nd < nDim; nd++)
136  {
137  Vmath::Zero(nPts, tmpfield[nd], 1);
138 
139  tracelist->MultiplyByQuadratureMetric(Fwdflux[nd][nv],
140  Fwdflux[nd][nv]);
141  tracelist->MultiplyByQuadratureMetric(Bwdflux[nd][nv],
142  Bwdflux[nd][nv]);
143 
144  fields[nv]->AddTraceQuadPhysToOffDiag(
145  Fwdflux[nd][nv], Bwdflux[nd][nv], tmpfield[nd]);
146  fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
147  }
148  fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
149  Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
150  }
151 }
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
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492

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

◆ ConsVarAveJump()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Get the average and jump value of conservative variables on trace.

Definition at line 384 of file Diffusion.h.

390  {
391  v_ConsVarAveJump(nConvectiveFields, npnts, vFwd, vBwd, aver, jump);
392  }
virtual SOLVER_UTILS_EXPORT void v_ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble >> &vFwd, const Array< OneD, const Array< OneD, NekDouble >> &vBwd, Array< OneD, Array< OneD, NekDouble >> &aver, Array< OneD, Array< OneD, NekDouble >> &jump)
Definition: Diffusion.cpp:210

References v_ConsVarAveJump().

Referenced by Nektar::SolverUtils::DiffusionIP::CalcTraceNumFlux().

◆ Diffuse() [1/2]

void Diffusion::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 
)

Definition at line 59 of file Diffusion.cpp.

66 {
67  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
68 }
virtual SOLVER_UTILS_EXPORT void v_Diffuse(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Definition: Diffusion.cpp:153

References v_Diffuse().

◆ Diffuse() [2/2]

void Diffusion::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 
)

Definition at line 87 of file Diffusion.cpp.

94 {
95  m_time = time;
96  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
97 }

References m_time, and v_Diffuse().

◆ DiffuseCalcDerivative()

void Diffusion::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 at line 194 of file Diffusion.cpp.

200 {
201  v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
202 }
virtual SOLVER_UTILS_EXPORT void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Diffusion Flux, calculate the physical derivatives.
Definition: Diffusion.cpp:221

References v_DiffuseCalcDerivative().

Referenced by Nektar::SolverUtils::DiffusionLDG::v_DiffuseCoeffs(), and Nektar::DiffusionLDGNS::v_DiffuseCoeffs().

◆ DiffuseCoeffs() [1/3]

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Definition at line 171 of file Diffusion.h.

179  {
180  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd,
181  pBwd, qfield, nonZeroIndex);
182  }
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)
Definition: Diffusion.cpp:166

References v_DiffuseCoeffs().

◆ DiffuseCoeffs() [2/3]

void Diffusion::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()

Definition at line 76 of file Diffusion.cpp.

83 {
84  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
85 }

References v_DiffuseCoeffs().

◆ DiffuseCoeffs() [3/3]

void Diffusion::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 
)

Definition at line 99 of file Diffusion.cpp.

106 {
107  m_time = time;
108  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
109 }

References m_time, and v_DiffuseCoeffs().

◆ DiffuseTraceFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Diffusion term Trace Flux.

Definition at line 217 of file Diffusion.h.

228  {
229  v_DiffuseTraceFlux(fields, inarray, qfields, VolumeFlux, TraceFlux,
230  pFwd, pBwd, nonZeroIndex);
231  }
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &Qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble >> &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
Diffusion term Trace Flux.
Definition: Diffusion.cpp:242

References v_DiffuseTraceFlux().

Referenced by Nektar::SolverUtils::DiffusionLDG::v_DiffuseCoeffs(), Nektar::DiffusionLDGNS::v_DiffuseCoeffs(), and Nektar::SolverUtils::DiffusionIP::v_DiffuseCoeffs().

◆ DiffuseTraceSymmFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Definition at line 233 of file Diffusion.h.

247  {
248  v_DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
249  VolumeFlux, SymmFlux, pFwd, pBwd, nonZeroIndex,
250  solution_Aver, solution_jump);
251  }
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)
Definition: Diffusion.h:477

References v_DiffuseTraceSymmFlux().

◆ DiffuseVolumeFlux()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inline

Diffusion Volume FLux.

Definition at line 206 of file Diffusion.h.

212  {
213  v_DiffuseVolumeFlux(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
214  }
virtual SOLVER_UTILS_EXPORT void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex)
Diffusion Volume Flux.
Definition: Diffusion.cpp:232

References v_DiffuseVolumeFlux().

Referenced by Nektar::SolverUtils::DiffusionLDG::v_DiffuseCoeffs(), Nektar::DiffusionLDGNS::v_DiffuseCoeffs(), and Nektar::SolverUtils::DiffusionIP::v_DiffuseCoeffs().

◆ GetDivCurl()

void Diffusion::GetDivCurl ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const TensorOfArray3D< NekDouble > &  pVarDer 
)
protected

Compute divergence and curl squared.

Compute and store scalars needed for shock sensor, i.e. divergence and curl squared.

Parameters
dimensions,points,primaryvariables derivative

Definition at line 281 of file Diffusion.cpp.

284 {
285  int nDim = fields[0]->GetCoordim(0);
286  int nPts = fields[0]->GetTotPoints();
287  // Compute and store scalars needed for shock sensor
288  Array<OneD, NekDouble> tmp3(nPts, 0.0);
289  Array<OneD, NekDouble> tmp4(nPts, 0.0);
290 
291  // Zero the variables for the current iteration step
292  Vmath::Zero(nPts, m_divVel, 1);
293  Vmath::Zero(nPts, m_divVelSquare, 1);
294  Vmath::Zero(nPts, m_curlVelSquare, 1);
295 
296  // div vel
297  for (int j = 0; j < nDim; ++j)
298  {
299  Vmath::Vadd(nPts, m_divVel, 1, pVarDer[j][j], 1, m_divVel, 1);
300  }
301  // (div vel)**2
302  Vmath::Vmul(nPts, m_divVel, 1, m_divVel, 1, m_divVelSquare, 1);
303 
304  // curl vel
305  if (nDim > 2)
306  {
307  // curl[0] 3/2
308  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][1], 1, tmp4, 1);
309  // curl[0]-2/3
310  Vmath::Vcopy(nPts, pVarDer[1][2], 1, tmp3, 1);
311  Vmath::Neg(nPts, tmp3, 1);
312  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
313  // square curl[0]
314  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
315  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
316 
317  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
318  // curl[1] 3/1
319  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][0], 1, tmp4, 1);
320  // curl[1]-1/3
321  Vmath::Vcopy(nPts, pVarDer[0][2], 1, tmp3, 1);
322  Vmath::Neg(nPts, tmp3, 1);
323  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
324  // square curl[1]
325  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
326  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
327 
328  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
329  // curl[2] 1/2
330  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
331  // curl[2]-2/1
332  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
333  Vmath::Neg(nPts, tmp3, 1);
334  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
335  // square curl[2]
336  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
337  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
338  }
339  else if (nDim > 1)
340  {
341  // curl[2] 1/2
342  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
343  // curl[2]-2/1
344  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
345  Vmath::Neg(nPts, tmp3, 1);
346  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
347  // square curl[2]
348  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
349  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
350  }
351  else
352  {
353  Vmath::Fill(nPts, 0.0, m_curlVelSquare, 1);
354  }
355 }
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:132
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:131
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:133
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 Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Vmath::Fill(), m_curlVelSquare, m_divVel, m_divVelSquare, Vmath::Neg(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Zero().

◆ GetFluxTensor()

virtual TensorOfArray3D<NekDouble>& Nektar::SolverUtils::Diffusion::GetFluxTensor ( )
inlinevirtual

Definition at line 378 of file Diffusion.h.

379  {
380  return v_GetFluxTensor();
381  }
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:528

References v_GetFluxTensor().

◆ GetTraceNormal()

SOLVER_UTILS_EXPORT const Array<OneD, const Array<OneD, NekDouble> >& Nektar::SolverUtils::Diffusion::GetTraceNormal ( )
inline

Get trace normal.

Definition at line 396 of file Diffusion.h.

397  {
398  return v_GetTraceNormal();
399  }
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:204

References v_GetTraceNormal().

◆ InitObject()

void Diffusion::InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)

Definition at line 47 of file Diffusion.cpp.

49 {
50  v_InitObject(pSession, pFields);
51 
52  // Div curl storage
53  int nPts = pFields[0]->GetTotPoints();
54  m_divVel = Array<OneD, NekDouble>(nPts, 0.0);
55  m_divVelSquare = Array<OneD, NekDouble>(nPts, 0.0);
56  m_curlVelSquare = Array<OneD, NekDouble>(nPts, 0.0);
57 }
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Definition: Diffusion.h:412

References m_curlVelSquare, m_divVel, m_divVelSquare, and v_InitObject().

◆ SetDiffusionFluxCons() [1/2]

void Nektar::SolverUtils::Diffusion::SetDiffusionFluxCons ( DiffusionFluxCons  flux)
inline

Definition at line 339 of file Diffusion.h.

340  {
342  }
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:405

References m_FunctorDiffusionfluxCons.

◆ SetDiffusionFluxCons() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetDiffusionFluxCons ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 331 of file Diffusion.h.

332  {
334  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
335  std::placeholders::_3, std::placeholders::_4,
336  std::placeholders::_5, std::placeholders::_6);
337  }

References m_FunctorDiffusionfluxCons.

◆ SetDiffusionFluxConsTrace() [1/2]

void Nektar::SolverUtils::Diffusion::SetDiffusionFluxConsTrace ( DiffusionFluxCons  flux)
inline

Definition at line 353 of file Diffusion.h.

354  {
356  }
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:406

References m_FunctorDiffusionfluxConsTrace.

◆ SetDiffusionFluxConsTrace() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetDiffusionFluxConsTrace ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 345 of file Diffusion.h.

346  {
348  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
349  std::placeholders::_3, std::placeholders::_4,
350  std::placeholders::_5, std::placeholders::_6);
351  }

References m_FunctorDiffusionfluxConsTrace.

◆ SetDiffusionSymmFluxCons()

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetDiffusionSymmFluxCons ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 365 of file Diffusion.h.

366  {
368  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
369  std::placeholders::_3, std::placeholders::_4,
370  std::placeholders::_5, std::placeholders::_6);
371  }
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:408

References m_FunctorSymmetricfluxCons.

◆ SetFluxPenaltyNS() [1/2]

void Nektar::SolverUtils::Diffusion::SetFluxPenaltyNS ( DiffusionFluxPenaltyNS  flux)
inline

Definition at line 325 of file Diffusion.h.

326  {
327  m_fluxPenaltyNS = flux;
328  }
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:404

References m_fluxPenaltyNS.

◆ SetFluxPenaltyNS() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetFluxPenaltyNS ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 318 of file Diffusion.h.

319  {
321  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
322  std::placeholders::_3);
323  }

References m_fluxPenaltyNS.

◆ SetFluxVector() [1/2]

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::SetFluxVector ( DiffusionFluxVecCB  fluxVector)
inline

Definition at line 299 of file Diffusion.h.

300  {
301  m_fluxVector = fluxVector;
302  }
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:402

References m_fluxVector.

◆ SetFluxVector() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetFluxVector ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 293 of file Diffusion.h.

294  {
295  m_fluxVector = std::bind(func, obj, std::placeholders::_1,
296  std::placeholders::_2, std::placeholders::_3);
297  }

References m_fluxVector.

◆ SetFluxVectorNS() [1/2]

void Nektar::SolverUtils::Diffusion::SetFluxVectorNS ( DiffusionFluxVecCBNS  fluxVector)
inline

Definition at line 312 of file Diffusion.h.

313  {
314  m_fluxVectorNS = fluxVector;
315  }
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:403

References m_fluxVectorNS.

◆ SetFluxVectorNS() [2/2]

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetFluxVectorNS ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 305 of file Diffusion.h.

306  {
308  std::bind(func, obj, std::placeholders::_1, std::placeholders::_2,
309  std::placeholders::_3);
310  }

References m_fluxVectorNS.

◆ SetHomoDerivs()

void Nektar::SolverUtils::Diffusion::SetHomoDerivs ( Array< OneD, Array< OneD, NekDouble >> &  deriv)
inline

Definition at line 373 of file Diffusion.h.

374  {
375  v_SetHomoDerivs(deriv);
376  }
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble >> &deriv)
Definition: Diffusion.h:523

References v_SetHomoDerivs().

◆ SetSpecialBndTreat()

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::SolverUtils::Diffusion::SetSpecialBndTreat ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 359 of file Diffusion.h.

360  {
361  m_SpecialBndTreat = std::bind(func, obj, std::placeholders::_1);
362  }
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:407

References m_SpecialBndTreat.

◆ v_AddDiffusionSymmFluxToCoeff()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 496 of file Diffusion.h.

505  {
506  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
507  VolumeFlux, outarray, pFwd, pBwd);
508  }

Referenced by AddDiffusionSymmFluxToCoeff().

◆ v_AddDiffusionSymmFluxToPhys()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 509 of file Diffusion.h.

518  {
519  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
520  VolumeFlux, outarray, pFwd, pBwd);
521  }

Referenced by AddDiffusionSymmFluxToPhys().

◆ v_ConsVarAveJump()

void Diffusion::v_ConsVarAveJump ( const std::size_t  nConvectiveFields,
const size_t  npnts,
const Array< OneD, const Array< OneD, NekDouble >> &  vFwd,
const Array< OneD, const Array< OneD, NekDouble >> &  vBwd,
Array< OneD, Array< OneD, NekDouble >> &  aver,
Array< OneD, Array< OneD, NekDouble >> &  jump 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 210 of file Diffusion.cpp.

216 {
217  boost::ignore_unused(nConvectiveFields, npnts, vFwd, vBwd, aver, jump);
218  NEKERROR(ErrorUtil::efatal, "v_ConsVarAveJump not defined");
219 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by ConsVarAveJump().

◆ v_Diffuse()

void Diffusion::v_Diffuse ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionLFR, Nektar::SolverUtils::DiffusionIP, Nektar::SolverUtils::Diffusion3DHomogeneous1D, Nektar::DiffusionLDGNS, Nektar::SolverUtils::DiffusionLFRNS, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 153 of file Diffusion.cpp.

160 {
161  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, pFwd,
162  pBwd);
163  NEKERROR(ErrorUtil::efatal, "v_Diffuse not defined");
164 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by Diffuse().

◆ v_DiffuseCalcDerivative()

void Diffusion::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 
)
protectedvirtual

Diffusion Flux, calculate the physical derivatives.

Reimplemented in Nektar::DiffusionLDGNS, Nektar::SolverUtils::DiffusionLDG, and Nektar::SolverUtils::DiffusionIP.

Definition at line 221 of file Diffusion.cpp.

227 {
228  boost::ignore_unused(fields, inarray, qfields, pFwd, pBwd);
229  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
230 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by DiffuseCalcDerivative().

◆ v_DiffuseCoeffs() [1/2]

void Diffusion::v_DiffuseCoeffs ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP, Nektar::DiffusionLDGNS, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 166 of file Diffusion.cpp.

173 {
174  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, pFwd,
175  pBwd);
176  NEKERROR(ErrorUtil::efatal, "v_DiffuseCoeffs not defined");
177 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by DiffuseCoeffs().

◆ v_DiffuseCoeffs() [2/2]

void Diffusion::v_DiffuseCoeffs ( const std::size_t  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  vFwd,
const Array< OneD, Array< OneD, NekDouble >> &  vBwd,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, int > &  nonZeroIndex 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 179 of file Diffusion.cpp.

187 {
188  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray, vFwd,
189  vBwd, qfield, nonZeroIndex);
190  NEKERROR(ErrorUtil::efatal, "v_DiffuseCoeffs not defined");
191 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

◆ v_DiffuseTraceFlux()

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

Diffusion term Trace Flux.

Reimplemented in Nektar::DiffusionLDGNS, Nektar::SolverUtils::DiffusionLDG, and Nektar::SolverUtils::DiffusionIP.

Definition at line 242 of file Diffusion.cpp.

250 {
251  boost::ignore_unused(fields, inarray, qfields, VolumeFlux, TraceFlux, pFwd,
252  pBwd, nonZeroIndex);
253  NEKERROR(ErrorUtil::efatal, "Not defined function DiffuseTraceFLux.");
254 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by DiffuseTraceFlux().

◆ v_DiffuseTraceSymmFlux()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::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 
)
inlineprotectedvirtual

Definition at line 477 of file Diffusion.h.

489  {
490  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
491  VolumeFlux, SymmFlux, pFwd, pBwd, nonZeroIndex,
492  solution_Aver, solution_jump);
493  ASSERTL0(false, "Not defined function v_DiffuseTraceSymmFlux.");
494  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

References ASSERTL0.

Referenced by DiffuseTraceSymmFlux().

◆ v_DiffuseVolumeFlux()

void Diffusion::v_DiffuseVolumeFlux ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  qfields,
TensorOfArray3D< NekDouble > &  VolumeFlux,
Array< OneD, int > &  nonZeroIndex 
)
protectedvirtual

Diffusion Volume Flux.

Reimplemented in Nektar::DiffusionLDGNS, Nektar::SolverUtils::DiffusionLDG, and Nektar::SolverUtils::DiffusionIP.

Definition at line 232 of file Diffusion.cpp.

237 {
238  boost::ignore_unused(fields, inarray, qfields, VolumeFlux, nonZeroIndex);
239  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
240 }

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by DiffuseVolumeFlux().

◆ v_GetFluxTensor()

virtual TensorOfArray3D<NekDouble>& Nektar::SolverUtils::Diffusion::v_GetFluxTensor ( )
inlineprotectedvirtual

Reimplemented in Nektar::DiffusionLDGNS, and Nektar::SolverUtils::DiffusionLFRNS.

Definition at line 528 of file Diffusion.h.

529  {
530  static TensorOfArray3D<NekDouble> tmp;
531  return tmp;
532  }

Referenced by GetFluxTensor().

◆ v_GetPrimVar()

void Diffusion::v_GetPrimVar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  primVar 
)
protectedvirtual

Compute primary derivatives.

Compute primary variables.

Reimplemented in Nektar::DiffusionLDGNS.

Definition at line 260 of file Diffusion.cpp.

264 {
265 
266  int nDim = fields[0]->GetCoordim(0);
267  int nPts = fields[0]->GetTotPoints();
268  for (int i = 0; i < nDim; ++i)
269  {
270  primVar[i] = Array<OneD, NekDouble>(nPts, 0.0);
271  Vmath::Vdiv(nPts, inarray[i + 1], 1, inarray[0], 1, primVar[i], 1);
272  }
273 }
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

References Vmath::Vdiv().

◆ v_GetTraceNormal()

const Array< OneD, const Array< OneD, NekDouble > > & Diffusion::v_GetTraceNormal ( )
protectedvirtual

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 204 of file Diffusion.cpp.

205 {
206  NEKERROR(ErrorUtil::efatal, "v_GetTraceNormal not defined");
208 }
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References Nektar::ErrorUtil::efatal, NEKERROR, and Nektar::NullNekDoubleArrayOfArray.

Referenced by GetTraceNormal().

◆ v_InitObject()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
inlineprotectedvirtual

◆ v_SetHomoDerivs()

virtual void Nektar::SolverUtils::Diffusion::v_SetHomoDerivs ( Array< OneD, Array< OneD, NekDouble >> &  deriv)
inlineprotectedvirtual

Reimplemented in Nektar::DiffusionLDGNS, and Nektar::SolverUtils::DiffusionLFRNS.

Definition at line 523 of file Diffusion.h.

524  {
525  boost::ignore_unused(deriv);
526  }

Referenced by SetHomoDerivs().

Member Data Documentation

◆ m_curlVelSquare

Array<OneD, NekDouble> Nektar::SolverUtils::Diffusion::m_curlVelSquare

Definition at line 133 of file Diffusion.h.

Referenced by GetDivCurl(), and InitObject().

◆ m_divVel

Array<OneD, NekDouble> Nektar::SolverUtils::Diffusion::m_divVel

Params for Ducros sensor.

Definition at line 131 of file Diffusion.h.

Referenced by GetDivCurl(), and InitObject().

◆ m_divVelSquare

Array<OneD, NekDouble> Nektar::SolverUtils::Diffusion::m_divVelSquare

Definition at line 132 of file Diffusion.h.

Referenced by GetDivCurl(), and InitObject().

◆ m_fluxPenaltyNS

DiffusionFluxPenaltyNS Nektar::SolverUtils::Diffusion::m_fluxPenaltyNS
protected

Definition at line 404 of file Diffusion.h.

Referenced by Nektar::DiffusionLDGNS::NumericalFluxO2(), and SetFluxPenaltyNS().

◆ m_fluxVector

DiffusionFluxVecCB Nektar::SolverUtils::Diffusion::m_fluxVector
protected

◆ m_fluxVectorNS

DiffusionFluxVecCBNS Nektar::SolverUtils::Diffusion::m_fluxVectorNS
protected

◆ m_FunctorDiffusionfluxCons

DiffusionFluxCons Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxCons
protected

◆ m_FunctorDiffusionfluxConsTrace

DiffusionFluxCons Nektar::SolverUtils::Diffusion::m_FunctorDiffusionfluxConsTrace
protected

◆ m_FunctorSymmetricfluxCons

DiffusionSymmFluxCons Nektar::SolverUtils::Diffusion::m_FunctorSymmetricfluxCons
protected

◆ m_SpecialBndTreat

SpecialBndTreat Nektar::SolverUtils::Diffusion::m_SpecialBndTreat
protected

◆ m_time

NekDouble Nektar::SolverUtils::Diffusion::m_time = 0.0
protected

Definition at line 410 of file Diffusion.h.

Referenced by Diffuse(), and DiffuseCoeffs().