Nektar++
PressureMachTemperatureBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PressureMachTemperatureBC.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: Boundary condition specified in terms of pressure, Mach number
32// and temperature
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
40namespace Nektar
41{
42
45 "PressureMachTemperature", PressureMachTemperatureBC::create,
46 "BC prescribed in terms of p, Ma and T.");
47
51 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
52 const int pSpaceDim, const int bcRegion, const int cnt)
53 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
54{
55 int nvariables = m_fields.size();
56 int numBCPts =
57 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
58
59 // Array for storing conserved variables on the boundary
61 for (int i = 0; i < nvariables; ++i)
62 {
63 m_bcStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
64 }
65
66 // We assume that the pressure is given in entry [0] of
67 // the BC ("rho" position) and the temperature in entry m_spacedim+1
68 // ("E" position)
70 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys();
71 const Array<OneD, const NekDouble> temperature =
72 m_fields[m_spacedim + 1]->GetBndCondExpansions()[m_bcRegion]->GetPhys();
73
74 // Calculate density
75 m_varConv->GetRhoFromPT(pressure, temperature, m_bcStorage[0]);
76 // Calculate the internal energy times density
77 m_varConv->GetEFromRhoP(m_bcStorage[0], pressure,
79 Vmath::Vmul(numBCPts, m_bcStorage[m_spacedim + 1], 1, m_bcStorage[0], 1,
80 m_bcStorage[m_spacedim + 1], 1);
81 // We can now obtain the sound speed at this (rho,e) condition
82 Array<OneD, NekDouble> soundSpeed(numBCPts);
83 m_varConv->GetSoundSpeed(m_bcStorage, soundSpeed);
84
85 // Now update momentum and add kinetic energy to E
86 Array<OneD, NekDouble> tmp(numBCPts);
87 for (int i = 0; i < m_spacedim; ++i)
88 {
89 // tmp = velocity in i direction
91 numBCPts,
92 m_fields[i + 1]->GetBndCondExpansions()[m_bcRegion]->GetPhys(), 1,
93 soundSpeed, 1, tmp, 1);
94 // rho*u
95 Vmath::Vmul(numBCPts, m_bcStorage[0], 1, tmp, 1, m_bcStorage[i + 1], 1);
96 // tmp = 0.5 * rho *(rhou) in vel
97 Vmath::Vmul(numBCPts, m_bcStorage[i + 1], 1, tmp, 1, tmp, 1);
98 Vmath::Smul(numBCPts, 0.5, tmp, 1, tmp, 1);
99 // Add to E
100 Vmath::Vadd(numBCPts, m_bcStorage[m_spacedim + 1], 1, tmp, 1,
101 m_bcStorage[m_spacedim + 1], 1);
102 }
103
104 // Copy to boundary condition
105 for (int i = 0; i < nvariables; ++i)
106 {
108 numBCPts, m_bcStorage[i], 1,
109 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys(), 1);
110 }
111}
112
114 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &Fwd,
115 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &physarray,
116 [[maybe_unused]] const NekDouble &time)
117{
118 int nvariables = m_fields.size();
119 int numBCPts =
120 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
121 // Copy conserved variables to boundary condition
122 for (int i = 0; i < nvariables; ++i)
123 {
125 numBCPts, m_bcStorage[i], 1,
126 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys(), 1);
127 }
128}
129
130} // 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
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:97
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
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
Array< OneD, Array< OneD, NekDouble > > m_bcStorage
PressureMachTemperatureBC(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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Vmul(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:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825