Nektar++
EnforceEntropyVelocity.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: EnforceEntropyVelocity.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 velocity 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 "EnforceEntropyVelocity", EnforceEntropyVelocity::create,
48 "Riemann boundary condition enforcing entropy and velocity.");
49
53 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
54 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
55 const int pSpaceDim, const int bcRegion, const int cnt)
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
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
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
102 Vmath::Vcopy(m_npts, BCvals[0], 1, m_rhoBC, 1);
103
105 // Evaluate velocity on boundary
106 for (int i = 0; i < m_spacedim; ++i)
107 {
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 }
113 m_varConv->GetPressure(BCvals, m_pBC);
114
115 // Computing the normal velocity for characteristics coming
116 // from outside the computational domain
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}
126
129 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &physarray,
130 [[maybe_unused]] const NekDouble &time)
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
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
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
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 (ghost state)
178 // R represents properties inside boundary (numerical state)
179 NekDouble rhoL, uL, pL, cL, sL;
180 NekDouble cR, uR; // rR, r = RiemannInvariant
181 NekDouble EBC; // cstar, pstar, rhostar, ustar, BC = star
182
183 NekDouble gamMinOne = m_gamma - 1.0;
184 NekDouble oneOverGamMinOne = 1.0 / gamMinOne;
185 NekDouble twoOverGamMinOne = 2.0 * oneOverGamMinOne;
186 NekDouble gamInv = 1.0 / m_gamma;
187
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 state variables
198 cR = sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]);
199 uR = -Vn[pnt];
200
201 uL = 2 * (-m_VnInf[pnt]) - uR;
202 cL = cR;
203 sL = m_pBC[pnt] / (pow(m_rhoBC[pnt], m_gamma));
204 rhoL = pow((cL * cL) / (m_gamma * sL), oneOverGamMinOne);
205 pL = cL * cL * rhoL * gamInv;
206 }
207 else // Supersonic inflow
208 {
209 // all characteristics are from left so just impose
210 // star state to left values
211 // Note: m_vnInf is the negative of the normal velocity
212 // across boundary
213 rhoL = m_rhoBC[pnt];
214 uL = -m_VnInf[pnt];
215 pL = m_pBC[pnt];
216 }
217
218 // Boundary energy
219 EBC = pL * twoOverGamMinOne * 0.5;
220
221 // evaluate the different between the left state normal
222 // velocity and that from the desired condition (note
223 // m_VnInf is using an outwards normal definition.
224 NekDouble VnDiff = uL + m_VnInf[pnt];
225
226 // Note: Can just use the BC values directly!!
227 for (j = 0; j < nDimensions; ++j)
228 {
229 // Set velocity to the desired conditions modified to
230 // take account of the normal state for Riemann
231 // problem. (Negative accounts for outwards normal definition)
232 // velBC[j] = m_velBC[j][pnt];
233 velBC[j] = m_velBC[j][pnt] -
234 VnDiff * m_traceNormals[j][m_bndToTraceMap[pnt]];
235
236 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
237 }
238
239 //-------------------------------------------------------------------------
240 // Impose Left hand Riemann Invariant boundary conditions
241 m_bndPhys[0][pnt] = rhoL;
242 for (j = 0; j < nDimensions; ++j)
243 {
244 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
245 }
246 m_bndPhys[nDimensions + 1][pnt] = EBC;
247 }
248 else // Outflow
249 {
250
251 // Note: Allowing the switch can cause worse convergence in this
252 // type BC.
253 // So improve it later.
254 if (Mach[pnt] < 1.00)
255 {
256 rhoL = m_rhoBC[pnt];
257 uL = -m_VnInf[pnt];
258 pL = m_pBC[pnt];
259 }
260 else
261 {
262 // supersonic outflow
263 // Just set to imposed state and let Riemann BC dictate values
264 rhoL = m_rhoBC[pnt];
265 uL = -m_VnInf[pnt];
266 pL = m_pBC[pnt];
267 }
268
269 // Boundary energy
270 EBC = pL * twoOverGamMinOne * 0.5;
271
272 // Boundary velocities & Kinite energy
273 // Note: normals are negated since they point outwards in
274 // the domain
275 for (j = 0; j < nDimensions; ++j)
276 {
277 velBC[j] = -1.0 * uL * m_traceNormals[j][m_bndToTraceMap[pnt]];
278 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
279 }
280
281 // Impose Left hand Riemann Invariant boundary conditions
282 m_bndPhys[0][pnt] = rhoL;
283 for (j = 0; j < nDimensions; ++j)
284 {
285 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
286 }
287 m_bndPhys[nDimensions + 1][pnt] = EBC;
288 }
289 }
290}
291
292} // 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
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, NekDouble > m_rhoBC
EnforceEntropyVelocity(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)
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
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_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