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

Outflow characteristic boundary conditions for compressible flow problems. More...

#include <EnforceEntropyTotalEnthalpy.h>

Inheritance diagram for Nektar::EnforceEntropyTotalEnthalpy:
[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 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
 
void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 
void FromToRotation (Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat)
 
SOLVER_UTILS_EXPORT void rotateToNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void rotateFromNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
- 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 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 ()
 

Protected Attributes

int m_npts
 
Array< OneD, NekDoublem_rhoBC
 
Array< OneD, Array< OneD, NekDouble > > m_velBC
 
Array< OneD, NekDoublem_pBC
 
Array< OneD, NekDoublem_VnInf
 Reference normal velocity. More...
 
Array< OneD, int > m_bndToTraceMap
 
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
 
- 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...
 
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
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Private Member Functions

 EnforceEntropyTotalEnthalpy (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const int pSpaceDim, const int bcRegion, const int cnt)
 
 ~EnforceEntropyTotalEnthalpy (void) override
 

Friends

class MemoryManager< EnforceEntropyTotalEnthalpy >
 

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...
 

Detailed Description

Outflow characteristic boundary conditions for compressible flow problems.

Definition at line 51 of file EnforceEntropyTotalEnthalpy.h.

Constructor & Destructor Documentation

◆ EnforceEntropyTotalEnthalpy()

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

Definition at line 51 of file EnforceEntropyTotalEnthalpy.cpp.

56 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
57{
58
60 m_fields[0]->GetBndCondExpansions()[m_bcRegion];
61
62 //-> Gather a list of index from trace to this boundary
63 m_npts = bndexp->GetTotPoints();
64
65 m_bndToTraceMap = Array<OneD, int>(m_npts, -1);
66
67 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
68
69 // Construct a map for the boundary to trace map for easy acess to
70 // phys space points
71 int cnt1 = 0;
72 for (int e = 0; e < bndexp->GetNumElmts(); ++e)
73 {
74 int nTracePts = bndexp->GetExp(e)->GetTotPoints();
75
76 int id =
77 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
78
79 // Loop on the points of the m_bcRegion
80 for (int i = 0; i < nTracePts; i++)
81 {
82 // the ith point in region e
83 m_bndToTraceMap[cnt1++] = id + i;
84 }
85 }
86
87 Array<OneD, Array<OneD, NekDouble>> BCvals(m_fields.size());
88 m_bndPhys = Array<OneD, Array<OneD, NekDouble>>(m_fields.size());
89
90 for (int i = 0; i < m_fields.size(); ++i)
91 {
92 m_bndPhys[i] =
93 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys();
94
95 BCvals[i] = Array<OneD, NekDouble>(m_npts);
96 Vmath::Vcopy(m_npts, m_bndPhys[i], 1, BCvals[i], 1);
97 }
98
99 // Set up boudnary required BCs
100 m_rhoBC = Array<OneD, NekDouble>(m_npts);
101 Vmath::Vcopy(m_npts, BCvals[0], 1, m_rhoBC, 1);
102
103 m_velBC = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
104 // Evaluate velocity on boundary
105 for (int i = 0; i < m_spacedim; ++i)
106 {
107 m_velBC[i] = Array<OneD, NekDouble>(m_npts);
108 Vmath::Vcopy(m_npts, BCvals[i + 1], 1, m_velBC[i], 1);
109 Vmath::Vdiv(m_npts, m_velBC[i], 1, m_rhoBC, 1, m_velBC[i], 1);
110 }
111 m_pBC = Array<OneD, NekDouble>(m_npts);
112 m_varConv->GetPressure(BCvals, m_pBC);
113
114 // Computing the normal velocity for characteristics coming
115 // from outside the computational domain
116 m_VnInf = Array<OneD, NekDouble>(m_npts, 0.0);
117 for (int i = 0; i < m_spacedim; i++)
118 {
119 for (int j = 0; j < m_npts; ++j)
120 {
121 m_VnInf[j] += m_traceNormals[i][m_bndToTraceMap[j]] * m_velBC[i][j];
122 }
123 }
124}
CFSBndCond(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const int pSpaceDim, const int bcRegion, const int cnt)
Constructor.
Definition: CFSBndCond.cpp:47
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:97
int m_offset
Offset.
Definition: CFSBndCond.h:111
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:91
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
Array< OneD, Array< OneD, NekDouble > > m_velBC
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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.hpp:126
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::CFSBndCond::m_bcRegion, m_bndPhys, m_bndToTraceMap, Nektar::CFSBndCond::m_fields, m_npts, Nektar::CFSBndCond::m_offset, m_pBC, m_rhoBC, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, m_velBC, m_VnInf, Vmath::Vcopy(), and Vmath::Vdiv().

◆ ~EnforceEntropyTotalEnthalpy()

Nektar::EnforceEntropyTotalEnthalpy::~EnforceEntropyTotalEnthalpy ( void  )
inlineoverrideprivate

Definition at line 118 of file EnforceEntropyTotalEnthalpy.h.

118{};

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 57 of file EnforceEntropyTotalEnthalpy.h.

62 {
65 pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt);
66 return p;
67 }
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.

