46 "Ringleb flow boundary condition.");
53 const int pSpaceDim,
const int bcRegion,
const int cnt)
54 :
CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
57 m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
60 if (
m_session->DefinesSolverInfo(
"HOMOGENEOUS"))
62 std::string HomoStr =
m_session->GetSolverInfo(
"HOMOGENEOUS");
63 if ((HomoStr ==
"HOMOGENEOUS1D") || (HomoStr ==
"Homogeneous1D") ||
64 (HomoStr ==
"1D") || (HomoStr ==
"Homo1D"))
75 int nvariables = physarray.size();
81 int nPointsTot =
m_fields[0]->GetTotPoints();
82 int nPointsTot_plane =
m_fields[0]->GetPlane(0)->GetTotPoints();
83 n_planes = nPointsTot / nPointsTot_plane;
86 int id2, id2_plane, e_max;
90 for (
int e = 0; e < e_max; ++e)
102 int m_offset_plane =
m_offset / n_planes;
104 int e_max_plane = e_max / n_planes;
105 int nTracePts_plane =
m_fields[0]->GetTrace()->GetNpoints();
107 int planeID = floor((e + 0.5) / e_max_plane);
108 e_plane = e - e_max_plane * planeID;
110 id2_plane =
m_fields[0]->GetTrace()->GetPhys_Offset(
111 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
112 m_offset_plane + e_plane));
113 id2 = id2_plane + planeID * nTracePts_plane;
117 id2 =
m_fields[0]->GetTrace()->GetPhys_Offset(
118 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
130 NekDouble c, k, phi, r, J, VV, pp, sint,
P, ss;
144 NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
147 for (
int j = 0; j < npoints; j++)
150 while ((
abs(errV) > toll) || (
abs(errTheta) > toll))
154 c =
sqrt(1.0 - gamma_1_2 * VV);
158 J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
159 1.0 / (5.0 * c * c * c * c * c) -
160 0.5 *
log((1.0 + c) / (1.0 - c));
162 r = pow(c, 1.0 / gamma_1_2);
163 xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
164 yi = phi / (r * V) *
sqrt(1.0 - VV * pp);
165 par1 = 25.0 - 5.0 * VV;
172 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) * V +
173 1562.5 / pow(par1, 2.5) *
174 (-2.0 / (VV * V) + 4.0 / (VV * V) * ss) +
175 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
176 7812.5 / pow(par1, 3.5) * V -
178 (-1.0 / pow(par1, 0.5) * V /
179 (1.0 - 0.2 * pow(par1, 0.5)) -
180 (1.0 + 0.2 * pow(par1, 0.5)) /
181 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
182 pow(par1, 0.5) * V) /
183 (1.0 + 0.2 * pow(par1, 0.5)) *
184 (1.0 - 0.2 * pow(par1, 0.5));
186 J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
188 -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
189 pow((1.0 - ss), 0.5) +
190 78125.0 / V * sint / pow(par1, 3.5) * pow((1.0 - ss), 0.5);
193 if (
abs(x1[j]) < toll &&
abs(cos(theta)) < toll)
195 J22 = -39062.5 / pow(par1, 3.5) / V +
196 3125 / pow(par1, 2.5) / (VV * V) +
197 12.5 / pow(par1, 1.5) * V +
198 312.5 / pow(par1, 2.5) * V +
199 7812.5 / pow(par1, 3.5) * V -
201 (-1.0 / pow(par1, 0.5) * V /
202 (1.0 - 0.2 * pow(par1, 0.5)) -
203 (1.0 + 0.2 * pow(par1, 0.5)) /
204 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
205 pow(par1, 0.5) * V) /
206 (1.0 + 0.2 * pow(par1, 0.5)) *
207 (1.0 - 0.2 * pow(par1, 0.5));
210 dV = -1.0 / J22 * Fx;
216 J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
217 pow((1.0 - ss), 0.5) -
218 3125.0 / VV * ss / pow(par1, 2.5) /
219 pow((1.0 - ss), 0.5) * cos(theta);
221 det = -1.0 / (J11 * J22 - J12 * J21);
224 dV = det * (J22 * Fx - J12 * Fy);
225 dtheta = det * (-J21 * Fx + J11 * Fy);
229 theta = theta + dtheta;
232 errTheta =
abs(dtheta);
235 c =
sqrt(1.0 - gamma_1_2 * VV);
239 if (time < timeramp &&
240 !(
m_session->DefinesFunction(
"InitialConditions") &&
241 m_session->GetFunctionType(
"InitialConditions", 0) ==
245 pow(c, 1.0 / gamma_1_2) * exp(-1.0 + time / timeramp);
248 Fwd[0][kk] * V * cos(theta) * exp(-1 + time / timeramp);
251 Fwd[0][kk] * V * sin(theta) * exp(-1 + time / timeramp);
255 Fwd[0][kk] = pow(c, 1.0 / gamma_1_2);
256 Fwd[1][kk] = Fwd[0][kk] * V * cos(theta);
257 Fwd[2][kk] = Fwd[0][kk] * V * sin(theta);
260 P = (c * c) * Fwd[0][kk] / gamma;
261 Fwd[3][kk] =
P / (gamma - 1.0) +
262 0.5 * (Fwd[1][kk] * Fwd[1][kk] / Fwd[0][kk] +
263 Fwd[2][kk] * Fwd[2][kk] / Fwd[0][kk]);
268 V = kExt * sin(theta);
271 for (
int i = 0; i < nvariables; ++i)
276 ->UpdatePhys())[id1],
Encapsulates the user-defined boundary conditions for compressible flow solver.
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
NekDouble m_gamma
Parameters of the flow.
int m_bcRegion
Id of the boundary region.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
RinglebFlowBC(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Creates an instance of this class.
static std::string className
Name of the class.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ P
Monomial polynomials .
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
scalarT< T > abs(scalarT< T > in)
scalarT< T > log(scalarT< T > in)
scalarT< T > sqrt(scalarT< T > in)