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

#include <DiffusionLDG.h>

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

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionLDG ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
virtual void v_Diffuse (const std::size_t nConvective, 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 void v_DiffuseCoeffs (const std::size_t nConvective, 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 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 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 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...
 
void NumFluxforScalar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
void ApplyScalarBCs (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
 
void NumFluxforVector (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to penalty term. More...
 
void ApplyVectorBCs (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
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_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

std::string m_shockCaptureType
 
NekDouble m_C11
 Coefficient of penalty term. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
LibUtilities::SessionReaderSharedPtr m_session
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
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
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void 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 inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 

Detailed Description

Definition at line 46 of file DiffusionLDG.h.

Constructor & Destructor Documentation

◆ DiffusionLDG()

Nektar::SolverUtils::DiffusionLDG::DiffusionLDG ( )
protected

Definition at line 49 of file DiffusionLDG.cpp.

50 {
51 }

Referenced by create().

Member Function Documentation

◆ ApplyScalarBCs()

void Nektar::SolverUtils::DiffusionLDG::ApplyScalarBCs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const std::size_t  var,
const Array< OneD, const NekDouble > &  ufield,
const Array< OneD, const NekDouble > &  Fwd,
const Array< OneD, const NekDouble > &  Bwd,
Array< OneD, NekDouble > &  penaltyflux 
)
protected

Definition at line 276 of file DiffusionLDG.cpp.

282 {
283  boost::ignore_unused(ufield, Bwd);
284  // Number of boundary regions
285  std::size_t nBndRegions =
286  fields[var]->GetBndCondExpansions().size();
287  std::size_t cnt = 0;
288  for (std::size_t i = 0; i < nBndRegions; ++i)
289  {
290  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
292  {
293  continue;
294  }
295 
296  // Number of boundary expansion related to that region
297  std::size_t nBndEdges =
298  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
299 
300  // Weakly impose boundary conditions by modifying flux values
301  for (std::size_t e = 0; e < nBndEdges; ++e)
302  {
303  std::size_t nBndEdgePts = fields[var]
304  ->GetBndCondExpansions()[i]
305  ->GetExp(e)
306  ->GetTotPoints();
307 
308  std::size_t id1 =
309  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
310 
311  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
312  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
313  cnt++));
314 
315  // AV boundary conditions
316  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
317  GetUserDefined(),"Wall") ||
318  boost::iequals(fields[var]->GetBndConditions()[i]->
319  GetUserDefined(),"Symmetry") ||
320  boost::iequals(fields[var]->GetBndConditions()[i]->
321  GetUserDefined(),"WallViscous") ||
322  boost::iequals(fields[var]->GetBndConditions()[i]->
323  GetUserDefined(),"WallAdiabatic"))
324  {
325  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
326  }
327  // For Dirichlet boundary condition: uflux = g_D
328  else if (fields[var]
329  ->GetBndConditions()[i]
330  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
331  {
332  Vmath::Vcopy(
333  nBndEdgePts,
334  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
335  1, &penaltyflux[id2], 1);
336  }
337  // For Neumann boundary condition: uflux = u+
338  else if ((fields[var]->GetBndConditions()[i])
339  ->GetBoundaryConditionType() ==
341  {
342  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
343  }
344  }
345  }
346 }
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, and Vmath::Vcopy().

Referenced by NumFluxforScalar().

◆ ApplyVectorBCs()

void Nektar::SolverUtils::DiffusionLDG::ApplyVectorBCs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const std::size_t  var,
const std::size_t  dir,
const Array< OneD, const NekDouble > &  qfield,
const Array< OneD, const NekDouble > &  qFwd,
const Array< OneD, const NekDouble > &  qBwd,
Array< OneD, NekDouble > &  penaltyflux 
)
protected

Diffusion: Imposing weak boundary condition for q with flux uflux = g_D on Dirichlet boundary condition uflux = u_Fwd on Neumann boundary condition

Definition at line 413 of file DiffusionLDG.cpp.

420 {
421  boost::ignore_unused(qfield, qBwd);
422 
423  std::size_t nBndRegions =
424  fields[var]->GetBndCondExpansions().size();
425  std::size_t cnt = 0;
426 
427  for (std::size_t i = 0; i < nBndRegions; ++i)
428  {
429  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
431  {
432  continue;
433  }
434  std::size_t nBndEdges =
435  fields[var]->GetBndCondExpansions()[i]->GetExpSize();
436 
437  // Weakly impose boundary conditions by modifying flux values
438  for (std::size_t e = 0; e < nBndEdges; ++e)
439  {
440  std::size_t nBndEdgePts = fields[var]
441  ->GetBndCondExpansions()[i]
442  ->GetExp(e)
443  ->GetTotPoints();
444 
445  std::size_t id1 =
446  fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
447 
448  std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
449  fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
450  cnt++));
451 
452  // AV boundary conditions
453  if ( boost::iequals(fields[var]->GetBndConditions()[i]->
454  GetUserDefined(),"Wall") ||
455  boost::iequals(fields[var]->GetBndConditions()[i]->
456  GetUserDefined(),"Symmetry") ||
457  boost::iequals(fields[var]->GetBndConditions()[i]->
458  GetUserDefined(),"WallViscous") ||
459  boost::iequals(fields[var]->GetBndConditions()[i]->
460  GetUserDefined(),"WallAdiabatic"))
461  {
462  Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
463  }
464  // For Dirichlet boundary condition:
465  // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
466  else if (fields[var]
467  ->GetBndConditions()[i]
468  ->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
469  {
470  Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
471  &qFwd[id2], 1, &penaltyflux[id2], 1);
472  }
473  // For Neumann boundary condition: qflux = g_N
474  else if ((fields[var]->GetBndConditions()[i])
475  ->GetBoundaryConditionType() ==
477  {
478  Vmath::Vmul(
479  nBndEdgePts, &m_traceNormals[dir][id2], 1,
480  &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
481  1, &penaltyflux[id2], 1);
482  }
483  }
484  }
485 }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:65
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, m_traceNormals, Vmath::Vmul(), and Vmath::Zero().

