Nektar++
NonlinearSWE.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: NonlinearSWE.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 Shallow water equations in conservative variables
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/algorithm/string.hpp>
36 #include <iomanip>
37 #include <iostream>
38 
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 string NonlinearSWE::className =
48  "NonlinearSWE", NonlinearSWE::create,
49  "Nonlinear shallow water equation in conservative variables.");
50 
51 NonlinearSWE::NonlinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession,
53  : ShallowWaterSystem(pSession, pGraph)
54 {
55 }
56 
57 void NonlinearSWE::v_InitObject(bool DeclareFields)
58 {
59  ShallowWaterSystem::v_InitObject(DeclareFields);
60 
62  {
65  }
66  else
67  {
68  ASSERTL0(false, "Implicit SWE not set up.");
69  }
70 
71  // Type of advection class to be used
72  switch (m_projectionType)
73  {
74  // Continuous field
76  {
77  // Do nothing
78  break;
79  }
80  // Discontinuous field
82  {
83  string advName;
84  string diffName;
85  string riemName;
86 
87  //---------------------------------------------------------------
88  // Setting up advection and diffusion operators
89  // NB: diffusion not set up for SWE at the moment
90  // but kept here for future use ...
91  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
92  // m_session->LoadSolverInfo("DiffusionType", diffName, "LDGEddy");
94  advName, advName);
95  // m_diffusion = SolverUtils::GetDiffusionFactory()
96  // .CreateInstance(diffName, diffName);
97 
98  m_advection->SetFluxVector(&NonlinearSWE::GetFluxVector, this);
99  // m_diffusion->SetFluxVectorNS(&ShallowWaterSystem::
100  // GetEddyViscosityFluxVector,
101  // this);
102 
103  // Setting up Riemann solver for advection operator
104  m_session->LoadSolverInfo("UpwindType", riemName, "Average");
107  riemName, m_session);
108 
109  // Setting up upwind solver for diffusion operator
110  // m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
111  // .CreateInstance("UpwindLDG");
112 
113  // Setting up parameters for advection operator Riemann solver
114  m_riemannSolver->SetParam("gravity", &NonlinearSWE::GetGravity,
115  this);
116  m_riemannSolver->SetAuxVec("vecLocs", &NonlinearSWE::GetVecLocs,
117  this);
118  m_riemannSolver->SetVector("N", &NonlinearSWE::GetNormals, this);
119  m_riemannSolver->SetScalar("depth", &NonlinearSWE::GetDepth, this);
120 
121  // Setting up parameters for diffusion operator Riemann solver
122  // m_riemannSolverLDG->AddParam (
123  // "gravity",
124  // &NonlinearSWE::GetGravity, this);
125  // m_riemannSolverLDG->SetAuxVec(
126  // "vecLocs",
127  // &NonlinearSWE::GetVecLocs, this);
128  // m_riemannSolverLDG->AddVector(
129  // "N",
130  // &NonlinearSWE::GetNormals, this);
131 
132  // Concluding initialisation of advection / diffusion operators
133  m_advection->SetRiemannSolver(m_riemannSolver);
134  // m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
135  m_advection->InitObject(m_session, m_fields);
136  // m_diffusion->InitObject (m_session, m_fields);
137  break;
138  }
139  default:
140  {
141  ASSERTL0(false, "Unsupported projection type.");
142  break;
143  }
144  }
145 }
146 
148 {
149 }
150 
151 // physarray contains the conservative variables
153  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
154  Array<OneD, Array<OneD, NekDouble>> &outarray)
155 {
156 
157  int ncoeffs = GetNcoeffs();
158  int nq = GetTotPoints();
159 
160  Array<OneD, NekDouble> tmp(nq);
161  Array<OneD, NekDouble> mod(ncoeffs);
162 
163  switch (m_projectionType)
164  {
166  {
167  // add to hu equation
168  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
169  m_fields[0]->IProductWRTBase(tmp, mod);
170  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
171  m_fields[0]->BwdTrans(mod, tmp);
172  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
173 
174  // add to hv equation
175  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
176  Vmath::Neg(nq, tmp, 1);
177  m_fields[0]->IProductWRTBase(tmp, mod);
178  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
179  m_fields[0]->BwdTrans(mod, tmp);
180  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
181  }
182  break;
185  {
186  // add to hu equation
187  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
188  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
189 
190  // add to hv equation
191  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
192  Vmath::Neg(nq, tmp, 1);
193  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
194  }
195  break;
196  default:
197  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
198  break;
199  }
200 }
201 
202 // physarray contains the conservative variables
204  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
205  Array<OneD, Array<OneD, NekDouble>> &outarray)
206 {
207 
208  int ncoeffs = GetNcoeffs();
209  int nq = GetTotPoints();
210 
211  Array<OneD, NekDouble> tmp(nq);
212  Array<OneD, NekDouble> mod(ncoeffs);
213 
214  switch (m_projectionType)
215  {
217  {
218  for (int i = 0; i < m_spacedim; ++i)
219  {
220  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
221  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
222  m_fields[0]->IProductWRTBase(tmp, mod);
223  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
224  m_fields[0]->BwdTrans(mod, tmp);
225  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
226  }
227  }
228  break;
231  {
232  for (int i = 0; i < m_spacedim; ++i)
233  {
234  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
235  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
236  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
237  }
238  }
239  break;
240  default:
241  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
242  break;
243  }
244 }
245 
247  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
248  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
249 {
250  int i, j;
251  int ndim = m_spacedim;
252  int nvariables = inarray.size();
253  int nq = GetTotPoints();
254 
255  switch (m_projectionType)
256  {
258  {
259 
260  //-------------------------------------------------------
261  // Compute the DG advection including the numerical flux
262  // by using SolverUtils/Advection
263  // Input and output in physical space
265 
266  m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
267  time);
268  //-------------------------------------------------------
269 
270  //-------------------------------------------------------
271  // negate the outarray since moving terms to the rhs
272  for (i = 0; i < nvariables; ++i)
273  {
274  Vmath::Neg(nq, outarray[i], 1);
275  }
276  //-------------------------------------------------------
277 
278  //-------------------------------------------------
279  // Add "source terms"
280  // Input and output in physical space
281 
282  // Coriolis forcing
283  if (m_coriolis.size() != 0)
284  {
285  AddCoriolis(inarray, outarray);
286  }
287 
288  // Variable Depth
289  if (m_constantDepth != true)
290  {
291  AddVariableDepth(inarray, outarray);
292  }
293  //-------------------------------------------------
294  }
295  break;
298  {
299 
300  //-------------------------------------------------------
301  // Compute the fluxvector in physical space
303  nvariables);
304 
305  for (i = 0; i < nvariables; ++i)
306  {
307  fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
308  for (j = 0; j < ndim; ++j)
309  {
310  fluxvector[i][j] = Array<OneD, NekDouble>(nq);
311  }
312  }
313 
314  NonlinearSWE::GetFluxVector(inarray, fluxvector);
315  //-------------------------------------------------------
316 
317  //-------------------------------------------------------
318  // Take the derivative of the flux terms
319  // and negate the outarray since moving terms to the rhs
320  Array<OneD, NekDouble> tmp(nq);
321  Array<OneD, NekDouble> tmp1(nq);
322 
323  for (i = 0; i < nvariables; ++i)
324  {
325  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],
326  fluxvector[i][0], tmp);
327  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],
328  fluxvector[i][1], tmp1);
329  Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
330  Vmath::Neg(nq, outarray[i], 1);
331  }
332 
333  //-------------------------------------------------
334  // Add "source terms"
335  // Input and output in physical space
336 
337  // Coriolis forcing
338  if (m_coriolis.size() != 0)
339  {
340  AddCoriolis(inarray, outarray);
341  }
342 
343  // Variable Depth
344  if (m_constantDepth != true)
345  {
346  AddVariableDepth(inarray, outarray);
347  }
348  //-------------------------------------------------
349  }
350  break;
351  default:
352  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
353  break;
354  }
355 }
356 
358  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
359  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
360 {
361  int i;
362  int nvariables = inarray.size();
363 
364  switch (m_projectionType)
365  {
367  {
368 
369  // Just copy over array
370  int npoints = GetNpoints();
371 
372  for (i = 0; i < nvariables; ++i)
373  {
374  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
375  }
376  SetBoundaryConditions(outarray, time);
377  break;
378  }
381  {
382 
384  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
385 
386  for (i = 0; i < nvariables; ++i)
387  {
388  m_fields[i]->FwdTrans(inarray[i], coeffs);
389  m_fields[i]->BwdTrans(coeffs, outarray[i]);
390  }
391  break;
392  }
393  default:
394  ASSERTL0(false, "Unknown projection scheme");
395  break;
396  }
397 }
398 
399 //----------------------------------------------------
401  Array<OneD, Array<OneD, NekDouble>> &inarray, NekDouble time)
402 {
403  std::string varName;
404  int nvariables = m_fields.size();
405  int cnt = 0;
406  int nTracePts = GetTraceTotPoints();
407 
408  // Extract trace for boundaries. Needs to be done on all processors to avoid
409  // deadlock.
410  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
411  for (int i = 0; i < nvariables; ++i)
412  {
413  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
414  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
415  }
416 
417  // Loop over Boundary Regions
418  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
419  {
420 
421  if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
423  {
424  continue;
425  }
426 
427  // Wall Boundary Condition
428  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
429  "Wall"))
430  {
431  WallBoundary2D(n, cnt, Fwd, inarray);
432  }
433 
434  // Time Dependent Boundary Condition (specified in meshfile)
435  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
436  {
437  for (int i = 0; i < nvariables; ++i)
438  {
439  varName = m_session->GetVariable(i);
440  m_fields[i]->EvaluateBoundaryConditions(time, varName);
441  }
442  }
443  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
444  }
445 }
446 
447 //----------------------------------------------------
448 /**
449  * @brief Wall boundary condition.
450  */
451 void NonlinearSWE::WallBoundary(int bcRegion, int cnt,
453  Array<OneD, Array<OneD, NekDouble>> &physarray)
454 {
455  int i;
456  int nvariables = physarray.size();
457 
458  // Adjust the physical values of the trace to take
459  // user defined boundaries into account
460  int e, id1, id2, npts;
461 
462  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
463  ++e)
464  {
465  npts = m_fields[0]
466  ->GetBndCondExpansions()[bcRegion]
467  ->GetExp(e)
468  ->GetTotPoints();
469  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
470  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
471  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
472 
473  // For 2D/3D, define: v* = v - 2(v.n)n
474  Array<OneD, NekDouble> tmp(npts, 0.0);
475 
476  // Calculate (v.n)
477  for (i = 0; i < m_spacedim; ++i)
478  {
479  Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
480  &tmp[0], 1, &tmp[0], 1);
481  }
482 
483  // Calculate 2.0(v.n)
484  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
485 
486  // Calculate v* = v - 2.0(v.n)n
487  for (i = 0; i < m_spacedim; ++i)
488  {
489  Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
490  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
491  }
492 
493  // copy boundary adjusted values into the boundary expansion
494  for (i = 0; i < nvariables; ++i)
495  {
496  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
497  &(m_fields[i]
498  ->GetBndCondExpansions()[bcRegion]
499  ->UpdatePhys())[id1],
500  1);
501  }
502  }
503 }
504 
506  int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd,
507  Array<OneD, Array<OneD, NekDouble>> &physarray)
508 {
509 
510  int i;
511  int nvariables = physarray.size();
512 
513  // Adjust the physical values of the trace to take
514  // user defined boundaries into account
515  int e, id1, id2, npts;
516 
517  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
518  ++e)
519  {
520  npts = m_fields[0]
521  ->GetBndCondExpansions()[bcRegion]
522  ->GetExp(e)
523  ->GetNumPoints(0);
524  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
525  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
526  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
527 
528  switch (m_expdim)
529  {
530  case 1:
531  {
532  // negate the forward flux
533  Vmath::Neg(npts, &Fwd[1][id2], 1);
534  }
535  break;
536  case 2:
537  {
538  Array<OneD, NekDouble> tmp_n(npts);
539  Array<OneD, NekDouble> tmp_t(npts);
540 
541  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
542  &tmp_n[0], 1);
543  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
544  &tmp_n[0], 1, &tmp_n[0], 1);
545 
546  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
547  &tmp_t[0], 1);
548  Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
549  &tmp_t[0], 1, &tmp_t[0], 1);
550 
551  // negate the normal flux
552  Vmath::Neg(npts, tmp_n, 1);
553 
554  // rotate back to Cartesian
555  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
556  &Fwd[1][id2], 1);
557  Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
558  &Fwd[1][id2], 1, &Fwd[1][id2], 1);
559 
560  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
561  &Fwd[2][id2], 1);
562  Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
563  &Fwd[2][id2], 1, &Fwd[2][id2], 1);
564  }
565  break;
566  case 3:
567  ASSERTL0(false,
568  "3D not implemented for Shallow Water Equations");
569  break;
570  default:
571  ASSERTL0(false, "Illegal expansion dimension");
572  }
573 
574  // copy boundary adjusted values into the boundary expansion
575  for (i = 0; i < nvariables; ++i)
576  {
577  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
578  &(m_fields[i]
579  ->GetBndCondExpansions()[bcRegion]
580  ->UpdatePhys())[id1],
581  1);
582  }
583  }
584 }
585 
586 // Physfield in conservative Form
588  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
590 {
591  int i, j;
592  int nq = m_fields[0]->GetTotPoints();
593 
594  NekDouble g = m_g;
596 
597  // Flux vector for the mass equation
598  for (i = 0; i < m_spacedim; ++i)
599  {
600  velocity[i] = Array<OneD, NekDouble>(nq);
601  Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
602  }
603 
604  GetVelocityVector(physfield, velocity);
605 
606  // Put (0.5 g h h) in tmp
607  Array<OneD, NekDouble> tmp(nq);
608  Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
609  Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
610 
611  // Flux vector for the momentum equations
612  for (i = 0; i < m_spacedim; ++i)
613  {
614  for (j = 0; j < m_spacedim; ++j)
615  {
616  Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
617  1);
618  }
619 
620  // Add (0.5 g h h) to appropriate field
621  Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
622  }
623 }
624 
626  const Array<OneD, const Array<OneD, NekDouble>> &physin,
627  Array<OneD, Array<OneD, NekDouble>> &physout)
628 {
629  int nq = GetTotPoints();
630 
631  if (physin.get() == physout.get())
632  {
633  // copy indata and work with tmp array
635  for (int i = 0; i < 3; ++i)
636  {
637  // deep copy
638  tmp[i] = Array<OneD, NekDouble>(nq);
639  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
640  }
641 
642  // \eta = h - d
643  Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
644 
645  // u = hu/h
646  Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
647 
648  // v = hv/ v
649  Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
650  }
651  else
652  {
653  // \eta = h - d
654  Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
655 
656  // u = hu/h
657  Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
658 
659  // v = hv/ v
660  Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
661  }
662 }
663 
665 {
666  int nq = GetTotPoints();
667 
668  // u = hu/h
669  Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
670  m_fields[1]->UpdatePhys(), 1);
671 
672  // v = hv/ v
673  Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
674  m_fields[2]->UpdatePhys(), 1);
675 
676  // \eta = h - d
677  Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
678  m_fields[0]->UpdatePhys(), 1);
679 }
680 
682  const Array<OneD, const Array<OneD, NekDouble>> &physin,
683  Array<OneD, Array<OneD, NekDouble>> &physout)
684 {
685 
686  int nq = GetTotPoints();
687 
688  if (physin.get() == physout.get())
689  {
690  // copy indata and work with tmp array
692  for (int i = 0; i < 3; ++i)
693  {
694  // deep copy
695  tmp[i] = Array<OneD, NekDouble>(nq);
696  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
697  }
698 
699  // h = \eta + d
700  Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
701 
702  // hu = h * u
703  Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
704 
705  // hv = h * v
706  Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
707  }
708  else
709  {
710  // h = \eta + d
711  Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
712 
713  // hu = h * u
714  Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
715 
716  // hv = h * v
717  Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
718  }
719 }
720 
722 {
723  int nq = GetTotPoints();
724 
725  // h = \eta + d
726  Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
727  m_fields[0]->UpdatePhys(), 1);
728 
729  // hu = h * u
730  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
731  m_fields[1]->UpdatePhys(), 1);
732 
733  // hv = h * v
734  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
735  m_fields[2]->UpdatePhys(), 1);
736 }
737 
738 /**
739  * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
740  * \f$ h\mathbf{v} \f$.
741  *
742  * @param physfield Momentum field.
743  * @param velocity Velocity field.
744  */
746  const Array<OneD, Array<OneD, NekDouble>> &physfield,
747  Array<OneD, Array<OneD, NekDouble>> &velocity)
748 {
749  const int npts = physfield[0].size();
750 
751  for (int i = 0; i < m_spacedim; ++i)
752  {
753  Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
754  }
755 }
756 
758 {
760  SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
761  SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
762  SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
763 }
764 
765 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
virtual void v_PrimitiveToConservative() override
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)
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
virtual void v_ConservativeToPrimitive() override
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble >> &physarray, NekDouble time)
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
Wall boundary condition.
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
Compute the velocity field given the momentum .
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
const Array< OneD, NekDouble > & GetDepth()
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, NekDouble > m_depth
Still water depth.
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
int m_spacedim
Spatial dimension (>= expansion dim).
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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:574
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:359
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.cpp:598
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
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:284
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
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:419