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

Wall boundary conditions for compressible flow problems. More...

#include <IsentropicVortexBC.h>

Inheritance diagram for Nektar::IsentropicVortexBC:
[legend]

Static Public Member Functions

static CFSBndCondSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 Creates an instance of this class. More...
 

Static Public Attributes

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

Protected Member Functions

void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
 
- Protected Member Functions inherited from Nektar::CFSBndCond
 CFSBndCond (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 Constructor. More...
 
virtual void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)=0
 
virtual void v_ApplyBwdWeight ()
 

Private Member Functions

 IsentropicVortexBC (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 
void EvaluateIsentropicVortex (const Array< OneD, NekDouble > &x, const Array< OneD, NekDouble > &y, const Array< OneD, NekDouble > &z, Array< OneD, Array< OneD, NekDouble > > &u, NekDouble time, const int o=0)
 
 ~IsentropicVortexBC (void) override
 

Friends

class MemoryManager< IsentropicVortexBC >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::CFSBndCond
virtual ~CFSBndCond ()
 
void Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time=0)
 Apply the boundary condition. More...
 
void ApplyBwdWeight ()
 Apply the Weight of boundary condition. More...
 
- Protected Attributes inherited from Nektar::CFSBndCond
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array of fields. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Trace normals. More...
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 Grid Velocity. More...
 
int m_spacedim
 Space dimension. More...
 
VariableConverterSharedPtr m_varConv
 Auxiliary object to convert variables. More...
 
NekDouble m_diffusionAveWeight
 Weight for average calculation of diffusion term. More...
 
NekDouble m_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
NekDouble m_pOut
 
Array< OneD, NekDoublem_velInf
 
NekDouble m_angVel
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Detailed Description

Wall boundary conditions for compressible flow problems.

Definition at line 46 of file IsentropicVortexBC.h.

Constructor & Destructor Documentation

◆ IsentropicVortexBC()

Nektar::IsentropicVortexBC::IsentropicVortexBC ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, Array< OneD, NekDouble > > &  pTraceNormals,
const Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
private

Definition at line 47 of file IsentropicVortexBC.cpp.

53 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
54 bcRegion, cnt)
55{
56}
CFSBndCond(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Constructor.
Definition: CFSBndCond.cpp:47

◆ ~IsentropicVortexBC()

Nektar::IsentropicVortexBC::~IsentropicVortexBC ( void  )
inlineoverrideprivate

Definition at line 88 of file IsentropicVortexBC.h.

88{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::IsentropicVortexBC::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, Array< OneD, NekDouble > > &  pTraceNormals,
const Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file IsentropicVortexBC.h.

58 {
61 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
62 bcRegion, cnt);
63 return p;
64 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:51

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

◆ EvaluateIsentropicVortex()

void Nektar::IsentropicVortexBC::EvaluateIsentropicVortex ( const Array< OneD, NekDouble > &  x,
const Array< OneD, NekDouble > &  y,
const Array< OneD, NekDouble > &  z,
Array< OneD, Array< OneD, NekDouble > > &  u,
NekDouble  time,
const int  o = 0 
)
private

Definition at line 101 of file IsentropicVortexBC.cpp.

105{
106 int nq = x.size();
107
108 // Flow parameters
109 const NekDouble x0 = 5.0;
110 const NekDouble y0 = 0.0;
111 const NekDouble beta = 5.0;
112 const NekDouble u0 = 1.0;
113 const NekDouble v0 = 0.5;
114 const NekDouble gamma = m_gamma;
115 NekDouble r, xbar, ybar, tmp;
116 NekDouble fac = 1.0 / (16.0 * gamma * M_PI * M_PI);
117
118 // In 3D zero rhow field.
119 if (m_spacedim == 3)
120 {
121 Vmath::Zero(nq, &u[3][o], 1);
122 }
123
124 // Fill storage
125 for (int i = 0; i < nq; ++i)
126 {
127 xbar = x[i] - u0 * time - x0;
128 ybar = y[i] - v0 * time - y0;
129 r = sqrt(xbar * xbar + ybar * ybar);
130 tmp = beta * exp(1 - r * r);
131 u[0][i + o] =
132 pow(1.0 - (gamma - 1.0) * tmp * tmp * fac, 1.0 / (gamma - 1.0));
133 u[1][i + o] = u[0][i + o] * (u0 - tmp * ybar / (2 * M_PI));
134 u[2][i + o] = u[0][i + o] * (v0 + tmp * xbar / (2 * M_PI));
135 u[m_spacedim + 1][i + o] =
136 pow(u[0][i + o], gamma) / (gamma - 1.0) +
137 0.5 * (u[1][i + o] * u[1][i + o] + u[2][i + o] * u[2][i + o]) /
138 u[0][i + o];
139 }
140}
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:105
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::beta, Nektar::CFSBndCond::m_gamma, Nektar::CFSBndCond::m_spacedim, tinysimd::sqrt(), and Vmath::Zero().

Referenced by v_Apply().

◆ v_Apply()

void Nektar::IsentropicVortexBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 58 of file IsentropicVortexBC.cpp.

61{
62 int nvariables = physarray.size();
63
64 const Array<OneD, const int> &bndTraceMap = m_fields[0]->GetTraceBndMap();
65 // loop over Boundary Regions
66 int npoints, id1, id2, e_max;
67
68 e_max = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
69
70 for (int e = 0; e < e_max; ++e)
71 {
72 npoints = m_fields[0]
73 ->GetBndCondExpansions()[m_bcRegion]
74 ->GetExp(e)
75 ->GetTotPoints();
76 id1 =
77 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
78 id2 =
79 m_fields[0]->GetTrace()->GetPhys_Offset(bndTraceMap[m_offset + e]);
80
81 Array<OneD, NekDouble> x(npoints, 0.0);
82 Array<OneD, NekDouble> y(npoints, 0.0);
83 Array<OneD, NekDouble> z(npoints, 0.0);
84
85 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
86 x, y, z);
87
88 EvaluateIsentropicVortex(x, y, z, Fwd, time, id2);
89
90 for (int i = 0; i < nvariables; ++i)
91 {
92 Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
93 &(m_fields[i]
94 ->GetBndCondExpansions()[m_bcRegion]
95 ->UpdatePhys())[id1],
96 1);
97 }
98 }
99}
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
int m_offset
Offset.
Definition: CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
void EvaluateIsentropicVortex(const Array< OneD, NekDouble > &x, const Array< OneD, NekDouble > &y, const Array< OneD, NekDouble > &z, Array< OneD, Array< OneD, NekDouble > > &u, NekDouble time, const int o=0)
std::vector< double > z(NPUPPER)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References EvaluateIsentropicVortex(), Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_offset, Vmath::Vcopy(), and Nektar::UnitTests::z().

Friends And Related Function Documentation

◆ MemoryManager< IsentropicVortexBC >

friend class MemoryManager< IsentropicVortexBC >
friend

Definition at line 1 of file IsentropicVortexBC.h.

Member Data Documentation

◆ className

std::string Nektar::IsentropicVortexBC::className
static
Initial value:
=
"IsentropicVortex", IsentropicVortexBC::create,
"Isentropic vortex boundary condition.")
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Name of the class.

Definition at line 67 of file IsentropicVortexBC.h.