Nektar++
ShallowWaterSystem.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ShallowWaterSystem.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: Generic timestepping for shallow water solvers
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iostream>
36
38
39using namespace std;
40
41namespace Nektar
42{
43/**
44 * @class ShallowWaterSystem
45 *
46 * Provides the underlying timestepping framework for shallow water flow solvers
47 * including the general timestepping routines. This class is not intended
48 * to be directly instantiated, but rather is a base class on which to
49 * define shallow water solvers, e.g. SWE, Boussinesq, linear and nonlinear
50 * versions.
51 *
52 * For details on implementing unsteady solvers see
53 * \ref sectionADRSolverModuleImplementation
54 */
55
56/**
57 * Processes SolverInfo parameters from the session file and sets up
58 * timestepping-specific code.
59 * @param pSession Session object to read parameters from.
60 */
61
64 "ShallowWaterSystem", ShallowWaterSystem::create,
65 "Auxiliary functions for the shallow water system.");
66
70 : UnsteadySystem(pSession, pGraph)
71{
72}
73
74void ShallowWaterSystem::v_InitObject(bool DeclareFields)
75{
76 UnsteadySystem::v_InitObject(DeclareFields);
77
78 // if discontinuous Galerkin determine numerical flux to use
80 {
81 ASSERTL0(m_session->DefinesSolverInfo("UPWINDTYPE"),
82 "No UPWINDTYPE defined in session.");
83 }
84
85 // Set up locations of velocity vector.
88 for (int i = 0; i < m_spacedim; ++i)
89 {
90 m_vecLocs[0][i] = 1 + i;
91 }
92
93 // Load acceleration of gravity
94 m_session->LoadParameter("Gravity", m_g, 9.81);
95
96 // input/output in primitive variables
97 m_primitive = true;
98
100
101 m_constantDepth = true;
102 NekDouble depth = m_depth[0];
103 for (int i = 0; i < GetTotPoints(); ++i)
104 {
105 if (m_depth[i] != depth)
106 {
107 m_constantDepth = false;
108 break;
109 }
110 }
111
112 // Compute the bottom slopes
113 int nq = GetTotPoints();
114 if (m_constantDepth != true)
115 {
117 for (int i = 0; i < m_spacedim; ++i)
118 {
121 m_bottomSlope[i]);
122 Vmath::Neg(nq, m_bottomSlope[i], 1);
123 }
124 }
125
127}
128
130{
132 if (m_constantDepth == true)
133 {
134 SolverUtils::AddSummaryItem(s, "Depth", "constant");
135 }
136 else
137 {
138 SolverUtils::AddSummaryItem(s, "Depth", "variable");
139 }
140}
141
143 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
144 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
145{
146 int i;
147 int nvariables = inarray.size();
148
149 switch (m_projectionType)
150 {
152 {
153
154 // Just copy over array
155 if (inarray != outarray)
156 {
157 int npoints = GetNpoints();
158
159 for (i = 0; i < nvariables; ++i)
160 {
161 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
162 }
163 }
164
165 SetBoundaryConditions(outarray, time);
166 break;
167 }
170 {
171
173 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
174
175 for (i = 0; i < nvariables; ++i)
176 {
177 m_fields[i]->FwdTrans(inarray[i], coeffs);
178 m_fields[i]->BwdTrans(coeffs, outarray[i]);
179 }
180 break;
181 }
182 default:
183 ASSERTL0(false, "Unknown projection scheme");
184 break;
185 }
186}
187
188//----------------------------------------------------
191{
192 std::string varName;
193 int nvariables = 3;
194 int cnt = 0;
195 int nTracePts = GetTraceTotPoints();
196
197 // Extract trace for boundaries. Needs to be done on all processors to avoid
198 // deadlock.
199 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
200 for (int i = 0; i < nvariables; ++i)
201 {
202 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
203 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
204 }
205
206 // Loop over Boundary Regions
207 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
208 {
209 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
211 {
212 continue;
213 }
214
215 // Wall Boundary Condition
216 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
217 "Wall"))
218 {
219 WallBoundary2D(n, cnt, Fwd);
220 }
221
222 // Time Dependent Boundary Condition (specified in meshfile)
223 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
224 {
225 for (int i = 0; i < nvariables; ++i)
226 {
227 varName = m_session->GetVariable(i);
228 m_fields[i]->EvaluateBoundaryConditions(time, varName);
229 }
230 }
231 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
232 }
233}
234
236 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd)
237{
238 int i;
239 int nvariables = 3;
240
241 // Adjust the physical values of the trace to take
242 // user defined boundaries into account
243 int e, id1, id2, npts;
244
245 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
246 ++e)
247 {
248 npts = m_fields[0]
249 ->GetBndCondExpansions()[bcRegion]
250 ->GetExp(e)
251 ->GetNumPoints(0);
252 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
253 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
254 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
255
256 switch (m_expdim)
257 {
258 case 1:
259 {
260 // negate the forward flux
261 Vmath::Neg(npts, &Fwd[1][id2], 1);
262 }
263 break;
264 case 2:
265 {
266 Array<OneD, NekDouble> tmp_n(npts);
267 Array<OneD, NekDouble> tmp_t(npts);
268
269 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
270 &tmp_n[0], 1);
271 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
272 &tmp_n[0], 1, &tmp_n[0], 1);
273
274 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
275 &tmp_t[0], 1);
276 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
277 &tmp_t[0], 1, &tmp_t[0], 1);
278
279 // negate the normal flux
280 Vmath::Neg(npts, tmp_n, 1);
281
282 // rotate back to Cartesian
283 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
284 &Fwd[1][id2], 1);
285 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
286 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
287
288 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
289 &Fwd[2][id2], 1);
290 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
291 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
292 }
293 break;
294 case 3:
295 ASSERTL0(false,
296 "3D not implemented for Shallow Water Equations");
297 break;
298 default:
299 ASSERTL0(false, "Illegal expansion dimension");
300 }
301
302 // copy boundary adjusted values into the boundary expansion
303 for (i = 0; i < nvariables; ++i)
304 {
305 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
306 &(m_fields[i]
307 ->GetBndCondExpansions()[bcRegion]
308 ->UpdatePhys())[id1],
309 1);
310 }
311 }
312}
313
315 int bcRegion, int cnt, Array<OneD, Array<OneD, NekDouble>> &Fwd,
316 Array<OneD, Array<OneD, NekDouble>> &physarray)
317{
318 int i;
319 int nvariables = physarray.size();
320
321 // Adjust the physical values of the trace to take
322 // user defined boundaries into account
323 int e, id1, id2, npts;
324
325 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
326 ++e)
327 {
328 npts = m_fields[0]
329 ->GetBndCondExpansions()[bcRegion]
330 ->GetExp(e)
331 ->GetTotPoints();
332 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
333 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
334 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
335
336 // For 2D/3D, define: v* = v - 2(v.n)n
337 Array<OneD, NekDouble> tmp(npts, 0.0);
338
339 // Calculate (v.n)
340 for (i = 0; i < m_spacedim; ++i)
341 {
342 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
343 &tmp[0], 1, &tmp[0], 1);
344 }
345
346 // Calculate 2.0(v.n)
347 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
348
349 // Calculate v* = v - 2.0(v.n)n
350 for (i = 0; i < m_spacedim; ++i)
351 {
352 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
353 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
354 }
355
356 // copy boundary adjusted values into the boundary expansion
357 for (i = 0; i < nvariables; ++i)
358 {
359 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
360 &(m_fields[i]
361 ->GetBndCondExpansions()[bcRegion]
362 ->UpdatePhys())[id1],
363 1);
364 }
365 }
366}
367
368// physarray contains the conservative variables
370 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
372{
373
374 int ncoeffs = GetNcoeffs();
375 int nq = GetTotPoints();
376
378 Array<OneD, NekDouble> mod(ncoeffs);
379
380 switch (m_projectionType)
381 {
383 {
384 // add to u equation
385 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
386 m_fields[0]->IProductWRTBase(tmp, mod);
387 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
388 m_fields[0]->BwdTrans(mod, tmp);
389 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
390
391 // add to v equation
392 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
393 Vmath::Neg(nq, tmp, 1);
394 m_fields[0]->IProductWRTBase(tmp, mod);
395 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
396 m_fields[0]->BwdTrans(mod, tmp);
397 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
398 }
399 break;
402 {
403 // add to u equation
404 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
405 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
406
407 // add to v equation
408 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
409 Vmath::Neg(nq, tmp, 1);
410 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
411 }
412 break;
413 default:
414 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
415 break;
416 }
417}
418
420{
421 int nq = GetTotPoints();
422
423 // u = hu/h
424 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
425 m_fields[1]->UpdatePhys(), 1);
426
427 // v = hv/ v
428 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
429 m_fields[2]->UpdatePhys(), 1);
430
431 // \eta = h - d
432 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
433 m_fields[0]->UpdatePhys(), 1);
434}
435
437{
438 int nq = GetTotPoints();
439
440 // h = \eta + d
441 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
442 m_fields[0]->UpdatePhys(), 1);
443
444 // hu = h * u
445 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
446 m_fields[1]->UpdatePhys(), 1);
447
448 // hv = h * v
449 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
450 m_fields[2]->UpdatePhys(), 1);
451}
452
454{
455 GetFunction("WaterDepth")->Evaluate("d", m_depth);
456}
457
459{
460 GetFunction("Coriolis")->Evaluate("f", m_coriolis);
461}
462
465{
466
467 int cnt = 0;
468 // loop over Boundary Regions
469 for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
470 ++bcRegion)
471 {
472 if (m_fields[0]
473 ->GetBndConditions()[bcRegion]
474 ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
475 {
476 continue;
477 }
478
479 // Copy the forward trace of the field to the backward trace
480 int e, id2, npts;
481
482 for (e = 0;
483 e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
484 ++e)
485 {
486 npts = m_fields[0]
487 ->GetBndCondExpansions()[bcRegion]
488 ->GetExp(e)
489 ->GetTotPoints();
490 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
491 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
492 e));
493
494 Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
495 }
496
497 cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
498 }
499}
500
501} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
NekDouble m_g
Acceleration of gravity.
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd)
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
void WallBoundary(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
bool m_constantDepth
Indicates if constant depth case.
bool m_primitive
Indicates if variables are primitive or conservative.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
ShallowWaterSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
Array< OneD, NekDouble > m_depth
Still water depth.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
static std::string className
Name of class.
int m_spacedim
Spatial dimension (>= expansion dim).
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT int GetTotPoints()
Base class for unsteady solvers.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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.hpp:366
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.hpp:180
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.hpp:381
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.hpp:100
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.hpp:126
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.hpp:220
STL namespace.