[petsc-users] do SNESVI objects support multigrid?

Ed Bueler elbueler at alaska.edu
Fri Feb 27 22:28:58 CST 2015


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)
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150227/81a6022c/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: out.ex55err
Type: application/octet-stream
Size: 310391 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150227/81a6022c/attachment-0001.obj>


More information about the petsc-users mailing list