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 #if defined(PETSC_HAVE_NVSHMEM) 29 PetscInt rootbuflen_rmax; /* max rootbuflen[REMOTE] over comm */ 30 PetscInt nRemoteLeafRanks; /* niranks - ndiranks */ 31 PetscInt nRemoteLeafRanksMax; /* max nRemoteLeafRanks over comm */ 32 33 PetscInt *leafbufdisp; /* [nRemoteLeafRanks]. For my i-th remote leaf rank, I will put to its leafbuf_shmem[] at offset leafbufdisp[i], in <unit> to be set */ 34 PetscInt *leafsigdisp; /* [nRemoteLeafRanks]. For my i-th remote leaf rank, I am its leafsigdisp[i]-th root rank */ 35 36 PetscInt *leafbufdisp_d; 37 PetscInt *leafsigdisp_d; /* Copy of leafsigdisp[] on device */ 38 PetscMPIInt *iranks_d; /* Copy of the remote part of (leaf) iranks[] on device */ 39 PetscInt *ioffset_d; /* Copy of the remote part of ioffset[] on device */ 40 #endif 41 } PetscSF_Basic; 42 43 static inline PetscErrorCode PetscSFGetRootInfo_Basic(PetscSF sf, PetscInt *nrootranks, PetscInt *ndrootranks, const PetscMPIInt **rootranks, const PetscInt **rootoffset, const PetscInt **rootloc) { 44 PetscSF_Basic *bas = (PetscSF_Basic *)sf->data; 45 46 PetscFunctionBegin; 47 if (nrootranks) *nrootranks = bas->niranks; 48 if (ndrootranks) *ndrootranks = bas->ndiranks; 49 if (rootranks) *rootranks = bas->iranks; 50 if (rootoffset) *rootoffset = bas->ioffset; 51 if (rootloc) *rootloc = bas->irootloc; 52 PetscFunctionReturn(0); 53 } 54 55 static inline PetscErrorCode PetscSFGetLeafInfo_Basic(PetscSF sf, PetscInt *nleafranks, PetscInt *ndleafranks, const PetscMPIInt **leafranks, const PetscInt **leafoffset, const PetscInt **leafloc, const PetscInt **leafrremote) { 56 PetscFunctionBegin; 57 if (nleafranks) *nleafranks = sf->nranks; 58 if (ndleafranks) *ndleafranks = sf->ndranks; 59 if (leafranks) *leafranks = sf->ranks; 60 if (leafoffset) *leafoffset = sf->roffset; 61 if (leafloc) *leafloc = sf->rmine; 62 if (leafrremote) *leafrremote = sf->rremote; 63 PetscFunctionReturn(0); 64 } 65 66 PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF); 67 PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF, PetscViewer); 68 PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF); 69 PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF); 70 PETSC_INTERN PetscErrorCode PetscSFBcastEnd_Basic(PetscSF, MPI_Datatype, const void *, void *, MPI_Op); 71 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Basic(PetscSF, MPI_Datatype, const void *, void *, MPI_Op); 72 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF, MPI_Datatype, PetscMemType, void *, PetscMemType, const void *, void *, MPI_Op); 73 PETSC_INTERN PetscErrorCode PetscSFCreateEmbeddedRootSF_Basic(PetscSF, PetscInt, const PetscInt *, PetscSF *); 74 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF, PetscInt *, const PetscMPIInt **, const PetscInt **, const PetscInt **); 75 76 #if defined(PETSC_HAVE_NVSHMEM) 77 PETSC_INTERN PetscErrorCode PetscSFReset_Basic_NVSHMEM(PetscSF); 78 #endif 79 80 #endif 81