xref: /petsc/src/sys/utils/mpitr.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
17d0a6c19SBarry 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 
7c6db04a5SJed Brown #include <petscsys.h>           /*I "petscsys.h" I*/
80f8e0872SSatish Balay 
9994fe344SLisandro Dalcin #if defined(PETSC_USE_LOG) && !defined(PETSC_HAVE_MPIUNI)
100f8e0872SSatish Balay 
110f8e0872SSatish Balay /*@C
120f8e0872SSatish Balay    PetscMPIDump - Dumps a listing of incomplete MPI operations, such as sends that
130f8e0872SSatish Balay    have never been received, etc.
140f8e0872SSatish Balay 
150f8e0872SSatish Balay    Collective on PETSC_COMM_WORLD
160f8e0872SSatish Balay 
170f8e0872SSatish Balay    Input Parameter:
180f8e0872SSatish Balay .  fp - file pointer.  If fp is NULL, stdout is assumed.
190f8e0872SSatish Balay 
200f8e0872SSatish Balay    Options Database Key:
210f8e0872SSatish Balay .  -mpidump - Dumps MPI incompleteness during call to PetscFinalize()
220f8e0872SSatish Balay 
230f8e0872SSatish Balay     Level: developer
240f8e0872SSatish Balay 
250f8e0872SSatish Balay .seealso:  PetscMallocDump()
260f8e0872SSatish Balay  @*/
277087cfbeSBarry Smith PetscErrorCode  PetscMPIDump(FILE *fd)
280f8e0872SSatish Balay {
290f8e0872SSatish Balay   PetscMPIInt    rank;
300f8e0872SSatish Balay   double         tsends,trecvs,work;
31f56c2debSBarry Smith   int            err;
320f8e0872SSatish Balay 
330f8e0872SSatish Balay   PetscFunctionBegin;
345f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
35da9f1d6bSBarry Smith   if (!fd) fd = PETSC_STDOUT;
360f8e0872SSatish Balay 
370f8e0872SSatish Balay   /* Did we wait on all the non-blocking sends and receives? */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1));
39ad39c06fSJed Brown   if (petsc_irecv_ct + petsc_isend_ct != petsc_sum_of_waits_ct) {
405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_SELF,fd,"[%d]You have not waited on all non-blocking sends and receives",rank));
415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_SELF,fd,"[%d]Number non-blocking sends %g receives %g number of waits %g\n",rank,petsc_isend_ct,petsc_irecv_ct,petsc_sum_of_waits_ct));
42f56c2debSBarry Smith     err  = fflush(fd);
43*28b400f6SJacob Faibussowitsch     PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
440f8e0872SSatish Balay   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1));
460f8e0872SSatish Balay   /* Did we receive all the messages that we sent? */
47ad39c06fSJed Brown   work = petsc_irecv_ct + petsc_recv_ct;
485f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Reduce(&work,&trecvs,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD));
49ad39c06fSJed Brown   work = petsc_isend_ct + petsc_send_ct;
505f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Reduce(&work,&tsends,1,MPI_DOUBLE,MPI_SUM,0,PETSC_COMM_WORLD));
51dd400576SPatrick Sanan   if (rank == 0 && tsends != trecvs) {
525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFPrintf(PETSC_COMM_SELF,fd,"Total number sends %g not equal receives %g\n",tsends,trecvs));
53f56c2debSBarry Smith     err  = fflush(fd);
54*28b400f6SJacob Faibussowitsch     PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file");
550f8e0872SSatish Balay   }
560f8e0872SSatish Balay   PetscFunctionReturn(0);
570f8e0872SSatish Balay }
580f8e0872SSatish Balay 
590f8e0872SSatish Balay #else
600f8e0872SSatish Balay 
617087cfbeSBarry Smith PetscErrorCode  PetscMPIDump(FILE *fd)
620f8e0872SSatish Balay {
630f8e0872SSatish Balay   PetscFunctionBegin;
640f8e0872SSatish Balay   PetscFunctionReturn(0);
650f8e0872SSatish Balay }
660f8e0872SSatish Balay 
670f8e0872SSatish Balay #endif
680f8e0872SSatish Balay 
69b674149eSJunchao Zhang #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
708198064fSBarry Smith /*
718198064fSBarry Smith     OpenMPI version of MPI_Win_allocate_shared() does not provide __float128 alignment so we provide
728198064fSBarry Smith     a utility that insures alignment up to data item size.
738198064fSBarry Smith */
748198064fSBarry Smith PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint sz,PetscMPIInt szind,MPI_Info info,MPI_Comm comm,void *ptr,MPI_Win *win)
758198064fSBarry Smith {
768198064fSBarry Smith   float          *tmp;
770f8e0872SSatish Balay 
788198064fSBarry Smith   PetscFunctionBegin;
795f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Win_allocate_shared(16+sz,szind,info,comm,&tmp,win));
808198064fSBarry Smith   tmp += ((size_t)tmp) % szind ? szind/4 - ((((size_t)tmp) % szind)/4) : 0;
818198064fSBarry Smith   *(void**)ptr = (void*)tmp;
828198064fSBarry Smith   PetscFunctionReturn(0);
838198064fSBarry Smith }
848198064fSBarry Smith 
858198064fSBarry Smith PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win win,PetscMPIInt rank,MPI_Aint *sz,PetscMPIInt *szind,void *ptr)
868198064fSBarry Smith {
878198064fSBarry Smith   float          *tmp;
888198064fSBarry Smith 
898198064fSBarry Smith   PetscFunctionBegin;
905f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Win_shared_query(win,rank,sz,szind,&tmp));
912c71b3e2SJacob Faibussowitsch   PetscCheckFalse(*szind <= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"szkind %d must be positive",*szind);
928198064fSBarry Smith   tmp += ((size_t)tmp) % *szind ? *szind/4 - ((((size_t)tmp) % *szind)/4) : 0;
938198064fSBarry Smith   *(void**)ptr = (void*)tmp;
948198064fSBarry Smith   PetscFunctionReturn(0);
958198064fSBarry Smith }
968198064fSBarry Smith 
978198064fSBarry Smith #endif
98