[petsc-dev] Segmentation fault GAMG + MatView

Hong hzhang at mcs.anl.gov
Thu Feb 9 16:53:17 CST 2017


I'll check it tomorrow and fix it.
Hong

On Thu, Feb 9, 2017 at 4:27 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

>
>    When called with a NULL set of numerical values
> MatCreateSeqAIJWithArrays() should turn off the function pointers for
> things that don't make sense for matrices without numerical values, for
> example MatMult() and MatView(), in fact most function pointers don't make
> sense and should be turned off.
> This could be done with a
> if (!a) {
>    ...
> }
> after the call to MatSetType().
>
>    Barry
>
> > On Feb 9, 2017, at 4:08 PM, Matthew Knepley <knepley at gmail.com> wrote:
> >
> > On Thu, Feb 9, 2017 at 4:00 PM, Pierre Jolivet <
> pierre.jolivet at enseeiht.fr> wrote:
> > Oh my, this example also segfaults (at least on my Mac) when running:
> > mpirun -np 1 ./ex1 -pc_type gamg -mat_view
> >
> > Hong, in
> >
> >   MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ() on line 1015
> >
> > you create a matrix without values, not surprising since this is the
> symbolic phase. However,
> > it calls MatAssemblyEnd(), which calls ViewFromOptions(), which gets an
> SEGV because
> > it does not think to check if a->a is NULL.
> >
> > a) Should we make View robust to missing matrix values?
> >
> > b) I don't think this internal matrix should really be viewed like this,
> so maybe change its prefix?
> >
> >   Thanks,
> >
> >     Matt
> >
> > [0]PETSC ERROR: [0] MatView_SeqAIJ_ASCII line 560
> /Volumes/Data/Repositories/petsc/src/mat/impls/aij/seq/aij.c
> > [0]PETSC ERROR: [0] MatView_SeqAIJ line 915 /Volumes/Data/Repositories/
> petsc/src/mat/impls/aij/seq/aij.c
> > [0]PETSC ERROR: [0] MatView line 920 /Volumes/Data/Repositories/
> petsc/src/mat/interface/matrix.c
> > [0]PETSC ERROR: [0] PetscObjectView line 92 /Volumes/Data/Repositories/
> petsc/src/sys/objects/destroy.c
> > [0]PETSC ERROR: [0] PetscObjectViewFromOptions line 2702
> /Volumes/Data/Repositories/petsc/src/sys/objects/options.c
> > [0]PETSC ERROR: [0] MatAssemblyEnd line 5151 /Volumes/Data/Repositories/
> petsc/src/mat/interface/matrix.c
> > [0]PETSC ERROR: [0] MatCreateSeqAIJWithArrays line 4137
> /Volumes/Data/Repositories/petsc/src/mat/impls/aij/seq/aij.c
> > [0]PETSC ERROR: [0] MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ line 1012
> /Volumes/Data/Repositories/petsc/src/mat/impls/aij/seq/matmatmult.c
> > [0]PETSC ERROR: [0] MatTransposeMatMult_SeqAIJ_SeqAIJ line 994
> /Volumes/Data/Repositories/petsc/src/mat/impls/aij/seq/matmatmult.c
> > [0]PETSC ERROR: [0] MatTransposeMatMult line 9607
> /Volumes/Data/Repositories/petsc/src/mat/interface/matrix.c
> > [0]PETSC ERROR: [0] PCGAMGCoarsen_AGG line 914
> /Volumes/Data/Repositories/petsc/src/ksp/pc/impls/gamg/agg.c
> > [0]PETSC ERROR: [0] PCSetUp_GAMG line 425 /Volumes/Data/Repositories/
> petsc/src/ksp/pc/impls/gamg/gamg.c
> > [0]PETSC ERROR: [0] PCSetUp line 886 /Volumes/Data/Repositories/
> petsc/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: [0] KSPSetUp line 291 /Volumes/Data/Repositories/
> petsc/src/ksp/ksp/interface/itfunc.c
> >
> > [0]PETSC ERROR: Petsc Development GIT revision: v3.7.5-3046-gde967d2
> GIT Date: 2017-02-09 18:59:37 +0000
> > [0]PETSC ERROR: Configure options --LDFLAGS=-L/opt/intel/lib
> --download-hypre --download-metis --download-mumps --download-parmetis
> --download-pastix --download-ptscotch --with-blacs-include=/opt/intel/mkl/include
> --with-blacs-lib=/opt/intel/mkl/lib/libmkl_blacs_mpich_lp64.so
> --with-blaslapack-dir=/opt/intel/mkl --with-chtml=0
> --with-mkl_pardiso-include=/opt/intel/mkl/include --with-mkl_pardiso=1
> --with-scalapack-include=/opt/intel/mkl/include
> --with-scalapack-lib="[/opt/intel/mkl/lib/libmkl_
> scalapack_lp64.so,/opt/intel/mkl/lib/libmkl_blacs_mpich_lp64.so]"
> --with-shared-libraries=1 --with-x=1 PETSC_ARCH=arch-darwin-c-debug
> >
> >> On Feb 9, 2017, at 10:52 PM, Pierre Jolivet <pierre.jolivet at enseeiht.fr>
> wrote:
> >>
> >> There, attached.
> >>
> >> Thanks.
> >>
> >> PS: here is the patch also, just in case
> >> diff --git a/src/ksp/ksp/examples/tutorials/ex1.c
> b/src/ksp/ksp/examples/tutorials/ex1.c
> >> index 6382c09..d9dc69a 100644
> >> --- a/src/ksp/ksp/examples/tutorials/ex1.c
> >> +++ b/src/ksp/ksp/examples/tutorials/ex1.c
> >> @@ -123,6 +123,11 @@ int main(int argc,char **args)
> >>      routines.
> >>    */
> >>    ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
> >> +  ierr = KSPSetUp(ksp);CHKERRQ(ierr);
> >> +  char right_options[] = "-ksp_norm_type unpreconditioned -ksp_pc_side
> right";
> >> +  PetscOptionsInsertString(NULL,right_options);
> >> +  ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
> >> +  ierr = KSPSetUp(ksp);CHKERRQ(ierr);
> >>
> >>    if (nonzeroguess) {
> >>      PetscScalar p = .5;
> >>
> >> <ex1.c>
> >>
> >>> On Feb 9, 2017, at 2:56 PM, Matthew Knepley <knepley at gmail.com> wrote:
> >>>
> >>> On Thu, Feb 9, 2017 at 7:48 AM, Pierre Jolivet <
> pierre.jolivet at enseeiht.fr> wrote:
> >>>> On Feb 9, 2017, at 2:31 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
> >>>>
> >>>> On Thu, Feb 9, 2017 at 7:29 AM, Pierre Jolivet <
> pierre.jolivet at enseeiht.fr> wrote:
> >>>>
> >>>>> On Feb 9, 2017, at 2:17 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
> >>>>>
> >>>>> On Thu, Feb 9, 2017 at 7:06 AM, Pierre Jolivet <
> pierre.jolivet at enseeiht.fr> wrote:
> >>>>>
> >>>>>> On Feb 9, 2017, at 1:37 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
> >>>>>>
> >>>>>> On Thu, Feb 9, 2017 at 4:56 AM, Pierre Jolivet <
> pierre.jolivet at enseeiht.fr> wrote:
> >>>>>> Hello,
> >>>>>> I _might_ be doing something wrong when calling
> KSPSetFromOptions/KSPSetUp, but still, am I supposed to get this kind of
> error?
> >>>>>> [0]PETSC ERROR: No support for this operation for this object type
> >>>>>> [0]PETSC ERROR: KSP gmres does not support PRECONDITIONED with RIGHT
> >>>>>>
> >>>>>> If not, I’ll try to send a MWE (it basically depends on where I set
> the option -ksp_pc_side right).
> >>>>>>
> >>>>>> I believe its telling you that it cannot display the preconditioned
> residual with right preconditioning, which is true.
> >>>>>> Are you requesting both?
> >>>>>
> >>>>> No, I’m just using -ksp_view and -ksp_pc_side right. In the
> “standard” case, it obviously work.
> >>>>> However, the error is triggered if I set -ksp_pc_side right after a
> first call to KSPSetOperators()/KSPSetUp() without destroying the KSP.
> >>>>> Should I destroy my KSP before changing from left- (default) to
> right-preconditioning if KSPSetOperators()/KSPSetUp has already been called?
> >>>>>
> >>>>> Its likely you have to use
> >>>>>
> >>>>>   http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/
> KSPSetNormType.html
> >>>>>
> >>>>> or
> >>>>>
> >>>>>   -ksp_norm_type unpreconditioned
> >>>>
> >>>> If I use this after my first call to KSPSetUp, I now end up with:
> >>>> [0]PETSC ERROR: KSP gmres does not support UNPRECONDITIONED with LEFT
> >>>>
> >>>> Yes, you would need to use it at the same time that you change the
> preconditioning side.
> >>>
> >>> Yes, I’m using PetscOptionsInsert, and whether I let PETSc parse
> "-ksp_pc_side right -ksp_norm_type unpreconditioned” or "-ksp_norm_type
> unpreconditioned -ksp_pc_side right” I end up with the same error
> afterwards.
> >>>
> >>> Can you make an example do it (SNES ex5)? Or can you send a small
> example that does it?
> >>>
> >>>   Thanks,
> >>>
> >>>     Matt
> >>>
> >>> Thanks.
> >>>
> >>>>   Thanks,
> >>>>
> >>>>      Matt
> >>>>> since we would not have a chance to change it since its been set
> already.
> >>>>>
> >>>>>   Thanks,
> >>>>>
> >>>>>      Matt
> >>>>>
> >>>>> Thanks.
> >>>>>
> >>>>>>   Thanks,
> >>>>>>
> >>>>>>      Matt
> >>>>>>
> >>>>>> Thanks,
> >>>>>> Pierre
> >>>>>>
> >>>>>>
> >>>>>>
> >>>>>> --
> >>>>>> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> >>>>>> -- Norbert Wiener
> >>>>>
> >>>>>
> >>>>>
> >>>>>
> >>>>> --
> >>>>> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> >>>>> -- Norbert Wiener
> >>>>
> >>>>
> >>>>
> >>>>
> >>>> --
> >>>> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> >>>> -- Norbert Wiener
> >>>
> >>>
> >>>
> >>>
> >>> --
> >>> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> >>> -- Norbert Wiener
> >>
> >
> >
> >
> >
> > --
> > What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> > -- Norbert Wiener
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20170209/59b3c376/attachment.html>


More information about the petsc-dev mailing list