Nektar++
NonlinearPeregrine.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NonlinearPeregrine.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 Boussinesq equations of Peregrine in
32// conservative variables (constant depth case)
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <iomanip>
37#include <iostream>
38
39#include <boost/algorithm/string.hpp>
40#include <boost/core/ignore_unused.hpp>
41
44
45using namespace std;
46
47namespace Nektar
48{
49
52 "NonlinearPeregrine", NonlinearPeregrine::create,
53 "Nonlinear Peregrine equations in conservative variables.");
54
58 : ShallowWaterSystem(pSession, pGraph), m_factors()
59{
62 // note: eFactorTau = 1.0 becomes unstable...
63 // we need to investigate the behaviuor w.r.t. tau
64}
65
66void NonlinearPeregrine::v_InitObject(bool DeclareFields)
67{
69
70 if (m_session->DefinesSolverInfo("PROBLEMTYPE"))
71 {
72 int i;
73 std::string ProblemTypeStr = m_session->GetSolverInfo("PROBLEMTYPE");
74 for (i = 0; i < (int)SIZE_ProblemType; ++i)
75 {
76 if (boost::iequals(ProblemTypeMap[i], ProblemTypeStr))
77 {
79 break;
80 }
81 }
82 }
83 else
84 {
86 }
87
89 {
92 }
93 else
94 {
95 ASSERTL0(false, "Implicit Peregrine not set up.");
96 }
97
98 // NB! At the moment only the constant depth case is
99 // supported for the Peregrine eq.
100 if (m_session->DefinesParameter("ConstDepth"))
101 {
102 m_const_depth = m_session->GetParameter("ConstDepth");
103 }
104 else
105 {
106 ASSERTL0(false, "Constant Depth not specified");
107 }
108
109 // Type of advection class to be used
110 switch (m_projectionType)
111 {
112 // Continuous field
114 {
115 ASSERTL0(false,
116 "Continuous projection type not supported for Peregrine.");
117 break;
118 }
119 // Discontinuous field
121 {
122 string advName;
123 string diffName;
124 string riemName;
125
126 //---------------------------------------------------------------
127 // Setting up advection and diffusion operators
128 // NB: diffusion not set up for SWE at the moment
129 // but kept here for future use ...
130 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
131 // m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
133 advName, advName);
134
136 this);
137
138 // Setting up Riemann solver for advection operator
139 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
140
143 riemName, m_session);
144
145 // Setting up parameters for advection operator Riemann solver
146 m_riemannSolver->SetParam("gravity",
148 m_riemannSolver->SetAuxVec("vecLocs",
151 this);
153 this);
154
155 // Concluding initialisation of advection / diffusion operators
156 m_advection->SetRiemannSolver(m_riemannSolver);
157 m_advection->InitObject(m_session, m_fields);
158 break;
159 }
160 default:
161 {
162 ASSERTL0(false, "Unsupported projection type.");
163 break;
164 }
165 }
166}
167
169{
170}
171
172// physarray contains the conservative variables
174 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
176{
177
178 int ncoeffs = GetNcoeffs();
179 int nq = GetTotPoints();
180
182 Array<OneD, NekDouble> mod(ncoeffs);
183
184 switch (m_projectionType)
185 {
187 {
188 // add to hu equation
189 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
190 m_fields[0]->IProductWRTBase(tmp, mod);
191 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
192 m_fields[0]->BwdTrans(mod, tmp);
193 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
194
195 // add to hv equation
196 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
197 Vmath::Neg(nq, tmp, 1);
198 m_fields[0]->IProductWRTBase(tmp, mod);
199 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
200 m_fields[0]->BwdTrans(mod, tmp);
201 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
202 break;
203 }
206 {
207 // add to hu equation
208 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
209 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
210
211 // add to hv equation
212 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
213 Vmath::Neg(nq, tmp, 1);
214 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
215 break;
216 }
217 default:
218 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
219 break;
220 }
221}
222
223// physarray contains the conservative variables
225 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
227{
228
229 int ncoeffs = GetNcoeffs();
230 int nq = GetTotPoints();
231
233 Array<OneD, NekDouble> mod(ncoeffs);
234
235 switch (m_projectionType)
236 {
238 {
239 for (int i = 0; i < m_spacedim; ++i)
240 {
241 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
242 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
243 m_fields[0]->IProductWRTBase(tmp, mod);
244 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
245 m_fields[0]->BwdTrans(mod, tmp);
246 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
247 }
248 break;
249 }
252 {
253 for (int i = 0; i < m_spacedim; ++i)
254 {
255 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
256 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
257 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
258 }
259 break;
260 }
261 default:
262 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
263 break;
264 }
265}
266
268 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
269 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
270{
271 int i;
272 int nvariables = inarray.size();
273 int ncoeffs = GetNcoeffs();
274 int nq = GetTotPoints();
275
276 switch (m_projectionType)
277 {
279 {
280
281 //-------------------------------------------------------
282 // inarray in physical space
283
284 Array<OneD, Array<OneD, NekDouble>> modarray(nvariables);
285 for (i = 0; i < nvariables; ++i)
286 {
287 modarray[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
288 }
289 //-------------------------------------------------------
290
291 //-------------------------------------------------------
292 // Compute the DG advection including the numerical flux
293 // by using SolverUtils/Advection
294 // Input and output in physical space
296
297 m_advection->Advect(nvariables - 1, m_fields, advVel, inarray,
298 outarray, time);
299 //-------------------------------------------------------
300
301 //-------------------------------------------------------
302 // negate the outarray since moving terms to the rhs
303 for (i = 0; i < nvariables - 1; ++i)
304 {
305 Vmath::Neg(nq, outarray[i], 1);
306 }
307 //-------------------------------------------------------
308
309 //-------------------------------------------------
310 // Add "source terms"
311 // Input and output in physical space
312
313 // Coriolis forcing
314 if (m_coriolis.size() != 0)
315 {
316 AddCoriolis(inarray, outarray);
317 }
318
319 // Variable Depth
320 if (m_constantDepth != true)
321 {
322 ASSERTL0(false,
323 "Variable depth not supported for the Peregrine "
324 "equations");
325 }
326
327 //-------------------------------------------------
328
329 //---------------------------------------
330 // As no more terms is required for the
331 // continuity equation and we have aleady evaluated
332 // the values for h_t we are done for h
333 //---------------------------------------
334
335 //-------------------------------------------------
336 // go to modal space
337 m_fields[0]->IProductWRTBase(outarray[1], modarray[1]);
338 m_fields[0]->IProductWRTBase(outarray[2], modarray[2]);
339
340 // store f1 and f2 for later use (modal space)
341 Array<OneD, NekDouble> f1(ncoeffs);
342 Array<OneD, NekDouble> f2(ncoeffs);
343
344 Vmath::Vcopy(ncoeffs, modarray[1], 1, f1, 1); // f1
345 Vmath::Vcopy(ncoeffs, modarray[2], 1, f2, 1); // f2
346
347 // Solve the remaining block-diagonal systems
348 m_fields[0]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
349 m_fields[0]->MultiplyByElmtInvMass(modarray[2], modarray[2]);
350 //---------------------------------------------
351
352 //---------------------------------------------
353
354 //-------------------------------------------------
355 // create tmp fields to be used during
356 // the dispersive section
357
360
361 for (i = 0; i < 2; ++i)
362 {
363 coeffsfield[i] = Array<OneD, NekDouble>(ncoeffs);
364 physfield[i] = Array<OneD, NekDouble>(nq);
365 }
366 //---------------------------------------------
367
368 //---------------------------------------------
369 // Go from modal to physical space
370 Vmath::Vcopy(nq, outarray[1], 1, physfield[0], 1);
371 Vmath::Vcopy(nq, outarray[2], 1, physfield[1], 1);
372 //---------------------------------------
373
374 //---------------------------------------
375 // Start for solve of mixed dispersive terms
376 // using the 'WCE method'
377 // (Eskilsson & Sherwin, JCP 2006)
378
379 // constant depth case
380 // \nabla \cdot (\nabla z) - invgamma z
381 // = - invgamma (\nabla \cdot {\bf f}_(2,3)
382
383 NekDouble gamma = (m_const_depth * m_const_depth) * (1.0 / 3.0);
384 NekDouble invgamma = 1.0 / gamma;
385
386 int nTraceNumPoints = GetTraceTotPoints();
387 Array<OneD, NekDouble> upwindX(nTraceNumPoints);
388 Array<OneD, NekDouble> upwindY(nTraceNumPoints);
389 //--------------------------------------------
390
391 //--------------------------------------------
392 // Compute the forcing function for the
393 // wave continuity equation
394
395 // Set boundary condidtions for z
396 SetBoundaryConditionsForcing(physfield, time);
397
398 // \nabla \phi \cdot f_{2,3}
399 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
400 m_fields[0]->IProductWRTDerivBase(1, physfield[1], coeffsfield[1]);
401 Vmath::Vadd(ncoeffs, coeffsfield[0], 1, coeffsfield[1], 1,
402 coeffsfield[0], 1);
403 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
404
405 // Evaluate upwind numerical flux (physical space)
406 NumericalFluxForcing(physfield, upwindX, upwindY);
407 Array<OneD, NekDouble> normflux(nTraceNumPoints);
408 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0], 1,
409 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
410 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
411 m_fields[0]->MultiplyByElmtInvMass(coeffsfield[0], coeffsfield[0]);
412 m_fields[0]->BwdTrans(coeffsfield[0], physfield[0]);
413
414 Vmath::Smul(nq, -invgamma, physfield[0], 1, physfield[0], 1);
415
416 // ok: forcing function for HelmSolve... done!
417 //--------------------------------------
418
419 //--------------------------------------
420 // Solve the Helmhotz-type equation
421 // for the wave continuity equation
422 // (missing slope terms...)
423
424 // note: this is just valid for the constant depth case:
425
426 // as of now we need not to specify any
427 // BC routine for the WCE: periodic
428 // and zero Neumann (for walls)
429
430 WCESolve(physfield[0], invgamma);
431
432 Vmath::Vcopy(nq, physfield[0], 1, outarray[3], 1); // store z
433
434 // ok: Wave Continuity Equation... done!
435 //------------------------------------
436
437 //------------------------------------
438 // Return to the primary variables
439
440 // (h {\bf u})_t = gamma \nabla z + {\bf f}_{2,3}
441
442 Vmath::Smul(nq, gamma, physfield[0], 1, physfield[0], 1);
443
444 // Set boundary conditions
445 SetBoundaryConditionsContVariables(physfield[0], time);
446
447 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
448 m_fields[1]->IProductWRTDerivBase(1, physfield[0], coeffsfield[1]);
449
450 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
451 Vmath::Neg(ncoeffs, coeffsfield[1], 1);
452
453 // Evaluate upwind numerical flux (physical space)
454 NumericalFluxConsVariables(physfield[0], upwindX, upwindY);
455
456 {
457 Array<OneD, NekDouble> uptemp(nTraceNumPoints, 0.0);
458 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0],
459 1, uptemp, 1, m_traceNormals[1], 1, normflux, 1);
460 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
461 Vmath::Vvtvvtp(nTraceNumPoints, uptemp, 1, m_traceNormals[0], 1,
462 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
463 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[1]);
464 }
465
466 Vmath::Vadd(ncoeffs, f1, 1, coeffsfield[0], 1, modarray[1], 1);
467 Vmath::Vadd(ncoeffs, f2, 1, coeffsfield[1], 1, modarray[2], 1);
468
469 m_fields[1]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
470 m_fields[2]->MultiplyByElmtInvMass(modarray[2], modarray[2]);
471
472 m_fields[1]->BwdTrans(modarray[1], outarray[1]);
473 m_fields[2]->BwdTrans(modarray[2], outarray[2]);
474
475 // ok: returned to conservative variables... done!
476 //---------------------
477
478 break;
479 }
482 ASSERTL0(false, "Unknown projection scheme for the Peregrine");
483 break;
484 default:
485 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
486 break;
487 }
488}
489
491 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
492 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
493{
494 int i;
495 int nvariables = inarray.size();
496
497 switch (m_projectionType)
498 {
500 {
501
502 // Just copy over array
503 if (inarray != outarray)
504 {
505 int npoints = GetNpoints();
506
507 for (i = 0; i < nvariables; ++i)
508 {
509 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
510 }
511 }
512
513 SetBoundaryConditions(outarray, time);
514 break;
515 }
518 {
519
521 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
522
523 for (i = 0; i < nvariables; ++i)
524 {
525 m_fields[i]->FwdTrans(inarray[i], coeffs);
526 m_fields[i]->BwdTrans(coeffs, outarray[i]);
527 }
528 break;
529 }
530 default:
531 ASSERTL0(false, "Unknown projection scheme");
532 break;
533 }
534}
535
536//----------------------------------------------------
539{
540
541 int nvariables = m_fields.size();
542 int cnt = 0;
543 int nTracePts = GetTraceTotPoints();
544
545 // Extract trace for boundaries. Needs to be done on all processors to avoid
546 // deadlock.
547 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
548 for (int i = 0; i < nvariables; ++i)
549 {
550 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
551 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
552 }
553
554 // loop over Boundary Regions
555 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
556 {
557
558 // Wall Boundary Condition
559 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
560 "Wall"))
561 {
562 WallBoundary2D(n, cnt, Fwd, inarray);
563 }
564
565 // Time Dependent Boundary Condition (specified in meshfile)
566 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
567 {
568 for (int i = 0; i < nvariables; ++i)
569 {
570 m_fields[i]->EvaluateBoundaryConditions(time);
571 }
572 }
573 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
574 }
575}
576
577//----------------------------------------------------
578/**
579 * @brief Wall boundary condition.
580 */
582 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd,
583 Array<OneD, Array<OneD, NekDouble>> &physarray)
584{
585 int i;
586 int nvariables = physarray.size();
587
588 // Adjust the physical values of the trace to take
589 // user defined boundaries into account
590 int e, id1, id2, npts;
592 m_fields[0]->GetBndCondExpansions()[bcRegion];
593 for (e = 0; e < bcexp->GetExpSize(); ++e)
594 {
595 npts = bcexp->GetExp(e)->GetTotPoints();
596 id1 = bcexp->GetPhys_Offset(e);
597 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
598 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
599
600 // For 2D/3D, define: v* = v - 2(v.n)n
601 Array<OneD, NekDouble> tmp(npts, 0.0);
602
603 // Calculate (v.n)
604 for (i = 0; i < m_spacedim; ++i)
605 {
606 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
607 &tmp[0], 1, &tmp[0], 1);
608 }
609
610 // Calculate 2.0(v.n)
611 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
612
613 // Calculate v* = v - 2.0(v.n)n
614 for (i = 0; i < m_spacedim; ++i)
615 {
616 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
617 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
618 }
619
620 // copy boundary adjusted values into the boundary expansion
621 for (i = 0; i < nvariables; ++i)
622 {
623 bcexp = m_fields[i]->GetBndCondExpansions()[bcRegion];
624 Vmath::Vcopy(npts, &Fwd[i][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
625 }
626 }
627}
628
630 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd,
631 Array<OneD, Array<OneD, NekDouble>> &physarray)
632{
633 boost::ignore_unused(physarray);
634
635 int i;
636 int nvariables = 3;
637
638 // Adjust the physical values of the trace to take
639 // user defined boundaries into account
640 int e, id1, id2, npts;
642 m_fields[0]->GetBndCondExpansions()[bcRegion];
643
644 for (e = 0; e < bcexp->GetExpSize(); ++e)
645 {
646 npts = bcexp->GetExp(e)->GetNumPoints(0);
647 id1 = bcexp->GetPhys_Offset(e);
648 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
649 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
650
651 switch (m_expdim)
652 {
653 case 1:
654 {
655 // negate the forward flux
656 Vmath::Neg(npts, &Fwd[1][id2], 1);
657 break;
658 }
659 case 2:
660 {
661 Array<OneD, NekDouble> tmp_n(npts);
662 Array<OneD, NekDouble> tmp_t(npts);
663
664 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
665 &tmp_n[0], 1);
666 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
667 &tmp_n[0], 1, &tmp_n[0], 1);
668
669 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
670 &tmp_t[0], 1);
671 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
672 &tmp_t[0], 1, &tmp_t[0], 1);
673
674 // negate the normal flux
675 Vmath::Neg(npts, tmp_n, 1);
676
677 // rotate back to Cartesian
678 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
679 &Fwd[1][id2], 1);
680 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
681 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
682
683 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
684 &Fwd[2][id2], 1);
685 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
686 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
687 break;
688 }
689 case 3:
690 ASSERTL0(false,
691 "3D not implemented for Shallow Water Equations");
692 break;
693 default:
694 ASSERTL0(false, "Illegal expansion dimension");
695 }
696
697 // copy boundary adjusted values into the boundary expansion
698 for (i = 0; i < nvariables; ++i)
699 {
700 bcexp = m_fields[i]->GetBndCondExpansions()[bcRegion];
701 Vmath::Vcopy(npts, &Fwd[i][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
702 }
703 }
704}
705
706// Physfield in conservative Form
708 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
710{
711 int i, j;
712 int nq = m_fields[0]->GetTotPoints();
713
714 NekDouble g = m_g;
716
717 // Flux vector for the mass equation
718 for (i = 0; i < m_spacedim; ++i)
719 {
720 velocity[i] = Array<OneD, NekDouble>(nq);
721 Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
722 }
723
724 GetVelocityVector(physfield, velocity);
725
726 // Put (0.5 g h h) in tmp
728 Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
729 Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
730
731 // Flux vector for the momentum equations
732 for (i = 0; i < m_spacedim; ++i)
733 {
734 for (j = 0; j < m_spacedim; ++j)
735 {
736 Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
737 1);
738 }
739
740 // Add (0.5 g h h) to appropriate field
741 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
742 }
743}
744
746 const Array<OneD, const Array<OneD, NekDouble>> &physin,
748{
749 int nq = GetTotPoints();
750
751 if (physin.get() == physout.get())
752 {
753 // copy indata and work with tmp array
755 for (int i = 0; i < 3; ++i)
756 {
757 // deep copy
758 tmp[i] = Array<OneD, NekDouble>(nq);
759 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
760 }
761
762 // \eta = h - d
763 Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
764
765 // u = hu/h
766 Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
767
768 // v = hv/ v
769 Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
770 }
771 else
772 {
773 // \eta = h - d
774 Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
775
776 // u = hu/h
777 Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
778
779 // v = hv/ v
780 Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
781 }
782}
783
785{
786 int nq = GetTotPoints();
787
788 // u = hu/h
789 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
790 m_fields[1]->UpdatePhys(), 1);
791
792 // v = hv/ v
793 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
794 m_fields[2]->UpdatePhys(), 1);
795
796 // \eta = h - d
797 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
798 m_fields[0]->UpdatePhys(), 1);
799}
800
802 const Array<OneD, const Array<OneD, NekDouble>> &physin,
804{
805
806 int nq = GetTotPoints();
807
808 if (physin.get() == physout.get())
809 {
810 // copy indata and work with tmp array
812 for (int i = 0; i < 3; ++i)
813 {
814 // deep copy
815 tmp[i] = Array<OneD, NekDouble>(nq);
816 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
817 }
818
819 // h = \eta + d
820 Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
821
822 // hu = h * u
823 Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
824
825 // hv = h * v
826 Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
827 }
828 else
829 {
830 // h = \eta + d
831 Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
832
833 // hu = h * u
834 Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
835
836 // hv = h * v
837 Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
838 }
839}
840
842{
843 int nq = GetTotPoints();
844
845 // h = \eta + d
846 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
847 m_fields[0]->UpdatePhys(), 1);
848
849 // hu = h * u
850 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
851 m_fields[1]->UpdatePhys(), 1);
852
853 // hv = h * v
854 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
855 m_fields[2]->UpdatePhys(), 1);
856}
857
858/**
859 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
860 * \f$ h\mathbf{v} \f$.
861 *
862 * @param physfield Momentum field.
863 * @param velocity Velocity field.
864 */
866 const Array<OneD, Array<OneD, NekDouble>> &physfield,
868{
869 const int npts = physfield[0].size();
870
871 for (int i = 0; i < m_spacedim; ++i)
872 {
873 Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
874 }
875}
876
878{
880 SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
881 SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
882 SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
883 SolverUtils::AddSummaryItem(s, "", "z should be in field[3]");
884}
885
887{
888 int nq = GetTotPoints();
889
891
892 for (int j = 0; j < nq; j++)
893 {
894 (m_fields[3]->UpdatePhys())[j] = fce[j];
895 }
896
897 m_fields[3]->SetPhysState(true);
898
899 m_fields[3]->HelmSolve(m_fields[3]->GetPhys(), m_fields[3]->UpdateCoeffs(),
900 m_factors);
901
902 m_fields[3]->BwdTrans(m_fields[3]->GetCoeffs(), m_fields[3]->UpdatePhys());
903
904 m_fields[3]->SetPhysState(true);
905
906 Vmath::Vcopy(nq, m_fields[3]->GetPhys(), 1, fce, 1);
907}
908
910 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
912{
913 int i;
914 int nTraceNumPoints = GetTraceTotPoints();
915
916 //-----------------------------------------------------
917 // get temporary arrays
920
921 for (i = 0; i < 2; ++i)
922 {
923 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
924 Bwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
925 }
926 //-----------------------------------------------------
927
928 //-----------------------------------------------------
929 // get the physical values at the trace
930 // (any time-dependent BC previuosly put in fields[1] and [2]
931
932 m_fields[1]->GetFwdBwdTracePhys(inarray[0], Fwd[0], Bwd[0]);
933 m_fields[2]->GetFwdBwdTracePhys(inarray[1], Fwd[1], Bwd[1]);
934 //-----------------------------------------------------
935
936 //-----------------------------------------------------
937 // use centred fluxes for the numerical flux
938 for (i = 0; i < nTraceNumPoints; ++i)
939 {
940 numfluxX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
941 numfluxY[i] = 0.5 * (Fwd[1][i] + Bwd[1][i]);
942 }
943 //-----------------------------------------------------
944}
945
948{
949 boost::ignore_unused(time);
950
951 int cnt = 0;
952
953 // loop over Boundary Regions
954 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
955 {
956 // Use wall for all BC...
957 // Wall Boundary Condition
958 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
959 "Wall"))
960 {
961 WallBoundaryForcing(n, cnt, inarray);
962 }
963
964 // Timedependent Boundary Condition
965 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
966 {
967 ASSERTL0(false, "time-dependent BC not implemented for Boussinesq");
968 }
969 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
970 }
971}
972
973// fills up boundary expansion for field[1] and [2]
975 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &inarray)
976{
977
978 // std::cout << " WallBoundaryForcing" << std::endl;
979
980 int nTraceNumPoints = GetTraceTotPoints();
981 int nvariables = 2;
982
983 // get physical values of f1 and f2 for the forward trace
984 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
985 for (int i = 0; i < nvariables; ++i)
986 {
987 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
988 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
989 }
990
991 // Adjust the physical values of the trace to take
992 // user defined boundaries into account
993 int e, id1, id2, npts;
995 m_fields[0]->GetBndCondExpansions()[bcRegion];
996 for (e = 0; e < bcexp->GetExpSize(); ++e)
997 {
998 npts = bcexp->GetExp(e)->GetTotPoints();
999 id1 = bcexp->GetPhys_Offset(e);
1000 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
1001 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
1002
1003 switch (m_expdim)
1004 {
1005 case 1:
1006 {
1007 ASSERTL0(false, "1D not yet implemented for Boussinesq");
1008 break;
1009 }
1010 case 2:
1011 {
1012 Array<OneD, NekDouble> tmp_n(npts);
1013 Array<OneD, NekDouble> tmp_t(npts);
1014
1015 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[0][id2], 1,
1016 &tmp_n[0], 1);
1017 Vmath::Vvtvp(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
1018 &tmp_n[0], 1, &tmp_n[0], 1);
1019
1020 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[1][id2], 1,
1021 &tmp_t[0], 1);
1022 Vmath::Vvtvm(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
1023 &tmp_t[0], 1, &tmp_t[0], 1);
1024
1025 // negate the normal flux
1026 Vmath::Neg(npts, tmp_n, 1);
1027
1028 // rotate back to Cartesian
1029 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
1030 &Fwd[0][id2], 1);
1031 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
1032 &Fwd[0][id2], 1, &Fwd[0][id2], 1);
1033
1034 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
1035 &Fwd[1][id2], 1);
1036 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
1037 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
1038 break;
1039 }
1040 case 3:
1041 ASSERTL0(false, "3D not implemented for Boussinesq equations");
1042 break;
1043 default:
1044 ASSERTL0(false, "Illegal expansion dimension");
1045 }
1046
1047 // copy boundary adjusted values into the boundary expansion
1048 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
1049 Vmath::Vcopy(npts, &Fwd[0][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1050
1051 bcexp = m_fields[2]->GetBndCondExpansions()[bcRegion];
1052 Vmath::Vcopy(npts, &Fwd[1][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1053 }
1054}
1055
1057 Array<OneD, NekDouble> &inarray, NekDouble time)
1058{
1059 boost::ignore_unused(time);
1060
1061 int cnt = 0;
1062
1063 // loop over Boundary Regions
1064 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
1065 {
1066 // Use wall for all
1067 // Wall Boundary Condition
1068 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
1069 "Wall"))
1070 {
1071 WallBoundaryContVariables(n, cnt, inarray);
1072 }
1073
1074 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
1075 {
1076 WallBoundaryContVariables(n, cnt, inarray);
1077 }
1078
1079 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize() - 1;
1080 }
1081}
1082
1084 int bcRegion, int cnt, Array<OneD, NekDouble> &inarray)
1085{
1086 int nTraceNumPoints = GetTraceTotPoints();
1087
1088 // get physical values of z for the forward trace
1089 Array<OneD, NekDouble> z(nTraceNumPoints);
1090 m_fields[0]->ExtractTracePhys(inarray, z);
1091
1092 // Adjust the physical values of the trace to take
1093 // user defined boundaries into account
1094 int e, id1, id2, npts;
1096 m_fields[0]->GetBndCondExpansions()[bcRegion];
1097
1098 for (e = 0; e < bcexp->GetExpSize(); ++e)
1099 {
1100 npts = bcexp->GetExp(e)->GetTotPoints();
1101 id1 = bcexp->GetPhys_Offset(e);
1102 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
1103 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
1104
1105 // copy boundary adjusted values into the boundary expansion
1106 // field[1] and field[2]
1107 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
1108 Vmath::Vcopy(npts, &z[id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1109 }
1110}
1111
1115{
1116 int i;
1117 int nTraceNumPoints = GetTraceTotPoints();
1118
1119 //-----------------------------------------------------
1120 // get temporary arrays
1123
1124 Fwd[0] = Array<OneD, NekDouble>(nTraceNumPoints);
1125 Bwd[0] = Array<OneD, NekDouble>(nTraceNumPoints);
1126 //-----------------------------------------------------
1127
1128 //-----------------------------------------------------
1129 // get the physical values at the trace
1130 // (we have put any time-dependent BC in field[1])
1131
1132 m_fields[1]->GetFwdBwdTracePhys(physfield, Fwd[0], Bwd[0]);
1133 //-----------------------------------------------------
1134
1135 //-----------------------------------------------------
1136 // use centred fluxes for the numerical flux
1137 for (i = 0; i < nTraceNumPoints; ++i)
1138 {
1139 outX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
1140 outY[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
1141 }
1142 //-----------------------------------------------------
1143}
1144
1145// initial condition Laitone's first order solitary wave
1147 NekDouble time, NekDouble x_offset)
1148{
1149 int nq = GetTotPoints();
1150
1151 NekDouble A = 1.0;
1152 NekDouble C = sqrt(m_g * d) * (1.0 + 0.5 * (amp / d));
1153
1156 Array<OneD, NekDouble> zeros(nq, 0.0);
1157
1158 // get the coordinates (assuming all fields have the same
1159 // discretisation)
1160 m_fields[0]->GetCoords(x0, x1);
1161
1162 for (int i = 0; i < nq; i++)
1163 {
1164 (m_fields[0]->UpdatePhys())[i] =
1165 amp * pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
1166 (A * (x0[i] + x_offset) - C * time))),
1167 2.0);
1168 (m_fields[1]->UpdatePhys())[i] =
1169 (amp / d) *
1170 pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
1171 (A * (x0[i] + x_offset) - C * time))),
1172 2.0) *
1173 sqrt(m_g * d);
1174 }
1175
1176 Vmath::Sadd(nq, d, m_fields[0]->GetPhys(), 1, m_fields[0]->UpdatePhys(), 1);
1177 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
1178 m_fields[1]->UpdatePhys(), 1);
1179 Vmath::Vcopy(nq, zeros, 1, m_fields[2]->UpdatePhys(), 1);
1180 Vmath::Vcopy(nq, zeros, 1, m_fields[3]->UpdatePhys(), 1);
1181
1182 // Forward transform to fill the coefficient space
1183 for (int i = 0; i < 4; ++i)
1184 {
1185 m_fields[i]->SetPhysState(true);
1186 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1187 m_fields[i]->UpdateCoeffs());
1188 }
1189}
1190
1191/**
1192 * @brief Set the initial conditions.
1193 */
1195 bool dumpInitialConditions,
1196 const int domain)
1197{
1198 boost::ignore_unused(domain);
1199
1200 switch (m_problemType)
1201 {
1202 case eSolitaryWave:
1203 {
1204 LaitoneSolitaryWave(0.1, m_const_depth, 0.0, 0.0);
1205 m_nchk++;
1206 break;
1207 }
1208 default:
1209 {
1210 EquationSystem::v_SetInitialConditions(initialtime, false);
1211 break;
1212 }
1213 }
1214
1215 if (dumpInitialConditions)
1216 {
1217 // Dump initial conditions to file
1219 }
1220}
1221
1222} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
virtual void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
Set the initial conditions.
virtual void v_PrimitiveToConservative() override
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
virtual void v_ConservativeToPrimitive() override
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void SetBoundaryConditionsForcing(Array< OneD, Array< OneD, NekDouble > > &inarray, NekDouble time)
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void WallBoundaryForcing(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &inarray)
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void LaitoneSolitaryWave(NekDouble amp, NekDouble d, NekDouble time, NekDouble x_offset)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
static std::string className
Name of class.
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Wall boundary condition.
void NumericalFluxForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &numfluxX, Array< OneD, NekDouble > &numfluxY)
void SetBoundaryConditionsContVariables(Array< OneD, NekDouble > &inarray, NekDouble time)
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
void WallBoundaryContVariables(int bcRegion, int cnt, Array< OneD, NekDouble > &inarray)
void NumericalFluxConsVariables(Array< OneD, NekDouble > &physfield, Array< OneD, NekDouble > &outX, Array< OneD, NekDouble > &outY)
StdRegions::ConstFactorMap m_factors
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
NonlinearPeregrine(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
virtual ~NonlinearPeregrine()
problem type selector
void WCESolve(Array< OneD, NekDouble > &fce, NekDouble lambda)
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.
virtual 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.
virtual 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.
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
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()
int m_nchk
Number of checkpoints written so far.
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
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
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:49
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
const char *const ProblemTypeMap[]
@ SIZE_ProblemType
Length of enum list.
@ eSolitaryWave
First order Laitone solitary wave.
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.cpp:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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.cpp:569
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.cpp:354
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.cpp:593
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.cpp:245
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.cpp:280
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:379
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.cpp:687
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
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.cpp:414
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294