Referenced by NumFluxforVector().

◆ create()

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

Definition at line 49 of file DiffusionLDG.h.

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

References DiffusionLDG().

◆ NumFluxforScalar()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforScalar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
TensorOfArray3D< NekDouble > &  uflux,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
protected

Definition at line 227 of file DiffusionLDG.cpp.

233 {
234  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
235  std::size_t nvariables = fields.size();
236  std::size_t nDim = fields[0]->GetCoordim(0);
237 
238  Array<OneD, NekDouble> Fwd{nTracePts};
239  Array<OneD, NekDouble> Bwd{nTracePts};
240  Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
241 
242  // Get the sign of (v \cdot n), v = an arbitrary vector
243  // Evaluate upwind flux:
244  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
245  for (std::size_t i = 0; i < nvariables; ++i)
246  {
247  // Compute Fwd and Bwd value of ufield of i direction
248  if (pFwd == NullNekDoubleArrayOfArray ||
250  {
251  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
252  }
253  else
254  {
255  Fwd = pFwd[i];
256  Bwd = pBwd[i];
257  }
258 
259  // Upwind
260  Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
261 
262  // Imposing weak boundary condition with flux
263  if (fields[0]->GetBndCondExpansions().size())
264  {
265  ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
266  }
267 
268  for (std::size_t j = 0; j < nDim; ++j)
269  {
270  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
271  uflux[j][i], 1);
272  }
273  }
274 }
void ApplyScalarBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References ApplyScalarBCs(), m_traceNormals, Nektar::NullNekDoubleArrayOfArray, Vmath::Vcopy(), and Vmath::Vmul().

Referenced by v_DiffuseCalcDerivative().

◆ NumFluxforVector()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforVector ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
protected

Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to penalty term.

Definition at line 352 of file DiffusionLDG.cpp.

