Nektar++
LinearSWE.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: LinearSWE.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: Linear Shallow water equations in primitive 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 LinearSWE::className =
48  "LinearSWE", LinearSWE::create,
49  "Linear shallow water equation in primitive variables.");
50 
51 LinearSWE::LinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession,
53  : ShallowWaterSystem(pSession, pGraph)
54 {
55 }
56 
57 void LinearSWE::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(&LinearSWE::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, "NoSolver");
105  if ((riemName == "LinearHLL") && (m_constantDepth != true))
106  {
107  ASSERTL0(false, "LinearHLL only valid for constant depth");
108  }
111  riemName, m_session);
112 
113  // Setting up upwind solver for diffusion operator
114  // m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
115  // .CreateInstance("UpwindLDG");
116 
117  // Setting up parameters for advection operator Riemann solver
118  m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
119  m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
120  m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
121 
122  // The numerical flux for linear SWE requires depth information
123  int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
124  m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
125  m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
126  m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
128  m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
129  this);
130  m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
131  this);
132 
133  // Setting up parameters for diffusion operator Riemann solver
134  // m_riemannSolverLDG->AddParam (
135  // "gravity",
136  // &NonlinearSWE::GetGravity, this);
137  // m_riemannSolverLDG->SetAuxVec(
138  // "vecLocs",
139  // &NonlinearSWE::GetVecLocs, this);
140  // m_riemannSolverLDG->AddVector(
141  // "N",
142  // &NonlinearSWE::GetNormals, this);
143 
144  // Concluding initialisation of advection / diffusion operators
145  m_advection->SetRiemannSolver(m_riemannSolver);
146  // m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
147  m_advection->InitObject(m_session, m_fields);
148  // m_diffusion->InitObject (m_session, m_fields);
149  break;
150  }
151  default:
152  {
153  ASSERTL0(false, "Unsupported projection type.");
154  break;
155  }
156  }
157 }
158 
160 {
161 }
162 
163 // physarray contains the conservative variables
165  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
166  Array<OneD, Array<OneD, NekDouble>> &outarray)
167 {
168 
169  int ncoeffs = GetNcoeffs();
170  int nq = GetTotPoints();
171 
172  Array<OneD, NekDouble> tmp(nq);
173  Array<OneD, NekDouble> mod(ncoeffs);
174 
175  switch (m_projectionType)
176  {
178  {
179  // add to u equation
180  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
181  m_fields[0]->IProductWRTBase(tmp, mod);
182  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
183  m_fields[0]->BwdTrans(mod, tmp);
184  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
185 
186  // add to v equation
187  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
188  Vmath::Neg(nq, tmp, 1);
189  m_fields[0]->IProductWRTBase(tmp, mod);
190  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
191  m_fields[0]->BwdTrans(mod, tmp);
192  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
193  }
194  break;
197  {
198  // add to u equation
199  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
200  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
201 
202  // add to v equation
203  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
204  Vmath::Neg(nq, tmp, 1);
205  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
206  }
207  break;
208  default:
209  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
210  break;
211  }
212 }
213 
215  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
216  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
217 {
218  int i, j;
219  int ndim = m_spacedim;
220  int nvariables = inarray.size();
221  int nq = GetTotPoints();
222 
223  switch (m_projectionType)
224  {
226  {
227 
228  //-------------------------------------------------------
229  // Compute the DG advection including the numerical flux
230  // by using SolverUtils/Advection
231  // Input and output in physical space
233 
234  m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
235  time);
236  //-------------------------------------------------------
237 
238  //-------------------------------------------------------
239  // negate the outarray since moving terms to the rhs
240  for (i = 0; i < nvariables; ++i)
241  {
242  Vmath::Neg(nq, outarray[i], 1);
243  }
244  //-------------------------------------------------------
245 
246  //-------------------------------------------------
247  // Add "source terms"
248  // Input and output in physical space
249 
250  // Coriolis forcing
251  if (m_coriolis.size() != 0)
252  {
253  AddCoriolis(inarray, outarray);
254  }
255  //-------------------------------------------------
256  }
257  break;
260  {
261 
262  //-------------------------------------------------------
263  // Compute the fluxvector in physical space
265  nvariables);
266 
267  for (i = 0; i < nvariables; ++i)
268  {
269  fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
270  for (j = 0; j < ndim; ++j)
271  {
272  fluxvector[i][j] = Array<OneD, NekDouble>(nq);
273  }
274  }
275 
276  LinearSWE::GetFluxVector(inarray, fluxvector);
277  //-------------------------------------------------------
278 
279  //-------------------------------------------------------
280  // Take the derivative of the flux terms
281  // and negate the outarray since moving terms to the rhs
282  Array<OneD, NekDouble> tmp(nq);
283  Array<OneD, NekDouble> tmp1(nq);
284 
285  for (i = 0; i < nvariables; ++i)
286  {
287  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],
288  fluxvector[i][0], tmp);
289  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],
290  fluxvector[i][1], tmp1);
291  Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
292  Vmath::Neg(nq, outarray[i], 1);
293  }
294 
295  //-------------------------------------------------
296  // Add "source terms"
297  // Input and output in physical space
298 
299  // Coriolis forcing
300  if (m_coriolis.size() != 0)
301  {
302  AddCoriolis(inarray, outarray);
303  }
304  //-------------------------------------------------
305  }
306  break;
307  default:
308  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
309  break;
310  }
311 }
312 
314  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
315  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
316 {
317  int i;
318  int nvariables = inarray.size();
319 
320  switch (m_projectionType)
321  {
323  {
324 
325  // Just copy over array
326  int npoints = GetNpoints();
327 
328  for (i = 0; i < nvariables; ++i)
329  {
330  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
331  }
332  SetBoundaryConditions(outarray, time);
333  break;
334  }
337  {
338 
340  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
341 
342  for (i = 0; i < nvariables; ++i)
343  {
344  m_fields[i]->FwdTrans(inarray[i], coeffs);
345  m_fields[i]->BwdTrans(coeffs, outarray[i]);
346  }
347  break;
348  }
349  default:
350  ASSERTL0(false, "Unknown projection scheme");
351  break;
352  }
353 }
354 
355 //----------------------------------------------------
357  Array<OneD, Array<OneD, NekDouble>> &inarray, NekDouble time)
358 {
359  std::string varName;
360  int nvariables = m_fields.size();
361  int cnt = 0;
362  int nTracePts = GetTraceTotPoints();
363 
364  // Extract trace for boundaries. Needs to be done on all processors to avoid
365  // deadlock.
366  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
367  for (int i = 0; i < nvariables; ++i)
368  {
369  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
370  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
371  }
372 
373  // loop over Boundary Regions
374  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
375  {
376  if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
378  {
379  continue;
380  }
381 
382  // Wall Boundary Condition
383  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
384  "Wall"))
385  {
386  WallBoundary2D(n, cnt, Fwd, inarray);
387  }
388 
389  // Time Dependent Boundary Condition (specified in meshfile)
390  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
391  {
392  for (int i = 0; i < nvariables; ++i)
393  {
394  varName = m_session->GetVariable(i);
395  m_fields[i]->EvaluateBoundaryConditions(time, varName);
396  }
397  }
398  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
399  }
400 }
401 
402 //----------------------------------------------------
403 /**
404  * @brief Wall boundary condition.
405  */
406 void LinearSWE::WallBoundary(int bcRegion, int cnt,
408  Array<OneD, Array<OneD, NekDouble>> &physarray)
409 {
410  int i;
411  int nvariables = physarray.size();
412 
413  // Adjust the physical values of the trace to take
414  // user defined boundaries into account
415  int e, id1, id2, npts;
416 
417  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
418  ++e)
419  {
420  npts = m_fields[0]
421  ->GetBndCondExpansions()[bcRegion]
422  ->GetExp(e)
423  ->GetTotPoints();
424  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
425  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
426  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
427 
428  // For 2D/3D, define: v* = v - 2(v.n)n
429  Array<OneD, NekDouble> tmp(npts, 0.0);
430 
431  // Calculate (v.n)
432  for (i = 0; i < m_spacedim; ++i)
433  {
434  Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
435  &tmp[0], 1, &tmp[0], 1);
436  }
437 
438  // Calculate 2.0(v.n)
439  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
440 
441  // Calculate v* = v - 2.0(v.n)n
442  for (i = 0; i < m_spacedim; ++i)
443  {
444  Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
445  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
446  }
447 
448  // copy boundary adjusted values into the boundary expansion
449  for (i = 0; i < nvariables; ++i)
450  {
451  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
452  &(m_fields[i]
453  ->GetBndCondExpansions()[bcRegion]
454  ->UpdatePhys())[id1],
455  1);
456  }
457  }
458 }
459 
460 void LinearSWE::WallBoundary2D(int bcRegion, int cnt,
462  Array<OneD, Array<OneD, NekDouble>> &physarray)
463 {
464 
465  int i;
466  int nvariables = physarray.size();
467 
468  // Adjust the physical values of the trace to take
469  // user defined boundaries into account
470  int e, id1, id2, npts;
471 
472  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
473  ++e)
474  {
475  npts = m_fields[0]
476  ->GetBndCondExpansions()[bcRegion]
477  ->GetExp(e)
478  ->GetNumPoints(0);
479  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
480  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
481  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
482 
483  switch (m_expdim)
484  {
485  case 1:
486  {
487  // negate the forward flux
488  Vmath::Neg(npts, &Fwd[1][id2], 1);
489  }
490  break;
491  case 2:
492  {
493  Array<OneD, NekDouble> tmp_n(npts);
494  Array<OneD, NekDouble> tmp_t(npts);
495 
496  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
497  &tmp_n[0], 1);
498  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
499  &tmp_n[0], 1, &tmp_n[0], 1);
500 
501  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
502  &tmp_t[0], 1);
503  Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
504  &tmp_t[0], 1, &tmp_t[0], 1);
505 
506  // negate the normal flux
507  Vmath::Neg(npts, tmp_n, 1);
508 
509  // rotate back to Cartesian
510  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
511  &Fwd[1][id2], 1);
512  Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
513  &Fwd[1][id2], 1, &Fwd[1][id2], 1);
514 
515  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
516  &Fwd[2][id2], 1);
517  Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
518  &Fwd[2][id2], 1, &Fwd[2][id2], 1);
519  }
520  break;
521  case 3:
522  ASSERTL0(false,
523  "3D not implemented for Shallow Water Equations");
524  break;
525  default:
526  ASSERTL0(false, "Illegal expansion dimension");
527  }
528 
529  // copy boundary adjusted values into the boundary expansion
530  for (i = 0; i < nvariables; ++i)
531  {
532  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
533  &(m_fields[i]
534  ->GetBndCondExpansions()[bcRegion]
535  ->UpdatePhys())[id1],
536  1);
537  }
538  }
539 }
540 
541 // Physfield in primitive Form
543  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
545 {
546  int i, j;
547  int nq = m_fields[0]->GetTotPoints();
548 
549  NekDouble g = m_g;
550 
551  // Flux vector for the mass equation
552  for (i = 0; i < m_spacedim; ++i)
553  {
554  Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
555  }
556 
557  // Put (g eta) in tmp
558  Array<OneD, NekDouble> tmp(nq);
559  Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
560 
561  // Flux vector for the momentum equations
562  for (i = 0; i < m_spacedim; ++i)
563  {
564  for (j = 0; j < m_spacedim; ++j)
565  {
566  // must zero fluxes as not initialised to zero in AdvectionWeakDG
567  // ...
568  Vmath::Zero(nq, flux[i + 1][j], 1);
569  }
570 
571  // Add (g eta) to appropriate field
572  Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
573  }
574 }
575 
577  const Array<OneD, const Array<OneD, NekDouble>> &physin,
578  Array<OneD, Array<OneD, NekDouble>> &physout)
579 {
580  int nq = GetTotPoints();
581 
582  if (physin.get() == physout.get())
583  {
584  // copy indata and work with tmp array
586  for (int i = 0; i < 3; ++i)
587  {
588  // deep copy
589  tmp[i] = Array<OneD, NekDouble>(nq);
590  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
591  }
592 
593  // \eta = h - d
594  Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
595 
596  // u = hu/h
597  Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
598 
599  // v = hv/ v
600  Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
601  }
602  else
603  {
604  // \eta = h - d
605  Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
606 
607  // u = hu/h
608  Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
609 
610  // v = hv/ v
611  Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
612  }
613 }
614 
616 {
617  int nq = GetTotPoints();
618 
619  // u = hu/h
620  Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
621  m_fields[1]->UpdatePhys(), 1);
622 
623  // v = hv/ v
624  Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
625  m_fields[2]->UpdatePhys(), 1);
626 
627  // \eta = h - d
628  Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
629  m_fields[0]->UpdatePhys(), 1);
630 }
631 
633  const Array<OneD, const Array<OneD, NekDouble>> &physin,
634  Array<OneD, Array<OneD, NekDouble>> &physout)
635 {
636 
637  int nq = GetTotPoints();
638 
639  if (physin.get() == physout.get())
640  {
641  // copy indata and work with tmp array
643  for (int i = 0; i < 3; ++i)
644  {
645  // deep copy
646  tmp[i] = Array<OneD, NekDouble>(nq);
647  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
648  }
649 
650  // h = \eta + d
651  Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
652 
653  // hu = h * u
654  Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
655 
656  // hv = h * v
657  Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
658  }
659  else
660  {
661  // h = \eta + d
662  Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
663 
664  // hu = h * u
665  Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
666 
667  // hv = h * v
668  Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
669  }
670 }
671 
673 {
674  int nq = GetTotPoints();
675 
676  // h = \eta + d
677  Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
678  m_fields[0]->UpdatePhys(), 1);
679 
680  // hu = h * u
681  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
682  m_fields[1]->UpdatePhys(), 1);
683 
684  // hv = h * v
685  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
686  m_fields[2]->UpdatePhys(), 1);
687 }
688 
689 /**
690  * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
691  * \f$ h\mathbf{v} \f$.
692  *
693  * @param physfield Velocity field.
694  * @param velocity Velocity field.
695  */
697  const Array<OneD, Array<OneD, NekDouble>> &physfield,
698  Array<OneD, Array<OneD, NekDouble>> &velocity)
699 {
700  const int npts = physfield[0].size();
701 
702  for (int i = 0; i < m_spacedim; ++i)
703  {
704  Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
705  }
706 }
707 
709 {
711  if (m_session->DefinesSolverInfo("UpwindType"))
712  {
713  std::string UpwindType;
714  UpwindType = m_session->GetSolverInfo("UpwindType");
715  if (UpwindType == "LinearAverage")
716  {
717  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
718  }
719  if (UpwindType == "LinearHLL")
720  {
721  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
722  }
723  }
724  SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
725  SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
726  SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
727 }
728 
729 } // 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)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Definition: LinearSWE.cpp:313
virtual ~LinearSWE()
Definition: LinearSWE.cpp:159
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Definition: LinearSWE.cpp:708
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
Compute the velocity field given the momentum .
Definition: LinearSWE.cpp:696
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: LinearSWE.cpp:164
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:96
virtual void v_ConservativeToPrimitive() override
Definition: LinearSWE.cpp:615
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
Definition: LinearSWE.cpp:460
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble >> &physarray, NekDouble time)
Definition: LinearSWE.cpp:356
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
Definition: LinearSWE.cpp:542
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
Definition: LinearSWE.cpp:57
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Definition: LinearSWE.cpp:214
virtual void v_PrimitiveToConservative() override
Definition: LinearSWE.cpp:672
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:75
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
Wall boundary condition.
Definition: LinearSWE.cpp:406
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:76
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:100
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
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