Nektar++
RinglebFlowBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RinglebFlowBC.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Ringleb flow boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "RinglebFlowBC.h"
37
38using namespace std;
39
40namespace Nektar
41{
42
43std::string RinglebFlowBC::className =
45 "RinglebFlow", RinglebFlowBC::create,
46 "Ringleb flow boundary condition.");
47
51 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
52 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
53 const int pSpaceDim, const int bcRegion, const int cnt)
54 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
55 bcRegion, cnt)
56{
57 m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
58
59 m_homo1D = false;
60 if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
61 {
62 std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
63 if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D") ||
64 (HomoStr == "1D") || (HomoStr == "Homo1D"))
65 {
66 m_homo1D = true;
67 }
68 }
69}
70
73 const NekDouble &time)
74{
75 int nvariables = physarray.size();
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]
94 ->GetExp(e)
95 ->GetTotPoints();
96 int id1 =
97 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
98
99 // For 3DHomogenoeus1D
100 if (m_expdim == 2 && m_homo1D)
101 {
102 int m_offset_plane = m_offset / n_planes;
103 int e_plane;
104 int e_max_plane = e_max / n_planes;
105 int nTracePts_plane = m_fields[0]->GetTrace()->GetNpoints();
106
107 int planeID = floor((e + 0.5) / e_max_plane);
108 e_plane = e - e_max_plane * planeID;
109
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;
114 }
115 else // For general case
116 {
117 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
118 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(
119 m_offset + e));
120 }
121
122 Array<OneD, NekDouble> x0(npoints, 0.0);
123 Array<OneD, NekDouble> x1(npoints, 0.0);
124 Array<OneD, NekDouble> x2(npoints, 0.0);
125
126 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
127 x0, x1, x2);
128
129 // Flow parameters
130 NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
131 NekDouble J11, J12, J21, J22, det;
132 NekDouble Fx, Fy;
133 NekDouble xi, yi;
134 NekDouble dV;
135 NekDouble dtheta;
136 NekDouble par1;
137 NekDouble theta = M_PI / 4.0;
138 NekDouble kExt = 0.7;
139 NekDouble V = kExt * sin(theta);
140 NekDouble toll = 1.0e-8;
141 NekDouble errV = 1.0;
142 NekDouble errTheta = 1.0;
143 NekDouble gamma = m_gamma;
144 NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
145
146 // Loop on all the points of that edge
147 for (int j = 0; j < npoints; j++)
148 {
149
150 while ((abs(errV) > toll) || (abs(errTheta) > toll))
151 {
152 VV = V * V;
153 sint = sin(theta);
154 c = sqrt(1.0 - gamma_1_2 * VV);
155 k = V / sint;
156 phi = 1.0 / k;
157 pp = phi * phi;
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));
161
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;
166 ss = sint * sint;
167
168 Fx = xi - x0[j];
169 Fy = yi - x1[j];
170
171 J11 =
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 -
177 0.25 *
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));
185
186 J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
187 J21 =
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);
191
192 // the matrix is singular when theta = pi/2
193 if (abs(x1[j]) < toll && abs(cos(theta)) < toll)
194 {
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 -
200 0.25 *
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));
208
209 // dV = -dV/dx * Fx
210 dV = -1.0 / J22 * Fx;
211 dtheta = 0.0;
212 theta = M_PI / 2.0;
213 }
214 else
215 {
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);
220
221 det = -1.0 / (J11 * J22 - J12 * J21);
222
223 // [dV dtheta]' = -[invJ]*[Fx Fy]'
224 dV = det * (J22 * Fx - J12 * Fy);
225 dtheta = det * (-J21 * Fx + J11 * Fy);
226 }
227
228 V = V + dV;
229 theta = theta + dtheta;
230
231 errV = abs(dV);
232 errTheta = abs(dtheta);
233 }
234
235 c = sqrt(1.0 - gamma_1_2 * VV);
236 int kk = id2 + j;
237 NekDouble timeramp = 200.0;
238 ;
239 if (time < timeramp &&
240 !(m_session->DefinesFunction("InitialConditions") &&
241 m_session->GetFunctionType("InitialConditions", 0) ==
243 {
244 Fwd[0][kk] =
245 pow(c, 1.0 / gamma_1_2) * exp(-1.0 + time / timeramp);
246
247 Fwd[1][kk] =
248 Fwd[0][kk] * V * cos(theta) * exp(-1 + time / timeramp);
249
250 Fwd[2][kk] =
251 Fwd[0][kk] * V * sin(theta) * exp(-1 + time / timeramp);
252 }
253 else
254 {
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);
258 }
259
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]);
264
265 errV = 1.0;
266 errTheta = 1.0;
267 theta = M_PI / 4.0;
268 V = kExt * sin(theta);
269 }
270
271 for (int i = 0; i < nvariables; ++i)
272 {
273 Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
274 &(m_fields[i]
275 ->GetBndCondExpansions()[m_bcRegion]
276 ->UpdatePhys())[id1],
277 1);
278 }
279 }
280}
281
282} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: CFSBndCond.h:90
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:105
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
int m_offset
Offset.
Definition: CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
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.
Definition: RinglebFlowBC.h:52
static std::string className
Name of the class.
Definition: RinglebFlowBC.h:66
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ P
Monomial polynomials .
Definition: BasisType.h:62
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294