46 "Ringleb flow boundary condition.");
52 const int pSpaceDim,
const int bcRegion,
const int cnt)
53 :
CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
58 if (
m_session->DefinesSolverInfo(
"HOMOGENEOUS"))
60 std::string HomoStr =
m_session->GetSolverInfo(
"HOMOGENEOUS");
61 if ((HomoStr ==
"HOMOGENEOUS1D") || (HomoStr ==
"Homogeneous1D") ||
62 (HomoStr ==
"1D") || (HomoStr ==
"Homo1D"))
73 int nvariables = physarray.size();
79 int nPointsTot =
m_fields[0]->GetTotPoints();
80 int nPointsTot_plane =
m_fields[0]->GetPlane(0)->GetTotPoints();
81 n_planes = nPointsTot / nPointsTot_plane;
84 int id2, id2_plane, e_max;
88 for (
int e = 0; e < e_max; ++e)
100 int m_offset_plane =
m_offset / n_planes;
102 int e_max_plane = e_max / n_planes;
103 int nTracePts_plane =
m_fields[0]->GetTrace()->GetNpoints();
105 int planeID = floor((e + 0.5) / e_max_plane);
106 e_plane = e - e_max_plane * planeID;
108 id2_plane =
m_fields[0]->GetTrace()->GetPhys_Offset(
109 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
110 m_offset_plane + e_plane));
111 id2 = id2_plane + planeID * nTracePts_plane;
115 id2 =
m_fields[0]->GetTrace()->GetPhys_Offset(
116 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
128 NekDouble c, k, phi, r, J, VV, pp, sint,
P, ss;
142 NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
145 for (
int j = 0; j < npoints; j++)
148 while ((
abs(errV) > toll) || (
abs(errTheta) > toll))
152 c =
sqrt(1.0 - gamma_1_2 * VV);
156 J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
157 1.0 / (5.0 * c * c * c * c * c) -
158 0.5 *
log((1.0 + c) / (1.0 - c));
160 r = pow(c, 1.0 / gamma_1_2);
161 xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
162 yi = phi / (r * V) *
sqrt(1.0 - VV * pp);
163 par1 = 25.0 - 5.0 * VV;
170 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) * V +
171 1562.5 / pow(par1, 2.5) *
172 (-2.0 / (VV * V) + 4.0 / (VV * V) * ss) +
173 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
174 7812.5 / pow(par1, 3.5) * V -
176 (-1.0 / pow(par1, 0.5) * V /
177 (1.0 - 0.2 * pow(par1, 0.5)) -
178 (1.0 + 0.2 * pow(par1, 0.5)) /
179 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
180 pow(par1, 0.5) * V) /
181 (1.0 + 0.2 * pow(par1, 0.5)) *
182 (1.0 - 0.2 * pow(par1, 0.5));
184 J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
186 -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
187 pow((1.0 - ss), 0.5) +
188 78125.0 / V * sint / pow(par1, 3.5) * pow((1.0 - ss), 0.5);
191 if (
abs(x1[j]) < toll &&
abs(cos(theta)) < toll)
193 J22 = -39062.5 / pow(par1, 3.5) / V +
194 3125 / pow(par1, 2.5) / (VV * V) +
195 12.5 / pow(par1, 1.5) * V +
196 312.5 / pow(par1, 2.5) * V +
197 7812.5 / pow(par1, 3.5) * V -
199 (-1.0 / pow(par1, 0.5) * V /
200 (1.0 - 0.2 * pow(par1, 0.5)) -
201 (1.0 + 0.2 * pow(par1, 0.5)) /
202 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
203 pow(par1, 0.5) * V) /
204 (1.0 + 0.2 * pow(par1, 0.5)) *
205 (1.0 - 0.2 * pow(par1, 0.5));
208 dV = -1.0 / J22 * Fx;
214 J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
215 pow((1.0 - ss), 0.5) -
216 3125.0 / VV * ss / pow(par1, 2.5) /
217 pow((1.0 - ss), 0.5) * cos(theta);
219 det = -1.0 / (J11 * J22 - J12 * J21);
222 dV = det * (J22 * Fx - J12 * Fy);
223 dtheta = det * (-J21 * Fx + J11 * Fy);
227 theta = theta + dtheta;
230 errTheta =
abs(dtheta);
233 c =
sqrt(1.0 - gamma_1_2 * VV);
237 if (time < timeramp &&
238 !(
m_session->DefinesFunction(
"InitialConditions") &&
239 m_session->GetFunctionType(
"InitialConditions", 0) ==
243 pow(c, 1.0 / gamma_1_2) * exp(-1.0 + time / timeramp);
246 Fwd[0][kk] * V * cos(theta) * exp(-1 + time / timeramp);
249 Fwd[0][kk] * V * sin(theta) * exp(-1 + time / timeramp);
253 Fwd[0][kk] = pow(c, 1.0 / gamma_1_2);
254 Fwd[1][kk] = Fwd[0][kk] * V * cos(theta);
255 Fwd[2][kk] = Fwd[0][kk] * V * sin(theta);
258 P = (c * c) * Fwd[0][kk] / gamma;
259 Fwd[3][kk] =
P / (gamma - 1.0) +
260 0.5 * (Fwd[1][kk] * Fwd[1][kk] / Fwd[0][kk] +
261 Fwd[2][kk] * Fwd[2][kk] / Fwd[0][kk]);
266 V = kExt * sin(theta);
269 for (
int i = 0; i < nvariables; ++i)
274 ->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.
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.
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
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)
static std::string className
Name of the class.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ P
Monomial polynomials .
The above copyright notice and this permission notice shall be included.
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)