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