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 Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, 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 Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, 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 Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, 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...
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 Grid Velocity. 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
 
NekDouble m_angVel
 
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 Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
private

Definition at line 47 of file PressureOutflowBC.cpp.

53 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
54 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:98
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
CFSBndCond(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Constructor.
Definition: CFSBndCond.cpp:47
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
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 83 of file PressureOutflowBC.h.

83{};

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 Array< OneD, Array< OneD, NekDouble > > &  pGridVelocity,
const int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file PressureOutflowBC.h.

59 {
62 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
63 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: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 67 of file PressureOutflowBC.cpp.

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

◆ m_pressureStorage

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

Definition at line 86 of file PressureOutflowBC.h.

Referenced by PressureOutflowBC(), and v_Apply().