Nektar++
NonlinearSWE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NonlinearSWE.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: Nonlinear Shallow water equations in conservative 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 "NonlinearSWE", NonlinearSWE::create,
49 "Nonlinear shallow water equation in conservative variables.");
50
53 : ShallowWaterSystem(pSession, pGraph)
54{
55}
56
57void NonlinearSWE::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(&NonlinearSWE::GetFluxVector, this);
93
94 // Setting up Riemann solver for advection operator
95 m_session->LoadSolverInfo("UpwindType", riemName, "Average");
98 riemName, m_session);
99
100 // Setting up parameters for advection operator Riemann solver
101 m_riemannSolver->SetParam("gravity", &NonlinearSWE::GetGravity,
102 this);
103 m_riemannSolver->SetAuxVec("vecLocs", &NonlinearSWE::GetVecLocs,
104 this);
105 m_riemannSolver->SetVector("N", &NonlinearSWE::GetNormals, this);
106 m_riemannSolver->SetScalar("depth", &NonlinearSWE::GetDepth, this);
107
108 // Concluding initialisation of advection operators
109 m_advection->SetRiemannSolver(m_riemannSolver);
110 m_advection->InitObject(m_session, m_fields);
111 break;
112 }
113 default:
114 {
115 ASSERTL0(false, "Unsupported projection type.");
116 break;
117 }
118 }
119}
120
122{
123}
124
125// physarray contains the conservative variables
127 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
129{
130
131 int ncoeffs = GetNcoeffs();
132 int nq = GetTotPoints();
133
135 Array<OneD, NekDouble> mod(ncoeffs);
136
137 switch (m_projectionType)
138 {
140 {
141 // add to hu equation
142 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
143 m_fields[0]->IProductWRTBase(tmp, mod);
144 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
145 m_fields[0]->BwdTrans(mod, tmp);
146 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
147
148 // add to hv equation
149 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
150 Vmath::Neg(nq, tmp, 1);
151 m_fields[0]->IProductWRTBase(tmp, mod);
152 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
153 m_fields[0]->BwdTrans(mod, tmp);
154 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
155 }
156 break;
159 {
160 // add to hu equation
161 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
162 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
163
164 // add to hv equation
165 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
166 Vmath::Neg(nq, tmp, 1);
167 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
168 }
169 break;
170 default:
171 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
172 break;
173 }
174}
175
176// physarray contains the conservative variables
178 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
180{
181
182 int ncoeffs = GetNcoeffs();
183 int nq = GetTotPoints();
184
186 Array<OneD, NekDouble> mod(ncoeffs);
187
188 switch (m_projectionType)
189 {
191 {
192 for (int i = 0; i < m_spacedim; ++i)
193 {
194 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
195 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
196 m_fields[0]->IProductWRTBase(tmp, mod);
197 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
198 m_fields[0]->BwdTrans(mod, tmp);
199 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
200 }
201 }
202 break;
205 {
206 for (int i = 0; i < m_spacedim; ++i)
207 {
208 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
209 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
210 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
211 }
212 }
213 break;
214 default:
215 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
216 break;
217 }
218}
219
221 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
222 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
223{
224 int i, j;
225 int ndim = m_spacedim;
226 int nvariables = inarray.size();
227 int nq = GetTotPoints();
228
229 switch (m_projectionType)
230 {
232 {
233
234 //-------------------------------------------------------
235 // Compute the DG advection including the numerical flux
236 // by using SolverUtils/Advection
237 // Input and output in physical space
239
240 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
241 time);
242 //-------------------------------------------------------
243
244 //-------------------------------------------------------
245 // negate the outarray since moving terms to the rhs
246 for (i = 0; i < nvariables; ++i)
247 {
248 Vmath::Neg(nq, outarray[i], 1);
249 }
250 //-------------------------------------------------------
251
252 //-------------------------------------------------
253 // Add "source terms"
254 // Input and output in physical space
255
256 // Coriolis forcing
257 if (m_coriolis.size() != 0)
258 {
259 AddCoriolis(inarray, outarray);
260 }
261
262 // Variable Depth
263 if (m_constantDepth != true)
264 {
265 AddVariableDepth(inarray, outarray);
266 }
267 //-------------------------------------------------
268 }
269 break;
272 {
273
274 //-------------------------------------------------------
275 // Compute the fluxvector in physical space
277 nvariables);
278
279 for (i = 0; i < nvariables; ++i)
280 {
281 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
282 for (j = 0; j < ndim; ++j)
283 {
284 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
285 }
286 }
287
288 NonlinearSWE::GetFluxVector(inarray, fluxvector);
289 //-------------------------------------------------------
290
291 //-------------------------------------------------------
292 // Take the derivative of the flux terms
293 // and negate the outarray since moving terms to the rhs
295 Array<OneD, NekDouble> tmp1(nq);
296
297 for (i = 0; i < nvariables; ++i)
298 {
300 fluxvector[i][0], tmp);
302 fluxvector[i][1], tmp1);
303 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
304 Vmath::Neg(nq, outarray[i], 1);
305 }
306
307 //-------------------------------------------------
308 // Add "source terms"
309 // Input and output in physical space
310
311 // Coriolis forcing
312 if (m_coriolis.size() != 0)
313 {
314 AddCoriolis(inarray, outarray);
315 }
316
317 // Variable Depth
318 if (m_constantDepth != true)
319 {
320 AddVariableDepth(inarray, outarray);
321 }
322 //-------------------------------------------------
323 }
324 break;
325 default:
326 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
327 break;
328 }
329}
330
332 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
333 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
334{
335 int i;
336 int nvariables = inarray.size();
337
338 switch (m_projectionType)
339 {
341 {
342
343 // Just copy over array
344 if (inarray != outarray)
345 {
346 int npoints = GetNpoints();
347
348 for (i = 0; i < nvariables; ++i)
349 {
350 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
351 }
352 }
353 SetBoundaryConditions(outarray, time);
354 break;
355 }
358 {
359
361 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
362
363 for (i = 0; i < nvariables; ++i)
364 {
365 m_fields[i]->FwdTrans(inarray[i], coeffs);
366 m_fields[i]->BwdTrans(coeffs, outarray[i]);
367 }
368 break;
369 }
370 default:
371 ASSERTL0(false, "Unknown projection scheme");
372 break;
373 }
374}
375
376//----------------------------------------------------
379{
380 std::string varName;
381 int nvariables = m_fields.size();
382 int cnt = 0;
383 int nTracePts = GetTraceTotPoints();
384
385 // Extract trace for boundaries. Needs to be done on all processors to avoid
386 // deadlock.
387 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
388 for (int i = 0; i < nvariables; ++i)
389 {
390 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
391 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
392 }
393
394 // Loop over Boundary Regions
395 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
396 {
397
398 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
400 {
401 continue;
402 }
403
404 // Wall Boundary Condition
405 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
406 "Wall"))
407 {
408 WallBoundary2D(n, cnt, Fwd, inarray);
409 }
410
411 // Time Dependent Boundary Condition (specified in meshfile)
412 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
413 {
414 for (int i = 0; i < nvariables; ++i)
415 {
416 varName = m_session->GetVariable(i);
417 m_fields[i]->EvaluateBoundaryConditions(time, varName);
418 }
419 }
420 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
421 }
422}
423
424//----------------------------------------------------
425/**
426 * @brief Wall boundary condition.
427 */
428void NonlinearSWE::WallBoundary(int bcRegion, int cnt,
430 Array<OneD, Array<OneD, NekDouble>> &physarray)
431{
432 int i;
433 int nvariables = physarray.size();
434
435 // Adjust the physical values of the trace to take
436 // user defined boundaries into account
437 int e, id1, id2, npts;
438
439 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
440 ++e)
441 {
442 npts = m_fields[0]
443 ->GetBndCondExpansions()[bcRegion]
444 ->GetExp(e)
445 ->GetTotPoints();
446 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
447 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
448 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
449
450 // For 2D/3D, define: v* = v - 2(v.n)n
451 Array<OneD, NekDouble> tmp(npts, 0.0);
452
453 // Calculate (v.n)
454 for (i = 0; i < m_spacedim; ++i)
455 {
456 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
457 &tmp[0], 1, &tmp[0], 1);
458 }
459
460 // Calculate 2.0(v.n)
461 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
462
463 // Calculate v* = v - 2.0(v.n)n
464 for (i = 0; i < m_spacedim; ++i)
465 {
466 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
467 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
468 }
469
470 // copy boundary adjusted values into the boundary expansion
471 for (i = 0; i < nvariables; ++i)
472 {
473 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
474 &(m_fields[i]
475 ->GetBndCondExpansions()[bcRegion]
476 ->UpdatePhys())[id1],
477 1);
478 }
479 }
480}
481
483 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd,
484 Array<OneD, Array<OneD, NekDouble>> &physarray)
485{
486
487 int i;
488 int nvariables = physarray.size();
489
490 // Adjust the physical values of the trace to take
491 // user defined boundaries into account
492 int e, id1, id2, npts;
493
494 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
495 ++e)
496 {
497 npts = m_fields[0]
498 ->GetBndCondExpansions()[bcRegion]
499 ->GetExp(e)
500 ->GetNumPoints(0);
501 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
502 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
503 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
504
505 switch (m_expdim)
506 {
507 case 1:
508 {
509 // negate the forward flux
510 Vmath::Neg(npts, &Fwd[1][id2], 1);
511 }
512 break;
513 case 2:
514 {
515 Array<OneD, NekDouble> tmp_n(npts);
516 Array<OneD, NekDouble> tmp_t(npts);
517
518 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
519 &tmp_n[0], 1);
520 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
521 &tmp_n[0], 1, &tmp_n[0], 1);
522
523 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
524 &tmp_t[0], 1);
525 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
526 &tmp_t[0], 1, &tmp_t[0], 1);
527
528 // negate the normal flux
529 Vmath::Neg(npts, tmp_n, 1);
530
531 // rotate back to Cartesian
532 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
533 &Fwd[1][id2], 1);
534 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
535 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
536
537 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
538 &Fwd[2][id2], 1);
539 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
540 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
541 }
542 break;
543 case 3:
544 ASSERTL0(false,
545 "3D not implemented for Shallow Water Equations");
546 break;
547 default:
548 ASSERTL0(false, "Illegal expansion dimension");
549 }
550
551 // copy boundary adjusted values into the boundary expansion
552 for (i = 0; i < nvariables; ++i)
553 {
554 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
555 &(m_fields[i]
556 ->GetBndCondExpansions()[bcRegion]
557 ->UpdatePhys())[id1],
558 1);
559 }
560 }
561}
562
563// Physfield in conservative Form
565 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
567{
568 int i, j;
569 int nq = m_fields[0]->GetTotPoints();
570
571 NekDouble g = m_g;
573
574 // Flux vector for the mass equation
575 for (i = 0; i < m_spacedim; ++i)
576 {
578 Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
579 }
580
581 GetVelocityVector(physfield, velocity);
582
583 // Put (0.5 g h h) in tmp
585 Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
586 Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
587
588 // Flux vector for the momentum equations
589 for (i = 0; i < m_spacedim; ++i)
590 {
591 for (j = 0; j < m_spacedim; ++j)
592 {
593 Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
594 1);
595 }
596
597 // Add (0.5 g h h) to appropriate field
598 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
599 }
600}
601
603 const Array<OneD, const Array<OneD, NekDouble>> &physin,
605{
606 int nq = GetTotPoints();
607
608 if (physin.get() == physout.get())
609 {
610 // copy indata and work with tmp array
612 for (int i = 0; i < 3; ++i)
613 {
614 // deep copy
615 tmp[i] = Array<OneD, NekDouble>(nq);
616 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
617 }
618
619 // \eta = h - d
620 Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
621
622 // u = hu/h
623 Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
624
625 // v = hv/ v
626 Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
627 }
628 else
629 {
630 // \eta = h - d
631 Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
632
633 // u = hu/h
634 Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
635
636 // v = hv/ v
637 Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
638 }
639}
640
642{
643 int nq = GetTotPoints();
644
645 // u = hu/h
646 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
647 m_fields[1]->UpdatePhys(), 1);
648
649 // v = hv/ v
650 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
651 m_fields[2]->UpdatePhys(), 1);
652
653 // \eta = h - d
654 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
655 m_fields[0]->UpdatePhys(), 1);
656}
657
659 const Array<OneD, const Array<OneD, NekDouble>> &physin,
661{
662 int nq = GetTotPoints();
663
664 if (physin.get() == physout.get())
665 {
666 // copy indata and work with tmp array
668 for (int i = 0; i < 3; ++i)
669 {
670 // deep copy
671 tmp[i] = Array<OneD, NekDouble>(nq);
672 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
673 }
674
675 // h = \eta + d
676 Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
677
678 // hu = h * u
679 Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
680
681 // hv = h * v
682 Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
683 }
684 else
685 {
686 // h = \eta + d
687 Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
688
689 // hu = h * u
690 Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
691
692 // hv = h * v
693 Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
694 }
695}
696
698{
699 int nq = GetTotPoints();
700
701 // h = \eta + d
702 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
703 m_fields[0]->UpdatePhys(), 1);
704
705 // hu = h * u
706 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
707 m_fields[1]->UpdatePhys(), 1);
708
709 // hv = h * v
710 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
711 m_fields[2]->UpdatePhys(), 1);
712}
713
714/**
715 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
716 * \f$ h\mathbf{v} \f$.
717 *
718 * @param physfield Momentum field.
719 * @param velocity Velocity field.
720 */
722 const Array<OneD, Array<OneD, NekDouble>> &physfield,
724{
725 const int npts = physfield[0].size();
726
727 for (int i = 0; i < m_spacedim; ++i)
728 {
729 Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
730 }
731}
732
734{
736 SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
737 SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
738 SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
739}
740
741} // 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_PrimitiveToConservative() override
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
static std::string className
Name of class.
Definition: NonlinearSWE.h:63
void v_ConservativeToPrimitive() override
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Wall boundary condition.
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: NonlinearSWE.h:53
NonlinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
const Array< OneD, NekDouble > & GetDepth()
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 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