1 #if !defined(__SFBASIC_H) 2 #define __SFBASIC_H 3 4 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 5 6 typedef struct _n_PetscSFLink* PetscSFLink; 7 8 #define SFBASICHEADER \ 9 PetscMPIInt niranks; /* Number of incoming ranks (ranks accessing my roots) */ \ 10 PetscMPIInt ndiranks; /* Number of incoming ranks (ranks accessing my roots) in distinguished set */ \ 11 PetscMPIInt *iranks; /* Array of ranks that reference my roots */ \ 12 PetscInt itotal; /* Total number of graph edges referencing my roots */ \ 13 PetscInt *ioffset; /* Array of length niranks+1 holding offset in irootloc[] for each rank */ \ 14 PetscInt *irootloc; /* Incoming roots referenced by ranks starting at ioffset[rank] */ \ 15 PetscInt *irootloc_d[2]; /* A copy of irootloc[local/remote] in device memory if needed */ \ 16 PetscInt rootbuflen[2]; /* Length (in unit) of root buffers, in layout of [PETSCSF_LOCAL/REMOTE] */ \ 17 PetscBool rootcontig[2]; /* True means the local/remote segments of indices in irootloc[] are contiguous ... */ \ 18 PetscInt rootstart[2]; /* ... and start from rootstart[0] and rootstart[1] respectively */ \ 19 PetscSFPackOpt rootpackopt[2]; /* Pack optimization plans based on patterns in irootloc[]. NULL for no optimizations */ \ 20 PetscSFPackOpt rootpackopt_d[2];/* Copy of rootpackopt[] on device if needed */ \ 21 PetscBool rootdups[2]; /* Indices of roots in irootloc[local/remote] have dups. Used for data-race test */ \ 22 PetscInt nrootreqs; /* Number of MPI reqests */ \ 23 PetscSFLink avail; /* One or more entries per MPI Datatype, lazily constructed */ \ 24 PetscSFLink inuse /* Buffers being used for transactions that have not yet completed */ 25 26 typedef struct { 27 SFBASICHEADER; 28 } PetscSF_Basic; 29 30 PETSC_STATIC_INLINE PetscErrorCode PetscSFGetRootInfo_Basic(PetscSF sf,PetscInt *nrootranks,PetscInt *ndrootranks,const PetscMPIInt **rootranks,const PetscInt **rootoffset,const PetscInt **rootloc) 31 { 32 PetscSF_Basic *bas = (PetscSF_Basic*)sf->data; 33 34 PetscFunctionBegin; 35 if (nrootranks) *nrootranks = bas->niranks; 36 if (ndrootranks) *ndrootranks = bas->ndiranks; 37 if (rootranks) *rootranks = bas->iranks; 38 if (rootoffset) *rootoffset = bas->ioffset; 39 if (rootloc) *rootloc = bas->irootloc; 40 PetscFunctionReturn(0); 41 } 42 43 PETSC_STATIC_INLINE PetscErrorCode PetscSFGetLeafInfo_Basic(PetscSF sf,PetscInt *nleafranks,PetscInt *ndleafranks,const PetscMPIInt **leafranks,const PetscInt **leafoffset,const PetscInt **leafloc,const PetscInt **leafrremote) 44 { 45 PetscFunctionBegin; 46 if (nleafranks) *nleafranks = sf->nranks; 47 if (ndleafranks) *ndleafranks = sf->ndranks; 48 if (leafranks) *leafranks = sf->ranks; 49 if (leafoffset) *leafoffset = sf->roffset; 50 if (leafloc) *leafloc = sf->rmine; 51 if (leafrremote) *leafrremote = sf->rremote; 52 PetscFunctionReturn(0); 53 } 54 55 PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF); 56 PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF,PetscViewer); 57 PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF); 58 PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF); 59 PETSC_INTERN PetscErrorCode PetscSFBcastAndOpEnd_Basic (PetscSF,MPI_Datatype,const void*,void*,MPI_Op); 60 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Basic (PetscSF,MPI_Datatype,const void*,void*,MPI_Op); 61 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF,MPI_Datatype,PetscMemType,void*,PetscMemType,const void*,void*,MPI_Op); 62 PETSC_INTERN PetscErrorCode PetscSFCreateEmbeddedSF_Basic(PetscSF,PetscInt,const PetscInt*,PetscSF*); 63 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF,PetscInt*,const PetscMPIInt**,const PetscInt**,const PetscInt**); 64 #endif 65