VecDestroy and memory leak

nicolas aunai nicolas.aunai at gmail.com
Sat Nov 28 08:44:19 CST 2009


Hi,

ah ok I thought the natural vector was destroyed with VecDestroy each
time it is called. I could indeed check what you said in the memory
log I'm writing.

Here is what -malloc_dump is writing at the end of the execution :

[0]Total space allocated 26648 bytes
[ 0]8 bytes PetscStrallocpy() line 82 in src/sys/utils/str.c
      [0]  PetscObjectChangeTypeName() line 107 in src/sys/objects/pname.c
      [0]  VecCreate_Seq_Private() line 714 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecCreate_Seq() line 804 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecSetType() line 39 in src/vec/vec/interface/vecreg.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]8 bytes PetscMapSetUp() line 140 in src/vec/vec/impls/mpi/pmap.c
      [0]  VecCreate_Seq_Private() line 714 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecCreate_Seq() line 804 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecSetType() line 39 in src/vec/vec/interface/vecreg.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]24 bytes VecCreate_Seq_Private() line 715 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecCreate_Seq() line 804 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecSetType() line 39 in src/vec/vec/interface/vecreg.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]25344 bytes VecCreate_Seq() line 809 in src/vec/vec/impls/seq/bvec2.c
      [0]  VecSetType() line 39 in src/vec/vec/interface/vecreg.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]40 bytes VecCreate() line 42 in src/vec/vec/interface/veccreate.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]496 bytes VecCreate() line 39 in src/vec/vec/interface/veccreate.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]64 bytes VecCreate() line 39 in src/vec/vec/interface/veccreate.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]656 bytes VecCreate() line 39 in src/vec/vec/interface/veccreate.c
      [0]  VecCreateSeq() line 37 in src/vec/vec/impls/seq/vseqcr.c
[ 0]8 bytes PetscCommDuplicate() line 221 in src/sys/objects/tagm.c
      [0]  PetscHeaderCreate_Private() line 26 in src/sys/objects/inherit.c
      [0]  VecCreate() line 32 in src/vec/vec/interface/veccreate.c
      [0]  VecCreateSeqWithArray() line 771 in src/vec/vec/impls/seq/bvec2.c
      [0]  DACreate2d() line 338 in src/dm/da/src/da2.c


All this disapears when I comment the call to the function I suspect
to be responsible for my leak.
However I still can't see where I don't free what I've created.


here is the function :


void getseqsol(void)
{
DA da;
Vec b[3], nsol, bx, by, bz, ssol;
PetscInt mx, my, i, j, ij, dof, ija;
VecScatter ctx;
PetscScalar *bxa, *bya, *bza;
Vec sol;
struct UserCtx   *uc = (struct UserCtx *) (*(solv.dmmg))->user;

sol = DMMGGetx(solv.dmmg);

PetscObjectQuery((PetscObject) sol, "DA", (PetscObject *) &da);


DAGetInfo( da, PETSC_IGNORE, &mx, &my,
           PETSC_IGNORE,PETSC_IGNORE,
           PETSC_IGNORE, PETSC_IGNORE,
           &dof, PETSC_IGNORE, PETSC_IGNORE ,
           PETSC_IGNORE);


DACreateNaturalVector(da, &nsol);
DAGlobalToNaturalBegin(da, sol, INSERT_VALUES, nsol);
DAGlobalToNaturalEnd(da, sol, INSERT_VALUES, nsol);


VecCreateSeq(PETSC_COMM_SELF, mx*my*dof, &ssol);
VecScatterCreateToAll(nsol, &ctx, &ssol);
VecScatterBegin(ctx, nsol, ssol, INSERT_VALUES, SCATTER_FORWARD);
VecScatterEnd(ctx, nsol, ssol, INSERT_VALUES,SCATTER_FORWARD);

VecDestroy(nsol);



VecCreateSeq(PETSC_COMM_SELF, mx*my, &bx);
VecCreateSeq(PETSC_COMM_SELF, mx*my, &by);
VecCreateSeq(PETSC_COMM_SELF, mx*my, &bz);


b[0] = bx;
b[1] = by;
b[2] = bz;

VecSetBlockSize(ssol, 3);
VecStrideGatherAll(ssol, b, INSERT_VALUES);

VecGetArray(bx, &bxa);
VecGetArray(by, &bya);
VecGetArray(bz, &bza);



if (uc->ipc == 0)
   {
   for (i=0; i < uc->nx; i++)
       {
       for (j=0; j < uc->ny+1; j++)
           {
           ij = i + j*(uc->nx+1);
           ija = i + j*uc->nx;

           uc->s1[ij].c[0] = bxa[ija];
           uc->s1[ij].c[1] = bya[ija];
           uc->s1[ij].c[2] = bza[ija];
           }
       }

   for (j = 0; j < uc->ny+1; j++)
       {
       ij  = uc->nx + j*(uc->nx+1);
       ija = 0 + j*(uc->nx+1);

       uc->s1[ij].c[0] = uc->s1[ija].c[0];
       uc->s1[ij].c[1] = uc->s1[ija].c[1];
       uc->s1[ij].c[2] = uc->s1[ija].c[2];
       }
   }


else if (uc->ipc == 1)
   {
   for (i=0; i < uc->nx; i++)
       {
       for (j=0; j < uc->ny+1; j++)
           {
           ij = i + j*(uc->nx+1);
           ija = i + j*uc->nx;

           uc->s1[ij].b[0] = bxa[ija];
           uc->s1[ij].b[1] = bya[ija];
           uc->s1[ij].b[2] = bza[ija];
           }
       }

   for (j = 0; j < uc->ny+1; j++)
       {
       ij  = uc->nx + j*(uc->nx+1);
       ija = 0 + j*(uc->nx+1);

       uc->s1[ij].b[0] = uc->s1[ija].b[0];
       uc->s1[ij].b[1] = uc->s1[ija].b[1];
       uc->s1[ij].b[2] = uc->s1[ija].b[2];
       }
   }


VecRestoreArray(bx, &bxa);
VecRestoreArray(by, &bya);
VecRestoreArray(bz, &bza);

VecDestroy(bx);
VecDestroy(by);
VecDestroy(bz);
VecScatterDestroy(ctx);
VecDestroy(ssol);

}



 The malloc_dump seems to refer to sequential vectors... well I think
I'm freeing everything I've created unless I misunderstood something ?


Nico



2009/11/28 Jed Brown <jed at 59a2.org>:
> On Sat, 28 Nov 2009 10:18:04 +0100, nicolas aunai <nicolas.aunai at gmail.com> wrote:
>> I've looked at the function VecDestroy() definition, it seems to check
>> a PetscObject member called 'refct' before actually free the memory,
>> I've printed this 'refct' for my natural vector and its value is '2',
>> while it is '1' for a vector correctly freed.
>
> The DA holds a reference to this vector so it can give it away the next
> time you call DACreateNaturalVector().  It's not a leak because the DA
> will destroy it's reference in DADestroy().  You can run with
> -malloc_dump to confirm that all PETSc objects have been destroyed.
>
> Jed
>


More information about the petsc-users mailing list