Nektar++
LEE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LEE.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2017 Kilian Lackhove
10// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11// Department of Aeronautics, Imperial College London (UK), and Scientific
12// Computing and Imaging Institute, University of Utah (USA).
13//
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Linearized Euler Equations
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <iostream>
37
39
40using namespace std;
41
42namespace Nektar
43{
45 "LEE", LEE::create, "Linearized Euler Equations");
46
49 : UnsteadySystem(pSession, pGraph), AcousticSystem(pSession, pGraph)
50{
51 m_ip = 0;
52 m_irho = 1;
53 m_iu = 2;
54
55 m_conservative = true;
56}
57
58/**
59 * @brief Initialization object for the LEE class.
60 */
61void LEE::v_InitObject(bool DeclareFields)
62{
63 AcousticSystem::v_InitObject(DeclareFields);
64
65 m_bfNames.push_back("gamma");
66
67 // Initialize basefield again
69 for (int i = 0; i < m_bf.size(); ++i)
70 {
72 }
73 GetFunction("Baseflow", m_fields[0], true)
74 ->Evaluate(m_bfNames, m_bf, m_time);
75
76 // Define the normal velocity fields
78 for (int i = 0; i < m_bfFwdBwd.size(); i++)
79 {
81 }
82
83 string riemName;
84 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
85 if (boost::to_lower_copy(riemName) == "characteristics" ||
86 boost::to_lower_copy(riemName) == "leeupwind" ||
87 boost::to_lower_copy(riemName) == "upwind")
88 {
89 riemName = "LEEUpwind";
90 }
91 if (boost::to_lower_copy(riemName) == "laxfriedrichs")
92 {
93 riemName = "LEELaxFriedrichs";
94 }
96 riemName, m_session);
97 m_riemannSolver->SetVector("N", &LEE::GetNormals, this);
98 m_riemannSolver->SetVector("basefieldFwdBwd", &LEE::GetBasefieldFwdBwd,
99 this);
100 m_riemannSolver->SetAuxVec("vecLocs", &LEE::GetVecLocs, this);
101
102 // Set up advection operator
103 string advName;
104 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
107 m_advection->SetFluxVector(&LEE::v_GetFluxVector, this);
108 m_advection->SetRiemannSolver(m_riemannSolver);
109 m_advection->InitObject(m_session, m_fields);
110
112 {
115 }
116 else
117 {
118 ASSERTL0(false, "Implicit LEE not set up.");
119 }
120}
121
122/**
123 * @brief Destructor for LEE class.
124 */
126{
127}
128
129/**
130 * @brief Return the flux vector for the LEE equations.
131 *
132 * @param physfield Fields.
133 * @param flux Resulting flux. flux[eq][dir][pt]
134 */
136 const Array<OneD, Array<OneD, NekDouble>> &physfield,
138{
139 int nq = physfield[0].size();
140
141 ASSERTL1(flux[0].size() == m_spacedim,
142 "Dimension of flux array and velocity array do not match");
143
146 for (int i = 0; i < m_spacedim; ++i)
147 {
148 u0[i] = m_bf[2 + i];
149 }
150
152 Array<OneD, const NekDouble> rho = physfield[m_irho];
154 for (int i = 0; i < m_spacedim; ++i)
155 {
156 ru[i] = physfield[m_iu + i];
157 }
158
159 // F_{adv,p',j} = c0^2 * ru_j + u0_j * p
160 for (int j = 0; j < m_spacedim; ++j)
161 {
162 int i = 0;
163 Vmath::Vvtvvtp(nq, c0sq, 1, ru[j], 1, u0[j], 1, p, 1, flux[i][j], 1);
164 }
165
166 // F_{adv,rho',j} = u0_j * rho' + ru_j
167 for (int j = 0; j < m_spacedim; ++j)
168 {
169 int i = 1;
170 // u0_j * rho' + ru_j
171 Vmath::Vvtvp(nq, u0[j], 1, rho, 1, ru[j], 1, flux[i][j], 1);
172 }
173
174 for (int i = 0; i < m_spacedim; ++i)
175 {
176 // F_{adv,u'_i,j} = ru_i * u0_j + delta_ij * p
177 for (int j = 0; j < m_spacedim; ++j)
178 {
179 // ru_i * u0_j
180 Vmath::Vmul(nq, ru[i], 1, u0[j], 1, flux[m_iu + i][j], 1);
181
182 // kronecker delta
183 if (i == j)
184 {
185 // delta_ij + p
186 Vmath::Vadd(nq, p, 1, flux[m_iu + i][j], 1, flux[m_iu + i][j],
187 1);
188 }
189 }
190 }
191}
192
195{
196 int nq = GetTotPoints();
197
199 for (int i = 0; i < m_spacedim + 2; ++i)
200 {
201 if (i == 1)
202 {
203 // skip rho
204 continue;
205 }
206
207 linTerm[i] = Array<OneD, NekDouble>(nq);
208 }
209
213
214 Array<OneD, NekDouble> gammaMinOne(nq);
215 Vmath::Sadd(nq, -1.0, gamma, 1, gammaMinOne, 1);
216
218 Vmath::Vmul(nq, c0sq, 1, rho0, 1, p0, 1);
219 Vmath::Vdiv(nq, p0, 1, gamma, 1, p0, 1);
220
222 for (int i = 0; i < m_spacedim; ++i)
223 {
224 u0[i] = m_bf[2 + i];
225 }
226
227 Array<OneD, const NekDouble> p = inarray[0];
228 Array<OneD, const NekDouble> rho = inarray[1];
229
231 for (int i = 0; i < m_spacedim; ++i)
232 {
233 ru[i] = inarray[2 + i];
234 }
235
236 Array<OneD, NekDouble> grad(nq);
237 Array<OneD, NekDouble> tmp1(nq);
238
239 // p
240 {
241 Vmath::Zero(nq, linTerm[m_ip], 1);
242 // (1-gamma) ( ru_j / rho0 * dp0/dx_j - p * du0_j/dx_j )
243 for (int j = 0; j < m_spacedim; ++j)
244 {
245 // ru_j / rho0 * dp0/dx_j
246 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], p0, grad);
247 Vmath::Vmul(nq, grad, 1, ru[j], 1, tmp1, 1);
248 Vmath::Vdiv(nq, tmp1, 1, rho0, 1, tmp1, 1);
249 // p * du0_j/dx_j - ru_j / rho0 * dp0/dx_j
250 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], u0[j],
251 grad);
252 Vmath::Vvtvm(nq, grad, 1, p, 1, tmp1, 1, tmp1, 1);
253 // (gamma-1) (p * du0_j/dx_j - ru_j / rho0 * dp0/dx_j)
254 Vmath::Vvtvp(nq, gammaMinOne, 1, tmp1, 1, linTerm[m_ip], 1,
255 linTerm[m_ip], 1);
256 }
257 }
258
259 // rho has no linTerm
260
261 // ru_i
262 for (int i = 0; i < m_spacedim; ++i)
263 {
264 Vmath::Zero(nq, linTerm[m_iu + i], 1);
265 // du0_i/dx_j * (u0_j * rho + ru_j)
266 for (int j = 0; j < m_spacedim; ++j)
267 {
268 // d u0_i / d x_j
269 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], u0[i],
270 grad);
271 // u0_j * rho + ru_j
272 Vmath::Vvtvp(nq, u0[j], 1, rho, 1, ru[j], 1, tmp1, 1);
273 // du0_i/dx_j * (u0_j * rho + ru_j)
274 Vmath::Vvtvp(nq, grad, 1, tmp1, 1, linTerm[m_iu + i], 1,
275 linTerm[m_iu + i], 1);
276 }
277 }
278
280 for (int i = 0; i < m_spacedim + 2; ++i)
281 {
282 if (i == 1)
283 {
284 // skip rho
285 continue;
286 }
287
288 m_fields[0]->FwdTrans(linTerm[i], tmpC);
289 m_fields[0]->BwdTrans(tmpC, linTerm[i]);
290
291 Vmath::Vadd(nq, outarray[i], 1, linTerm[i], 1, outarray[i], 1);
292 }
293}
294
295/**
296 * @brief Outflow characteristic boundary conditions for compressible
297 * flow problems.
298 */
299void LEE::v_RiemannInvariantBC(int bcRegion, int cnt,
302 Array<OneD, Array<OneD, NekDouble>> &physarray)
303{
304 int id1, id2, nBCEdgePts;
305 int nVariables = physarray.size();
306
307 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
308
309 int eMax = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
310
311 for (int e = 0; e < eMax; ++e)
312 {
313 nBCEdgePts = m_fields[0]
314 ->GetBndCondExpansions()[bcRegion]
315 ->GetExp(e)
316 ->GetTotPoints();
317 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
318 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[cnt + e]);
319
320 // Calculate (v.n)
321 Array<OneD, NekDouble> RVn(nBCEdgePts, 0.0);
322 for (int i = 0; i < m_spacedim; ++i)
323 {
324 Vmath::Vvtvp(nBCEdgePts, &Fwd[m_iu + i][id2], 1,
325 &m_traceNormals[i][id2], 1, &RVn[0], 1, &RVn[0], 1);
326 }
327
328 // Calculate (v0.n)
329 Array<OneD, NekDouble> RVn0(nBCEdgePts, 0.0);
330 for (int i = 0; i < m_spacedim; ++i)
331 {
332 Vmath::Vvtvp(nBCEdgePts, &BfFwd[2 + i][id2], 1,
333 &m_traceNormals[i][id2], 1, &RVn0[0], 1, &RVn0[0], 1);
334 }
335
336 for (int i = 0; i < nBCEdgePts; ++i)
337 {
338 NekDouble c = sqrt(BfFwd[0][id2 + i]);
339
340 NekDouble h1, h4, h5;
341
342 if (RVn0[i] > 0)
343 {
344 // rho - p / c^2
345 h1 = Fwd[m_irho][id2 + i] - Fwd[m_ip][id2 + i] / (c * c);
346 }
347 else
348 {
349 h1 = 0.0;
350 }
351
352 if (RVn0[i] - c > 0)
353 {
354 // ru / 2 - p / (2*c)
355 h4 = RVn[i] / 2 - Fwd[m_ip][id2 + i] / (2 * c);
356 }
357 else
358 {
359 h4 = 0.0;
360 }
361
362 if (RVn0[i] + c > 0)
363 {
364 // ru / 2 + p / (2*c)
365 h5 = RVn[i] / 2 + Fwd[m_ip][id2 + i] / (2 * c);
366 }
367 else
368 {
369 h5 = 0.0;
370 }
371
372 // compute conservative variables
373 // p = c0*(h5-h4)
374 // rho = h1 + (h5-h4)/c0
375 // ru = h4+h5
376 Fwd[m_ip][id2 + i] = c * (h5 - h4);
377 Fwd[m_irho][id2 + i] = h1 + (h5 - h4) / c;
378 NekDouble RVnNew = h4 + h5;
379
380 // adjust velocity pert. according to new value
381 // here we just omit the wall parallel velocity components, i.e
382 // setting them to zero. This is equivalent to setting the two
383 // vorticity characteristics h2 and h3 to zero. Mathematically,
384 // this is only legitimate for incoming characteristics. However,
385 // as h2 and h3 are convected by the flow, the value we precribe at
386 // an the boundary for putgoing characteristics does not matter.
387 // This implementation saves a few operations and is more robust
388 // for mixed in/outflow boundaries and at the boundaries edges.
389 for (int j = 0; j < m_spacedim; ++j)
390 {
391 Fwd[m_iu + j][id2 + i] = RVnNew * m_traceNormals[j][id2 + i];
392 }
393 }
394
395 // Copy boundary adjusted values into the boundary expansion
396 for (int i = 0; i < nVariables; ++i)
397 {
398 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
399 &(m_fields[i]
400 ->GetBndCondExpansions()[bcRegion]
401 ->UpdatePhys())[id1],
402 1);
403 }
404 }
405}
406
407} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
Get the locations of the components of the directed fields within the fields array.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
Get the normal vectors.
std::vector< std::string > m_bfNames
bool m_conservative
we are dealing with a conservative formualtion
SolverUtils::AdvectionSharedPtr m_advection
int m_ip
indices of the fields
Array< OneD, Array< OneD, NekDouble > > m_bfFwdBwd
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
Array< OneD, Array< OneD, NekDouble > > m_bf
void v_InitObject(bool DeclareFields=true) override
Initialization object for the AcousticSystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetBasefieldFwdBwd()
Get the baseflow field.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: LEE.h:52
static std::string className
Name of class.
Definition: LEE.h:62
void v_GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux) override
Return the flux vector for the LEE equations.
Definition: LEE.cpp:135
void v_InitObject(bool DeclareFields=true) override
Initialization object for the LEE class.
Definition: LEE.cpp:61
~LEE() override
Destructor.
Definition: LEE.cpp:125
LEE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
Definition: LEE.cpp:47
void v_AddLinTerm(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray) override
Definition: LEE.cpp:193
void v_RiemannInvariantBC(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &BfFwd, Array< OneD, Array< OneD, NekDouble > > &physarray) override
Outflow characteristic boundary conditions for compressible flow problems.
Definition: LEE.cpp:299
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNcoeffs()
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT int GetTotPoints()
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
double NekDouble
void Vmul(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.hpp:72
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.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Vvtvm(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)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381
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.hpp:126
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.hpp:439
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294