Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
SOLVER_UTILS_EXPORT void FluxVec (TensorOfArray3D< NekDouble > &fluxvector)
 
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 SetArtificialDiffusionVector (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetSpecialBndTreat (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetCalcViscosity (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 GetAVmu (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &muvar, Array< OneD, NekDouble > &MuVarTrace)
 Get the mu of artifical viscosity(AV) More...
 
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
 
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
CalcViscosity m_CalcViscosity
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
NekDouble m_time =0.0
 

Detailed Description

Definition at line 155 of file Diffusion.h.

Constructor & Destructor Documentation

◆ ~Diffusion()

Diffusion::~Diffusion ( )
inlinevirtual

Definition at line 163 of file Diffusion.h.

164  {};

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 285 of file Diffusion.h.

294  {
295  v_AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields,
296  inarray, qfield, VolumeFlux, outarray, pFwd, pBwd);
297  }
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:568

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 300 of file Diffusion.h.

309  {
310  v_AddDiffusionSymmFluxToPhys(nConvectiveFields, fields, inarray,
311  qfield, VolumeFlux, outarray, pFwd, pBwd);
312  }
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:582

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 137 of file Diffusion.cpp.

147  {
148  boost::ignore_unused(nTracePts);
149 
150  int nCoeffs = outarray[nConvectiveFields-1].size();
151  Array<OneD, NekDouble> tmpCoeff(nCoeffs,0.0);
152  Array<OneD, Array<OneD, NekDouble> > tmpfield(nDim);
153 
154  for(int i = 0;i<nDim;i++)
155  {
156  tmpfield[i] = Array<OneD, NekDouble>(nPts,0.0);
157  }
158  int nv = 0;
159 
160  for(int j=0;j<nonZeroIndex.size();j++)
161  {
162  nv = nonZeroIndex[j];
163  MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
164  for(int nd=0;nd<nDim;nd++)
165  {
166  Vmath::Zero(nPts,tmpfield[nd],1);
167 
168  tracelist->MultiplyByQuadratureMetric(Fwdflux[nd][nv],
169  Fwdflux[nd][nv]);
170  tracelist->MultiplyByQuadratureMetric(Bwdflux[nd][nv],
171  Bwdflux[nd][nv]);
172 
173  fields[nv]->AddTraceQuadPhysToOffDiag(Fwdflux[nd][nv],
174  Bwdflux[nd][nv],tmpfield[nd]);
175  fields[nv]->DivideByQuadratureMetric(tmpfield[nd],tmpfield[nd]);
176  }
177  fields[nv]->IProductWRTDerivBase(tmpfield,tmpCoeff);
178  Vmath::Vadd(nCoeffs,tmpCoeff,1,outarray[nv],1,outarray[nv],1);
179  }
180  }
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:322
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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 446 of file Diffusion.h.

453  {
454  v_ConsVarAveJump(nConvectiveFields, npnts, vFwd, vBwd, aver,
455  jump);
456  }
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:242

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 60 of file Diffusion.cpp.

67  {
68  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
69  }
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:182

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 89 of file Diffusion.cpp.

97  {
98  m_time = time;
99  v_Diffuse(nConvectiveFields, fields, inarray, outarray, pFwd, pBwd);
100  }

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 224 of file Diffusion.cpp.

230  {
231  v_DiffuseCalcDerivative(fields, inarray, qfields, pFwd, pBwd);
232  }
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:255

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 199 of file Diffusion.h.

208  {
209  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray,
210  pFwd, pBwd,qfield,nonZeroIndex);
211  }
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:195

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 77 of file Diffusion.cpp.

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

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 102 of file Diffusion.cpp.

110  {
111  m_time = time;
112  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, outarray,
113  pFwd, pBwd);
114  }

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 249 of file Diffusion.h.

261  {
262  v_DiffuseTraceFlux(fields, inarray, qfields,
263  VolumeFlux, TraceFlux, pFwd, pBwd, nonZeroIndex);
264  }
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:279

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 266 of file Diffusion.h.

278  {
279  v_DiffuseTraceSymmFlux(nConvectiveFields,fields,inarray,qfield,
280  VolumeFlux,SymmFlux,pFwd,pBwd,
281  nonZeroIndex,solution_Aver,solution_jump);
282  }
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:549

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 236 of file Diffusion.h.

243  {
244  v_DiffuseVolumeFlux(fields, inarray, qfields,
245  VolumeFlux, nonZeroIndex);
246  }
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:267

References v_DiffuseVolumeFlux().

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

◆ FluxVec()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Diffusion::FluxVec ( TensorOfArray3D< NekDouble > &  fluxvector)

