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