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 
39 using namespace std;
40 
41 namespace 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 
62 string ShallowWaterSystem::className =
64  "ShallowWaterSystem", ShallowWaterSystem::create,
65  "Auxiliary functions for the shallow water system.");
66 
67 ShallowWaterSystem::ShallowWaterSystem(
70  : UnsteadySystem(pSession, pGraph)
71 {
72 }
73 
74 void 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 generic input parameters
94  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
95 
96  // Load acceleration of gravity
97  m_session->LoadParameter("Gravity", m_g, 9.81);
98 
99  // input/output in primitive variables
100  m_primitive = true;
101 
103 
104  m_constantDepth = true;
105  NekDouble depth = m_depth[0];
106  for (int i = 0; i < GetTotPoints(); ++i)
107  {
108  if (m_depth[i] != depth)
109  {
110  m_constantDepth = false;
111  break;
112  }
113  }
114 
115  // Compute the bottom slopes
116  int nq = GetTotPoints();
117  if (m_constantDepth != true)
118  {
120  for (int i = 0; i < m_spacedim; ++i)
121  {
124  m_bottomSlope[i]);
125  Vmath::Neg(nq, m_bottomSlope[i], 1);
126  }
127  }
128 
130 }
131 
132 /**
133  *
134  */
136 {
137 }
138 
140 {
142  if (m_constantDepth == true)
143  {
144  SolverUtils::AddSummaryItem(s, "Depth", "constant");
145  }
146  else
147  {
148  SolverUtils::AddSummaryItem(s, "Depth", "variable");
149  }
150 }
151 
153 {
154  ASSERTL0(false, "This function is not implemented for this equation.");
155 }
156 
158 {
159  ASSERTL0(false, "This function is not implemented for this equation.");
160 }
161 
163 {
164  GetFunction("WaterDepth")->Evaluate("d", m_depth);
165 }
166 
168 {
169  GetFunction("Coriolis")->Evaluate("f", m_coriolis);
170 }
171 
174 {
175 
176  int cnt = 0;
177  // loop over Boundary Regions
178  for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
179  ++bcRegion)
180  {
181  if (m_fields[0]
182  ->GetBndConditions()[bcRegion]
183  ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
184  {
185  continue;
186  }
187 
188  // Copy the forward trace of the field to the backward trace
189  int e, id2, npts;
190 
191  for (e = 0;
192  e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
193  ++e)
194  {
195  npts = m_fields[0]
196  ->GetBndCondExpansions()[bcRegion]
197  ->GetExp(e)
198  ->GetTotPoints();
199  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
200  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
201  e));
202 
203  Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
204  }
205 
206  cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
207  }
208 }
209 
210 } // 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
NekDouble m_g
Acceleration of gravity.
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
bool m_constantDepth
Indicates if constant depth case.
bool m_primitive
Indicates if variables are primitive or conservative.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
virtual ~ShallowWaterSystem()
Destructor.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
Array< OneD, NekDouble > m_depth
Still water depth.
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
int m_infosteps
Number of time steps between outputting status information.
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 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.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
virtual 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:91
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
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: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
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255