Actual source code: sfbasic.h

  1: #pragma once

  3: #include <petsc/private/sfimpl.h>

  5: #define SFBASICHEADER \
  6:   PetscMPIInt    niranks;          /* Number of incoming ranks (ranks accessing my roots) */ \
  7:   PetscMPIInt    ndiranks;         /* Number of incoming ranks (ranks accessing my roots) in distinguished set */ \
  8:   PetscMPIInt   *iranks;           /* Array of ranks that reference my roots */ \
  9:   PetscInt       itotal;           /* Total number of graph edges referencing my roots */ \
 10:   PetscInt      *ioffset;          /* Array of length niranks+1 holding offset in irootloc[] for each rank */ \
 11:   PetscInt      *irootloc;         /* Incoming roots referenced by ranks starting at ioffset[rank] */ \
 12:   PetscInt      *irootloc_d[2];    /* A copy of irootloc[local/remote] in device memory if needed */ \
 13:   PetscInt       rootbuflen[2];    /* Length (in unit) of root buffers, in layout of [PETSCSF_LOCAL/REMOTE] */ \
 14:   PetscBool      rootcontig[2];    /* True means the local/remote segments of indices in irootloc[] are contiguous ... */ \
 15:   PetscInt       rootstart[2];     /* ... and start from rootstart[0] and rootstart[1] respectively */ \
 16:   PetscSFPackOpt rootpackopt[2];   /* Pack optimization plans based on patterns in irootloc[]. NULL for no optimizations */ \
 17:   PetscSFPackOpt rootpackopt_d[2]; /* Copy of rootpackopt[] on device if needed */ \
 18:   PetscBool      rootdups[2];      /* Indices of roots in irootloc[local/remote] have dups. Used for data-race test */ \
 19:   PetscMPIInt    nrootreqs;        /* Number of MPI requests */ \
 20:   PetscSFLink    avail;            /* One or more entries per MPI Datatype, lazily constructed */ \
 21:   PetscSFLink    inuse             /* Buffers being used for transactions that have not yet completed */

 23: typedef struct {
 24:   SFBASICHEADER;
 25: #if defined(PETSC_HAVE_NVSHMEM)
 26:   PetscInt    rootbuflen_rmax;     /* max rootbuflen[REMOTE] over comm */
 27:   PetscMPIInt nRemoteLeafRanks;    /* niranks - ndiranks */
 28:   PetscMPIInt nRemoteLeafRanksMax; /* max nRemoteLeafRanks over comm */

 30:   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 */
 31:   PetscInt *leafsigdisp; /* [nRemoteLeafRanks]. For my i-th remote leaf rank, I am its leafsigdisp[i]-th root rank */

 33:   PetscInt    *leafbufdisp_d;
 34:   PetscInt    *leafsigdisp_d; /* Copy of leafsigdisp[] on device */
 35:   PetscMPIInt *iranks_d;      /* Copy of the remote part of (leaf) iranks[] on device */
 36:   PetscInt    *ioffset_d;     /* Copy of the remote part of ioffset[] on device */
 37: #endif
 38: } PetscSF_Basic;

 40: static inline PetscErrorCode PetscSFGetRootInfo_Basic(PetscSF sf, PetscMPIInt *nrootranks, PetscMPIInt *ndrootranks, const PetscMPIInt **rootranks, const PetscInt **rootoffset, const PetscInt **rootloc)
 41: {
 42:   PetscSF_Basic *bas = (PetscSF_Basic *)sf->data;

 44:   PetscFunctionBegin;
 45:   if (nrootranks) *nrootranks = bas->niranks;
 46:   if (ndrootranks) *ndrootranks = bas->ndiranks;
 47:   if (rootranks) *rootranks = bas->iranks;
 48:   if (rootoffset) *rootoffset = bas->ioffset;
 49:   if (rootloc) *rootloc = bas->irootloc;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 53: static inline PetscErrorCode PetscSFGetLeafInfo_Basic(PetscSF sf, PetscMPIInt *nleafranks, PetscMPIInt *ndleafranks, const PetscMPIInt **leafranks, const PetscInt **leafoffset, const PetscInt **leafloc, const PetscInt **leafrremote)
 54: {
 55:   PetscFunctionBegin;
 56:   if (nleafranks) *nleafranks = sf->nranks;
 57:   if (ndleafranks) *ndleafranks = sf->ndranks;
 58:   if (leafranks) *leafranks = sf->ranks;
 59:   if (leafoffset) *leafoffset = sf->roffset;
 60:   if (leafloc) *leafloc = sf->rmine;
 61:   if (leafrremote) *leafrremote = sf->rremote;
 62:   PetscFunctionReturn(PETSC_SUCCESS);
 63: }

 65: PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF);
 66: PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF, PetscViewer);
 67: PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF);
 68: PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF);
 69: PETSC_INTERN PetscErrorCode PetscSFBcastBegin_Basic(PetscSF, MPI_Datatype, PetscMemType, const void *, PetscMemType, void *, MPI_Op);
 70: PETSC_INTERN PetscErrorCode PetscSFBcastEnd_Basic(PetscSF, MPI_Datatype, const void *, void *, MPI_Op);
 71: PETSC_INTERN PetscErrorCode PetscSFReduceBegin_Basic(PetscSF, MPI_Datatype, PetscMemType, const void *, PetscMemType, void *, MPI_Op);
 72: PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Basic(PetscSF, MPI_Datatype, const void *, void *, MPI_Op);
 73: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF, MPI_Datatype, PetscMemType, void *, PetscMemType, const void *, void *, MPI_Op);
 74: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF, MPI_Datatype, void *, const void *, void *, MPI_Op);
 75: PETSC_INTERN PetscErrorCode PetscSFCreateEmbeddedRootSF_Basic(PetscSF, PetscInt, const PetscInt *, PetscSF *);
 76: PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF, PetscMPIInt *, const PetscMPIInt **, const PetscInt **, const PetscInt **);

 78: #if defined(PETSC_HAVE_NVSHMEM)
 79: PETSC_INTERN PetscErrorCode PetscSFReset_Basic_NVSHMEM(PetscSF);
 80: #endif