Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::RinglebFlowBC Class Reference

Wall boundary conditions for compressible flow problems. More...

#include <RinglebFlowBC.h>

Inheritance diagram for Nektar::RinglebFlowBC:
[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

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

Private Attributes

int m_expdim
 
bool m_homo1D
 

Friends

class MemoryManager< RinglebFlowBC >
 

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

Wall boundary conditions for compressible flow problems.

Definition at line 47 of file RinglebFlowBC.h.

Constructor & Destructor Documentation

◆ RinglebFlowBC()

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

References m_expdim, m_homo1D, and Nektar::CFSBndCond::m_session.

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56  m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
57 
58  m_homo1D = false;
59  if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
60  {
61  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
62  if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D")
63  || (HomoStr == "1D") || (HomoStr == "Homo1D"))
64  {
65  m_homo1D = true;
66  }
67  }
68 }
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
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: CFSBndCond.h:83

◆ ~RinglebFlowBC()

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

Definition at line 84 of file RinglebFlowBC.h.

References m_expdim.

84 {};

Member Function Documentation

◆ create()

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

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

59  {
61  AllocateSharedPtr(pSession, pFields,
62  pTraceNormals, pSpaceDim, bcRegion, cnt);
63  return p;
64  }
std::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:48
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ v_Apply()

void Nektar::RinglebFlowBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
protectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 70 of file RinglebFlowBC.cpp.

References Nektar::LibUtilities::eFunctionTypeFile, Nektar::CFSBndCond::m_bcRegion, m_expdim, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_gamma, m_homo1D, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_session, class_topology::P, and Vmath::Vcopy().