357 {
358  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
359  std::size_t nvariables = fields.size();
360  std::size_t nDim = qfield.size();
361 
362  Array<OneD, NekDouble> Fwd{nTracePts};
363  Array<OneD, NekDouble> Bwd{nTracePts};
364  Array<OneD, NekDouble> qFwd{nTracePts};
365  Array<OneD, NekDouble> qBwd{nTracePts};
366  Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
367  Array<OneD, NekDouble> uterm{nTracePts};
368 
369  // Evaulate upwind flux:
370  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
371  for (std::size_t i = 0; i < nvariables; ++i)
372  {
373  // Generate Stability term = - C11 ( u- - u+ )
374  fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
375  Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
376  Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
377 
378  qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
379  for (std::size_t j = 0; j < nDim; ++j)
380  {
381  // Compute Fwd and Bwd value of ufield of jth direction
382  fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
383 
384  // Downwind
385  Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
386 
387  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
388  qfluxtemp, 1);
389 
390  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
391  Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
392 
393  // Imposing weak boundary condition with flux
394  if (fields[0]->GetBndCondExpansions().size())
395  {
396  ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
397  qfluxtemp);
398  }
399 
400  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
401  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
402  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
403  Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
404  }
405  }
406 }
void ApplyVectorBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:63
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:372

References ApplyVectorBCs(), m_C11, m_traceNormals, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_DiffuseTraceFlux().

◆ v_Diffuse()

void Nektar::SolverUtils::DiffusionLDG::v_Diffuse ( const std::size_t  nConvective,
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 from Nektar::SolverUtils::Diffusion.

Definition at line 76 of file DiffusionLDG.cpp.

83 {
84  std::size_t nCoeffs = fields[0]->GetNcoeffs();
85 
86  Array<OneD, Array<OneD, NekDouble> > tmp{nConvectiveFields};
87  for (std::size_t i=0; i < nConvectiveFields; ++i)
88  {
89  tmp[i] = Array<OneD, NekDouble> {nCoeffs, 0.0};
90  }
91 
92  DiffusionLDG::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, tmp,
93  pFwd, pBwd);
94 
95  for (std::size_t i = 0; i < nConvectiveFields; ++i)
96  {
97  fields[i]->BwdTrans (tmp[i], outarray[i]);
98  }
99 }
virtual void v_DiffuseCoeffs(const std::size_t nConvective, 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)

References v_DiffuseCoeffs().

◆ v_DiffuseCalcDerivative()

void Nektar::SolverUtils::DiffusionLDG::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 from Nektar::SolverUtils::Diffusion.

Definition at line 164 of file DiffusionLDG.cpp.

170 {
171  std::size_t nConvectiveFields = fields.size();
172  std::size_t nDim = fields[0]->GetCoordim(0);
173  std::size_t nCoeffs = fields[0]->GetNcoeffs();
174  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
175 
176  Array<OneD, NekDouble> tmp{nCoeffs};
177  TensorOfArray3D<NekDouble> flux {nDim};
178  for (std::size_t j = 0; j < nDim; ++j)
179  {
180  flux[j] = Array<OneD, Array<OneD, NekDouble> > {nConvectiveFields};
181  for (std::size_t i = 0; i < nConvectiveFields; ++i)
182  {
183  flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
184  }
185  }
186 
187  NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
188 
189  for (std::size_t j = 0; j < nDim; ++j)
190  {
191  for (std::size_t i = 0; i < nConvectiveFields; ++i)
192  {
193  fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
194  Vmath::Neg (nCoeffs, tmp, 1);
195  fields[i]->AddTraceIntegral (flux[j][i], tmp);
196  fields[i]->SetPhysState (false);
197  fields[i]->MultiplyByElmtInvMass(tmp, tmp);
198  fields[i]->BwdTrans (tmp, qfield[j][i]);
199  }
200  }
201 }
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461

References Vmath::Neg(), and NumFluxforScalar().

◆ v_DiffuseCoeffs()

