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

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

#include <EnforceEntropyPressure.h>

Inheritance diagram for Nektar::EnforceEntropyPressure:
[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
 
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 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 ()
 

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

Private Member Functions

 EnforceEntropyPressure (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)
 
 ~EnforceEntropyPressure (void) override
 

Friends

class MemoryManager< EnforceEntropyPressure >
 

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

Constructor & Destructor Documentation

◆ EnforceEntropyPressure()

Nektar::EnforceEntropyPressure::EnforceEntropyPressure ( 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 50 of file EnforceEntropyPressure.cpp.

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

◆ ~EnforceEntropyPressure()

Nektar::EnforceEntropyPressure::~EnforceEntropyPressure ( void  )
inlineoverrideprivate

Definition at line 121 of file EnforceEntropyPressure.h.

121{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::EnforceEntropyPressure::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 57 of file EnforceEntropyPressure.h.

63 {
66 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
67 bcRegion, cnt);
68 return p;
69 }
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::EnforceEntropyPressure::FromToRotation ( Array< OneD, const NekDouble > &  from,
Array< OneD, const NekDouble > &  to,
NekDouble mat 
)
protected

◆ GenerateRotationMatrices()

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

◆ rotateFromNormal()

SOLVER_UTILS_EXPORT void Nektar::EnforceEntropyPressure::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::EnforceEntropyPressure::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::EnforceEntropyPressure::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 127 of file EnforceEntropyPressure.cpp.

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

friend class MemoryManager< EnforceEntropyPressure >
friend

Definition at line 1 of file EnforceEntropyPressure.h.

Member Data Documentation

◆ className

std::string Nektar::EnforceEntropyPressure::className
static
Initial value:
=
"EnforceEntropyPressure", EnforceEntropyPressure::create,
"Riemann boundary condition enforcing entropy and pressure.")
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 72 of file EnforceEntropyPressure.h.

◆ m_bndPhys

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

Definition at line 89 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_bndToTraceMap

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

Definition at line 85 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_npts

int Nektar::EnforceEntropyPressure::m_npts
protected

Definition at line 75 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_pBC

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

Definition at line 81 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_rhoBC

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

Definition at line 77 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_velBC

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

Definition at line 79 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().

◆ m_VnInf

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

Reference normal velocity.

Definition at line 83 of file EnforceEntropyPressure.h.

Referenced by EnforceEntropyPressure(), and v_Apply().