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