Nektar++
NonlinearSWE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NonlinearSWE.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// 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: Nonlinear Shallow water equations in conservative variables
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36#include <iostream>
37
40
41namespace Nektar
42{
43std::string NonlinearSWE::className =
45 "NonlinearSWE", NonlinearSWE::create,
46 "Nonlinear shallow water equation in conservative variables.");
47
50 : ShallowWaterSystem(pSession, pGraph)
51{
52}
53
54void NonlinearSWE::v_InitObject(bool DeclareFields)
55{
57
59 {
62 }
63 else
64 {
65 ASSERTL0(false, "Implicit SWE not set up.");
66 }
67
68 // Type of advection class to be used
69 switch (m_projectionType)
70 {
71 // Continuous field
73 {
74 // Do nothing
75 break;
76 }
77 // Discontinuous field
79 {
80 std::string advName;
81 std::string diffName;
82 std::string riemName;
83
84 //---------------------------------------------------------------
85 // Setting up advection and diffusion operators
86 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
88 advName, advName);
89 m_advection->SetFluxVector(&NonlinearSWE::GetFluxVector, this);
90
91 // Setting up Riemann solver for advection operator
92 m_session->LoadSolverInfo("UpwindType", riemName, "Average");
95 riemName, m_session);
96
97 // Setting up parameters for advection operator Riemann solver
98 m_riemannSolver->SetParam("gravity", &NonlinearSWE::GetGravity,
99 this);
100 m_riemannSolver->SetAuxVec("vecLocs", &NonlinearSWE::GetVecLocs,
101 this);
102 m_riemannSolver->SetVector("N", &NonlinearSWE::GetNormals, this);
103 m_riemannSolver->SetScalar("depth", &NonlinearSWE::GetDepth, this);
104
105 // Concluding initialisation of advection operators
106 m_advection->SetRiemannSolver(m_riemannSolver);
107 m_advection->InitObject(m_session, m_fields);
108 break;
109 }
110 default:
111 {
112 ASSERTL0(false, "Unsupported projection type.");
113 break;
114 }
115 }
116}
117
119{
121 SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
122 SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
123 SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
124}
125
127 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
128 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
129{
130 int i, j;
131 int ndim = m_spacedim;
132 int nvariables = inarray.size();
133 int nq = GetTotPoints();
134
135 switch (m_projectionType)
136 {
138 {
139
140 //-------------------------------------------------------
141 // Compute the DG advection including the numerical flux
142 // by using SolverUtils/Advection
143 // Input and output in physical space
145
146 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
147 time);
148 //-------------------------------------------------------
149
150 //-------------------------------------------------------
151 // negate the outarray since moving terms to the rhs
152 for (i = 0; i < nvariables; ++i)
153 {
154 Vmath::Neg(nq, outarray[i], 1);
155 }
156 //-------------------------------------------------------
157
158 //-------------------------------------------------
159 // Add "source terms"
160 // Input and output in physical space
161
162 // Coriolis forcing
163 if (m_coriolis.size() != 0)
164 {
165 AddCoriolis(inarray, outarray);
166 }
167
168 // Variable Depth
169 if (m_constantDepth != true)
170 {
171 AddVariableDepth(inarray, outarray);
172 }
173 //-------------------------------------------------
174 }
175 break;
178 {
179
180 //-------------------------------------------------------
181 // Compute the fluxvector in physical space
183 nvariables);
184
185 for (i = 0; i < nvariables; ++i)
186 {
187 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
188 for (j = 0; j < ndim; ++j)
189 {
190 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
191 }
192 }
193
194 NonlinearSWE::GetFluxVector(inarray, fluxvector);
195 //-------------------------------------------------------
196
197 //-------------------------------------------------------
198 // Take the derivative of the flux terms
199 // and negate the outarray since moving terms to the rhs
201 Array<OneD, NekDouble> tmp1(nq);
202
203 for (i = 0; i < nvariables; ++i)
204 {
206 fluxvector[i][0], tmp);
208 fluxvector[i][1], tmp1);
209 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
210 Vmath::Neg(nq, outarray[i], 1);
211 }
212
213 //-------------------------------------------------
214 // Add "source terms"
215 // Input and output in physical space
216
217 // Coriolis forcing
218 if (m_coriolis.size() != 0)
219 {
220 AddCoriolis(inarray, outarray);
221 }
222
223 // Variable Depth
224 if (m_constantDepth != true)
225 {
226 AddVariableDepth(inarray, outarray);
227 }
228 //-------------------------------------------------
229 }
230 break;
231 default:
232 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
233 break;
234 }
235}
236
237// Physfield in conservative Form
239 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
241{
242 int i, j;
243 int nq = m_fields[0]->GetTotPoints();
244
245 NekDouble g = m_g;
247
248 // Flux vector for the mass equation
249 for (i = 0; i < m_spacedim; ++i)
250 {
251 velocity[i] = Array<OneD, NekDouble>(nq);
252 Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
253 }
254
255 GetVelocityVector(physfield, velocity);
256
257 // Put (0.5 g h h) in tmp
259 Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
260 Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
261
262 // Flux vector for the momentum equations
263 for (i = 0; i < m_spacedim; ++i)
264 {
265 for (j = 0; j < m_spacedim; ++j)
266 {
267 Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
268 1);
269 }
270
271 // Add (0.5 g h h) to appropriate field
272 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
273 }
274}
275
276/**
277 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
278 * \f$ h\mathbf{v} \f$.
279 *
280 * @param physfield Momentum field.
281 * @param velocity Velocity field.
282 */
284 const Array<OneD, Array<OneD, NekDouble>> &physfield,
286{
287 const int npts = physfield[0].size();
288
289 for (int i = 0; i < m_spacedim; ++i)
290 {
291 Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
292 }
293}
294
295// physarray contains the conservative variables
297 const Array<OneD, const Array<OneD, NekDouble>> &physarray,
299{
300 int ncoeffs = GetNcoeffs();
301 int nq = GetTotPoints();
302
304 Array<OneD, NekDouble> mod(ncoeffs);
305
306 switch (m_projectionType)
307 {
309 {
310 for (int i = 0; i < m_spacedim; ++i)
311 {
312 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
313 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
314 m_fields[0]->IProductWRTBase(tmp, mod);
315 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
316 m_fields[0]->BwdTrans(mod, tmp);
317 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
318 }
319 }
320 break;
323 {
324 for (int i = 0; i < m_spacedim; ++i)
325 {
326 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
327 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
328 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
329 }
330 }
331 break;
332 default:
333 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
334 break;
335 }
336}
337
338} // 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.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
static std::string className
Name of class.
Definition: NonlinearSWE.h:59
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: NonlinearSWE.h:49
NonlinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
const Array< OneD, NekDouble > & GetDepth()
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)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
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
RiemannSolverFactory & GetRiemannSolverFactory()
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 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 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