File: HandsOnExercise.html

package info (click to toggle)
petsc 3.7.5%2Bdfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,864 kB
  • ctags: 618,438
  • sloc: ansic: 515,133; python: 29,793; makefile: 20,458; fortran: 18,998; cpp: 6,515; f90: 3,914; sh: 1,012; xml: 621; objc: 445; csh: 240; java: 13
file content (270 lines) | stat: -rw-r--r-- 12,082 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
<html><head> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/tutorials/HandsOnExercise.html" />
    <title>PETSc Hands On</title>
</head><body>
   <div id="version" align=right><b>petsc-3.7.5 2017-01-01</b></div>
   <div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.7.5 v3.7.5 tutorials/HandsOnExercise.html "><small>Report Typos and Errors</small></a></div>

<center>
<h2> PETSc Hands On</h2>
</center>
<p>
PETSc comes with a large number of example codes to illustrate usage. The three main sets are in:

<ul>
<li> ODEs - <a href="../src/ts/examples/tutorials/index.html">src/ts/examples/tutorials</a></li>
<li> Nonlinear problems - <a href="../src/snes/examples/tutorials/index.html">src/snes/examples/tutorials</a></li>
<li> Linear problems - <a href="../src/ksp/ksp/examples/tutorials/index.html">src/ksp/ksp/examples/tutorials</a></li>
</ul>

<!--------------------------------------------------------------------------->
<h2>Example 1: Linear Poisson equation on a 2D grid</h2>

<p>WHAT THIS EXAMPLE DEMONSTRATES:</p>
<ul>
  <li>Using command line options</li>
  <li>Using Linear Solvers</li>
  <li>Handling a simple structured grid</li>
</ul>

<p>FURTHER DETAILS:</p>
<ul>
  <li><a href="../src/ksp/ksp/examples/tutorials/ex50.c.html#line1">Mathematical description of the problem</a></li>
  <li><a href="../src/ksp/ksp/examples/tutorials/ex50.c.html#line21">the source code</a></li>
</ul>

<p>DO THE FOLLOWING:</p>
<ul>
  <li> Compile src/ksp/ksp/examples/tutorials/ex50.c
    <pre>
      cd petsc/src/ksp/ksp/examples/tutorials
      make ex50
  </pre>
  </li>
  <li> Run a 1 processor example with a 3x3 mesh and view the matrix assembled
    <pre>
     mpiexec -n 1 ./ex50  -da_grid_x 4 -da_grid_y 4 -mat_view     <a href="../src/ksp/ksp/examples/tutorials/output/ex50_tut_1.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Run with a 120x120 mesh on 16 processors using superlu_dist and view the solver options used
    <pre>
     mpiexec -n 16 ./ex50  -da_grid_x 120 -da_grid_y 120 -pc_type lu -pc_factor_mat_solver_package superlu_dist -ksp_monitor -ksp_view     <a href="../src/ksp/ksp/examples/tutorials/output/ex50_tut_2.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Run with a 2049x2049 grid using multigrid solver on 16 processors with 10 multigrid levels
    <pre>
     mpiexec -n 16 ./ex50 -da_grid_x 2049 -da_grid_y 2049 -pc_type mg -pc_mg_levels 10 -ksp_monitor     <a href="../src/ksp/ksp/examples/tutorials/output/ex50_tut_3.out.html">[Expected output]</a>
    </pre>
  </li>
</ul>


<!--------------------------------------------------------------------------->

<h2>Example 2: Nonlinear ODE arising from a time-dependent one dimensional PDE</h2>

<p>WHAT THIS EXAMPLE DEMONSTRATES:</p>
<ul>
  <li>Using command line options</li>
  <li>Handling a simple structured grid</li>
  <li>Using the ODE integrator</li>
  <li>Using call-back functions</li>
</ul>

<p>FURTHER DETAILS:</p>
<ul>
  <li><a href="../src/ts/examples/tutorials/ex2.c.html#line13">Mathematical description of the problem</a></li>
  <li><a href="../src/ts/examples/tutorials/ex2.c.html#line36">the source code</a></li>
</ul>

