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
58 // Type of advection class to be used
59 switch (m_projectionType)
60 {
61 // Continuous field
63 {
64 // Do nothing
65 break;
66 }
67 // Discontinuous field
69 {
70 std::string advName;
71 std::string diffName;
72 std::string riemName;
73
74 //---------------------------------------------------------------
75 // Setting up advection and diffusion operators
76 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
78 advName, advName);
79 m_advection->SetFluxVector(&LinearSWE::GetFluxVector, this);
80
81 // Setting up Riemann solver for advection operator
82 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
83 if ((riemName == "LinearHLL") && (m_constantDepth != true))
84 {
85 ASSERTL0(false, "LinearHLL only valid for constant depth");
86 }
89 riemName, m_session);
90
91 // Setting up parameters for advection operator Riemann solver
92 m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
93 m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
94 m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
95
96 // The numerical flux for linear SWE requires depth information
97 int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
98 m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
99 m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
100 m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
102 m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
103 this);
104 m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
105 this);
106
107 // Concluding initialisation of advection operators
108 m_advection->SetRiemannSolver(m_riemannSolver);
109 m_advection->InitObject(m_session, m_fields);
110 break;
111 }
112 default:
113 {
114 ASSERTL0(false, "Unsupported projection type.");
115 break;
116 }
117 }
118
122}
123
125 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
126 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
127{
128 int ndim = m_spacedim;
129 int nvariables = inarray.size();
130 int nq = GetTotPoints();
131
132 switch (m_projectionType)
133 {
135 {
136 //-------------------------------------------------------
137 // Compute the DG advection including the numerical flux
138 // by using SolverUtils/Advection
139 // Input and output in physical space
141 inarray, outarray, time);
142 //-------------------------------------------------------
143
144 //-------------------------------------------------------
145 // negate the outarray since moving terms to the rhs
146 for (int i = 0; i < nvariables; ++i)
147 {
148 Vmath::Neg(nq, outarray[i], 1);
149 }
150 //-------------------------------------------------------
151
152 //-------------------------------------------------
153 // Add "source terms"
154 // Input and output in physical space
155
156 // Coriolis forcing
157 if (m_coriolis.size() != 0)
158 {
159 AddCoriolis(inarray, outarray);
160 }
161 //-------------------------------------------------
162 }
163 break;
165 {
166 //-------------------------------------------------------
167 // Compute the fluxvector in physical space
169 nvariables);
170
171 for (int i = 0; i < nvariables; ++i)
172 {
173 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
174 for (int j = 0; j < ndim; ++j)
175 {
176 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
177 }
178 }
179
180 LinearSWE::GetFluxVector(inarray, fluxvector);
181 //-------------------------------------------------------
182
183 //-------------------------------------------------------
184 // Take the derivative of the flux terms
185 // and negate the outarray since moving terms to the rhs
186 Array<OneD, NekDouble> tmp0(nq);
187 Array<OneD, NekDouble> tmp1(nq);
188
189 for (int i = 0; i < nvariables; ++i)
190 {
192 fluxvector[i][0], tmp0);
194 fluxvector[i][1], tmp1);
195 Vmath::Vadd(nq, tmp0, 1, tmp1, 1, outarray[i], 1);
196 Vmath::Neg(nq, outarray[i], 1);
197 }
198
199 //-------------------------------------------------
200 // Add "source terms"
201 // Input and output in physical space
202
203 // Coriolis forcing
204 if (m_coriolis.size() != 0)
205 {
206 AddCoriolis(inarray, outarray);
207 }
208 //-------------------------------------------------
209 }
210 break;
211 default:
212 ASSERTL0(false, "Unknown projection scheme for the LinearSWE");
213 break;
214 }
215}
216
218{
220 if (m_session->DefinesSolverInfo("UpwindType"))
221 {
222 std::string UpwindType;
223 UpwindType = m_session->GetSolverInfo("UpwindType");
224 if (UpwindType == "LinearAverage")
225 {
226 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
227 }
228 else if (UpwindType == "LinearHLL")
229 {
230 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
231 }
232 }
233 SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
234 SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
235 SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
236}
237
238// Physfield in primitive Form
240 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
242{
243 int nq = m_fields[0]->GetTotPoints();
244
245 NekDouble g = m_g;
246
247 // Flux vector for the mass equation
248 for (int i = 0; i < m_spacedim; ++i)
249 {
250 Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
251 }
252
253 // Flux vector for the momentum equations
254 for (int i = 0; i < m_spacedim; ++i)
255 {
256 for (int j = 0; j < m_spacedim; ++j)
257 {
258 // must zero fluxes as not initialised to zero in AdvectionWeakDG
259 Vmath::Zero(nq, flux[i + 1][j], 1);
260 }
261
262 // Add (g eta) to appropriate field
263 Vmath::Smul(nq, g, physfield[0], 1, flux[i + 1][i], 1);
264 }
265}
266
267/**
268 * @brief Compute the velocity field \f$ \mathbf{v} \f$ given the momentum
269 * \f$ h\mathbf{v} \f$.
270 *
271 * @param physfield Velocity field.
272 * @param velocity Velocity field.
273 */
275 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
277{
278 const int npts = physfield[0].size();
279
280 for (int i = 0; i < m_spacedim; ++i)
281 {
282 Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
283 }
284}
285
288{
289 int cnt = 0;
290
291 // loop over Boundary Regions
292 for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
293 ++bcRegion)
294 {
295 if (m_fields[0]
296 ->GetBndConditions()[bcRegion]
297 ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
298 {
299 continue;
300 }
301
302 // Copy the forward trace of the field to the backward trace
303 int id2, npts;
304
305 for (int e = 0;
306 e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
307 ++e)
308 {
309 npts = m_fields[0]
310 ->GetBndCondExpansions()[bcRegion]
311 ->GetExp(e)
312 ->GetTotPoints();
313 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
314 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
315 e));
316
317 Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
318 }
319
320 cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
321 }
322}
323
324} // 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)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
static std::string className
Name of class.
Definition: LinearSWE.h:60
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Definition: LinearSWE.cpp:217
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:92
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:87
void v_DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time) override
Definition: LinearSWE.cpp:124
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
Definition: LinearSWE.cpp:54
void GetVelocityVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
Definition: LinearSWE.cpp:274
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:239
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:99
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:100
void CopyBoundaryTrace(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: LinearSWE.cpp:286
Base class for unsteady solvers.
NekDouble m_g
Acceleration of gravity.
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.
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, Array< OneD, Array< OneD, NekDouble > > &outpnt, const NekDouble time, const NekDouble lambda)
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.
SOLVER_UTILS_EXPORT int GetExpSize()
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.
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
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
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