Nektar++
RiemannInvariantBC.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: RiemannInvariantBC.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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Riemann invariant boundary condition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
37 #include "RiemannInvariantBC.h"
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string RiemannInvariantBC::className =
46  "RiemannInvariant", RiemannInvariantBC::create,
47  "Riemann invariant boundary condition.");
48 
49 RiemannInvariantBC::RiemannInvariantBC(
52  const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
53  const int pSpaceDim, const int bcRegion, const int cnt)
54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56  // Calculate VnInf
57  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
58  m_VnInf = Array<OneD, NekDouble>(nTracePts, 0.0);
59 
60  // Computing the normal velocity for characteristics coming
61  // from outside the computational domain
62  for (int i = 0; i < m_spacedim; i++)
63  {
64  Vmath::Svtvp(nTracePts, m_velInf[i], m_traceNormals[i], 1, m_VnInf, 1,
65  m_VnInf, 1);
66  }
67 }
68 
70  Array<OneD, Array<OneD, NekDouble>> &physarray,
71  const NekDouble &time)
72 {
73  boost::ignore_unused(physarray, time);
74 
75  int i, j;
76  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
77  int nDimensions = m_spacedim;
78 
79  const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
80 
81  NekDouble gammaInv = 1.0 / m_gamma;
82  NekDouble gammaMinusOne = m_gamma - 1.0;
83  NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
84 
85  // Computing the normal velocity for characteristics coming
86  // from inside the computational domain
87  Array<OneD, NekDouble> Vn(nTracePts, 0.0);
88  Array<OneD, NekDouble> Vel(nTracePts, 0.0);
89  for (i = 0; i < nDimensions; ++i)
90  {
91  Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
92  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
93  }
94 
95  // Computing the absolute value of the velocity in order to compute the
96  // Mach number to decide whether supersonic or subsonic
97  Array<OneD, NekDouble> absVel(nTracePts, 0.0);
98  m_varConv->GetAbsoluteVelocity(Fwd, absVel);
99 
100  // Get speed of sound
101  Array<OneD, NekDouble> pressure(nTracePts);
102  Array<OneD, NekDouble> soundSpeed(nTracePts);
103 
104  m_varConv->GetPressure(Fwd, pressure);
105  m_varConv->GetSoundSpeed(Fwd, soundSpeed);
106 
107  // Get Mach
108  Array<OneD, NekDouble> Mach(nTracePts, 0.0);
109  Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
110  Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
111 
112  // Auxiliary variables
113  int eMax;
114  int e, id1, id2, nBCEdgePts, pnt;
115  NekDouble cPlus, rPlus, cMinus, rMinus, VDBC, VNBC;
116  Array<OneD, NekDouble> velBC(nDimensions, 0.0);
117  Array<OneD, NekDouble> rhoVelBC(nDimensions, 0.0);
118  NekDouble rhoBC, EBC, cBC, sBC, pBC;
119 
120  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
121 
122  // Loop on m_bcRegions
123  for (e = 0; e < eMax; ++e)
124  {
125  nBCEdgePts = m_fields[0]
126  ->GetBndCondExpansions()[m_bcRegion]
127  ->GetExp(e)
128  ->GetTotPoints();
129 
130  id1 =
131  m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
132  id2 =
133  m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
134 
135  // Loop on the points of the m_bcRegion
136  for (i = 0; i < nBCEdgePts; i++)
137  {
138  pnt = id2 + i;
139 
140  // Impose inflow Riemann invariant
141  if (Vn[pnt] <= 0.0)
142  {
143  // Subsonic flows
144  if (Mach[pnt] < 1.00)
145  {
146  // + Characteristic from inside
147  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
148  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
149 
150  // - Characteristic from boundary
151  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
152  rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
153  }
154  else
155  {
156  // + Characteristic from inside
157  cPlus = sqrt(m_gamma * m_pInf / m_rhoInf);
158  rPlus = m_VnInf[pnt] + 2.0 * cPlus * gammaMinusOneInv;
159 
160  // + Characteristic from inside
161  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
162  rMinus = m_VnInf[pnt] - 2.0 * cPlus * gammaMinusOneInv;
163  }
164 
165  // Riemann boundary variables
166  VNBC = 0.5 * (rPlus + rMinus);
167  cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
168  VDBC = VNBC - m_VnInf[pnt];
169 
170  // Thermodynamic boundary variables
171  sBC = m_pInf / (pow(m_rhoInf, m_gamma));
172  rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
173  pBC = rhoBC * cBC * cBC * gammaInv;
174 
175  // Kinetic energy initialiasation
176  NekDouble EkBC = 0.0;
177 
178  // Boundary velocities
179  for (j = 0; j < nDimensions; ++j)
180  {
181  velBC[j] = m_velInf[j] + VDBC * m_traceNormals[j][pnt];
182  rhoVelBC[j] = rhoBC * velBC[j];
183  EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
184  }
185 
186  // Boundary energy
187  EBC = pBC * gammaMinusOneInv + EkBC;
188 
189  // Imposing Riemann Invariant boundary conditions
190  (m_fields[0]
191  ->GetBndCondExpansions()[m_bcRegion]
192  ->UpdatePhys())[id1 + i] = rhoBC;
193  for (j = 0; j < nDimensions; ++j)
194  {
195  (m_fields[j + 1]
196  ->GetBndCondExpansions()[m_bcRegion]
197  ->UpdatePhys())[id1 + i] = rhoVelBC[j];
198  }
199  (m_fields[nDimensions + 1]
200  ->GetBndCondExpansions()[m_bcRegion]
201  ->UpdatePhys())[id1 + i] = EBC;
202  }
203  else // Impose outflow Riemann invariant
204  {
205  // Subsonic flows
206  if (Mach[pnt] < 1.00)
207  {
208  // + Characteristic from inside
209  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
210  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
211 
212  // - Characteristic from boundary
213  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
214  rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
215  }
216  else
217  {
218  // + Characteristic from inside
219  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
220  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
221 
222  // + Characteristic from inside
223  cMinus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
224  rMinus = Vn[pnt] - 2.0 * cPlus * gammaMinusOneInv;
225  }
226 
227  // Riemann boundary variables
228  VNBC = 0.5 * (rPlus + rMinus);
229  cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
230  VDBC = VNBC - Vn[pnt];
231 
232  // Thermodynamic boundary variables
233  sBC = pressure[pnt] / (pow(Fwd[0][pnt], m_gamma));
234  rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
235  pBC = rhoBC * cBC * cBC * gammaInv;
236 
237  // Kinetic energy initialiasation
238  NekDouble EkBC = 0.0;
239 
240  // Boundary velocities
241  for (j = 0; j < nDimensions; ++j)
242  {
243  velBC[j] = Fwd[j + 1][pnt] / Fwd[0][pnt] +
244  VDBC * m_traceNormals[j][pnt];
245  rhoVelBC[j] = rhoBC * velBC[j];
246  EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
247  }
248 
249  // Boundary energy
250  EBC = pBC * gammaMinusOneInv + EkBC;
251 
252  // Imposing Riemann Invariant boundary conditions
253  (m_fields[0]
254  ->GetBndCondExpansions()[m_bcRegion]
255  ->UpdatePhys())[id1 + i] = rhoBC;
256  for (j = 0; j < nDimensions; ++j)
257  {
258  (m_fields[j + 1]
259  ->GetBndCondExpansions()[m_bcRegion]
260  ->UpdatePhys())[id1 + i] = rhoVelBC[j];
261  }
262  (m_fields[nDimensions + 1]
263  ->GetBndCondExpansions()[m_bcRegion]
264  ->UpdatePhys())[id1 + i] = EBC;
265  }
266  }
267  }
268 }
269 
270 } // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
NekDouble m_rhoInf
Definition: CFSBndCond.h:103
NekDouble m_pInf
Definition: CFSBndCond.h:104
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:106
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray, const NekDouble &time) override
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:553
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
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.cpp:284
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294