◆ GetAVmu()

void Diffusion::GetAVmu ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  muvar,
Array< OneD, NekDouble > &  MuVarTrace 
)

Get the mu of artifical viscosity(AV)

Definition at line 115 of file Diffusion.cpp.

120  {
121  size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
122 
123  Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
124  Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
125 
126  m_ArtificialDiffusionVector(inarray, muvar);
127 
128  // BwdMuvar is left to be 0.0 according to DiffusionLDG.cpp
129  fields[0]->GetFwdBwdTracePhys(muvar, Fwd, Bwd, false);
130 
131  for (int k = 0; k < nTracePts; ++k)
132  {
133  MuVarTrace[k] = 0.5 * (Fwd[k] + Bwd[k]) ;
134  }
135  }
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
Definition: Diffusion.h:469

References m_ArtificialDiffusionVector.

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

◆ 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 319 of file Diffusion.cpp.

322  {
323  int nDim = fields[0]->GetCoordim(0);
324  int nPts = fields[0]->GetTotPoints();
325  // Compute and store scalars needed for shock sensor
326  Array<OneD, NekDouble> tmp3(nPts, 0.0);
327  Array<OneD, NekDouble> tmp4(nPts, 0.0);
328 
329  // Zero the variables for the current iteration step
330  Vmath::Zero(nPts, m_divVel, 1);
331  Vmath::Zero(nPts, m_divVelSquare, 1);
332  Vmath::Zero(nPts, m_curlVelSquare, 1);
333 
334  // div vel
335  for (int j = 0; j < nDim; ++j)
336  {
337  Vmath::Vadd(nPts, m_divVel, 1, pVarDer[j][j], 1, m_divVel, 1);
338  }
339  // (div vel)**2
340  Vmath::Vmul(nPts, m_divVel, 1, m_divVel, 1, m_divVelSquare, 1);
341 
342  // curl vel
343  if ( nDim > 2 )
344  {
345  // curl[0] 3/2
346  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][1], 1, tmp4, 1);
347  // curl[0]-2/3
348  Vmath::Vcopy(nPts, pVarDer[1][2], 1, tmp3, 1);
349  Vmath::Neg(nPts, tmp3, 1);
350  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
351  // square curl[0]
352  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
353  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
354 
355  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
356  // curl[1] 3/1
357  Vmath::Vadd(nPts, tmp4, 1, pVarDer[2][0], 1, tmp4, 1);
358  // curl[1]-1/3
359  Vmath::Vcopy(nPts, pVarDer[0][2], 1, tmp3, 1);
360  Vmath::Neg(nPts, tmp3, 1);
361  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
362  // square curl[1]
363  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
364  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
365 
366  tmp4 = Array<OneD, NekDouble>(nPts, 0.0);
367  // curl[2] 1/2
368  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
369  // curl[2]-2/1
370  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
371  Vmath::Neg(nPts, tmp3, 1);
372  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
373  // square curl[2]
374  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
375  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
376 
377  }
378  else if ( nDim > 1 )
379  {
380  // curl[2] 1/2
381  Vmath::Vadd(nPts, tmp4, 1, pVarDer[0][1], 1, tmp4, 1);
382  // curl[2]-2/1
383  Vmath::Vcopy(nPts, pVarDer[1][0], 1, tmp3, 1);
384  Vmath::Neg(nPts, tmp3, 1);
385  Vmath::Vadd(nPts, tmp4, 1, tmp3, 1, tmp4, 1);
386  // square curl[2]
387  Vmath::Vmul(nPts, tmp4, 1, tmp4, 1, tmp3, 1);
388  Vmath::Vadd(nPts, m_curlVelSquare, 1, tmp3, 1, m_curlVelSquare, 1);
389  }
390  else
391  {
392  Vmath::Fill(nPts, 0.0, m_curlVelSquare, 1);
393  }
394  }
Array< OneD, NekDouble > m_divVelSquare
Definition: Diffusion.h:160
Array< OneD, NekDouble > m_divVel
Params for Ducros sensor.
Definition: Diffusion.h:159
Array< OneD, NekDouble > m_curlVelSquare
Definition: Diffusion.h:161
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:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:1199

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 434 of file Diffusion.h.

435  {
436  return v_GetFluxTensor();
437  }
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor()
Definition: Diffusion.h:603

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 460 of file Diffusion.h.

461  {
462  return v_GetTraceNormal();
463  }
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal()
Definition: Diffusion.cpp:236

References v_GetTraceNormal().

◆ InitObject()

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

