Nektar++
MMFSWE.cpp
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////////
2 //
3 // File: MMFSWE.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: MMF solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iomanip>
36 #include <iostream>
37 
38 #include <boost/algorithm/string/predicate.hpp>
39 #include <boost/core/ignore_unused.hpp>
40 
45 
46 namespace Nektar
47 {
48 std::string MMFSWE::className =
50  "MMFSWE", MMFSWE::create, "MMFSWE equation.");
51 
54  : UnsteadySystem(pSession, pGraph), MMFSystem(pSession, pGraph)
55 {
56  m_planeNumber = 0;
57 }
58 
59 /**
60  * @brief Initialisation object for the unsteady linear advection equation.
61  */
62 void MMFSWE::v_InitObject(bool DeclareFields)
63 {
64  // Call to the initialisation object
65  UnsteadySystem::v_InitObject(DeclareFields);
66 
67  int nq = m_fields[0]->GetNpoints();
68  int shapedim = m_fields[0]->GetShapeDimension();
69  Array<OneD, Array<OneD, NekDouble>> Anisotropy(shapedim);
70  for (int j = 0; j < shapedim; ++j)
71  {
72  Anisotropy[j] = Array<OneD, NekDouble>(nq, 1.0);
73  }
74 
75  MMFSystem::MMFInitObject(Anisotropy);
76 
77  // Load acceleration of gravity
78  m_session->LoadParameter("Gravity", m_g, 9.81);
79 
80  // Add Coriolois effects
81  m_session->LoadParameter("AddCoriolis", m_AddCoriolis, 1);
82 
83  // Add Rotation of the sphere along the pole
84  m_session->LoadParameter("AddRotation", m_AddRotation, 1);
85 
86  m_session->LoadParameter("AddRossbyDisturbance", m_RossbyDisturbance, 0);
87  m_session->LoadParameter("PurturbedJet", m_PurturbedJet, 0);
88 
89  // Define TestType
90  ASSERTL0(m_session->DefinesSolverInfo("TESTTYPE"),
91  "No TESTTYPE defined in session.");
92  std::string TestTypeStr = m_session->GetSolverInfo("TESTTYPE");
93  for (int i = 0; i < (int)SIZE_TestType; ++i)
94  {
95  if (boost::iequals(TestTypeMap[i], TestTypeStr))
96  {
97  m_TestType = (TestType)i;
98  break;
99  }
100  }
101 
102  // Variable Setting for each test problem
103  NekDouble gms = 9.80616;
104 
105  switch (m_TestType)
106  {
107  case eTestSteadyZonal:
108  {
109  NekDouble SecondToDay = 60.0 * 60.0 * 24.0;
110  NekDouble rad_earth = 6.37122 * 1000000;
111  NekDouble Omegams;
112 
113  // Nondimensionalized coeffs.
114  m_g = (gms * SecondToDay * SecondToDay) / rad_earth;
115 
116  m_session->LoadParameter("RotationAngle", m_alpha, 0.0);
117  m_session->LoadParameter("u0", m_u0, 2.0 * m_pi / 12.0);
118  m_session->LoadParameter("Omega", Omegams, 7.292 * 0.00001);
119  m_Omega = Omegams * SecondToDay;
120 
121  m_session->LoadParameter("H0", m_H0, 2.94 * 10000);
122  m_H0 = m_H0 / (rad_earth * gms);
123 
124  m_Hvar = (1.0 / m_g) * (m_Omega * m_u0 + 0.5 * m_u0 * m_u0);
125  }
126  break;
127 
128  case eTestUnsteadyZonal:
129  {
130  NekDouble SecondToDay = 60.0 * 60.0 * 24.0;
131  NekDouble rad_earth = 6.37122 * 1000000;
132  NekDouble Omegams;
133 
134  // Nondimensionalized coeffs.
135  m_g = (gms * SecondToDay * SecondToDay) / rad_earth;
136 
137  m_session->LoadParameter("RotationAngle", m_alpha, 0.0);
138  m_session->LoadParameter("u0", m_u0, 2.0 * m_pi / 12.0);
139  m_session->LoadParameter("Omega", Omegams, 7.292 * 0.00001);
140  m_Omega = Omegams * SecondToDay;
141 
142  m_H0 = 133681.0 / (rad_earth * gms); // m^2 / s^2
143  m_k2 = 10.0 / (rad_earth * gms); // m^2 / s^2
144  }
145  break;
146 
148  {
149  NekDouble SecondToDay = 60.0 * 60.0 * 24.0;
150  NekDouble rad_earth = 6.37122 * 1000000;
151  NekDouble Omegams;
152 
153  // Nondimensionalized coeffs.
154  m_g = (gms * SecondToDay * SecondToDay) / rad_earth;
155 
156  m_session->LoadParameter("RotationAngle", m_alpha, 0.0);
157 
158  m_session->LoadParameter("u0", m_u0, 20.0);
159  m_u0 = m_u0 * SecondToDay / rad_earth;
160 
161  m_session->LoadParameter("Omega", Omegams, 7.292 * 0.00001);
162  m_Omega = Omegams * SecondToDay;
163 
164  m_H0 = 5960.0 / rad_earth;
165  m_hs0 = 2000.0 / rad_earth;
166  }
167  break;
168 
169  case eTestUnstableJet:
170  {
171  NekDouble SecondToDay = 60.0 * 60.0 * 24.0;
172  NekDouble rad_earth = 6.37122 * 1000000;
173  NekDouble Omegams;
174 
175  // Nondimensionalized coeffs.
176  m_g = (gms * SecondToDay * SecondToDay) / rad_earth;
177 
178  m_session->LoadParameter("u0", m_u0, 80.0);
179  m_u0 = m_u0 * SecondToDay / rad_earth;
180 
181  m_session->LoadParameter("Omega", Omegams, 7.292 * 0.00001);
182  m_Omega = Omegams * SecondToDay;
183 
184  m_session->LoadParameter("H0", m_H0, 10000.0);
185  m_H0 = m_H0 / rad_earth;
186 
187  m_uthetamax = 80 * SecondToDay / rad_earth;
188  m_theta0 = m_pi / 7.0;
189  m_theta1 = m_pi / 2.0 - m_theta0;
190  m_en = exp(-4.0 / (m_theta1 - m_theta0) / (m_theta1 - m_theta0));
191  m_hbar = 120.0 / rad_earth;
192 
193  std::cout << "m_theta0 = " << m_theta0
194  << ", m_theta1 = " << m_theta1 << ", m_en = " << m_en
195  << ", m_hbar = " << m_hbar << std::endl;
196  }
197  break;
198 
199  case eTestRossbyWave:
200  {
201  NekDouble SecondToDay = 60.0 * 60.0 * 24.0;
202  NekDouble rad_earth = 6.37122 * 1000000;
203  NekDouble Omegams;
204 
205  // Nondimensionalized coeffs.
206  m_g = (gms * SecondToDay * SecondToDay) / rad_earth;
207 
208  m_session->LoadParameter("Omega", Omegams, 7.292 * 0.00001);
209  m_Omega = Omegams * SecondToDay;
210 
211  m_session->LoadParameter("H0", m_H0, 8000.0);
212  m_H0 = m_H0 / rad_earth;
213 
214  m_angfreq = 7.848 * 0.000001 * SecondToDay;
215  m_K = 7.848 * 0.000001 * SecondToDay;
216  }
217  break;
218 
219  default:
220  break;
221  }
222 
223  // TestVorticityComputation
225  {
227  }
228 
229  // If explicit it computes RHS and PROJECTION for the time integration
231  {
234  }
235  // Otherwise it gives an error (no implicit integration)
236  else
237  {
238  ASSERTL0(false, "Implicit unsteady Advection not set up.");
239  }
240 }
241 
242 /**
243  * @brief Unsteady linear advection equation destructor.
244  */
246 {
247 }
248 
250 {
251  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
252 
253  int i, nchk = 1;
254  int nvariables = 0;
255  int nfields = m_fields.size();
256  int nq = m_fields[0]->GetNpoints();
257 
258  if (m_intVariables.empty())
259  {
260  for (i = 0; i < nfields; ++i)
261  {
262  m_intVariables.push_back(i);
263  }
264  nvariables = nfields;
265  }
266  else
267  {
268  nvariables = m_intVariables.size();
269  }
270 
271  // Set up wrapper to fields data storage.
272  Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
273  Array<OneD, Array<OneD, NekDouble>> tmp(nvariables);
274 
275  // Order storage to list time-integrated fields first.
276  for (i = 0; i < nvariables; ++i)
277  {
278  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
279  m_fields[m_intVariables[i]]->SetPhysState(false);
280  }
281 
282  // Initialise time integration scheme
283  m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
284 
285  // Check uniqueness of checkpoint output
286  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
287  (m_checktime > 0.0 && m_checksteps == 0) ||
288  (m_checktime == 0.0 && m_checksteps > 0),
289  "Only one of IO_CheckTime and IO_CheckSteps "
290  "should be set!");
291 
292  LibUtilities::Timer timer;
293  bool doCheckTime = false;
294  int step = 0;
295  NekDouble intTime = 0.0;
296  NekDouble cpuTime = 0.0;
297  NekDouble elapsed = 0.0;
298 
299  NekDouble Mass = 0.0, Energy = 0.0, Enstrophy = 0.0, Vorticity = 0.0;
300  Array<OneD, NekDouble> zeta(nq);
301  Array<OneD, Array<OneD, NekDouble>> fieldsprimitive(nvariables);
302  for (int i = 0; i < nvariables; ++i)
303  {
304  fieldsprimitive[i] = Array<OneD, NekDouble>(nq);
305  }
306 
307  while (step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol)
308  {
309  timer.Start();
310  fields = m_intScheme->TimeIntegrate(step, m_timestep, m_ode);
311  timer.Stop();
312 
313  m_time += m_timestep;
314  elapsed = timer.TimePerTest(1);
315  intTime += elapsed;
316  cpuTime += elapsed;
317 
318  // Write out status information
319  if (m_infosteps && !((step + 1) % m_infosteps) &&
320  m_session->GetComm()->GetRank() == 0)
321  {
322  std::cout << "Steps: " << std::setw(8) << std::left << step + 1
323  << " "
324  << "Time: " << std::setw(12) << std::left << m_time;
325 
326  std::stringstream ss;
327  ss << cpuTime << "s";
328  std::cout << " CPU Time: " << std::setw(8) << std::left << ss.str()
329  << std::endl;
330 
331  // Printout Mass, Energy, Enstrophy
332  ConservativeToPrimitive(fields, fieldsprimitive);
333 
334  // Vorticity zeta
335  ComputeVorticity(fieldsprimitive[1], fieldsprimitive[2], zeta);
336  Vorticity = std::abs(m_fields[0]->Integral(zeta) - m_Vorticity0);
337 
338  // Masss = h^*
339  Mass = (ComputeMass(fieldsprimitive[0]) - m_Mass0) / m_Mass0;
340 
341  // Energy = 0.5*( h^*(u^2 + v^2) + g ( h^2 - h_s^s ) )
342  Energy = (ComputeEnergy(fieldsprimitive[0], fieldsprimitive[1],
343  fieldsprimitive[2]) -
344  m_Energy0) /
345  m_Energy0;
346 
347  // Enstrophy = 0.5/h^* ( \mathbf{k} \cdot (\nabla \times \mathbf{v}
348  // ) + f )^2
349  Enstrophy =
350  (ComputeEnstrophy(fieldsprimitive[0], fieldsprimitive[1],
351  fieldsprimitive[2]) -
352  m_Enstrophy0) /
353  m_Enstrophy0;
354 
355  std::cout << "dMass = " << std::setw(8) << std::left << Mass << " "
356  << ", dEnergy = " << std::setw(8) << std::left << Energy
357  << " "
358  << ", dEnstrophy = " << std::setw(8) << std::left
359  << Enstrophy << " "
360  << ", dVorticity = " << std::setw(8) << std::left
361  << Vorticity << std::endl
362  << std::endl;
363 
364  cpuTime = 0.0;
365  }
366 
367  // (h, hu, hv) -> (\eta, u, v)
369 
370  // Transform data into coefficient space
371  for (i = 0; i < nvariables; ++i)
372  {
373  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
374  m_fields[m_intVariables[i]]->FwdTransLocalElmt(
375  fields[i], m_fields[m_intVariables[i]]->UpdateCoeffs());
376  m_fields[m_intVariables[i]]->SetPhysState(false);
377  }
378 
379  // Write out checkpoint files
380  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
381  doCheckTime)
382  {
383  Checkpoint_Output(nchk++);
384  doCheckTime = false;
385 
386  // (\eta, u, v) -> (h, hu, hv)
388  }
389 
390  // Step advance
391  ++step;
392  }
393 
394  // Print out summary statistics
395  if (m_session->GetComm()->GetRank() == 0)
396  {
397  if (m_cflSafetyFactor > 0.0)
398  {
399  std::cout << "CFL safety factor : " << m_cflSafetyFactor
400  << std::endl
401  << "CFL time-step : " << m_timestep << std::endl;
402  }
403 
404  if (m_session->GetSolverInfo("Driver") != "SteadyState")
405  {
406  std::cout << "Time-integration : " << intTime << "s" << std::endl;
407  }
408  }
409 
410  for (i = 0; i < nvariables; ++i)
411  {
412  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
413  m_fields[m_intVariables[i]]->SetPhysState(true);
414  }
415 
416  // (h, hu, hv) -> (\eta, u, v)
418 
419  for (i = 0; i < nvariables; ++i)
420  {
421  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
422  m_fields[i]->UpdateCoeffs());
423  }
424 }
425 
426 void MMFSWE::DoOdeRhs(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
427  Array<OneD, Array<OneD, NekDouble>> &outarray,
428  const NekDouble time)
429 {
430  boost::ignore_unused(time);
431 
432  int i;
433  int nvariables = inarray.size();
434  int ncoeffs = GetNcoeffs();
435  int nq = GetTotPoints();
436 
437  // inarray in physical space
438  Array<OneD, Array<OneD, NekDouble>> physarray(nvariables);
439  Array<OneD, Array<OneD, NekDouble>> modarray(nvariables);
440 
441  for (i = 0; i < nvariables; ++i)
442  {
443  physarray[i] = Array<OneD, NekDouble>(nq);
444  modarray[i] = Array<OneD, NekDouble>(ncoeffs);
445  }
446 
447  // (h, hu, hv) -> (\eta, u, v)
448  ConservativeToPrimitive(inarray, physarray);
449 
450  // Weak Directional Derivative
451  WeakDGSWEDirDeriv(physarray, modarray);
452  AddDivForGradient(physarray, modarray);
453 
454  for (i = 0; i < nvariables; ++i)
455  {
456  Vmath::Neg(ncoeffs, modarray[i], 1);
457  }
458 
459  // coriolis forcing
460  if (m_AddCoriolis)
461  {
462  AddCoriolis(physarray, modarray);
463  }
464 
465  // Bottom Elevation Effect
466  // Add g H \nabla m_depth
467  AddElevationEffect(physarray, modarray);
468 
469  // Add terms concerning the rotation of the moving frame
470  if (m_AddRotation)
471  {
472  AddRotation(physarray, modarray);
473  }
474 
475  for (i = 0; i < nvariables; ++i)
476  {
477  m_fields[i]->MultiplyByElmtInvMass(modarray[i], modarray[i]);
478  m_fields[i]->BwdTrans(modarray[i], outarray[i]);
479  }
480 }
481 
483  const Array<OneD, Array<OneD, NekDouble>> &InField,
484  Array<OneD, Array<OneD, NekDouble>> &OutField)
485 {
486  int i;
487  int nq = GetNpoints();
488  int ncoeffs = GetNcoeffs();
489  int nTracePointsTot = GetTraceNpoints();
490  int nvariables = m_fields.size();
491 
493  Array<OneD, Array<OneD, NekDouble>> physfield(nvariables);
494 
495  for (i = 0; i < m_shapedim; ++i)
496  {
497  fluxvector[i] = Array<OneD, NekDouble>(nq);
498  }
499 
500  // InField is Primitive
501  for (i = 0; i < nvariables; ++i)
502  {
503  physfield[i] = InField[i];
504  }
505 
506  // Get the ith component of the flux vector in (physical space)
507  // fluxvector[0] = component for e^1 cdot \nabla \varphi
508  // fluxvector[1] = component for e^2 cdot \nabla \varphi
509  Array<OneD, NekDouble> fluxtmp(nq);
510  Array<OneD, NekDouble> tmp(ncoeffs);
511 
512  // Compute Divergence Components
513  for (i = 0; i < nvariables; ++i)
514  {
515  GetSWEFluxVector(i, physfield, fluxvector);
516 
517  OutField[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
518  for (int j = 0; j < m_shapedim; ++j)
519  {
520  // Directional derivation with respect to the j'th moving frame
521  // tmp_j = ( \nabla \phi, fluxvector[j] \mathbf{e}^j )
522  m_fields[i]->IProductWRTDirectionalDerivBase(m_movingframes[j],
523  fluxvector[j], tmp);
524  Vmath::Vadd(ncoeffs, &tmp[0], 1, &OutField[i][0], 1,
525  &OutField[i][0], 1);
526  }
527  }
528 
529  // Numerical Flux
530  Array<OneD, Array<OneD, NekDouble>> numfluxFwd(nvariables);
531  Array<OneD, Array<OneD, NekDouble>> numfluxBwd(nvariables);
532 
533  for (i = 0; i < nvariables; ++i)
534  {
535  numfluxFwd[i] = Array<OneD, NekDouble>(nTracePointsTot);
536  numfluxBwd[i] = Array<OneD, NekDouble>(nTracePointsTot);
537  }
538 
539  NumericalSWEFlux(physfield, numfluxFwd, numfluxBwd);
540 
541  // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
542  for (i = 0; i < nvariables; ++i)
543  {
544  Vmath::Neg(ncoeffs, OutField[i], 1);
545  m_fields[i]->AddFwdBwdTraceIntegral(numfluxFwd[i], numfluxBwd[i],
546  OutField[i]);
547  m_fields[i]->SetPhysState(false);
548  }
549 }
550 
551 // Substract 0.5 * g * H * H / || e^m ||^2 \nalba \cdot e^m
553  Array<OneD, Array<OneD, NekDouble>> &outarray)
554 {
555  // routine works for both primitive and conservative formulations
556  int ncoeffs = outarray[0].size();
557  int nq = physarray[0].size();
558 
560  Array<OneD, NekDouble> tmp(nq);
562  for (int i = 0; i < m_shapedim; ++i)
563  {
564  fluxvector[i] = Array<OneD, NekDouble>(nq);
565  }
566 
567  // Get 0.5 g H*H / || e^m ||^2
568  GetSWEFluxVector(3, physarray, fluxvector);
569 
570  Array<OneD, NekDouble> tmpc(ncoeffs);
571  for (int j = 0; j < m_shapedim; ++j)
572  {
573  Vmath::Vmul(nq, &fluxvector[0][0], 1, &m_DivMF[j][0], 1, &tmp[0], 1);
574 
575  Vmath::Neg(nq, &tmp[0], 1);
576  m_fields[0]->IProductWRTBase(tmp, tmpc);
577 
578  Vmath::Vadd(ncoeffs, outarray[j + 1], 1, tmpc, 1, outarray[j + 1], 1);
579  }
580 }
581 
583  const int i, const Array<OneD, const Array<OneD, NekDouble>> &physfield,
585 {
586  int nq = m_fields[0]->GetTotPoints();
587 
588  switch (i)
589  {
590  // flux function for the h equation = [(\eta + d ) u, (\eta + d) v ]
591  case 0:
592  {
593  // h in flux 1
594  Vmath::Vadd(nq, physfield[0], 1, m_depth, 1, flux[1], 1);
595 
596  // hu in flux 0
597  Vmath::Vmul(nq, flux[1], 1, physfield[1], 1, flux[0], 1);
598 
599  // hv in flux 1
600  Vmath::Vmul(nq, flux[1], 1, physfield[2], 1, flux[1], 1);
601  }
602  break;
603 
604  // flux function for the hu equation = [Hu^2 + 0.5 g H^2, Huv]
605  case 1:
606  {
607  Array<OneD, NekDouble> tmp(nq);
608 
609  // h in tmp
610  Vmath::Vadd(nq, physfield[0], 1, m_depth, 1, tmp, 1);
611 
612  // hu in flux 1
613  Vmath::Vmul(nq, tmp, 1, physfield[1], 1, flux[1], 1);
614 
615  // huu in flux 0
616  Vmath::Vmul(nq, flux[1], 1, physfield[1], 1, flux[0], 1);
617 
618  // hh in tmp
619  Vmath::Vmul(nq, tmp, 1, tmp, 1, tmp, 1);
620 
621  // huu + 0.5 g hh in flux 0
622  // Daxpy overwrites flux[0] on exit
623  Blas::Daxpy(nq, 0.5 * m_g, tmp, 1, flux[0], 1);
624 
625  // huv in flux 1
626  Vmath::Vmul(nq, flux[1], 1, physfield[2], 1, flux[1], 1);
627  }
628  break;
629 
630  // flux function for the hv equation = [Huv, Hv^2]
631  case 2:
632  {
633  Array<OneD, NekDouble> tmp(nq);
634 
635  // h in tmp
636  Vmath::Vadd(nq, physfield[0], 1, m_depth, 1, tmp, 1);
637 
638  // hv in flux 0
639  Vmath::Vmul(nq, tmp, 1, physfield[2], 1, flux[0], 1);
640 
641  // hvv in flux 1
642  Vmath::Vmul(nq, flux[0], 1, physfield[2], 1, flux[1], 1);
643 
644  // huv in flux 0
645  Vmath::Vmul(nq, flux[0], 1, physfield[1], 1, flux[0], 1);
646 
647  // hh in tmp
648  Vmath::Vmul(nq, tmp, 1, tmp, 1, tmp, 1);
649 
650  // hvv + 0.5 g hh in flux 1
651  Blas::Daxpy(nq, 0.5 * m_g, tmp, 1, flux[1], 1);
652  }
653  break;
654 
655  // flux function 0.5 g h * h
656  case 3:
657  {
659  flux[0] = Array<OneD, NekDouble>(nq, 0.0);
660 
661  // h in tmp
662  Vmath::Vadd(nq, physfield[0], 1, m_depth, 1, h, 1);
663 
664  // hh in tmp
665  Vmath::Vmul(nq, h, 1, h, 1, h, 1);
666 
667  // 0.5 g hh in flux 0
668  Blas::Daxpy(nq, 0.5 * m_g, h, 1, flux[0], 1);
669  }
670  break;
671 
672  default:
673  break;
674  }
675 }
676 
678  Array<OneD, Array<OneD, NekDouble>> &numfluxFwd,
679  Array<OneD, Array<OneD, NekDouble>> &numfluxBwd)
680 {
681 
682  int i, k;
683  int nTraceNumPoints = GetTraceTotPoints();
684  int nvariables = 3; // only the dependent variables
685 
686  // get temporary arrays
687  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
688  Array<OneD, Array<OneD, NekDouble>> Bwd(nvariables);
689  Array<OneD, NekDouble> DepthFwd(nTraceNumPoints);
690  Array<OneD, NekDouble> DepthBwd(nTraceNumPoints);
691 
692  for (i = 0; i < nvariables; ++i)
693  {
694  Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
695  Bwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
696  }
697 
698  // get the physical values at the trace
699  for (i = 0; i < nvariables; ++i)
700  {
701  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
702 
703  // Copy Fwd to Bwd at boundaries
704  CopyBoundaryTrace(Fwd[i], Bwd[i], SolverUtils::eFwdEQBwd);
705  }
706  m_fields[0]->GetFwdBwdTracePhys(m_depth, DepthFwd, DepthBwd);
707  CopyBoundaryTrace(DepthFwd, DepthBwd, SolverUtils::eFwdEQBwd);
708 
709  // note that we are using the same depth - i.e. the depth is assumed
710  // continuous...
711  switch (m_upwindType)
712  {
713  case SolverUtils::eHLLC:
714  {
715  // rotate the values to the normal direction
716  NekDouble tmpX, tmpY;
717 
718  // Fwd[1] = hu^+ = ( h^1^+ (e^1 \cdot n) + h^2^+ ( e^2 \cdot n) )
719  // Fwd[2] = hv^+ = ( h^1^+ (e^1 \cdot n^{\perp}) + h^2^+ ( e^2 \cdot
720  // n^{\perp}) )
721  // Bwd[1] = hu^- = ( h^1^- (e^1 \cdot n) + h^2^- ( e^2 \cdot n) )
722  // Bwd[2] = hv^- = ( h^1^- (e^1 \cdot n^{\perp}) + h^2^- ( e^2 \cdot
723  // n^{\perp}) )
724 
725  for (k = 0; k < nTraceNumPoints; ++k)
726  {
727  tmpX = Fwd[1][k] * m_ncdotMFFwd[0][k] +
728  Fwd[2][k] * m_ncdotMFFwd[1][k];
729  tmpY = Fwd[1][k] * m_nperpcdotMFFwd[0][k] +
730  Fwd[2][k] * m_nperpcdotMFFwd[1][k];
731  Fwd[1][k] = tmpX;
732  Fwd[2][k] = tmpY;
733 
734  tmpX = Bwd[1][k] * m_ncdotMFBwd[0][k] +
735  Bwd[2][k] * m_ncdotMFBwd[1][k];
736  tmpY = Bwd[1][k] * m_nperpcdotMFBwd[0][k] +
737  Bwd[2][k] * m_nperpcdotMFBwd[1][k];
738  Bwd[1][k] = tmpX;
739  Bwd[2][k] = tmpY;
740  }
741 
742  // Solve the Riemann problem
743  NekDouble denomFwd, denomBwd;
744  Array<OneD, NekDouble> numfluxF(nvariables);
745  Array<OneD, NekDouble> numfluxB(nvariables);
746 
747  NekDouble eF1n, eF2n, eB1n, eB2n;
748  NekDouble eF1t, eF2t, eB1t, eB2t;
749  for (k = 0; k < nTraceNumPoints; ++k)
750  {
751  RiemannSolverHLLC(k, Fwd[0][k] + DepthFwd[k], Fwd[1][k],
752  Fwd[2][k], Bwd[0][k] + DepthFwd[k], Bwd[1][k],
753  Bwd[2][k], numfluxF, numfluxB);
754 
755  // uflux = 1/[ ( e^1 \cdot n) ( e^2 \cdot n^{\perp} ) - (e^2
756  // \cdot n)(e^1 \cdot n^{\perp} ) ] ( e^2 \cdot n^{\perp} huflux
757  // - e^2 \cdot n hvflux
758  // vflux = -1/[ ( e^1 \cdot n) ( e^2 \cdot n^{\perp} ) - (e^2
759  // \cdot n)(e^1 \cdot n^{\perp} ) ] ( -e^1 \cdot n^{\perp}
760  // huflux + e^1 \cdot n hvflux
761  eF1n = m_ncdotMFFwd[0][k];
762  eF2n = m_ncdotMFFwd[1][k];
763  eB1n = m_ncdotMFBwd[0][k];
764  eB2n = m_ncdotMFBwd[1][k];
765 
766  eF1t = m_nperpcdotMFFwd[0][k];
767  eF2t = m_nperpcdotMFFwd[1][k];
768  eB1t = m_nperpcdotMFBwd[0][k];
769  eB2t = m_nperpcdotMFBwd[1][k];
770 
771  denomFwd = eF1n * eF2t - eF2n * eF1t;
772  denomBwd = eB1n * eB2t - eB2n * eB1t;
773 
774  numfluxFwd[0][k] = numfluxF[0];
775  numfluxFwd[1][k] = (1.0 / denomFwd) *
776  (eF2t * numfluxF[1] - eF2n * numfluxF[2]);
777  numfluxFwd[2][k] =
778  (1.0 / denomFwd) *
779  (-1.0 * eF1t * numfluxF[1] + eF1n * numfluxF[2]);
780 
781  numfluxBwd[0][k] = 1.0 * numfluxB[0];
782  numfluxBwd[1][k] = (1.0 / denomBwd) *
783  (eB2t * numfluxB[1] - eB2n * numfluxB[2]);
784  numfluxBwd[2][k] =
785  (1.0 / denomBwd) *
786  (-1.0 * eB1t * numfluxB[1] + eB1n * numfluxB[2]);
787  }
788  }
789  break;
790 
794  {
795  Array<OneD, NekDouble> numfluxF(nvariables * (m_shapedim + 1));
796  Array<OneD, NekDouble> numfluxB(nvariables * (m_shapedim + 1));
797 
798  for (k = 0; k < nTraceNumPoints; ++k)
799  {
801  {
802  AverageFlux(k, Fwd[0][k] + DepthFwd[k], Fwd[1][k],
803  Fwd[2][k], Bwd[0][k] + DepthFwd[k], Bwd[1][k],
804  Bwd[2][k], numfluxF, numfluxB);
805  }
806 
808  {
809  LaxFriedrichFlux(k, Fwd[0][k] + DepthFwd[k], Fwd[1][k],
810  Fwd[2][k], Bwd[0][k] + DepthFwd[k],
811  Bwd[1][k], Bwd[2][k], numfluxF, numfluxB);
812  }
813 
815  {
816  RusanovFlux(k, Fwd[0][k] + DepthFwd[k], Fwd[1][k],
817  Fwd[2][k], Bwd[0][k] + DepthFwd[k], Bwd[1][k],
818  Bwd[2][k], numfluxF, numfluxB);
819  }
820 
821  int indx;
822  NekDouble tmpF0, tmpF1, tmpB0, tmpB1;
823  for (i = 0; i < nvariables; ++i)
824  {
825  indx = i * (m_shapedim + 1);
826 
827  tmpF0 = numfluxF[indx] * m_ncdotMFFwd[0][k];
828  tmpF1 = numfluxF[indx + 1] * m_ncdotMFFwd[1][k];
829 
830  tmpB0 = numfluxB[indx] * m_ncdotMFBwd[0][k];
831  tmpB1 = numfluxB[indx + 1] * m_ncdotMFBwd[1][k];
832 
833  numfluxFwd[i][k] = tmpF0 + tmpF1 + numfluxF[indx + 2];
834  numfluxBwd[i][k] = tmpB0 + tmpB1 + numfluxB[indx + 2];
835  }
836  }
837  }
838  break;
839 
840  default:
841  {
842  ASSERTL0(false, "populate switch statement for upwind flux");
843  }
844  break;
845  }
846 }
847 
848 void MMFSWE::RiemannSolverHLLC(const int index, NekDouble hL, NekDouble uL,
849  NekDouble vL, NekDouble hR, NekDouble uR,
850  NekDouble vR, Array<OneD, NekDouble> &numfluxF,
851  Array<OneD, NekDouble> &numfluxB)
852 {
853  boost::ignore_unused(index);
854 
855  NekDouble g = m_g;
856 
857  NekDouble cL = sqrt(g * hL);
858  NekDouble cR = sqrt(g * hR);
859 
860  NekDouble uRF, vRF, uLB, vLB;
861  NekDouble hstarF, hstarB;
862 
863  // Temporary assignments
864  uRF = uR;
865  vRF = vR;
866  uLB = uL;
867  vLB = vL;
868 
869  // the two-rarefaction wave assumption
870  hstarF = 0.5 * (cL + cR) + 0.25 * (uL - uRF);
871  hstarF *= hstarF;
872  hstarF *= (1.0 / g);
873 
874  hstarB = 0.5 * (cL + cR) + 0.25 * (uLB - uR);
875  hstarB *= hstarB;
876  hstarB *= (1.0 / g);
877 
878  NekDouble hfluxF, hufluxF, hvfluxF;
879  NekDouble hfluxB, hufluxB, hvfluxB;
880  Computehhuhvflux(hL, uL, vL, hR, uRF, vRF, hstarF, hfluxF, hufluxF,
881  hvfluxF);
882  Computehhuhvflux(hL, uLB, vLB, hR, uR, vR, hstarB, hfluxB, hufluxB,
883  hvfluxB);
884 
885  numfluxF[0] = hfluxF;
886  numfluxF[1] = hufluxF;
887  numfluxF[2] = hvfluxF;
888 
889  numfluxB[0] = hfluxB;
890  numfluxB[1] = hufluxB;
891  numfluxB[2] = hvfluxB;
892 }
893 
895  NekDouble hR, NekDouble uR, NekDouble vR,
896  NekDouble hstar, NekDouble &hflux,
897  NekDouble &huflux, NekDouble &hvflux)
898 {
899  NekDouble g = m_g;
900 
901  NekDouble hC, huC, hvC, SL, SR, Sstar;
902  NekDouble cL = sqrt(g * hL);
903  NekDouble cR = sqrt(g * hR);
904 
905  // Compute SL
906  if (hstar > hL)
907  SL = uL - cL * sqrt(0.5 * ((hstar * hstar + hstar * hL) / (hL * hL)));
908  else
909  SL = uL - cL;
910 
911  // Compute SR
912  if (hstar > hR)
913  SR = uR + cR * sqrt(0.5 * ((hstar * hstar + hstar * hR) / (hR * hR)));
914  else
915  SR = uR + cR;
916 
917  if (fabs(hR * (uR - SR) - hL * (uL - SL)) <= 1.0e-15)
918  Sstar = 0.0;
919  else
920  Sstar = (SL * hR * (uR - SR) - SR * hL * (uL - SL)) /
921  (hR * (uR - SR) - hL * (uL - SL));
922 
923  if (SL >= 0)
924  {
925  hflux = hL * uL;
926  huflux = uL * uL * hL + 0.5 * g * hL * hL;
927  hvflux = hL * uL * vL;
928  }
929  else if (SR <= 0)
930  {
931  hflux = hR * uR;
932  huflux = uR * uR * hR + 0.5 * g * hR * hR;
933  hvflux = hR * uR * vR;
934  }
935  else
936  {
937  if ((SL < 0) && (Sstar >= 0))
938  {
939  hC = hL * ((SL - uL) / (SL - Sstar));
940  huC = hC * Sstar;
941  hvC = hC * vL;
942 
943  hflux = hL * uL + SL * (hC - hL);
944  huflux = (uL * uL * hL + 0.5 * g * hL * hL) + SL * (huC - hL * uL);
945  hvflux = (uL * vL * hL) + SL * (hvC - hL * vL);
946  }
947  else
948  {
949  hC = hR * ((SR - uR) / (SR - Sstar));
950  huC = hC * Sstar;
951  hvC = hC * vR;
952 
953  hflux = hR * uR + SR * (hC - hR);
954  huflux = (uR * uR * hR + 0.5 * g * hR * hR) + SR * (huC - hR * uR);
955  hvflux = (uR * vR * hR) + SR * (hvC - hR * vR);
956  }
957  }
958 }
959 
960 void MMFSWE::AverageFlux(const int index, NekDouble hL, NekDouble uL,
961  NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR,
962  Array<OneD, NekDouble> &numfluxF,
963  Array<OneD, NekDouble> &numfluxB)
964 {
965  NekDouble MageF1, MageF2, MageB1, MageB2;
966  NekDouble eF1_cdot_eB1, eF1_cdot_eB2;
967  NekDouble eF2_cdot_eB1, eF2_cdot_eB2;
968 
969  NekDouble g = m_g;
970  NekDouble uRF, vRF, uLB, vLB;
971 
972  ComputeMagAndDot(index, MageF1, MageF2, MageB1, MageB2, eF1_cdot_eB1,
973  eF1_cdot_eB2, eF2_cdot_eB1, eF2_cdot_eB2);
974 
975  // uRF = uR component in moving frames e^{Fwd}
976  // vRF = vR component in moving frames e^{Fwd}
977  uRF = (uR * eF1_cdot_eB1 + vR * eF1_cdot_eB2) / MageF1;
978  vRF = (uR * eF2_cdot_eB1 + vR * eF2_cdot_eB2) / MageF2;
979 
980  numfluxF[0] = 0.5 * (hL * uL + hR * uRF);
981  numfluxF[1] = 0.5 * (hL * vL + hR * vRF);
982  numfluxF[2] = 0.0;
983 
984  numfluxF[3] =
985  0.5 * (hL * uL * uL + hR * uRF * uRF + 0.5 * g * (hL * hL + hR * hR));
986  numfluxF[4] = 0.5 * (hL * uL * vL + hR * uRF * vRF);
987  numfluxF[5] = 0.0;
988 
989  numfluxF[6] = 0.5 * (hL * uL * vL + hR * uRF * vRF);
990  numfluxF[7] =
991  0.5 * (hL * vL * vL + hR * vRF * vRF + 0.5 * g * (hL * hL + hR * hR));
992  numfluxF[8] = 0.0;
993 
994  // uLB = uL component in moving frames e^{Bwd}
995  // vLB = vL component in moving frames e^{Bwd}
996  uLB = (uL * eF1_cdot_eB1 + vL * eF2_cdot_eB1) / MageB1;
997  vLB = (uL * eF1_cdot_eB2 + vL * eF2_cdot_eB2) / MageB2;
998 
999  numfluxB[0] = 0.5 * (hR * uR + hR * uLB);
1000  numfluxB[1] = 0.5 * (hR * vR + hR * vLB);
1001  numfluxB[2] = 0.0;
1002 
1003  numfluxB[3] =
1004  0.5 * (hR * uR * uR + hR * uLB * uLB + 0.5 * g * (hR * hR + hL * hL));
1005  numfluxB[4] = 0.5 * (hR * uR * vR + hR * uLB * vLB);
1006  numfluxB[5] = 0.0;
1007 
1008  numfluxB[6] = 0.5 * (hR * uR * vR + hR * uLB * vLB);
1009  numfluxB[7] =
1010  0.5 * (hR * vR * vR + hR * vLB * vLB + 0.5 * g * (hR * hR + hL * hL));
1011  numfluxB[8] = 0.0;
1012 }
1013 
1014 void MMFSWE::LaxFriedrichFlux(const int index, NekDouble hL, NekDouble uL,
1015  NekDouble vL, NekDouble hR, NekDouble uR,
1016  NekDouble vR, Array<OneD, NekDouble> &numfluxF,
1017  Array<OneD, NekDouble> &numfluxB)
1018 {
1019  int nvariables = 3;
1020  NekDouble MageF1, MageF2, MageB1, MageB2;
1021  NekDouble eF1_cdot_eB1, eF1_cdot_eB2;
1022  NekDouble eF2_cdot_eB1, eF2_cdot_eB2;
1023 
1024  NekDouble g = m_g;
1025  NekDouble uRF, vRF, uLB, vLB;
1026  NekDouble velL, velR, lambdaF, lambdaB;
1027 
1028  Array<OneD, NekDouble> EigF(nvariables);
1029  Array<OneD, NekDouble> EigB(nvariables);
1030 
1031  // Compute Magnitude and Dot product of moving frames for the index
1032  ComputeMagAndDot(index, MageF1, MageF2, MageB1, MageB2, eF1_cdot_eB1,
1033  eF1_cdot_eB2, eF2_cdot_eB1, eF2_cdot_eB2);
1034 
1035  // Get the velocity in the normal to the edge
1036  velL = uL * m_ncdotMFFwd[0][index] + vL * m_ncdotMFFwd[1][index];
1037  velR = -1.0 * (uR * m_ncdotMFBwd[0][index] + vR * m_ncdotMFBwd[1][index]);
1038 
1039  EigF[0] = velL - sqrt(g * hL);
1040  EigF[1] = velL;
1041  EigF[2] = velL + sqrt(g * hL);
1042 
1043  EigB[0] = velR - sqrt(g * hR);
1044  EigB[1] = velR;
1045  EigB[2] = velR + sqrt(g * hR);
1046 
1047  lambdaF = Vmath::Vamax(nvariables, EigF, 1);
1048  lambdaB = Vmath::Vamax(nvariables, EigB, 1);
1049 
1050  // uRF = uR component in moving frames e^{Fwd}
1051  // vRF = vR component in moving frames e^{Fwd}
1052  uRF = (uR * eF1_cdot_eB1 + vR * eF1_cdot_eB2) / MageF1;
1053  vRF = (uR * eF2_cdot_eB1 + vR * eF2_cdot_eB2) / MageF2;
1054 
1055  numfluxF[0] = 0.5 * (hL * uL + hR * uRF);
1056  numfluxF[1] = 0.5 * (hL * vL + hR * vRF);
1057  numfluxF[2] = 0.5 * lambdaF * (hL - hR);
1058 
1059  numfluxF[3] = 0.5 * (hL * uL * uL * MageF1 + hR * uRF * uRF * MageB1 +
1060  0.5 * g * (hL * hL + hR * hR));
1061  numfluxF[4] = 0.5 * (hL * uL * vL * MageF1 + hR * uRF * vRF * MageB1);
1062  numfluxF[5] = 0.5 * lambdaF * (uL * hL - uRF * hR);
1063 
1064  numfluxF[6] = 0.5 * (hL * uL * vL * MageF2 + hR * uRF * vRF * MageB2);
1065  numfluxF[7] = 0.5 * (hL * vL * vL * MageF2 + hR * vRF * vRF * MageB2 +
1066  0.5 * g * (hL * hL + hR * hR));
1067  numfluxF[8] = 0.5 * lambdaF * (vL * hL - vRF * hR);
1068 
1069  // uLB = uL component in moving frames e^{Bwd}
1070  // vLB = vL component in moving frames e^{Bwd}
1071  uLB = (uL * eF1_cdot_eB1 + vL * eF2_cdot_eB1) / MageB1;
1072  vLB = (uL * eF1_cdot_eB2 + vL * eF2_cdot_eB2) / MageB2;
1073 
1074  numfluxB[0] = 0.5 * (hR * uR + hR * uLB);
1075  numfluxB[1] = 0.5 * (hR * vR + hR * vLB);
1076  numfluxB[2] = 0.5 * lambdaB * (hL - hR);
1077 
1078  numfluxB[3] = 0.5 * (hR * uR * uR * MageB1 + hR * uLB * uLB * MageF1 +
1079  0.5 * g * (hR * hR + hL * hL));
1080  numfluxB[4] = 0.5 * (hR * uR * vR * MageB1 + hR * uLB * vLB * MageF1);
1081  numfluxB[5] = 0.5 * lambdaB * (uLB * hL - uR * hR);
1082 
1083  numfluxB[6] = 0.5 * (hR * uR * vR * MageB2 + hR * uLB * vLB * MageF2);
1084  numfluxB[7] = 0.5 * (hR * vR * vR * MageB2 + hR * vLB * vLB * MageF2 +
1085  0.5 * g * (hR * hR + hL * hL));
1086  numfluxB[8] = 0.5 * lambdaB * (vLB * hL - vR * hR);
1087 }
1088 
1089 void MMFSWE::RusanovFlux(const int index, NekDouble hL, NekDouble uL,
1090  NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR,
1091  Array<OneD, NekDouble> &numfluxF,
1092  Array<OneD, NekDouble> &numfluxB)
1093 {
1094  int nvariables = 3;
1095  NekDouble MageF1, MageF2, MageB1, MageB2;
1096  NekDouble eF1_cdot_eB1, eF1_cdot_eB2;
1097  NekDouble eF2_cdot_eB1, eF2_cdot_eB2;
1098 
1099  NekDouble g = m_g;
1100  NekDouble uRF, vRF, uLB, vLB;
1101  NekDouble velL, velR;
1102 
1103  Array<OneD, NekDouble> EigF(nvariables);
1104  Array<OneD, NekDouble> EigB(nvariables);
1105 
1106  // Compute Magnitude and Dot product of moving frames for the index
1107  ComputeMagAndDot(index, MageF1, MageF2, MageB1, MageB2, eF1_cdot_eB1,
1108  eF1_cdot_eB2, eF2_cdot_eB1, eF2_cdot_eB2);
1109 
1110  // Get the velocity in the normal to the edge
1111  velL = uL * m_ncdotMFFwd[0][index] + vL * m_ncdotMFFwd[1][index];
1112  velR = -1.0 * (uR * m_ncdotMFBwd[0][index] + vR * m_ncdotMFBwd[1][index]);
1113 
1114  NekDouble SL, SR;
1115  SL = fabs(velL) + sqrt(g * hL);
1116  SR = fabs(velR) + sqrt(g * hR);
1117 
1118  NekDouble S;
1119  if (SL > SR)
1120  S = SL;
1121  else
1122  S = SR;
1123 
1124  // uRF = uR component in moving frames e^{Fwd}
1125  // vRF = vR component in moving frames e^{Fwd}
1126  uRF = (uR * eF1_cdot_eB1 + vR * eF1_cdot_eB2) / MageF1;
1127  vRF = (uR * eF2_cdot_eB1 + vR * eF2_cdot_eB2) / MageF2;
1128 
1129  numfluxF[0] = 0.5 * (hL * uL + hR * uRF);
1130  numfluxF[1] = 0.5 * (hL * vL + hR * vRF);
1131  numfluxF[2] = 0.5 * S * (hL - hR);
1132 
1133  numfluxF[3] =
1134  0.5 * (hL * uL * uL + hR * uRF * uRF + 0.5 * g * (hL * hL + hR * hR));
1135  numfluxF[4] = 0.5 * (hL * uL * vL + hR * uRF * vRF);
1136  numfluxF[5] = 0.5 * S * (uL * hL - uRF * hR);
1137 
1138  numfluxF[6] = 0.5 * (hL * uL * vL + hR * uRF * vRF);
1139  numfluxF[7] =
1140  0.5 * (hL * vL * vL + hR * vRF * vRF + 0.5 * g * (hL * hL + hR * hR));
1141  numfluxF[8] = 0.5 * S * (vL * hL - vRF * hR);
1142 
1143  // uLB = uL component in moving frames e^{Bwd}
1144  // vLB = vL component in moving frames e^{Bwd}
1145  uLB = (uL * eF1_cdot_eB1 + vL * eF2_cdot_eB1) / MageB1;
1146  vLB = (uL * eF1_cdot_eB2 + vL * eF2_cdot_eB2) / MageB2;
1147 
1148  numfluxB[0] = 0.5 * (hR * uR + hR * uLB);
1149  numfluxB[1] = 0.5 * (hR * vR + hR * vLB);
1150  numfluxB[2] = 0.5 * S * (hL - hR);
1151 
1152  numfluxB[3] =
1153  0.5 * (hR * uR * uR + hR * uLB * uLB + 0.5 * g * (hR * hR + hL * hL));
1154  numfluxB[4] = 0.5 * (hR * uR * vR + hR * uLB * vLB);
1155  numfluxB[5] = 0.5 * S * (uLB * hL - uR * hR);
1156 
1157  numfluxB[6] = 0.5 * (hR * uR * vR + hR * uLB * vLB);
1158  numfluxB[7] =
1159  0.5 * (hR * vR * vR + hR * vLB * vLB + 0.5 * g * (hR * hR + hL * hL));
1160  numfluxB[8] = 0.5 * S * (vLB * hL - vR * hR);
1161 }
1162 
1163 void MMFSWE::ComputeMagAndDot(const int index, NekDouble &MageF1,
1164  NekDouble &MageF2, NekDouble &MageB1,
1165  NekDouble &MageB2, NekDouble &eF1_cdot_eB1,
1166  NekDouble &eF1_cdot_eB2, NekDouble &eF2_cdot_eB1,
1167  NekDouble &eF2_cdot_eB2)
1168 {
1169  NekDouble MF1x, MF1y, MF1z, MF2x, MF2y, MF2z;
1170  NekDouble MB1x, MB1y, MB1z, MB2x, MB2y, MB2z;
1171 
1172  MF1x = m_MFtraceFwd[0][0][index];
1173  MF1y = m_MFtraceFwd[0][1][index];
1174  MF1z = m_MFtraceFwd[0][2][index];
1175 
1176  MF2x = m_MFtraceFwd[1][0][index];
1177  MF2y = m_MFtraceFwd[1][1][index];
1178  MF2z = m_MFtraceFwd[1][2][index];
1179 
1180  MB1x = m_MFtraceBwd[0][0][index];
1181  MB1y = m_MFtraceBwd[0][1][index];
1182  MB1z = m_MFtraceBwd[0][2][index];
1183 
1184  MB2x = m_MFtraceBwd[1][0][index];
1185  MB2y = m_MFtraceBwd[1][1][index];
1186  MB2z = m_MFtraceBwd[1][2][index];
1187 
1188  // MFtrace = MFtrace [ j*spacedim + k ], j = shape, k = sapce
1189  MageF1 = MF1x * MF1x + MF1y * MF1y + MF1z * MF1z;
1190  MageF2 = MF2x * MF2x + MF2y * MF2y + MF2z * MF2z;
1191  MageB1 = MB1x * MB1x + MB1y * MB1y + MB1z * MB1z;
1192  MageB2 = MB2x * MB2x + MB2y * MB2y + MB2z * MB2z;
1193 
1194  eF1_cdot_eB1 = MF1x * MB1x + MF1y * MB1y + MF1z * MB1z;
1195  eF1_cdot_eB2 = MF1x * MB2x + MF1y * MB2y + MF1z * MB2z;
1196  eF2_cdot_eB1 = MF2x * MB1x + MF2y * MB1y + MF2z * MB1z;
1197  eF2_cdot_eB2 = MF2x * MB2x + MF2y * MB2y + MF2z * MB2z;
1198 }
1199 
1200 // Add Coriolis factors
1202  Array<OneD, Array<OneD, NekDouble>> &outarray)
1203 {
1204  int ncoeffs = outarray[0].size();
1205  int nq = physarray[0].size();
1206 
1207  Array<OneD, NekDouble> h(nq);
1208  Array<OneD, NekDouble> tmp(nq);
1209  Array<OneD, NekDouble> tmpc(ncoeffs);
1210 
1211  // physarray is primitive
1212  // conservative formulation compute h
1213  // h = \eta + d
1214  Vmath::Vadd(nq, physarray[0], 1, m_depth, 1, h, 1);
1215 
1216  int indx = 0;
1217  for (int j = 0; j < m_shapedim; ++j)
1218  {
1219  if (j == 0)
1220  {
1221  indx = 2;
1222  }
1223 
1224  else if (j == 1)
1225  {
1226  indx = 1;
1227  }
1228 
1229  // add to hu equation
1230  Vmath::Vmul(nq, m_coriolis, 1, physarray[indx], 1, tmp, 1);
1231  Vmath::Vmul(nq, h, 1, tmp, 1, tmp, 1);
1232 
1233  if (j == 1)
1234  {
1235  Vmath::Neg(nq, tmp, 1);
1236  }
1237 
1238  // N \cdot (e^1 \times e^2 )
1239  // Vmath::Vmul(nq, &m_MF1crossMF2dotSN[0], 1, &tmp[0], 1, &tmp[0], 1);
1240  m_fields[0]->IProductWRTBase(tmp, tmpc);
1241  Vmath::Vadd(ncoeffs, tmpc, 1, outarray[j + 1], 1, outarray[j + 1], 1);
1242  }
1243 }
1244 
1245 // Compuate g H \nabla m_depth
1247  Array<OneD, Array<OneD, NekDouble>> &outarray)
1248 {
1249  int ncoeffs = outarray[0].size();
1250  int nq = physarray[0].size();
1251 
1252  Array<OneD, NekDouble> h(nq);
1253  Array<OneD, NekDouble> tmp(nq);
1254  Array<OneD, NekDouble> tmpc(ncoeffs);
1255 
1256  // physarray is primitive
1257  // conservative formulation compute h
1258  // h = \eta + d
1259  Vmath::Vadd(nq, physarray[0], 1, m_depth, 1, h, 1);
1260 
1261  for (int j = 0; j < m_shapedim; ++j)
1262  {
1263  Vmath::Vmul(nq, &h[0], 1, &m_Derivdepth[j][0], 1, &tmp[0], 1);
1264  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
1265 
1266  m_fields[0]->IProductWRTBase(tmp, tmpc);
1267 
1268  Vmath::Vadd(ncoeffs, tmpc, 1, outarray[j + 1], 1, outarray[j + 1], 1);
1269  }
1270 }
1271 
1272 // =================================================
1273 // Add rotational factors
1274 // =================================================
1276  Array<OneD, Array<OneD, NekDouble>> &outarray)
1277 {
1278  // routine works for both primitive and conservative formulations
1279  int ncoeffs = outarray[0].size();
1280  int nq = physarray[0].size();
1281 
1282  // Compute h
1283  Array<OneD, NekDouble> h(nq);
1284  Vmath::Vadd(nq, &physarray[0][0], 1, &m_depth[0], 1, &h[0], 1);
1285 
1286  Array<OneD, NekDouble> de0dt_cdot_e0;
1287  Array<OneD, NekDouble> de0dt_cdot_e1;
1288  Array<OneD, NekDouble> de1dt_cdot_e0;
1289  Array<OneD, NekDouble> de1dt_cdot_e1;
1290  Compute_demdt_cdot_ek(0, 0, physarray, de0dt_cdot_e0);
1291  Compute_demdt_cdot_ek(1, 0, physarray, de1dt_cdot_e0);
1292  Compute_demdt_cdot_ek(0, 1, physarray, de0dt_cdot_e1);
1293  Compute_demdt_cdot_ek(1, 1, physarray, de1dt_cdot_e1);
1294 
1295  Array<OneD, NekDouble> Rott1(nq);
1296  Array<OneD, NekDouble> Rott2(nq);
1297  Vmath::Vmul(nq, physarray[1], 1, de0dt_cdot_e0, 1, Rott1, 1);
1298  Vmath::Vmul(nq, physarray[1], 1, de0dt_cdot_e1, 1, Rott2, 1);
1299  Vmath::Vvtvp(nq, physarray[2], 1, de1dt_cdot_e0, 1, Rott1, 1, Rott1, 1);
1300  Vmath::Vvtvp(nq, physarray[2], 1, de1dt_cdot_e1, 1, Rott2, 1, Rott2, 1);
1301 
1302  // Multiply H and \partial \phi / \partial t which is assumed to be u_{\phi}
1303  Vmath::Vmul(nq, &h[0], 1, &Rott1[0], 1, &Rott1[0], 1);
1304  Vmath::Vmul(nq, &h[0], 1, &Rott2[0], 1, &Rott2[0], 1);
1305 
1306  Vmath::Neg(nq, Rott1, 1);
1307  Vmath::Neg(nq, Rott2, 1);
1308 
1309  Array<OneD, NekDouble> tmpc1(ncoeffs);
1310  Array<OneD, NekDouble> tmpc2(ncoeffs);
1311  m_fields[0]->IProductWRTBase(Rott1, tmpc1);
1312  m_fields[0]->IProductWRTBase(Rott2, tmpc2);
1313 
1314  Vmath::Vadd(ncoeffs, tmpc1, 1, outarray[1], 1, outarray[1], 1);
1315  Vmath::Vadd(ncoeffs, tmpc2, 1, outarray[2], 1, outarray[2], 1);
1316 }
1317 
1318 // =====================================================================
1319 // Compute \frac{d e^{indm}}{dt} \cdot e^{indk} / \| e^{indk} \|
1320 // Equivalently, [ \frac{d e^{indm}}{d \xi_1} \frac{ d \xi_1}{d \phi} + \frac{d
1321 // e^{indm}}{d \xi_2} \frac{ d \xi_2}{d \phi} ] \frac{d \phi}{dt}
1323  const int indm, const int indk,
1324  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
1325  Array<OneD, NekDouble> &outarray)
1326 {
1327  int j, k;
1328  int nq = m_fields[0]->GetNpoints();
1329 
1330  Array<OneD, NekDouble> tmp(nq, 0.0);
1331  Array<OneD, NekDouble> tmpderiv(nq);
1332 
1333  outarray = Array<OneD, NekDouble>(nq, 0.0);
1334  for (j = 0; j < m_shapedim; ++j)
1335  {
1336  for (k = 0; k < m_spacedim; ++k)
1337  {
1338  // Compute d e^m / d \xi_1 and d e^m / d \xi_2
1339  Vmath::Vcopy(nq, &m_movingframes[indm][k * nq], 1, &tmp[0], 1);
1340  m_fields[0]->PhysDirectionalDeriv(m_movingframes[j], tmp, tmpderiv);
1341 
1342  Vmath::Vmul(nq, &physarray[j + 1][0], 1, &tmpderiv[0], 1,
1343  &tmpderiv[0], 1);
1344 
1345  Vmath::Vvtvp(nq, &tmpderiv[0], 1, &m_movingframes[indk][k * nq], 1,
1346  &outarray[0], 1, &outarray[0], 1);
1347  }
1348  }
1349 }
1350 
1351 /**
1352  * @brief Compute the projection for the linear advection equation.
1353  *
1354  * @param inarray Given fields.
1355  * @param outarray Calculated solution.
1356  * @param time Time.
1357  */
1359  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1360  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
1361 {
1362  switch (m_projectionType)
1363  {
1365  {
1366  ConservativeToPrimitive(inarray, outarray);
1367  SetBoundaryConditions(outarray, time);
1368  PrimitiveToConservative(outarray, outarray);
1369  }
1370  break;
1371  default:
1372  ASSERTL0(false, "Unknown projection scheme");
1373  break;
1374  }
1375 }
1376 
1378  NekDouble time)
1379 {
1380 
1381  int nvariables = m_fields.size();
1382  int cnt = 0;
1383 
1384  // loop over Boundary Regions
1385  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
1386  {
1387 
1388  // Zonal Boundary Condition
1389  if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() == "eMG")
1390  {
1391  if (m_expdim == 1)
1392  {
1393  ASSERTL0(false, "Illegal dimension");
1394  }
1395  else if (m_expdim == 2)
1396  {
1397  // ZonalBoundary2D(n,cnt,inarray);
1398  }
1399  }
1400 
1401  // Wall Boundary Condition
1402  if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() == "eWall")
1403  {
1404  if (m_expdim == 1)
1405  {
1406  ASSERTL0(false, "Illegal dimension");
1407  }
1408  else if (m_expdim == 2)
1409  {
1410  WallBoundary2D(n, cnt, inarray);
1411  }
1412  }
1413 
1414  // Time Dependent Boundary Condition (specified in meshfile)
1415  if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
1416  "eTimeDependent")
1417  {
1418  for (int i = 0; i < nvariables; ++i)
1419  {
1420  m_fields[i]->EvaluateBoundaryConditions(time);
1421  }
1422  }
1423  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
1424  }
1425 }
1426 
1427 void MMFSWE::WallBoundary2D(int bcRegion, int cnt,
1428  Array<OneD, Array<OneD, NekDouble>> &physarray)
1429 {
1430 
1431  int i;
1432  int nTraceNumPoints = GetTraceTotPoints();
1433  int nvariables = physarray.size();
1434 
1435  // get physical values of the forward trace
1436  Array<OneD, Array<OneD, NekDouble>> Fwd0(nvariables);
1437  for (i = 0; i < nvariables; ++i)
1438  {
1439  Fwd0[i] = Array<OneD, NekDouble>(nTraceNumPoints);
1440  m_fields[i]->ExtractTracePhys(physarray[i], Fwd0[i]);
1441  }
1442 
1443  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
1444  for (i = 0; i < nvariables; ++i)
1445  {
1446  Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
1447  m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
1448  }
1449 
1450  // Adjust the physical values of the trace to take
1451  // user defined boundaries into account
1452  int e, id1, id2, npts;
1453 
1454  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
1455  ++e)
1456  {
1457  npts = m_fields[0]
1458  ->GetBndCondExpansions()[bcRegion]
1459  ->GetExp(e)
1460  ->GetTotPoints();
1461  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
1462  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
1463  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
1464 
1465  switch (m_expdim)
1466  {
1467  case 1:
1468  {
1469  // negate the forward flux
1470  Vmath::Neg(npts, &Fwd[1][id2], 1);
1471  }
1472  break;
1473  case 2:
1474  {
1475  Array<OneD, NekDouble> tmp_n(npts);
1476  Array<OneD, NekDouble> tmp_t(npts);
1477 
1478  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_ncdotMFFwd[0][id2], 1,
1479  &tmp_n[0], 1);
1480  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_ncdotMFFwd[1][id2], 1,
1481  &tmp_n[0], 1, &tmp_n[0], 1);
1482 
1483  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_nperpcdotMFFwd[0][id2], 1,
1484  &tmp_t[0], 1);
1485  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_nperpcdotMFFwd[1][id2],
1486  1, &tmp_t[0], 1, &tmp_t[0], 1);
1487 
1488  // negate the normal flux
1489  Vmath::Neg(npts, tmp_n, 1);
1490 
1491  Array<OneD, NekDouble> denom(npts);
1492  Array<OneD, NekDouble> tmp_u(npts);
1493  Array<OneD, NekDouble> tmp_v(npts);
1494 
1495  // denom = (e^1 \cdot n ) (e^2 \cdot t) - (e^2 \cdot n ) (e^1
1496  // \cdot t)
1497  Vmath::Vmul(npts, &m_ncdotMFFwd[1][id2], 1,
1498  &m_nperpcdotMFFwd[0][id2], 1, &denom[0], 1);
1499  Vmath::Vvtvm(npts, &m_ncdotMFFwd[0][id2], 1,
1500  &m_nperpcdotMFFwd[1][id2], 1, &denom[0], 1,
1501  &denom[0], 1);
1502 
1503  Vmath::Vmul(npts, &m_ncdotMFFwd[1][id2], 1, &tmp_t[0], 1,
1504  &tmp_u[0], 1);
1505  Vmath::Vvtvm(npts, &m_nperpcdotMFFwd[1][id2], 1, &tmp_n[0], 1,
1506  &tmp_u[0], 1, &tmp_u[0], 1);
1507  Vmath::Vdiv(npts, &tmp_u[0], 1, &denom[0], 1, &tmp_u[0], 1);
1508 
1509  Vmath::Vcopy(npts, &tmp_u[0], 1, &Fwd[1][id2], 1);
1510 
1511  Vmath::Vmul(npts, &m_nperpcdotMFFwd[0][id2], 1, &tmp_n[0], 1,
1512  &tmp_v[0], 1);
1513  Vmath::Vvtvm(npts, &m_ncdotMFFwd[0][id2], 1, &tmp_t[0], 1,
1514  &tmp_v[0], 1, &tmp_v[0], 1);
1515  Vmath::Vdiv(npts, &tmp_v[0], 1, &denom[0], 1, &tmp_v[0], 1);
1516 
1517  Vmath::Vcopy(npts, &tmp_v[0], 1, &Fwd[2][id2], 1);
1518  }
1519  break;
1520 
1521  default:
1522  ASSERTL0(false, "Illegal expansion dimension");
1523  }
1524 
1525  // copy boundary adjusted values into the boundary expansion
1526  for (i = 0; i < nvariables; ++i)
1527  {
1528  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
1529  &(m_fields[i]
1530  ->GetBndCondExpansions()[bcRegion]
1531  ->UpdatePhys())[id1],
1532  1);
1533  }
1534  }
1535 }
1536 
1537 /**
1538  *
1539  */
1541 {
1542  // Compute m_depth and m_Derivdepth
1544  EvaluateCoriolis();
1547 
1548  // transfer the initial conditions to modal values
1549  for (int i = 0; i < m_fields.size(); ++i)
1550  {
1551  m_fields[i]->SetPhysState(true);
1552  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1553  m_fields[i]->UpdateCoeffs());
1554  }
1555 }
1556 
1558 {
1559  int nq = GetTotPoints();
1560 
1561  m_depth = Array<OneD, NekDouble>(nq, 0.0);
1562 
1563  switch (m_TestType)
1564  {
1565  case eTestPlane:
1566  {
1567  Vmath::Fill(nq, 1.0, m_depth, 1);
1568  }
1569  break;
1570 
1571  case eTestUnsteadyZonal:
1572  {
1573  // H_0 = k_1 - k_2 - (1/2/g) (Omega sin \phi)
1574  Array<OneD, NekDouble> x(nq);
1575  Array<OneD, NekDouble> y(nq);
1576  Array<OneD, NekDouble> z(nq);
1577 
1578  m_fields[0]->GetCoords(x, y, z);
1579 
1580  NekDouble x0j, x1j, x2j;
1581  NekDouble sin_varphi, cos_varphi, sin_theta, cos_theta;
1582 
1583  for (int j = 0; j < nq; ++j)
1584  {
1585  x0j = x[j];
1586  x1j = y[j];
1587  x2j = z[j];
1588 
1589  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi,
1590  sin_theta, cos_theta);
1591  m_depth[j] =
1592  m_H0 - m_k2 -
1593  (0.5 / m_g) * (m_Omega * sin_theta) * (m_Omega * sin_theta);
1594  }
1595  }
1596  break;
1597 
1598  case eTestIsolatedMountain:
1599  {
1600  // H_0 = k_1 - k_2 - (1/2/g) (Omega sin \phi)
1601  Array<OneD, NekDouble> x(nq);
1602  Array<OneD, NekDouble> y(nq);
1603  Array<OneD, NekDouble> z(nq);
1604 
1605  m_fields[0]->GetCoords(x, y, z);
1606 
1607  int indx = 0;
1608  NekDouble x0j, x1j, x2j, dist2;
1609  NekDouble phi, theta, sin_varphi, cos_varphi, sin_theta, cos_theta;
1610  NekDouble hRad, phic, thetac;
1611  NekDouble Tol = 0.000001;
1612 
1613  hRad = m_pi / 9.0;
1614  phic = -m_pi / 2.0;
1615  thetac = m_pi / 6.0;
1616 
1617  for (int j = 0; j < nq; ++j)
1618  {
1619  x0j = x[j];
1620  x1j = y[j];
1621  x2j = z[j];
1622 
1623  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi,
1624  sin_theta, cos_theta);
1625 
1626  if ((std::abs(sin(phic) - sin_varphi) +
1627  std::abs(sin(thetac) - sin_theta)) < Tol)
1628  {
1629  std::cout << "A point " << j
1630  << " is coincient with the singularity "
1631  << std::endl;
1632  indx = 1;
1633  }
1634 
1635  phi = atan2(sin_varphi, cos_varphi);
1636  theta = atan2(sin_theta, cos_theta);
1637 
1638  // Compute r
1639  dist2 = (phi - phic) * (phi - phic) +
1640  (theta - thetac) * (theta - thetac);
1641 
1642  if (dist2 > hRad * hRad)
1643  {
1644  dist2 = hRad * hRad;
1645  }
1646 
1647  m_depth[j] = m_H0 - m_hs0 * (1.0 - sqrt(dist2) / hRad);
1648  }
1649 
1650  if (!indx)
1651  {
1652  std::cout << "No point is coincident with the singularity point"
1653  << std::endl;
1654  }
1655  }
1656  break;
1657 
1658  case eTestUnstableJet:
1659  {
1660  for (int j = 0; j < nq; ++j)
1661  {
1662  m_depth[j] = m_H0;
1663  }
1664  }
1665  break;
1666 
1667  case eTestSteadyZonal:
1668  case eTestRossbyWave:
1669  {
1670  Vmath::Zero(nq, m_depth, 1);
1671  }
1672  break;
1673 
1674  default:
1675  {
1676  Vmath::Zero(nq, m_depth, 1);
1677  }
1678  break;
1679  }
1680 
1681  // Comptue \nabla m_depth \cdot e^i
1683 
1684  for (int j = 0; j < m_shapedim; j++)
1685  {
1687  m_fields[0]->PhysDirectionalDeriv(m_movingframes[j], m_depth,
1688  m_Derivdepth[j]);
1689  }
1690 
1691  std::cout << "Water Depth (m_depth) was generated with mag = "
1692  << AvgAbsInt(m_depth) << " with max. deriv = ( "
1693  << Vmath::Vamax(nq, m_Derivdepth[0], 1) << " , "
1694  << Vmath::Vamax(nq, m_Derivdepth[1], 1) << " ) " << std::endl;
1695 }
1696 
1698 {
1699  switch (m_TestType)
1700  {
1701  case eTestPlane:
1702  {
1703  GetFunction("Coriolis")->Evaluate("f", m_coriolis);
1704  }
1705  break;
1706 
1707  case eTestSteadyZonal:
1708  {
1710  }
1711  break;
1712 
1713  case eTestUnsteadyZonal:
1714  case eTestIsolatedMountain:
1715  case eTestUnstableJet:
1716  case eTestRossbyWave:
1717  {
1719  }
1720  break;
1721 
1722  default:
1723  break;
1724  }
1725 }
1726 
1728 {
1729  int nq = GetTotPoints();
1730  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
1731 
1732  Array<OneD, NekDouble> x(nq);
1733  Array<OneD, NekDouble> y(nq);
1734  Array<OneD, NekDouble> z(nq);
1735 
1736  m_fields[0]->GetCoords(x, y, z);
1737 
1738  NekDouble x0j, x1j, x2j;
1739  NekDouble tmp;
1740 
1741  outarray = Array<OneD, NekDouble>(nq, 0.0);
1742  for (int j = 0; j < nq; ++j)
1743  {
1744  x0j = x[j];
1745  x1j = y[j];
1746  x2j = z[j];
1747 
1748  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
1749  cos_theta);
1750 
1751  // H = 2 \Omega *(- \cos \phi \cos \theta \sin \alpha + \sin \theta \cos
1752  // \alpha )
1753  tmp = -1.0 * cos_varphi * cos_theta * sin(m_alpha) +
1754  sin_theta * cos(m_alpha);
1755  outarray[j] = 2.0 * m_Omega * tmp;
1756  }
1757 }
1758 
1760 {
1761  int nq = GetTotPoints();
1762 
1763  NekDouble x0j, x1j, x2j;
1764  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
1765 
1766  Array<OneD, NekDouble> x(nq);
1767  Array<OneD, NekDouble> y(nq);
1768  Array<OneD, NekDouble> z(nq);
1769 
1770  m_fields[0]->GetCoords(x, y, z);
1771 
1772  outarray = Array<OneD, NekDouble>(nq, 0.0);
1773  for (int j = 0; j < nq; ++j)
1774  {
1775  x0j = x[j];
1776  x1j = y[j];
1777  x2j = z[j];
1778 
1779  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
1780  cos_theta);
1781 
1782  outarray[j] = 2.0 * m_Omega * sin_theta;
1783  }
1784 }
1785 
1787  bool dumpInitialConditions,
1788  const int domain)
1789 {
1790  boost::ignore_unused(domain);
1791 
1792  int nq = GetTotPoints();
1793 
1794  switch (m_TestType)
1795  {
1796  case eTestPlane:
1797  {
1798  Array<OneD, NekDouble> eta0(nq);
1799  Array<OneD, NekDouble> u0(nq);
1800  Array<OneD, NekDouble> v0(nq);
1801 
1802  TestSWE2Dproblem(initialtime, 0, eta0);
1803  m_fields[0]->SetPhys(eta0);
1804 
1805  TestSWE2Dproblem(initialtime, 1, u0);
1806  m_fields[1]->SetPhys(u0);
1807 
1808  TestSWE2Dproblem(initialtime, 2, v0);
1809  m_fields[2]->SetPhys(v0);
1810 
1811  // forward transform to fill the modal coeffs
1812  for (int i = 0; i < m_fields.size(); ++i)
1813  {
1814  m_fields[i]->SetPhysState(true);
1815  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1816  m_fields[i]->UpdateCoeffs());
1817  }
1818  }
1819  break;
1820 
1821  case eTestSteadyZonal:
1822  {
1823  Array<OneD, NekDouble> eta0(nq);
1824  Array<OneD, NekDouble> u0(nq);
1825  Array<OneD, NekDouble> v0(nq);
1826  Array<OneD, NekDouble> zeta0(nq);
1827 
1828  SteadyZonalFlow(0, eta0);
1829  m_fields[0]->SetPhys(eta0);
1830 
1831  SteadyZonalFlow(1, u0);
1832  m_fields[1]->SetPhys(u0);
1833 
1834  SteadyZonalFlow(2, v0);
1835  m_fields[2]->SetPhys(v0);
1836 
1837  // ComputeVorticity(u0, v0, zeta0);
1838  m_Vorticity0 = m_fields[0]->Integral(zeta0);
1839 
1840  m_Mass0 = ComputeMass(eta0);
1841  m_Energy0 = ComputeEnergy(eta0, u0, v0);
1842  m_Enstrophy0 = ComputeEnstrophy(eta0, u0, v0);
1843 
1844  // forward transform to fill the modal coeffs
1845  for (int i = 0; i < m_fields.size(); ++i)
1846  {
1847  m_fields[i]->SetPhysState(true);
1848  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1849  m_fields[i]->UpdateCoeffs());
1850  }
1851  }
1852  break;
1853 
1854  case eTestUnsteadyZonal:
1855  {
1856  Array<OneD, NekDouble> eta0(nq);
1857  Array<OneD, NekDouble> u0(nq);
1858  Array<OneD, NekDouble> v0(nq);
1859 
1860  UnsteadyZonalFlow(0, initialtime, eta0);
1861  m_fields[0]->SetPhys(eta0);
1862 
1863  UnsteadyZonalFlow(1, initialtime, u0);
1864  m_fields[1]->SetPhys(u0);
1865 
1866  UnsteadyZonalFlow(2, initialtime, v0);
1867  m_fields[2]->SetPhys(v0);
1868 
1869  m_Mass0 = ComputeMass(eta0);
1870  m_Energy0 = ComputeEnergy(eta0, u0, v0);
1871  m_Enstrophy0 = ComputeEnstrophy(eta0, u0, v0);
1872 
1873  // forward transform to fill the modal coeffs
1874  for (int i = 0; i < m_fields.size(); ++i)
1875  {
1876  m_fields[i]->SetPhysState(true);
1877  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1878  m_fields[i]->UpdateCoeffs());
1879  }
1880  }
1881  break;
1882 
1883  case eTestIsolatedMountain:
1884  {
1885  Array<OneD, NekDouble> eta0(nq);
1886  Array<OneD, NekDouble> u0(nq);
1887  Array<OneD, NekDouble> v0(nq);
1888 
1889  IsolatedMountainFlow(0, initialtime, eta0);
1890  m_fields[0]->SetPhys(eta0);
1891 
1892  IsolatedMountainFlow(1, initialtime, u0);
1893  m_fields[1]->SetPhys(u0);
1894 
1895  IsolatedMountainFlow(2, initialtime, v0);
1896  m_fields[2]->SetPhys(v0);
1897 
1898  m_Mass0 = ComputeMass(eta0);
1899  m_Energy0 = ComputeEnergy(eta0, u0, v0);
1900  m_Enstrophy0 = ComputeEnstrophy(eta0, u0, v0);
1901 
1902  // forward transform to fill the modal coeffs
1903  for (int i = 0; i < m_fields.size(); ++i)
1904  {
1905  m_fields[i]->SetPhysState(true);
1906  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1907  m_fields[i]->UpdateCoeffs());
1908  }
1909  }
1910  break;
1911 
1912  case eTestUnstableJet:
1913  {
1914  Array<OneD, NekDouble> eta0(nq);
1915  Array<OneD, NekDouble> u0(nq);
1916  Array<OneD, NekDouble> v0(nq);
1917 
1918  UnstableJetFlow(0, initialtime, eta0);
1919  m_fields[0]->SetPhys(eta0);
1920 
1921  UnstableJetFlow(1, initialtime, u0);
1922  m_fields[1]->SetPhys(u0);
1923 
1924  UnstableJetFlow(2, initialtime, v0);
1925  m_fields[2]->SetPhys(v0);
1926 
1927  m_Mass0 = ComputeMass(eta0);
1928  m_Energy0 = ComputeEnergy(eta0, u0, v0);
1929  m_Enstrophy0 = ComputeEnstrophy(eta0, u0, v0);
1930 
1931  // forward transform to fill the modal coeffs
1932  for (int i = 0; i < m_fields.size(); ++i)
1933  {
1934  m_fields[i]->SetPhysState(true);
1935  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1936  m_fields[i]->UpdateCoeffs());
1937  }
1938  }
1939  break;
1940 
1941  case eTestRossbyWave:
1942  {
1943  Array<OneD, NekDouble> eta0(nq);
1944  Array<OneD, NekDouble> u0(nq);
1945  Array<OneD, NekDouble> v0(nq);
1946 
1947  RossbyWave(0, eta0);
1948  m_fields[0]->SetPhys(eta0);
1949 
1950  RossbyWave(1, u0);
1951  m_fields[1]->SetPhys(u0);
1952 
1953  RossbyWave(2, v0);
1954  m_fields[2]->SetPhys(v0);
1955 
1956  m_Mass0 = ComputeMass(eta0);
1957  m_Energy0 = ComputeEnergy(eta0, u0, v0);
1958  m_Enstrophy0 = ComputeEnstrophy(eta0, u0, v0);
1959 
1960  // forward transform to fill the modal coeffs
1961  for (int i = 0; i < m_fields.size(); ++i)
1962  {
1963  m_fields[i]->SetPhysState(true);
1964  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1965  m_fields[i]->UpdateCoeffs());
1966  }
1967  }
1968  break;
1969 
1970  default:
1971  break;
1972  }
1973 
1974  if (dumpInitialConditions)
1975  {
1976  // dump initial conditions to file
1977  std::string outname = m_sessionName + "_initial.chk";
1978  WriteFld(outname);
1979 
1980  outname = m_sessionName + "_initialCART.chk";
1981  Checkpoint_Output_Cartesian(outname);
1982  }
1983 }
1984 
1985 void MMFSWE::TestSWE2Dproblem(const NekDouble time, unsigned int field,
1986  Array<OneD, NekDouble> &outfield)
1987 {
1988  boost::ignore_unused(time);
1989 
1990  int nq = m_fields[0]->GetNpoints();
1991 
1992  Array<OneD, NekDouble> x0(nq);
1993  Array<OneD, NekDouble> x1(nq);
1994  Array<OneD, NekDouble> x2(nq);
1995 
1996  m_fields[0]->GetCoords(x0, x1, x2);
1997 
1998  Array<OneD, NekDouble> eta0(nq);
1999  Array<OneD, NekDouble> u0(nq);
2000  Array<OneD, NekDouble> v0(nq);
2001 
2003 
2004  for (int i = 0; i < m_spacedim; ++i)
2005  {
2006  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2007  }
2008 
2009  for (int i = 0; i < nq; ++i)
2010  {
2011  eta0[i] = (0.771 * 0.395 * 0.395 * (1.0 / cosh(0.395 * x0[i])) *
2012  (1.0 / cosh(0.395 * x0[i]))) *
2013  (3.0 + 6.0 * x1[i] * x1[i]) / (4.0) *
2014  exp(-0.5 * x1[i] * x1[i]);
2015  uvec[0][i] = (0.771 * 0.395 * 0.395 * (1.0 / cosh(0.395 * x0[i])) *
2016  (1.0 / cosh(0.395 * x0[i]))) *
2017  (-9.0 + 6.0 * x1[i] * x1[i]) / (4.0) *
2018  exp(-0.5 * x1[i] * x1[i]);
2019  uvec[1][i] = (-2.0 * 0.395 * tanh(0.395 * x0[i])) *
2020  (0.771 * 0.395 * 0.395 * (1.0 / cosh(0.395 * x0[i])) *
2021  (1.0 / cosh(0.395 * x0[i]))) *
2022  (2.0 * x1[i]) * exp(-0.5 * x1[i] * x1[i]);
2023  }
2024 
2025  u0 = CartesianToMovingframes(uvec, 0);
2026  v0 = CartesianToMovingframes(uvec, 1);
2027 
2028  switch (field)
2029  {
2030  case (0):
2031  {
2032  outfield = eta0;
2033  }
2034  break;
2035 
2036  case (1):
2037  {
2038  outfield = u0;
2039  }
2040  break;
2041 
2042  case (2):
2043  {
2044  outfield = v0;
2045  }
2046  break;
2047  }
2048 }
2049 
2050 void MMFSWE::SteadyZonalFlow(unsigned int field,
2051  Array<OneD, NekDouble> &outfield)
2052 {
2053  int nq = GetTotPoints();
2054  NekDouble uhat, vhat;
2055  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
2056  NekDouble x0j, x1j, x2j, tmp;
2057 
2058  Array<OneD, NekDouble> eta(nq, 0.0);
2059  Array<OneD, NekDouble> u(nq, 0.0);
2060  Array<OneD, NekDouble> v(nq, 0.0);
2061 
2063  for (int i = 0; i < m_spacedim; ++i)
2064  {
2065  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2066  }
2067 
2068  Array<OneD, NekDouble> x(nq);
2069  Array<OneD, NekDouble> y(nq);
2070  Array<OneD, NekDouble> z(nq);
2071 
2072  m_fields[0]->GetCoords(x, y, z);
2073 
2074  for (int j = 0; j < nq; ++j)
2075  {
2076  x0j = x[j];
2077  x1j = y[j];
2078  x2j = z[j];
2079 
2080  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
2081  cos_theta);
2082 
2083  // H = H_0 - (1/g)*(a \Omega u_0 + 0.5*u_0^2 )*(- \cos \phi \cos \theta
2084  // \sin \alpha + \sin \theta \cos \alpha )^2
2085  tmp = -1.0 * cos_varphi * cos_theta * sin(m_alpha) +
2086  sin_theta * cos(m_alpha);
2087  eta[j] = m_H0 - m_Hvar * tmp * tmp;
2088 
2089  // u = (\vec{u} \cdot e^1 )/ || e^1 ||^2 , v = (\vec{u} \cdot e^2 )/
2090  // || e^2 ||^2
2091  uhat = m_u0 * (cos_theta * cos(m_alpha) +
2092  sin_theta * cos_varphi * sin(m_alpha));
2093  vhat = -1.0 * m_u0 * sin_varphi * sin(m_alpha);
2094 
2095  uvec[0][j] = -1.0 * uhat * sin_varphi - vhat * sin_theta * cos_varphi;
2096  uvec[1][j] = uhat * cos_varphi - vhat * sin_theta * sin_varphi;
2097  uvec[2][j] = vhat * cos_theta;
2098  }
2099 
2100  // Projection of u onto the tangent plane with conserving the mag. of the
2101  // velocity.
2103 
2104  for (int i = 0; i < m_spacedim; ++i)
2105  {
2106  uvecproj[i] = Array<OneD, NekDouble>(nq, 0.0);
2107  }
2108 
2109  // u is projected on the tangent plane with preserving its length
2110  // GramSchumitz(m_surfaceNormal, uvec, uvecproj, true);
2111 
2112  // Change it to the coordinate of moving frames
2113  // CartesianToMovingframes(0,uvecproj,u);
2114  // CartesianToMovingframes(1,uvecproj,v);
2115 
2116  u = CartesianToMovingframes(uvec, 0);
2117  v = CartesianToMovingframes(uvec, 1);
2118 
2119  switch (field)
2120  {
2121  case (0):
2122  {
2123  outfield = eta;
2124  }
2125  break;
2126 
2127  case (1):
2128  {
2129  outfield = u;
2130  }
2131  break;
2132 
2133  case (2):
2134  {
2135  outfield = v;
2136  }
2137  break;
2138  }
2139 }
2140 
2142 {
2143  int nq = m_fields[0]->GetTotPoints();
2144 
2145  Array<OneD, NekDouble> tmp(nq);
2146  Vmath::Vadd(nq, eta, 1, m_depth, 1, tmp, 1);
2147 
2148  return m_fields[0]->Integral(tmp);
2149 }
2150 
2154 {
2155  int nq = m_fields[0]->GetTotPoints();
2156 
2157  Array<OneD, NekDouble> tmp(nq);
2158  Array<OneD, NekDouble> htmp(nq);
2159  Array<OneD, NekDouble> hstmp(nq);
2160 
2161  Vmath::Vmul(nq, u, 1, u, 1, tmp, 1);
2162  Vmath::Vvtvp(nq, v, 1, v, 1, tmp, 1, tmp, 1);
2163  Vmath::Vmul(nq, eta, 1, tmp, 1, tmp, 1);
2164 
2165  Vmath::Sadd(nq, m_H0, eta, 1, htmp, 1);
2166  Vmath::Vmul(nq, htmp, 1, htmp, 1, htmp, 1);
2167 
2168  Vmath::Sadd(nq, -1.0 * m_H0, m_depth, 1, hstmp, 1);
2169  Vmath::Vmul(nq, hstmp, 1, hstmp, 1, hstmp, 1);
2170 
2171  Vmath::Vsub(nq, htmp, 1, hstmp, 1, htmp, 1);
2172  Vmath::Smul(nq, m_g, htmp, 1, htmp, 1);
2173 
2174  Vmath::Vadd(nq, htmp, 1, tmp, 1, tmp, 1);
2175  Vmath::Smul(nq, 0.5, tmp, 1, tmp, 1);
2176 
2177  return m_fields[0]->Integral(tmp);
2178 }
2179 
2183 {
2184  int nq = m_fields[0]->GetTotPoints();
2185 
2186  Array<OneD, NekDouble> hstartmp(nq);
2187  Array<OneD, NekDouble> tmp(nq);
2188 
2189  Vmath::Vadd(nq, eta, 1, m_depth, 1, hstartmp, 1);
2190 
2193  for (int i = 0; i < m_spacedim; ++i)
2194  {
2195  uvec[i] = Array<OneD, NekDouble>(nq);
2196  Curlu[i] = Array<OneD, NekDouble>(nq);
2197  }
2198 
2199  ComputeVorticity(u, v, tmp);
2200 
2201  Vmath::Vadd(nq, m_coriolis, 1, tmp, 1, tmp, 1);
2202  Vmath::Vmul(nq, tmp, 1, tmp, 1, tmp, 1);
2203  Vmath::Vdiv(nq, tmp, 1, hstartmp, 1, tmp, 1);
2204  Vmath::Smul(nq, 0.5, tmp, 1, tmp, 1);
2205 
2206  return m_fields[0]->Integral(tmp);
2207 }
2208 
2209 // Vorticity = \nabla v \cdot e^1 + v \nabla \cdot e^1 - ( \nabla u \cdot e^2 +
2210 // u \nabla \cdot e^2 )
2213  Array<OneD, NekDouble> &Vorticity)
2214 {
2215  int nq = m_fields[0]->GetTotPoints();
2216 
2217  Array<OneD, NekDouble> tmp(nq);
2218 
2219  Vorticity = Array<OneD, NekDouble>(nq, 0.0);
2220 
2221  m_fields[0]->PhysDirectionalDeriv(m_movingframes[0], v, Vorticity);
2222  Vmath::Vvtvp(nq, &v[0], 1, &m_CurlMF[1][2][0], 1, &Vorticity[0], 1,
2223  &Vorticity[0], 1);
2224 
2225  m_fields[0]->PhysDirectionalDeriv(m_movingframes[1], u, tmp);
2226  Vmath::Neg(nq, tmp, 1);
2227  Vmath::Vvtvp(nq, &u[0], 1, &m_CurlMF[0][2][0], 1, &tmp[0], 1, &tmp[0], 1);
2228 
2229  Vmath::Vadd(nq, tmp, 1, Vorticity, 1, Vorticity, 1);
2230 }
2231 
2233 {
2234  int nq = m_fields[0]->GetNpoints();
2235 
2236  Array<OneD, NekDouble> velcoeff(nq, 0.0);
2237 
2238  Array<OneD, NekDouble> Dtmp0(nq);
2239  Array<OneD, NekDouble> Dtmp1(nq);
2240  Array<OneD, NekDouble> Dtmp2(nq);
2241  Array<OneD, NekDouble> Drv(nq);
2242 
2243  vellc = Array<OneD, NekDouble>(nq, 0.0);
2244 
2245  // m_vellc = \nabla m_vel \cdot tan_i
2246  Array<OneD, NekDouble> tmp(nq);
2247  Array<OneD, NekDouble> vessel(nq);
2248 
2249  for (int j = 0; j < m_shapedim; ++j)
2250  {
2251  Vmath::Zero(nq, velcoeff, 1);
2252  for (int k = 0; k < m_spacedim; ++k)
2253  {
2254  // a_j = tan_j cdot m_vel
2255  Vmath::Vvtvp(nq, &m_movingframes[j][k * nq], 1, &m_velocity[k][0],
2256  1, &velcoeff[0], 1, &velcoeff[0], 1);
2257  }
2258 
2259  // d a_j / d x^k
2260  m_fields[0]->PhysDeriv(velcoeff, Dtmp0, Dtmp1, Dtmp2);
2261 
2262  for (int k = 0; k < m_spacedim; ++k)
2263  {
2264  // tan_j^k ( d a_j / d x^k )
2265  switch (k)
2266  {
2267  case (0):
2268  {
2269  Vmath::Vvtvp(nq, &Dtmp0[0], 1, &m_movingframes[j][k * nq],
2270  1, &vellc[0], 1, &vellc[0], 1);
2271  }
2272  break;
2273 
2274  case (1):
2275  {
2276  Vmath::Vvtvp(nq, &Dtmp1[0], 1, &m_movingframes[j][k * nq],
2277  1, &vellc[0], 1, &vellc[0], 1);
2278  }
2279  break;
2280 
2281  case (2):
2282  {
2283  Vmath::Vvtvp(nq, &Dtmp2[0], 1, &m_movingframes[j][k * nq],
2284  1, &vellc[0], 1, &vellc[0], 1);
2285  }
2286  break;
2287  }
2288  }
2289  }
2290 }
2291 
2292 void MMFSWE::UnsteadyZonalFlow(unsigned int field, const NekDouble time,
2293  Array<OneD, NekDouble> &outfield)
2294 {
2295  int nq = GetTotPoints();
2296  NekDouble uhat, vhat;
2297  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
2298  NekDouble x0j, x1j, x2j, tmp;
2299 
2300  NekDouble TR, Ttheta;
2301 
2302  Array<OneD, NekDouble> eta(nq, 0.0);
2303  Array<OneD, NekDouble> u(nq, 0.0);
2304  Array<OneD, NekDouble> v(nq, 0.0);
2305 
2307  for (int i = 0; i < m_spacedim; ++i)
2308  {
2309  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2310  }
2311 
2312  Array<OneD, NekDouble> x(nq);
2313  Array<OneD, NekDouble> y(nq);
2314  Array<OneD, NekDouble> z(nq);
2315 
2316  m_fields[0]->GetCoords(x, y, z);
2317 
2318  for (int j = 0; j < nq; ++j)
2319  {
2320  x0j = x[j];
2321  x1j = y[j];
2322  x2j = z[j];
2323 
2324  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
2325  cos_theta);
2326 
2327  // \eta = ( - ( u_0 ( - T_R sin \alpha \cos \theta + \cos \alpha \sin
2328  // \theta ) + \Omega \sin \theta )^2 + \Omega \sin \theta )^2 + (\Omega
2329  // \sin \theta )^2
2330  TR =
2331  cos_varphi * cos(m_Omega * time) - sin_varphi * sin(m_Omega * time);
2332  Ttheta =
2333  sin_varphi * cos(m_Omega * time) + cos_varphi * sin(m_Omega * time);
2334  tmp = -1.0 * TR * sin(m_alpha) * cos_theta + cos(m_alpha) * sin_theta;
2335 
2336  eta[j] = -1.0 * (m_u0 * tmp + m_Omega * sin_theta) *
2337  (m_u0 * tmp + m_Omega * sin_theta) +
2338  m_Omega * m_Omega * sin_theta * sin_theta;
2339  eta[j] = 0.5 * eta[j] / m_g;
2340 
2341  // u = u_0*(TR*\sin \alpha * \sin \theta + \cos \alpha * \cos \theta
2342  // v = - u_0 Ttheta * \sin \alpha
2343  uhat =
2344  m_u0 * (TR * sin(m_alpha) * sin_theta + cos(m_alpha) * cos_theta);
2345  vhat = -1.0 * m_u0 * Ttheta * sin(m_alpha);
2346 
2347  uvec[0][j] = -1.0 * uhat * sin_varphi - vhat * sin_theta * cos_varphi;
2348  uvec[1][j] = uhat * cos_varphi - vhat * sin_theta * sin_varphi;
2349  uvec[2][j] = vhat * cos_theta;
2350  }
2351 
2352  // Projection of u onto the tangent plane with conserving the mag. of the
2353  // velocity.
2355 
2356  for (int i = 0; i < m_spacedim; ++i)
2357  {
2358  uvecproj[i] = Array<OneD, NekDouble>(nq, 0.0);
2359  }
2360 
2361  // u is projected on the tangent plane with preserving its length
2362  // GramSchumitz(m_surfaceNormal, uvec, uvecproj, true);
2363 
2364  // Change it to the coordinate of moving frames
2365  // CartesianToMovingframes(0,uvecproj,u);
2366  // CartesianToMovingframes(1,uvecproj,v);
2367 
2368  u = CartesianToMovingframes(uvec, 0);
2369  v = CartesianToMovingframes(uvec, 1);
2370 
2371  switch (field)
2372  {
2373  case (0):
2374  {
2375  outfield = eta;
2376  }
2377  break;
2378 
2379  case (1):
2380  {
2381  outfield = u;
2382  }
2383  break;
2384 
2385  case (2):
2386  {
2387  outfield = v;
2388  }
2389  break;
2390  }
2391 }
2392 
2393 void MMFSWE::IsolatedMountainFlow(unsigned int field, const NekDouble time,
2394  Array<OneD, NekDouble> &outfield)
2395 {
2396  boost::ignore_unused(time);
2397 
2398  int nq = GetTotPoints();
2399 
2400  NekDouble uhat, vhat;
2401  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
2402  NekDouble x0j, x1j, x2j;
2403 
2404  Array<OneD, NekDouble> eta(nq, 0.0);
2405  Array<OneD, NekDouble> u(nq, 0.0);
2406  Array<OneD, NekDouble> v(nq, 0.0);
2407 
2409  for (int i = 0; i < m_spacedim; ++i)
2410  {
2411  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2412  }
2413 
2414  Array<OneD, NekDouble> x(nq);
2415  Array<OneD, NekDouble> y(nq);
2416  Array<OneD, NekDouble> z(nq);
2417 
2418  m_fields[0]->GetCoords(x, y, z);
2419 
2420  for (int j = 0; j < nq; ++j)
2421  {
2422  x0j = x[j];
2423  x1j = y[j];
2424  x2j = z[j];
2425 
2426  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
2427  cos_theta);
2428 
2429  // eta = - (1/g) (\Omega u_0 + 0.4 u^2 ) \sin^2 \theta
2430  eta[j] = (-1.0 / m_g) * (m_Omega * m_u0 + 0.5 * m_u0 * m_u0) *
2431  sin_theta * sin_theta;
2432 
2433  uhat = m_u0 * cos_theta;
2434  vhat = 0.0;
2435 
2436  uvec[0][j] = -1.0 * uhat * sin_varphi - vhat * sin_theta * cos_varphi;
2437  uvec[1][j] = uhat * cos_varphi - vhat * sin_theta * sin_varphi;
2438  uvec[2][j] = vhat * cos_theta;
2439  }
2440 
2441  // Projection of u onto the tangent plane with conserving the mag. of the
2442  // velocity.
2444 
2445  for (int i = 0; i < m_spacedim; ++i)
2446  {
2447  uvecproj[i] = Array<OneD, NekDouble>(nq, 0.0);
2448  }
2449 
2450  // u is projected on the tangent plane with preserving its length
2451  // GramSchumitz(m_surfaceNormal, uvec, uvecproj, true);
2452 
2453  // Change it to the coordinate of moving frames
2454  // CartesianToMovingframes(0,uvecproj,u);
2455  // CartesianToMovingframes(1,uvecproj,v);
2456 
2457  u = CartesianToMovingframes(uvec, 0);
2458  v = CartesianToMovingframes(uvec, 1);
2459 
2460  switch (field)
2461  {
2462  case (0):
2463  {
2464  outfield = eta;
2465  }
2466  break;
2467 
2468  case (1):
2469  {
2470  outfield = u;
2471  }
2472  break;
2473 
2474  case (2):
2475  {
2476  outfield = v;
2477  }
2478  break;
2479  }
2480 }
2481 
2482 void MMFSWE::UnstableJetFlow(unsigned int field, const NekDouble time,
2483  Array<OneD, NekDouble> &outfield)
2484 {
2485  boost::ignore_unused(time);
2486 
2487  int nq = GetTotPoints();
2488 
2489  NekDouble uhat, vhat;
2490  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
2491  NekDouble x0j, x1j, x2j;
2492  NekDouble Ttheta, Tphi;
2493 
2494  Array<OneD, NekDouble> eta(nq, 0.0);
2495  Array<OneD, NekDouble> u(nq, 0.0);
2496  Array<OneD, NekDouble> v(nq, 0.0);
2497 
2499  for (int i = 0; i < m_spacedim; ++i)
2500  {
2501  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2502  }
2503 
2504  Array<OneD, NekDouble> x(nq);
2505  Array<OneD, NekDouble> y(nq);
2506  Array<OneD, NekDouble> z(nq);
2507 
2508  m_fields[0]->GetCoords(x, y, z);
2509 
2510  int Nint = 1000;
2511  NekDouble dth, intj;
2512  for (int j = 0; j < nq; ++j)
2513  {
2514  x0j = x[j];
2515  x1j = y[j];
2516  x2j = z[j];
2517 
2518  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
2519  cos_theta);
2520 
2521  Ttheta = atan2(sin_theta, cos_theta);
2522  Tphi = atan2(sin_varphi, cos_varphi);
2523 
2524  uhat = ComputeUnstableJetuphi(Ttheta);
2525  vhat = 0.0;
2526 
2527  // eta = - (1/g) (\Omega u_0 + 0.4 u^2 ) \sin^2 \theta
2528  dth = Ttheta / Nint;
2529  eta[j] = dth * 0.5 *
2531  for (int i = 1; i < Nint - 1; i++)
2532  {
2533  intj = i * dth;
2534  eta[j] = eta[j] + dth * ComputeUnstableJetEta(intj);
2535  }
2536  eta[j] = (-1.0 / m_g) * eta[j];
2537 
2538  // Add perturbation
2539  if (m_PurturbedJet)
2540  {
2541  eta[j] = eta[j] + m_hbar * cos_theta * exp(-9.0 * Tphi * Tphi) *
2542  exp(-225.0 * (m_pi / 4.0 - Ttheta) *
2543  (m_pi / 4.0 - Ttheta));
2544  }
2545 
2546  uvec[0][j] = -1.0 * uhat * sin_varphi - vhat * sin_theta * cos_varphi;
2547  uvec[1][j] = uhat * cos_varphi - vhat * sin_theta * sin_varphi;
2548  uvec[2][j] = vhat * cos_theta;
2549  }
2550 
2551  // Projection of u onto the tangent plane with conserving the mag. of the
2552  // velocity.
2554 
2555  for (int i = 0; i < m_spacedim; ++i)
2556  {
2557  uvecproj[i] = Array<OneD, NekDouble>(nq, 0.0);
2558  }
2559 
2560  // u is projected on the tangent plane with preserving its length
2561  // GramSchumitz(m_surfaceNormal, uvec, uvecproj, true);
2562 
2563  // Change it to the coordinate of moving frames
2564  // CartesianToMovingframes(0,uvecproj,u);
2565  // CartesianToMovingframes(1,uvecproj,v);
2566 
2567  u = CartesianToMovingframes(uvec, 0);
2568  v = CartesianToMovingframes(uvec, 1);
2569 
2570  switch (field)
2571  {
2572  case (0):
2573  {
2574  outfield = eta;
2575  }
2576  break;
2577 
2578  case (1):
2579  {
2580  outfield = u;
2581  }
2582  break;
2583 
2584  case (2):
2585  {
2586  outfield = v;
2587  }
2588  break;
2589  }
2590 }
2591 
2592 void MMFSWE::RossbyWave(unsigned int field, Array<OneD, NekDouble> &outfield)
2593 {
2594  int nq = GetTotPoints();
2595  NekDouble uhat, vhat;
2596  NekDouble sin_theta, cos_theta, sin_varphi, cos_varphi;
2597  NekDouble x0j, x1j, x2j;
2598  NekDouble Ath, Bth, Cth, tmp;
2599 
2600  Array<OneD, NekDouble> eta(nq, 0.0);
2601  Array<OneD, NekDouble> u(nq, 0.0);
2602  Array<OneD, NekDouble> v(nq, 0.0);
2603 
2605  for (int i = 0; i < m_spacedim; ++i)
2606  {
2607  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2608  }
2609 
2610  Array<OneD, NekDouble> x(nq);
2611  Array<OneD, NekDouble> y(nq);
2612  Array<OneD, NekDouble> z(nq);
2613 
2614  m_fields[0]->GetCoords(x, y, z);
2615 
2616  NekDouble R = 4.0;
2617  NekDouble cos2theta, cosRtheta, cos6theta, cos2Rtheta, cosRm1theta;
2618  NekDouble cos2phi, cos4phi, sin4phi, cos8phi;
2619 
2620  // disturbancees of Rossby-Haurwitz Wave
2621  NekDouble x0d, y0d, z0d, phi0, theta0;
2622  // NekDouble rad_earth = 6.37122 * 1000000;
2623 
2624  phi0 = 40.0 * m_pi / 180.0;
2625  theta0 = 50.0 * m_pi / 180.0;
2626 
2627  x0d = cos(phi0) * cos(theta0);
2628  y0d = sin(phi0) * cos(theta0);
2629  z0d = sin(theta0);
2630 
2631  for (int j = 0; j < nq; ++j)
2632  {
2633  x0j = x[j];
2634  x1j = y[j];
2635  x2j = z[j];
2636 
2637  CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
2638  cos_theta);
2639 
2640  // H = H_0 - (1/g)*(a \Omega u_0 + 0.5*u_0^2 )*(- \cos \phi \cos \theta
2641  // \sin \alpha + \sin \theta \cos \alpha )^2
2642 
2643  // tmp = cos^{2R} \theta, R = 4;
2644  cos2theta = cos_theta * cos_theta;
2645  cosRm1theta = cos_theta * cos2theta;
2646  cosRtheta = cos2theta * cos2theta;
2647  cos6theta = cos2theta * cosRtheta;
2648  cos2Rtheta = cosRtheta * cosRtheta;
2649 
2650  tmp = (0.5 * m_angfreq) * (2.0 * m_Omega + m_angfreq);
2651  Ath = tmp * cos2theta +
2652  0.25 * m_K * m_K * cos6theta *
2653  ((R + 1.0) * cosRtheta + (2 * R * R - R - 2.0) * cos2theta -
2654  2 * R * R);
2655 
2656  tmp = (2.0 * m_K) * (m_Omega + m_angfreq) / ((R + 1.0) * (R + 2.0));
2657  Bth = tmp * cosRtheta *
2658  ((R * R + 2 * R + 2) - (R + 1.0) * (R + 1.0) * cos2theta);
2659 
2660  Cth =
2661  0.25 * m_K * m_K * cos2Rtheta * ((R + 1.0) * cos2theta - (R + 2.0));
2662 
2663  // cos (2 \phi) = 2 * cos \phi * cos \phi - 1.0
2664  cos2phi = 2.0 * cos_varphi * cos_varphi - 1.0;
2665  cos4phi = 2.0 * cos2phi * cos2phi - 1.0;
2666  cos8phi = 2.0 * cos4phi * cos4phi - 1.0;
2667 
2668  // sin (2 \phi) = 2 * cos \phi * sin \phi
2669  sin4phi = 4.0 * sin_varphi * cos_varphi * cos2phi;
2670 
2671  eta[j] = m_H0 + (1.0 / m_g) * (Ath + Bth * cos4phi + Cth * cos8phi);
2672 
2673  // disturbances is added
2674  if (m_RossbyDisturbance)
2675  {
2676  eta[j] = eta[j] *
2677  (1.0 + (1.0 / 40.0) * (x0j * x0d + x1j * y0d + x2j * z0d));
2678  }
2679 
2680  // u = (\vec{u} \cdot e^1 )/ || e^1 ||^2 , v = (\vec{u} \cdot e^2 )/
2681  // || e^2 ||^2
2682  uhat = m_angfreq * cos_theta +
2683  m_K * cosRm1theta * (R * sin_theta * sin_theta - cos2theta) *
2684  cos4phi;
2685  vhat = -1.0 * m_K * R * cosRm1theta * sin_theta * sin4phi;
2686 
2687  uvec[0][j] = -1.0 * uhat * sin_varphi - vhat * sin_theta * cos_varphi;
2688  uvec[1][j] = uhat * cos_varphi - vhat * sin_theta * sin_varphi;
2689  uvec[2][j] = vhat * cos_theta;
2690  }
2691 
2692  // NekDouble etamin, etaaver;
2693  // etamin = Vmath::Vmin(nq, eta, 1)*rad_earth;
2694  // etaaver = Vmath::Vsum(nq, eta, 1)/nq*rad_earth;
2695 
2696  // Projection of u onto the tangent plane with conserving the mag. of the
2697  // velocity.
2699 
2700  for (int i = 0; i < m_spacedim; ++i)
2701  {
2702  uvecproj[i] = Array<OneD, NekDouble>(nq, 0.0);
2703  }
2704 
2705  // u is projected on the tangent plane with preserving its length
2706  // GramSchumitz(m_surfaceNormal, uvec, uvecproj, true);
2707 
2708  // Change it to the coordinate of moving frames
2709  // CartesianToMovingframes(0,uvecproj,u);
2710  // CartesianToMovingframes(1,uvecproj,v);
2711 
2712  u = CartesianToMovingframes(uvec, 0);
2713  v = CartesianToMovingframes(uvec, 1);
2714 
2715  switch (field)
2716  {
2717  case (0):
2718  {
2719  outfield = eta;
2720  }
2721  break;
2722 
2723  case (1):
2724  {
2725  outfield = u;
2726  }
2727  break;
2728 
2729  case (2):
2730  {
2731  outfield = v;
2732  }
2733  break;
2734  }
2735 }
2736 
2738 {
2739  NekDouble uphi, f, dh;
2740 
2741  uphi = ComputeUnstableJetuphi(theta);
2742  f = 2.0 * m_Omega * sin(theta);
2743 
2744  dh = f * uphi + tan(theta) * uphi * uphi;
2745 
2746  return dh;
2747 }
2748 
2750 {
2751  NekDouble uphi;
2752 
2753  if ((theta > m_theta0) && (theta < m_theta1))
2754  {
2755  uphi = (m_uthetamax / m_en) *
2756  exp(1.0 / (theta - m_theta0) / (theta - m_theta1));
2757  }
2758 
2759  else
2760  {
2761  uphi = 0.0;
2762  }
2763 
2764  return uphi;
2765 }
2766 
2767 void MMFSWE::Checkpoint_Output_Cartesian(std::string outname)
2768 {
2769  int nq = m_fields[0]->GetTotPoints();
2770  int ncoeffs = m_fields[0]->GetNcoeffs();
2771 
2772  NekDouble rad_earth = 6.37122 * 1000000;
2773 
2774  int nvariables = 7;
2775 
2776  // vector u in Cartesian coordinates
2777  std::vector<std::string> variables(nvariables);
2778 
2779  variables[0] = "eta";
2780  variables[1] = "hstar";
2781  variables[2] = "vorticity";
2782  variables[3] = "ux";
2783  variables[4] = "uy";
2784  variables[5] = "uz";
2785  variables[6] = "null";
2786 
2787  // Obtain \vec{u} in cartesian coordinate
2788  Array<OneD, Array<OneD, NekDouble>> fieldphys(nvariables);
2789  std::vector<Array<OneD, NekDouble>> fieldcoeffs(nvariables);
2790  for (int i = 0; i < nvariables; ++i)
2791  {
2792  fieldphys[i] = Array<OneD, NekDouble>(nq, 0.0);
2793  fieldcoeffs[i] = Array<OneD, NekDouble>(ncoeffs);
2794  }
2795 
2796  Vmath::Smul(nq, rad_earth, &(m_fields[0]->GetPhys())[0], 1,
2797  &fieldphys[0][0], 1);
2798  // Vmath::Vadd(nq, &(m_fields[0]->GetPhys())[0], 1, &m_depth[0], 1,
2799  // &fieldphys[1][0], 1);
2800 
2801  Vmath::Vcopy(nq, &m_depth[0], 1, &fieldphys[1][0], 1);
2802  Vmath::Neg(nq, &fieldphys[1][0], 1);
2803  Vmath::Sadd(nq, m_H0, &fieldphys[1][0], 1, &fieldphys[1][0], 1);
2804  Vmath::Smul(nq, rad_earth, &fieldphys[1][0], 1, &fieldphys[1][0], 1);
2805 
2806  Array<OneD, NekDouble> utmp(nq);
2807  Array<OneD, NekDouble> vtmp(nq);
2808 
2809  Vmath::Vcopy(nq, &(m_fields[1]->GetPhys())[0], 1, &utmp[0], 1);
2810  Vmath::Vcopy(nq, &(m_fields[2]->GetPhys())[0], 1, &vtmp[0], 1);
2811 
2812  ComputeVorticity(utmp, vtmp, fieldphys[2]);
2813  // u_x = u e^1_x + v e^2_x
2814  int indx = 3;
2815  for (int k = 0; k < m_spacedim; ++k)
2816  {
2817  Vmath::Vmul(nq, &utmp[0], 1, &m_movingframes[0][k * nq], 1,
2818  &fieldphys[k + indx][0], 1);
2819  Vmath::Vvtvp(nq, &vtmp[0], 1, &m_movingframes[1][k * nq], 1,
2820  &fieldphys[k + indx][0], 1, &fieldphys[k + indx][0], 1);
2821  }
2822 
2823  for (int i = 0; i < nvariables; ++i)
2824  {
2825  m_fields[0]->FwdTrans(fieldphys[i], fieldcoeffs[i]);
2826  }
2827 
2828  WriteFld(outname, m_fields[0], fieldcoeffs, variables);
2829 }
2830 
2831 // (h, hu, hv) -> (\eta, u, v)
2833  const Array<OneD, const Array<OneD, NekDouble>> &physin,
2834  Array<OneD, Array<OneD, NekDouble>> &physout)
2835 {
2836  int nq = GetTotPoints();
2837 
2838  if (physin[0].get() == physout[0].get())
2839  {
2840  // copy indata and work with tmp array
2842  for (int i = 0; i < 3; ++i)
2843  {
2844  // deep copy
2845  tmp[i] = Array<OneD, NekDouble>(nq);
2846  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
2847  }
2848 
2849  // \eta = h - d
2850  Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
2851 
2852  // u = hu/h
2853  Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
2854 
2855  // v = hv/h
2856  Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
2857  }
2858  else
2859  {
2860  // \eta = h - d
2861  Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
2862 
2863  // u = hu/h
2864  Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
2865 
2866  // v = hv/h
2867  Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
2868  }
2869 }
2870 
2872  const Array<OneD, const Array<OneD, NekDouble>> &physin,
2873  Array<OneD, Array<OneD, NekDouble>> &physout)
2874 {
2875 
2876  int nq = GetTotPoints();
2877 
2878  if (physin[0].get() == physout[0].get())
2879  {
2880  // copy indata and work with tmp array
2882  for (int i = 0; i < 3; ++i)
2883  {
2884  // deep copy
2885  tmp[i] = Array<OneD, NekDouble>(nq);
2886  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
2887  }
2888 
2889  // h = \eta + d
2890  Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
2891 
2892  // hu = h * u
2893  Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
2894 
2895  // hv = h * v
2896  Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
2897  }
2898  else
2899  {
2900  // h = \eta + d
2901  Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
2902 
2903  // hu = h * u
2904  Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
2905 
2906  // hv = h * v
2907  Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
2908  }
2909 }
2910 
2911 // To obtain [eta, u, v ] from [H, Hu, Hv]
2913 {
2914  int nq = GetTotPoints();
2915 
2916  // u = hu/h
2917  Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
2918  m_fields[1]->UpdatePhys(), 1);
2919 
2920  // v = hv/ v
2921  Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
2922  m_fields[2]->UpdatePhys(), 1);
2923 
2924  // \eta = h - d
2925  Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
2926  m_fields[0]->UpdatePhys(), 1);
2927 }
2928 
2930 {
2931  int nq = GetTotPoints();
2932 
2933  // h = \eta + d
2934  Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
2935  m_fields[0]->UpdatePhys(), 1);
2936 
2937  // hu = h * u
2938  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
2939  m_fields[1]->UpdatePhys(), 1);
2940 
2941  // hv = h * v
2942  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
2943  m_fields[2]->UpdatePhys(), 1);
2944 }
2945 
2947 {
2948  // Construct beta
2949  int i, k;
2950  int n = 1, m = 1;
2951  int nq = m_fields[0]->GetTotPoints();
2952 
2953  NekDouble alpha, beta_theta, beta_phi;
2954 
2955  NekDouble xp, yp, zp, Re;
2956  NekDouble theta, phi, sin_theta, cos_theta, sin_varphi, cos_varphi;
2957  NekDouble cosntheta3;
2958 
2959  NekDouble thetax, thetay, thetaz;
2960  NekDouble phix, phiy, phiz;
2961 
2962  Array<OneD, NekDouble> x(nq);
2963  Array<OneD, NekDouble> y(nq);
2964  Array<OneD, NekDouble> z(nq);
2965 
2966  Array<OneD, NekDouble> u(nq);
2967  Array<OneD, NekDouble> v(nq);
2968 
2969  Array<OneD, NekDouble> vorticitycompt(nq);
2970  Array<OneD, NekDouble> vorticityexact(nq);
2971 
2972  m_fields[0]->GetCoords(x, y, z);
2973 
2975  for (i = 0; i < m_spacedim; ++i)
2976  {
2977  uvec[i] = Array<OneD, NekDouble>(nq, 0.0);
2978  }
2979 
2980  // Generate SphericalCoords
2981  for (k = 0; k < nq; ++k)
2982  {
2983  xp = x[k];
2984  yp = y[k];
2985  zp = z[k];
2986 
2987  Re = sqrt(xp * xp + yp * yp + zp * zp);
2988 
2989  CartesianToSpherical(xp, yp, zp, sin_varphi, cos_varphi, sin_theta,
2990  cos_theta);
2991 
2992  alpha = sin_varphi;
2993 
2994  theta = atan2(sin_theta, cos_theta);
2995  phi = atan2(sin_varphi, cos_varphi);
2996 
2997  cosntheta3 = cos(n * theta) * cos(n * theta) * cos(n * theta);
2998 
2999  beta_theta = -4.0 * n * cosntheta3 * cos(m * phi) * sin(n * theta) / Re;
3000  beta_phi = -m * cosntheta3 * sin(m * phi) / Re;
3001 
3002  thetax = -1.0 * cos_varphi * sin_theta;
3003  thetay = -1.0 * sin_varphi * sin_theta;
3004  thetaz = cos_theta;
3005 
3006  phix = -1.0 * sin_varphi;
3007  phiy = cos_varphi;
3008  phiz = 0.0;
3009 
3010  uvec[0][k] = alpha * (beta_theta * thetax + beta_phi * phix);
3011  uvec[1][k] = alpha * (beta_theta * thetay + beta_phi * phiy);
3012  uvec[2][k] = alpha * (beta_theta * thetaz + beta_phi * phiz);
3013 
3014  vorticityexact[k] = -4.0 * n / Re / Re * cos_theta * cos_theta *
3015  cos_varphi * cos(m * phi) * sin(n * theta);
3016  }
3017 
3018  u = CartesianToMovingframes(uvec, 0);
3019  v = CartesianToMovingframes(uvec, 1);
3020 
3021  std::cout << "chi migi1" << std::endl;
3022 
3023  ComputeVorticity(u, v, vorticitycompt);
3024  /*for (int k=0; k < nq; k++)
3025  {
3026 
3027  std::cout << "vorticitycompt[ " << k << "]"<< "\t"<<vorticitycompt[k]<<
3028  std::endl;
3029 
3030  }*/
3031 
3032  Vmath::Vsub(nq, vorticityexact, 1, vorticitycompt, 1, vorticitycompt, 1);
3033 
3034  std::cout << "Vorticity: L2 error = " << AvgAbsInt(vorticitycompt)
3035  << ", Linf error = " << Vmath::Vamax(nq, vorticitycompt, 1)
3036  << std::endl;
3037 }
3038 
3039 NekDouble MMFSWE::v_L2Error(unsigned int field,
3040  const Array<OneD, NekDouble> &exactsoln,
3041  bool Normalised)
3042 {
3043  boost::ignore_unused(exactsoln);
3044 
3045  int nq = m_fields[field]->GetNpoints();
3046  NekDouble L2error = -1.0;
3047 
3048  if (m_NumQuadPointsError == 0)
3049  {
3050  if (m_fields[field]->GetPhysState() == false)
3051  {
3052  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
3053  m_fields[field]->UpdatePhys());
3054  }
3055 
3056  switch (field)
3057  {
3058  case (0):
3059  {
3060  // I(h) = I (h - h_exact ) / I (h_exact)
3061  Array<OneD, NekDouble> exactsolution(nq);
3062  v_EvaluateExactSolution(0, exactsolution, m_time);
3063 
3064  // exactsoln = u - u_T so that L2 compute u_T
3065  NekDouble L2exact = m_fields[0]->Integral(exactsolution);
3066 
3067  Vmath::Vsub(nq, &(m_fields[0]->GetPhys())[0], 1,
3068  &exactsolution[0], 1, &exactsolution[0], 1);
3069  Vmath::Vabs(nq, exactsolution, 1, exactsolution, 1);
3070 
3071  L2error = (m_fields[0]->Integral(exactsolution)) / L2exact;
3072  }
3073  break;
3074 
3075  case (1):
3076  {
3077  // I2 (u) = I( (u - u_ext)^2 + (v - v_ext)^2 )^{1/2} / I(
3078  // u_ext^2 + v_ext^2 )^{1/2}
3079  Array<OneD, NekDouble> exactu(nq);
3080  Array<OneD, NekDouble> exactv(nq);
3081  Array<OneD, NekDouble> tmp(nq);
3082 
3083  // L2exact = \int (\sqrt{exactu*exactu+exactv*exactv})
3084  NekDouble L2exact;
3085  v_EvaluateExactSolution(1, exactu, m_time);
3086  v_EvaluateExactSolution(2, exactv, m_time);
3087  Vmath::Vmul(nq, exactu, 1, exactu, 1, tmp, 1);
3088  Vmath::Vvtvp(nq, exactv, 1, exactv, 1, tmp, 1, tmp, 1);
3089  Vmath::Vsqrt(nq, tmp, 1, tmp, 1);
3090 
3091  L2exact = m_fields[1]->Integral(tmp);
3092 
3093  // L2exact = \int
3094  // (\sqrt{(u-exactu)*(u-exactu)+(v-exactv)*(v-exactv)})
3095  Vmath::Vsub(nq, &(m_fields[1]->GetPhys())[0], 1, &exactu[0], 1,
3096  &exactu[0], 1);
3097  Vmath::Vsub(nq, &(m_fields[2]->GetPhys())[0], 1, &exactv[0], 1,
3098  &exactv[0], 1);
3099  Vmath::Vmul(nq, exactu, 1, exactu, 1, tmp, 1);
3100  Vmath::Vvtvp(nq, exactv, 1, exactv, 1, tmp, 1, tmp, 1);
3101  Vmath::Vsqrt(nq, tmp, 1, tmp, 1);
3102 
3103  L2error = (m_fields[1]->Integral(tmp)) / L2exact;
3104  }
3105  break;
3106 
3107  case (2):
3108  {
3109  L2error = 0.0;
3110  }
3111  break;
3112 
3113  default:
3114  break;
3115  }
3116 
3117  if (Normalised == true)
3118  {
3119  Array<OneD, NekDouble> one(m_fields[field]->GetNpoints(), 1.0);
3120 
3121  NekDouble Vol = m_fields[field]->Integral(one);
3122  m_comm->AllReduce(Vol, LibUtilities::ReduceSum);
3123 
3124  L2error = sqrt(L2error * L2error / Vol);
3125  }
3126  }
3127  else
3128  {
3129  Array<OneD, NekDouble> L2INF(2);
3130  L2INF = ErrorExtraPoints(field);
3131  L2error = L2INF[0];
3132  }
3133 
3134  return L2error;
3135 }
3136 
3137 NekDouble MMFSWE::v_LinfError(unsigned int field,
3138  const Array<OneD, NekDouble> &exactsoln)
3139 {
3140  boost::ignore_unused(exactsoln);
3141 
3142  NekDouble LinfError = -1.0;
3143 
3144  if (m_fields[field]->GetPhysState() == false)
3145  {
3146  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
3147  m_fields[field]->UpdatePhys());
3148  }
3149 
3150  int nq = m_fields[field]->GetNpoints();
3151 
3152  // Obtain \vec{u} in cartesian coordinate
3153  Array<OneD, NekDouble> x(nq);
3154  Array<OneD, NekDouble> y(nq);
3155  Array<OneD, NekDouble> z(nq);
3156 
3157  m_fields[0]->GetCoords(x, y, z);
3158 
3159  switch (field)
3160  {
3161  case (0):
3162  {
3163  NekDouble Letaint;
3164 
3165  Array<OneD, NekDouble> exactsolution(nq);
3166 
3167  EvaluateExactSolution(field, exactsolution, m_time);
3168  LinfError = m_fields[field]->Linf(m_fields[field]->GetPhys(),
3169  exactsolution);
3170 
3171  Letaint = Vmath::Vamax(nq, exactsolution, 1);
3172 
3173  Vmath::Vsub(nq, &(m_fields[0]->GetPhys())[0], 1, &exactsolution[0],
3174  1, &exactsolution[0], 1);
3175 
3176  LinfError = fabs(LinfError / Letaint);
3177  }
3178  break;
3179 
3180  case (1):
3181  {
3182  Array<OneD, NekDouble> exactu(nq);
3183  Array<OneD, NekDouble> exactv(nq);
3184  Array<OneD, NekDouble> tmpu(nq);
3185  Array<OneD, NekDouble> tmpv(nq);
3186  Array<OneD, NekDouble> Lerr(nq);
3187  Array<OneD, NekDouble> uT(nq);
3188 
3189  EvaluateExactSolution(1, exactu, m_time);
3190  EvaluateExactSolution(2, exactv, m_time);
3191 
3192  // Compute max[sqrt{(u-uex)^2 + (v-vex)^2}]
3193  Vmath::Vcopy(nq, &(m_fields[1]->UpdatePhys())[0], 1, &tmpu[0], 1);
3194  Vmath::Vcopy(nq, &(m_fields[2]->UpdatePhys())[0], 1, &tmpv[0], 1);
3195 
3196  Vmath::Vsub(nq, &exactu[0], 1, &tmpu[0], 1, &tmpu[0], 1);
3197  Vmath::Vsub(nq, &exactv[0], 1, &tmpv[0], 1, &tmpv[0], 1);
3198 
3199  Vmath::Vmul(nq, &tmpu[0], 1, &tmpu[0], 1, &tmpu[0], 1);
3200  Vmath::Vmul(nq, &tmpv[0], 1, &tmpv[0], 1, &tmpv[0], 1);
3201 
3202  Vmath::Vadd(nq, &tmpu[0], 1, &tmpv[0], 1, &Lerr[0], 1);
3203  Vmath::Vsqrt(nq, &Lerr[0], 1, &Lerr[0], 1);
3204 
3205  // uT = max[sqrt( u_T^2 + v_T^2 ) ]
3206  Vmath::Vmul(nq, &exactu[0], 1, &exactu[0], 1, &tmpu[0], 1);
3207  Vmath::Vmul(nq, &exactv[0], 1, &exactv[0], 1, &tmpv[0], 1);
3208  Vmath::Vadd(nq, &tmpu[0], 1, &tmpv[0], 1, &uT[0], 1);
3209  Vmath::Vsqrt(nq, &uT[0], 1, &uT[0], 1);
3210 
3211  LinfError = Vmath::Vamax(nq, Lerr, 1) / Vmath::Vamax(nq, uT, 1);
3212  }
3213  break;
3214 
3215  case (2):
3216  {
3217  LinfError = 0.0;
3218  }
3219  break;
3220 
3221  default:
3222  break;
3223  }
3224 
3225  return LinfError;
3226 }
3227 
3228 void MMFSWE::v_EvaluateExactSolution(unsigned int field,
3229  Array<OneD, NekDouble> &outfield,
3230  const NekDouble time)
3231 {
3232  switch (m_TestType)
3233  {
3234  case eTestSteadyZonal:
3235  {
3236  SteadyZonalFlow(field, outfield);
3237  }
3238  break;
3239 
3240  case eTestUnsteadyZonal:
3241  {
3242  UnsteadyZonalFlow(field, time, outfield);
3243  }
3244  break;
3245 
3246  case eTestIsolatedMountain:
3247  {
3248  IsolatedMountainFlow(field, time, outfield);
3249  }
3250  break;
3251 
3252  case eTestUnstableJet:
3253  {
3254  UnstableJetFlow(field, time, outfield);
3255  }
3256  break;
3257 
3258  case eTestRossbyWave:
3259  {
3260  RossbyWave(field, outfield);
3261  }
3262  break;
3263 
3264  case eTestPlane:
3265  {
3266  TestSWE2Dproblem(time, field, outfield);
3267  }
3268  break;
3269 
3270  default:
3271  break;
3272  }
3273 }
3274 
3276 {
3279 
3280  switch (m_TestType)
3281  {
3282  case eTestSteadyZonal:
3283  {
3284  SolverUtils::AddSummaryItem(s, "Rotation Angle", m_alpha);
3285  }
3286  break;
3287 
3288  case eTestRossbyWave:
3289  {
3290  SolverUtils::AddSummaryItem(s, "RossbyDistrubance",
3292  }
3293  break;
3294 
3295  case eTestUnstableJet:
3296  {
3297  SolverUtils::AddSummaryItem(s, "PurturbedJet", m_PurturbedJet);
3298  }
3299  break;
3300 
3301  default:
3302  break;
3303  }
3304 }
3305 
3306 } // end 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
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:68
void ComputeVorticity(const Array< OneD, const NekDouble > &u, const Array< OneD, const NekDouble > &v, Array< OneD, NekDouble > &Vorticity)
Definition: MMFSWE.cpp:2211
void AddCoriolis(Array< OneD, Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: MMFSWE.cpp:1201
void UnstableJetFlow(unsigned int field, const NekDouble time, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:2482
virtual void v_InitObject(bool DeclareFields=true) override
Initialise the object.
Definition: MMFSWE.cpp:62
int m_PurturbedJet
Definition: MMFSWE.h:276
NekDouble m_k2
Definition: MMFSWE.h:96
NekDouble ComputeUnstableJetuphi(const NekDouble theta)
Definition: MMFSWE.cpp:2749
MMFSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
Definition: MMFSWE.cpp:52
NekDouble m_alpha
Definition: MMFSWE.h:95
void AddElevationEffect(Array< OneD, Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: MMFSWE.cpp:1246
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: MMFSWE.h:67
virtual ~MMFSWE()
Destructor.
Definition: MMFSWE.cpp:245
NekDouble m_en
Definition: MMFSWE.h:98
void RossbyWave(unsigned int field, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:2592
void ComputeNablaCdotVelocity(Array< OneD, NekDouble > &vellc)
Definition: MMFSWE.cpp:2232
virtual NekDouble v_L2Error(unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised) override
Virtual function for the L_2 error computation between fields and a given exact solution.
Definition: MMFSWE.cpp:3039
NekDouble m_Energy0
Definition: MMFSWE.h:101
void GetSWEFluxVector(const int i, const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &flux)
Definition: MMFSWE.cpp:582
void AddDivForGradient(Array< OneD, Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: MMFSWE.cpp:552
void WeakDGSWEDirDeriv(const Array< OneD, Array< OneD, NekDouble >> &InField, Array< OneD, Array< OneD, NekDouble >> &OutField)
Definition: MMFSWE.cpp:482
void EvaluateStandardCoriolis(Array< OneD, NekDouble > &outarray)
Definition: MMFSWE.cpp:1759
static std::string className
Name of class.
Definition: MMFSWE.h:77
NekDouble m_theta0
Definition: MMFSWE.h:98
void PrimitiveToConservative()
Definition: MMFSWE.cpp:2929
virtual void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time) override
Definition: MMFSWE.cpp:3228
void IsolatedMountainFlow(unsigned int field, const NekDouble time, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:2393
NekDouble m_hbar
Definition: MMFSWE.h:98
NekDouble m_K
Definition: MMFSWE.h:99
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Definition: MMFSWE.h:90
void ComputeMagAndDot(const int index, NekDouble &MageF1, NekDouble &MageF2, NekDouble &MageB1, NekDouble &MageB2, NekDouble &eF1_cdot_eB1, NekDouble &eF1_cdot_eB2, NekDouble &eF2_cdot_eB1, NekDouble &eF2_cdot_eB2)
Definition: MMFSWE.cpp:1163
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the RHS.
Definition: MMFSWE.cpp:426
NekDouble m_Hvar
Definition: MMFSWE.h:95
void Compute_demdt_cdot_ek(const int indm, const int indk, const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, NekDouble > &outarray)
Definition: MMFSWE.cpp:1322
TestType m_TestType
Definition: MMFSWE.h:79
void AverageFlux(const int index, NekDouble hL, NekDouble uL, NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR, Array< OneD, NekDouble > &numfluxF, Array< OneD, NekDouble > &numfluxB)
Definition: MMFSWE.cpp:960
NekDouble m_Vorticity0
Definition: MMFSWE.h:101
void Checkpoint_Output_Cartesian(std::string outname)
Definition: MMFSWE.cpp:2767
void Computehhuhvflux(NekDouble hL, NekDouble uL, NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR, NekDouble hstar, NekDouble &hflux, NekDouble &huflux, NekDouble &hvflux)
Definition: MMFSWE.cpp:894
virtual NekDouble v_LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln) override
Virtual function for the L_inf error computation between fields and a given exact solution.
Definition: MMFSWE.cpp:3137
void TestVorticityComputation()
Definition: MMFSWE.cpp:2946
void RiemannSolverHLLC(const int index, NekDouble hL, NekDouble uL, NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR, Array< OneD, NekDouble > &numfluxF, Array< OneD, NekDouble > &numfluxB)
Definition: MMFSWE.cpp:848
virtual void v_DoInitialise() override
Sets up initial conditions.
Definition: MMFSWE.cpp:1540
NekDouble ComputeEnstrophy(const Array< OneD, const NekDouble > &eta, const Array< OneD, const NekDouble > &u, const Array< OneD, const NekDouble > &v)
Definition: MMFSWE.cpp:2180
void NumericalSWEFlux(Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &numfluxFwd, Array< OneD, Array< OneD, NekDouble >> &numfluxBwd)
Definition: MMFSWE.cpp:677
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
Definition: MMFSWE.cpp:3275
virtual void v_SetInitialConditions(const NekDouble initialtime, bool dumpInitialConditions, const int domain) override
Definition: MMFSWE.cpp:1786
NekDouble m_uthetamax
Definition: MMFSWE.h:98
NekDouble ComputeMass(const Array< OneD, const NekDouble > &eta)
Definition: MMFSWE.cpp:2141
void EvaluateWaterDepth(void)
Definition: MMFSWE.cpp:1557
void UnsteadyZonalFlow(unsigned int field, const NekDouble time, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:2292
NekDouble ComputeEnergy(const Array< OneD, const NekDouble > &eta, const Array< OneD, const NekDouble > &u, const Array< OneD, const NekDouble > &v)
Definition: MMFSWE.cpp:2151
void LaxFriedrichFlux(const int index, NekDouble hL, NekDouble uL, NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR, Array< OneD, NekDouble > &numfluxF, Array< OneD, NekDouble > &numfluxB)
Definition: MMFSWE.cpp:1014
NekDouble ComputeUnstableJetEta(const NekDouble theta)
Definition: MMFSWE.cpp:2737
NekDouble m_Mass0
Definition: MMFSWE.h:101
void EvaluateCoriolis(void)
Definition: MMFSWE.cpp:1697
void EvaluateCoriolisForZonalFlow(Array< OneD, NekDouble > &outarray)
Definition: MMFSWE.cpp:1727
NekDouble m_theta1
Definition: MMFSWE.h:98
int m_AddRotation
Definition: MMFSWE.h:92
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
Definition: MMFSWE.h:107
int m_planeNumber
Definition: MMFSWE.h:115
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &physarray)
Definition: MMFSWE.cpp:1427
virtual void v_DoSolve() override
Solves an unsteady problem.
Definition: MMFSWE.cpp:249
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble >> &inarray, NekDouble time)
Definition: MMFSWE.cpp:1377
NekDouble m_angfreq
Definition: MMFSWE.h:99
void TestSWE2Dproblem(const NekDouble time, unsigned int field, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:1985
NekDouble m_Omega
Definition: MMFSWE.h:95
int m_RossbyDisturbance
Definition: MMFSWE.h:275
NekDouble m_H0
Definition: MMFSWE.h:95
int m_AddCoriolis
Definition: MMFSWE.h:92
void AddRotation(Array< OneD, Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
Definition: MMFSWE.cpp:1275
Array< OneD, Array< OneD, NekDouble > > m_Derivdepth
Definition: MMFSWE.h:87
NekDouble m_u0
Definition: MMFSWE.h:95
void ConservativeToPrimitive()
Definition: MMFSWE.cpp:2912
NekDouble m_Enstrophy0
Definition: MMFSWE.h:101
void SteadyZonalFlow(unsigned int field, Array< OneD, NekDouble > &outfield)
Definition: MMFSWE.cpp:2050
void RusanovFlux(const int index, NekDouble hL, NekDouble uL, NekDouble vL, NekDouble hR, NekDouble uR, NekDouble vR, Array< OneD, NekDouble > &numfluxF, Array< OneD, NekDouble > &numfluxB)
Definition: MMFSWE.cpp:1089
Array< OneD, NekDouble > m_depth
Still water depth.
Definition: MMFSWE.h:86
NekDouble m_g
Definition: MMFSWE.h:94
NekDouble m_hs0
Definition: MMFSWE.h:97
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the projection.
Definition: MMFSWE.cpp:1358
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
int m_expdim
Expansion dimension.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_time
Current time of simulation.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT NekDouble LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
Linf error computation.
SOLVER_UTILS_EXPORT void EvaluateExactSolution(int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
Evaluates an exact solution.
NekDouble m_fintime
Finish time of the simulation.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
NekDouble m_checktime
Time between checkpoints.
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT int GetExpSize()
std::string m_sessionName
Name of the session.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > ErrorExtraPoints(unsigned int field)
Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
int m_steps
Number of steps to take.
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT int GetTotPoints()
A base class for PDEs which include an advection component.
Definition: MMFSystem.h:147
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > CartesianToMovingframes(const Array< OneD, const Array< OneD, NekDouble >> &uvec, unsigned int field)
Definition: MMFSystem.cpp:774
Array< OneD, Array< OneD, NekDouble > > m_DivMF
Definition: MMFSystem.h:195
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFFwd
Definition: MMFSystem.h:192
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
Definition: MMFSystem.cpp:2469
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFBwd
Definition: MMFSystem.h:193
SOLVER_UTILS_EXPORT NekDouble AvgAbsInt(const Array< OneD, const NekDouble > &inarray)
Definition: MMFSystem.cpp:2313
Array< OneD, Array< OneD, NekDouble > > m_movingframes
Definition: MMFSystem.h:186
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceFwd
Definition: MMFSystem.h:199
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceBwd
Definition: MMFSystem.h:200
SOLVER_UTILS_EXPORT void MMFInitObject(const Array< OneD, const Array< OneD, NekDouble >> &Anisotropy, const int TangentXelem=-1)
Definition: MMFSystem.cpp:53
SOLVER_UTILS_EXPORT void CopyBoundaryTrace(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, const BoundaryCopyType BDCopyType, const int var=0, const std::string btype="NoUserDefined")
Definition: MMFSystem.cpp:838
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFBwd
Definition: MMFSystem.h:190
SOLVER_UTILS_EXPORT void CartesianToSpherical(const NekDouble x0j, const NekDouble x1j, const NekDouble x2j, NekDouble &sin_varphi, NekDouble &cos_varphi, NekDouble &sin_theta, NekDouble &cos_theta)
Definition: MMFSystem.cpp:795
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFFwd
Definition: MMFSystem.h:189
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_CurlMF
Definition: MMFSystem.h:196
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static const NekDouble kNekZeroTol
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
EquationSystemFactory & GetEquationSystemFactory()
@ eAverage
averaged (or centred) flux
Definition: MMFSystem.h:80
@ eLaxFriedrich
Lax-Friedrich flux.
Definition: MMFSystem.h:81
@ eHLLC
Harten-Lax-Leer Contact wave flux.
Definition: MMFSystem.h:85
@ eRusanov
Upwind.
Definition: MMFSystem.h:83
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
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
@ eTestUnstableJet
Definition: MMFSWE.h:50
@ eTestSteadyZonal
Definition: MMFSWE.h:47
@ eTestPlane
Definition: MMFDiffusion.h:48
@ eTestUnsteadyZonal
Definition: MMFSWE.h:48
@ SIZE_TestType
Length of enum list.
Definition: MMFDiffusion.h:55
@ eTestIsolatedMountain
Definition: MMFSWE.h:49
@ eTestRossbyWave
Definition: MMFSWE.h:51
const char *const TestTypeMap[]
Definition: MMFDiffusion.h:58
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:534
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 Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:553
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 Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
T Vamax(int n, const T *x, const int incx)
Return the maximum absolute element in x called vamax to avoid conflict with max.
Definition: Vmath.cpp:999
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384
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
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294