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

#include <DiffusionLDGNS.h>

Inheritance diagram for Nektar::DiffusionLDGNS:
[legend]

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionLDGNS ()
 
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)
 Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations: More...
 
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 NumericalFluxO1 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &numericalFluxO1, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Builds the numerical flux for the 1st order derivatives. More...
 
void ApplyBCsO1 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, Array< OneD, NekDouble > > &flux01)
 Imposes appropriate bcs for the 1st order derivatives. More...
 
void NumericalFluxO2 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Build the numerical flux for the 2nd order derivatives. More...
 
void ApplyBCsO2 (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)
 Imposes appropriate bcs for the 2nd order derivatives. More...
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual void v_GetPrimVar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
 Compute primary variables. More...
 
- 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 SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal ()
 
SOLVER_UTILS_EXPORT void GetDivCurl (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
 Compute divergence and curl squared. More...
 

Protected Attributes

NekDouble m_C11
 Penalty coefficient for LDGNS. More...
 
Array< OneD, NekDoublem_traceOneOverH
 h scaling for penalty term More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceVel
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
LibUtilities::SessionReaderSharedPtr m_session
 
NekDouble m_Twall
 
EquationOfStateSharedPtr m_eos
 Equation of system for computing temperature. More...
 
TensorOfArray3D< NekDoublem_viscTensor
 
Array< OneD, Array< OneD, NekDouble > > m_homoDerivs
 
std::size_t m_spaceDim
 
std::size_t m_diffDim
 
- 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 50 of file DiffusionLDGNS.h.

Constructor & Destructor Documentation

◆ DiffusionLDGNS()

Nektar::DiffusionLDGNS::DiffusionLDGNS ( )
protected

Definition at line 49 of file DiffusionLDGNS.cpp.

50 {
51 }

Member Function Documentation

◆ ApplyBCsO1()

void Nektar::DiffusionLDGNS::ApplyBCsO1 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd,
Array< OneD, Array< OneD, NekDouble > > &  flux01 
)
protected

Imposes appropriate bcs for the 1st order derivatives.

Definition at line 426 of file DiffusionLDGNS.cpp.

