Memory corruption error in SNESSolve

Dave May dave.mayhem23 at gmail.com
Tue Jan 22 18:55:08 CST 2008


Hey Vijay,
    I can spot something wrong. VecGetSize gives you the global length of
the vector. However VecGetArray only gives you access to the LOCAL portion
of the vector data for that processor. Thus you need to use a local index if
you are going to set the values directly into the array.

In parallel, xx[n-1] will always be outside the bounds of the processors
local portion of the vector if n is the global vector length.

Cheers,
    Dave.


Here's the function.
>
> PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
> {
>  PetscScalar    *xx,*ff,d;
>  PetscErrorCode ierr;
>  PetscInt       i,n;
>
>  PetscFunctionBegin;
>  ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
>  ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
>  ierr = VecGetSize(x,&n);CHKERRQ(ierr);
>
>  double h = 1.0/(n-1) ;
>  d = (PetscReal)(n - 1); d = d*d;
>  double pos ;
>
>  ff[0]   = xx[0];      // left boundary with Dirichlet = 0
>  // residual : (un-1 -2*un + un+1)/h^2 + Reaction*u - f
>  for (i=1; i<n-1; i++) {
>    pos = i*h ;
>    ff[i] = -Diffusion(pos, 0, 0, 0, xx[i])*d*(xx[i-1] - 2.0*xx[i] +
> xx[i+1]) + Reaction(pos, 0, 0, 0, xx[i])*xx[i] - RHS(pos, 0, 0, 0,
> xx[i]) ;
>  }
>  ff[n-1] = xx[n-1] - 1.0;            // right boundary with Dirichlet = 1
>
>  ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
>  ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
>  PetscFunctionReturn(0);
> }
>
> One question though. If i use GetArray and RestoreArray to change the
> values directly in memory, do i still have to make the
> VecAssemblyBegin-End calls ?
>
> Thanks again for the help Matt.
>
> Vijay
>
> On Jan 22, 2008 5:02 PM, Matthew Knepley <knepley at gmail.com> wrote:
> > On Jan 22, 2008 3:54 PM, Vijay M <vijay.m at gmail.com> wrote:
> > > Hi all,
> > >
> > > I am currently trying to use LibMesh with PETSc SNES solver to solve a
> > > nonlinear Diffusion-reaction problem. The code works fine with 1
> > > processor but when i use 2 processors for the exact same problem, it
> > > throws me this error regarding memory corruption. I am not sure what
> > > this error is about and was wondering if someone can help me out on
> > > this. Note: My residual function is being called successfully and
> > > returns the residual back without any problem since i can see some
> > > print message at the end of the call ! Hence, this problem is not in
> > > my custom FormFunction.
> >
> > Actually, it is very likely that the problem is in your FormFunction.
> The
> > error is being reported because after every call to FormFunction, we
> > walk the PETSc heap and validate all the pieces of memory we allocated.
> > Something wrote over a piece of memory we already own, so the OS
> > did not complain and you saw your printout, however there is still
> corruption.
> > I would advise you to run 'valgrind' on the code. It is the best tool I
> have
> > used for finding these kind of errors.
> >
> >   Matt
> >
> >
> > > [1]PETSC ERROR: PetscMallocValidate: error detected at
> > > SNESComputeFunction() line 950 in src/snes/interface/snes.c
> > > [1]PETSC ERROR: Memory at address 0xa5cf40 is corrupted
> > > [1]PETSC ERROR: Probably write past beginning or end of array
> > > [1]PETSC ERROR: Last intact block allocated in PetscStrallocpy() line
> > > 82 in src/sys/utils/str.c
> > > [0]PETSC ERROR: PetscMallocValidate: error detected at
> > > SNESComputeFunction() line 950 in src/snes/interface/snes.c
> > > [0]PETSC ERROR: Memory at address 0xa5cf30 is corrupted
> > > [0]PETSC ERROR: Probably write past beginning or end of array
> > > [0]PETSC ERROR: Last intact block allocated in PetscStrallocpy() line
> > > 82 in src/sys/utils/str.c
> > > [0]PETSC ERROR: --------------------- Error Message
> > > ------------------------------------
> > > [0]PETSC ERROR: Memory corruption!
> > > [0]PETSC ERROR:  !
> > > [0]PETSC ERROR:
> > >
> ------------------------------------------------------------------------
> > > [0]PETSC ERROR: Petsc Release Version 2.3.3, Patch 7, Fri Oct 26
> > > 14:21:35 CDT 2007 HG revision:
> > > 2e223033ba960114833e1f9713ab393ec78c056f
> > > [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> > > [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> > > [0]PETSC ERROR: See docs/index.html for manual pages.
> > > [0]PETSC ERROR:
> > >
> ------------------------------------------------------------------------
> > > [0]PETSC ERROR: ./myProg-opt on a linux-gnu named [1]PETSC ERROR:
> > > --------------------- Error Message
> > > ------------------------------------
> > > [1]PETSC ERROR: Memory corruption!
> > > [1]PETSC ERROR:  !
> > > [1]PETSC ERROR:
> > >
> ------------------------------------------------------------------------
> > > [1]PETSC ERROR: Petsc Release Version 2.3.3, Patch 7, Fri Oct 26
> > > 14:21:35 CDT 2007 HG revision:
> > > 2e223033ba960114833e1f9713ab393ec78c056f
> > > [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> > > [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> > > [1]PETSC ERROR: See docs/index.html for manual pages.
> > > [1]PETSC ERROR: ----------------------------[cli_0]: aborting job:
> > > application called MPI_Abort(comm=0x84000000, 78) - process 0
> > > [cli_1]: aborting job:
> > > application called MPI_Abort(comm=0x84000000, 78) - process 1
> > > grove.ne.tamu.edu by vijaysm Tue Jan 22 15:47:58 2008
> > > [0]PETSC ERROR: Libraries linked from
> > > /state/partition1/local/petsc-2.3.3-p7/lib/linux-gnu-c-debug
> > > [0]PETSC ERROR: Configure run at Fri Nov  9 15:02:46 2007
> > > [0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=ifort
> > > --with-cxx=g++ --download-f-blas-lapack=1 --download-mpich=1
> > > --with-shared=1
> > > [0]PETSC ERROR:
> > >
> ------------------------------------------------------------------------
> > > [0]PETSC ERROR: PetscMallocValidate() line 139 in src/sys/memory/mtr.c
> > > [0]PETSC ERROR: SNESComputeFunction() line 950 in
> src/snes/interface/snes.c
> > > [0]PETSC ERROR: SNESSolve_LS() line 158 in src/snes/impls/ls/ls.c
> > > [0]PETSC ERROR: SNESSolve() line 1871 in src/snes/interface/snes.c
> > > [0]PETSC ERROR: main() line 395 in unknowndirectory/myProg.C
> > >
> > > If you need some more information, feel free to let me know. Thanks in
> advance !
> > >
> > > Vijay
> > >
> > >
> >
> >
> >
> > --
> > 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-users/attachments/20080123/558c8b00/attachment.htm>


More information about the petsc-users mailing list