[petsc-users] [petsc-maint] GAMG with PETSc

Mark Adams mfadams at lbl.gov
Fri Aug 12 10:03:30 CDT 2016


Yes, we should check the type.

#define MATAIJ             "aij"
#define MATSEQAIJ          "seqaij"
#define MATMPIAIJ          "mpiaij"

I assume that MATAIJ is mapped to SEQ or MPI, and so just need to check for
the latter two.  I'll do that now.

Thanks,
Mark

On Fri, Aug 12, 2016 at 7:02 AM, Lawrence Mitchell <
lawrence.mitchell at imperial.ac.uk> wrote:

> [Added petsc-maint to cc, since I think this is an actual bug]
>
> > On 12 Aug 2016, at 01:16, Gaetan Kenway <gaetank at gmail.com> wrote:
> >
> > Hi
> >
> > I'm attempting to try out using GAMG for a preconditioner for my
> compressible CFD problem. However, I'm getting segfaults when trying to run
> the code. The code is based on ksp ex23.c which is attached. It just reads
> in two precomputed matrices (the actual jacobian and an approximate
> jacobian used to build the PC) and solves with a RHS of ones.
> >
> > My normal approach to solving the system is with ASM+ILU. With the
> following options, everything works fine.
>
> This appears to be a problem that GAMG doesn't work with BAIJ matrices.
>
> But there is no checking of the input matrix type anywhere.
>
> For example, with a debug PETSc:
>
> cd src/ksp/ksp/examples/tutorials
> make ex23
> ./ex23 -pc_type gamg -mat_type baij
>
> [0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> [0]PETSC ERROR: Object is in wrong state
> [0]PETSC ERROR: Must call MatXXXSetPreallocation() or MatSetUp() on
> argument 1 "mat" before MatSetValues()
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
> for trouble shooting.
> [0]PETSC ERROR: Petsc Development GIT revision: v3.7.3-1123-g21143a8  GIT
> Date: 2016-08-08 17:24:17 -0700
> [0]PETSC ERROR: ./ex23 on a arch-darwin-c-dbg named yam-laptop.local by
> lmitche1 Fri Aug 12 11:51:17 2016
> [0]PETSC ERROR: Configure options --download-chaco=1 --download-ctetgen=1
> --download-exodusii=1 --download-hdf5=1 --download-hypre=1
> --download-metis=1 --download-ml=1 --download-mumps=1 --download-netcdf=1
> --download-parmetis=1 --download-ptscotch=1 --download-scalapack=1
> --download-triangle=1 --with-c2html=0 --with-debugging=1
> --with-shared-libraries=1 PETSC_ARCH=arch-darwin-c-dbg
> [0]PETSC ERROR: #1 MatSetValues() line 1195 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/mat/interface/matrix.c
> [0]PETSC ERROR: #2 PCGAMGFilterGraph() line 342 in
> /Users/lmitche1/Documents/work/src/deps/petsc/src/ksp/pc/impls/gamg/util.c
> [0]PETSC ERROR: #3 PCGAMGGraph_AGG() line 908 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/pc/impls/gamg/agg.c
> [0]PETSC ERROR: #4 PCSetUp_GAMG() line 525 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/pc/impls/gamg/gamg.c
> [0]PETSC ERROR: #5 PCSetUp() line 968 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #6 KSPSetUp() line 393 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #7 KSPSolve() line 602 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #8 main() line 158 in /Users/lmitche1/Documents/
> work/src/deps/petsc/src/ksp/ksp/examples/tutorials/ex23.c
> [0]PETSC ERROR: PETSc Option Table entries:
> [0]PETSC ERROR: -mat_type baij
> [0]PETSC ERROR: -pc_type gamg
> [0]PETSC ERROR: ----------------End of Error Message -------send entire
> error message to petsc-maint at mcs.anl.gov----------
>
> Looking in PCGAMGFilterGraph I see:
>
>   ierr = MatGetType(Gmat,&mtype);CHKERRQ(ierr);
>   ierr = MatCreate(comm, &tGmat);CHKERRQ(ierr);
>   ierr = MatSetSizes(tGmat,nloc,nloc,MM,MM);CHKERRQ(ierr);
>   ierr = MatSetBlockSizes(tGmat, 1, 1);CHKERRQ(ierr);
>   ierr = MatSetType(tGmat, mtype);CHKERRQ(ierr);
>   ierr = MatSeqAIJSetPreallocation(tGmat,0,d_nnz);CHKERRQ(ierr);
>   ierr = MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
>
>
> ...
>           ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
> CHKERRQ(ierr);
>
>
> So if Gmat is neither SEQAIJ nor MPIAIJ, then no preallocation has
> happened (and MatSetUp is not called).
>
> Fixing the few instances here by just changing the type of these matrices
> to AIJ.  One runs into to the problem that creating the coarse grid
> operators doesn't work, since MatMatMult and friends do not exist for BAIJ
> matrices.
>
> I guess GAMG could MatConvert from BAIJ to AIJ (but this uses extra
> memory).
>
> But it should probably barf with a comprehensible error message.
>
> Thoughts?
>
> Lawrence
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160812/f9cbe78c/attachment.html>


More information about the petsc-users mailing list