Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // 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: Pressure inflow boundary condition
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include "PressureInflowFileBC.h"
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42 
43 std::string PressureInflowFileBC::className = GetCFSBndCondFactory().
44  RegisterCreatorFunction("PressureInflowFile",
45  PressureInflowFileBC::create,
46  "Pressure inflow (file) boundary condition.");
47 
48 PressureInflowFileBC::PressureInflowFileBC(const LibUtilities::SessionReaderSharedPtr& pSession,
50  const Array<OneD, Array<OneD, NekDouble> >& pTraceNormals,
51  const int pSpaceDim,
52  const int bcRegion,
53  const int cnt)
54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56  int nvariables = m_fields.num_elements();
57  // Loop over Boundary Regions for PressureInflowFileBC
59 
60  int numBCPts = m_fields[0]->
61  GetBndCondExpansions()[m_bcRegion]->GetNpoints();
62  for (int i = 0; i < nvariables; ++i)
63  {
64  m_fieldStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
66  numBCPts,
67  m_fields[i]->GetBndCondExpansions()[m_bcRegion]->GetPhys(), 1,
68  m_fieldStorage[i], 1);
69  }
70 }
71 
74  Array<OneD, Array<OneD, NekDouble> > &physarray,
75  const NekDouble &time)
76 {
77  int i, j;
78  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
79  int nVariables = physarray.num_elements();
80  int nDimensions = m_spacedim;
81 
82  const Array<OneD, const int> &traceBndMap
83  = m_fields[0]->GetTraceBndMap();
84 
85  NekDouble gammaMinusOne = m_gamma - 1.0;
86  NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
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, pressure, 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  // Loop on points of m_bcRegion 'e'
130  for (i = 0; i < npts; i++)
131  {
132  pnt = id2+i;
133 
134  // Subsonic flows
135  if (Mach[pnt] < 0.99)
136  {
137  // Partial extrapolation for subsonic cases
138  for (j = 0; j < nVariables-1; ++j)
139  {
140  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
141  UpdatePhys())[id1+i] = m_fieldStorage[j][id1+i];
142  }
143 
144  // Kinetic energy calculation
145  NekDouble Ek = 0.0;
146  for (j = 1; j < nVariables-1; ++j)
147  {
148  Ek += 0.5 * (m_fieldStorage[j][id1+i] *
149  m_fieldStorage[j][id1+i]) /
150  m_fieldStorage[0][id1+i];
151  }
152  rhoeb = gammaMinusOneInv * pressure[pnt] + Ek;
153 
154  (m_fields[nVariables-1]->GetBndCondExpansions()[m_bcRegion]->
155  UpdatePhys())[id1+i] = rhoeb;
156  }
157  // Supersonic flows
158  else
159  {
160  for (j = 0; j < nVariables; ++j)
161  {
162  // Extrapolation for supersonic cases
163  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
164  UpdatePhys())[id1+i] = Fwd[j][pnt];
165  }
166  }
167  }
168  }
169 }
170 
171 }
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:90
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)
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:442
STL namespace.
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:241
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:86
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:424
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:42
static std::string npts
Definition: InputFld.cpp:43
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
double NekDouble
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:95
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:92
int m_offset
Offset.
Definition: CFSBndCond.h:103
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:88
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:101
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage