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 = GetCFSBndCondFactory().
45  RegisterCreatorFunction("RiemannInvariant",
46  RiemannInvariantBC::create,
47  "Riemann invariant boundary condition.");
48 
49 RiemannInvariantBC::RiemannInvariantBC(
52  const Array<OneD, Array<OneD, NekDouble> >& pTraceNormals,
53  const int pSpaceDim,
54  const int bcRegion,
55  const int cnt)
56  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
57 {
58  // Calculate VnInf
59  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
60  m_VnInf = Array<OneD, NekDouble> (nTracePts, 0.0);
61 
62  // Computing the normal velocity for characteristics coming
63  // from outside the computational domain
64  for( int i =0; i < m_spacedim; i++)
65  {
66  Vmath::Svtvp(nTracePts, m_velInf[i],
67  m_traceNormals[i], 1,
68  m_VnInf, 1,
69  m_VnInf, 1);
70  }
71 }
72 
75  Array<OneD, Array<OneD, NekDouble> > &physarray,
76  const NekDouble &time)
77 {
78  boost::ignore_unused(physarray, time);
79 
80  int i, j;
81  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
82  int nDimensions = m_spacedim;
83 
84  const Array<OneD, const int> &traceBndMap
85  = m_fields[0]->GetTraceBndMap();
86 
87  NekDouble gammaInv = 1.0 / m_gamma;
88  NekDouble gammaMinusOne = m_gamma - 1.0;
89  NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
90 
91  // Computing the normal velocity for characteristics coming
92  // from inside the computational domain
93  Array<OneD, NekDouble > Vn (nTracePts, 0.0);
94  Array<OneD, NekDouble > Vel(nTracePts, 0.0);
95  for (i = 0; i < nDimensions; ++i)
96  {
97  Vmath::Vdiv(nTracePts, Fwd[i+1], 1, Fwd[0], 1, Vel, 1);
98  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
99  }
100 
101  // Computing the absolute value of the velocity in order to compute the
102  // Mach number to decide whether supersonic or subsonic
103  Array<OneD, NekDouble > absVel(nTracePts, 0.0);
104  m_varConv->GetAbsoluteVelocity(Fwd, absVel);
105 
106  // Get speed of sound
107  Array<OneD, NekDouble > pressure (nTracePts);
108  Array<OneD, NekDouble > soundSpeed(nTracePts);
109 
110  m_varConv->GetPressure(Fwd, pressure);
111  m_varConv->GetSoundSpeed(Fwd, soundSpeed);
112 
113  // Get Mach
114  Array<OneD, NekDouble > Mach(nTracePts, 0.0);
115  Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
116  Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
117 
118  // Auxiliary variables
119  int eMax;
120  int e, id1, id2, nBCEdgePts, pnt;
121  NekDouble cPlus, rPlus, cMinus, rMinus, VDBC, VNBC;
122  Array<OneD, NekDouble> velBC(nDimensions, 0.0);
123  Array<OneD, NekDouble> rhoVelBC(nDimensions, 0.0);
124  NekDouble rhoBC, EBC, cBC, sBC, pBC;
125 
126  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
127 
128  // Loop on m_bcRegions
129  for (e = 0; e < eMax; ++e)
130  {
131  nBCEdgePts = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
132  GetExp(e)->GetTotPoints();
133 
134  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
135  GetPhys_Offset(e);
136  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset+e]);
137 
138  // Loop on the points of the m_bcRegion
139  for (i = 0; i < nBCEdgePts; i++)
140  {
141  pnt = id2+i;
142 
143  // Impose inflow Riemann invariant
144  if (Vn[pnt] <= 0.0)
145  {
146  // Subsonic flows
147  if (Mach[pnt] < 1.00)
148  {
149  // + Characteristic from inside
150  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
151  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
152 
153  // - Characteristic from boundary
154  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
155  rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
156  }
157  else
158  {
159  // + Characteristic from inside
160  cPlus = sqrt(m_gamma * m_pInf / m_rhoInf);
161  rPlus = m_VnInf[pnt] + 2.0 * cPlus * gammaMinusOneInv;
162 
163  // + Characteristic from inside
164  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
165  rMinus = m_VnInf[pnt] - 2.0 * cPlus * gammaMinusOneInv;
166  }
167 
168  // Riemann boundary variables
169  VNBC = 0.5 * (rPlus + rMinus);
170  cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
171  VDBC = VNBC - m_VnInf[pnt];
172 
173  // Thermodynamic boundary variables
174  sBC = m_pInf / (pow(m_rhoInf, m_gamma));
175  rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
176  pBC = rhoBC * cBC * cBC * gammaInv;
177 
178  // Kinetic energy initialiasation
179  NekDouble EkBC = 0.0;
180 
181  // Boundary velocities
182  for ( j = 0; j < nDimensions; ++j)
183  {
184  velBC[j] = m_velInf[j] + VDBC * m_traceNormals[j][pnt];
185  rhoVelBC[j] = rhoBC * velBC[j];
186  EkBC += 0.5 * rhoBC * velBC[j]*velBC[j];
187  }
188 
189  // Boundary energy
190  EBC = pBC * gammaMinusOneInv + EkBC;
191 
192  // Imposing Riemann Invariant boundary conditions
193  (m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
194  UpdatePhys())[id1+i] = rhoBC;
195  for (j = 0; j < nDimensions; ++j)
196  {
197  (m_fields[j+1]->GetBndCondExpansions()[m_bcRegion]->
198  UpdatePhys())[id1+i] = rhoVelBC[j];
199  }
200  (m_fields[nDimensions+1]->GetBndCondExpansions()[m_bcRegion]->
201  UpdatePhys())[id1+i] = EBC;
202 
203  }
204  else // Impose outflow Riemann invariant
205  {
206  // Subsonic flows
207  if (Mach[pnt] < 1.00)
208  {
209  // + Characteristic from inside
210  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
211  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
212 
213  // - Characteristic from boundary
214  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
215  rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
216  }
217  else
218  {
219  // + Characteristic from inside
220  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
221  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
222 
223  // + Characteristic from inside
224  cMinus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
225  rMinus = Vn[pnt] - 2.0 * cPlus * gammaMinusOneInv;
226  }
227 
228  // Riemann boundary variables
229  VNBC = 0.5 * (rPlus + rMinus);
230  cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
231  VDBC = VNBC - Vn[pnt];
232 
233  // Thermodynamic boundary variables
234  sBC = pressure[pnt] / (pow(Fwd[0][pnt], m_gamma));
235  rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
236  pBC = rhoBC * cBC * cBC * gammaInv;
237 
238  // Kinetic energy initialiasation
239  NekDouble EkBC = 0.0;
240 
241  // Boundary velocities
242  for ( j = 0; j < nDimensions; ++j)
243  {
244  velBC[j] = Fwd[j+1][pnt] / Fwd[0][pnt] +
245  VDBC * m_traceNormals[j][pnt];
246  rhoVelBC[j] = rhoBC * velBC[j];
247  EkBC += 0.5 * rhoBC * velBC[j]*velBC[j];
248  }
249 
250  // Boundary energy
251  EBC = pBC * gammaMinusOneInv + EkBC;
252 
253  // Imposing Riemann Invariant boundary conditions
254  (m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
255  UpdatePhys())[id1+i] = rhoBC;
256  for (j = 0; j < nDimensions; ++j)
257  {
258  (m_fields[j+1]->GetBndCondExpansions()[m_bcRegion]->
259  UpdatePhys())[id1+i] = rhoVelBC[j];
260  }
261  (m_fields[nDimensions+1]->GetBndCondExpansions()[m_bcRegion]->
262  UpdatePhys())[id1+i] = EBC;
263  }
264  }
265  }
266 }
267 
268 }
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
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
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)
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:1
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:565
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:493
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:513
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:257
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267