Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | 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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::PressureOutflowBC:
Collaboration graph
[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)
 

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::ExpListSharedPtr
m_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 49 of file PressureOutflowBC.h.

Constructor & Destructor Documentation

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 48 of file PressureOutflowBC.cpp.

55  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56 {
57 }
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:51
virtual Nektar::PressureOutflowBC::~PressureOutflowBC ( void  )
inlineprivatevirtual

Definition at line 86 of file PressureOutflowBC.h.

86 {};

Member Function Documentation

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 56 of file PressureOutflowBC.h.

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

61  {
63  AllocateSharedPtr(pSession, pFields,
64  pTraceNormals, pSpaceDim, bcRegion, cnt);
65  return p;
66  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:49
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 59 of file PressureOutflowBC.cpp.

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_gamma, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_pInf, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, npts, Vmath::Vabs(), Vmath::Vdiv(), and Vmath::Vvtvp().

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

Friends And Related Function Documentation

friend class MemoryManager< PressureOutflowBC >
friend

Definition at line 53 of file PressureOutflowBC.h.

Member Data Documentation

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

Name of the class.

Definition at line 69 of file PressureOutflowBC.h.