<p>DO THE FOLLOWING:</p>
<ul>
  <li> Compile src/ts/examples/tutorials/ex2.c
    <pre>
      cd petsc/src/ts/examples/tutorials
      make ex2
  </pre>
  </li>
  <li> Run a 1 processor example on the default grid with all the default solver options
    <pre>
     mpiexec -n 1 ./ex2 -ts_max_steps 10 -ts_monitor     <a href="../src/ts/examples/tutorials/output/ex2_tut_1.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Run with the same options on 4 processors plus monitor convergence of the nonlinear and linear solvers
    <pre>
     mpiexec -n 4 ./ex2 -ts_max_steps 10 -ts_monitor -snes_monitor -ksp_monitor     <a href="../src/ts/examples/tutorials/output/ex2_tut_2.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Run with the same options on 16 processors with 128 grid points
    <pre>
     mpiexec -n 16 ./ex2 -ts_max_steps 10 -ts_monitor -M 128     <a href="../src/ts/examples/tutorials/output/ex2_tut_3.out.html">[Expected output]</a>
    </pre>
  </li>
</ul>


<!--------------------------------------------------------------------------->

<h2>Example 3: Nonlinear PDE on a structured grid</h2>

<p>WHAT THIS EXAMPLE DEMONSTRATES:</p>
<ul>
  <li>Handling a 2d structured grid</li>
  <li>Using the nonlinear solvers</li>
  <li>Changing the default linear solver</li>
</ul>

<p>FURTHER DETAILS:</p>
<ul>
  <li><a href="../src/snes/examples/tutorials/ex19.c.html#line19">Mathematical description of the problem</a></li>
  <li><a href="../src/snes/examples/tutorials/ex19.c.html#line94">main program source code</a></li>
  <li><a href="../src/snes/examples/tutorials/ex19.c.html#line246">physics source code</a></li>
</ul>

