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

Stagnation conditions inflow boundary conditions for compressible flow problems where the energy and density are prescribed. More...

#include <StagnationInflowBC.h>

Inheritance diagram for Nektar::StagnationInflowBC:
[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

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

Private Attributes

Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
 

Friends

class MemoryManager< StagnationInflowBC >
 

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

Stagnation conditions inflow boundary conditions for compressible flow problems where the energy and density are prescribed.

Definition at line 47 of file StagnationInflowBC.h.

Constructor & Destructor Documentation

◆ StagnationInflowBC()

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

52 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
53{
54 const size_t nvariables = m_fields.size();
55 const int expdim = m_fields[0]->GetGraph()->GetMeshDimension();
56 const int spacedim = m_fields[0]->GetGraph()->GetSpaceDimension();
57 const int numBCPts =
58 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetNpoints();
59
60 // Allocate internal storage for boundary condition values
61 m_fieldStorage = Array<OneD, Array<OneD, NekDouble>>(nvariables);
62 for (int i = 0; i < nvariables; ++i)
63 {
64 m_fieldStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
65 Vmath::Vcopy(numBCPts,
66 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->GetPhys(),
67 1, m_fieldStorage[i], 1);
68 }
69
70 // Compute the flow direction that will be imposed
71 for (int j = 0; j < numBCPts; ++j)
72 {
73 // Compute norm of user-specified flow direction
74 NekDouble dirNorm = 0.;
75 for (int i = 0; i < spacedim; ++i)
76 {
77 dirNorm += std::pow(m_fieldStorage[i + 1][j], 2);
78 }
79 dirNorm = std::sqrt(dirNorm);
80
81 // Use the user-specified flow direction if it's nonzero
82 if (dirNorm > 1.E-8)
83 {
84 for (int i = 0; i < spacedim; ++i)
85 {
86 m_fieldStorage[i + 1][j] /= dirNorm;
87 }
88 }
89 else
90 {
91 for (int i = 0; i < expdim; ++i)
92 {
93 m_fieldStorage[i + 1][j] = -m_traceNormals[i][j];
94 }
95 for (int i = expdim; i < spacedim; ++i)
96 {
97 m_fieldStorage[i + 1][j] = 0.;
98 }
99 }
100 }
101}
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
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
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, Array< OneD, NekDouble > > m_fieldStorage
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, Nektar::CFSBndCond::m_traceNormals, tinysimd::sqrt(), and Vmath::Vcopy().

◆ ~StagnationInflowBC()

Nektar::StagnationInflowBC::~StagnationInflowBC ( void  )
inlineoverrideprivate

Definition at line 80 of file StagnationInflowBC.h.

80{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::StagnationInflowBC::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 StagnationInflowBC.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::StagnationInflowBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 103 of file StagnationInflowBC.cpp.

106{
107 const size_t nTracePts = Fwd[0].size();
108 const size_t nVariables = physarray.size();
109
110 ASSERTL0(nTracePts == m_fields[0]->GetTrace()->GetNpoints(),
111 "Number of trace points does not match in "
112 "StagnationInflowBC::v_Apply()");
113
114 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
115
116 // Get velocity magnitude from Fwd
117 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
118 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
119
120 // Loop over all elements on the boundary
121 for (int e = 0;
122 e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize(); ++e)
123 {
124 // Number of quadrature points in this element
125 const int npts = m_fields[0]
126 ->GetBndCondExpansions()[m_bcRegion]
127 ->GetExp(e)
128 ->GetTotPoints();
129 const int id1 =
130 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
131 const int id2 =
132 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
133
134 // Loop over all quadrature points on this element
135 for (int i = 0; i < npts; ++i)
136 {
137 // Copy conserved variables at stagnation state from BC
138 NekDouble rhoStag = m_fieldStorage[0][id1 + i];
139 NekDouble EStag = m_fieldStorage[nVariables - 1][id1 + i];
140
141 // Compute stagnation enthalpy
142 NekDouble hStag = m_gamma * EStag / rhoStag;
143
144 // Compute static enthalpy by solving 1D energy balance
145 NekDouble hStat = hStag - 0.5 * pow(absVel[id2 + i], 2);
146
147 // Compute density from isentropic flow relation
148 NekDouble rho = rhoStag * pow(rhoStag * hStat / (EStag * m_gamma),
149 1. / (m_gamma - 1.));
150
151 // Update density for BC
152 (m_fields[0]
153 ->GetBndCondExpansions()[m_bcRegion]
154 ->UpdatePhys())[id1 + i] = rho;
155
156 // Update momentum for BC
157 for (int n = 0; n < m_spacedim; ++n)
158 {
159 (m_fields[1 + n]
160 ->GetBndCondExpansions()[m_bcRegion]
161 ->UpdatePhys())[id1 + i] =
162 rho * absVel[id2 + i] * m_fieldStorage[1 + n][id1 + i];
163 }
164
165 // Update total energy for BC
166 (m_fields[nVariables - 1]
167 ->GetBndCondExpansions()[m_bcRegion]
168 ->UpdatePhys())[id1 + i] =
169 rho * (hStat / m_gamma + 0.5 * pow(absVel[id2 + i], 2));
170 }
171 }
172}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:97
int m_offset
Offset.
Definition: CFSBndCond.h:111

References ASSERTL0, Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, Nektar::CFSBndCond::m_gamma, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_spacedim, and Nektar::CFSBndCond::m_varConv.

Friends And Related Function Documentation

◆ MemoryManager< StagnationInflowBC >

friend class MemoryManager< StagnationInflowBC >
friend

Definition at line 1 of file StagnationInflowBC.h.

Member Data Documentation

◆ className

std::string Nektar::StagnationInflowBC::className
static
Initial value:
=
"StagnationInflow", StagnationInflowBC::create,
"Stagnation conditions inflow 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 StagnationInflowBC.h.

◆ m_fieldStorage

Array<OneD, Array<OneD, NekDouble> > Nektar::StagnationInflowBC::m_fieldStorage
private

Definition at line 83 of file StagnationInflowBC.h.

Referenced by StagnationInflowBC(), and v_Apply().