[petsc-users] do SNESVI objects support multigrid?

Barry Smith bsmith at mcs.anl.gov
Sun Mar 1 16:26:30 CST 2015


   Ed,

   I have fixed ex9 so it works super nicely with multigrid in the branch pr276/elbueler/fix-ex9/master and next

$ ./ex9 -da_refine 5 -pc_type mg -snes_monitor  -ksp_monitor_true_residual -snes_vi_monitor
setup done: grid  Mx,My = 321,321  with spacing  dx,dy = 0.013,0.013
  0 SNES Function norm 1.691947292888e+01 
  0 SNES VI Function norm 1.691947292888e+01 Active lower constraints 19981/20079 upper constraints 0/0 Percent of total 0.193913 Percent of bounded 0
    0 KSP preconditioned resid norm 2.107314279885e+01 true resid norm 1.691947292888e+01 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP preconditioned resid norm 2.771357241740e-01 true resid norm 2.720618543577e-01 ||r(i)||/||b|| 1.607980671155e-02
    2 KSP preconditioned resid norm 4.702355133797e-03 true resid norm 8.162965549257e-03 ||r(i)||/||b|| 4.824598014118e-04
    3 KSP preconditioned resid norm 1.104885629534e-04 true resid norm 2.674390325875e-04 ||r(i)||/||b|| 1.580658178370e-05

  The problem was your scaling of the Dirichlet boundary conditions made geometric multigrid not work. I first turned off the VI part of the example and got terrible multigrid convergence and this led me to the fix. So it shows while debugging with constraints it is always good to debug without them first.


  Barry


