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 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 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 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...
 
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
 
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 int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
private

Definition at line 47 of file PressureInflowFileBC.cpp.

52 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
53{
54 int nvariables = m_fields.size();
55 // Loop over Boundary Regions for PressureInflowFileBC
56 m_fieldStorage = Array<OneD, Array<OneD, NekDouble>>(nvariables);
57
58 int numBCPts =
59 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
60 for (int i = 0; i < nvariables; ++i)
61 {
62 m_fieldStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
63 Vmath::Vcopy(numBCPts,
64 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->GetPhys(),
65 1, m_fieldStorage[i], 1);
66 }
67}
CFSBndCond(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)
Constructor.
Definition: CFSBndCond.cpp:47
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:91
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 82 of file PressureInflowFileBC.h.

82{};

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 int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
inlinestatic

Creates an instance of this class.

Definition at line 55 of file PressureInflowFileBC.h.

60 {
63 pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt);
64 return p;
65 }
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 69 of file PressureInflowFileBC.cpp.

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

Name of the class.

Definition at line 68 of file PressureInflowFileBC.h.

◆ m_fieldStorage

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

Definition at line 85 of file PressureInflowFileBC.h.

Referenced by PressureInflowFileBC(), and v_Apply().