[petsc-users] Using -malloc_dump to examine memory leak

Yuyun Yang yyang85 at stanford.edu
Thu Apr 18 10:37:58 CDT 2019


Oh thank you, sorry I probably missed that!

-----Original Message-----
From: Smith, Barry F. <bsmith at mcs.anl.gov> 
Sent: Thursday, April 18, 2019 8:36 AM
To: Yuyun Yang <yyang85 at stanford.edu>
Cc: petsc-users at mcs.anl.gov
Subject: Re: [petsc-users] Using -malloc_dump to examine memory leak


   Please see my mail from yesterday; perhaps it did not get through

 On Apr 16, 2019, at 11:56 PM, Yuyun Yang <yyang85 at stanford.edu> wrote:
> 
> 
> So using -objects_dump showed nothing below the line:
> The following objects were never freed
> -----------------------------------------
> Does that mean I’ve freed all the objects? If so, why would -malloc_dump give me a bunch of messages?

 No idea. If you have messages from -malloc_dump then memory has not been freed and almost always this is due to undestroyed objects. Presumably you didn't turn off logging in ./configure (nobody does that). If you use the option -objects_dump this is processed in PetscOptionsCheckInitial_Private() which is called from PetscInitialize() with #if defined(PETSC_USE_LOG)  ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&PetscObjectsLog);CHKERRQ(ierr);
#endif
then in inherit.c each object creation gets logged with

PetscErrorCode  PetscHeaderCreate_Private(PetscObject h,PetscClassId classid,const char class_name[],const char descr[],const char mansec[],
                                         MPI_Comm comm,PetscObjectDestroyFunction destroy,PetscObjectViewFunction view) {  static PetscInt idcnt = 1;  PetscErrorCode  ierr; #if defined(PETSC_USE_LOG)
 PetscObject     *newPetscObjects;
 PetscInt         newPetscObjectsMaxCounts,i;
#endif

 PetscFunctionBegin;
 h->classid               = classid;
 h->type                  = 0;
 h->class_name            = (char*)class_name;
 h->description           = (char*)descr;
 h->mansec                = (char*)mansec;
 h->prefix                = 0;
 h->refct                 = 1;
#if defined(PETSC_HAVE_SAWS)
 h->amsmem                = PETSC_FALSE;
#endif
 h->id                    = idcnt++;
 h->parentid              = 0;
 h->qlist                 = 0;
 h->olist                 = 0;
 h->bops->destroy         = destroy;
 h->bops->view            = view;
 h->bops->getcomm         = PetscObjectGetComm_Petsc;
 h->bops->compose         = PetscObjectCompose_Petsc;
 h->bops->query           = PetscObjectQuery_Petsc;
 h->bops->composefunction = PetscObjectComposeFunction_Petsc;
 h->bops->queryfunction   = PetscObjectQueryFunction_Petsc;

 ierr = PetscCommDuplicate(comm,&h->comm,&h->tag);CHKERRQ(ierr);

#if defined(PETSC_USE_LOG)
 /* Keep a record of object created */
 if (PetscObjectsLog) {
   PetscObjectsCounts++;
   for (i=0; i<PetscObjectsMaxCounts; i++) {
     if (!PetscObjects[i]) {
       PetscObjects[i] = h;
       PetscFunctionReturn(0);
     }
   }
   /* Need to increase the space for storing PETSc objects */
   if (!PetscObjectsMaxCounts) newPetscObjectsMaxCounts = 100;
   else                        newPetscObjectsMaxCounts = 2*PetscObjectsMaxCounts;
   ierr = PetscMalloc1(newPetscObjectsMaxCounts,&newPetscObjects);CHKERRQ(ierr);
   ierr = PetscMemcpy(newPetscObjects,PetscObjects,PetscObjectsMaxCounts*sizeof(PetscObject));CHKERRQ(ierr);
   ierr = PetscMemzero(newPetscObjects+PetscObjectsMaxCounts,(newPetscObjectsMaxCounts - PetscObjectsMaxCounts)*sizeof(PetscObject));CHKERRQ(ierr);
   ierr = PetscFree(PetscObjects);CHKERRQ(ierr);

   PetscObjects                        = newPetscObjects;
   PetscObjects[PetscObjectsMaxCounts] = h;
   PetscObjectsMaxCounts               = newPetscObjectsMaxCounts;
 }
#endif
 PetscFunctionReturn(0);
}