Definition at line 47 of file Diffusion.cpp.

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

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

◆ SetArtificialDiffusionVector()

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

Definition at line 399 of file Diffusion.h.

401  {
402  m_ArtificialDiffusionVector = std::bind(
403  func, obj, std::placeholders::_1, std::placeholders::_2);
404  }

References m_ArtificialDiffusionVector.

◆ SetCalcViscosity()

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

Definition at line 414 of file Diffusion.h.

415  {
416  m_CalcViscosity = std::bind(
417  func, obj, std::placeholders::_1, std::placeholders::_2);
418  }

References m_CalcViscosity.

◆ SetDiffusionFluxCons() [1/2]

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

Definition at line 378 of file Diffusion.h.

379  {
381  }
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:470

References m_FunctorDiffusionfluxCons.

◆ SetDiffusionFluxCons() [2/2]

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

Definition at line 369 of file Diffusion.h.

370  {
371  m_FunctorDiffusionfluxCons = std::bind(
372  func, obj, std::placeholders::_1, std::placeholders::_2,
373  std::placeholders::_3, std::placeholders::_4,
374  std::placeholders::_5, std::placeholders::_6,
375  std::placeholders::_7);
376  }

References m_FunctorDiffusionfluxCons.

◆ SetDiffusionFluxConsTrace() [1/2]

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

Definition at line 393 of file Diffusion.h.

394  {
396  }
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:471

References m_FunctorDiffusionfluxConsTrace.

◆ SetDiffusionFluxConsTrace() [2/2]

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

Definition at line 384 of file Diffusion.h.

385  {
387  func, obj, std::placeholders::_1, std::placeholders::_2,
388  std::placeholders::_3, std::placeholders::_4,
389  std::placeholders::_5, std::placeholders::_6,
390  std::placeholders::_7);
391  }

References m_FunctorDiffusionfluxConsTrace.

◆ SetDiffusionSymmFluxCons()

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

Definition at line 421 of file Diffusion.h.

422  {
423  m_FunctorSymmetricfluxCons = std::bind(
424  func, obj, std::placeholders::_1, std::placeholders::_2,
425  std::placeholders::_3, std::placeholders::_4,
426  std::placeholders::_5, std::placeholders::_6);
427  }
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:474

References m_FunctorSymmetricfluxCons.

◆ SetFluxPenaltyNS() [1/2]

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

Definition at line 363 of file Diffusion.h.

364  {
365  m_fluxPenaltyNS = flux;
366  }
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:468

References m_fluxPenaltyNS.

◆ SetFluxPenaltyNS() [2/2]

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

Definition at line 356 of file Diffusion.h.

357  {
358  m_fluxPenaltyNS = std::bind(
359  func, obj, std::placeholders::_1, std::placeholders::_2,
360  std::placeholders::_3);
361  }

References m_fluxPenaltyNS.

◆ SetFluxVector() [1/2]

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

Definition at line 336 of file Diffusion.h.

338  {
339  m_fluxVector = fluxVector;
340  }
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:466

References m_fluxVector.

◆ SetFluxVector() [2/2]

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

Definition at line 329 of file Diffusion.h.

330  {
331  m_fluxVector = std::bind(
332  func, obj, std::placeholders::_1, std::placeholders::_2,
333  std::placeholders::_3);
334  }

References m_fluxVector.

◆ SetFluxVectorNS() [1/2]

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

Definition at line 350 of file Diffusion.h.

351  {
352  m_fluxVectorNS = fluxVector;
353  }
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:467

References m_fluxVectorNS.

◆ SetFluxVectorNS() [2/2]

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

Definition at line 343 of file Diffusion.h.

344  {
345  m_fluxVectorNS = std::bind(
346  func, obj, std::placeholders::_1, std::placeholders::_2,
347  std::placeholders::_3);
348  }

References m_fluxVectorNS.

◆ SetHomoDerivs()

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

Definition at line 429 of file Diffusion.h.

430  {
431  v_SetHomoDerivs(deriv);
432  }
virtual void v_SetHomoDerivs(Array< OneD, Array< OneD, NekDouble > > &deriv)
Definition: Diffusion.h:597

References v_SetHomoDerivs().

◆ SetSpecialBndTreat()

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

Definition at line 407 of file Diffusion.h.

408  {
409  m_SpecialBndTreat = std::bind(
410  func, obj, std::placeholders::_1);
411  }
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:472

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

Definition at line 568 of file Diffusion.h.

577  {
578  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
579  VolumeFlux, outarray, pFwd, pBwd);
580 
581  }

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

Definition at line 582 of file Diffusion.h.

