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 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

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

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

Definition at line 47 of file PressureOutflowNonReflectiveBC.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
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().

◆ ~PressureOutflowNonReflectiveBC()

Nektar::PressureOutflowNonReflectiveBC::~PressureOutflowNonReflectiveBC ( void  )
inlineoverrideprivate

Definition at line 84 of file PressureOutflowNonReflectiveBC.h.

84{};

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 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 54 of file PressureOutflowNonReflectiveBC.h.

60 {
63 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
64 bcRegion, cnt);
65 return p;
66 }
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::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.

71{
72 int i, j;
73 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
74 int nVariables = physarray.size();
75 int nDimensions = m_spacedim;
76
77 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
78
79 // Computing the normal velocity for characteristics coming
80 // from inside the computational domain
81 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
82 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
83 for (i = 0; i < nDimensions; ++i)
84 {
85 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
86 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
87 }
88
89 // Computing the absolute value of the velocity in order to compute the
90 // Mach number to decide whether supersonic or subsonic
91 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
92 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
93
94 // Get speed of sound
95 Array<OneD, NekDouble> soundSpeed(nTracePts);
96 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
97
98 // Get Mach
99 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
100 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
101 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
102
103 // Auxiliary variables
104 int e, id1, id2, npts, pnt;
105 NekDouble rhoeb;
106
107 // Loop on the m_bcRegions
108 for (e = 0;
109 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
110 {
111 npts = m_fields[0]
112 ->GetBndCondExpansions()[m_bcRegion]
113 ->GetExp(e)
114 ->GetTotPoints();
115 id1 =
116 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
117 id2 =
118 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]
147 ->GetBndCondExpansions()[m_bcRegion]
148 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
149 }
150
151 (m_fields[nVariables - 1]
152 ->GetBndCondExpansions()[m_bcRegion]
153 ->UpdatePhys())[id1 + i] =
154 2.0 * rhoeb - Fwd[nVariables - 1][pnt];
155 }
156 // Supersonic flows
157 else
158 {
159 for (j = 0; j < nVariables; ++j)
160 {
161 // Extrapolation for supersonic cases
162 (m_fields[j]
163 ->GetBndCondExpansions()[m_bcRegion]
164 ->UpdatePhys())[id1 + i] = Fwd[j][pnt];
165 }
166 }
167 }
168 }
169}
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< 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.
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 69 of file PressureOutflowNonReflectiveBC.h.

◆ m_pressureStorage

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

Definition at line 87 of file PressureOutflowNonReflectiveBC.h.

Referenced by PressureOutflowNonReflectiveBC(), and v_Apply().