Nektar++
PressureInflowFileBC.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PressureInflowFileBC.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: Pressure inflow boundary condition
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
37 #include "PressureInflowFileBC.h"
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string PressureInflowFileBC::className = GetCFSBndCondFactory().
45  RegisterCreatorFunction("PressureInflowFile",
46  PressureInflowFileBC::create,
47  "Pressure inflow (file) boundary condition.");
48 
49 PressureInflowFileBC::PressureInflowFileBC(const LibUtilities::SessionReaderSharedPtr& pSession,
51  const Array<OneD, Array<OneD, NekDouble> >& pTraceNormals,
52  const int pSpaceDim,
53  const int bcRegion,
54  const int cnt)
55  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56 {
57  int nvariables = m_fields.size();
58  // Loop over Boundary Regions for PressureInflowFileBC
60 
61  int numBCPts = m_fields[0]->
62  GetBndCondExpansions()[m_bcRegion]->GetNpoints();
63  for (int i = 0; i < nvariables; ++i)
64  {
65  m_fieldStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
67  numBCPts,
68  m_fields[i]->GetBndCondExpansions()[m_bcRegion]->GetPhys(), 1,
69  m_fieldStorage[i], 1);
70  }
71 }
72 
75  Array<OneD, Array<OneD, NekDouble> > &physarray,
76  const NekDouble &time)
77 {
78  boost::ignore_unused(time);
79 
80  int i, j;
81  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
82  int nVariables = physarray.size();
83  int nDimensions = m_spacedim;
84 
85  const Array<OneD, const int> &traceBndMap
86  = m_fields[0]->GetTraceBndMap();
87 
88  // Computing the normal velocity for characteristics coming
89  // from inside the computational domain
90  Array<OneD, NekDouble > Vn (nTracePts, 0.0);
91  Array<OneD, NekDouble > Vel(nTracePts, 0.0);
92  for (i = 0; i < nDimensions; ++i)
93  {
94  Vmath::Vdiv(nTracePts, Fwd[i+1], 1, Fwd[0], 1, Vel, 1);
95  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
96  }
97 
98  // Computing the absolute value of the velocity in order to compute the
99  // Mach number to decide whether supersonic or subsonic
100  Array<OneD, NekDouble > absVel(nTracePts, 0.0);
101  m_varConv->GetAbsoluteVelocity(Fwd, absVel);
102 
103  // Get speed of sound
104  Array<OneD, NekDouble > pressure (nTracePts);
105  Array<OneD, NekDouble > soundSpeed(nTracePts);
106 
107  m_varConv->GetPressure(Fwd, pressure);
108  m_varConv->GetSoundSpeed(Fwd, soundSpeed);
109 
110  // Get Mach
111  Array<OneD, NekDouble > Mach(nTracePts, 0.0);
112  Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
113  Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
114 
115  // Auxiliary variables
116  int e, id1, id2, npts, pnt;
117  NekDouble rhoeb;
118 
119  // Loop on the m_bcRegions
120  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
121  GetExpSize(); ++e)
122  {
123  npts = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
124  GetExp(e)->GetTotPoints();
125  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
126  GetPhys_Offset(e);
127  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset+e]);
128 
129  // Get internal energy
130  Array<OneD, NekDouble> tmpPressure (npts, pressure+id2);
131  Array<OneD, NekDouble> rho (npts, Fwd[0]+id2);
132  Array<OneD, NekDouble> e(npts);
133  m_varConv->GetEFromRhoP(rho, tmpPressure, e);
134 
135  // Loop on points of m_bcRegion 'e'
136  for (i = 0; i < npts; i++)
137  {
138  pnt = id2+i;
139 
140  // Subsonic flows
141  if (Mach[pnt] < 0.99)
142  {
143  // Partial extrapolation for subsonic cases
144  for (j = 0; j < nVariables-1; ++j)
145  {
146  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
147  UpdatePhys())[id1+i] = m_fieldStorage[j][id1+i];
148  }
149 
150  // Kinetic energy calculation
151  NekDouble Ek = 0.0;
152  for (j = 1; j < nVariables-1; ++j)
153  {
154  Ek += 0.5 * (m_fieldStorage[j][id1+i] *
155  m_fieldStorage[j][id1+i]) /
156  m_fieldStorage[0][id1+i];
157  }
158 
159  rhoeb = Fwd[0][pnt] * e[i] + Ek;
160 
161  (m_fields[nVariables-1]->GetBndCondExpansions()[m_bcRegion]->
162  UpdatePhys())[id1+i] = rhoeb;
163  }
164  // Supersonic flows
165  else
166  {
167  for (j = 0; j < nVariables; ++j)
168  {
169  // Extrapolation for supersonic cases
170  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
171  UpdatePhys())[id1+i] = Fwd[j][pnt];
172  }
173  }
174  }
175  }
176 }
177 
178 }
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
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
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)
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 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
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199