591  {
592  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
593  VolumeFlux, outarray, pFwd, pBwd);
594 
595  }

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

Definition at line 242 of file Diffusion.cpp.

249  {
250  boost::ignore_unused(nConvectiveFields, npnts, vFwd, vBwd,
251  aver, jump);
252  NEKERROR(ErrorUtil::efatal, "v_ConsVarAveJump not defined");
253  }
#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::Diffusion3DHomogeneous1D, Nektar::DiffusionLDGNS, Nektar::SolverUtils::DiffusionLFRNS, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 182 of file Diffusion.cpp.

189  {
190  boost::ignore_unused(nConvectiveFields, fields, inarray, outarray,
191  pFwd, pBwd);
192  NEKERROR(ErrorUtil::efatal, "v_Diffuse not defined");
193  }

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, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 255 of file Diffusion.cpp.

261  {
262  boost::ignore_unused(fields, inarray, qfields,
263  pFwd, pBwd);
264  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
265  }

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::DiffusionLDGNS, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 195 of file Diffusion.cpp.

202  {
203  boost::ignore_unused(nConvectiveFields, fields, inarray,
204  outarray, pFwd, pBwd);
205  NEKERROR(ErrorUtil::efatal, "v_DiffuseCoeffs not defined");
206  }

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

Definition at line 208 of file Diffusion.cpp.

217  {
218  boost::ignore_unused(nConvectiveFields, fields, inarray,
219  outarray, vFwd, vBwd,qfield,nonZeroIndex);
220  NEKERROR(ErrorUtil::efatal,"v_DiffuseCoeffs not defined");
221  }

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, and Nektar::SolverUtils::DiffusionLDG.

Definition at line 279 of file Diffusion.cpp.

288  {
289  boost::ignore_unused(fields, inarray, qfields, VolumeFlux,
290  TraceFlux, pFwd, pBwd, nonZeroIndex);
291  NEKERROR(ErrorUtil::efatal, "Not defined function DiffuseTraceFLux.");
292  }

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 549 of file Diffusion.h.

561  {
562  boost::ignore_unused(nConvectiveFields, fields, inarray, qfield,
563  VolumeFlux, SymmFlux, pFwd, pBwd,
564  nonZeroIndex,solution_Aver,solution_jump);
565  ASSERTL0(false, "Not defined function v_DiffuseTraceSymmFlux.");
566  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216

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 267 of file Diffusion.cpp.

273  {
274  boost::ignore_unused(fields, inarray, qfields, VolumeFlux,
275  nonZeroIndex);
276  NEKERROR(ErrorUtil::efatal, "Not defined for function DiffuseVolumeFLux.");
277  }

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 603 of file Diffusion.h.

604  {
605  static TensorOfArray3D<NekDouble> tmp;
606  return tmp;
607  }

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 298 of file Diffusion.cpp.

302  {
303 
304  int nDim = fields[0]->GetCoordim(0);
305  int nPts = fields[0]->GetTotPoints();
306  for(int i = 0; i < nDim; ++i)
307  {
308  primVar[i] = Array<OneD, NekDouble>(nPts, 0.0);
309  Vmath::Vdiv(nPts, inarray[i+1], 1, inarray[0], 1, primVar[i], 1);
310  }
311  }
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:257

References Vmath::Vdiv().

◆ v_GetTraceNormal()

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

Reimplemented in Nektar::SolverUtils::DiffusionIP.

Definition at line 236 of file Diffusion.cpp.

237  {
238  NEKERROR(ErrorUtil::efatal,"v_GetTraceNormal not defined");
240  }
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 597 of file Diffusion.h.

599  {
600  boost::ignore_unused(deriv);
601  }

Referenced by SetHomoDerivs().

Member Data Documentation

◆ m_ArtificialDiffusionVector

DiffusionArtificialDiffusion Nektar::SolverUtils::Diffusion::m_ArtificialDiffusionVector
protected

◆ m_CalcViscosity

CalcViscosity Nektar::SolverUtils::Diffusion::m_CalcViscosity
protected

Definition at line 473 of file Diffusion.h.

Referenced by SetCalcViscosity().

◆ m_curlVelSquare

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

Definition at line 161 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 159 of file Diffusion.h.

Referenced by GetDivCurl(), and InitObject().

◆ m_divVelSquare

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

Definition at line 160 of file Diffusion.h.

Referenced by GetDivCurl(), and InitObject().

◆ m_fluxPenaltyNS

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

Definition at line 468 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 476 of file Diffusion.h.

Referenced by Diffuse(), and DiffuseCoeffs().