Actual source code: mpitr.c

  1: #define PETSC_DLL
  2: /*
  3:     Code for tracing mistakes in MPI usage. For example, sends that are never received,
  4:   nonblocking messages that are not correctly waited for, etc.
  5: */

 7:  #include petsc.h

  9: #if defined(PETSC_USE_LOG) && !defined(_petsc_mpi_uni)

 13: /*@C
 14:    PetscMPIDump - Dumps a listing of incomplete MPI operations, such as sends that
 15:    have never been received, etc.

 17:    Collective on PETSC_COMM_WORLD

 19:    Input Parameter:
 20: .  fp - file pointer.  If fp is NULL, stdout is assumed.

 22:    Options Database Key:
 23: .  -mpidump - Dumps MPI incompleteness during call to PetscFinalize()

 25:     Level: developer

 27: .seealso:  PetscMallocDump()
 28:  @*/
 29: PetscErrorCode  PetscMPIDump(FILE *fd)
 30: {
 32:   PetscMPIInt    rank;
 33:   double         tsends,trecvs,work;

 36:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 37:   if (!fd) fd = stdout;
 38: 
 39:   /* Did we wait on all the non-blocking sends and receives? */
 40:   PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);
 41:   if (irecv_ct + isend_ct != sum_of_waits_ct) {
 42:     PetscFPrintf(PETSC_COMM_SELF,fd,"[%d]You have not waited on all non-blocking sends and receives",rank);
 43:     PetscFPrintf(PETSC_COMM_SELF,fd,"[%d]Number non-blocking sends %g receives %g number of waits %g\n",rank,isend_ct,irecv_ct,sum_of_waits_ct);
 44:     fflush(fd);
 45:   }
 46:   PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);
 47:   /* Did we receive all the messages that we sent? */
 48:   work = irecv_ct + recv_ct;
 49:   MPI_Reduce(&work,&trecvs,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
 50:   work = isend_ct + send_ct;
 51:   MPI_Reduce(&work,&tsends,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);
 52:   if (!rank && tsends != trecvs) {
 53:     PetscFPrintf(PETSC_COMM_SELF,fd,"Total number sends %g not equal receives %g\n",tsends,trecvs);
 54:     fflush(fd);
 55:   }
 56:   return(0);
 57: }

 59: #else

 63: PetscErrorCode  PetscMPIDump(FILE *fd)
 64: {
 66:   return(0);
 67: }

 69: #endif