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