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

56  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
57 {
58  // Calculate VnInf
59  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
60  m_VnInf = Array<OneD, NekDouble> (nTracePts, 0.0);
61 
62  // Computing the normal velocity for characteristics coming
63  // from outside the computational domain
64  for( int i =0; i < m_spacedim; i++)
65  {
66  Vmath::Svtvp(nTracePts, m_velInf[i],
67  m_traceNormals[i], 1,
68  m_VnInf, 1,
69  m_VnInf, 1);
70  }
71 }
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.cpp:565

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 88 of file RiemannInvariantBC.h.

88 {};

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 55 of file RiemannInvariantBC.h.

60  {
62  AllocateSharedPtr(pSession, pFields,
63  pTraceNormals, pSpaceDim, 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: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 
)
protectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 73 of file RiemannInvariantBC.cpp.

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

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:
RegisterCreatorFunction("RiemannInvariant",
"Riemann invariant boundary condition.")
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 68 of file RiemannInvariantBC.h.

◆ m_VnInf

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

Reference normal velocity.

Definition at line 78 of file RiemannInvariantBC.h.

Referenced by RiemannInvariantBC(), and v_Apply().