◆ FromToRotation()

void Nektar::EnforceEntropyTotalEnthalpy::FromToRotation ( Array< OneD, const NekDouble > &  from,
Array< OneD, const NekDouble > &  to,
NekDouble mat 
)
protected

◆ GenerateRotationMatrices()

void Nektar::EnforceEntropyTotalEnthalpy::GenerateRotationMatrices ( const Array< OneD, const Array< OneD, NekDouble > > &  normals)
protected

◆ rotateFromNormal()

SOLVER_UTILS_EXPORT void Nektar::EnforceEntropyTotalEnthalpy::rotateFromNormal ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
const Array< OneD, const Array< OneD, NekDouble > > &  normals,
const Array< OneD, const Array< OneD, NekDouble > > &  vecLocs,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protected

◆ rotateToNormal()

SOLVER_UTILS_EXPORT void Nektar::EnforceEntropyTotalEnthalpy::rotateToNormal ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
const Array< OneD, const Array< OneD, NekDouble > > &  normals,
const Array< OneD, const Array< OneD, NekDouble > > &  vecLocs,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protected

◆ v_Apply()

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

Implements Nektar::CFSBndCond.

Definition at line 126 of file EnforceEntropyTotalEnthalpy.cpp.

130{
131 int i, j;
132 int nDimensions = m_spacedim;
133
134 Array<OneD, Array<OneD, NekDouble>> FwdBnd(Fwd.size());
135 Array<OneD, Array<OneD, NekDouble>> bndPhys(Fwd.size());
136
137 // make a local copy of Fwd along boundary of interest
138 for (i = 0; i < Fwd.size(); ++i)
139 {
140 FwdBnd[i] = Array<OneD, NekDouble>(m_npts);
141 for (j = 0; j < m_npts; ++j)
142 {
143 FwdBnd[i][j] = Fwd[i][m_bndToTraceMap[j]];
144 }
145 }
146
147 // Computing the normal velocity for characteristics coming
148 // from inside the computational domain
149 Array<OneD, NekDouble> Vn(m_npts, 0.0);
150
151 for (i = 0; i < nDimensions; ++i)
152 {
153 for (j = 0; j < m_npts; ++j)
154 {
155 Vn[j] += m_traceNormals[i][m_bndToTraceMap[j]] * FwdBnd[i + 1][j];
156 }
157 }
158 // divide by density.
159 Vmath::Vdiv(m_npts, Vn, 1, FwdBnd[0], 1, Vn, 1);
160
161 // Get speed of sound
162 Array<OneD, NekDouble> pressure(m_npts);
163 Array<OneD, NekDouble> soundSpeed(m_npts);
164
165 m_varConv->GetPressure(FwdBnd, pressure);
166 m_varConv->GetSoundSpeed(FwdBnd, soundSpeed);
167
168 // Get Mach. Note: it is computed by Vn/c
169 Array<OneD, NekDouble> Mach(m_npts, 0.0);
170 Vmath::Vdiv(m_npts, Vn, 1, soundSpeed, 1, Mach, 1);
171 Vmath::Vabs(m_npts, Mach, 1, Mach, 1);
172
173 // Auxiliary variables
174 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
175
176 // L represents properties outside boundary
177 // R represents properties inside boundary (numerical state)
178 NekDouble rhoL, uL, pL;
179 NekDouble EBC, rR, cstar, pstar, rhostar, ustar, Sstar, vn;
180
181 NekDouble gamMinOne = m_gamma - 1.0;
182 NekDouble twoOverGamMinOne = 2.0 / gamMinOne;
183 NekDouble gamInv = 1.0 / m_gamma;
184
185 NekDouble tmp1, tmp2;
186 // Loop on m_bcRegions
187 for (int pnt = 0; pnt < m_npts; ++pnt)
188 {
189 // Impose inflow Riemann invariant
190 if (Vn[pnt] <= 0.0)
191 {
192 // Subsonic flows
193 if (Mach[pnt] < 1.00)
194 {
195 // right characteristic
196 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
197 twoOverGamMinOne;
198 vn = -m_VnInf[pnt]; // vn BC
199
200 // fix total entropy and entropy to be input values
201 // compute cstar and ustar using left-pointint characteristic
202 // IR^-
203 tmp1 = twoOverGamMinOne * rR;
204 tmp1 = tmp1 * tmp1;
205 tmp2 = rR * rR - vn * vn -
206 twoOverGamMinOne * m_gamma * m_pBC[pnt] / m_rhoBC[pnt];
207 tmp2 =
208 (twoOverGamMinOne * twoOverGamMinOne + twoOverGamMinOne) *
209 tmp2;
210 cstar = -twoOverGamMinOne * rR + sqrt(tmp1 - tmp2);
211 cstar = cstar / (twoOverGamMinOne * twoOverGamMinOne +
212 twoOverGamMinOne);
213 ustar = rR + twoOverGamMinOne * cstar;
214 Sstar = m_pBC[pnt] / pow(m_rhoBC[pnt], m_gamma);
215 rhostar = pow(cstar * cstar / (m_gamma * Sstar),
216 0.5 * twoOverGamMinOne);
217 pstar = rhostar * cstar * cstar * gamInv;
218
219 // add supplement equation that rhoL=rhostar
220 // then pL=pstar, according to IL^0
221 // and uL=ustar, according to IL^+
222 rhoL = rhostar;
223 pL = pstar;
224 uL = ustar;
225 }
226 else // Supersonic inflow
227 {
228 // all characteristics are from left so just impose
229 // star state to left values
230 // Note: m_vnInf is the negative of the normal velocity
231 // across boundary
232 rhoL = m_rhoBC[pnt];
233 uL = -m_VnInf[pnt];
234 pL = m_pBC[pnt];
235 }
236
237 // Boundary energy
238 EBC = pL * twoOverGamMinOne * 0.5;
239
240 // evaluate the different between the left state normal
241 // velocity and that from the desired condition (note
242 // m_VnInf is using an outwards normal definition.
243 NekDouble VnDiff = uL + m_VnInf[pnt];
244
245 // Boundary velocities & Kinite energy
246 // Note: normals are negated since they point outwards in
247 // the domain
248
249 // Note: Can just use the BC values directly!!
250 for (j = 0; j < nDimensions; ++j)
251 {
252 // Set velocity to the desired conditions modified to
253 // take account of the normal state for Riemann
254 // problem. (Negative accounts for outwards normal definition)
255 // velBC[j] = m_velBC[j][pnt];
256 velBC[j] = m_velBC[j][pnt] -
257 VnDiff * m_traceNormals[j][m_bndToTraceMap[pnt]];
258
259 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
260 }
261
262 //-------------------------------------------------------------------------
263 // Impose Left hand Riemann Invariant boundary conditions
264 m_bndPhys[0][pnt] = rhoL;
265 for (j = 0; j < nDimensions; ++j)
266 {
267 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
268 }
269 m_bndPhys[nDimensions + 1][pnt] = EBC;
270 }
271 else // Outflow
272 {
273
274 // Note: Allowing the switch can cause worse convergence in this
275 // type BC.
276 // So improve it later.
277 if (Mach[pnt] < 1.00)
278 {
279 // subsonic outflow: fix pstar
280 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
281 twoOverGamMinOne;
282 // vn = -m_VnInf[pnt];
283
284 pstar = m_pBC[pnt];
285 rhostar = FwdBnd[0][pnt] * pow((pstar / pressure[pnt]), gamInv);
286 cstar = sqrt(m_gamma * pstar / rhostar);
287 ustar = rR + cstar * twoOverGamMinOne;
288
289 rhoL = rhostar;
290 uL = ustar;
291 pL = pstar;
292 }
293 else
294 {
295 // supersonic outflow
296 // Just set to imposed state and let Riemann BC dictate values
297 rhoL = m_rhoBC[pnt];
298 uL = -m_VnInf[pnt];
299 pL = m_pBC[pnt];
300 }
301
302 // Boundary energy
303 EBC = pL * twoOverGamMinOne * 0.5;
304
305 // Boundary velocities & Kinite energy
306 // Note: normals are negated since they point outwards in
307 // the domain
308 for (j = 0; j < nDimensions; ++j)
309 {
310 velBC[j] = -1.0 * uL * m_traceNormals[j][m_bndToTraceMap[pnt]];
311 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
312 }
313
314 // Impose Left hand Riemann Invariant boundary conditions
315 m_bndPhys[0][pnt] = rhoL;
316 for (j = 0; j < nDimensions; ++j)
317 {
318 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
319 }
320 m_bndPhys[nDimensions + 1][pnt] = EBC;
321 }
322 }
323}
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
double NekDouble
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References m_bndPhys, m_bndToTraceMap, Nektar::CFSBndCond::m_gamma, m_npts, m_pBC, m_rhoBC, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, m_velBC, m_VnInf, CG_Iterations::pressure, tinysimd::sqrt(), Vmath::Vabs(), and Vmath::Vdiv().

