[petsc-users] pc gamg did not converge in sinv where it used to

Mark Adams mfadams at lbl.gov
Tue Feb 7 16:57:50 CST 2017


This is failing on a test for a data cache being NULL. It is not. It is in
the reset routine, which might not be used by many people. This code is
pretty stable as far as I know, but maybe someone changed something in
here. This is probably a false positive -- you could just comment out that
line ("this should not happen, cleaned up in SetUp") in gamg.c:38

This data is freed in a few places with something like:

ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);

But, not set to NULL.

Perhaps PetscFree should be called with the address of the pointer so that
it can NULL it or explicitly zeroed out.

PETSc: what is the proper way to do this?

Mark


On Tue, Feb 7, 2017 at 4:51 PM, Denis Davydov <davydden at gmail.com> wrote:

> Dear all,
>
> I rerun calculations (unit tests) which used to work with slightly older
> versions of PETSc/SLEPc (a year ago, or so) and
> see the  "KSPSolve has not converged” error for shift and invert
> transformation with gamg preconditioner (below).
> Shifted matrix is SPD, could have bad condition number but it worked
> before.
>
> Judging from a pre-last line in the log below, something is wrong with
> GAMG or its interaction with SLEPc.
>
> Regards,
> Denis.
>
> [0]PETSC ERROR: Configure options --prefix=/Users/davydden/spack
> /opt/spack/darwin-sierra-x86_64/clang-8.0.0-apple/petsc-3.
> 7.5-e35uxj5pf7ndtedbnzvq7dot5oie3psc --with-ssl=0 --with-x=0
> --download-c2html=0 --download-sowing=0 --download-hwloc=0 --with-mpi=1
> --with-mpi-dir=/Users/davydden/spack/opt/spack/darwin-
> sierra-x86_64/clang-8.0.0-apple/openmpi-2.0.2-oos4cxprn2bislzhc3rbc3lo4dadginw
> --with-precision=double --with-scalar-type=real --with-shared-libraries=1
> --with-debugging=0 --with-64-bit-indices=0 --with-blas-lapack-lib=/Users/
> davydden/spack/opt/spack/darwin-sierra-x86_64/clang-8.0.0-
> apple/openblas-0.2.19-xg7etkjyo7xjnvuojsdc2xoixesxoerh/lib/libopenblas.dylib
> --with-metis=1 --with-metis-dir=/Users/davydden/spack/opt/spack/darwin-
> sierra-x86_64/clang-8.0.0-apple/metis-5.1.0-c2muo3b6k5aea2743g34s35u2ve23yxq
> --with-boost=0 --with-hdf5=1 --with-hdf5-dir=/Users/davydde
> n/spack/opt/spack/darwin-sierra-x86_64/clang-8.0.0-
> apple/hdf5-1.10.0-patch1-kwkdkp5mo4olqpkitojab3xg3xhl6gfz --with-hypre=1
> --with-hypre-dir=/Users/davydden/spack/opt/spack/darwin-
> sierra-x86_64/clang-8.0.0-apple/hypre-2.11.1-dprfz6t4t24wesj74dcfdtajwm5n5qzu
> --with-parmetis=1 --with-parmetis-dir=/Users/dav
> ydden/spack/opt/spack/darwin-sierra-x86_64/clang-8.0.0-
> apple/parmetis-4.0.3-tm5dd6fqnbw6hiepajhkvmqk42xrrgos --with-mumps=1
> --with-mumps-dir=/Users/davydden/spack/opt/spack/darwin-
> sierra-x86_64/clang-8.0.0-apple/mumps-5.0.2-4pwiz7bplhnmk2lfvg746u4pq4z5wuc3
> --with-scalapack=1 --with-scalapack-dir=/Users/da
> vydden/spack/opt/spack/darwin-sierra-x86_64/clang-8.0.0-
> apple/netlib-scalapack-2.0.2-gu4tjf4bjvtia3ohzjqbickedmywqslk
> --with-superlu_dist-include=/Users/davydden/spack/opt/spack/
> darwin-sierra-x86_64/clang-8.0.0-apple/superlu-dist-5.1.1-h
> ut4gapp5v6qzvjlqahattdvr3tyynoy/include --with-superlu_dist-lib=/Users
> /davydden/spack/opt/spack/darwin-sierra-x86_64/clang-8.
> 0.0-apple/superlu-dist-5.1.1-hut4gapp5v6qzvjlqahattdvr3tyynoy/lib/libsuperlu_dist.a
> --with-superlu_dist=1
> [0]PETSC ERROR: #1 KSPSolve() line 847 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-zUUkMX/petsc-3.7.5/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #2 PCGAMGOptProlongator_AGG() line 1227 in
> /private/var/folders/5k/sqpp24tx3ylds4fgm13pfht00000gn/T/
> davydden/spack-stage/spack-stage-zUUkMX/petsc-3.7.5/src/ksp/
> pc/impls/gamg/agg.c
> [0]PETSC ERROR: #3 PCSetUp_GAMG() line 581 in
> /private/var/folders/5k/sqpp24tx3ylds4fgm13pfht00000gn/T/
> davydden/spack-stage/spack-stage-zUUkMX/petsc-3.7.5/src/ksp/
> pc/impls/gamg/gamg.c
> [0]PETSC ERROR: #4 PCSetUp() line 968 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-zUUkMX/petsc-3.7.5/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #5 KSPSetUp() line 390 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-zUUkMX/petsc-3.7.5/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #6 STSetUp_Sinvert() line 153 in
> /private/var/folders/5k/sqpp24tx3ylds4fgm13pfht00000gn/T/
> davydden/spack-stage/spack-stage-5WvuFP/slepc-3.7.3/src/sys/
> classes/st/impls/sinvert/sinvert.c
> [0]PETSC ERROR: #7 STSetUp() line 301 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-5WvuFP/slepc-3.7.3/src/sys/classes/st/interface/stsolve.c
> [0]PETSC ERROR: #8 EPSSetUp() line 205 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-5WvuFP/slepc-3.7.3/src/eps/interface/epssetup.c
> [0]PETSC ERROR: #9 EPSSolve() line 89 in /private/var/folders/5k/sqpp24
> tx3ylds4fgm13pfht00000gn/T/davydden/spack-stage/spack-sta
> ge-5WvuFP/slepc-3.7.3/src/eps/interface/epssolve.c
> ***[0]PCReset_GAMG this should not happen, cleaned up in SetUp
> ERROR: Uncaught exception in MPI_InitFinalize on proc 0. Skipping
> MPI_Finalize() to avoid a deadlock.
>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20170207/a4cd8391/attachment.html>


More information about the petsc-users mailing list