Nektar++
LinearSWE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LinearSWE.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: Linear Shallow water equations in primitive variables
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/algorithm/string.hpp>
36#include <iomanip>
37#include <iostream>
38
41
42using namespace std;
43
44namespace Nektar
45{
48 "LinearSWE", LinearSWE::create,
49 "Linear shallow water equation in primitive variables.");
50
53 : ShallowWaterSystem(pSession, pGraph)
54{
55}
56
57void LinearSWE::v_InitObject(bool DeclareFields)
58{
60
62 {
65 }
66 else
67 {
68 ASSERTL0(false, "Implicit SWE not set up.");
69 }
70
71 // Type of advection class to be used
72 switch (m_projectionType)
73 {
74 // Continuous field
76 {
77 // Do nothing
78 break;
79 }
80 // Discontinuous field
82 {
83 string advName;
84 string diffName;
85 string riemName;
86
87 //---------------------------------------------------------------
88 // Setting up advection and diffusion operators
89 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
91 advName, advName);
92 m_advection->SetFluxVector(&LinearSWE::GetFluxVector, this);
93
94 // Setting up Riemann solver for advection operator
95 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
96 if ((riemName == "LinearHLL") && (m_constantDepth != true))
97 {
98 ASSERTL0(false, "LinearHLL only valid for constant depth");
99 }
102 riemName, m_session);
103
104 // Setting up parameters for advection operator Riemann solver
105 m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
106 m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
107 m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
108
109 // The numerical flux for linear SWE requires depth information
110 int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
111 m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
112 m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
113 m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
115 m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
116 this);
117 m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
118 this);
119
120 // Concluding initialisation of advection operators
121 m_advection->SetRiemannSolver(m_riemannSolver);
122 m_advection->InitObject(m_session, m_fields);
123 break;
124 }
125 default:
126 {
127 ASSERTL0(false, "Unsupported projection type.");
128 break;
129 }
130 }
131}
132
134{
135}
136
137// physarray contains the conservative variables
139 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
141{
142
143 int ncoeffs = GetNcoeffs();
144 int nq = GetTotPoints();
145
147 Array<OneD, NekDouble> mod(ncoeffs);
148
149 switch (m_projectionType)
150 {
152 {
153 // add to u equation
154 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
155 m_fields[0]->IProductWRTBase(tmp, mod);
156 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
157 m_fields[0]->BwdTrans(mod, tmp);
158 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
159
160 // add to v equation
161 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
162 Vmath::Neg(nq, tmp, 1);
163 m_fields[0]->IProductWRTBase(tmp, mod);
164 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
165 m_fields[0]->BwdTrans(mod, tmp);
166 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
167 }
168 break;
171 {
172 // add to u equation
173 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
174 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
175
176 // add to v equation
177 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
178 Vmath::Neg(nq, tmp, 1);
179 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
180 }
181 break;
182 default:
183 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
184 break;
185 }
186}
187
189 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
190 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
191{
192 int i, j;
193 int ndim = m_spacedim;
194 int nvariables = inarray.size();
195 int nq = GetTotPoints();
196
197 switch (m_projectionType)
198 {
200 {
201
202 //-------------------------------------------------------
203 // Compute the DG advection including the numerical flux
204 // by using SolverUtils/Advection
205 // Input and output in physical space
207
208 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
209 time);
210 //-------------------------------------------------------
211
212 //-------------------------------------------------------
213 // negate the outarray since moving terms to the rhs
214 for (i = 0; i < nvariables; ++i)
215 {
216 Vmath::Neg(nq, outarray[i], 1);
217 }
218 //-------------------------------------------------------
219
220 //-------------------------------------------------
221 // Add "source terms"
222 // Input and output in physical space
223
224 // Coriolis forcing
225 if (m_coriolis.size() != 0)
226 {
227 AddCoriolis(inarray, outarray);
228 }
229 //-------------------------------------------------
230 }
231 break;
234 {
235
236 //-------------------------------------------------------
237 // Compute the fluxvector in physical space
239 nvariables);
240
241 for (i = 0; i < nvariables; ++i)
242 {
243 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
244 for (j = 0; j < ndim; ++j)
245 {
246 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
247 }
248 }
249
250 LinearSWE::GetFluxVector(inarray, fluxvector);
251 //-------------------------------------------------------
252
253 //-------------------------------------------------------
254 // Take the derivative of the flux terms
255 // and negate the outarray since moving terms to the rhs
257 Array<OneD, NekDouble> tmp1(nq);
258
259 for (i = 0; i < nvariables; ++i)
260 {
262 fluxvector[i][0], tmp);
264 fluxvector[i][1], tmp1);
265 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
266 Vmath::Neg(nq, outarray[i], 1);
267 }
268
269 //-------------------------------------------------
270 // Add "source terms"
271 // Input and output in physical space
272
273 // Coriolis forcing
274 if (m_coriolis.size() != 0)
275 {
276 AddCoriolis(inarray, outarray);
277 }
278 //-------------------------------------------------
279 }
280 break;
281 default:
282 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
283 break;
284 }
285}
286
288 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
289 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
290{
291 int i;
292 int nvariables = inarray.size();
293
294 switch (m_projectionType)
295 {
297 {
298
299 // Just copy over array
300 if (inarray != outarray)
301 {
302 int npoints = GetNpoints();
303
304 for (i = 0; i < nvariables; ++i)
305 {
306 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
307 }
308 }
309 SetBoundaryConditions(outarray, time);
310 break;
311 }
314 {
315
317 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
318
319 for (i = 0; i < nvariables; ++i)
320 {
321 m_fields[i]->FwdTrans(inarray[i], coeffs);
322 m_fields[i]->BwdTrans(coeffs, outarray[i]);
323 }
324 break;
325 }
326 default:
327 ASSERTL0(false, "Unknown projection scheme");
328 break;
329 }
330}
331
332//----------------------------------------------------
335{
336 std::string varName;
337 int nvariables = m_fields.size();
338 int cnt = 0;
339 int nTracePts = GetTraceTotPoints();
340
341 // Extract trace for boundaries. Needs to be done on all processors to avoid
342 // deadlock.
343 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
344 for (int i = 0; i < nvariables; ++i)
345 {
346 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
347 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
348 }
349
350 // loop over Boundary Regions
351 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
352 {
353 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
355 {
356 continue;
357 }
358
359 // Wall Boundary Condition
360 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
361 "Wall"))
362 {
363 WallBoundary2D(n, cnt, Fwd, inarray);
364 }
365
366 // Time Dependent Boundary Condition (specified in meshfile)
367 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
368 {
369 for (int i = 0; i < nvariables; ++i)
370 {
371 varName = m_session->GetVariable(i);
372 m_fields[i]->EvaluateBoundaryConditions(time, varName);
373 }
374 }
375 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
376 }
377}
378
379//----------------------------------------------------
380/**
381 * @brief Wall boundary condition.
382 */
383void LinearSWE::WallBoundary(int bcRegion, int cnt,
385 Array<OneD, Array<OneD, NekDouble>> &physarray)
386{
387 int i;
388 int nvariables = physarray.size();
389
390 // Adjust the physical values of the trace to take
391 // user defined boundaries into account
392 int e, id1, id2, npts;
393
394 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
395 ++e)
396 {
397 npts = m_fields[0]
398 ->GetBndCondExpansions()[bcRegion]
399 ->GetExp(e)
400 ->GetTotPoints();
401 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
402 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
403 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
404
405 // For 2D/3D, define: v* = v - 2(v.n)n
406 Array<OneD, NekDouble> tmp(npts, 0.0);
407
408 // Calculate (v.n)
409 for (i = 0; i < m_spacedim; ++i)
410 {
411 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
412 &tmp[0], 1, &tmp[0], 1);
413 }
414
415 // Calculate 2.0(v.n)
416 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
417
418 // Calculate v* = v - 2.0(v.n)n
419 for (i = 0; i < m_spacedim; ++i)
420 {
421 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
422 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
423 }
424
425 // copy boundary adjusted values into the boundary expansion
426 for (i = 0; i < nvariables; ++i)
427 {
428 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
429 &(m_fields[i]
430 ->GetBndCondExpansions()[bcRegion]
431 ->UpdatePhys())[id1],
432 1);
433 }
434 }
435}
436
437void LinearSWE::WallBoundary2D(int bcRegion, int cnt,
439 Array<OneD, Array<OneD, NekDouble>> &physarray)
440{
441
442 int i;
443 int nvariables = physarray.size();
444
445 // Adjust the physical values of the trace to take
446 // user defined boundaries into account
447 int e, id1, id2, npts;
448
449 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
450 ++e)
451 {
452 npts = m_fields[0]
453 ->GetBndCondExpansions()[bcRegion]
454 ->GetExp(e)
455 ->GetNumPoints(0);
456 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
457 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
458 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
459
460 switch (m_expdim)
461 {
462 case 1:
463 {
464 // negate the forward flux
465 Vmath::Neg(npts, &Fwd[1][id2], 1);
466 }
467 break;
468 case 2:
469 {
470 Array<OneD, NekDouble> tmp_n(npts);
471 Array<OneD, NekDouble> tmp_t(npts);
472
473 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
474 &tmp_n[0], 1);
475 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
476 &tmp_n[0], 1, &tmp_n[0], 1);
477
478 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
479 &tmp_t[0], 1);
480 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
481 &tmp_t[0], 1, &tmp_t[0], 1);
482
483 // negate the normal flux
484 Vmath::Neg(npts, tmp_n, 1);
485
486 // rotate back to Cartesian
487 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
488 &Fwd[1][id2], 1);
489 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
490 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
491
492 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
493 &Fwd[2][id2], 1);
494 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
495 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
496 }
497 break;
498 case 3:
499 ASSERTL0(false,
500 "3D not implemented for Shallow Water Equations");
501 break;
502 default:
503 ASSERTL0(false, "Illegal expansion dimension");
504 }
505
506 // copy boundary adjusted values into the boundary expansion
507 for (i = 0; i < nvariables; ++i)
508 {
509 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
510 &(m_fields[i]
511 ->GetBndCondExpansions()[bcRegion]
512 ->UpdatePhys())[id1],
513 1);
514 }
515 }
516}
517
518// Physfield in primitive Form
520 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
522{
523 int i, j;
524 int nq = m_fields[0]->GetTotPoints();
525
526 NekDouble g = m_g;
527
528 // Flux vector for the mass equation
529 for (i = 0; i < m_spacedim; ++i)
530 {
531 Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
532 }
533
534 // Put (g eta) in tmp
536 Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
537
538 // Flux vector for the momentum equations
539 for (i = 0; i < m_spacedim; ++i)
540 {
541 for (j = 0; j < m_spacedim; ++j)
542 {
543 // must zero fluxes as not initialised to zero in AdvectionWeakDG
544 // ...
545 Vmath::Zero(nq, flux[i + 1][j], 1);
546 }
547
548 // Add (g eta) to appropriate field
549 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
550 }
551}
552
554 const Array<OneD, const Array<OneD, NekDouble>> &physin,
556{
557 int nq = GetTotPoints();
558
559 if (physin.get() == physout.get())
560 {
561 // copy indata and work with tmp array
563 for (int i = 0; i < 3; ++i)
564 {
565 // deep copy
566 tmp[i] = Array<OneD, NekDouble>(nq);
567 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
568 }
569
570 // \eta = h - d
571 Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
572
573 // u = hu/h
574 Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
575
576 // v = hv/ v
577 Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
578 }
579 else
580 {
581 // \eta = h - d
582 Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
583
584 // u = hu/h
585 Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
586
587 // v = hv/ v
588 Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
589 }
590}
591
593{
594 int nq = GetTotPoints();
595
596 // u = hu/h
597 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
598 m_fields[1]->UpdatePhys(), 1);
599
600 // v = hv/ v
601 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
602 m_fields[2]->UpdatePhys(), 1);
603
604 // \eta = h - d
605 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
606 m_fields[0]->UpdatePhys(), 1);
607}
608
610 const Array<OneD, const Array<OneD, NekDouble>> &physin,
612{
613
614 int nq = GetTotPoints();
615
616 if (physin.get() == physout.get())
617 {
618 // copy indata and work with tmp array
620 for (int i = 0; i < 3; ++i)
621 {
622 // deep copy
623 tmp[i] = Array<OneD, NekDouble>(nq);
624 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
625 }
626
627 // h = \eta + d
628 Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
629
630 // hu = h * u
631 Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
632
633 // hv = h * v
634 Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
635 }
636 else
637 {
638 // h = \eta + d
639 Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
640
641 // hu = h * u
642 Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
643
644 // hv = h * v
645 Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
646 }
647}
648
650{
651 int nq = GetTotPoints();
652
653 // h = \eta + d
654 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
655 m_fields[0]->UpdatePhys(), 1);
656
657 // hu = h * u
658 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
659 m_fields[1]->UpdatePhys(), 1);
660
661 // hv = h * v
662 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
663 m_fields[2]->UpdatePhys(), 1);
664}
665
666/**
667 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
668 * \f$ h\mathbf{v} \f$.
669 *
670 * @param physfield Velocity field.
671 * @param velocity Velocity field.
672 */
674 const Array<OneD, Array<OneD, NekDouble>> &physfield,
676{
677 const int npts = physfield[0].size();
678
679 for (int i = 0; i < m_spacedim; ++i)
680 {
681 Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
682 }
683}
684
686{
688 if (m_session->DefinesSolverInfo("UpwindType"))
689 {
690 std::string UpwindType;
691 UpwindType = m_session->GetSolverInfo("UpwindType");
692 if (UpwindType == "LinearAverage")
693 {
694 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
695 }
696 if (UpwindType == "LinearHLL")
697 {
698 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
699 }
700 }
701 SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
702 SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
703 SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
704}
705
706} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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)
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Definition: LinearSWE.cpp:685
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
Definition: LinearSWE.cpp:673
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:188
LinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: LinearSWE.cpp:51
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: LinearSWE.h:54
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:100
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:287
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:96
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
Definition: LinearSWE.cpp:333
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Definition: LinearSWE.cpp:437
void v_ConservativeToPrimitive() override
Definition: LinearSWE.cpp:592
~LinearSWE() override
Definition: LinearSWE.cpp:133
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
Definition: LinearSWE.cpp:57
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:519
static std::string className
Name of class.
Definition: LinearSWE.h:64
void v_PrimitiveToConservative() override
Definition: LinearSWE.cpp:649
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:75
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:76
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Wall boundary condition.
Definition: LinearSWE.cpp:383
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: LinearSWE.cpp:138
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, NekDouble > m_depth
Still water depth.
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
int m_spacedim
Spatial dimension (>= expansion dim).
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
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 GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetTotPoints()
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
const std::vector< NekDouble > velocity
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
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 Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220