Friends And Related Function Documentation

◆ MemoryManager< EnforceEntropyTotalEnthalpy >

Definition at line 1 of file EnforceEntropyTotalEnthalpy.h.

Member Data Documentation

◆ className

std::string Nektar::EnforceEntropyTotalEnthalpy::className
static
Initial value:
=
"EnforceEntropyTotalEnthalpy", EnforceEntropyTotalEnthalpy::create,
"Riemann invariant boundary condition, \
fixing H and S.")
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, 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.
Definition: NekFactory.hpp:197
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Name of the class.

Definition at line 70 of file EnforceEntropyTotalEnthalpy.h.

◆ m_bndPhys

Array<OneD, Array<OneD, NekDouble> > Nektar::EnforceEntropyTotalEnthalpy::m_bndPhys
protected

Definition at line 87 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_bndToTraceMap

Array<OneD, int> Nektar::EnforceEntropyTotalEnthalpy::m_bndToTraceMap
protected

Definition at line 83 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_npts

int Nektar::EnforceEntropyTotalEnthalpy::m_npts
protected

Definition at line 73 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_pBC

Array<OneD, NekDouble> Nektar::EnforceEntropyTotalEnthalpy::m_pBC
protected

Definition at line 79 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_rhoBC

Array<OneD, NekDouble> Nektar::EnforceEntropyTotalEnthalpy::m_rhoBC
protected

Definition at line 75 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_velBC

Array<OneD, Array<OneD, NekDouble> > Nektar::EnforceEntropyTotalEnthalpy::m_velBC
protected

Definition at line 77 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().

◆ m_VnInf

Array<OneD, NekDouble> Nektar::EnforceEntropyTotalEnthalpy::m_VnInf
protected

Reference normal velocity.

Definition at line 81 of file EnforceEntropyTotalEnthalpy.h.

Referenced by EnforceEntropyTotalEnthalpy(), and v_Apply().