File: handson.md.txt

package info (click to toggle)
petsc 3.23.1%2Bdfsg1-1exp1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 515,576 kB
  • sloc: ansic: 751,607; cpp: 51,542; python: 38,598; f90: 17,352; javascript: 3,493; makefile: 3,157; sh: 1,502; xml: 619; objc: 445; java: 13; csh: 1
file content (362 lines) | stat: -rw-r--r-- 9,052 bytes parent folder | download
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
(handson)=

# Tutorials, by Mathematical Problem

TODO: Add link to Python example here

(handson-example-1)=

## Linear elliptic PDE on a 2D grid

WHAT THIS EXAMPLE DEMONSTRATES:

- Using command line options
- Using Linear Solvers
- Handling a simple structured grid

FURTHER DETAILS:

- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex50.c.html#line1">Mathematical description of the problem</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex50.c.html#line21">the source code</a>

DO THE FOLLOWING:

- Compile `src/ksp/ksp/tutorials/ex50.c`

  ```console
  $ cd petsc/src/ksp/ksp/tutorials
  $ make ex50
  ```

- Run a 1 processor example with a 3x3 mesh and view the matrix
  assembled

  ```console
  $ mpiexec -n 1 ./ex50  -da_grid_x 4 -da_grid_y 4 -mat_view
  ```

  Expected output:

  ```{literalinclude} /../src/ksp/ksp/tutorials/output/ex50_tut_1.out
  :language: none
  ```

- Run with a 120x120 mesh on 4 processors using superlu_dist and
  view the solver options used

  ```console
  $ mpiexec -n 4 ./ex50  -da_grid_x 120 -da_grid_y 120 -pc_type lu -pc_factor_mat_solver_type superlu_dist -ksp_monitor -ksp_view
  ```

  Expected output:

  ```{literalinclude} /../src/ksp/ksp/tutorials/output/ex50_tut_2.out
  :language: none
  ```

- Run with a 1025x1025 grid using multigrid solver on 4
  processors with 9 multigrid levels

  ```console
  $ mpiexec -n 4 ./ex50 -da_grid_x 1025 -da_grid_y 1025 -pc_type mg -pc_mg_levels 9 -ksp_monitor
  ```

  Expected output:

  ```{literalinclude} /../src/ksp/ksp/tutorials/output/ex50_tut_3.out
  :language: none
  ```

(handson-example-2)=

## Nonlinear ODE arising from a time-dependent one-dimensional PDE

WHAT THIS EXAMPLE DEMONSTRATES:

- Using command line options
- Handling a simple structured grid
- Using the ODE integrator
- Using call-back functions

FURTHER DETAILS:

- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex2.c.html#line13">Mathematical description of the problem</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex2.c.html#line36">the source code</a>

DO THE FOLLOWING:

- Compile `src/ts/tutorials/ex2.c`

  ```console
  $ cd petsc/src/ts/tutorials
  $ make ex2
  ```

- Run a 1 processor example on the default grid with all the
  default solver options

  ```console
  $ mpiexec -n 1 ./ex2 -ts_max_steps 10 -ts_monitor
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex2_tut_1.out
  :language: none
  ```

- Run with the same options on 4 processors plus monitor
  convergence of the nonlinear and linear solvers

  ```console
  $ mpiexec -n 4 ./ex2 -ts_max_steps 10 -ts_monitor -snes_monitor -ksp_monitor
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex2_tut_2.out
  :language: none
  ```

- Run with the same options on 4 processors with 128 grid points

  ```console
  $ mpiexec -n 16 ./ex2 -ts_max_steps 10 -ts_monitor -M 128
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex2_tut_3.out
  :language: none
  ```

(handson-example-3)=

## Nonlinear PDE on a structured grid

WHAT THIS EXAMPLE DEMONSTRATES:

- Handling a 2d structured grid
- Using the nonlinear solvers
- Changing the default linear solver

FURTHER DETAILS:

- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex19.c.html#line19">Mathematical description of the problem</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex19.c.html#line94">main program source code</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex19.c.html#line246">physics source code</a>

DO THE FOLLOWING:

- Compile `src/snes/tutorials/ex19.c`

  ```console
  $ cd petsc/src/snes/tutorials/
  $ make ex19
  ```

- Run a 4 processor example with 5 levels of grid refinement,
  monitor the convergence of the nonlinear and linear solver and
  examine the exact solver used

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -snes_monitor -ksp_monitor -snes_view
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_1.out
  :language: none
  ```

- Run with the same options but use geometric multigrid as the
  linear solver

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type mg
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_2.out
  :language: none
  ```

  Note this requires many fewer iterations than the default
  solver

- Run with the same options but use algebraic multigrid (hypre's
  BoomerAMG) as the linear solver

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type hypre
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_3.out
  :language: none
  ```

  Note this requires many fewer iterations than the default
  solver but requires more linear solver iterations than
  geometric multigrid.

- Run with the same options but use the ML preconditioner from
  Trilinos

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type ml
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_8.out
  :language: none
  ```

- Run on 1 processor with the default linear solver and profile
  the run

  ```console
  $ mpiexec -n 1 ./ex19 -da_refine 5 -log_view
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_4.out
  :language: none
  ```

  Search for the line beginning with SNESSolve, the fourth column
  gives the time for the nonlinear solve.

- Run on 1 processor with the geometric multigrid linear solver
  and profile the run

  ```console
  $ mpiexec -n 1 ./ex19 -da_refine 5 -log_view -pc_type mg
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_5.out
  :language: none
  ```

  Compare the runtime for SNESSolve to the case with the default
  solver

- Run on 4 processors with the default linear solver and profile
  the run

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -log_view
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_6.out
  :language: none
  ```

  Compare the runtime for `SNESSolve` to the 1 processor case with
  the default solver. What is the speedup?

- Run on 4 processors with the geometric multigrid linear solver
  and profile the run

  ```console
  $ mpiexec -n 4 ./ex19 -da_refine 5 -log_view -pc_type mg
  ```

  Expected output:

  ```{literalinclude} /../src/snes/tutorials/output/ex19_tut_7.out
  :language: none
  ```

  Compare the runtime for SNESSolve to the 1 processor case with
  multigrid. What is the speedup? Why is the speedup for
  multigrid lower than the speedup for the default solver?

(handson-example-4)=

## Nonlinear time dependent PDE on unstructured grid

WHAT THIS EXAMPLE DEMONSTRATES:

- Changing the default ODE integrator
- Handling unstructured grids
- Registering your own interchangeable physics and algorithm
  modules

FURTHER DETAILS:

- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex11.c.html">Mathematical description of the problem</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex11.c.html#line1403">main program source code</a>
- <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex11.c.html#line186">source code of physics modules</a>

DO THE FOLLOWING:

- Compile `src/ts/tutorials/ex11.c`

  ```console
  $ cd petsc/src/ts/tutorials
  $ make ex11
  ```

- Run simple advection through a tiny hybrid mesh

  ```console
  $ mpiexec -n 1 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex11_tut_1.out
  :language: none
  ```

- Run simple advection through a small mesh with a Rosenbrock-W
  solver

  ```console
  $ mpiexec -n 1 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex11_tut_2.out
  :language: none
  ```

- Run simple advection through a larger quadrilateral mesh of an
  annulus with least squares reconstruction and no limiting,
  monitoring the error

  ```console
  $ mpiexec -n 4 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex11_tut_3.out
  :language: none
  ```

  Compare turning to the error after turning off reconstruction.

- Run shallow water on the larger mesh with least squares
  reconstruction and minmod limiting, monitoring water Height
  (integral is conserved) and Energy (not conserved)

  ```console
  $ mpiexec -n 4 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod
  ```

  Expected output:

  ```{literalinclude} /../src/ts/tutorials/output/ex11_tut_4.out
  :language: none
  ```