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

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

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

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

Definition at line 47 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 47 of file PressureOutflowBC.cpp.

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

◆ ~PressureOutflowBC()

Nektar::PressureOutflowBC::~PressureOutflowBC ( void  )
inlineoverrideprivate

Definition at line 80 of file PressureOutflowBC.h.

80{};

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

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

Implements Nektar::CFSBndCond.

Definition at line 65 of file PressureOutflowBC.cpp.

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

friend class MemoryManager< PressureOutflowBC >
friend

Definition at line 1 of file PressureOutflowBC.h.

Member Data Documentation

◆ className

std::string Nektar::PressureOutflowBC::className
static
Initial value:
=
"PressureOutflow", PressureOutflowBC::create,
"Pressure outflow 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 66 of file PressureOutflowBC.h.

◆ m_pressureStorage

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

Definition at line 83 of file PressureOutflowBC.h.

Referenced by PressureOutflowBC(), and v_Apply().