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

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

Private Member Functions

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

Definition at line 47 of file RiemannInvariantBC.cpp.

53 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
54 bcRegion, cnt)
55{
56 // Calculate VnInf
57 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
58 m_VnInf = Array<OneD, NekDouble>(nTracePts, 0.0);
59
60 // Computing the normal velocity for characteristics coming
61 // from outside the computational domain
62 for (int i = 0; i < m_spacedim; i++)
63 {
64 Vmath::Svtvp(nTracePts, m_velInf[i], m_traceNormals[i], 1, m_VnInf, 1,
65 m_VnInf, 1);
66 }
67}
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:94
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:109
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_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 86 of file RiemannInvariantBC.h.

86{};

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 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 RiemannInvariantBC.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::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 69 of file RiemannInvariantBC.cpp.

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

◆ m_VnInf

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

Reference normal velocity.

Definition at line 76 of file RiemannInvariantBC.h.

Referenced by RiemannInvariantBC(), and v_Apply().