Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::PressureInflowFileBC Class Reference

Pressure inflow boundary conditions for compressible flow problems where either the density and the velocities are assigned from a file or the full state is assigned from a file (depending on the problem type, either subsonic or supersonic). More...

#include <PressureInflowFileBC.h>

Inheritance diagram for Nektar::PressureInflowFileBC:
[legend]

Static Public Member Functions

static CFSBndCondSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
 
- Protected Member Functions inherited from Nektar::CFSBndCond
 CFSBndCond (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 Constructor. More...
 
virtual void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)=0
 
virtual void v_ApplyBwdWeight ()
 

Private Member Functions

 PressureInflowFileBC (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
 
 ~PressureInflowFileBC (void) override
 

Private Attributes

Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
 

Friends

class MemoryManager< PressureInflowFileBC >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::CFSBndCond
virtual ~CFSBndCond ()
 
void Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time=0)
 Apply the boundary condition. More...
 
void ApplyBwdWeight ()
 Apply the Weight of boundary condition. More...
 
- Protected Attributes inherited from Nektar::CFSBndCond
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array of fields. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Trace normals. More...
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 Grid Velocity. More...
 
int m_spacedim
 Space dimension. More...
 
VariableConverterSharedPtr m_varConv
 Auxiliary object to convert variables. More...
 
NekDouble m_diffusionAveWeight
 Weight for average calculation of diffusion term. More...
 
NekDouble m_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
NekDouble m_pOut
 
Array< OneD, NekDoublem_velInf
 
NekDouble m_angVel
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Detailed Description

Pressure inflow boundary conditions for compressible flow problems where either the density and the velocities are assigned from a file or the full state is assigned from a file (depending on the problem type, either subsonic or supersonic).

Definition at line 49 of file PressureInflowFileBC.h.

Constructor & Destructor Documentation

◆ PressureInflowFileBC()

Nektar::PressureInflowFileBC::PressureInflowFileBC ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, Array< OneD, NekDouble > > &  pTraceNormals,
const Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
private

Definition at line 47 of file PressureInflowFileBC.cpp.

53 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
54 bcRegion, cnt)
55{
56 int nvariables = m_fields.size();
57 // Loop over Boundary Regions for PressureInflowFileBC
58 m_fieldStorage = Array<OneD, Array<OneD, NekDouble>>(nvariables);
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}
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
CFSBndCond(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Constructor.
Definition: CFSBndCond.cpp:47
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, and Vmath::Vcopy().

◆ ~PressureInflowFileBC()

Nektar::PressureInflowFileBC::~PressureInflowFileBC ( void  )
inlineoverrideprivate

Definition at line 85 of file PressureInflowFileBC.h.

85{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::PressureInflowFileBC::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const Array< OneD, Array< OneD, NekDouble > > &  pTraceNormals,
const Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
inlinestatic

Creates an instance of this class.

Definition at line 55 of file PressureInflowFileBC.h.

61 {
64 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
65 bcRegion, cnt);
66 return p;
67 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:51

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_Apply()

void Nektar::PressureInflowFileBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 71 of file PressureInflowFileBC.cpp.

75{
76 int i, j;
77 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
78 int nVariables = physarray.size();
79 int nDimensions = m_spacedim;
80
81 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
82
83 // Computing the normal velocity for characteristics coming
84 // from inside the computational domain
85 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
86 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
87 for (i = 0; i < nDimensions; ++i)
88 {
89 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
90 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
91 }
92
93 // Computing the absolute value of the velocity in order to compute the
94 // Mach number to decide whether supersonic or subsonic
95 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
96 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
97
98 // Get speed of sound
99 Array<OneD, NekDouble> pressure(nTracePts);
100 Array<OneD, NekDouble> soundSpeed(nTracePts);
101
102 m_varConv->GetPressure(Fwd, pressure);
103 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
104
105 // Get Mach
106 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
107 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
108 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
109
110 // Auxiliary variables
111 int e, id1, id2, npts, pnt;
112 NekDouble rhoeb;
113
114 // Loop on the m_bcRegions
115 for (e = 0;
116 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
117 {
118 npts = m_fields[0]
119 ->GetBndCondExpansions()[m_bcRegion]
120 ->GetExp(e)
121 ->GetTotPoints();
122 id1 =
123 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
124 id2 =
125 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
126
127 // Get internal energy
128 Array<OneD, NekDouble> tmpPressure(npts, pressure + id2);
129 Array<OneD, NekDouble> rho(npts, Fwd[0] + id2);
130 Array<OneD, NekDouble> Ei(npts);
131 m_varConv->GetEFromRhoP(rho, tmpPressure, Ei);
132
133 // Loop on points of m_bcRegion 'e'
134 for (i = 0; i < npts; i++)
135 {
136 pnt = id2 + i;
137
138 // Subsonic flows
139 if (Mach[pnt] < 0.99)
140 {
141 // Partial extrapolation for subsonic cases
142 for (j = 0; j < nVariables - 1; ++j)
143 {
144 (m_fields[j]
145 ->GetBndCondExpansions()[m_bcRegion]
146 ->UpdatePhys())[id1 + i] = m_fieldStorage[j][id1 + i];
147 }
148
149 // Kinetic energy calculation
150 NekDouble Ek = 0.0;
151 for (j = 1; j < nVariables - 1; ++j)
152 {
153 Ek += 0.5 *
154 (m_fieldStorage[j][id1 + i] *
155 m_fieldStorage[j][id1 + i]) /
156 m_fieldStorage[0][id1 + i];
157 }
158
159 rhoeb = Fwd[0][pnt] * Ei[i] + Ek;
160
161 (m_fields[nVariables - 1]
162 ->GetBndCondExpansions()[m_bcRegion]
163 ->UpdatePhys())[id1 + i] = rhoeb;
164 }
165 // Supersonic flows
166 else
167 {
168 for (j = 0; j < nVariables; ++j)
169 {
170 // Extrapolation for supersonic cases
171 (m_fields[j]
172 ->GetBndCondExpansions()[m_bcRegion]
173 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
174 }
175 }
176 }
177 }
178}
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:94
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:100
int m_offset
Offset.
Definition: CFSBndCond.h:115
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

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, CG_Iterations::pressure, Vmath::Vabs(), Vmath::Vdiv(), and Vmath::Vvtvp().

Friends And Related Function Documentation

◆ MemoryManager< PressureInflowFileBC >

friend class MemoryManager< PressureInflowFileBC >
friend

Definition at line 1 of file PressureInflowFileBC.h.

Member Data Documentation

◆ className

std::string Nektar::PressureInflowFileBC::className
static
Initial value:
=
"PressureInflowFile", PressureInflowFileBC::create,
"Pressure inflow (file) boundary condition.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Creates an instance of this class.
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Name of the class.

Definition at line 70 of file PressureInflowFileBC.h.

◆ m_fieldStorage

Array<OneD, Array<OneD, NekDouble> > Nektar::PressureInflowFileBC::m_fieldStorage
private

Definition at line 88 of file PressureInflowFileBC.h.

Referenced by PressureInflowFileBC(), and v_Apply().