[petsc-dev] programming with PETSc options

Barry Smith bsmith at mcs.anl.gov
Sat Jul 9 22:44:38 CDT 2011


  Jungho and I have implemented a simple geometric multigrid method for a saddle point problem variational inequalities that works (so far) surprisingly well using the PETSc options programming language. Because it is a fun example of programming a compossible hierarchical nested solver I thought I would share the command line.  This is essentially multigrid on the entire matrix with a Schur complement smoother on each level. Changing the options a bit would give you a Schur complement solver on the entire matrix with multigrid on the A00 block.

./ex55 -ksp_type fgmres -pc_type mg -mg_levels_ksp_type fgmres -mg_levels_pc_type fieldsplit -mg_levels_pc_fieldsplit_detect_saddle_point -mg_levels_pc_fieldsplit_type schur -mg_levels_pc_fieldsplit_factorization_type full -mg_levels_pc_fieldsplit_schur_precondition user -mg_levels_fieldsplit_1_ksp_type gmres -mg_levels_fieldsplit_1_pc_type none -mg_levels_fieldsplit_0_ksp_type preonly -mg_levels_fieldsplit_0_pc_type sor -mg_levels_fieldsplit_0_pc_sor_forward -snes_vi_monitor -ksp_monitor_true_residual -pc_mg_levels 5 -pc_mg_galerkin -mg_levels_ksp_monitor -mg_levels_fieldsplit_ksp_monitor -mg_levels_ksp_max_it 2 -mg_levels_fieldsplit_ksp_max_it 5 -snes_atol 1.e-11 -mg_coarse_ksp_type preonly -mg_coarse_pc_type svd -da_grid_x 65 -da_grid_y 65


Run flexible GMRES with 5 levels of multigrid as the preconditioner
--------------
./ex55 -ksp_type fgmres -pc_type mg -pc_mg_levels 5 -da_grid_x 65 -da_grid_y 65   


Use the Galerkin process to compute the coarse grid operators
----------------
-pc_mg_galerkin


Use SVD as the coarse grid saddle point solver
---------------
-mg_coarse_ksp_type preonly -mg_coarse_pc_type svd


Use 2 iterations of flexible GMRES with a Schur complement preconditioner as the smoother
---------------
-mg_levels_ksp_type fgmres -mg_levels_ksp_max_it 2  -mg_levels_pc_type fieldsplit -mg_levels_pc_fieldsplit_detect_saddle_point -mg_levels_pc_fieldsplit_type schur -mg_levels_pc_fieldsplit_factorization_type full -mg_levels_pc_fieldsplit_schur_precondition diag


Iterate on the Shur complement with GMRES with no preconditioner for 5 iterations in the Schur complement smoother
--------------
-mg_levels_fieldsplit_1_ksp_type gmres -mg_levels_fieldsplit_1_pc_type none  -mg_levels_fieldsplit_ksp_max_it 5


Approximate the action of the Shur complement by using only the lower diagonal part of A00 
--------------
 -mg_levels_fieldsplit_0_ksp_type preonly -mg_levels_fieldsplit_0_pc_type sor -mg_levels_fieldsplit_0_pc_sor_forward 


The bulk of the computation is just forward Gauss-Seidel on the non-saddle point part of the matrix. (I've removed the monitoring of the inner linear iterations) 

 0 SNES VI Function norm 3.976439267245e-05 Active constraints 2315
    0 KSP unpreconditioned resid norm 3.976439267245e-05 true resid norm 3.976439267245e-05 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP unpreconditioned resid norm 3.543930653675e-05 true resid norm 3.543930653676e-05 ||r(i)||/||b|| 8.912321842479e-01
    2 KSP unpreconditioned resid norm 4.005676188871e-08 true resid norm 4.005676188723e-08 ||r(i)||/||b|| 1.007352538166e-03
    3 KSP unpreconditioned resid norm 1.073296617426e-10 true resid norm 1.073296621498e-10 ||r(i)||/||b|| 2.699139982696e-06
  1 SNES VI Function norm 8.882771573091e-03 Active constraints 2201
    0 KSP unpreconditioned resid norm 8.882771573091e-03 true resid norm 8.882771573091e-03 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP unpreconditioned resid norm 2.591449737029e-07 true resid norm 2.591449737028e-07 ||r(i)||/||b|| 2.917388695301e-05
    2 KSP unpreconditioned resid norm 6.018436991709e-10 true resid norm 6.018436991512e-10 ||r(i)||/||b|| 6.775404435417e-08
  2 SNES VI Function norm 6.604821607805e-10 Active constraints 2203
    0 KSP unpreconditioned resid norm 6.604820670139e-10 true resid norm 6.604820670139e-10 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP unpreconditioned resid norm 5.516796916907e-10 true resid norm 5.516796916907e-10 ||r(i)||/||b|| 8.352682370089e-01
    2 KSP unpreconditioned resid norm 2.370052538357e-12 true resid norm 2.370052538361e-12 ||r(i)||/||b|| 3.588367734308e-03
    3 KSP unpreconditioned resid norm 4.953732708209e-14 true resid norm 4.953732708254e-14 ||r(i)||/||b|| 7.500177454704e-05
    4 KSP unpreconditioned resid norm 2.777240109831e-16 true resid norm 2.777239994389e-16 ||r(i)||/||b|| 4.204868130554e-07
  3 SNES VI Function norm 3.519585974269e-13 Active constraints 2204

  

   Eventually we'll program all our solvers with command line options :-)


   Barry




More information about the petsc-dev mailing list