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