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
41
42namespace Nektar
43{
44
47 "NonlinearPeregrine", NonlinearPeregrine::create,
48 "Nonlinear Peregrine equations in conservative variables.");
49
53 : NonlinearSWE(pSession, pGraph), m_factors()
54{
57 // note: eFactorTau = 1.0 becomes unstable...
58 // we need to investigate the behaviuor w.r.t. tau
59}
60
61void NonlinearPeregrine::v_InitObject(bool DeclareFields)
62{
64
65 if (m_session->DefinesSolverInfo("PROBLEMTYPE"))
66 {
67 int i;
68 std::string ProblemTypeStr = m_session->GetSolverInfo("PROBLEMTYPE");
69 for (i = 0; i < (int)SIZE_ProblemType; ++i)
70 {
71 if (boost::iequals(ProblemTypeMap[i], ProblemTypeStr))
72 {
74 break;
75 }
76 }
77 }
78 else
79 {
81 }
82
84 {
87 }
88 else
89 {
90 ASSERTL0(false, "Implicit Peregrine not set up.");
91 }
92
93 // NB! At the moment only the constant depth case is
94 // supported for the Peregrine eq.
95 if (m_session->DefinesParameter("ConstDepth"))
96 {
97 m_const_depth = m_session->GetParameter("ConstDepth");
98 }
99 else
100 {
101 ASSERTL0(false, "Constant Depth not specified");
102 }
103
104 // Type of advection class to be used
105 switch (m_projectionType)
106 {
107 // Continuous field
109 {
110 ASSERTL0(false,
111 "Continuous projection type not supported for Peregrine.");
112 break;
113 }
114 // Discontinuous field
116 {
117 std::string advName;
118 std::string diffName;
119 std::string riemName;
120
121 //---------------------------------------------------------------
122 // Setting up advection and diffusion operators
123 // NB: diffusion not set up for SWE at the moment
124 // but kept here for future use ...
125 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
126 // m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
128 advName, advName);
129
131 this);
132
133 // Setting up Riemann solver for advection operator
134 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
135
138 riemName, m_session);
139
140 // Setting up parameters for advection operator Riemann solver
141 m_riemannSolver->SetParam("gravity",
143 m_riemannSolver->SetAuxVec("vecLocs",
146 this);
148 this);
149
150 // Concluding initialisation of advection / diffusion operators
151 m_advection->SetRiemannSolver(m_riemannSolver);
152 m_advection->InitObject(m_session, m_fields);
153 break;
154 }
155 default:
156 {
157 ASSERTL0(false, "Unsupported projection type.");
158 break;
159 }
160 }
161}
162
164{
166 SolverUtils::AddSummaryItem(s, "", "z should be in field[3]");
167}
168
169/**
170 * @brief Set the initial conditions.
171 */
173 NekDouble initialtime, bool dumpInitialConditions,
174 [[maybe_unused]] const int domain)
175{
176 switch (m_problemType)
177 {
178 case eSolitaryWave:
179 {
180 LaitoneSolitaryWave(0.1, m_const_depth, 0.0, 0.0);
181 break;
182 }
183 default:
184 {
185 EquationSystem::v_SetInitialConditions(initialtime, false);
186 m_nchk--; // Note: m_nchk has been incremented in EquationSystem.
187 break;
188 }
189 }
190
191 if (dumpInitialConditions && m_checksteps && m_nchk == 0 &&
192 !m_comm->IsParallelInTime())
193 {
195 }
196 else if (dumpInitialConditions && m_nchk == 0 && m_comm->IsParallelInTime())
197 {
198 std::string newdir = m_sessionName + ".pit";
199 if (!fs::is_directory(newdir))
200 {
201 fs::create_directory(newdir);
202 }
203 if (m_comm->GetTimeComm()->GetRank() == 0)
204 {
205 WriteFld(newdir + "/" + m_sessionName + "_0.fld");
206 }
207 }
208 m_nchk++;
209}
210
212 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
213 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
214{
215 int i;
216 int nvariables = inarray.size();
217 int ncoeffs = GetNcoeffs();
218 int nq = GetTotPoints();
219
220 switch (m_projectionType)
221 {
223 {
224
225 //-------------------------------------------------------
226 // inarray in physical space
228 for (i = 0; i < 2; ++i)
229 {
230 modarray[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
231 }
232 //-------------------------------------------------------
233
234 //-------------------------------------------------------
235 // Compute the DG advection including the numerical flux
236 // by using SolverUtils/Advection
237 // Input and output in physical space
238 m_advection->Advect(nvariables - 1, m_fields,
239 NullNekDoubleArrayOfArray, inarray, outarray,
240 time);
241 //-------------------------------------------------------
242
243 //-------------------------------------------------------
244 // negate the outarray since moving terms to the rhs
245 for (i = 0; i < nvariables - 1; ++i)
246 {
247 Vmath::Neg(nq, outarray[i], 1);
248 }
249 //-------------------------------------------------------
250
251 //-------------------------------------------------
252 // Add "source terms"
253 // Input and output in physical space
254
255 // Coriolis forcing
256 if (m_coriolis.size() != 0)
257 {
258 AddCoriolis(inarray, outarray);
259 }
260
261 // Variable Depth
262 if (m_constantDepth != true)
263 {
264 ASSERTL0(false,
265 "Variable depth not supported for the Peregrine "
266 "equations");
267 }
268
269 //-------------------------------------------------
270
271 //---------------------------------------
272 // As no more terms is required for the
273 // continuity equation and we have aleady evaluated
274 // the values for h_t we are done for h
275 //---------------------------------------
276
277 //-------------------------------------------------
278 // go to modal space
279 m_fields[0]->IProductWRTBase(outarray[1], modarray[0]);
280 m_fields[0]->IProductWRTBase(outarray[2], modarray[1]);
281
282 // store f1 and f2 for later use (modal space)
283 Array<OneD, NekDouble> f1(ncoeffs);
284 Array<OneD, NekDouble> f2(ncoeffs);
285
286 Vmath::Vcopy(ncoeffs, modarray[0], 1, f1, 1); // f1
287 Vmath::Vcopy(ncoeffs, modarray[1], 1, f2, 1); // f2
288
289 // Solve the remaining block-diagonal systems
290 m_fields[0]->MultiplyByElmtInvMass(modarray[0], modarray[0]);
291 m_fields[0]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
292 //---------------------------------------------
293
294 //---------------------------------------------
295
296 //-------------------------------------------------
297 // create tmp fields to be used during
298 // the dispersive section
299
302
303 for (i = 0; i < 2; ++i)
304 {
305 coeffsfield[i] = Array<OneD, NekDouble>(ncoeffs);
306 physfield[i] = Array<OneD, NekDouble>(nq);
307 }
308 //---------------------------------------------
309
310 //---------------------------------------------
311 // Go from modal to physical space
312 Vmath::Vcopy(nq, outarray[1], 1, physfield[0], 1);
313 Vmath::Vcopy(nq, outarray[2], 1, physfield[1], 1);
314 //---------------------------------------
315
316 //---------------------------------------
317 // Start for solve of mixed dispersive terms
318 // using the 'WCE method'
319 // (Eskilsson & Sherwin, JCP 2006)
320
321 // constant depth case
322 // \nabla \cdot (\nabla z) - invgamma z
323 // = - invgamma (\nabla \cdot {\bf f}_(2,3)
324
325 NekDouble gamma = (m_const_depth * m_const_depth) * (1.0 / 3.0);
326 NekDouble invgamma = 1.0 / gamma;
327
328 int nTraceNumPoints = GetTraceTotPoints();
329 Array<OneD, NekDouble> upwindX(nTraceNumPoints);
330 Array<OneD, NekDouble> upwindY(nTraceNumPoints);
331 //--------------------------------------------
332
333 //--------------------------------------------
334 // Compute the forcing function for the
335 // wave continuity equation
336
337 // Set boundary condidtions for z
338 SetBoundaryConditionsForcing(physfield, time);
339
340 // \nabla \phi \cdot f_{2,3}
341 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
342 m_fields[0]->IProductWRTDerivBase(1, physfield[1], coeffsfield[1]);
343 Vmath::Vadd(ncoeffs, coeffsfield[0], 1, coeffsfield[1], 1,
344 coeffsfield[0], 1);
345 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
346
347 // Evaluate upwind numerical flux (physical space)
348 NumericalFluxForcing(physfield, upwindX, upwindY);
349 Array<OneD, NekDouble> normflux(nTraceNumPoints);
350 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0], 1,
351 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
352 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
353 m_fields[0]->MultiplyByElmtInvMass(coeffsfield[0], coeffsfield[0]);
354 m_fields[0]->BwdTrans(coeffsfield[0], physfield[0]);
355
356 Vmath::Smul(nq, -invgamma, physfield[0], 1, physfield[0], 1);
357
358 // ok: forcing function for HelmSolve... done!
359 //--------------------------------------
360
361 //--------------------------------------
362 // Solve the Helmhotz-type equation
363 // for the wave continuity equation
364 // (missing slope terms...)
365
366 // note: this is just valid for the constant depth case:
367
368 // as of now we need not to specify any
369 // BC routine for the WCE: periodic
370 // and zero Neumann (for walls)
371
372 WCESolve(physfield[0], invgamma);
373
374 Vmath::Vcopy(nq, physfield[0], 1, outarray[3], 1); // store z
375
376 // ok: Wave Continuity Equation... done!
377 //------------------------------------
378
379 //------------------------------------
380 // Return to the primary variables
381
382 // (h {\bf u})_t = gamma \nabla z + {\bf f}_{2,3}
383
384 Vmath::Smul(nq, gamma, physfield[0], 1, physfield[0], 1);
385
386 // Set boundary conditions
387 SetBoundaryConditionsContVariables(physfield[0], time);
388
389 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
390 m_fields[0]->IProductWRTDerivBase(1, physfield[0], coeffsfield[1]);
391
392 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
393 Vmath::Neg(ncoeffs, coeffsfield[1], 1);
394
395 // Evaluate upwind numerical flux (physical space)
396 NumericalFluxConsVariables(physfield[0], upwindX, upwindY);
397
398 {
399 Array<OneD, NekDouble> uptemp(nTraceNumPoints, 0.0);
400 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0],
401 1, uptemp, 1, m_traceNormals[1], 1, normflux, 1);
402 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
403 Vmath::Vvtvvtp(nTraceNumPoints, uptemp, 1, m_traceNormals[0], 1,
404 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
405 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[1]);
406 }
407
408 Vmath::Vadd(ncoeffs, f1, 1, coeffsfield[0], 1, modarray[0], 1);
409 Vmath::Vadd(ncoeffs, f2, 1, coeffsfield[1], 1, modarray[1], 1);
410
411 m_fields[1]->MultiplyByElmtInvMass(modarray[0], modarray[0]);
412 m_fields[2]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
413
414 m_fields[1]->BwdTrans(modarray[0], outarray[1]);
415 m_fields[2]->BwdTrans(modarray[1], outarray[2]);
416
417 // ok: returned to conservative variables... done!
418 //---------------------
419
420 break;
421 }
424 ASSERTL0(false, "Unknown projection scheme for the Peregrine");
425 break;
426 default:
427 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
428 break;
429 }
430}
431
432// initial condition Laitone's first order solitary wave
434 NekDouble time, NekDouble x_offset)
435{
436 int nq = GetTotPoints();
437
438 NekDouble A = 1.0;
439 NekDouble C = sqrt(m_g * d) * (1.0 + 0.5 * (amp / d));
440
443 Array<OneD, NekDouble> zeros(nq, 0.0);
444
445 // get the coordinates (assuming all fields have the same
446 // discretisation)
447 m_fields[0]->GetCoords(x0, x1);
448
449 for (int i = 0; i < nq; i++)
450 {
451 (m_fields[0]->UpdatePhys())[i] =
452 amp * pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
453 (A * (x0[i] + x_offset) - C * time))),
454 2.0);
455 (m_fields[1]->UpdatePhys())[i] =
456 (amp / d) *
457 pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
458 (A * (x0[i] + x_offset) - C * time))),
459 2.0) *
460 sqrt(m_g * d);
461 }
462
463 Vmath::Sadd(nq, d, m_fields[0]->GetPhys(), 1, m_fields[0]->UpdatePhys(), 1);
464 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
465 m_fields[1]->UpdatePhys(), 1);
466 Vmath::Vcopy(nq, zeros, 1, m_fields[2]->UpdatePhys(), 1);
467 Vmath::Vcopy(nq, zeros, 1, m_fields[3]->UpdatePhys(), 1);
468
469 // Forward transform to fill the coefficient space
470 for (int i = 0; i < 4; ++i)
471 {
472 m_fields[i]->SetPhysState(true);
473 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
474 m_fields[i]->UpdateCoeffs());
475 }
476}
477
479{
480 int nq = GetTotPoints();
481
483
484 for (int j = 0; j < nq; j++)
485 {
486 (m_fields[3]->UpdatePhys())[j] = fce[j];
487 }
488
489 m_fields[3]->SetPhysState(true);
490
491 m_fields[3]->HelmSolve(m_fields[3]->GetPhys(), m_fields[3]->UpdateCoeffs(),
492 m_factors);
493
494 m_fields[3]->BwdTrans(m_fields[3]->GetCoeffs(), m_fields[3]->UpdatePhys());
495
496 m_fields[3]->SetPhysState(true);
497
498 Vmath::Vcopy(nq, m_fields[3]->GetPhys(), 1, fce, 1);
499}
500
502 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
504{
505 int i;
506 int nTraceNumPoints = GetTraceTotPoints();
507
508 //-----------------------------------------------------
509 // get temporary arrays
512
513 for (i = 0; i < 2; ++i)
514 {
515 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
516 Bwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
517 }
518 //-----------------------------------------------------
519
520 //-----------------------------------------------------
521 // get the physical values at the trace
522 // (any time-dependent BC previously put in fields[1] and [2]
523
524 m_fields[1]->GetFwdBwdTracePhys(inarray[0], Fwd[0], Bwd[0]);
525 m_fields[2]->GetFwdBwdTracePhys(inarray[1], Fwd[1], Bwd[1]);
526 //-----------------------------------------------------
527
528 //-----------------------------------------------------
529 // use centred fluxes for the numerical flux
530 for (i = 0; i < nTraceNumPoints; ++i)
531 {
532 numfluxX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
533 numfluxY[i] = 0.5 * (Fwd[1][i] + Bwd[1][i]);
534 }
535 //-----------------------------------------------------
536}
537
540 [[maybe_unused]] NekDouble time)
541{
542 int cnt = 0;
543
544 // loop over Boundary Regions
545 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
546 {
547 // Use wall for all BC...
548 // Wall Boundary Condition
549 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
550 "Wall"))
551 {
552 WallBoundaryForcing(n, cnt, inarray);
553 }
554
555 // Timedependent Boundary Condition
556 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
557 {
558 ASSERTL0(false, "time-dependent BC not implemented for Boussinesq");
559 }
560 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
561 }
562}
563
564// fills up boundary expansion for field[1] and [2]
566 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &inarray)
567{
568 int nTraceNumPoints = GetTraceTotPoints();
569 int nvariables = 2;
570
571 // get physical values of f1 and f2 for the forward trace
572 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
573 for (int i = 0; i < nvariables; ++i)
574 {
575 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
576 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
577 }
578
579 // Adjust the physical values of the trace to take
580 // user defined boundaries into account
581 int e, id1, id2, npts;
583 m_fields[0]->GetBndCondExpansions()[bcRegion];
584 for (e = 0; e < bcexp->GetExpSize(); ++e)
585 {
586 npts = bcexp->GetExp(e)->GetTotPoints();
587 id1 = bcexp->GetPhys_Offset(e);
588 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
589 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
590
591 switch (m_expdim)
592 {
593 case 1:
594 {
595 ASSERTL0(false, "1D not yet implemented for Boussinesq");
596 break;
597 }
598 case 2:
599 {
600 Array<OneD, NekDouble> tmp_n(npts);
601 Array<OneD, NekDouble> tmp_t(npts);
602
603 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[0][id2], 1,
604 &tmp_n[0], 1);
605 Vmath::Vvtvp(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
606 &tmp_n[0], 1, &tmp_n[0], 1);
607
608 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[1][id2], 1,
609 &tmp_t[0], 1);
610 Vmath::Vvtvm(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
611 &tmp_t[0], 1, &tmp_t[0], 1);
612
613 // negate the normal flux
614 Vmath::Neg(npts, tmp_n, 1);
615
616 // rotate back to Cartesian
617 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
618 &Fwd[0][id2], 1);
619 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
620 &Fwd[0][id2], 1, &Fwd[0][id2], 1);
621
622 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
623 &Fwd[1][id2], 1);
624 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
625 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
626 break;
627 }
628 case 3:
629 ASSERTL0(false, "3D not implemented for Boussinesq equations");
630 break;
631 default:
632 ASSERTL0(false, "Illegal expansion dimension");
633 }
634
635 // copy boundary adjusted values into the boundary expansion
636 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
637 Vmath::Vcopy(npts, &Fwd[0][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
638
639 bcexp = m_fields[2]->GetBndCondExpansions()[bcRegion];
640 Vmath::Vcopy(npts, &Fwd[1][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
641 }
642}
643
645 Array<OneD, NekDouble> &inarray, [[maybe_unused]] NekDouble time)
646{
647 int cnt = 0;
648
649 // loop over Boundary Regions
650 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
651 {
652 // Use wall for all
653 // Wall Boundary Condition
654 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
655 "Wall"))
656 {
657 WallBoundaryContVariables(n, cnt, inarray);
658 }
659
660 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
661 {
662 WallBoundaryContVariables(n, cnt, inarray);
663 }
664
665 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize() - 1;
666 }
667}
668
670 int bcRegion, int cnt, Array<OneD, NekDouble> &inarray)
671{
672 int nTraceNumPoints = GetTraceTotPoints();
673
674 // get physical values of z for the forward trace
675 Array<OneD, NekDouble> z(nTraceNumPoints);
676 m_fields[0]->ExtractTracePhys(inarray, z);
677
678 // Adjust the physical values of the trace to take
679 // user defined boundaries into account
680 int e, id1, id2, npts;
682 m_fields[0]->GetBndCondExpansions()[bcRegion];
683
684 for (e = 0; e < bcexp->GetExpSize(); ++e)
685 {
686 npts = bcexp->GetExp(e)->GetTotPoints();
687 id1 = bcexp->GetPhys_Offset(e);
688 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
689 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
690
691 // copy boundary adjusted values into the boundary expansion
692 // field[1] and field[2]
693 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
694 Vmath::Vcopy(npts, &z[id2], 1, &(bcexp->UpdatePhys())[id1], 1);
695 }
696}
697
701{
702 int i;
703 int nTraceNumPoints = GetTraceTotPoints();
704
705 //-----------------------------------------------------
706 // get temporary arrays
707 Array<OneD, NekDouble> Fwd(nTraceNumPoints);
708 Array<OneD, NekDouble> Bwd(nTraceNumPoints);
709 //-----------------------------------------------------
710
711 //-----------------------------------------------------
712 // get the physical values at the trace
713 // (we have put any time-dependent BC in field[1])
714
715 m_fields[1]->GetFwdBwdTracePhys(physfield, Fwd, Bwd);
716 //-----------------------------------------------------
717
718 //-----------------------------------------------------
719 // use centred fluxes for the numerical flux
720 for (i = 0; i < nTraceNumPoints; ++i)
721 {
722 outX[i] = 0.5 * (Fwd[i] + Bwd[i]);
723 outY[i] = 0.5 * (Fwd[i] + Bwd[i]);
724 }
725 //-----------------------------------------------------
726}
727
728} // 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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
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.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
void SetBoundaryConditionsForcing(Array< OneD, Array< OneD, NekDouble > > &inarray, 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)
static std::string className
Name of class.
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 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
problem type selector
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)
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
NekDouble m_g
Acceleration of gravity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
const Array< OneD, NekDouble > & GetDepth()
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
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_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 GetNcoeffs()
int m_nchk
Number of checkpoints written so far.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
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
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 Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.hpp:439
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294