<p>DO THE FOLLOWING:</p>
<ul>
  <li> Compile src/snes/examples/tutorials/ex19.c
    <pre>
      cd petsc/src/snes/examples/tutorials/
      make ex19
  </pre>
  </li>
  <li> Run a 16 processor example with six levels of grid refinement, monitor the convergence of the nonlinear and linear solver
      and examine the exact solver used
    <pre>
     mpiexec -n 16 ./ex19 -da_refine 6 -snes_monitor -ksp_monitor -snes_view     <a href="../src/snes/examples/tutorials/output/ex19_tut_1.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Run with the same options but use geometric multigrid as the linear solver
    <pre>
     mpiexec -n 16 ./ex19 -da_refine 6 -snes_monitor -ksp_monitor -snes_view -pc_type mg     <a href="../src/snes/examples/tutorials/output/ex19_tut_2.out.html">[Expected output]</a>
    </pre>
    Note this requires many fewer iterations than the default solver
  </li>
  <li> Run with the same options but use algebraic multigrid (hypre's BoomerAMG) as the linear solver
    <pre>
     mpiexec -n 16 ./ex19 -da_refine 6 -snes_monitor -ksp_monitor -snes_view -pc_type hypre     <a href="../src/snes/examples/tutorials/output/ex19_tut_3.out.html">[Expected output]</a>
    </pre>
    Note this requires many fewer iterations than the default solver but requires more iterations than geometric multigrid
  </li>
  <li> Run on 4 processors with the default linear solver and profile the run
    <pre>
     mpiexec -n 4 ./ex19 -da_refine 6 -log_summary     <a href="../src/snes/examples/tutorials/output/ex19_tut_4.out.html">[Expected output]</a>
    </pre>
    Search for the line beginning with SNESSolve, the fourth column gives the time for the nonlinear solve.
  </li>
  <li> Run on 4 processors with the geometric multigrid linear solver and profile the run
    <pre>
     mpiexec -n 4 ./ex19 -da_refine 6 -log_summary -pc_type mg     <a href="../src/snes/examples/tutorials/output/ex19_tut_5.out.html">[Expected output]</a>
    </pre>
    Compare the runtime for SNESSolve to the case with the default solver
  </li>
  <li> Run on 16 processors with the default linear solver and profile the run
    <pre>
     mpiexec -n 16 ./ex19 -da_refine 6 -log_summary     <a href="../src/snes/examples/tutorials/output/ex19_tut_6.out.html">[Expected output]</a>
    </pre>
    Compare the runtime for SNESSolve to the 4 processor case with the default solver, what is the speedup?
  </li>
  <li> Run on 16 processors with the geometric multigrid linear solver and profile the run
    <pre>
     mpiexec -n 16 ./ex19 -da_refine 6 -log_summary -pc_type mg     <a href="../src/snes/examples/tutorials/output/ex19_tut_7.out.html">[Expected output]</a>
    </pre>
    Compare the runtime for the SNESSolve to the 4 processor case with multigrid, what is the speedup? 
    Why is the speedup for multigrid lower than the speedup for the default solver? 
  </li>
</ul>

<!--------------------------------------------------------------------------->

<h2>Example 4: Linear Stokes-type PDE on a structured grid</h2>

<p>WHAT THIS EXAMPLE DEMONSTRATES:</p>
<ul>
  <li>Handling a 3d structured grid</li>
  <li>Controlling linear solver options</li>
  <li>Selecting composible preconditioners</li>
  <li>Solving a Stokes problem</li>
  <li>Adding your own problem specific visualization</li>
</ul>

<p>FURTHER DETAILS:</p>
<ul>
  <li><a href="../src/ksp/ksp/examples/tutorials/ex42.c.html">Mathematical description of the problem</a></li>
  <li><a href="../src/ksp/ksp/examples/tutorials/ex42.c.html#line2059">main program source code</a></li>
  <li><a href="../src/ksp/ksp/examples/tutorials/ex42.c.html#line819">physics source code</a></li>
</ul>

<p>DO THE FOLLOWING:</p>
<ul>
  <li> Compile src/ksp/ksp/examples/tutorials/ex42.c
    <pre>
      cd petsc/src/ksp/ksp/examples/tutorials
      make ex42
  </pre>
  </li>
  <li> Solve with the default solver
    <pre>
     mpiexec -n 4 ./ex42  -stokes_ksp_monitor   <a href="../src/ksp/ksp/examples/tutorials/output/ex42_tut_1.out.html">[Expected output]</a>
    </pre>
    Note the poor convergence for even a very small problem
  </li>
  <li> Solve with a solver appropriate for Stoke's problems -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur
    <pre>
     mpiexec -n 4 ./ex42  -stokes_ksp_monitor -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur  <a href="../src/ksp/ksp/examples/tutorials/output/ex42_tut_2.out.html">[Expected output]</a>
    </pre>
  </li>
  <li> Solve with a finer mesh
    <pre>
     mpiexec -n 16 ./ex42  -mx 40 -stokes_ksp_monitor  -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur  <a href="../src/ksp/ksp/examples/tutorials/output/ex42_tut_3.out.html">[Expected output]</a>
    </pre>
  </li>
</ul>

<!--------------------------------------------------------------------------->

<h2>Example 5: Nonlinear time dependent PDE on Unstructured Grid</h2>

<p>WHAT THIS EXAMPLE DEMONSTRATES:</p>
<ul>
  <li>Changing the default ODE integrator</li>
  <li>Handling unstructured grids</li>
  <li>Registering your own interchangeable physics and algorithm modules</li>
</ul>

<p>FURTHER DETAILS:</p>
<ul>
  <li><a href="../src/ts/examples/tutorials/ex11.c.html">Mathematical description of the problem</a></li>
  <li><a href="../src/ts/examples/tutorials/ex11.c.html#line1403">main program source code</a></li>
  <li><a href="../src/ts/examples/tutorials/ex11.c.html#line186">source code of physics modules</a></li>
</ul>

<p>DO THE FOLLOWING:</p>
<ul>
  <li> Compile src/ts/examples/tutorials/ex11.c
    <pre>
      cd petsc/src/ts/examples/tutorials
      make ex11
  </pre>
  <li>Run simple advection through a tiny hybrid mesh
    <pre>
     mpiexec -n 1 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo      <a href="../src/ts/examples/tutorials/output/ex11_tut_1.out.html">[Expected output]</a>
    </pre>
  </li>
  <li>Run simple advection through a small mesh with a Rosenbrock-W solver
    <pre>
     mpiexec -n 1 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw      <a href="../src/ts/examples/tutorials/output/ex11_tut_2.out.html">[Expected output]</a>
    </pre>
  </li>
  <li>Run simple advection through a larger quadrilateral mesh of an annulus with least squares reconstruction and no limiting, monitoring the error
    <pre>
     mpiexec -n 16 ./ex11 -f ${PETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin <a href="../src/ts/examples/tutorials/output/ex11_tut_4.out.html">[Expected output]</a>
    </pre>
    Compare turning to the error after turning off reconstruction.
  </li>
  <li>Run shallow water on the larger mesh with least squares reconstruction and minmod limiting, monitoring water Height (integral is conserved) and Energy (not conserved)
    <pre>
     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 <a href="../src/ts/examples/tutorials/output/ex11_tut_5.out.html">[Expected output]</a>
    </pre>
  </li>
</ul>

</body></html>