Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // 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 Shallow water equations in conservative variables
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 #include <iomanip>
38 #include <boost/algorithm/string.hpp>
39 
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47  string NonlinearSWE::className =
49  "NonlinearSWE", NonlinearSWE::create,
50  "Nonlinear shallow water equation in conservative variables.");
51 
52  NonlinearSWE::NonlinearSWE(
54  : ShallowWaterSystem(pSession)
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);
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  // Wall Boundary Condition
442  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),"Wall"))
443  {
444  WallBoundary2D(n, cnt, Fwd, inarray);
445  }
446 
447  // Time Dependent Boundary Condition (specified in meshfile)
448  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
449  {
450  for (int i = 0; i < nvariables; ++i)
451  {
452  varName = m_session->GetVariable(i);
453  m_fields[i]->EvaluateBoundaryConditions(time, varName);
454  }
455  }
456  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
457  }
458  }
459 
460  //----------------------------------------------------
461  /**
462  * @brief Wall boundary condition.
463  */
465  int bcRegion,
466  int cnt,
468  Array<OneD, Array<OneD, NekDouble> > &physarray)
469  {
470  int i;
471  int nvariables = physarray.num_elements();
472 
473  // Adjust the physical values of the trace to take
474  // user defined boundaries into account
475  int e, id1, id2, npts;
476 
477  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]
478  ->GetExpSize(); ++e)
479  {
480  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->
481  GetExp(e)->GetTotPoints();
482  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->
483  GetPhys_Offset(e);
484  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
485  m_fields[0]->GetTraceMap()->
486  GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
487 
488  // For 2D/3D, define: v* = v - 2(v.n)n
489  Array<OneD, NekDouble> tmp(npts, 0.0);
490 
491  // Calculate (v.n)
492  for (i = 0; i < m_spacedim; ++i)
493  {
494  Vmath::Vvtvp(npts,
495  &Fwd[1+i][id2], 1,
496  &m_traceNormals[i][id2], 1,
497  &tmp[0], 1,
498  &tmp[0], 1);
499  }
500 
501  // Calculate 2.0(v.n)
502  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
503 
504  // Calculate v* = v - 2.0(v.n)n
505  for (i = 0; i < m_spacedim; ++i)
506  {
507  Vmath::Vvtvp(npts,
508  &tmp[0], 1,
509  &m_traceNormals[i][id2], 1,
510  &Fwd[1+i][id2], 1,
511  &Fwd[1+i][id2], 1);
512  }
513 
514  // copy boundary adjusted values into the boundary expansion
515  for (i = 0; i < nvariables; ++i)
516  {
517  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
518  &(m_fields[i]->GetBndCondExpansions()[bcRegion]->
519  UpdatePhys())[id1], 1);
520  }
521  }
522  }
523 
524 
526  {
527 
528  int i;
529  int nvariables = physarray.num_elements();
530 
531  // Adjust the physical values of the trace to take
532  // user defined boundaries into account
533  int e, id1, id2, npts;
534 
535  for(e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize(); ++e)
536  {
537  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExp(e)->GetNumPoints(0);
538  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e) ;
539  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
540 
541  switch(m_expdim)
542  {
543  case 1:
544  {
545  // negate the forward flux
546  Vmath::Neg(npts,&Fwd[1][id2],1);
547  }
548  break;
549  case 2:
550  {
551  Array<OneD, NekDouble> tmp_n(npts);
552  Array<OneD, NekDouble> tmp_t(npts);
553 
554  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[0][id2],1,&tmp_n[0],1);
555  Vmath::Vvtvp(npts,&Fwd[2][id2],1,&m_traceNormals[1][id2],1,&tmp_n[0],1,&tmp_n[0],1);
556 
557  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[1][id2],1,&tmp_t[0],1);
558  Vmath::Vvtvm(npts,&Fwd[2][id2],1,&m_traceNormals[0][id2],1,&tmp_t[0],1,&tmp_t[0],1);
559 
560  // negate the normal flux
561  Vmath::Neg(npts,tmp_n,1);
562 
563  // rotate back to Cartesian
564  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[1][id2],1,&Fwd[1][id2],1);
565  Vmath::Vvtvm(npts,&tmp_n[0],1,&m_traceNormals[0][id2],1,&Fwd[1][id2],1,&Fwd[1][id2],1);
566 
567  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[0][id2],1,&Fwd[2][id2],1);
568  Vmath::Vvtvp(npts,&tmp_n[0],1,&m_traceNormals[1][id2],1,&Fwd[2][id2],1,&Fwd[2][id2],1);
569  }
570  break;
571  case 3:
572  ASSERTL0(false,"3D not implemented for Shallow Water Equations");
573  break;
574  default:
575  ASSERTL0(false,"Illegal expansion dimension");
576  }
577 
578 
579 
580  // copy boundary adjusted values into the boundary expansion
581  for (i = 0; i < nvariables; ++i)
582  {
583  Vmath::Vcopy(npts,&Fwd[i][id2], 1,&(m_fields[i]->GetBndCondExpansions()[bcRegion]->UpdatePhys())[id1],1);
584  }
585  }
586  }
587 
588 
589  // Physfield in conservative Form
591  const Array<OneD, const Array<OneD, NekDouble> > &physfield,
593  {
594  int i, j;
595  int nq = m_fields[0]->GetTotPoints();
596 
597  NekDouble g = m_g;
599 
600  // Flux vector for the mass equation
601  for (i = 0; i < m_spacedim; ++i)
602  {
603  velocity[i] = Array<OneD, NekDouble>(nq);
604  Vmath::Vcopy(nq, physfield[i+1], 1, flux[0][i], 1);
605  }
606 
607  GetVelocityVector(physfield, velocity);
608 
609  // Put (0.5 g h h) in tmp
610  Array<OneD, NekDouble> tmp(nq);
611  Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
612  Vmath::Smul(nq, 0.5*g, tmp, 1, tmp, 1);
613 
614  // Flux vector for the momentum equations
615  for (i = 0; i < m_spacedim; ++i)
616  {
617  for (j = 0; j < m_spacedim; ++j)
618  {
619  Vmath::Vmul(nq, velocity[j], 1, physfield[i+1], 1,
620  flux[i+1][j], 1);
621  }
622 
623  // Add (0.5 g h h) to appropriate field
624  Vmath::Vadd(nq, flux[i+1][i], 1, tmp, 1, flux[i+1][i], 1);
625  }
626 
627  }
628 
630  Array<OneD, Array<OneD, NekDouble> >&physout)
631  {
632  int nq = GetTotPoints();
633 
634  if(physin.get() == physout.get())
635  {
636  // copy indata and work with tmp array
638  for (int i = 0; i < 3; ++i)
639  {
640  // deep copy
641  tmp[i] = Array<OneD, NekDouble>(nq);
642  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
643  }
644 
645  // \eta = h - d
646  Vmath::Vsub(nq,tmp[0],1,m_depth,1,physout[0],1);
647 
648  // u = hu/h
649  Vmath::Vdiv(nq,tmp[1],1,tmp[0],1,physout[1],1);
650 
651  // v = hv/ v
652  Vmath::Vdiv(nq,tmp[2],1,tmp[0],1,physout[2],1);
653  }
654  else
655  {
656  // \eta = h - d
657  Vmath::Vsub(nq,physin[0],1,m_depth,1,physout[0],1);
658 
659  // u = hu/h
660  Vmath::Vdiv(nq,physin[1],1,physin[0],1,physout[1],1);
661 
662  // v = hv/ v
663  Vmath::Vdiv(nq,physin[2],1,physin[0],1,physout[2],1);
664  }
665  }
666 
667 
669  {
670  int nq = GetTotPoints();
671 
672  // u = hu/h
673  Vmath::Vdiv(nq,m_fields[1]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
674 
675  // v = hv/ v
676  Vmath::Vdiv(nq,m_fields[2]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
677 
678  // \eta = h - d
679  Vmath::Vsub(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
680  }
681 
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  }
709  else
710  {
711  // h = \eta + d
712  Vmath::Vadd(nq,physin[0],1,m_depth,1,physout[0],1);
713 
714  // hu = h * u
715  Vmath::Vmul(nq,physout[0],1,physin[1],1,physout[1],1);
716 
717  // hv = h * v
718  Vmath::Vmul(nq,physout[0],1,physin[2],1,physout[2],1);
719 
720  }
721 
722  }
723 
725  {
726  int nq = GetTotPoints();
727 
728  // h = \eta + d
729  Vmath::Vadd(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
730 
731  // hu = h * u
732  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[1]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
733 
734  // hv = h * v
735  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[2]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
736  }
737 
738 
739  /**
740  * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
741  * \f$ h\mathbf{v} \f$.
742  *
743  * @param physfield Momentum field.
744  * @param velocity Velocity field.
745  */
747  const Array<OneD, Array<OneD, NekDouble> > &physfield,
748  Array<OneD, Array<OneD, NekDouble> > &velocity)
749  {
750  const int npts = physfield[0].num_elements();
751 
752  for (int i = 0; i < m_spacedim; ++i)
753  {
754  Vmath::Vdiv(npts, physfield[1+i], 1, physfield[0], 1,
755  velocity[i], 1);
756  }
757  }
758 
759 
761  {
763  SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
764  SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
765  SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
766  }
767 
768 } //end of namespace
769 
Array< OneD, NekDouble > m_coriolis
Coriolis force.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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.
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 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:47
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:442
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:241
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
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
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:213
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:50
RiemannSolverFactory & GetRiemannSolverFactory()
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 Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
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:343
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
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:465
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)
NekDouble m_g
Acceleration of gravity.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
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:299
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:183
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215