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

Pressure outflow non-reflective boundary conditions for compressible flow problems. More...

#include <PressureOutflowNonReflectiveBC.h>

Inheritance diagram for Nektar::PressureOutflowNonReflectiveBC:
[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

 PressureOutflowNonReflectiveBC (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 ~PressureOutflowNonReflectiveBC (void)
 

Private Attributes

Array< OneD, NekDoublem_pressureStorage
 

Friends

class MemoryManager< PressureOutflowNonReflectiveBC >
 

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 outflow non-reflective boundary conditions for compressible flow problems.

Definition at line 48 of file PressureOutflowNonReflectiveBC.h.

Constructor & Destructor Documentation

◆ PressureOutflowNonReflectiveBC()

Nektar::PressureOutflowNonReflectiveBC::PressureOutflowNonReflectiveBC ( 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 PressureOutflowNonReflectiveBC.cpp.

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56  int numBCPts =
57  m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
58  m_pressureStorage = Array<OneD, NekDouble>(numBCPts, 0.0);
59 
60  // Get Pressure
62  numBCPts,
63  m_fields[m_spacedim + 1]->GetBndCondExpansions()[m_bcRegion]->GetPhys(),
64  1, m_pressureStorage, 1);
65 }
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
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
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_pressureStorage, Nektar::CFSBndCond::m_spacedim, and Vmath::Vcopy().

◆ ~PressureOutflowNonReflectiveBC()

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

Definition at line 81 of file PressureOutflowNonReflectiveBC.h.

81 {};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::PressureOutflowNonReflectiveBC::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 54 of file PressureOutflowNonReflectiveBC.h.

59  {
62  pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt);
63  return p;
64  }
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::PressureOutflowNonReflectiveBC::v_Apply ( Array< OneD, Array< OneD, NekDouble >> &  Fwd,
Array< OneD, Array< OneD, NekDouble >> &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 67 of file PressureOutflowNonReflectiveBC.cpp.

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

Definition at line 1 of file PressureOutflowNonReflectiveBC.h.

Member Data Documentation

◆ className

std::string Nektar::PressureOutflowNonReflectiveBC::className
static
Initial value:
=
"PressureOutflowNonReflective", PressureOutflowNonReflectiveBC::create,
"Pressure outflow non-reflective 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 67 of file PressureOutflowNonReflectiveBC.h.

◆ m_pressureStorage

Array<OneD, NekDouble> Nektar::PressureOutflowNonReflectiveBC::m_pressureStorage
private

Definition at line 84 of file PressureOutflowNonReflectiveBC.h.

Referenced by PressureOutflowNonReflectiveBC(), and v_Apply().