-
Notifications
You must be signed in to change notification settings - Fork 8
/
Copy pathbouss2d.html
449 lines (436 loc) · 27.8 KB
/
bouss2d.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
<!DOCTYPE html>
<html lang="en">
<head>
<meta charset="utf-8" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" /><meta name="generator" content="Docutils 0.19: https://docutils.sourceforge.io/" />
<title>Boussinesq solvers in Two Space Dimensions — Clawpack 5.11.x documentation</title>
<link rel="stylesheet" href="_static/base.css" type="text/css" />
<link rel="stylesheet" href="_static/layout.css" type="text/css" />
<link rel="stylesheet" type="text/css" href="_static/pygments.css" />
<link rel="stylesheet" type="text/css" href="_static/flasky.css" />
<link rel="stylesheet" type="text/css" href="_static/graphviz.css" />
<script data-url_root="./" id="documentation_options" src="_static/documentation_options.js"></script>
<script src="_static/jquery.js"></script>
<script src="_static/underscore.js"></script>
<script src="_static/_sphinx_javascript_frameworks_compat.js"></script>
<script src="_static/doctools.js"></script>
<script src="_static/sphinx_highlight.js"></script>
<link rel="shortcut icon" href="_static/clawicon.ico"/>
<link rel="author" title="About these documents" href="about.html" />
<link rel="index" title="Index" href="genindex.html" />
<link rel="search" title="Search" href="search.html" />
<link rel="next" title="PyClaw" href="pyclaw/index.html" />
<link rel="prev" title="Boussinesq solvers in One Space Dimension" href="bouss1d.html" />
</head><body>
<div id="main-wrapper" class="sphinx">
<div id="header-wrapper">
<section id="header">
<!-- <h1><a href="http://clawpack.org/">Clawpack</a></h1> -->
<h1><a href="http://clawpack.org/">Clawpack-5</a></h1>
<nav>
<ul>
<li>
<a href="contents.html">Docs</a>
</li>
<li>
<a href="installing.html">Install</a>
</li>
<li>
<a class="" href="http://clawpack.org/gallery/index.html">Gallery</a>
</li>
<li>
<a href="about.html">Citation</a>
</li>
<li>
<a class="active" href="http://github.com/clawpack">GitHub</a>
</li>
<li>
<a class="" href="community.html">Community</a>
</li>
<li>
<a class="" href="developers.html">Contribute</a>
</li>
</ul>
</nav>
</section>
<div class="decoration"></div>
</div>
<div class="related" role="navigation" aria-label="related navigation">
<h3>Navigation</h3>
<ul>
<li class="right" style="margin-right: 10px">
<a href="genindex.html" title="General Index"
accesskey="I">index</a></li>
<li class="right" >
<a href="py-modindex.html" title="Python Module Index"
>modules</a> |</li>
<li class="right" >
<a href="pyclaw/index.html" title="PyClaw"
accesskey="N">next</a> |</li>
<li class="right" >
<a href="bouss1d.html" title="Boussinesq solvers in One Space Dimension"
accesskey="P">previous</a> |</li>
<li class="nav-item nav-item-0"><a href="contents.html">Clawpack 5.11.x documentation</a> »</li>
<li class="nav-item nav-item-1"><a href="geoclaw.html" accesskey="U">GeoClaw Description and Detailed Contents</a> »</li>
<li class="nav-item nav-item-this"><a href="">Boussinesq solvers in Two Space Dimensions</a></li>
</ul>
</div>
<div class="document">
<div class="documentwrapper">
<div class="bodywrapper">
<div class="body" role="main">
<section id="boussinesq-solvers-in-two-space-dimensions">
<span id="bouss2d"></span><h1>Boussinesq solvers in Two Space Dimensions<a class="headerlink" href="#boussinesq-solvers-in-two-space-dimensions" title="Permalink to this heading">¶</a></h1>
<p>As of Version 5.10.0, GeoClaw includes the option to solve a
dispersive Boussinesq-type equation known as Serre-Green-Naghdi (SGN)
instead of the usual shallow water equations (SWE).
This equation and related depth-averaged equations have been used
extensively in the literature
to better model wave propagation in situations where the wavelength is not
sufficiently long relative to the fluid depth for the SWE
approximation to be accurate. Applications include the study
of tsunamis generated by landslides, asteroid impacts, or other localized
phenomena. Including dispersive terms has also been found to give more
realistic results for earthquake-generated tsunamis in some situations.
See <a class="reference internal" href="biblio.html#bergerleveque2023" id="id1"><span>[BergerLeVeque2023]</span></a> for references to some of this literature along
with more discussion of the GeoClaw implementation and test problems.</p>
<p>The one-dimensional version of these capabilities are
described in <a class="reference internal" href="bouss1d.html#bouss1d"><span class="std std-ref">Boussinesq solvers in One Space Dimension</span></a>.</p>
<p>The SGN equations are still depth-averaged equations, with the same
conserved quantities (fluid depth <cite>h</cite> and momenta <cite>hu</cite> and <cite>hv</cite> in 2D) as the
shallow water equations (SWE), but the
equations contain higher order derivative terms and so they are no longer
hyperbolic and require implicit methods for efficient solution with a
physically reasonable time step. This adds considerable complexity to the
code since adaptive mesh refinement (AMR) is still supported.
The implementation proceeds by alternating time
steps on the shallow water equations (SWE) with the solution of elliptic
equations where the operator involves second order derivatives in <cite>x</cite> and <cite>y</cite>
of a new set of variables used to modify the momenta each time step.
The right hand side also involves third order derivatives of the topography.</p>
<p>In two space dimensions, solving this
elliptic equation requires setting up and solving a sparse
linear system of equations in each time step, at each refinement level when
AMR is being used. All grid cells from all patches at
the same refinement level
are included in the linear system. Boundary conditions at the edge of
patches must be interpolated from coarser level solutions, in much the same
way that the boundary conditions for <cite>h</cite>, <cite>hu</cite>, and <cite>hv</cite> are interpolated
when solving the SWE with AMR. Because the solution of the elliptic system
yields correction terms to the momenta (denoted here by <cite>huc</cite> and <cite>hvc</cite>),
when solving the Boussinesq equations the array <cite>q</cite> of state variables
is expanded to include these correction terms as well, and so the number of
equations and variables is 5 instead of 3. This change must be made in
<cite>setrun.py</cite>, along with other changes discussed below, in order to use
the Boussinesq solvers.</p>
<p>Currently the only linear system solver supported for solving the large
sparse elliptic systems is <a class="reference external" href="https://petsc.org/release/">PETSc</a>,
which can use MPI to solve these these systems. Using the Boussinesq solvers
requires these prerequesites, as discussed further in <a class="reference internal" href="#bouss2d-prereqs"><span class="std std-ref">Prerequisites for the 2d Boussinesq code</span></a>.</p>
<p>See <a class="reference internal" href="biblio.html#bergerleveque2023" id="id2"><span>[BergerLeVeque2023]</span></a> for more discussion of the equations solved and the
AMR algorithms developed for these equations.</p>
<section id="boussinesq-type-dispersive-equations">
<span id="bouss2d-eqns"></span><h2>Boussinesq-type dispersive equations<a class="headerlink" href="#boussinesq-type-dispersive-equations" title="Permalink to this heading">¶</a></h2>
<p>The equations we solve are not the original depth-averaged dispersive
equations derived by Boussinesq, but for simplicity
in this documentation and the code, we often refer to the
equations simply as “Boussinesq equations”, following common practice.
Many variants of these equations have been derived and fine-tuned to
better match the dispersion relation of the linearized
<a class="reference external" href="https://en.wikipedia.org/wiki/Airy_wave_theory">Airy wave theory</a>
and to incorporate variable bottom topography.</p>
<p>Two variants are currently implement in GeoClaw, described below.
In practice we recommend using only the SGN equations, which we have found
to be more stable.</p>
</section>
<section id="the-sgn-equations">
<span id="bouss2d-sgn"></span><h2>The SGN equations<a class="headerlink" href="#the-sgn-equations" title="Permalink to this heading">¶</a></h2>
<p>The Serre-Green-Naghdi (SGN) equations implemented in GeoClaw
are generalized to include a parameter <cite>alpha</cite>
suggested by Bonneton et al. Both the 1D and 2D versions of these equations
and their GeoClaw implementation are discussed in <a class="reference internal" href="biblio.html#bergerleveque2023" id="id3"><span>[BergerLeVeque2023]</span></a>.
The value <cite>alpha = 1.153</cite> is
recommended since it gives a better approximation to the linear dispersion
relation of the Airy solution to the full 3d problem.
This value is
hardwired into <cite>$CLAW/geoclaw/src/2d/bouss/bouss_module.f90</cite>. To change
this value, you must modify this module. (See <a class="reference internal" href="makefiles_library.html#makefiles-library"><span class="std std-ref">Library routines in Makefiles</span></a>
for tips on modifying a library routine.)
Setting <cite>alpha = 1</cite> gives the original SGN equations.</p>
</section>
<section id="the-madsen-sorensen-ms-equations">
<span id="bouss2d-ms"></span><h2>The Madsen-Sorensen (MS) equations<a class="headerlink" href="#the-madsen-sorensen-ms-equations" title="Permalink to this heading">¶</a></h2>
<p>Primarily for historical reasons, GeoClaw also includes an implementation of
another Boussinesq-type dispersive system, the Madsen-Sorensen (MS) equations.
These equations also give a good approximation to the linear dispersion
relation of the Airy solution when the parameter <cite>Bparam = 1/15</cite> is used,
which is hardwired into <cite>$CLAW/geoclaw/src/2d/bouss/bouss_module.f90</cite>.
These equations were used in an earlier GeoClaw implementation
by Jihwan Kim, known as BoussClaw <a class="reference internal" href="biblio.html#kimetal2017" id="id4"><span>[KimEtAl2017]</span></a>.
This implementation was successfully used in a number of studies
(see <a class="reference internal" href="biblio.html#bergerleveque2023" id="id5"><span>[BergerLeVeque2023]</span></a> for some citations).
However, extensive tests with these equations have revealed stability issues,
particularly with the use of AMR, which have also been reported by others.
Implementations of MS in both 1D and 2D are included in GeoClaw,
but are generally not
recommended except for those interested in comparing different
formulations for small numbers of time steps,
and perhaps further investigating the stability issues.</p>
</section>
<section id="using-the-2d-boussinesq-code">
<span id="bouss2d-usage"></span><h2>Using the 2d Boussinesq code<a class="headerlink" href="#using-the-2d-boussinesq-code" title="Permalink to this heading">¶</a></h2>
<p>Provided the <a class="reference internal" href="#bouss2d-prereqs"><span class="std std-ref">Prerequisites for the 2d Boussinesq code</span></a> have been installed, switching from the
SWE to a Boussinesq solver in GeoClaw requires only minor changes to
<cite>setrun.py</cite> and the <cite>Makefile</cite>.</p>
<p>See the files in <cite>$CLAW/geoclaw/examples/bouss/radial_flat</cite> for an example.</p>
<section id="setrun-py">
<span id="bouss2d-setrun"></span><h3>setrun.py<a class="headerlink" href="#setrun-py" title="Permalink to this heading">¶</a></h3>
<p>As mentioned above, it is necessary to set:</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span><span class="n">clawdata</span><span class="o">.</span><span class="n">num_eqn</span> <span class="o">=</span> <span class="mi">5</span>
</pre></div>
</div>
<p>instead of the usual value 3 used for SWE, since correction terms for the
momenta are also stored in the <cite>q</cite> arrays to facilitate interpolation to
ghost cells of finer level patches each time step.</p>
<p>It is also necessary to set some parameters that are specific to the
Boussinesq solvers. Somewhere in the <cite>setrun</cite> function you must include</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span><span class="kn">from</span> <span class="nn">clawpack.geoclaw.data</span> <span class="kn">import</span> <span class="n">BoussData</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">add_data</span><span class="p">(</span><span class="n">BoussData</span><span class="p">(),</span><span class="s1">'bouss_data'</span><span class="p">)</span>
</pre></div>
</div>
<p>and then the following parameters can be adjusted (the values shown here
are the default values that will be used if you do not specify a value
directly):</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span><span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_equations</span> <span class="o">=</span> <span class="mi">2</span> <span class="c1"># 0=SWE, 1=MS, 2=SGN</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_min_level</span> <span class="o">=</span> <span class="mi">1</span> <span class="c1"># coarsest level to apply bouss</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_max_level</span> <span class="o">=</span> <span class="mi">10</span> <span class="c1"># finest level to apply bouss</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_min_depth</span> <span class="o">=</span> <span class="mf">10.</span> <span class="c1"># depth (meters) to switch to SWE</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_solver</span> <span class="o">=</span> <span class="mi">3</span> <span class="c1"># 3=PETSc</span>
<span class="n">rundata</span><span class="o">.</span><span class="n">bouss_data</span><span class="o">.</span><span class="n">bouss_tstart</span> <span class="o">=</span> <span class="mf">0.</span> <span class="c1"># time to switch from SWE</span>
</pre></div>
</div>
<p>These parameters are described below:</p>
<ul>
<li><p><cite>bouss_equations</cite>: The system of equations being solved. Setting this to 2
gives the recommended SGN equations, while 1 gives Madsen-Sorensen.</p>
<p>Setting <cite>bouss_equations = 0</cite> causes the code to revert to the shallow
water equations, useful for comparing dispersive and nondispersive results.
(But if <cite>bouss_data</cite> is being set, it still requires <cite>clawdata.num_eqn = 5</cite>
and the two new components in q are always 0 in this case, so this is
slightly less efficient than using the standard GeoClaw.)</p>
</li>
<li><p><cite>bouss_min_level</cite>: The minimum AMR level on which Boussinesq correction
terms should be applied. In some cases it may be desirable to use the SWE
on the coarsest grids in the ocean while Boussinesq corrections are only
applied on fine levels near shore, for example.</p></li>
<li><p><cite>bouss_max_level</cite>: The finest AMR level on which Boussinesq correction
terms should be applied. In some cases it may be desirable to use the SWE
only on coarser grids if the finest level grid only exists in very shallow
regions or onshore, where the the equations switch to SWE for inundation
modeling. Since much of the computational work is often on the finest level,
avoiding the Boussinesq terms altogether on these levels may be advantageous
in some situations.</p></li>
<li><p><cite>bouss_min_depth</cite>: The criterion used for switching from Boussinesq to SWE
in shallow water and onshore. If the original water depth <cite>h</cite> at time <cite>t0</cite>
is less than <cite>bouss_min_depth</cite> in a cell or any of its nearest
neighbors in a 3-by-3 neighborhood,
then this cell is omitted from set of unknowns in the elliptic equation
solve and no dispersive correction terms are calculated for this cell.
This is discussed further below in <a class="reference internal" href="#bouss2d-switch"><span class="std std-ref">Wave breaking and switching to SWE</span></a>.</p></li>
<li><p><cite>bouss_solver</cite>: What linear system solver to use. Currently only the value
3 for <a class="reference external" href="https://petsc.org/release/">PETSc</a> is recognized.</p></li>
<li><p><cite>bouss_tstart</cite>: The time <cite>t</cite> at which to start applying Boussinesq terms.
Normally you will want this to be less than or equal to <cite>t0</cite>, the starting
time of the calculation (which is not always 0). However,
there are some cases in which the initial data results in extreme
motion in the first few time steps and it is necessary to get things going
with the SWE. For most applications this is not necessary and you need
only change this parameter if you are solving a problem for which <cite>t0 < 0</cite>.</p></li>
</ul>
</section>
<section id="makefile">
<span id="bouss2d-makefile"></span><h3>Makefile<a class="headerlink" href="#makefile" title="Permalink to this heading">¶</a></h3>
<p>You can copy the <cite>Makefile</cite> from
<cite>$CLAW/geoclaw/examples/bouss/radial_flat/Makefile</cite> and make any adjustments
needed.</p>
<p>This <cite>Makefile</cite> reads in the standard Boussinesq solver file
<cite>$CLAW/geoclaw/src/2d/bouss/Makefile.bouss</cite>, which lists the Fortran modules
and source code files that are used by default from the library
<cite>$CLAW/geoclaw/src/2d/bouss</cite>, or from <cite>$CLAW/amrclaw/src/2d</cite> or
<cite>$CLAW/geoclaw/src/2d/shallow</cite> in the case of files that did not need to
be modified for the Boussinesq code.</p>
<p>Two <cite>Makefile</cite> variables <cite>PETSC_DIR</cite> and <cite>PETSC_ARCH</cite> must be set (perhaps as
environment variables in the shell from which <cite>make</cite> is invoked). These are
described further below in <a class="reference internal" href="#bouss2d-prereqs"><span class="std std-ref">Prerequisites for the 2d Boussinesq code</span></a>.</p>
<p>The <cite>FFLAGS</cite> specified in the <cite>Makefile</cite> should include <cite>-DHAVE_PETSC</cite>
to indicate that <cite>PETSc</cite> is being used, necessary when compiling the
source code for Boussinesq solvers.</p>
<p>The <cite>Makefile</cite> should also include a line of the form:</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span>PETSC_OPTIONS=-options_file $(CLAW)/geoclaw/examples/bouss/petscMPIoptions
</pre></div>
</div>
<p>with a pointer to the file that sets various <cite>PETSc</cite> options. The file
<cite>$CLAW/geoclaw/examples/bouss/petscMPIoptions</cite> gives the options used in
the examples, which may be adequate for other problems too.
This file includes some comments briefly explaining the options set.
We use a GMRES Krylov space method as the main solver
and algebraic multigrid as the preconditioner.
For more about the options for these methods, see:</p>
<blockquote>
<div><ul class="simple">
<li><p><a class="reference external" href="https://petsc.org/release/manualpages/KSP/KSPSetFromOptions">https://petsc.org/release/manualpages/KSP/KSPSetFromOptions</a></p></li>
<li><p><a class="reference external" href="https://petsc.org/release/manualpages/PC/PCSetFromOptions/">https://petsc.org/release/manualpages/PC/PCSetFromOptions/</a></p></li>
</ul>
</div></blockquote>
<p>In addition to a line of the form</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span><span class="n">EXE</span> <span class="o">=</span> <span class="n">xgeoclaw</span>
</pre></div>
</div>
<p>that specifies the name and location of the executable to be generated, the
<cite>Makefile</cite> should also contain a line of the form:</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span><span class="n">RUNEXE</span><span class="o">=</span><span class="s2">"$</span><span class="si">{PETSC_DIR}</span><span class="s2">/$</span><span class="si">{PETSC_ARCH}</span><span class="s2">/bin/mpiexec -n 6"</span>
</pre></div>
</div>
<p>This is the command that should be used in order to run the executable.
In other words, if you set <cite>PETSC_DIR</cite> and <cite>PETSC_ARCH</cite> as environment
variables, and the executable is named <cite>xgeoclaw</cite> as usual, then the command</p>
<div class="highlight-default notranslate"><div class="highlight"><pre><span></span>$PETSC_DIR/$PETSC_ARCH/bin/mpiexec -n 6 xgeoclaw
</pre></div>
</div>
<p>given in the shell should run the executable (invoking MPI with 6 processes in
this example). If this does not work then one of the environment variables
may be set incorrectly to find the <cite>mpiexec</cite> command.</p>
</section>
</section>
<section id="prerequisites-for-the-2d-boussinesq-code">
<span id="bouss2d-prereqs"></span><h2>Prerequisites for the 2d Boussinesq code<a class="headerlink" href="#prerequisites-for-the-2d-boussinesq-code" title="Permalink to this heading">¶</a></h2>
<p>Currently the only linear solver supported is <cite>PETSc</cite>, so this must be
installed, see <a class="reference external" href="https://petsc.org/release/install/">https://petsc.org/release/install/</a> for instructions
and also note the <a class="reference external" href="https://petsc.org/release/install/install_tutorial/#prerequisites">PETSc prerequisites</a>.
Note that MPI, LAPACK, and the BLAS are required and will be installed as
part of installing PETSc. If you already have some of the prerequisites
installed, be sure to read <a class="reference external" href="https://petsc.org/release/install/install/#configuring-petsc">Configuring PETSc</a>
before installing.</p>
<p>The environment variables <cite>$PETSC_DIR</cite> and <cite>$PETSC_ARCH</cite> must be set
appropriately based on your PETSc installation, either as environment
variables or directly in the <cite>Makefile</cite>.
See the PETSc documentation page
<a class="reference external" href="https://petsc.org/release/install/multibuild/#environmental-variables-petsc-dir-and-petsc-arch">Environmental Variables $PETSC_DIR And $PETSC_ARCH</a>.</p>
</section>
<section id="wave-breaking-and-switching-to-swe">
<span id="bouss2d-switch"></span><h2>Wave breaking and switching to SWE<a class="headerlink" href="#wave-breaking-and-switching-to-swe" title="Permalink to this heading">¶</a></h2>
<p>The <cite>bouss_min_depth</cite> parameter is needed because in very shallow water, and for
modeling onshore inundation, the Boussinesq equations are not suitable.
So some criterion is needed to drop these correction terms and revert to
solving SWE near shore. Many different approaches have been used in the
literature. So far we have only implemented the simplest commonly used approach,
which is to revert to SWE in any grid cell where the initial water depth (at
the initial time) is less than <cite>bouss_min_depth</cite>.</p>
</section>
<section id="examples">
<h2>Examples<a class="headerlink" href="#examples" title="Permalink to this heading">¶</a></h2>
<p>In addition to one example application included in GeoClaw, found in the
directory <cite>$CLAW/geoclaw/examples/bouss/radial_flat</cite>, several other examples
of usage can be found in the code repository
<a class="reference external" href="https://github.com/rjleveque/ImplicitAMR-paper">https://github.com/rjleveque/ImplicitAMR-paper</a>, which was
developed to accompany the paper <a class="reference internal" href="biblio.html#bergerleveque2023" id="id6"><span>[BergerLeVeque2023]</span></a>.</p>
</section>
</section>
<div class="clearer"></div>
</div>
</div>
</div>
<div class="sphinxsidebar" role="navigation" aria-label="main navigation">
<div class="sphinxsidebarwrapper">
<p><a href="http://clawpack.org/">
<img class="logo" src= "_static/clawlogo.jpg" alt="Logo"/>
</a>
<h2>Version 5.11.x</h2>
</p>
<div>
<h3><a href="contents.html">Table of Contents</a></h3>
<ul>
<li><a class="reference internal" href="#">Boussinesq solvers in Two Space Dimensions</a><ul>
<li><a class="reference internal" href="#boussinesq-type-dispersive-equations">Boussinesq-type dispersive equations</a></li>
<li><a class="reference internal" href="#the-sgn-equations">The SGN equations</a></li>
<li><a class="reference internal" href="#the-madsen-sorensen-ms-equations">The Madsen-Sorensen (MS) equations</a></li>
<li><a class="reference internal" href="#using-the-2d-boussinesq-code">Using the 2d Boussinesq code</a><ul>
<li><a class="reference internal" href="#setrun-py">setrun.py</a></li>
<li><a class="reference internal" href="#makefile">Makefile</a></li>
</ul>
</li>
<li><a class="reference internal" href="#prerequisites-for-the-2d-boussinesq-code">Prerequisites for the 2d Boussinesq code</a></li>
<li><a class="reference internal" href="#wave-breaking-and-switching-to-swe">Wave breaking and switching to SWE</a></li>
<li><a class="reference internal" href="#examples">Examples</a></li>
</ul>
</li>
</ul>
</div><h3>Related Topics</h3>
<ul>
<li><a href="contents.html">Documentation overview</a><ul>
<li><a href="geoclaw.html">GeoClaw Description and Detailed Contents</a><ul>
<li>Previous: <a href="bouss1d.html" title="previous chapter">Boussinesq solvers in One Space Dimension</a></li>
<li>Next: <a href="pyclaw/index.html" title="next chapter">PyClaw</a></li>
</ul></li>
</ul></li>
</ul>
<div class="widget navlinks">
<h3>This Page</h3>
<ul class="this-page-menu">
<li><a href="_sources/bouss2d.rst.txt"
rel="nofollow"
target="_blank">Source .rst</a></li>
<li><a href="https://github.com/clawpack/doc/blob/dev/doc/bouss2d.rst"
rel="nofollow"
target="_blank">Source on GitHub</a></li>
<li><a href="https://github.com/clawpack/doc/commits/dev/doc/bouss2d.rst"
rel="nofollow"
target="_blank">History</a></li>
<li><a href="https://github.com/clawpack/doc/edit/dev/doc/bouss2d.rst"
rel="nofollow"
target="_blank">Suggest Edits</a></li>
<li><a href="https://github.com/clawpack/doc/issues/new/choose"
rel="nofollow"
target="_blank">Raise an Issue</a></li>
</ul>
</div>
<div id="searchbox" style="display: none" role="search">
<h3 id="searchlabel">Quick search</h3>
<div class="searchformwrapper">
<form class="search" action="search.html" method="get">
<input type="text" name="q" aria-labelledby="searchlabel" autocomplete="off" autocorrect="off" autocapitalize="off" spellcheck="false"/>
<input type="submit" value="Go" />
</form>
</div>
</div>
<script>document.getElementById('searchbox').style.display = "block"</script>
<h4>Latest Version</h4>
<ul>
<li><a href="./dev/bouss2d.html">dev</a></li>
<li><a href="bouss2d.html">v5.11.x</a></li>
</ul>
<h4>Older Versions</h4>
<ul>
<li><a href="./v5.10.x/bouss2d.html">v5.10.x</a></li>
<li><a href="./v5.7.x/index.html">v5.7.x</a></li>
<li><a href="./v5.8.x/index.html">v5.8.x</a></li>
<li><a href="./v5.9.x/bouss2d.html">v5.9.x</a></li>
</ul>
</div>
</div>
<div class="clearer"></div>
</div>
<div class="footer">
© Copyright CC-BY 2024, The Clawpack Development Team.
Created using <a href="http://sphinx.pocoo.org/">Sphinx</a>.
</div>
<script>
(function(i,s,o,g,r,a,m){i['GoogleAnalyticsObject']=r;i[r]=i[r]||function(){
(i[r].q=i[r].q||[]).push(arguments)},i[r].l=1*new Date();a=s.createElement(o),
m=s.getElementsByTagName(o)[0];a.async=1;a.src=g;m.parentNode.insertBefore(a,m)
})(window,document,'script','//www.google-analytics.com/analytics.js','ga');
ga('create', 'UA-44811544-1', 'auto');
ga('send', 'pageview');
</script>
</body>
</html>