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

Pressure outflow boundary conditions for compressible flow problems. More...

#include <PressureOutflowBC.h>

Inheritance diagram for Nektar::PressureOutflowBC:
[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)
 
- 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...
 

Private Member Functions

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

Private Attributes

Array< OneD, NekDoublem_pressureStorage
 

Friends

class MemoryManager< PressureOutflowBC >
 

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...
 
- 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_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
Array< OneD, NekDoublem_velInf
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Detailed Description

Pressure outflow boundary conditions for compressible flow problems.

Definition at line 48 of file PressureOutflowBC.h.

Constructor & Destructor Documentation

◆ PressureOutflowBC()

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

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

56  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
57 {
58  int numBCPts = m_fields[0]->
59  GetBndCondExpansions()[m_bcRegion]->GetNpoints();
60  m_pressureStorage = Array<OneD, NekDouble>(numBCPts, 0.0);
61 
62  // Get Pressure
63  Vmath::Vcopy(numBCPts,
64  m_fields[m_spacedim+1]->GetBndCondExpansions()[m_bcRegion]->GetPhys(), 1,
66 }
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_spacedim
Space dimension.
Definition: CFSBndCond.h:89
Array< OneD, NekDouble > m_pressureStorage
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:85
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:100

◆ ~PressureOutflowBC()

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

Definition at line 85 of file PressureOutflowBC.h.

References m_pressureStorage.

85 {};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::PressureOutflowBC::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 PressureOutflowBC.h.

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

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

◆ v_Apply()

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

Implements Nektar::CFSBndCond.

Definition at line 68 of file PressureOutflowBC.cpp.

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

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

Friends And Related Function Documentation

◆ MemoryManager< PressureOutflowBC >

friend class MemoryManager< PressureOutflowBC >
friend

Definition at line 52 of file PressureOutflowBC.h.

Member Data Documentation

◆ className

std::string Nektar::PressureOutflowBC::className
static
Initial value:
RegisterCreatorFunction("PressureOutflow",
"Pressure outflow boundary condition.")

Name of the class.

Definition at line 68 of file PressureOutflowBC.h.

◆ m_pressureStorage

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

Definition at line 85 of file PressureOutflowBC.h.

Referenced by PressureOutflowBC(), v_Apply(), and ~PressureOutflowBC().