Nektar++
RiemannInvariantBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RiemannInvariantBC.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: Riemann invariant boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "RiemannInvariantBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
44 "RiemannInvariant", RiemannInvariantBC::create,
45 "Riemann invariant 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 // Calculate VnInf
55 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
56 m_VnInf = Array<OneD, NekDouble>(nTracePts, 0.0);
57
58 // Computing the normal velocity for characteristics coming
59 // from outside the computational domain
60 for (int i = 0; i < m_spacedim; i++)
61 {
62 Vmath::Svtvp(nTracePts, m_velInf[i], m_traceNormals[i], 1, m_VnInf, 1,
63 m_VnInf, 1);
64 }
65}
66
69 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &physarray,
70 [[maybe_unused]] const NekDouble &time)
71{
72 int i, j;
73 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
74 int nDimensions = m_spacedim;
75
76 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
77
78 NekDouble gammaInv = 1.0 / m_gamma;
79 NekDouble gammaMinusOne = m_gamma - 1.0;
80 NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
81
82 // Computing the normal velocity for characteristics coming
83 // from inside the computational domain
84 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
85 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
86 for (i = 0; i < nDimensions; ++i)
87 {
88 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
89 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
90 }
91
92 // Computing the absolute value of the velocity in order to compute the
93 // Mach number to decide whether supersonic or subsonic
94 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
95 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
96
97 // Get speed of sound
99 Array<OneD, NekDouble> soundSpeed(nTracePts);
100
101 m_varConv->GetPressure(Fwd, pressure);
102 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
103
104 // Get Mach
105 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
106 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
107 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
108
109 // Auxiliary variables
110 int eMax;
111 int e, id1, id2, nBCEdgePts, pnt;
112 NekDouble cPlus, rPlus, cMinus, rMinus, VDBC, VNBC;
113 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
114 Array<OneD, NekDouble> rhoVelBC(nDimensions, 0.0);
115 NekDouble rhoBC, EBC, cBC, sBC, pBC;
116
117 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
118
119 // Loop on m_bcRegions
120 for (e = 0; e < eMax; ++e)
121 {
122 nBCEdgePts = m_fields[0]
123 ->GetBndCondExpansions()[m_bcRegion]
124 ->GetExp(e)
125 ->GetTotPoints();
126
127 id1 =
128 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
129 id2 =
130 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
131
132 // Loop on the points of the m_bcRegion
133 for (i = 0; i < nBCEdgePts; i++)
134 {
135 pnt = id2 + i;
136
137 // Impose inflow Riemann invariant
138 if (Vn[pnt] <= 0.0)
139 {
140 // Subsonic flows
141 if (Mach[pnt] < 1.00)
142 {
143 // + Characteristic from inside
144 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
145 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
146
147 // - Characteristic from boundary
148 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
149 rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
150 }
151 else
152 {
153 // + Characteristic from inside
154 cPlus = sqrt(m_gamma * m_pInf / m_rhoInf);
155 rPlus = m_VnInf[pnt] + 2.0 * cPlus * gammaMinusOneInv;
156
157 // + Characteristic from inside
158 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
159 rMinus = m_VnInf[pnt] - 2.0 * cPlus * gammaMinusOneInv;
160 }
161
162 // Riemann boundary variables
163 VNBC = 0.5 * (rPlus + rMinus);
164 cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
165 VDBC = VNBC - m_VnInf[pnt];
166
167 // Thermodynamic boundary variables
168 sBC = m_pInf / (pow(m_rhoInf, m_gamma));
169 rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
170 pBC = rhoBC * cBC * cBC * gammaInv;
171
172 // Kinetic energy initialiasation
173 NekDouble EkBC = 0.0;
174
175 // Boundary velocities
176 for (j = 0; j < nDimensions; ++j)
177 {
178 velBC[j] = m_velInf[j] + VDBC * m_traceNormals[j][pnt];
179 rhoVelBC[j] = rhoBC * velBC[j];
180 EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
181 }
182
183 // Boundary energy
184 EBC = pBC * gammaMinusOneInv + EkBC;
185
186 // Imposing Riemann Invariant boundary conditions
187 (m_fields[0]
188 ->GetBndCondExpansions()[m_bcRegion]
189 ->UpdatePhys())[id1 + i] = rhoBC;
190 for (j = 0; j < nDimensions; ++j)
191 {
192 (m_fields[j + 1]
193 ->GetBndCondExpansions()[m_bcRegion]
194 ->UpdatePhys())[id1 + i] = rhoVelBC[j];
195 }
196 (m_fields[nDimensions + 1]
197 ->GetBndCondExpansions()[m_bcRegion]
198 ->UpdatePhys())[id1 + i] = EBC;
199 }
200 else // Impose outflow Riemann invariant
201 {
202 // Subsonic flows
203 if (Mach[pnt] < 1.00)
204 {
205 // + Characteristic from inside
206 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
207 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
208
209 // - Characteristic from boundary
210 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
211 rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
212 }
213 else
214 {
215 // + Characteristic from inside
216 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
217 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
218
219 // + Characteristic from inside
220 cMinus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
221 rMinus = Vn[pnt] - 2.0 * cPlus * gammaMinusOneInv;
222 }
223
224 // Riemann boundary variables
225 VNBC = 0.5 * (rPlus + rMinus);
226 cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
227 VDBC = VNBC - Vn[pnt];
228
229 // Thermodynamic boundary variables
230 sBC = pressure[pnt] / (pow(Fwd[0][pnt], m_gamma));
231 rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
232 pBC = rhoBC * cBC * cBC * gammaInv;
233
234 // Kinetic energy initialiasation
235 NekDouble EkBC = 0.0;
236
237 // Boundary velocities
238 for (j = 0; j < nDimensions; ++j)
239 {
240 velBC[j] = Fwd[j + 1][pnt] / Fwd[0][pnt] +
241 VDBC * m_traceNormals[j][pnt];
242 rhoVelBC[j] = rhoBC * velBC[j];
243 EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
244 }
245
246 // Boundary energy
247 EBC = pBC * gammaMinusOneInv + EkBC;
248
249 // Imposing Riemann Invariant boundary conditions
250 (m_fields[0]
251 ->GetBndCondExpansions()[m_bcRegion]
252 ->UpdatePhys())[id1 + i] = rhoBC;
253 for (j = 0; j < nDimensions; ++j)
254 {
255 (m_fields[j + 1]
256 ->GetBndCondExpansions()[m_bcRegion]
257 ->UpdatePhys())[id1 + i] = rhoVelBC[j];
258 }
259 (m_fields[nDimensions + 1]
260 ->GetBndCondExpansions()[m_bcRegion]
261 ->UpdatePhys())[id1 + i] = EBC;
262 }
263 }
264 }
265}
266
267} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
NekDouble m_rhoInf
Definition: CFSBndCond.h:103
NekDouble m_pInf
Definition: CFSBndCond.h:104
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:106
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
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
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
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.
RiemannInvariantBC(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)
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294