xref: /petsc/src/sys/utils/mpitr.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
1*7d0a6c19SBarry Smith 
20f8e0872SSatish Balay /*
30f8e0872SSatish Balay     Code for tracing mistakes in MPI usage. For example, sends that are never received,
40f8e0872SSatish Balay   nonblocking messages that are not correctly waited for, etc.
50f8e0872SSatish Balay */
60f8e0872SSatish Balay 
7d382aafbSBarry Smith #include "petscsys.h"           /*I "petscsys.h" I*/
80f8e0872SSatish Balay 
9c8217ed5SSatish Balay #if defined(PETSC_USE_LOG) && !defined(__MPIUNI_H)
100f8e0872SSatish Balay 
110f8e0872SSatish Balay #undef __FUNCT__
120f8e0872SSatish Balay #define __FUNCT__ "PetscMPIDump"
130f8e0872SSatish Balay /*@C
140f8e0872SSatish Balay    PetscMPIDump - Dumps a listing of incomplete MPI operations, such as sends that
150f8e0872SSatish Balay    have never been received, etc.
160f8e0872SSatish Balay 
170f8e0872SSatish Balay    Collective on PETSC_COMM_WORLD
180f8e0872SSatish Balay 
190f8e0872SSatish Balay    Input Parameter:
200f8e0872SSatish Balay .  fp - file pointer.  If fp is NULL, stdout is assumed.
210f8e0872SSatish Balay 
220f8e0872SSatish Balay    Options Database Key:
230f8e0872SSatish Balay .  -mpidump - Dumps MPI incompleteness during call to PetscFinalize()
240f8e0872SSatish Balay 
250f8e0872SSatish Balay     Level: developer
260f8e0872SSatish Balay 
270f8e0872SSatish Balay .seealso:  PetscMallocDump()
280f8e0872SSatish Balay  @*/
297087cfbeSBarry Smith PetscErrorCode  PetscMPIDump(FILE *fd)
300f8e0872SSatish Balay {
310f8e0872SSatish Balay   PetscErrorCode ierr;
320f8e0872SSatish Balay   PetscMPIInt    rank;
330f8e0872SSatish Balay   double         tsends,trecvs,work;
34f56c2debSBarry Smith   int            err;
350f8e0872SSatish Balay 
360f8e0872SSatish Balay   PetscFunctionBegin;
370f8e0872SSatish Balay   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
38da9f1d6bSBarry Smith   if (!fd) fd = PETSC_STDOUT;
390f8e0872SSatish Balay 
400f8e0872SSatish Balay   /* Did we wait on all the non-blocking sends and receives? */
410f8e0872SSatish Balay   ierr = PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);CHKERRQ(ierr);
420f8e0872SSatish Balay   if (irecv_ct + isend_ct != sum_of_waits_ct) {
430f8e0872SSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,fd,"[%d]You have not waited on all non-blocking sends and receives",rank);CHKERRQ(ierr);
440f8e0872SSatish Balay     ierr = 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);CHKERRQ(ierr);
45f56c2debSBarry Smith     err = fflush(fd);
46e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
470f8e0872SSatish Balay   }
480f8e0872SSatish Balay   ierr = PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);CHKERRQ(ierr);
490f8e0872SSatish Balay   /* Did we receive all the messages that we sent? */
500f8e0872SSatish Balay   work = irecv_ct + recv_ct;
510f8e0872SSatish Balay   ierr = MPI_Reduce(&work,&trecvs,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
520f8e0872SSatish Balay   work = isend_ct + send_ct;
530f8e0872SSatish Balay   ierr = MPI_Reduce(&work,&tsends,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
540f8e0872SSatish Balay   if (!rank && tsends != trecvs) {
550f8e0872SSatish Balay     ierr = PetscFPrintf(PETSC_COMM_SELF,fd,"Total number sends %g not equal receives %g\n",tsends,trecvs);CHKERRQ(ierr);
56f56c2debSBarry Smith     err = fflush(fd);
57e32f2f54SBarry Smith     if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
580f8e0872SSatish Balay   }
590f8e0872SSatish Balay   PetscFunctionReturn(0);
600f8e0872SSatish Balay }
610f8e0872SSatish Balay 
620f8e0872SSatish Balay #else
630f8e0872SSatish Balay 
640f8e0872SSatish Balay #undef __FUNCT__
650f8e0872SSatish Balay #define __FUNCT__ "PetscMPIDump"
667087cfbeSBarry Smith PetscErrorCode  PetscMPIDump(FILE *fd)
670f8e0872SSatish Balay {
680f8e0872SSatish Balay   PetscFunctionBegin;
690f8e0872SSatish Balay   PetscFunctionReturn(0);
700f8e0872SSatish Balay }
710f8e0872SSatish Balay 
720f8e0872SSatish Balay #endif
730f8e0872SSatish Balay 
740f8e0872SSatish Balay 
750f8e0872SSatish Balay 
760f8e0872SSatish Balay 
770f8e0872SSatish Balay 
780f8e0872SSatish Balay 
790f8e0872SSatish Balay 
800f8e0872SSatish Balay 
810f8e0872SSatish Balay 
82