432 {
433  boost::ignore_unused(pBwd);
434 
435  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
436  std::size_t nScalars = inarray.size();
437 
438  Array<OneD, NekDouble> tmp1{nTracePts, 0.0};
439  Array<OneD, NekDouble> tmp2{nTracePts, 0.0};
440  Array<OneD, NekDouble> Tw{nTracePts, m_Twall};
441 
442  Array< OneD, Array<OneD, NekDouble > > scalarVariables{nScalars};
443 
444  // Extract internal values of the scalar variables for Neumann bcs
445  for (std::size_t i = 0; i < nScalars; ++i)
446  {
447  scalarVariables[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
448  }
449 
450  // Compute boundary conditions for velocities
451  for (std::size_t i = 0; i < nScalars-1; ++i)
452  {
453  // Note that cnt has to loop on nBndRegions and nBndEdges
454  // and has to be reset to zero per each equation
455  std::size_t cnt = 0;
456  std::size_t nBndRegions = fields[i+1]->
457  GetBndCondExpansions().size();
458  for (std::size_t j = 0; j < nBndRegions; ++j)
459  {
460  if (fields[i+1]->GetBndConditions()[j]->
461  GetBoundaryConditionType() == SpatialDomains::ePeriodic)
462  {
463  continue;
464  }
465 
466  std::size_t nBndEdges = fields[i+1]->
467  GetBndCondExpansions()[j]->GetExpSize();
468  for (std::size_t e = 0; e < nBndEdges; ++e)
469  {
470  std::size_t nBndEdgePts = fields[i+1]->
471  GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
472 
473  std::size_t id1 = fields[i+1]->
474  GetBndCondExpansions()[j]->GetPhys_Offset(e);
475 
476  std::size_t id2 = fields[0]->GetTrace()->
477  GetPhys_Offset(fields[0]->GetTraceMap()->
478  GetBndCondIDToGlobalTraceID(cnt++));
479 
480  if (boost::iequals(fields[i]->GetBndConditions()[j]->
481  GetUserDefined(),"WallViscous") ||
482  boost::iequals(fields[i]->GetBndConditions()[j]->
483  GetUserDefined(),"WallAdiabatic"))
484  {
485  // Reinforcing bcs for velocity in case of Wall bcs
486  Vmath::Zero(nBndEdgePts, &scalarVariables[i][id2], 1);
487 
488  }
489  else if (
490  boost::iequals(fields[i]->GetBndConditions()[j]->
491  GetUserDefined(),"Wall") ||
492  boost::iequals(fields[i]->GetBndConditions()[j]->
493  GetUserDefined(),"Symmetry"))
494  {
495  // Symmetry bc: normal velocity is zero
496  // get all velocities at once because we need u.n
497  if (i==0)
498  {
499  // tmp1 = -(u.n)
500  Vmath::Zero(nBndEdgePts, tmp1, 1);
501  for (std::size_t k = 0; k < nScalars-1; ++k)
502  {
503  Vmath::Vdiv(nBndEdgePts,
504  &(fields[k+1]->GetBndCondExpansions()[j]->
505  UpdatePhys())[id1], 1,
506  &(fields[0]->GetBndCondExpansions()[j]->
507  UpdatePhys())[id1], 1,
508  &scalarVariables[k][id2], 1);
509  Vmath::Vvtvp(nBndEdgePts,
510  &m_traceNormals[k][id2], 1,
511  &scalarVariables[k][id2], 1,
512  &tmp1[0], 1,
513  &tmp1[0], 1);
514  }
515  Vmath::Smul(nBndEdgePts, -1.0,
516  &tmp1[0], 1,
517  &tmp1[0], 1);
518 
519  // u_i - (u.n)n_i
520  for (std::size_t k = 0; k < nScalars-1; ++k)
521  {
522  Vmath::Vvtvp(nBndEdgePts,
523  &tmp1[0], 1,
524  &m_traceNormals[k][id2], 1,
525  &scalarVariables[k][id2], 1,
526  &scalarVariables[k][id2], 1);
527  }
528  }
529  }
530  else if (fields[i]->GetBndConditions()[j]->
531  GetBoundaryConditionType() ==
533  {
534  // Imposing velocity bcs if not Wall
535  Vmath::Vdiv(nBndEdgePts,
536  &(fields[i+1]->GetBndCondExpansions()[j]->
537  UpdatePhys())[id1], 1,
538  &(fields[0]->GetBndCondExpansions()[j]->
539  UpdatePhys())[id1], 1,
540  &scalarVariables[i][id2], 1);
541  }
542 
543  // For Dirichlet boundary condition: uflux = u_bcs
544  if (fields[i]->GetBndConditions()[j]->
545  GetBoundaryConditionType() ==
547  {
548  Vmath::Vcopy(nBndEdgePts,
549  &scalarVariables[i][id2], 1,
550  &fluxO1[i][id2], 1);
551  }
552 
553  // For Neumann boundary condition: uflux = u_+
554  else if ((fields[i]->GetBndConditions()[j])->
555  GetBoundaryConditionType() ==
557  {
558  Vmath::Vcopy(nBndEdgePts,
559  &pFwd[i][id2], 1,
560  &fluxO1[i][id2], 1);
561  }
562 
563  // Building kinetic energy to be used for T bcs
564  Vmath::Vmul(nBndEdgePts,
565  &scalarVariables[i][id2], 1,
566  &scalarVariables[i][id2], 1,
567  &tmp1[id2], 1);
568 
569  Vmath::Smul(nBndEdgePts, 0.5,
570  &tmp1[id2], 1,
571  &tmp1[id2], 1);
572 
573  Vmath::Vadd(nBndEdgePts,
574  &tmp2[id2], 1,
575  &tmp1[id2], 1,
576  &tmp2[id2], 1);
577  }
578  }
579  }
580 
581  // Compute boundary conditions for temperature
582  std::size_t cnt = 0;
583  std::size_t nBndRegions = fields[nScalars]->
584  GetBndCondExpansions().size();
585  for (std::size_t j = 0; j < nBndRegions; ++j)
586  {
587  if (fields[nScalars]->GetBndConditions()[j]->
588  GetBoundaryConditionType() ==
590  {
591  continue;
592  }
593 
594  std::size_t nBndEdges = fields[nScalars]->
595  GetBndCondExpansions()[j]->GetExpSize();
596  for (std::size_t e = 0; e < nBndEdges; ++e)
597  {
598  std::size_t nBndEdgePts = fields[nScalars]->
599  GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
600 
601  std::size_t id1 = fields[nScalars]->
602  GetBndCondExpansions()[j]->GetPhys_Offset(e);
603 
604  std::size_t id2 = fields[0]->GetTrace()->
605  GetPhys_Offset(fields[0]->GetTraceMap()->
606  GetBndCondIDToGlobalTraceID(cnt++));
607 
608  // Imposing Temperature Twall at the wall
609  if (boost::iequals(fields[nScalars]->GetBndConditions()[j]->
610  GetUserDefined(),"WallViscous"))
611  {
612  Vmath::Vcopy(nBndEdgePts,
613  &Tw[0], 1,
614  &scalarVariables[nScalars-1][id2], 1);
615  }
616  // Imposing Temperature through condition on the Energy
617  // for no wall boundaries (e.g. farfield)
618  else if (fields[nScalars]->GetBndConditions()[j]->
619  GetBoundaryConditionType() ==
621  {
622  // Use equation of state to evaluate temperature
623  NekDouble rho, ene;
624  for (std::size_t n = 0; n < nBndEdgePts; ++n)
625  {
626  rho = fields[0]->
627  GetBndCondExpansions()[j]->
628  GetPhys()[id1+n];
629  ene = fields[nScalars]->
630  GetBndCondExpansions()[j]->
631  GetPhys()[id1 +n] / rho - tmp2[id2+n];
632  scalarVariables[nScalars-1][id2+n] =
633  m_eos->GetTemperature(rho, ene);
634  }
635  }
636 
637  // For Dirichlet boundary condition: uflux = u_bcs
638  if (fields[nScalars]->GetBndConditions()[j]->
639  GetBoundaryConditionType() == SpatialDomains::eDirichlet &&
640  !boost::iequals(fields[nScalars]->GetBndConditions()[j]
641  ->GetUserDefined(), "WallAdiabatic"))
642  {
643  Vmath::Vcopy(nBndEdgePts,
644  &scalarVariables[nScalars-1][id2], 1,
645  &fluxO1[nScalars-1][id2], 1);
646 
647  }
648 
649  // For Neumann boundary condition: uflux = u_+
650  else if (((fields[nScalars]->GetBndConditions()[j])->
651  GetBoundaryConditionType() == SpatialDomains::eNeumann) ||
652  boost::iequals(fields[nScalars]->GetBndConditions()[j]->
653  GetUserDefined(), "WallAdiabatic"))
654  {
655  Vmath::Vcopy(nBndEdgePts,
656  &pFwd[nScalars-1][id2], 1,
657  &fluxO1[nScalars-1][id2], 1);
658 
659  }
660  }
661  }
662 }
EquationOfStateSharedPtr m_eos
Equation of system for computing temperature.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
double NekDouble
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 Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513
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 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
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
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, m_eos, m_traceNormals, m_Twall, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by NumericalFluxO1().

◆ ApplyBCsO2()

void Nektar::DiffusionLDGNS::ApplyBCsO2 ( 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

Imposes appropriate bcs for the 2nd order derivatives.

Definition at line 733 of file DiffusionLDGNS.cpp.

741 {
742  boost::ignore_unused(qfield, qBwd);
743 
744  std::size_t cnt = 0;
745  std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
746  // Loop on the boundary regions to apply appropriate bcs
747  for (std::size_t i = 0; i < nBndRegions; ++i)
748  {
749  // Number of boundary regions related to region 'i'
750  std::size_t nBndEdges = fields[var]->
751  GetBndCondExpansions()[i]->GetExpSize();
752 
753  if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType()
755  {
756  continue;
757  }
758 
759  // Weakly impose bcs by modifying flux values
760  for (std::size_t e = 0; e < nBndEdges; ++e)
761  {
762  std::size_t nBndEdgePts = fields[var]->
763  GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
764 
765  std::size_t id2 = fields[0]->GetTrace()->
766  GetPhys_Offset(fields[0]->GetTraceMap()->
767  GetBndCondIDToGlobalTraceID(cnt++));
768 
769  // In case of Dirichlet bcs:
770  // uflux = gD
771  if (fields[var]->GetBndConditions()[i]->
772  GetBoundaryConditionType() == SpatialDomains::eDirichlet
773  && !boost::iequals(fields[var]->GetBndConditions()[i]->
774  GetUserDefined(), "WallAdiabatic"))
775  {
776  Vmath::Vmul(nBndEdgePts,
777  &m_traceNormals[dir][id2], 1,
778  &qFwd[id2], 1,
779  &penaltyflux[id2], 1);
780  }
781  // 3.4) In case of Neumann bcs:
782  // uflux = u+
783  else if ((fields[var]->GetBndConditions()[i])->
784  GetBoundaryConditionType() == SpatialDomains::eNeumann)
785  {
786  ASSERTL0(false,
787  "Neumann bcs not implemented for LDGNS");
788 
789  /*
790  Vmath::Vmul(nBndEdgePts,
791  &m_traceNormals[dir][id2], 1,
792  &(fields[var]->
793  GetBndCondExpansions()[i]->
794  UpdatePhys())[id1], 1,
795  &penaltyflux[id2], 1);
796  */
797  }
798  else if (boost::iequals(fields[var]->GetBndConditions()[i]->
799  GetUserDefined(), "WallAdiabatic"))
800  {
801  if ((var == m_spaceDim + 1))
802  {
803  Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
804  }
805  else
806  {
807 
808  Vmath::Vmul(nBndEdgePts,
809  &m_traceNormals[dir][id2], 1,
810  &qFwd[id2], 1,
811  &penaltyflux[id2], 1);
812 
813  }
814  }
815  }
816  }
817 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216

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

Referenced by NumericalFluxO2().

◆ create()

static DiffusionSharedPtr Nektar::DiffusionLDGNS::create ( std::string  diffType)
inlinestatic

Definition at line 53 of file DiffusionLDGNS.h.

54  {
55  boost::ignore_unused(diffType);
56  return DiffusionSharedPtr(new DiffusionLDGNS());
57  }
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627

◆ NumericalFluxO1()

void Nektar::DiffusionLDGNS::NumericalFluxO1 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
TensorOfArray3D< NekDouble > &  numericalFluxO1,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
protected

Builds the numerical flux for the 1st order derivatives.

Definition at line 388 of file DiffusionLDGNS.cpp.

394 {
395  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
396  std::size_t nScalars = inarray.size();
397 
398  //Upwind
399  Array<OneD, Array<OneD, NekDouble> > numflux{nScalars};
400  for (std::size_t i = 0; i < nScalars; ++i)
401  {
402  numflux[i] = {pFwd[i]};
403  }
404 
405  // Modify the values in case of boundary interfaces
406  if (fields[0]->GetBndCondExpansions().size())
407  {
408  ApplyBCsO1(fields, inarray, pFwd, pBwd, numflux);
409  }
410 
411  // Splitting the numerical flux into the dimensions
412  for (std::size_t j = 0; j < m_spaceDim; ++j)
413  {
414  for (std::size_t i = 0; i < nScalars; ++i)
415  {
416  Vmath::Vmul(nTracePts, m_traceNormals[j], 1,numflux[i], 1,
417  numericalFluxO1[j][i], 1);
418  }
419  }
420 }
void ApplyBCsO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, Array< OneD, NekDouble > > &flux01)
Imposes appropriate bcs for the 1st order derivatives.

References ApplyBCsO1(), m_spaceDim, m_traceNormals, and Vmath::Vmul().

Referenced by v_DiffuseCalcDerivative().

◆ NumericalFluxO2()

void Nektar::DiffusionLDGNS::NumericalFluxO2 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
protected

Build the numerical flux for the 2nd order derivatives.

Definition at line 668 of file DiffusionLDGNS.cpp.

674 {
675  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
676  std::size_t nVariables = fields.size();
677 
678  // Initialize penalty flux
679  Array<OneD, Array<OneD, NekDouble> > fluxPen{nVariables-1};
680  for (std::size_t i = 0; i < nVariables-1; ++i)
681  {
682  fluxPen[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
683  }
684 
685  // Get penalty flux
686  m_fluxPenaltyNS(uFwd, uBwd, fluxPen);
687 
688  // Evaluate Riemann flux
689  // qflux = \hat{q} \cdot u = q \cdot n
690  // Notice: i = 1 (first row of the viscous tensor is zero)
691 
692  Array<OneD, NekDouble > qFwd{nTracePts};
693  Array<OneD, NekDouble > qBwd{nTracePts};
694  Array<OneD, NekDouble > qtemp{nTracePts};
695  Array<OneD, NekDouble > qfluxtemp{nTracePts, 0.0};
696  std::size_t nDim = fields[0]->GetCoordim(0);
697  for (std::size_t i = 1; i < nVariables; ++i)
698  {
699  qflux[i] = Array<OneD, NekDouble> {nTracePts, 0.0};
700  for (std::size_t j = 0; j < nDim; ++j)
701  {
702  // Compute qFwd and qBwd value of qfield in position 'ji'
703  fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
704 
705  // Downwind
706  Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
707 
708  // Multiply the Riemann flux by the trace normals
709  Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
710  qfluxtemp, 1);
711 
712  // Add penalty term
713  Vmath::Vvtvp(nTracePts, m_traceOneOverH, 1, fluxPen[i-1], 1,
714  qfluxtemp, 1, qfluxtemp, 1);
715 
716  // Impose weak boundary condition with flux
717  if (fields[0]->GetBndCondExpansions().size())
718  {
719  ApplyBCsO2(fields, i, j, qfield[j][i], qFwd, qBwd, qfluxtemp);
720  }
721 
722  // Store the final flux into qflux
723  Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
724  }
725  }
726 }
void ApplyBCsO2(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)
Imposes appropriate bcs for the 2nd order derivatives.
Array< OneD, NekDouble > m_traceOneOverH
h scaling for penalty term
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:468

References ApplyBCsO2(), Nektar::SolverUtils::Diffusion::m_fluxPenaltyNS, m_traceNormals, m_traceOneOverH, Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vvtvp().

Referenced by v_DiffuseTraceFlux().

◆ v_Diffuse()

void Nektar::DiffusionLDGNS::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

Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations:

\( \langle\psi, \hat{u}\cdot n\rangle - \langle\nabla\psi \cdot u\rangle \langle\phi, \hat{q}\cdot n\rangle - (\nabla \phi \cdot q) \rangle \)

The equations that need a diffusion operator are those related with the velocities and with the energy.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 210 of file DiffusionLDGNS.cpp.

217 {
218  std::size_t nCoeffs = fields[0]->GetNcoeffs();
219  Array<OneD, Array<OneD, NekDouble> > tmp2{nConvectiveFields};
220  for (std::size_t i = 0; i < nConvectiveFields; ++i)
221  {
222  tmp2[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
223  }
224  v_DiffuseCoeffs(nConvectiveFields, fields, inarray, tmp2, pFwd, pBwd);
225  for (std::size_t i = 0; i < nConvectiveFields; ++i)
226  {
227  fields[i]->BwdTrans (tmp2[i], outarray[i]);
228  }
229 }
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::DiffusionLDGNS::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 306 of file DiffusionLDGNS.cpp.

312 {
313  std::size_t nDim = fields[0]->GetCoordim(0);
314  std::size_t nCoeffs = fields[0]->GetNcoeffs();
315  std::size_t nScalars = inarray.size();
316  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
317  std::size_t nConvectiveFields = fields.size();
318 
319  Array<OneD, NekDouble> tmp1{nCoeffs};
320  Array<OneD, Array<OneD, NekDouble> > tmp2{nConvectiveFields};
321  TensorOfArray3D<NekDouble> numericalFluxO1{m_spaceDim};
322 
323  for (std::size_t j = 0; j < m_spaceDim; ++j)
324  {
325  numericalFluxO1[j] = Array<OneD, Array<OneD, NekDouble> >{nScalars};
326 
327  for (std::size_t i = 0; i < nScalars; ++i)
328  {
329  numericalFluxO1[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
330  }
331  }
332 
333  NumericalFluxO1(fields, inarray, numericalFluxO1, pFwd, pBwd);
334 
335  for (std::size_t j = 0; j < nDim; ++j)
336  {
337  for (std::size_t i = 0; i < nScalars; ++i)
338  {
339  fields[i]->IProductWRTDerivBase (j, inarray[i], tmp1);
340  Vmath::Neg (nCoeffs, tmp1, 1);
341  fields[i]->AddTraceIntegral (numericalFluxO1[j][i], tmp1);
342  fields[i]->SetPhysState (false);
343  fields[i]->MultiplyByElmtInvMass(tmp1, tmp1);
344  fields[i]->BwdTrans (tmp1, qfields[j][i]);
345  }
346  }
347  // For 3D Homogeneous 1D only take derivatives in 3rd direction
348  if (m_diffDim == 1)
349  {
350  for (std::size_t i = 0; i < nScalars; ++i)
351  {
352  qfields[2][i] = m_homoDerivs[i];
353  }
354  }
355 }
Array< OneD, Array< OneD, NekDouble > > m_homoDerivs
void NumericalFluxO1(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &numericalFluxO1, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Builds the numerical flux for the 1st order derivatives.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461

References m_diffDim, m_homoDerivs, m_spaceDim, Vmath::Neg(), and NumericalFluxO1().

◆ v_DiffuseCoeffs()

void Nektar::DiffusionLDGNS::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 231 of file DiffusionLDGNS.cpp.

238 {
239  std::size_t nDim = fields[0]->GetCoordim(0);
240  std::size_t nPts = fields[0]->GetTotPoints();
241  std::size_t nCoeffs = fields[0]->GetNcoeffs();
242  std::size_t nScalars = inarray.size();
243  std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
244 
245 
246  TensorOfArray3D<NekDouble> derivativesO1{m_spaceDim};
247 
248  for (std::size_t j = 0; j < m_spaceDim; ++j)
249  {
250  derivativesO1[j] = Array<OneD, Array<OneD, NekDouble> >{nScalars};
251 
252  for (std::size_t i = 0; i < nScalars; ++i)
253  {
254  derivativesO1[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
255  }
256  }
257 
258  DiffuseCalcDerivative(fields, inarray, derivativesO1, pFwd, pBwd);
259 
260  // Initialisation viscous tensor
261  m_viscTensor = TensorOfArray3D<NekDouble> {m_spaceDim};
262  Array<OneD, Array<OneD, NekDouble> > viscousFlux{nConvectiveFields};
263 
264  for (std::size_t j = 0; j < m_spaceDim; ++j)
265  {
266  m_viscTensor[j] = Array<OneD, Array<OneD, NekDouble> >{nScalars+1};
267  for (std::size_t i = 0; i < nScalars+1; ++i)
268  {
269  m_viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
270  }
271  }
272 
273  for (std::size_t i = 0; i < nConvectiveFields; ++i)
274  {
275  viscousFlux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
276  }
277 
278  DiffuseVolumeFlux(fields, inarray, derivativesO1, m_viscTensor);
279 
280  // Compute u from q_{\eta} and q_{\xi}
281  // Obtain numerical fluxes
282  DiffuseTraceFlux(fields, inarray, derivativesO1, m_viscTensor, viscousFlux,
283  pFwd, pBwd);
284 
285  Array<OneD, NekDouble> tmpOut{nCoeffs};
286  Array<OneD, Array<OneD, NekDouble>> tmpIn{nDim};
287 
288  for (std::size_t i = 0; i < nConvectiveFields; ++i)
289  {
290  // Temporary fix to call collection op
291  // we can reorder m_viscTensor later
292  for (std::size_t j = 0; j < nDim; ++j)
293  {
294  tmpIn[j] = m_viscTensor[j][i];
295  }
296  fields[i]->IProductWRTDerivBase(tmpIn, tmpOut);
297 
298  // Evaulate <\phi, \hat{F}\cdot n> - outarray[i]
299  Vmath::Neg (nCoeffs, tmpOut, 1);
300  fields[i]->AddTraceIntegral (viscousFlux[i], tmpOut);
301  fields[i]->SetPhysState (false);
302  fields[i]->MultiplyByElmtInvMass(tmpOut, outarray[i]);
303  }
304 }
TensorOfArray3D< NekDouble > m_viscTensor
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(), m_spaceDim, m_viscTensor, and Vmath::Neg().

Referenced by v_Diffuse().

◆ v_DiffuseTraceFlux()

void Nektar::DiffusionLDGNS::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 370 of file DiffusionLDGNS.cpp.

379 {
380  boost::ignore_unused(inarray, qfields, nonZeroIndex);
381  NumericalFluxO2(fields, VolumeFlux, TraceFlux, pFwd, pBwd);
382 }
void NumericalFluxO2(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Build the numerical flux for the 2nd order derivatives.

References NumericalFluxO2().

◆ v_DiffuseVolumeFlux()

void Nektar::DiffusionLDGNS::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 357 of file DiffusionLDGNS.cpp.

363 {
364 
365  boost::ignore_unused(fields, nonZeroIndex);
366  m_fluxVectorNS(inarray, qfields, VolumeFlux);
367 
368 }
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:467

References Nektar::SolverUtils::Diffusion::m_fluxVectorNS.

◆ v_GetFluxTensor()

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

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 162 of file DiffusionLDGNS.h.

163  {
164  return m_viscTensor;
165  }

◆ v_GetPrimVar()

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

Compute primary variables.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 823 of file DiffusionLDGNS.cpp.

827 {
828  int nDim = fields[0]->GetCoordim(0);
829  int nPts = fields[0]->GetTotPoints();
830  for(int i = 0; i < nDim; ++i)
831  {
832  primVar[i] = Array<OneD, NekDouble>(nPts, 0.0);
833  primVar[i] = inarray[i];
834  }
835 }

◆ v_InitObject()

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

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 53 of file DiffusionLDGNS.cpp.

56 {
57  m_session = pSession;
58  m_session->LoadParameter ("Twall", m_Twall, 300.15);
59 
60  // Setting up the normals
61  std::size_t nDim = pFields[0]->GetCoordim(0);
62  std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
63 
64  m_spaceDim = nDim;
65  if (pSession->DefinesSolverInfo("HOMOGENEOUS"))
66  {
67  m_spaceDim = 3;
68  }
69 
70  m_diffDim = m_spaceDim - nDim;
71 
72  m_traceVel = Array<OneD, Array<OneD, NekDouble> >{m_spaceDim};
73  m_traceNormals = Array<OneD, Array<OneD, NekDouble> >{m_spaceDim};
74  for (std::size_t i = 0; i < m_spaceDim; ++i)
75  {
76  m_traceVel[i] = Array<OneD, NekDouble> {nTracePts, 0.0};
77  m_traceNormals[i] = Array<OneD, NekDouble> {nTracePts};
78  }
79  pFields[0]->GetTrace()->GetNormals(m_traceNormals);
80 
81  // Create equation of state object
82  std::string eosType;
83  m_session->LoadSolverInfo("EquationOfState",
84  eosType, "IdealGas");
86  .CreateInstance(eosType, m_session);
87 
88 
89  // Set up {h} reference on the trace for penalty term
90  //
91  // Note, this shold be replaced with something smarter when merging
92  // LDG with IP
93 
94  // Get min h per element
95  std::size_t nElements = pFields[0]->GetExpSize();
96  Array<OneD, NekDouble> hEle{nElements, 1.0};
97  for (std::size_t e = 0; e < nElements; e++)
98  {
99  NekDouble h{1.0e+10};
100  std::size_t expDim = pFields[0]->GetShapeDimension();
101  switch(expDim)
102  {
103  case 3:
104  {
106  exp3D = pFields[0]->GetExp(e)->as<LocalRegions::Expansion3D>();
107  for(std::size_t i = 0; i < exp3D->GetNtraces(); ++i)
108  {
109  h = std::min(h, exp3D->GetGeom3D()->GetEdge(i)->GetVertex(0)->
110  dist(*(exp3D->GetGeom3D()->GetEdge(i)->GetVertex(1))));
111  }
112  break;
113  }
114 
115  case 2:
116  {
118  exp2D = pFields[0]->GetExp(e)->as<LocalRegions::Expansion2D>();
119  for(std::size_t i = 0; i < exp2D->GetNtraces(); ++i)
120  {
121  h = std::min(h, exp2D->GetGeom2D()->GetEdge(i)->GetVertex(0)->
122  dist(*(exp2D->GetGeom2D()->GetEdge(i)->GetVertex(1))));
123  }
124  break;
125  }
126  case 1:
127  {
129  exp1D = pFields[0]->GetExp(e)->as<LocalRegions::Expansion1D>();
130 
131  h = std::min(h, exp1D->GetGeom1D()->GetVertex(0)->
132  dist(*(exp1D->GetGeom1D()->GetVertex(1))));
133  break;
134  }
135  default:
136  {
137  ASSERTL0(false,"Dimension out of bound.")
138  }
139  }
140 
141  // Store scaling
142  hEle[e] = h;
143  }
144  // Expand h from elements to points
145  std::size_t nPts = pFields[0]->GetTotPoints();
146  Array<OneD, NekDouble> hElePts{nPts, 0.0};
147  Array<OneD, NekDouble> tmp;
148  for (std::size_t e = 0; e < pFields[0]->GetExpSize(); e++)
149  {
150  std::size_t nElmtPoints = pFields[0]->GetExp(e)->GetTotPoints();
151  std::size_t physOffset = pFields[0]->GetPhys_Offset(e);
152  Vmath::Fill(nElmtPoints, hEle[e],
153  tmp = hElePts + physOffset, 1);
154  }
155  // Get Fwd and Bwd traces
156  Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
157  Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
158  pFields[0]->GetFwdBwdTracePhys(hElePts, Fwd, Bwd);
159  // Fix boundaries
160  std::size_t cnt = 0;
161  std::size_t nBndRegions = pFields[0]->GetBndCondExpansions().size();
162  // Loop on the boundary regions
163  for (std::size_t i = 0; i < nBndRegions; ++i)
164  {
165  // Number of boundary regions related to region 'i'
166  std::size_t nBndEdges = pFields[0]->
167  GetBndCondExpansions()[i]->GetExpSize();
168 
169  if (pFields[0]->GetBndConditions()[i]->GetBoundaryConditionType()
171  {
172  continue;
173  }
174 
175  // Get value from interior
176  for (std::size_t e = 0; e < nBndEdges; ++e)
177  {
178  std::size_t nBndEdgePts = pFields[0]->
179  GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
180 
181  std::size_t id2 = pFields[0]->GetTrace()->
182  GetPhys_Offset(pFields[0]->GetTraceMap()->
183  GetBndCondIDToGlobalTraceID(cnt++));
184 
185  Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &Bwd[id2], 1);
186  }
187  }
188  // Get average of traces
189  Array<OneD, NekDouble> traceH{nTracePts, 1.0};
190  m_traceOneOverH = Array<OneD, NekDouble> {nTracePts, 1.0};
191  Vmath::Svtsvtp(nTracePts, 0.5, Fwd, 1, 0.5, Bwd, 1, traceH, 1);
192  // Multiply by coefficient = - C11 / h
193  m_session->LoadParameter ("LDGNSc11", m_C11, 1.0);
194  Vmath::Sdiv(nTracePts, -m_C11, traceH, 1, m_traceOneOverH, 1);
195 }
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, Array< OneD, NekDouble > > m_traceVel
NekDouble m_C11
Penalty coefficient for LDGNS.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:51
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:49
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:691
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
Definition: Vmath.cpp:291
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::SpatialDomains::ePeriodic, Vmath::Fill(), Nektar::GetEquationOfStateFactory(), Nektar::LocalRegions::Expansion1D::GetGeom1D(), m_C11, m_diffDim, m_eos, m_session, m_spaceDim, m_traceNormals, m_traceOneOverH, m_traceVel, m_Twall, Vmath::Sdiv(), Vmath::Svtsvtp(), and Vmath::Vcopy().

◆ v_SetHomoDerivs()

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

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 156 of file DiffusionLDGNS.h.

158  {
159  m_homoDerivs = deriv;
160  }

Member Data Documentation

◆ m_C11

NekDouble Nektar::DiffusionLDGNS::m_C11
protected

Penalty coefficient for LDGNS.

Definition at line 65 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

◆ m_diffDim

std::size_t Nektar::DiffusionLDGNS::m_diffDim
protected

Definition at line 82 of file DiffusionLDGNS.h.

Referenced by v_DiffuseCalcDerivative(), and v_InitObject().

◆ m_eos

EquationOfStateSharedPtr Nektar::DiffusionLDGNS::m_eos
protected

Equation of system for computing temperature.

Definition at line 75 of file DiffusionLDGNS.h.

Referenced by ApplyBCsO1(), and v_InitObject().

◆ m_homoDerivs

Array<OneD, Array<OneD, NekDouble> > Nektar::DiffusionLDGNS::m_homoDerivs
protected

Definition at line 79 of file DiffusionLDGNS.h.

Referenced by v_DiffuseCalcDerivative().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::DiffusionLDGNS::m_session
protected

Definition at line 72 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

◆ m_spaceDim

std::size_t Nektar::DiffusionLDGNS::m_spaceDim
protected

◆ m_traceNormals

Array<OneD, Array<OneD, NekDouble> > Nektar::DiffusionLDGNS::m_traceNormals
protected

◆ m_traceOneOverH

Array<OneD, NekDouble> Nektar::DiffusionLDGNS::m_traceOneOverH
protected

h scaling for penalty term

Definition at line 68 of file DiffusionLDGNS.h.

Referenced by NumericalFluxO2(), and v_InitObject().

◆ m_traceVel

Array<OneD, Array<OneD, NekDouble> > Nektar::DiffusionLDGNS::m_traceVel
protected

Definition at line 70 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

◆ m_Twall

NekDouble Nektar::DiffusionLDGNS::m_Twall
protected

Definition at line 73 of file DiffusionLDGNS.h.

Referenced by ApplyBCsO1(), and v_InitObject().

◆ m_viscTensor

TensorOfArray3D<NekDouble> Nektar::DiffusionLDGNS::m_viscTensor
protected

Definition at line 77 of file DiffusionLDGNS.h.

Referenced by v_DiffuseCoeffs().

◆ type

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

Definition at line 59 of file DiffusionLDGNS.h.