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
38
39using namespace std;
40
41namespace Nektar
42{
43
46 "PressureInflowFile", PressureInflowFileBC::create,
47 "Pressure inflow (file) boundary condition.");
48
52 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
53 const int pSpaceDim, const int bcRegion, const int cnt)
54 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, 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
73 Array<OneD, Array<OneD, NekDouble>> &physarray, const NekDouble &time)
74{
75 boost::ignore_unused(time);
76
77 int i, j;
78 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
79 int nVariables = physarray.size();
80 int nDimensions = m_spacedim;
81
82 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
83
84 // Computing the normal velocity for characteristics coming
85 // from inside the computational domain
86 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
87 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
88 for (i = 0; i < nDimensions; ++i)
89 {
90 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
91 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
92 }
93
94 // Computing the absolute value of the velocity in order to compute the
95 // Mach number to decide whether supersonic or subsonic
96 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
97 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
98
99 // Get speed of sound
101 Array<OneD, NekDouble> soundSpeed(nTracePts);
102
103 m_varConv->GetPressure(Fwd, pressure);
104 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
105
106 // Get Mach
107 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
108 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
109 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
110
111 // Auxiliary variables
112 int e, id1, id2, npts, pnt;
113 NekDouble rhoeb;
114
115 // Loop on the m_bcRegions
116 for (e = 0;
117 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
118 {
119 npts = m_fields[0]
120 ->GetBndCondExpansions()[m_bcRegion]
121 ->GetExp(e)
122 ->GetTotPoints();
123 id1 =
124 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
125 id2 =
126 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
127
128 // Get internal energy
129 Array<OneD, NekDouble> tmpPressure(npts, pressure + id2);
130 Array<OneD, NekDouble> rho(npts, Fwd[0] + id2);
131 Array<OneD, NekDouble> Ei(npts);
132 m_varConv->GetEFromRhoP(rho, tmpPressure, Ei);
133
134 // Loop on points of m_bcRegion 'e'
135 for (i = 0; i < npts; i++)
136 {
137 pnt = id2 + i;
138
139 // Subsonic flows
140 if (Mach[pnt] < 0.99)
141 {
142 // Partial extrapolation for subsonic cases
143 for (j = 0; j < nVariables - 1; ++j)
144 {
145 (m_fields[j]
146 ->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 *
155 (m_fieldStorage[j][id1 + i] *
156 m_fieldStorage[j][id1 + i]) /
157 m_fieldStorage[0][id1 + i];
158 }
159
160 rhoeb = Fwd[0][pnt] * Ei[i] + Ek;
161
162 (m_fields[nVariables - 1]
163 ->GetBndCondExpansions()[m_bcRegion]
164 ->UpdatePhys())[id1 + i] = rhoeb;
165 }
166 // Supersonic flows
167 else
168 {
169 for (j = 0; j < nVariables; ++j)
170 {
171 // Extrapolation for supersonic cases
172 (m_fields[j]
173 ->GetBndCondExpansions()[m_bcRegion]
174 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
175 }
176 }
177 }
178 }
179}
180
181} // 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:198
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) 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
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:548
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:569
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:280
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191