Nektar++
StagnationInflowBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StagnationInflowBC.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: Stagnation conditions inflow boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "StagnationInflowBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
44 "StagnationInflow", StagnationInflowBC::create,
45 "Stagnation conditions inflow 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 const size_t nvariables = m_fields.size();
55 const int expdim = m_fields[0]->GetGraph()->GetMeshDimension();
56 const int spacedim = m_fields[0]->GetGraph()->GetSpaceDimension();
57 const int numBCPts =
58 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
59
60 // Allocate internal storage for boundary condition values
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 // Compute the flow direction that will be imposed
71 for (int j = 0; j < numBCPts; ++j)
72 {
73 // Compute norm of user-specified flow direction
74 NekDouble dirNorm = 0.;
75 for (int i = 0; i < spacedim; ++i)
76 {
77 dirNorm += std::pow(m_fieldStorage[i + 1][j], 2);
78 }
79 dirNorm = std::sqrt(dirNorm);
80
81 // Use the user-specified flow direction if it's nonzero
82 if (dirNorm > 1.E-8)
83 {
84 for (int i = 0; i < spacedim; ++i)
85 {
86 m_fieldStorage[i + 1][j] /= dirNorm;
87 }
88 }
89 else
90 {
91 for (int i = 0; i < expdim; ++i)
92 {
93 m_fieldStorage[i + 1][j] = -m_traceNormals[i][j];
94 }
95 for (int i = expdim; i < spacedim; ++i)
96 {
97 m_fieldStorage[i + 1][j] = 0.;
98 }
99 }
100 }
101}
102
104 Array<OneD, Array<OneD, NekDouble>> &physarray,
105 [[maybe_unused]] const NekDouble &time)
106{
107 const size_t nTracePts = Fwd[0].size();
108 const size_t nVariables = physarray.size();
109
110 ASSERTL0(nTracePts == m_fields[0]->GetTrace()->GetNpoints(),
111 "Number of trace points does not match in "
112 "StagnationInflowBC::v_Apply()");
113
114 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
115
116 // Get velocity magnitude from Fwd
117 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
118 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
119
120 // Loop over all elements on the boundary
121 for (int e = 0;
122 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
123 {
124 // Number of quadrature points in this element
125 const int npts = m_fields[0]
126 ->GetBndCondExpansions()[m_bcRegion]
127 ->GetExp(e)
128 ->GetTotPoints();
129 const int id1 =
130 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
131 const int id2 =
132 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
133
134 // Loop over all quadrature points on this element
135 for (int i = 0; i < npts; ++i)
136 {
137 // Copy conserved variables at stagnation state from BC
138 NekDouble rhoStag = m_fieldStorage[0][id1 + i];
139 NekDouble EStag = m_fieldStorage[nVariables - 1][id1 + i];
140
141 // Compute stagnation enthalpy
142 NekDouble hStag = m_gamma * EStag / rhoStag;
143
144 // Compute static enthalpy by solving 1D energy balance
145 NekDouble hStat = hStag - 0.5 * pow(absVel[id2 + i], 2);
146
147 // Compute density from isentropic flow relation
148 NekDouble rho = rhoStag * pow(rhoStag * hStat / (EStag * m_gamma),
149 1. / (m_gamma - 1.));
150
151 // Update density for BC
152 (m_fields[0]
153 ->GetBndCondExpansions()[m_bcRegion]
154 ->UpdatePhys())[id1 + i] = rho;
155
156 // Update momentum for BC
157 for (int n = 0; n < m_spacedim; ++n)
158 {
159 (m_fields[1 + n]
160 ->GetBndCondExpansions()[m_bcRegion]
161 ->UpdatePhys())[id1 + i] =
162 rho * absVel[id2 + i] * m_fieldStorage[1 + n][id1 + i];
163 }
164
165 // Update total energy for BC
166 (m_fields[nVariables - 1]
167 ->GetBndCondExpansions()[m_bcRegion]
168 ->UpdatePhys())[id1 + i] =
169 rho * (hStat / m_gamma + 0.5 * pow(absVel[id2 + i], 2));
170 }
171 }
172}
173
174} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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
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
StagnationInflowBC(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)
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.
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.
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
std::shared_ptr< SessionReader > SessionReaderSharedPtr
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294