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

void v_GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu) override
 
- 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)
 
virtual void v_DoArtificialDiffusionCoeff (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual void v_GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu)=0
 
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)
 
 ~NonSmoothShockCapture (void) override
 

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 DoArtificialDiffusionCoeff (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Apply the artificial diffusion the outarray is in coeff space. More...
 
void GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu)
 Calculate the artificial viscosity. 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...
 

Detailed Description

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

Definition at line 47 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.

51 : ArtificialDiffusion(pSession, pFields, spacedim)
52{
53 m_session->LoadParameter("SensorOffset", m_offset, 1);
54
55 // init h/p
56 m_varConv->SetElmtMinHP(m_fields);
57}
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
ArtificialDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
Constructor.
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.

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

◆ ~NonSmoothShockCapture()

Nektar::NonSmoothShockCapture::~NonSmoothShockCapture ( void  )
inlineoverrideprivate

Definition at line 78 of file NonSmoothShockCapture.h.

78{};

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 53 of file NonSmoothShockCapture.h.

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

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

◆ v_GetArtificialViscosity()

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

Implements Nektar::ArtificialDiffusion.

Definition at line 62 of file NonSmoothShockCapture.cpp.

65{
66 int nPts = m_fields[0]->GetTotPoints();
67 int nElements = m_fields[0]->GetExpSize();
68
69 // Determine the maximum wavespeed
70 Array<OneD, NekDouble> Lambda(nPts, 0.0);
71 Array<OneD, NekDouble> soundspeed(nPts, 0.0);
72 Array<OneD, NekDouble> absVelocity(nPts, 0.0);
73 m_varConv->GetSoundSpeed(physfield, soundspeed);
74 m_varConv->GetAbsoluteVelocity(physfield, absVelocity);
75 Vmath::Vadd(nPts, absVelocity, 1, soundspeed, 1, Lambda, 1);
76
77 // Compute sensor based on rho
78 Array<OneD, NekDouble> Sensor(nPts, 0.0);
79 m_varConv->GetSensor(m_fields[0], physfield, Sensor, mu, m_offset);
80
81 // Scale AV
82 Array<OneD, NekDouble> tmp;
83 for (int e = 0; e < nElements; e++)
84 {
85 int physOffset = m_fields[0]->GetPhys_Offset(e);
86 int nElmtPoints = m_fields[0]->GetExp(e)->GetTotPoints();
87
88 // Get max wavespeed per element
89 NekDouble LambdaElmt = 0.0;
90 LambdaElmt = Vmath::Vmax(nElmtPoints, tmp = Lambda + physOffset, 1);
91
92 // Scale viscosity by the maximum wave speed and h/p
93 LambdaElmt *= m_mu0 * m_varConv->GetElmtMinHP()[e];
94 Vmath::Smul(nElmtPoints, LambdaElmt, tmp = mu + physOffset, 1,
95 tmp = mu + physOffset, 1);
96 }
97}
NekDouble m_mu0
Constant scaling.
double NekDouble
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.hpp:180
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.hpp:100
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.hpp:644

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

Friends And Related Function Documentation

◆ MemoryManager< NonSmoothShockCapture >

friend class MemoryManager< NonSmoothShockCapture >
friend

Definition at line 1 of file NonSmoothShockCapture.h.

Member Data Documentation

◆ className

std::string Nektar::NonSmoothShockCapture::className
static
Initial value:
=
"NonSmooth artificial diffusion for shock capture.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static ArtificialDiffusionSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
Creates an instance of this class.
ArtificialDiffusionFactory & GetArtificialDiffusionFactory()
Declaration of the artificial diffusion factory singleton.

Name of the class.

Definition at line 65 of file NonSmoothShockCapture.h.

◆ m_offset

int Nektar::NonSmoothShockCapture::m_offset
private

Parameters.

Definition at line 81 of file NonSmoothShockCapture.h.

Referenced by NonSmoothShockCapture(), and v_GetArtificialViscosity().