Nektar++
EnforceEntropyPressure.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: EnforceEntropyPressure.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 pressure at the inflow boundary;
34// Enforce the pressure 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 "EnforceEntropyPressure", EnforceEntropyPressure::create,
48 "Riemann boundary condition enforcing entropy and pressure.");
49
53 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
54 const int pSpaceDim, const int bcRegion, const int cnt)
55 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56{
57
59 m_fields[0]->GetBndCondExpansions()[m_bcRegion];
60
61 //-> Gather a list of index from trace to this boundary
62 m_npts = bndexp->GetTotPoints();
63
65
66 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
67
68 // Construct a map for the boundary to trace map for easy acess to
69 // phys space points
70 int cnt1 = 0;
71 for (int e = 0; e < bndexp->GetNumElmts(); ++e)
72 {
73 int nTracePts = bndexp->GetExp(e)->GetTotPoints();
74
75 int id =
76 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
77
78 // Loop on the points of the m_bcRegion
79 for (int i = 0; i < nTracePts; i++)
80 {
81 // the ith point in region e
82 m_bndToTraceMap[cnt1++] = id + i;
83 }
84 }
85
88
89 for (int i = 0; i < m_fields.size(); ++i)
90 {
91 m_bndPhys[i] =
92 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys();
93
94 BCvals[i] = Array<OneD, NekDouble>(m_npts);
95 Vmath::Vcopy(m_npts, m_bndPhys[i], 1, BCvals[i], 1);
96 }
97
98 // Set up boudnary required BCs
100 Vmath::Vcopy(m_npts, BCvals[0], 1, m_rhoBC, 1);
101
103 // Evaluate velocity on boundary
104 for (int i = 0; i < m_spacedim; ++i)
105 {
107 Vmath::Vcopy(m_npts, BCvals[i + 1], 1, m_velBC[i], 1);
108 Vmath::Vdiv(m_npts, m_velBC[i], 1, m_rhoBC, 1, m_velBC[i], 1);
109 }
111 m_varConv->GetPressure(BCvals, m_pBC);
112
113 // Computing the normal velocity for characteristics coming
114 // from outside the computational domain
116 for (int i = 0; i < m_spacedim; i++)
117 {
118 for (int j = 0; j < m_npts; ++j)
119 {
120 m_VnInf[j] += m_traceNormals[i][m_bndToTraceMap[j]] * m_velBC[i][j];
121 }
122 }
123}
124
127 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &physarray,
128 [[maybe_unused]] const NekDouble &time)
129{
130 int i, j;
131 int nDimensions = m_spacedim;
132
133 Array<OneD, Array<OneD, NekDouble>> FwdBnd(Fwd.size());
134 Array<OneD, Array<OneD, NekDouble>> bndPhys(Fwd.size());
135
136 // make a local copy of Fwd along boundary of interest
137 for (i = 0; i < Fwd.size(); ++i)
138 {
139 FwdBnd[i] = Array<OneD, NekDouble>(m_npts);
140 for (j = 0; j < m_npts; ++j)
141 {
142 FwdBnd[i][j] = Fwd[i][m_bndToTraceMap[j]];
143 }
144 }
145
146 // Computing the normal velocity for characteristics coming
147 // from inside the computational domain
149
150 for (i = 0; i < nDimensions; ++i)
151 {
152 for (j = 0; j < m_npts; ++j)
153 {
154 Vn[j] += m_traceNormals[i][m_bndToTraceMap[j]] * FwdBnd[i + 1][j];
155 }
156 }
157 // divide by density.
158 Vmath::Vdiv(m_npts, Vn, 1, FwdBnd[0], 1, Vn, 1);
159
160 // Get speed of sound
162 Array<OneD, NekDouble> soundSpeed(m_npts);
163
164 m_varConv->GetPressure(FwdBnd, pressure);
165 m_varConv->GetSoundSpeed(FwdBnd, soundSpeed);
166
167 // Get Mach. Note: it is computed by Vn/c
169 Vmath::Vdiv(m_npts, Vn, 1, soundSpeed, 1, Mach, 1);
170 Vmath::Vabs(m_npts, Mach, 1, Mach, 1);
171
172 // Auxiliary variables
173 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
174
175 // L represents properties outside boundary
176 // R represents properties inside boundary (numerical state)
177 NekDouble rhoL, uL, pL;
178 NekDouble EBC, rR, cstar, pstar, rhostar, ustar; // vn
179
180 NekDouble gamMinOne = m_gamma - 1.0;
181 NekDouble twoOverGamMinOne = 2.0 / gamMinOne;
182 NekDouble gamInv = 1.0 / m_gamma;
183
184 // Loop on m_bcRegions
185 for (int pnt = 0; pnt < m_npts; ++pnt)
186 {
187 // Impose inflow Riemann invariant
188 if (Vn[pnt] <= 0.0)
189 {
190 // Subsonic flows
191 if (Mach[pnt] < 1.00)
192 {
193 // right characteristic
194 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
195 twoOverGamMinOne;
196 // vn = -m_VnInf[pnt]; //vn BC
197
198 // fix rhostar and pstar to be the input values
199 // compute ustar using left-pointing characteristic line IR^-
200 pstar = m_pBC[pnt];
201 rhostar = m_rhoBC[pnt];
202 cstar = sqrt(m_gamma * pstar / rhostar);
203 ustar = rR + cstar * twoOverGamMinOne;
204
205 // add supplement equation that rhoL=rhostar
206 // then pL=pstar, according to IL^0
207 // and uL=ustar, according to IL^+
208 rhoL = rhostar;
209 pL = pstar;
210 uL = ustar;
211
212 // std subsnoic inflow
213 // rhoL = m_rhoBC[pnt];
214 // uL = -m_VnInf[pnt];
215 // pL = m_pBC[pnt];
216 }
217 else // Supersonic inflow
218 {
219 // all characteristics are from left so just impose
220 // star state to left values
221 // Note: m_vnInf is the negative of the normal velocity
222 // across boundary
223 rhoL = m_rhoBC[pnt];
224 uL = -m_VnInf[pnt];
225 pL = m_pBC[pnt];
226 }
227
228 // Boundary energy
229 EBC = pL * twoOverGamMinOne * 0.5;
230
231 // evaluate the different between the left state normal
232 // velocity and that from the desired condition (note
233 // m_VnInf is using an outwards normal definition.
234 NekDouble VnDiff = uL + m_VnInf[pnt];
235
236 // Boundary velocities & Kinite energy
237 // Note: normals are negated since they point outwards in
238 // the domain
239
240 // Note: Can just use the BC values directly!!
241 for (j = 0; j < nDimensions; ++j)
242 {
243 // Set velocity to the desired conditions modified to
244 // take account of the normal state for Riemann
245 // problem. (Negative accounts for outwards normal definition)
246 // velBC[j] = m_velBC[j][pnt];
247 velBC[j] = m_velBC[j][pnt] -
248 VnDiff * m_traceNormals[j][m_bndToTraceMap[pnt]];
249
250 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
251 }
252
253 //-------------------------------------------------------------------------
254 // Impose Left hand Riemann Invariant boundary conditions
255 m_bndPhys[0][pnt] = rhoL;
256 for (j = 0; j < nDimensions; ++j)
257 {
258 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
259 }
260 m_bndPhys[nDimensions + 1][pnt] = EBC;
261 }
262 else // Outflow
263 {
264
265 // Note: Allowing the switch can cause worse convergence in this
266 // type BC.
267 // So improve it later.
268 if (Mach[pnt] < 1.00)
269 {
270 // subsonic outflow: fix pstar
271 rR = -Vn[pnt] - sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]) *
272 twoOverGamMinOne;
273 // vn = -m_VnInf[pnt];
274
275 pstar = m_pBC[pnt];
276 rhostar = FwdBnd[0][pnt] * pow((pstar / pressure[pnt]), gamInv);
277 cstar = sqrt(m_gamma * pstar / rhostar);
278 ustar = rR + cstar * twoOverGamMinOne;
279
280 rhoL = rhostar;
281 uL = ustar;
282 pL = pstar;
283 }
284 else
285 {
286 // supersonic outflow
287 // Just set to imposed state and let Riemann BC dictate values
288 rhoL = m_rhoBC[pnt];
289 uL = -m_VnInf[pnt];
290 pL = m_pBC[pnt];
291 }
292
293 // Boundary energy
294 EBC = pL * twoOverGamMinOne * 0.5;
295
296 // Boundary velocities & Kinite energy
297 // Note: normals are negated since they point outwards in
298 // the domain
299 for (j = 0; j < nDimensions; ++j)
300 {
301 velBC[j] = -1.0 * uL * m_traceNormals[j][m_bndToTraceMap[pnt]];
302 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
303 }
304
305 // Impose Left hand Riemann Invariant boundary conditions
306 m_bndPhys[0][pnt] = rhoL;
307 for (j = 0; j < nDimensions; ++j)
308 {
309 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
310 }
311 m_bndPhys[nDimensions + 1][pnt] = EBC;
312 }
313 }
314}
315
316} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
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
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.
static std::string className
Name of the class.
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
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
EnforceEntropyPressure(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)
Array< OneD, NekDouble > m_rhoBC
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294