Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::RiemannInvariantBC:
Collaboration graph
[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...
 

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::ExpListSharedPtr
m_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_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
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...
 

Detailed Description

Outflow characteristic boundary conditions for compressible flow problems.

Definition at line 49 of file RiemannInvariantBC.h.

Constructor & Destructor Documentation

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

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

55  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56 {
57  // Calculate VnInf
58  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
59  m_VnInf = Array<OneD, NekDouble> (nTracePts, 0.0);
60 
61  // Computing the normal velocity for characteristics coming
62  // from outside the computational domain
63  for( int i =0; i < m_spacedim; i++)
64  {
65  Vmath::Svtvp(nTracePts, m_velInf[i],
66  m_traceNormals[i], 1,
67  m_VnInf, 1,
68  m_VnInf, 1);
69  }
70 }
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:51
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:90
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:485
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:86
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:88
virtual Nektar::RiemannInvariantBC::~RiemannInvariantBC ( void  )
inlineprivatevirtual

Definition at line 89 of file RiemannInvariantBC.h.

89 {};

Member Function Documentation

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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

61  {
63  AllocateSharedPtr(pSession, pFields,
64  pTraceNormals, pSpaceDim, bcRegion, cnt);
65  return p;
66  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:49
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 72 of file RiemannInvariantBC.cpp.

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, Vmath::Vabs(), Vmath::Vdiv(), and Vmath::Vvtvp().

76 {
77  int i, j;
78  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
79  int nDimensions = m_spacedim;
80 
81  const Array<OneD, const int> &traceBndMap
82  = m_fields[0]->GetTraceBndMap();
83 
84  NekDouble gammaInv = 1.0 / m_gamma;
85  NekDouble gammaMinusOne = m_gamma - 1.0;
86  NekDouble gammaMinusOneInv = 1.0 / gammaMinusOne;
87 
88  // Computing the normal velocity for characteristics coming
89  // from inside the computational domain
90  Array<OneD, NekDouble > Vn (nTracePts, 0.0);
91  Array<OneD, NekDouble > Vel(nTracePts, 0.0);
92  for (i = 0; i < nDimensions; ++i)
93  {
94  Vmath::Vdiv(nTracePts, Fwd[i+1], 1, Fwd[0], 1, Vel, 1);
95  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Vel, 1, Vn, 1, Vn, 1);
96  }
97 
98  // Computing the absolute value of the velocity in order to compute the
99  // Mach number to decide whether supersonic or subsonic
100  Array<OneD, NekDouble > absVel(nTracePts, 0.0);
101  m_varConv->GetAbsoluteVelocity(Fwd, absVel);
102 
103  // Get speed of sound
104  Array<OneD, NekDouble > pressure (nTracePts);
105  Array<OneD, NekDouble > soundSpeed(nTracePts);
106 
107  m_varConv->GetPressure(Fwd, pressure);
108  m_varConv->GetSoundSpeed(Fwd, pressure, soundSpeed);
109 
110  // Get Mach
111  Array<OneD, NekDouble > Mach(nTracePts, 0.0);
112  Vmath::Vdiv(nTracePts, Vn, 1, soundSpeed, 1, Mach, 1);
113  Vmath::Vabs(nTracePts, Mach, 1, Mach, 1);
114 
115  // Auxiliary variables
116  int eMax;
117  int e, id1, id2, nBCEdgePts, pnt;
118  NekDouble cPlus, rPlus, cMinus, rMinus, VDBC, VNBC;
119  Array<OneD, NekDouble> velBC(nDimensions, 0.0);
120  Array<OneD, NekDouble> rhoVelBC(nDimensions, 0.0);
121  NekDouble rhoBC, EBC, cBC, sBC, pBC;
122 
123  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
124 
125  // Loop on m_bcRegions
126  for (e = 0; e < eMax; ++e)
127  {
128  nBCEdgePts = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
129  GetExp(e)->GetTotPoints();
130 
131  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
132  GetPhys_Offset(e);
133  id2 = 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]->GetBndCondExpansions()[m_bcRegion]->
191  UpdatePhys())[id1+i] = rhoBC;
192  for (j = 0; j < nDimensions; ++j)
193  {
194  (m_fields[j+1]->GetBndCondExpansions()[m_bcRegion]->
195  UpdatePhys())[id1+i] = rhoVelBC[j];
196  }
197  (m_fields[nDimensions+1]->GetBndCondExpansions()[m_bcRegion]->
198  UpdatePhys())[id1+i] = EBC;
199 
200  }
201  else // Impose outflow Riemann invariant
202  {
203  // Subsonic flows
204  if (Mach[pnt] < 1.00)
205  {
206  // + Characteristic from inside
207  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
208  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
209 
210  // - Characteristic from boundary
211  cMinus = sqrt(m_gamma * m_pInf / m_rhoInf);
212  rMinus = m_VnInf[pnt] - 2.0 * cMinus * gammaMinusOneInv;
213  }
214  else
215  {
216  // + Characteristic from inside
217  cPlus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
218  rPlus = Vn[pnt] + 2.0 * cPlus * gammaMinusOneInv;
219 
220  // + Characteristic from inside
221  cMinus = sqrt(m_gamma * pressure[pnt] / Fwd[0][pnt]);
222  rMinus = Vn[pnt] - 2.0 * cPlus * gammaMinusOneInv;
223  }
224 
225  // Riemann boundary variables
226  VNBC = 0.5 * (rPlus + rMinus);
227  cBC = 0.25 * gammaMinusOne * (rPlus - rMinus);
228  VDBC = VNBC - Vn[pnt];
229 
230  // Thermodynamic boundary variables
231  sBC = pressure[pnt] / (pow(Fwd[0][pnt], m_gamma));
232  rhoBC = pow((cBC * cBC) / (m_gamma * sBC), gammaMinusOneInv);
233  pBC = rhoBC * cBC * cBC * gammaInv;
234 
235  // Kinetic energy initialiasation
236  NekDouble EkBC = 0.0;
237 
238  // Boundary velocities
239  for ( j = 0; j < nDimensions; ++j)
240  {
241  velBC[j] = Fwd[j+1][pnt] / Fwd[0][pnt] +
242  VDBC * m_traceNormals[j][pnt];
243  rhoVelBC[j] = rhoBC * velBC[j];
244  EkBC += 0.5 * rhoBC * velBC[j]*velBC[j];
245  }
246 
247  // Boundary energy
248  EBC = pBC * gammaMinusOneInv + EkBC;
249 
250  // Imposing Riemann Invariant boundary conditions
251  (m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
252  UpdatePhys())[id1+i] = rhoBC;
253  for (j = 0; j < nDimensions; ++j)
254  {
255  (m_fields[j+1]->GetBndCondExpansions()[m_bcRegion]->
256  UpdatePhys())[id1+i] = rhoVelBC[j];
257  }
258  (m_fields[nDimensions+1]->GetBndCondExpansions()[m_bcRegion]->
259  UpdatePhys())[id1+i] = EBC;
260  }
261  }
262  }
263 }
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:90
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:442
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:241
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:86
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:424
NekDouble m_pInf
Definition: CFSBndCond.h:97
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
double NekDouble
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:95
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:92
Array< OneD, NekDouble > m_velInf
Definition: CFSBndCond.h:98
int m_offset
Offset.
Definition: CFSBndCond.h:103
NekDouble m_rhoInf
Definition: CFSBndCond.h:96
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:88
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:101

Friends And Related Function Documentation

friend class MemoryManager< RiemannInvariantBC >
friend

Definition at line 53 of file RiemannInvariantBC.h.

Member Data Documentation

std::string Nektar::RiemannInvariantBC::className
static
Initial value:
RegisterCreatorFunction("RiemannInvariant",
"Riemann invariant boundary condition.")

Name of the class.

Definition at line 69 of file RiemannInvariantBC.h.

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

Reference normal velocity.

Definition at line 79 of file RiemannInvariantBC.h.

Referenced by RiemannInvariantBC(), and v_Apply().