#include /*I "petscsys.h" I*/ const char *const PetscBuildTwoSidedTypes[] = { "ALLREDUCE", "IBARRIER", "PetscBuildTwoSidedType", "PETSC_BUILDTWOSIDED_", 0 }; static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET; #undef __FUNCT__ #define __FUNCT__ "PetscCommBuildTwoSidedSetType" /*@ PetscCommBuildTwoSidedSetType - set algorithm to use when building two-sided communication Logically Collective Input Arguments: + comm - PETSC_COMM_WORLD - twosided - algorithm to use in subsequent calls to PetscCommBuildTwoSided() Level: developer Note: This option is currently global, but could be made per-communicator. .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedGetType() @*/ PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm comm,PetscBuildTwoSidedType twosided) { PetscFunctionBegin; #if defined(PETSC_USE_DEBUG) { /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */ PetscMPIInt ierr; PetscMPIInt b1[2],b2[2]; b1[0] = -(PetscMPIInt)twosided; b1[1] = (PetscMPIInt)twosided; ierr = MPI_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); if (-b2[0] != b2[1]) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes"); } #endif _twosided_type = twosided; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscCommBuildTwoSidedGetType" /*@ PetscCommBuildTwoSidedGetType - set algorithm to use when building two-sided communication Logically Collective Output Arguments: + comm - communicator on which to query algorithm - twosided - algorithm to use for PetscCommBuildTwoSided() Level: developer .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType() @*/ PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm comm,PetscBuildTwoSidedType *twosided) { PetscErrorCode ierr; PetscFunctionBegin; *twosided = PETSC_BUILDTWOSIDED_NOTSET; if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) { #if defined(PETSC_HAVE_MPI_IBARRIER) # if defined(PETSC_HAVE_MPICH_CH3_SOCK) && !defined(PETSC_HAVE_MPICH_CH3_SOCK_FIXED_NBC_PROGRESS) /* Deadlock in Ibarrier: http://trac.mpich.org/projects/mpich/ticket/1785 */ _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; # else _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER; # endif #else _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; #endif ierr = PetscOptionsGetEnum(NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL);CHKERRQ(ierr); } *twosided = _twosided_type; PetscFunctionReturn(0); } #if defined(PETSC_HAVE_MPI_IBARRIER) /* Segmented (extendable) array implementation */ typedef struct _SegArray *SegArray; struct _SegArray { PetscInt unitbytes; PetscInt alloc; PetscInt used; PetscInt tailused; SegArray tail; union { /* Dummy types to ensure alignment */ PetscReal dummy_real; PetscInt dummy_int; char array[1]; } u; }; #undef __FUNCT__ #define __FUNCT__ "SegArrayCreate" static PetscErrorCode SegArrayCreate(PetscInt unitbytes,PetscInt expected,SegArray *seg) { PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscMalloc(offsetof(struct _SegArray,u)+expected*unitbytes,seg);CHKERRQ(ierr); ierr = PetscMemzero(*seg,offsetof(struct _SegArray,u));CHKERRQ(ierr); (*seg)->unitbytes = unitbytes; (*seg)->alloc = expected; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SegArrayAlloc_Private" static PetscErrorCode SegArrayAlloc_Private(SegArray *seg,PetscInt count) { PetscErrorCode ierr; SegArray newseg,s; PetscInt alloc; PetscFunctionBegin; s = *seg; /* Grow at least fast enough to hold next item, like Fibonacci otherwise (up to 1MB chunks) */ alloc = PetscMax(s->used+count,PetscMin(1000000/s->unitbytes+1,s->alloc+s->tailused)); ierr = PetscMalloc(offsetof(struct _SegArray,u)+alloc*s->unitbytes,&newseg);CHKERRQ(ierr); ierr = PetscMemzero(newseg,offsetof(struct _SegArray,u));CHKERRQ(ierr); newseg->unitbytes = s->unitbytes; newseg->tailused = s->used + s->tailused; newseg->tail = s; newseg->alloc = alloc; *seg = newseg; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SegArrayGet" static PetscErrorCode SegArrayGet(SegArray *seg,PetscInt count,void *array) { PetscErrorCode ierr; SegArray s; PetscFunctionBegin; s = *seg; if (PetscUnlikely(s->used + count > s->alloc)) {ierr = SegArrayAlloc_Private(seg,count);CHKERRQ(ierr);} s = *seg; *(char**)array = &s->u.array[s->used*s->unitbytes]; s->used += count; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SegArrayDestroy" static PetscErrorCode SegArrayDestroy(SegArray *seg) { PetscErrorCode ierr; SegArray s; PetscFunctionBegin; for (s=*seg; s;) { SegArray tail = s->tail; ierr = PetscFree(s);CHKERRQ(ierr); s = tail; } *seg = NULL; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "SegArrayExtract" /* Extracts contiguous data and resets segarray */ static PetscErrorCode SegArrayExtract(SegArray *seg,void *contiguous) { PetscErrorCode ierr; PetscInt unitbytes; SegArray s,t; char *contig,*ptr; PetscFunctionBegin; s = *seg; unitbytes = s->unitbytes; ierr = PetscMalloc((s->used+s->tailused)*unitbytes,&contig);CHKERRQ(ierr); ptr = contig + s->tailused*unitbytes; ierr = PetscMemcpy(ptr,s->u.array,s->used*unitbytes);CHKERRQ(ierr); for (t=s->tail; t;) { SegArray tail = t->tail; ptr -= t->used*unitbytes; ierr = PetscMemcpy(ptr,t->u.array,t->used*unitbytes);CHKERRQ(ierr); ierr = PetscFree(t);CHKERRQ(ierr); t = tail; } if (ptr != contig) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Tail count does not match"); s->tailused = 0; s->tail = NULL; *(char**)contiguous = contig; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscCommBuildTwoSided_Ibarrier" static PetscErrorCode PetscCommBuildTwoSided_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata) { PetscErrorCode ierr; PetscMPIInt nrecvs,tag,unitbytes,done; PetscInt i; char *tdata; MPI_Request *sendreqs,barrier; SegArray segrank,segdata; PetscFunctionBegin; ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); tdata = (char*)todata; ierr = PetscMalloc(nto*sizeof(MPI_Request),&sendreqs);CHKERRQ(ierr); for (i=0; i