Actual source code: sfbasic.h
1: #pragma once
3: #include <petsc/private/sfimpl.h>
5: typedef struct _n_PetscSFLink *PetscSFLink;
7: #define SFBASICHEADER \
8: PetscMPIInt niranks; /* Number of incoming ranks (ranks accessing my roots) */ \
9: PetscMPIInt ndiranks; /* Number of incoming ranks (ranks accessing my roots) in distinguished set */ \
10: PetscMPIInt *iranks; /* Array of ranks that reference my roots */ \
11: PetscInt itotal; /* Total number of graph edges referencing my roots */ \
12: PetscInt *ioffset; /* Array of length niranks+1 holding offset in irootloc[] for each rank */ \
13: PetscInt *irootloc; /* Incoming roots referenced by ranks starting at ioffset[rank] */ \
14: PetscInt *irootloc_d[2]; /* A copy of irootloc[local/remote] in device memory if needed */ \
15: PetscInt rootbuflen[2]; /* Length (in unit) of root buffers, in layout of [PETSCSF_LOCAL/REMOTE] */ \
16: PetscBool rootcontig[2]; /* True means the local/remote segments of indices in irootloc[] are contiguous ... */ \
17: PetscInt rootstart[2]; /* ... and start from rootstart[0] and rootstart[1] respectively */ \
18: PetscSFPackOpt rootpackopt[2]; /* Pack optimization plans based on patterns in irootloc[]. NULL for no optimizations */ \
19: PetscSFPackOpt rootpackopt_d[2]; /* Copy of rootpackopt[] on device if needed */ \
20: PetscBool rootdups[2]; /* Indices of roots in irootloc[local/remote] have dups. Used for data-race test */ \
21: PetscInt nrootreqs; /* Number of MPI requests */ \
22: PetscSFLink avail; /* One or more entries per MPI Datatype, lazily constructed */ \
23: PetscSFLink inuse /* Buffers being used for transactions that have not yet completed */
25: typedef struct {
26: SFBASICHEADER;
27: #if defined(PETSC_HAVE_NVSHMEM)
28: PetscInt rootbuflen_rmax; /* max rootbuflen[REMOTE] over comm */
29: PetscInt nRemoteLeafRanks; /* niranks - ndiranks */
30: PetscInt nRemoteLeafRanksMax; /* max nRemoteLeafRanks over comm */
32: 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 */
33: PetscInt *leafsigdisp; /* [nRemoteLeafRanks]. For my i-th remote leaf rank, I am its leafsigdisp[i]-th root rank */
35: PetscInt *leafbufdisp_d;
36: PetscInt *leafsigdisp_d; /* Copy of leafsigdisp[] on device */
37: PetscMPIInt *iranks_d; /* Copy of the remote part of (leaf) iranks[] on device */
38: PetscInt *ioffset_d; /* Copy of the remote part of ioffset[] on device */
39: #endif
40: } PetscSF_Basic;
42: static inline PetscErrorCode PetscSFGetRootInfo_Basic(PetscSF sf, PetscInt *nrootranks, PetscInt *ndrootranks, const PetscMPIInt **rootranks, const PetscInt **rootoffset, const PetscInt **rootloc)
43: {
44: PetscSF_Basic *bas = (PetscSF_Basic *)sf->data;
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(PETSC_SUCCESS);
53: }
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: {
57: PetscFunctionBegin;
58: if (nleafranks) *nleafranks = sf->nranks;
59: if (ndleafranks) *ndleafranks = sf->ndranks;
60: if (leafranks) *leafranks = sf->ranks;
61: if (leafoffset) *leafoffset = sf->roffset;
62: if (leafloc) *leafloc = sf->rmine;
63: if (leafrremote) *leafrremote = sf->rremote;
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF);
68: PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF, PetscViewer);
69: PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF);
70: PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF);
71: PETSC_INTERN PetscErrorCode PetscSFBcastEnd_Basic(PetscSF, MPI_Datatype, const void *, 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 PetscSFCreateEmbeddedRootSF_Basic(PetscSF, PetscInt, const PetscInt *, PetscSF *);
75: PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF, PetscInt *, const PetscMPIInt **, const PetscInt **, const PetscInt **);
77: #if defined(PETSC_HAVE_NVSHMEM)
78: PETSC_INTERN PetscErrorCode PetscSFReset_Basic_NVSHMEM(PetscSF);
79: #endif