Nektar++
LinearSWE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LinearSWE.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Linear Shallow water equations in primitive variables
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36#include <iostream>
37
40
41namespace Nektar
42{
43std::string LinearSWE::className =
45 "LinearSWE", LinearSWE::create,
46 "Linear shallow water equation in primitive variables.");
47
50 : ShallowWaterSystem(pSession, pGraph)
51{
52}
53
54void LinearSWE::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(&LinearSWE::GetFluxVector, this);
90
91 // Setting up Riemann solver for advection operator
92 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
93 if ((riemName == "LinearHLL") && (m_constantDepth != true))
94 {
95 ASSERTL0(false, "LinearHLL only valid for constant depth");
96 }
99 riemName, m_session);
100
101 // Setting up parameters for advection operator Riemann solver
102 m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
103 m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
104 m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
105
106 // The numerical flux for linear SWE requires depth information
107 int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
108 m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
109 m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
110 m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
112 m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
113 this);
114 m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
115 this);
116
117 // Concluding initialisation of advection operators
118 m_advection->SetRiemannSolver(m_riemannSolver);
119 m_advection->InitObject(m_session, m_fields);
120 break;
121 }
122 default:
123 {
124 ASSERTL0(false, "Unsupported projection type.");
125 break;
126 }
127 }
128}
129
131{
133 if (m_session->DefinesSolverInfo("UpwindType"))
134 {
135 std::string UpwindType;
136 UpwindType = m_session->GetSolverInfo("UpwindType");
137 if (UpwindType == "LinearAverage")
138 {
139 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
140 }
141 if (UpwindType == "LinearHLL")
142 {
143 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
144 }
145 }
146 SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
147 SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
148 SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
149}
150
152 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
153 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
154{
155 int i, j;
156 int ndim = m_spacedim;
157 int nvariables = inarray.size();
158 int nq = GetTotPoints();
159
160 switch (m_projectionType)
161 {
163 {
164
165 //-------------------------------------------------------
166 // Compute the DG advection including the numerical flux
167 // by using SolverUtils/Advection
168 // Input and output in physical space
170
171 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
172 time);
173 //-------------------------------------------------------
174
175 //-------------------------------------------------------
176 // negate the outarray since moving terms to the rhs
177 for (i = 0; i < nvariables; ++i)
178 {
179 Vmath::Neg(nq, outarray[i], 1);
180 }
181 //-------------------------------------------------------
182
183 //-------------------------------------------------
184 // Add "source terms"
185 // Input and output in physical space
186
187 // Coriolis forcing
188 if (m_coriolis.size() != 0)
189 {
190 AddCoriolis(inarray, outarray);
191 }
192 //-------------------------------------------------
193 }
194 break;
197 {
198
199 //-------------------------------------------------------
200 // Compute the fluxvector in physical space
202 nvariables);
203
204 for (i = 0; i < nvariables; ++i)
205 {
206 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
207 for (j = 0; j < ndim; ++j)
208 {
209 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
210 }
211 }
212
213 LinearSWE::GetFluxVector(inarray, fluxvector);
214 //-------------------------------------------------------
215
216 //-------------------------------------------------------
217 // Take the derivative of the flux terms
218 // and negate the outarray since moving terms to the rhs
220 Array<OneD, NekDouble> tmp1(nq);
221
222 for (i = 0; i < nvariables; ++i)
223 {
225 fluxvector[i][0], tmp);
227 fluxvector[i][1], tmp1);
228 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
229 Vmath::Neg(nq, outarray[i], 1);
230 }
231
232 //-------------------------------------------------
233 // Add "source terms"
234 // Input and output in physical space
235
236 // Coriolis forcing
237 if (m_coriolis.size() != 0)
238 {
239 AddCoriolis(inarray, outarray);
240 }
241 //-------------------------------------------------
242 }
243 break;
244 default:
245 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
246 break;
247 }
248}
249
250// Physfield in primitive Form
252 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
254{
255 int i, j;
256 int nq = m_fields[0]->GetTotPoints();
257
258 NekDouble g = m_g;
259
260 // Flux vector for the mass equation
261 for (i = 0; i < m_spacedim; ++i)
262 {
263 Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
264 }
265
266 // Put (g eta) in tmp
268 Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
269
270 // Flux vector for the momentum equations
271 for (i = 0; i < m_spacedim; ++i)
272 {
273 for (j = 0; j < m_spacedim; ++j)
274 {
275 // must zero fluxes as not initialised to zero in AdvectionWeakDG
276 // ...
277 Vmath::Zero(nq, flux[i + 1][j], 1);
278 }
279
280 // Add (g eta) to appropriate field
281 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
282 }
283}
284
285/**
286 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
287 * \f$ h\mathbf{v} \f$.
288 *
289 * @param physfield Velocity field.
290 * @param velocity Velocity field.
291 */
293 const Array<OneD, Array<OneD, NekDouble>> &physfield,
295{
296 const int npts = physfield[0].size();
297
298 for (int i = 0; i < m_spacedim; ++i)
299 {
300 Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
301 }
302}
303
304} // 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: LinearSWE.h:59
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Definition: LinearSWE.cpp:130
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
Definition: LinearSWE.cpp:292
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:151
LinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: LinearSWE.cpp:48
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: LinearSWE.h:49
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:87
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:82
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
Definition: LinearSWE.cpp:54
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:251
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:94
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:95
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
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.
Array< OneD, NekDouble > m_depth
Still water depth.
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.
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825