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

Non Smooth artificial diffusion for shock capture for compressible flow problems. More...

#include <NonSmoothShockCapture.h>

Inheritance diagram for Nektar::NonSmoothShockCapture:
[legend]

Static Public Member Functions

static ArtificialDiffusionSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

virtual void v_GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu)
 
- Protected Member Functions inherited from Nektar::ArtificialDiffusion
 ArtificialDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
 Constructor. More...
 
virtual void v_DoArtificialDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
 Return the flux vector for the artificial viscosity operator. More...
 

Private Member Functions

 NonSmoothShockCapture (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
 
virtual ~NonSmoothShockCapture (void)
 

Private Attributes

int m_offset
 Parameters. More...
 

Friends

class MemoryManager< NonSmoothShockCapture >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::ArtificialDiffusion
virtual ~ArtificialDiffusion ()
 
void DoArtificialDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Apply the artificial diffusion. More...
 
void GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu)
 Calculate the artificial viscosity. More...
 
void SetElmtHP (const Array< OneD, NekDouble > &hOverP)
 Set h/p scaling. More...
 
- Protected Attributes inherited from Nektar::ArtificialDiffusion
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array of fields. More...
 
VariableConverterSharedPtr m_varConv
 Auxiliary object to convert variables. More...
 
SolverUtils::DiffusionSharedPtr m_diffusion
 LDG Diffusion operator. More...
 
NekDouble m_mu0
 Constant scaling. More...
 
Array< OneD, NekDoublem_hOverP
 h/p scaling More...
 

Detailed Description

Non Smooth artificial diffusion for shock capture for compressible flow problems.

Definition at line 48 of file NonSmoothShockCapture.h.

Constructor & Destructor Documentation

◆ NonSmoothShockCapture()

Nektar::NonSmoothShockCapture::NonSmoothShockCapture ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const int  spacedim 
)
private

Definition at line 47 of file NonSmoothShockCapture.cpp.

References m_offset, and Nektar::ArtificialDiffusion::m_session.

51  : ArtificialDiffusion(pSession, pFields, spacedim)
52 {
53  m_session->LoadParameter ("SensorOffset", m_offset, 1);
54 }
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
ArtificialDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
Constructor.

◆ ~NonSmoothShockCapture()

virtual Nektar::NonSmoothShockCapture::~NonSmoothShockCapture ( void  )
inlineprivatevirtual

Definition at line 81 of file NonSmoothShockCapture.h.

References m_offset.

81 {};

Member Function Documentation

◆ create()

static ArtificialDiffusionSharedPtr Nektar::NonSmoothShockCapture::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const int  spacedim 
)
inlinestatic

Creates an instance of this class.

Definition at line 55 of file NonSmoothShockCapture.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

59  {
62  AllocateSharedPtr(pSession, pFields, spacedim);
63  return p;
64  }
std::shared_ptr< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ v_GetArtificialViscosity()

void Nektar::NonSmoothShockCapture::v_GetArtificialViscosity ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, NekDouble > &  mu 
)
protectedvirtual

Implements Nektar::ArtificialDiffusion.

Definition at line 56 of file NonSmoothShockCapture.cpp.

References Nektar::ArtificialDiffusion::m_fields, Nektar::ArtificialDiffusion::m_hOverP, Nektar::ArtificialDiffusion::m_mu0, m_offset, Nektar::ArtificialDiffusion::m_varConv, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmax().

59 {
60  int nPts = m_fields[0]->GetTotPoints();
61  int nElements = m_fields[0]->GetExpSize();
62 
63  // Determine the maximum wavespeed
64  Array <OneD, NekDouble > Lambda(nPts, 0.0);
65  Array <OneD, NekDouble > soundspeed(nPts, 0.0);
66  Array <OneD, NekDouble > absVelocity(nPts, 0.0);
67  m_varConv->GetSoundSpeed(physfield, soundspeed);
68  m_varConv->GetAbsoluteVelocity(physfield, absVelocity);
69  Vmath::Vadd(nPts, absVelocity, 1, soundspeed, 1, Lambda, 1);
70 
71  // Compute sensor based on rho
72  Array<OneD, NekDouble> Sensor(nPts, 0.0);
73  m_varConv->GetSensor(m_fields[0], physfield, Sensor, mu, m_offset);
74 
75  // Scale AV
76  Array<OneD, NekDouble> tmp;
77  for (int e = 0; e < nElements; e++)
78  {
79  int physOffset = m_fields[0]->GetPhys_Offset(e);
80  int nElmtPoints = m_fields[0]->GetExp(e)->GetTotPoints();
81 
82  // Get max wavespeed per element
83  NekDouble LambdaElmt = 0.0;
84  LambdaElmt = Vmath::Vmax(nElmtPoints, tmp = Lambda + physOffset, 1);
85 
86  // Scale viscosity by the maximum wave speed and h/p
87  LambdaElmt *= m_mu0 * m_hOverP[e];
88  Vmath::Smul(nElmtPoints, LambdaElmt, tmp = mu + physOffset, 1,
89  tmp = mu + physOffset, 1);
90  }
91 }
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:782
Array< OneD, NekDouble > m_hOverP
h/p scaling
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
double NekDouble
NekDouble m_mu0
Constant scaling.
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
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:302

Friends And Related Function Documentation

◆ MemoryManager< NonSmoothShockCapture >

friend class MemoryManager< NonSmoothShockCapture >
friend

Definition at line 52 of file NonSmoothShockCapture.h.

Member Data Documentation

◆ className

std::string Nektar::NonSmoothShockCapture::className
static
Initial value:
RegisterCreatorFunction("NonSmooth",
"NonSmooth artificial diffusion for shock capture.")

Name of the class.

Definition at line 67 of file NonSmoothShockCapture.h.

◆ m_offset

int Nektar::NonSmoothShockCapture::m_offset
private

Parameters.

Definition at line 81 of file NonSmoothShockCapture.h.

Referenced by NonSmoothShockCapture(), v_GetArtificialViscosity(), and ~NonSmoothShockCapture().