Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::PressureInflowFileBC Class Reference

Pressure inflow boundary conditions for compressible flow problems where either the density and the velocities are assigned from a file or the full state is assigned from a file (depending on the problem type, either subsonic or supersonic). More...

#include <PressureInflowFileBC.h>

Inheritance diagram for Nektar::PressureInflowFileBC:
Inheritance graph
[legend]
Collaboration diagram for Nektar::PressureInflowFileBC:
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...
 

Private Member Functions

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

Private Attributes

Array< OneD, Array< OneD,
NekDouble > > 
m_fieldStorage
 

Friends

class MemoryManager< PressureInflowFileBC >
 

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

Detailed Description

Pressure inflow boundary conditions for compressible flow problems where either the density and the velocities are assigned from a file or the full state is assigned from a file (depending on the problem type, either subsonic or supersonic).

Definition at line 51 of file PressureInflowFileBC.h.

Constructor & Destructor Documentation

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

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, and Vmath::Vcopy().

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56  int nvariables = m_fields.num_elements();
57  // Loop over Boundary Regions for PressureInflowFileBC
58  m_fieldStorage = Array<OneD, Array<OneD, NekDouble> > (nvariables);
59 
60  int numBCPts = m_fields[0]->
61  GetBndCondExpansions()[m_bcRegion]->GetNpoints();
62  for (int i = 0; i < nvariables; ++i)
63  {
64  m_fieldStorage[i] = Array<OneD, NekDouble>(numBCPts, 0.0);
66  numBCPts,
67  m_fields[i]->GetBndCondExpansions()[m_bcRegion]->GetPhys(), 1,
68  m_fieldStorage[i], 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
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:86
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:101
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage
virtual Nektar::PressureInflowFileBC::~PressureInflowFileBC ( void  )
inlineprivatevirtual

Definition at line 88 of file PressureInflowFileBC.h.

88 {};

Member Function Documentation

static CFSBndCondSharedPtr Nektar::PressureInflowFileBC::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 58 of file PressureInflowFileBC.h.

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

63  {
65  AllocateSharedPtr(pSession, pFields,
66  pTraceNormals, pSpaceDim, bcRegion, cnt);
67  return p;
68  }
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::PressureInflowFileBC::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 PressureInflowFileBC.cpp.

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, m_fieldStorage, Nektar::CFSBndCond::m_gamma, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, npts, Vmath::Vabs(), Vmath::Vdiv(), and Vmath::Vvtvp().

76 {
77  int i, j;
78  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
79  int nVariables = physarray.num_elements();
80  int nDimensions = m_spacedim;
81 
82  const Array<OneD, const int> &traceBndMap
83  = m_fields[0]->GetTraceBndMap();
84 
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 e, id1, id2, npts, pnt;
117  NekDouble rhoeb;
118 
119  // Loop on the m_bcRegions
120  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
121  GetExpSize(); ++e)
122  {
123  npts = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
124  GetExp(e)->GetTotPoints();
125  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
126  GetPhys_Offset(e);
127  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset+e]);
128 
129  // Loop on points of m_bcRegion 'e'
130  for (i = 0; i < npts; i++)
131  {
132  pnt = id2+i;
133 
134  // Subsonic flows
135  if (Mach[pnt] < 0.99)
136  {
137  // Partial extrapolation for subsonic cases
138  for (j = 0; j < nVariables-1; ++j)
139  {
140  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
141  UpdatePhys())[id1+i] = m_fieldStorage[j][id1+i];
142  }
143 
144  // Kinetic energy calculation
145  NekDouble Ek = 0.0;
146  for (j = 1; j < nVariables-1; ++j)
147  {
148  Ek += 0.5 * (m_fieldStorage[j][id1+i] *
149  m_fieldStorage[j][id1+i]) /
150  m_fieldStorage[0][id1+i];
151  }
152  rhoeb = gammaMinusOneInv * pressure[pnt] + Ek;
153 
154  (m_fields[nVariables-1]->GetBndCondExpansions()[m_bcRegion]->
155  UpdatePhys())[id1+i] = rhoeb;
156  }
157  // Supersonic flows
158  else
159  {
160  for (j = 0; j < nVariables; ++j)
161  {
162  // Extrapolation for supersonic cases
163  (m_fields[j]->GetBndCondExpansions()[m_bcRegion]->
164  UpdatePhys())[id1+i] = Fwd[j][pnt];
165  }
166  }
167  }
168  }
169 }
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
static std::string npts
Definition: InputFld.cpp:43
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
int m_offset
Offset.
Definition: CFSBndCond.h:103
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
Array< OneD, Array< OneD, NekDouble > > m_fieldStorage

Friends And Related Function Documentation

friend class MemoryManager< PressureInflowFileBC >
friend

Definition at line 55 of file PressureInflowFileBC.h.

Member Data Documentation

std::string Nektar::PressureInflowFileBC::className
static
Initial value:
RegisterCreatorFunction("PressureInflowFile",
"Pressure inflow (file) boundary condition.")

Name of the class.

Definition at line 71 of file PressureInflowFileBC.h.

Array<OneD, Array<OneD, NekDouble> > Nektar::PressureInflowFileBC::m_fieldStorage
private

Definition at line 88 of file PressureInflowFileBC.h.

Referenced by PressureInflowFileBC(), and v_Apply().