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

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

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

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, 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:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:106
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, 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.cpp:622

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

◆ ~RiemannInvariantBC()

virtual Nektar::RiemannInvariantBC::~RiemannInvariantBC ( void  )
inlineprivatevirtual

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

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.

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