Nektar++
PressureOutflowBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PressureOutflowBC.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 outflow boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "PressureOutflowBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
44 "PressureOutflow", PressureOutflowBC::create,
45 "Pressure outflow 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 numBCPts =
55 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
57
58 // Get Pressure
60 numBCPts,
61 m_fields[m_spacedim + 1]->GetBndCondExpansions()[m_bcRegion]->GetPhys(),
62 1, m_pressureStorage, 1);
63}
64
67 [[maybe_unused]] const NekDouble &time)
68{
69 int i, j;
70 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
71 int nVariables = physarray.size();
72 int nDimensions = m_spacedim;
73
74 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
75
76 // Computing the normal velocity for characteristics coming
77 // from inside the computational domain
78 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
79 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
80 for (i = 0; i < nDimensions; ++i)
81 {
82 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
83 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
84 }
85
86 // Computing the absolute value of the velocity in order to compute the
87 // Mach number to decide whether supersonic or subsonic
88 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
89 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
90
91 // Get speed of sound
92 Array<OneD, NekDouble> soundSpeed(nTracePts);
93 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
94
95 // Get Mach
96 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
97 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
98 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
99
100 // Auxiliary variables
101 int e, id1, id2, npts, pnt;
102 NekDouble rhoeb;
103
104 // Loop on the m_bcRegions
105 for (e = 0;
106 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
107 {
108 npts = m_fields[0]
109 ->GetBndCondExpansions()[m_bcRegion]
110 ->GetExp(e)
111 ->GetTotPoints();
112 id1 =
113 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
114 id2 =
115 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
116
117 // Get internal energy
119 Array<OneD, NekDouble> rho(npts, Fwd[0] + id2);
120 Array<OneD, NekDouble> Ei(npts);
121 m_varConv->GetEFromRhoP(rho, pressure, Ei);
122
123 // Loop on points of m_bcRegion 'e'
124 for (i = 0; i < npts; i++)
125 {
126 pnt = id2 + i;
127
128 // Subsonic flows
129 if (Mach[pnt] < 0.99)
130 {
131 // Kinetic energy calculation
132 NekDouble Ek = 0.0;
133 for (j = 1; j < nVariables - 1; ++j)
134 {
135 Ek += 0.5 * (Fwd[j][pnt] * Fwd[j][pnt]) / Fwd[0][pnt];
136 }
137
138 rhoeb = Fwd[0][pnt] * Ei[i] + Ek;
139
140 // Partial extrapolation for subsonic cases
141 for (j = 0; j < nVariables - 1; ++j)
142 {
143 (m_fields[j]
144 ->GetBndCondExpansions()[m_bcRegion]
145 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
146 }
147
148 (m_fields[nVariables - 1]
149 ->GetBndCondExpansions()[m_bcRegion]
150 ->UpdatePhys())[id1 + i] = rhoeb;
151 }
152 // Supersonic flows
153 else
154 {
155 for (j = 0; j < nVariables; ++j)
156 {
157 // Extrapolation for supersonic cases
158 (m_fields[j]
159 ->GetBndCondExpansions()[m_bcRegion]
160 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
161 }
162 }
163 }
164 }
165}
166
167} // 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, NekDouble > m_pressureStorage
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.
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
PressureOutflowBC(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)
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