Then in PetscFinalize() 

#if defined(PETSC_USE_LOG)
 /*
      List all objects the user may have forgot to free  */  if (PetscObjectsLog) {
   ierr = PetscOptionsHasName(NULL,NULL,"-objects_dump",&flg1);CHKERRQ(ierr);
   if (flg1) {
     MPI_Comm local_comm;
     char     string[64];

     ierr = PetscOptionsGetString(NULL,NULL,"-objects_dump",string,64,NULL);CHKERRQ(ierr);
     ierr = MPI_Comm_dup(MPI_COMM_WORLD,&local_comm);CHKERRQ(ierr);
     ierr = PetscSequentialPhaseBegin_Private(local_comm,1);CHKERRQ(ierr);
     ierr = PetscObjectsDump(stdout,(string[0] == 'a') ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
     ierr = PetscSequentialPhaseEnd_Private(local_comm,1);CHKERRQ(ierr);
     ierr = MPI_Comm_free(&local_comm);CHKERRQ(ierr);
   }
 }
#endif

Very simple code. 

Try 

-objects_dump all 

and see if that changes things. Otherwise you'll need to use the debugger to track why it is not printing out the undestroyed objects.


> 
> I also ran valgrind, and am getting the following error when using --track-origins=yes:
> Uninitialised value was created by a stack allocation
> ==617==    at 0x417D66: ComputeVel_qd::computeVel(double*, double, int&, int) (fault.cpp:953)
> 
> The first argument of this function is a Petsc array (obtained from VecGetArray, and after calling this function, VecRestoreArray was called on that same object). I can’t think of a way why valgrind thinks there is an uninitialized value. Is this by any chance related to the -malloc_dump messages?

 Should be unrelated to the -malloc_dump message.

 A "stack allocation" is usually a local variable, while the value from VecGetArray() is from the heap. Check the use of all your local variables in this function.

 Barry


> On Apr 18, 2019, at 10:06 AM, Yuyun Yang <yyang85 at stanford.edu> wrote:
> 
> Just wanted to follow up on this question, thanks for your time!
> 
> Best,
> Yuyun
> 
> Get Outlook for iOS
> _____________________________
> From: Yuyun Yang <yyang85 at stanford.edu>
> Sent: Tuesday, April 16, 2019 21:56
> Subject: RE: [petsc-users] Using -malloc_dump to examine memory leak
> To: Yuyun Yang <yyang85 at stanford.edu>, Smith, Barry F. 
> <bsmith at mcs.anl.gov>
> 
> 
> So using -objects_dump showed nothing below the line:
> The following objects were never freed
> -----------------------------------------
> Does that mean I’ve freed all the objects? If so, why would -malloc_dump give me a bunch of messages?
>  
> I also ran valgrind, and am getting the following error when using --track-origins=yes:
> Uninitialised value was created by a stack allocation
> ==617==    at 0x417D66: ComputeVel_qd::computeVel(double*, double, int&, int) (fault.cpp:953)
>  
> The first argument of this function is a Petsc array (obtained from VecGetArray, and after calling this function, VecRestoreArray was called on that same object). I can’t think of a way why valgrind thinks there is an uninitialized value. Is this by any chance related to the -malloc_dump messages?
>  
> Thanks for your help on this.
>  
> Best regards,
> Yuyun
>  
> From: petsc-users <petsc-users-bounces at mcs.anl.gov> On Behalf Of Yuyun 
> Yang via petsc-users
> Sent: Tuesday, April 16, 2019 7:30 AM
> To: Smith, Barry F. <bsmith at mcs.anl.gov>
> Cc: petsc-users at mcs.anl.gov
> Subject: Re: [petsc-users] Using -malloc_dump to examine memory leak
>  
> Great, thank you for the advice!
>  
> Best regards,
> Yuyun
>  
> Get Outlook for iOS
> From: Smith, Barry F. <bsmith at mcs.anl.gov>
> Sent: Tuesday, April 16, 2019 5:54:15 AM
> To: Yuyun Yang
> Cc: petsc-users at mcs.anl.gov
> Subject: Re: [petsc-users] Using -malloc_dump to examine memory leak
>  
> 
>   Please try the flag -options_dump this tries to give a much more 
> concise view of what objects have not been freed. For example I commented out the last VecDestroy() in src/snes/examples/tutorials/ex19.c and then obtained:
> 
> ./ex19 -objects_dump
> lid velocity = 0.0625, prandtl # = 1., grashof # = 1.
> Number of SNES iterations = 2
> The following objects were never freed
> -----------------------------------------
> [0] DM da DM_0x84000000_0
>       [0]  DMDACreate2d() in /Users/barrysmith/Src/petsc/src/dm/impls/da/da2.c
>       [0]  main() in 
> /Users/barrysmith/Src/petsc/src/snes/examples/tutorials/ex19.c
> [0] Vec seq Vec_0x84000000_6
>       [0]  DMCreateGlobalVector() in /Users/barrysmith/Src/petsc/src/dm/interface/dm.c
>       [0]  main() in 
> /Users/barrysmith/Src/petsc/src/snes/examples/tutorials/ex19.c
> 
> Now I just need to look at the calls to DMCreateGlobalVector and DMDACreate2d in main to see what I did not free.
> 
> Note that since PETSc objects may hold references to other PETSc objects some items may not be freed for which you DID call destroy on.
> For example because the unfreed vector holds a reference to the DM the 
> DM is listed as not freed. Once you properly destroy the vector you'll not that the DM is no longer listed as non freed.
> 
> It can be a little overwhelming at first to figure out what objects 
> have not been freed. We recommending setting the environmental 
> variable export PETSC_OPTIONS=-malloc_test so that every run of your code reports memory issues and you can keep them under control from the beginning (when the code is small and growing) rather than wait until the code is large and there are many unfreed objects to chase down.
> 
>    Good luck
> 
> 
> 
>    Barry
> 
> 
> > On Apr 16, 2019, at 1:14 AM, Yuyun Yang via petsc-users <petsc-users at mcs.anl.gov> wrote:
> > 
> > Hello team,
> >  
> > I’m trying to use the options -malloc_dump and -malloc_debug to examine memory leaks. The messages however, are quite generic, and don’t really tell me where the problems occur, for example:
> >  
> > [ 0]1520 bytes VecCreate() line 35 in /home/yyy910805/petsc/src/vec/vec/interface/veccreate.c
> >       [0]  PetscMallocA() line 35 in /home/yyy910805/petsc/src/sys/memory/mal.c
> >       [0]  VecCreate() line 30 in /home/yyy910805/petsc/src/vec/vec/interface/veccreate.c
> >       [0]  VecDuplicate_Seq() line 804 in /home/yyy910805/petsc/src/vec/vec/impls/seq/bvec2.c
> >       [0]  VecDuplicate() line 375 in 
> > /home/yyy910805/petsc/src/vec/vec/interface/vector.c
> >  
> > The code is huge, so going through every single VecCreate/VecDuplicate and VecDestroy is going to be time-consuming. Meanwhile, running valgrind gave me some uninitialized values errors that don’t seem to be related to the above message (or maybe they are?).
> >  
> > How can I use this option to debug effectively?
> >  
> > Thanks a lot,
> > Yuyun
> 



More information about the petsc-users mailing list