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

virtual 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_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)
 
virtual ~PressureInflowFileBC (void)
 

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 49 of file PressureInflowFileBC.cpp.

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, 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 }
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.cpp:1255

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

◆ ~PressureInflowFileBC()

virtual Nektar::PressureInflowFileBC::~PressureInflowFileBC ( void  )
inlineprivatevirtual

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:48

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.

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

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:198
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().