#include /*I "petscsys.h" I*/ #include typedef enum {BUILDTWOSIDED_NOTSET = -1, #if defined(PETSC_HAVE_MPI_IBARRIER) BUILDTWOSIDED_IBARRIER, #endif BUILDTWOSIDED_ALLREDUCE} BuildTwoSidedType; static const char *const BuildTwoSidedTypes[] = { #if defined(PETSC_HAVE_MPI_IBARRIER) "IBARRIER", #endif "ALLREDUCE", "BuildTwoSidedType", "BUILDTWOSIDED_", 0 }; static BuildTwoSidedType _twosided_type = BUILDTWOSIDED_NOTSET; #undef __FUNCT__ #define __FUNCT__ "PetscBuildTwoSidedGetType" static PetscErrorCode PetscBuildTwoSidedGetType(BuildTwoSidedType *twosided) { PetscErrorCode ierr; PetscFunctionBegin; if (_twosided_type == BUILDTWOSIDED_NOTSET) { #if defined(PETSC_HAVE_MPI_IBARRIER) _twosided_type = BUILDTWOSIDED_IBARRIER; #else _twosided_type = BUILDTWOSIDED_ALLREDUCE; #endif ierr = PetscOptionsGetEnum(PETSC_NULL,"-build_twosided",BuildTwoSidedTypes,(PetscEnum*)&_twosided_type,PETSC_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 = PETSC_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 = PETSC_NULL; *(char**)contiguous = contig; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "BuildTwoSided_Ibarrier" static PetscErrorCode BuildTwoSided_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