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