74 {
75  int nvariables = physarray.num_elements();
76 
77  // For 3DHomogenoeus1D
78  int n_planes = 1;
79  if (m_expdim == 2 && m_homo1D)
80  {
81  int nPointsTot = m_fields[0]->GetTotPoints();
82  int nPointsTot_plane = m_fields[0]->GetPlane(0)->GetTotPoints();
83  n_planes = nPointsTot/nPointsTot_plane;
84  }
85 
86  int id2, id2_plane, e_max;
87 
88  e_max = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
89 
90  for(int e = 0; e < e_max; ++e)
91  {
92  int npoints = m_fields[0]->
93  GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetTotPoints();
94  int id1 = m_fields[0]->
95  GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
96 
97  // For 3DHomogenoeus1D
98  if (m_expdim == 2 && m_homo1D)
99  {
100  int m_offset_plane = m_offset/n_planes;
101  int e_plane;
102  int e_max_plane = e_max/n_planes;
103  int nTracePts_plane = m_fields[0]->GetTrace()->GetNpoints();
104 
105  int planeID = floor((e + 0.5 )/ e_max_plane );
106  e_plane = e - e_max_plane*planeID;
107 
108  id2_plane = m_fields[0]->GetTrace()->GetPhys_Offset(
109  m_fields[0]->GetTraceMap()->
110  GetBndCondCoeffsToGlobalCoeffsMap(
111  m_offset_plane + e_plane));
112  id2 = id2_plane + planeID*nTracePts_plane;
113  }
114  else // For general case
115  {
116  id2 = m_fields[0]->
117  GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->
118  GetBndCondTraceToGlobalTraceMap(m_offset+e));
119  }
120 
121  Array<OneD,NekDouble> x0(npoints, 0.0);
122  Array<OneD,NekDouble> x1(npoints, 0.0);
123  Array<OneD,NekDouble> x2(npoints, 0.0);
124 
125  m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
126  GetExp(e)->GetCoords(x0, x1, x2);
127 
128  // Flow parameters
129  NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
130  NekDouble J11, J12, J21, J22, det;
131  NekDouble Fx, Fy;
132  NekDouble xi, yi;
133  NekDouble dV;
134  NekDouble dtheta;
135  NekDouble par1;
136  NekDouble theta = M_PI / 4.0;
137  NekDouble kExt = 0.7;
138  NekDouble V = kExt * sin(theta);
139  NekDouble toll = 1.0e-8;
140  NekDouble errV = 1.0;
141  NekDouble errTheta = 1.0;
142  NekDouble gamma = m_gamma;
143  NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
144 
145  // Loop on all the points of that edge
146  for (int j = 0; j < npoints; j++)
147  {
148 
149  while ((abs(errV) > toll) || (abs(errTheta) > toll))
150  {
151  VV = V * V;
152  sint = sin(theta);
153  c = sqrt(1.0 - gamma_1_2 * VV);
154  k = V / sint;
155  phi = 1.0 / k;
156  pp = phi * phi;
157  J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
158  1.0 / (5.0 * c * c * c * c * c) -
159  0.5 * log((1.0 + c) / (1.0 - c));
160 
161  r = pow(c, 1.0 / gamma_1_2);
162  xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
163  yi = phi / (r * V) * sqrt(1.0 - VV * pp);
164  par1 = 25.0 - 5.0 * VV;
165  ss = sint * sint;
166 
167  Fx = xi - x0[j];
168  Fy = yi - x1[j];
169 
170  J11 = 39062.5 / pow(par1, 3.5) *
171  (1.0 / VV - 2.0 / VV * ss) * V + 1562.5 /
172  pow(par1, 2.5) * (-2.0 / (VV * V) + 4.0 /
173  (VV * V) * ss) + 12.5 / pow(par1, 1.5) * V +
174  312.5 / pow(par1, 2.5) * V + 7812.5 /
175  pow(par1, 3.5) * V - 0.25 *
176  (-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 *
177  pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
178  pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
179  pow(par1, 0.5) * V) / (1.0 + 0.2 * pow(par1, 0.5)) *
180  (1.0 - 0.2 * pow(par1, 0.5));
181 
182  J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
183  J21 = -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
184  pow((1.0 - ss), 0.5) + 78125.0 / V * sint /
185  pow(par1, 3.5) * pow((1.0 - ss), 0.5);
186 
187  // the matrix is singular when theta = pi/2
188  if (abs(x1[j]) < toll && abs(cos(theta)) < toll)
189  {
190  J22 = -39062.5 / pow(par1, 3.5) / V + 3125 /
191  pow(par1, 2.5) / (VV * V) + 12.5 /
192  pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) *
193  V + 7812.5 / pow(par1, 3.5) * V - 0.25 *
194  (-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 *
195  pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
196  pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
197  pow(par1, 0.5) * V) / (1.0 + 0.2 *
198  pow(par1, 0.5)) * (1.0 - 0.2 * pow(par1, 0.5));
199 
200  // dV = -dV/dx * Fx
201  dV = -1.0 / J22 * Fx;
202  dtheta = 0.0;
203  theta = M_PI / 2.0;
204  }
205  else
206  {
207  J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
208  pow((1.0 - ss), 0.5) - 3125.0 / VV * ss /
209  pow(par1, 2.5) / pow((1.0 - ss), 0.5) *
210  cos(theta);
211 
212  det = -1.0 / (J11 * J22 - J12 * J21);
213 
214  // [dV dtheta]' = -[invJ]*[Fx Fy]'
215  dV = det * ( J22 * Fx - J12 * Fy);
216  dtheta = det * (-J21 * Fx + J11 * Fy);
217  }
218 
219  V = V + dV;
220  theta = theta + dtheta;
221 
222  errV = abs(dV);
223  errTheta = abs(dtheta);
224  }
225 
226  c = sqrt(1.0 - gamma_1_2 * VV);
227  int kk = id2 + j;
228  NekDouble timeramp = 200.0;;
229  if (time<timeramp &&
230  !(m_session->DefinesFunction("InitialConditions") &&
231  m_session->GetFunctionType("InitialConditions", 0) ==
233  {
234  Fwd[0][kk] = pow(c, 1.0 / gamma_1_2) *
235  exp(-1.0 + time /timeramp);
236 
237  Fwd[1][kk] = Fwd[0][kk] * V * cos(theta) *
238  exp(-1 + time / timeramp);
239 
240  Fwd[2][kk] = Fwd[0][kk] * V * sin(theta) *
241  exp(-1 + time / timeramp);
242  }
243  else
244  {
245  Fwd[0][kk] = pow(c, 1.0 / gamma_1_2);
246  Fwd[1][kk] = Fwd[0][kk] * V * cos(theta);
247  Fwd[2][kk] = Fwd[0][kk] * V * sin(theta);
248  }
249 
250  P = (c * c) * Fwd[0][kk] / gamma;
251  Fwd[3][kk] = P / (gamma - 1.0) + 0.5 *
252  (Fwd[1][kk] * Fwd[1][kk] / Fwd[0][kk] +
253  Fwd[2][kk] * Fwd[2][kk] / Fwd[0][kk]);
254 
255  errV = 1.0;
256  errTheta = 1.0;
257  theta = M_PI / 4.0;
258  V = kExt * sin(theta);
259  }
260 
261  for (int i = 0; i < nvariables; ++i)
262  {
263  Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
264  &(m_fields[i]->GetBndCondExpansions()[m_bcRegion]->
265  UpdatePhys())[id1],1);
266  }
267  }
268 }
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: CFSBndCond.h:83
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:85
double NekDouble
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:94
int m_offset
Offset.
Definition: CFSBndCond.h:102
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:100

Friends And Related Function Documentation

◆ MemoryManager< RinglebFlowBC >

friend class MemoryManager< RinglebFlowBC >
friend

Definition at line 51 of file RinglebFlowBC.h.

Member Data Documentation

◆ className

std::string Nektar::RinglebFlowBC::className
static
Initial value:
RegisterCreatorFunction("RinglebFlow",
"Ringleb flow boundary condition.")

Name of the class.

Definition at line 67 of file RinglebFlowBC.h.

◆ m_expdim

int Nektar::RinglebFlowBC::m_expdim
private

Definition at line 84 of file RinglebFlowBC.h.

Referenced by RinglebFlowBC(), v_Apply(), and ~RinglebFlowBC().

◆ m_homo1D

bool Nektar::RinglebFlowBC::m_homo1D
private

Definition at line 87 of file RinglebFlowBC.h.

Referenced by RinglebFlowBC(), and v_Apply().