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