Nektar++
EnforceEntropyTotalEnthalpy.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: EnforceEntropyTotalEnthalpy.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Modified Riemann invariant boundary condition.
33// Enforce the entropy and total enthalpy at the inflow boundary;
34// Enforce the Riemann invariant at the outflow boundary.
35// The input can be either VALUE or FILE.
36//
37///////////////////////////////////////////////////////////////////////////////
38
40using namespace std;
41
42namespace Nektar
43{
44
47 "EnforceEntropyTotalEnthalpy", EnforceEntropyTotalEnthalpy::create,
48 "Riemann invariant boundary condition, \
49 fixing H and S.");
50
54 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
55 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
56 const int pSpaceDim, const int bcRegion, const int cnt)
57 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
58 bcRegion, cnt)
59{
60
62 m_fields[0]->GetBndCondExpansions()[m_bcRegion];
63
64 //-> Gather a list of index from trace to this boundary
65 m_npts = bndexp->GetTotPoints();
66
68
69 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
70
71 // Construct a map for the boundary to trace map for easy acess to
72 // phys space points
73 int cnt1 = 0;
74 for (int e = 0; e < bndexp->GetNumElmts(); ++e)
75 {
76 int nTracePts = bndexp->GetExp(e)->GetTotPoints();
77
78 int id =
79 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
80
81 // Loop on the points of the m_bcRegion
82 for (int i = 0; i < nTracePts; i++)
83 {
84 // the ith point in region e
85 m_bndToTraceMap[cnt1++] = id + i;
86 }
87 }
88
91
92 for (int i = 0; i < m_fields.size(); ++i)
93 {
94 m_bndPhys[i] =
95 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys();
96
97 BCvals[i] = Array<OneD, NekDouble>(m_npts);
98 Vmath::Vcopy(m_npts, m_bndPhys[i], 1, BCvals[i], 1);
99 }
100
101 // Set up boudnary required BCs
103 Vmath::Vcopy(m_npts, BCvals[0], 1, m_rhoBC, 1);
104
106 // Evaluate velocity on boundary
107 for (int i = 0; i < m_spacedim; ++i)
108 {
110 Vmath::Vcopy(m_npts, BCvals[i + 1], 1, m_velBC[i], 1);
111 Vmath::Vdiv(m_npts, m_velBC[i], 1, m_rhoBC, 1, m_velBC[i], 1);
112 }
114 m_varConv->GetPressure(BCvals, m_pBC);
115
116 // Computing the normal velocity for characteristics coming
117 // from outside the computational domain
119 for (int i = 0; i < m_spacedim; i++)
120 {
121 for (int j = 0; j < m_npts; ++j)
122 {
123 m_VnInf[j] += m_traceNormals[i][m_bndToTraceMap[j]] * m_velBC[i][j];
124 }
125 }
126}
127
130 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &physarray,
131 [[maybe_unused]] const NekDouble &time)
132{
133 int i, j;
134 int nDimensions = m_spacedim;
135
136 Array<OneD, Array<OneD, NekDouble>> FwdBnd(Fwd.size());
137 Array<OneD, Array<OneD, NekDouble>> bndPhys(Fwd.size());
138
139 // make a local copy of Fwd along boundary of interest
140 for (i = 0; i < Fwd.size(); ++i)
141 {
142 FwdBnd[i] = Array<OneD, NekDouble>(m_npts);
143 for (j = 0; j < m_npts; ++j)
144 {
145 FwdBnd[i][j] = Fwd[i][m_bndToTraceMap[j]];
146 }
147 }
148
149 // Computing the normal velocity for characteristics coming
150 // from inside the computational domain
152
153 for (i = 0; i < nDimensions; ++i)
154 {
155 for (j = 0; j < m_npts; ++j)
156 {
157 Vn[j] += m_traceNormals[i][m_bndToTraceMap[j]] * FwdBnd[i + 1][j];
158 }
159 }
160 // divide by density.
161 Vmath::Vdiv(m_npts, Vn, 1, FwdBnd[0], 1, Vn, 1);
162
163 // Get speed of sound
165 Array<OneD, NekDouble> soundSpeed(m_npts);
166
167 m_varConv->GetPressure(FwdBnd, pressure);
168 m_varConv->GetSoundSpeed(FwdBnd, soundSpeed);
169
170 // Get Mach. Note: it is computed by Vn/c
172 Vmath::Vdiv(m_npts, Vn, 1, soundSpeed, 1, Mach, 1);
173 Vmath::Vabs(m_npts, Mach, 1, Mach, 1);
174
175 // Auxiliary variables
176 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
177
178 // L represents properties outside boundary
179 // R represents properties inside boundary (numerical state)
180 NekDouble rhoL, uL, pL;
181 NekDouble EBC, rR, cstar, pstar, rhostar, ustar, Sstar, vn;
182
183 NekDouble gamMinOne = m_gamma - 1.0;
184 NekDouble twoOverGamMinOne = 2.0 / gamMinOne;
185 NekDouble gamInv = 1.0 / m_gamma;
186
187 NekDouble tmp1, tmp2;
188 // Loop on m_bcRegions
189 for (int pnt = 0; pnt < m_npts; ++pnt)
190 {
191 // Impose inflow Riemann invariant
192 if (Vn[pnt] <= 0.0)
193 {
194 // Subsonic flows
195 if (Mach[pnt] < 1.00)
196 {
197 // right characteristic
198 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
199 twoOverGamMinOne;
200 vn = -m_VnInf[pnt]; // vn BC
201
202 // fix total entropy and entropy to be input values
203 // compute cstar and ustar using left-pointint characteristic
204 // IR^-
205 tmp1 = twoOverGamMinOne * rR;
206 tmp1 = tmp1 * tmp1;
207 tmp2 = rR * rR - vn * vn -
208 twoOverGamMinOne * m_gamma * m_pBC[pnt] / m_rhoBC[pnt];
209 tmp2 =
210 (twoOverGamMinOne * twoOverGamMinOne + twoOverGamMinOne) *
211 tmp2;
212 cstar = -twoOverGamMinOne * rR + sqrt(tmp1 - tmp2);
213 cstar = cstar / (twoOverGamMinOne * twoOverGamMinOne +
214 twoOverGamMinOne);
215 ustar = rR + twoOverGamMinOne * cstar;
216 Sstar = m_pBC[pnt] / pow(m_rhoBC[pnt], m_gamma);
217 rhostar = pow(cstar * cstar / (m_gamma * Sstar),
218 0.5 * twoOverGamMinOne);
219 pstar = rhostar * cstar * cstar * gamInv;
220
221 // add supplement equation that rhoL=rhostar
222 // then pL=pstar, according to IL^0
223 // and uL=ustar, according to IL^+
224 rhoL = rhostar;
225 pL = pstar;
226 uL = ustar;
227 }
228 else // Supersonic inflow
229 {
230 // all characteristics are from left so just impose
231 // star state to left values
232 // Note: m_vnInf is the negative of the normal velocity
233 // across boundary
234 rhoL = m_rhoBC[pnt];
235 uL = -m_VnInf[pnt];
236 pL = m_pBC[pnt];
237 }
238
239 // Boundary energy
240 EBC = pL * twoOverGamMinOne * 0.5;
241
242 // evaluate the different between the left state normal
243 // velocity and that from the desired condition (note
244 // m_VnInf is using an outwards normal definition.
245 NekDouble VnDiff = uL + m_VnInf[pnt];
246
247 // Boundary velocities & Kinite energy
248 // Note: normals are negated since they point outwards in
249 // the domain
250
251 // Note: Can just use the BC values directly!!
252 for (j = 0; j < nDimensions; ++j)
253 {
254 // Set velocity to the desired conditions modified to
255 // take account of the normal state for Riemann
256 // problem. (Negative accounts for outwards normal definition)
257 // velBC[j] = m_velBC[j][pnt];
258 velBC[j] = m_velBC[j][pnt] -
259 VnDiff * m_traceNormals[j][m_bndToTraceMap[pnt]];
260
261 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
262 }
263
264 //-------------------------------------------------------------------------
265 // Impose Left hand Riemann Invariant boundary conditions
266 m_bndPhys[0][pnt] = rhoL;
267 for (j = 0; j < nDimensions; ++j)
268 {
269 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
270 }
271 m_bndPhys[nDimensions + 1][pnt] = EBC;
272 }
273 else // Outflow
274 {
275
276 // Note: Allowing the switch can cause worse convergence in this
277 // type BC.
278 // So improve it later.
279 if (Mach[pnt] < 1.00)
280 {
281 // subsonic outflow: fix pstar
282 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
283 twoOverGamMinOne;
284 // vn = -m_VnInf[pnt];
285
286 pstar = m_pBC[pnt];
287 rhostar = FwdBnd[0][pnt] * pow((pstar / pressure[pnt]), gamInv);
288 cstar = sqrt(m_gamma * pstar / rhostar);
289 ustar = rR + cstar * twoOverGamMinOne;
290
291 rhoL = rhostar;
292 uL = ustar;
293 pL = pstar;
294 }
295 else
296 {
297 // supersonic outflow
298 // Just set to imposed state and let Riemann BC dictate values
299 rhoL = m_rhoBC[pnt];
300 uL = -m_VnInf[pnt];
301 pL = m_pBC[pnt];
302 }
303
304 // Boundary energy
305 EBC = pL * twoOverGamMinOne * 0.5;
306
307 // Boundary velocities & Kinite energy
308 // Note: normals are negated since they point outwards in
309 // the domain
310 for (j = 0; j < nDimensions; ++j)
311 {
312 velBC[j] = -1.0 * uL * m_traceNormals[j][m_bndToTraceMap[pnt]];
313 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
314 }
315
316 // Impose Left hand Riemann Invariant boundary conditions
317 m_bndPhys[0][pnt] = rhoL;
318 for (j = 0; j < nDimensions; ++j)
319 {
320 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
321 }
322 m_bndPhys[nDimensions + 1][pnt] = EBC;
323 }
324 }
325}
326
327} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:94
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:105
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:100
int m_offset
Offset.
Definition: CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
EnforceEntropyTotalEnthalpy(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)
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.
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
Array< OneD, Array< OneD, NekDouble > > m_velBC
static std::string className
Name of the class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294