void Nektar::SolverUtils::DiffusionLDG::v_DiffuseCoeffs ( const std::size_t  nConvective,
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 from Nektar::SolverUtils::Diffusion.

Definition at line 101 of file DiffusionLDG.cpp.

108 {
109  std::size_t nDim = fields[0]->GetCoordim(0);
110  std::size_t nPts = fields[0]->GetTotPoints();
111  std::size_t nCoeffs = fields[0]->GetNcoeffs();
112  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
113 
114  Array<OneD, NekDouble> tmp{nCoeffs};
115 
116  TensorOfArray3D<NekDouble> qfield{nDim};
117  for (std::size_t j = 0; j < nDim; ++j)
118  {
119  qfield[j] = Array<OneD, Array<OneD, NekDouble> > {nConvectiveFields};
120  for (std::size_t i = 0; i < nConvectiveFields; ++i)
121  {
122  qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
123  }
124  }
125 
126  Array<OneD, Array<OneD, NekDouble > > traceflux{nConvectiveFields};
127  for (std::size_t i = 0; i < nConvectiveFields; ++i)
128  {
129  traceflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
130  }
131 
132  DiffuseCalcDerivative(fields, inarray, qfield, pFwd, pBwd);
133 
134  // Initialize viscous tensor
135  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> viscTensor{nDim};
136  for (std::size_t j = 0; j < nDim; ++j)
137  {
138  viscTensor[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
139  for (std::size_t i = 0; i < nConvectiveFields; ++i)
140  {
141  viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
142  }
143  }
144  DiffuseVolumeFlux(fields, inarray, qfield, viscTensor);
145  DiffuseTraceFlux(fields, inarray, qfield, viscTensor, traceflux, pFwd, pBwd);
146 
147  Array<OneD, Array<OneD, NekDouble> > qdbase{nDim};
148 
149  for (std::size_t i = 0; i < nConvectiveFields; ++i)
150  {
151  for (std::size_t j = 0; j < nDim; ++j)
152  {
153  qdbase[j] = viscTensor[j][i];
154  }
155  fields[i]->IProductWRTDerivBase(qdbase, tmp);
156 
157  Vmath::Neg (nCoeffs, tmp, 1);
158  fields[i]->AddTraceIntegral (traceflux[i], tmp);
159  fields[i]->SetPhysState (false);
160  fields[i]->MultiplyByElmtInvMass(tmp, outarray[i]);
161  }
162 }
SOLVER_UTILS_EXPORT void DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
Definition: Diffusion.h:249
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
Definition: Diffusion.h:236
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.cpp:224

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

Referenced by v_Diffuse().

◆ v_DiffuseTraceFlux()

void Nektar::SolverUtils::DiffusionLDG::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 from Nektar::SolverUtils::Diffusion.

Definition at line 213 of file DiffusionLDG.cpp.

222 {
223  boost::ignore_unused(qfield, pFwd, pBwd, nonZeroIndex);
224  NumFluxforVector(fields, inarray, viscTensor, TraceFlux);
225 }
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to p...

References NumFluxforVector().

◆ v_DiffuseVolumeFlux()

void Nektar::SolverUtils::DiffusionLDG::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 from Nektar::SolverUtils::Diffusion.

Definition at line 203 of file DiffusionLDG.cpp.

209 {
210  boost::ignore_unused(fields, nonZeroIndex);
211  m_fluxVector(inarray, qfield, viscTensor);
212 }
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:466

References Nektar::SolverUtils::Diffusion::m_fluxVector.

◆ v_InitObject()

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

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 53 of file DiffusionLDG.cpp.

56 {
57  m_session = pSession;
58 
59  m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
60 
61  // Set up penalty term for LDG
62  m_session->LoadParameter("LDGc11", m_C11, 1.0);
63 
64  // Setting up the normals
65  std::size_t nDim = pFields[0]->GetCoordim(0);
66  std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
67 
68  m_traceNormals = Array<OneD, Array<OneD, NekDouble>>{nDim};
69  for (std::size_t i = 0; i < nDim; ++i)
70  {
71  m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts};
72  }
73  pFields[0]->GetTrace()->GetNormals(m_traceNormals);
74 }
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:66

References m_C11, m_session, m_shockCaptureType, and m_traceNormals.

Member Data Documentation

◆ m_C11

NekDouble Nektar::SolverUtils::DiffusionLDG::m_C11
protected

Coefficient of penalty term.

Definition at line 63 of file DiffusionLDG.h.

Referenced by NumFluxforVector(), and v_InitObject().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionLDG::m_session
protected

Definition at line 66 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_shockCaptureType

std::string Nektar::SolverUtils::DiffusionLDG::m_shockCaptureType
protected

Definition at line 60 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_traceNormals

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

Definition at line 65 of file DiffusionLDG.h.

Referenced by ApplyVectorBCs(), NumFluxforScalar(), NumFluxforVector(), and v_InitObject().

◆ type

std::string Nektar::SolverUtils::DiffusionLDG::type
static
Initial value:
=
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionLDG.h:49
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

Definition at line 55 of file DiffusionLDG.h.