Actual source code: pbarrier.c
petsc-3.7.3 2016-08-01
2: #include <petsc/private/petscimpl.h> /*I "petscsys.h" I*/
4: /* Logging support */
5: PetscLogEvent PETSC_Barrier=0;
7: static int hash(const char *str)
8: {
9: int c,hash = 5381;
11: while ((c = *str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
12: return hash;
13: }
15: PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm comm,PetscMPIInt ctn,int line,const char *func,const char *file)
16: {
17: PetscMPIInt err;
18: PetscMPIInt b1[6],b2[6];
20: b1[0] = -(PetscMPIInt)line; b1[1] = -b1[0];
21: b1[2] = -(PetscMPIInt)hash(func); b1[3] = -b1[2];
22: b1[4] = -(PetscMPIInt)ctn; b1[5] = -b1[4];
23: err = MPI_Allreduce(b1,b2,6,MPI_INT,MPI_MAX,comm);
24: if (err) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"MPI_Allreduced() failed");
25: if (-b2[0] != b2[1]) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"MPI_Allreduce() called in different locations (code lines) on different processors");
26: if (-b2[2] != b2[3]) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"MPI_Allreduce() called in different locations (functions) on different processors");
27: if (-b2[4] != b2[5]) return PetscError(PETSC_COMM_SELF,line,func,file,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL,"MPI_Allreduce() called with different counts %d on different processors",ctn);
28: return 0;
29: }
33: /*@C
34: PetscBarrier - Blocks until this routine is executed by all
35: processors owning the object A.
37: Input Parameters:
38: . A - PETSc object (Mat, Vec, IS, SNES etc...)
39: Must be caste with a (PetscObject), can use NULL (for MPI_COMM_WORLD)
41: Level: intermediate
43: Notes:
44: This routine calls MPI_Barrier with the communicator of the PETSc Object "A".
46: With fortran Use NULL_OBJECT (instead of NULL)
48: Concepts: barrier
50: @*/
51: PetscErrorCode PetscBarrier(PetscObject obj)
52: {
54: MPI_Comm comm;
58: PetscLogEventBegin(PETSC_Barrier,obj,0,0,0);
59: if (obj) {
60: PetscObjectGetComm(obj,&comm);
61: } else comm = PETSC_COMM_WORLD;
62: MPI_Barrier(comm);
63: PetscLogEventEnd(PETSC_Barrier,obj,0,0,0);
64: return(0);
65: }