> On Feb 28, 2015, at 9:00 PM, Ed Bueler <elbueler at alaska.edu> wrote:
> 
> Barry --
> 
> >  If you would like I could take a look at trying to get it working again as it did before.
> 
> I would like that.
> 
> I think the availability of robust SNESVI is the key to hugely improving how we model conservation in thin layers (e.g. ice sheets, sea ice, surface hydrology, ...); see http://www2.gi.alaska.edu/snowice/glaciers/iceflow/bueler-liwg2015.pdf, a short and easy talk.  (View with acroread for the movies.)
> 
> Currently everyone in ice sheets (specifically) is doing explicit time stepping, in substantial part because that's how they think about the free boundary of the ice sheet, where the melt exceeds snowfall+flow.  By changing to thinking of each time-step as a VI you get to add "thickness of ice is nonnegative", as a constraint, to the problem, instead of ad hoc time-splitting discrete ways of moving the free boundary.
> 
> Ed
> 
> 
> On Sat, Feb 28, 2015 at 10:40 AM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> 
>   Ed,
> 
>    At one time it did. However it has suffered from neglect and bit rot over several years. Part of the reason is that we have had big plans to "do it much better" and hence neglected the old code; while we never "did it much better".
> 
>    So algorithmically what we implemented as an active set method for the VI with multigrid solving and I feel algorithmically it is a pretty good method (for problems where multigrid works well without the VI). My implementation was a bit ad hoc and didn't mesh great with the rest of PETSc which is why it hasn't been maintained properly. If you would like I could take a look at trying to get it working again as it did before.
> 
>   Barry
> 
> 
> 
> > On Feb 27, 2015, at 10:28 PM, Ed Bueler <elbueler at alaska.edu> wrote:
> >
> > Dear Petsc --
> >
> > I am confused on whether  -pc_type mg  is supported for SNESVI solvers.  The error messages sure aren't helping.
> >
> > I am using petsc-maint.
> >
> > First I looked at couple of simplest examples in src/snes/examples/tutorials, namely ex9 and ex54.  I got mixed results (below), namely crashes with 'nonconforming object sizes' or 'zero pivot' errors in 3 of 4 cases.  But I am not sure from the very limited docs what is supposed to work.
> >
> > The documented ex54 run itself (see line 12 of ex54.c) just hangs/stagnates for me:
> >
> > $ ./ex54 -pc_type mg -pc_mg_galerkin -T .01 -da_grid_x 65 -da_grid_y 65 -pc_mg_levels 4 -ksp_type fgmres -snes_atol 1.e-14 -mat_no_inode -snes_vi_monitor -snes_monitor
> >   0 SNES Function norm 6.177810506000e-03
> >   0 SNES VI Function norm 6.177810506000e-03 Active lower constraints 0/0 upper constraints 0/0 Percent of total 0 Percent of bounded 0
> >   0 SNES Function norm 6.177810506000e-03
> >   0 SNES VI Function norm 6.177810506000e-03 Active lower constraints 0/0 upper constraints 0/0 Percent of total 0 Percent of bounded 0
> >   0 SNES Function norm 6.177810506000e-03
> >   0 SNES VI Function norm 6.177810506000e-03 Active lower constraints 0/0 upper constraints 0/0 Percent of total 0 Percent of bounded 0
> >   0 SNES Function norm 6.177810506000e-03
> >   0 SNES VI Function norm 6.177810506000e-03 Active lower constraints 0/0 upper constraints 0/0 Percent of total 0 Percent of bounded 0
> >   0 SNES Function norm 6.177810506000e-03
> >   0 SNES VI Function norm 6.177810506000e-03 Active lower constraints 0/0 upper constraints 0/0 Percent of total 0 Percent of bounded 0
> > ...
> >
> > The wordy documented run of ex55.c stagnates and then seg faults for me (see attachment):
> >
> > ~/petsc-maint/src/snes/examples/tutorials[maint*]$ ./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 -ksp_rtol 1.e-8 &> out.ex55err
> >
> > These examples all seem to perform o.k. with default (and non-multigrid) options.
> >
> > Thanks for help!
> >
> > Ed
> >
> >
> > errors from simple ex9, ex54 runs below:
> >
> > ~/petsc-maint/src/snes/examples/tutorials[maint]$ ./ex9 -da_refine 2 -snes_type vinewtonrsls -pc_type mg
> > setup done: square       side length = 4.000
> >             grid               Mx,My = 41,41
> >             spacing            dx,dy = 0.100,0.100
> > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > [0]PETSC ERROR: Nonconforming object sizes
> > [0]PETSC ERROR: Mat mat,Vec x: global dim 1388 1681
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.5.3, unknown
> > [0]PETSC ERROR: ./ex9 on a linux-c-opt named bueler-leopard by ed Fri Feb 27 21:02:20 2015
> > [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-fblaslapack --download-mpich --with-debugging=0
> > [0]PETSC ERROR: #1 MatMultTranspose() line 2232 in /home/ed/petsc-maint/src/mat/interface/matrix.c
> > [0]PETSC ERROR: #2 MatRestrict() line 7529 in /home/ed/petsc-maint/src/mat/interface/matrix.c
> > [0]PETSC ERROR: #3 DMRestrictHook_SNESVecSol() line 480 in /home/ed/petsc-maint/src/snes/interface/snes.c
> > [0]PETSC ERROR: #4 DMRestrict() line 2022 in /home/ed/petsc-maint/src/dm/interface/dm.c
> > [0]PETSC ERROR: #5 PCSetUp_MG() line 699 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #6 PCSetUp() line 902 in /home/ed/petsc-maint/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #7 KSPSetUp() line 306 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #8 SNESSolve_VINEWTONRSLS() line 506 in /home/ed/petsc-maint/src/snes/impls/vi/rs/virs.c
> > [0]PETSC ERROR: #9 SNESSolve() line 3743 in /home/ed/petsc-maint/src/snes/interface/snes.c
> > [0]PETSC ERROR: #10 main() line 122 in /home/ed/petsc-maint/src/snes/examples/tutorials/ex9.c
> > [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> > application called MPI_Abort(MPI_COMM_WORLD, 60) - process 0
> > [unset]: aborting job:
> > application called MPI_Abort(MPI_COMM_WORLD, 60) - process 0
> >
> >
> > ~/petsc-maint/src/snes/examples/tutorials[maint]$ ./ex9 -da_refine 2 -snes_type vinewtonssls -pc_type mg
> > setup done: square       side length = 4.000
> >             grid               Mx,My = 41,41
> >             spacing            dx,dy = 0.100,0.100
> > number of Newton iterations = 8; result = CONVERGED_FNORM_RELATIVE
> > errors:    av |u-uexact|  = 2.909e-04
> >            |u-uexact|_inf = 1.896e-03
> >
> >
> > ~/petsc-maint/src/snes/examples/tutorials[maint]$ ./ex54 -snes_monitor -da_refine 2 -snes_type vinewtonssls -pc_type mg
> >   0 SNES Function norm 3.156635589354e-02
> > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > [0]PETSC ERROR: Zero pivot in LU factorization: http://www.mcs.anl.gov/petsc/documentation/faq.html#ZeroPivot
> > [0]PETSC ERROR: Zero pivot, row 0
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.5.3, unknown
> > [0]PETSC ERROR: ./ex54 on a linux-c-opt named bueler-leopard by ed Fri Feb 27 21:02:43 2015
> > [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-fblaslapack --download-mpich --with-debugging=0
> > [0]PETSC ERROR: #1 PetscKernel_A_gets_inverse_A_2() line 50 in /home/ed/petsc-maint/src/mat/impls/baij/seq/dgefa2.c
> > [0]PETSC ERROR: #2 MatSOR_SeqAIJ_Inode() line 2806 in /home/ed/petsc-maint/src/mat/impls/aij/seq/inode.c
> > [0]PETSC ERROR: #3 MatSOR() line 3643 in /home/ed/petsc-maint/src/mat/interface/matrix.c
> > [0]PETSC ERROR: #4 PCApply_SOR() line 35 in /home/ed/petsc-maint/src/ksp/pc/impls/sor/sor.c
> > [0]PETSC ERROR: #5 PCApply() line 440 in /home/ed/petsc-maint/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #6 KSP_PCApply() line 230 in /home/ed/petsc-maint/include/petsc-private/kspimpl.h
> > [0]PETSC ERROR: #7 KSPInitialResidual() line 56 in /home/ed/petsc-maint/src/ksp/ksp/interface/itres.c
> > [0]PETSC ERROR: #8 KSPSolve_GMRES() line 234 in /home/ed/petsc-maint/src/ksp/ksp/impls/gmres/gmres.c
> > [0]PETSC ERROR: #9 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #10 KSPSolve_Chebyshev() line 368 in /home/ed/petsc-maint/src/ksp/ksp/impls/cheby/cheby.c
> > [0]PETSC ERROR: #11 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #12 PCMGMCycle_Private() line 19 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #13 PCMGMCycle_Private() line 48 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #14 PCApply_MG() line 337 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #15 PCApply() line 440 in /home/ed/petsc-maint/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #16 KSP_PCApply() line 230 in /home/ed/petsc-maint/include/petsc-private/kspimpl.h
> > [0]PETSC ERROR: #17 KSPInitialResidual() line 63 in /home/ed/petsc-maint/src/ksp/ksp/interface/itres.c
> > [0]PETSC ERROR: #18 KSPSolve_GMRES() line 234 in /home/ed/petsc-maint/src/ksp/ksp/impls/gmres/gmres.c
> > [0]PETSC ERROR: #19 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #20 SNESSolve_VINEWTONSSLS() line 317 in /home/ed/petsc-maint/src/snes/impls/vi/ss/viss.c
> > [0]PETSC ERROR: #21 SNESSolve() line 3743 in /home/ed/petsc-maint/src/snes/interface/snes.c
> > [0]PETSC ERROR: #22 main() line 98 in /home/ed/petsc-maint/src/snes/examples/tutorials/ex54.c
> > [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> > application called MPI_Abort(MPI_COMM_WORLD, 71) - process 0
> > [unset]: aborting job:
> > application called MPI_Abort(MPI_COMM_WORLD, 71) - process 0
> >
> >
> > ~/petsc-maint/src/snes/examples/tutorials[maint]$ ./ex54 -snes_monitor -da_refine 2 -snes_type vinewtonrsls -pc_type mg
> >   0 SNES Function norm 3.160548858489e-02
> > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > [0]PETSC ERROR: Zero pivot in LU factorization: http://www.mcs.anl.gov/petsc/documentation/faq.html#ZeroPivot
> > [0]PETSC ERROR: Zero pivot, row 0
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.5.3, unknown
> > [0]PETSC ERROR: ./ex54 on a linux-c-opt named bueler-leopard by ed Fri Feb 27 21:02:48 2015
> > [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --download-fblaslapack --download-mpich --with-debugging=0
> > [0]PETSC ERROR: #1 PetscKernel_A_gets_inverse_A_2() line 50 in /home/ed/petsc-maint/src/mat/impls/baij/seq/dgefa2.c
> > [0]PETSC ERROR: #2 MatSOR_SeqAIJ_Inode() line 2806 in /home/ed/petsc-maint/src/mat/impls/aij/seq/inode.c
> > [0]PETSC ERROR: #3 MatSOR() line 3643 in /home/ed/petsc-maint/src/mat/interface/matrix.c
> > [0]PETSC ERROR: #4 PCApply_SOR() line 35 in /home/ed/petsc-maint/src/ksp/pc/impls/sor/sor.c
> > [0]PETSC ERROR: #5 PCApply() line 440 in /home/ed/petsc-maint/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #6 KSP_PCApply() line 230 in /home/ed/petsc-maint/include/petsc-private/kspimpl.h
> > [0]PETSC ERROR: #7 KSPInitialResidual() line 56 in /home/ed/petsc-maint/src/ksp/ksp/interface/itres.c
> > [0]PETSC ERROR: #8 KSPSolve_GMRES() line 234 in /home/ed/petsc-maint/src/ksp/ksp/impls/gmres/gmres.c
> > [0]PETSC ERROR: #9 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #10 KSPSolve_Chebyshev() line 368 in /home/ed/petsc-maint/src/ksp/ksp/impls/cheby/cheby.c
> > [0]PETSC ERROR: #11 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #12 PCMGMCycle_Private() line 19 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #13 PCMGMCycle_Private() line 48 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #14 PCApply_MG() line 337 in /home/ed/petsc-maint/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #15 PCApply() line 440 in /home/ed/petsc-maint/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #16 KSP_PCApply() line 230 in /home/ed/petsc-maint/include/petsc-private/kspimpl.h
> > [0]PETSC ERROR: #17 KSPInitialResidual() line 63 in /home/ed/petsc-maint/src/ksp/ksp/interface/itres.c
> > [0]PETSC ERROR: #18 KSPSolve_GMRES() line 234 in /home/ed/petsc-maint/src/ksp/ksp/impls/gmres/gmres.c
> > [0]PETSC ERROR: #19 KSPSolve() line 460 in /home/ed/petsc-maint/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #20 SNESSolve_VINEWTONRSLS() line 536 in /home/ed/petsc-maint/src/snes/impls/vi/rs/virs.c
> > [0]PETSC ERROR: #21 SNESSolve() line 3743 in /home/ed/petsc-maint/src/snes/interface/snes.c
> > [0]PETSC ERROR: #22 main() line 98 in /home/ed/petsc-maint/src/snes/examples/tutorials/ex54.c
> > [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> > application called MPI_Abort(MPI_COMM_WORLD, 71) - process 0
> > [unset]: aborting job:
> > application called MPI_Abort(MPI_COMM_WORLD, 71) - process 0
> >
> >
> > --
> > Ed Bueler
> > Dept of Math and Stat and Geophysical Institute
> > University of Alaska Fairbanks
> > Fairbanks, AK 99775-6660
> > 301C Chapman and 410D Elvey
> > 907 474-7693 and 907 474-7199  (fax 907 474-5394)
> > <out.ex55err>
> 
> 
> 
> 
> -- 
> Ed Bueler
> Dept of Math and Stat and Geophysical Institute
> University of Alaska Fairbanks
> Fairbanks, AK 99775-6660
> 301C Chapman and 410D Elvey
> 907 474-7693 and 907 474-7199  (fax 907 474-5394)



More information about the petsc-users mailing list