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

Outflow characteristic boundary conditions for compressible flow problems. More...

#include <RiemannInvariantBC.h>

Inheritance diagram for Nektar::RiemannInvariantBC:
[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 ()
 

Protected Attributes

Array< OneD, NekDoublem_VnInf
 Reference normal velocity. 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...
 

Private Member Functions

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

Friends

class MemoryManager< RiemannInvariantBC >
 

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

Detailed Description

Outflow characteristic boundary conditions for compressible flow problems.

Definition at line 47 of file RiemannInvariantBC.h.

Constructor & Destructor Documentation

◆ RiemannInvariantBC()

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

52 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
53{
54 // Calculate VnInf
55 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
56 m_VnInf = Array<OneD, NekDouble>(nTracePts, 0.0);
57
58 // Computing the normal velocity for characteristics coming
59 // from outside the computational domain
60 for (int i = 0; i < m_spacedim; i++)
61 {
62 Vmath::Svtvp(nTracePts, m_velInf[i], m_traceNormals[i], 1, m_VnInf, 1,
63 m_VnInf, 1);
64 }
65}
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
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:106
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:91
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396

References Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_velInf, m_VnInf, and Vmath::Svtvp().

◆ ~RiemannInvariantBC()

Nektar::RiemannInvariantBC::~RiemannInvariantBC ( void  )
inlineoverrideprivate

Definition at line 83 of file RiemannInvariantBC.h.

83{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::RiemannInvariantBC::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 RiemannInvariantBC.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::RiemannInvariantBC::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 RiemannInvariantBC.cpp.

71{
72 int i, j;
73 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
74 int nDimensions = m_spacedim;
75
76 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
77
78 NekDouble gammaInv = 1.0 / m_gamma;
79 NekDouble gammaMinusOne = m_gamma - 1.0;
80 NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
81
82 // Computing the normal velocity for characteristics coming
83 // from inside the computational domain
84 Array<OneD, NekDouble> Vn(nTracePts, 0.0);
85 Array<OneD, NekDouble> Vel(nTracePts, 0.0);
86 for (i = 0; i < nDimensions; ++i)
87 {
88 Vmath::Vdiv(nTracePts, Fwd[i + 1], 1, Fwd[0], 1, Vel, 1);
89 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
90 }
91
92 // Computing the absolute value of the velocity in order to compute the
93 // Mach number to decide whether supersonic or subsonic
94 Array<OneD, NekDouble> absVel(nTracePts, 0.0);
95 m_varConv->GetAbsoluteVelocity(Fwd, absVel);
96
97 // Get speed of sound
98 Array<OneD, NekDouble> pressure(nTracePts);
99 Array<OneD, NekDouble> soundSpeed(nTracePts);
100
101 m_varConv->GetPressure(Fwd, pressure);
102 m_varConv->GetSoundSpeed(Fwd, soundSpeed);
103
104 // Get Mach
105 Array<OneD, NekDouble> Mach(nTracePts, 0.0);
106 Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
107 Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
108
109 // Auxiliary variables
110 int eMax;
111 int e, id1, id2, nBCEdgePts, pnt;
112 NekDouble cPlus, rPlus, cMinus, rMinus, VDBC, VNBC;
113 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
114 Array<OneD, NekDouble> rhoVelBC(nDimensions, 0.0);
115 NekDouble rhoBC, EBC, cBC, sBC, pBC;
116
117 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
118
119 // Loop on m_bcRegions
120 for (e = 0; e < eMax; ++e)
121 {
122 nBCEdgePts = m_fields[0]
123 ->GetBndCondExpansions()[m_bcRegion]
124 ->GetExp(e)
125 ->GetTotPoints();
126
127 id1 =
128 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
129 id2 =
130 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
131
132 // Loop on the points of the m_bcRegion
133 for (i = 0; i < nBCEdgePts; i++)
134 {
135 pnt = id2 + i;
136
137 // Impose inflow Riemann invariant
138 if (Vn[pnt] <= 0.0)
139 {
140 // Subsonic flows
141 if (Mach[pnt] < 1.00)
142 {
143 // + Characteristic from inside
144 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
145 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
146
147 // - Characteristic from boundary
148 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
149 rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
150 }
151 else
152 {
153 // + Characteristic from inside
154 cPlus = sqrt(m_gamma * m_pInf / m_rhoInf);
155 rPlus = m_VnInf[pnt] + 2.0 * cPlus * gammaMinusOneInv;
156
157 // + Characteristic from inside
158 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
159 rMinus = m_VnInf[pnt] - 2.0 * cPlus * gammaMinusOneInv;
160 }
161
162 // Riemann boundary variables
163 VNBC = 0.5 * (rPlus + rMinus);
164 cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
165 VDBC = VNBC - m_VnInf[pnt];
166
167 // Thermodynamic boundary variables
168 sBC = m_pInf / (pow(m_rhoInf, m_gamma));
169 rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
170 pBC = rhoBC * cBC * cBC * gammaInv;
171
172 // Kinetic energy initialiasation
173 NekDouble EkBC = 0.0;
174
175 // Boundary velocities
176 for (j = 0; j < nDimensions; ++j)
177 {
178 velBC[j] = m_velInf[j] + VDBC * m_traceNormals[j][pnt];
179 rhoVelBC[j] = rhoBC * velBC[j];
180 EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
181 }
182
183 // Boundary energy
184 EBC = pBC * gammaMinusOneInv + EkBC;
185
186 // Imposing Riemann Invariant boundary conditions
187 (m_fields[0]
188 ->GetBndCondExpansions()[m_bcRegion]
189 ->UpdatePhys())[id1 + i] = rhoBC;
190 for (j = 0; j < nDimensions; ++j)
191 {
192 (m_fields[j + 1]
193 ->GetBndCondExpansions()[m_bcRegion]
194 ->UpdatePhys())[id1 + i] = rhoVelBC[j];
195 }
196 (m_fields[nDimensions + 1]
197 ->GetBndCondExpansions()[m_bcRegion]
198 ->UpdatePhys())[id1 + i] = EBC;
199 }
200 else // Impose outflow Riemann invariant
201 {
202 // Subsonic flows
203 if (Mach[pnt] < 1.00)
204 {
205 // + Characteristic from inside
206 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
207 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
208
209 // - Characteristic from boundary
210 cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
211 rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
212 }
213 else
214 {
215 // + Characteristic from inside
216 cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
217 rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
218
219 // + Characteristic from inside
220 cMinus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
221 rMinus = Vn[pnt] - 2.0 * cPlus * gammaMinusOneInv;
222 }
223
224 // Riemann boundary variables
225 VNBC = 0.5 * (rPlus + rMinus);
226 cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
227 VDBC = VNBC - Vn[pnt];
228
229 // Thermodynamic boundary variables
230 sBC = pressure[pnt] / (pow(Fwd[0][pnt], m_gamma));
231 rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
232 pBC = rhoBC * cBC * cBC * gammaInv;
233
234 // Kinetic energy initialiasation
235 NekDouble EkBC = 0.0;
236
237 // Boundary velocities
238 for (j = 0; j < nDimensions; ++j)
239 {
240 velBC[j] = Fwd[j + 1][pnt] / Fwd[0][pnt] +
241 VDBC * m_traceNormals[j][pnt];
242 rhoVelBC[j] = rhoBC * velBC[j];
243 EkBC += 0.5 * rhoBC * velBC[j] * velBC[j];
244 }
245
246 // Boundary energy
247 EBC = pBC * gammaMinusOneInv + EkBC;
248
249 // Imposing Riemann Invariant boundary conditions
250 (m_fields[0]
251 ->GetBndCondExpansions()[m_bcRegion]
252 ->UpdatePhys())[id1 + i] = rhoBC;
253 for (j = 0; j < nDimensions; ++j)
254 {
255 (m_fields[j + 1]
256 ->GetBndCondExpansions()[m_bcRegion]
257 ->UpdatePhys())[id1 + i] = rhoVelBC[j];
258 }
259 (m_fields[nDimensions + 1]
260 ->GetBndCondExpansions()[m_bcRegion]
261 ->UpdatePhys())[id1 + i] = EBC;
262 }
263 }
264 }
265}
NekDouble m_rhoInf
Definition: CFSBndCond.h:103
NekDouble m_pInf
Definition: CFSBndCond.h:104
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_gamma, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_pInf, Nektar::CFSBndCond::m_rhoInf, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, Nektar::CFSBndCond::m_velInf, m_VnInf, CG_Iterations::pressure, tinysimd::sqrt(), Vmath::Vabs(), Vmath::Vdiv(), and Vmath::Vvtvp().

Friends And Related Function Documentation

◆ MemoryManager< RiemannInvariantBC >

friend class MemoryManager< RiemannInvariantBC >
friend

Definition at line 1 of file RiemannInvariantBC.h.

Member Data Documentation

◆ className

std::string Nektar::RiemannInvariantBC::className
static
Initial value:
=
"RiemannInvariant", RiemannInvariantBC::create,
"Riemann invariant 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 RiemannInvariantBC.h.

◆ m_VnInf

Array<OneD, NekDouble> Nektar::RiemannInvariantBC::m_VnInf
protected

Reference normal velocity.

Definition at line 74 of file RiemannInvariantBC.h.

Referenced by RiemannInvariantBC(), and v_Apply().