Actual source code: sfneighbor.c
petsc-3.13.6 2020-09-29
1: #include <../src/vec/is/sf/impls/basic/sfpack.h>
2: #include <../src/vec/is/sf/impls/basic/sfbasic.h>
4: #if defined(PETSC_HAVE_MPI_NEIGHBORHOOD_COLLECTIVES)
6: typedef struct {
7: SFBASICHEADER;
8: MPI_Comm comms[2]; /* Communicators with distributed topology in both directions */
9: PetscBool initialized[2]; /* Are the two communicators initialized? */
10: PetscMPIInt *rootdispls,*rootcounts,*leafdispls,*leafcounts; /* displs/counts for non-distinguished ranks */
11: PetscInt rootdegree,leafdegree;
12: } PetscSF_Neighbor;
14: /*===================================================================================*/
15: /* Internal utility routines */
16: /*===================================================================================*/
18: /* Get the communicator with distributed graph topology, which is not cheap to build so we do it on demand (instead of at PetscSFSetUp time) */
19: static PetscErrorCode PetscSFGetDistComm_Neighbor(PetscSF sf,PetscSFDirection direction,MPI_Comm *distcomm)
20: {
21: PetscErrorCode ierr;
22: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
23: PetscInt nrootranks,ndrootranks,nleafranks,ndleafranks;
24: const PetscMPIInt *rootranks,*leafranks;
25: MPI_Comm comm;
28: PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,&rootranks,NULL,NULL); /* Which ranks will access my roots (I am a destination) */
29: PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,&leafranks,NULL,NULL,NULL); /* My leaves will access whose roots (I am a source) */
31: if (!dat->initialized[direction]) {
32: const PetscMPIInt indegree = nrootranks-ndrootranks,*sources = rootranks+ndrootranks;
33: const PetscMPIInt outdegree = nleafranks-ndleafranks,*destinations = leafranks+ndleafranks;
34: MPI_Comm *mycomm = &dat->comms[direction];
35: PetscObjectGetComm((PetscObject)sf,&comm);
36: if (direction == PETSCSF_LEAF2../../../../../..) {
37: MPI_Dist_graph_create_adjacent(comm,indegree,sources,dat->rootcounts/*src weights*/,outdegree,destinations,dat->leafcounts/*dest weights*/,MPI_INFO_NULL,1/*reorder*/,mycomm);
38: } else { /* PETSCSF_../../../../../..2LEAF, reverse src & dest */
39: MPI_Dist_graph_create_adjacent(comm,outdegree,destinations,dat->leafcounts/*src weights*/,indegree,sources,dat->rootcounts/*dest weights*/,MPI_INFO_NULL,1/*reorder*/,mycomm);
40: }
41: dat->initialized[direction] = PETSC_TRUE;
42: }
43: *distcomm = dat->comms[direction];
44: return(0);
45: }
47: /*===================================================================================*/
48: /* Implementations of SF public APIs */
49: /*===================================================================================*/
50: static PetscErrorCode PetscSFSetUp_Neighbor(PetscSF sf)
51: {
52: PetscErrorCode ierr;
53: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
54: PetscInt i,j,nrootranks,ndrootranks,nleafranks,ndleafranks;
55: const PetscInt *rootoffset,*leafoffset;
56: PetscMPIInt m,n;
59: /* SFNeighbor inherits from Basic */
60: PetscSFSetUp_Basic(sf);
61: /* SFNeighbor specific */
62: sf->persistent = PETSC_FALSE;
63: PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,NULL,&rootoffset,NULL);
64: PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,NULL,&leafoffset,NULL,NULL);
65: dat->rootdegree = nrootranks-ndrootranks;
66: dat->leafdegree = nleafranks-ndleafranks;
67: sf->nleafreqs = 0;
68: dat->nrootreqs = 1;
70: /* Only setup MPI displs/counts for non-distinguished ranks. Distinguished ranks use shared memory */
71: PetscMalloc4(dat->rootdegree,&dat->rootdispls,dat->rootdegree,&dat->rootcounts,dat->leafdegree,&dat->leafdispls,dat->leafdegree,&dat->leafcounts);
72: for (i=ndrootranks,j=0; i<nrootranks; i++,j++) {
73: PetscMPIIntCast(rootoffset[i]-rootoffset[ndrootranks],&m); dat->rootdispls[j] = m;
74: PetscMPIIntCast(rootoffset[i+1]-rootoffset[i], &n); dat->rootcounts[j] = n;
75: }
77: for (i=ndleafranks,j=0; i<nleafranks; i++,j++) {
78: PetscMPIIntCast(leafoffset[i]-leafoffset[ndleafranks],&m); dat->leafdispls[j] = m;
79: PetscMPIIntCast(leafoffset[i+1]-leafoffset[i], &n); dat->leafcounts[j] = n;
80: }
81: return(0);
82: }
84: static PetscErrorCode PetscSFReset_Neighbor(PetscSF sf)
85: {
86: PetscErrorCode ierr;
87: PetscInt i;
88: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
91: if (dat->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
92: PetscFree4(dat->rootdispls,dat->rootcounts,dat->leafdispls,dat->leafcounts);
93: for (i=0; i<2; i++) {
94: if (dat->initialized[i]) {
95: MPI_Comm_free(&dat->comms[i]);
96: dat->initialized[i] = PETSC_FALSE;
97: }
98: }
99: PetscSFReset_Basic(sf); /* Common part */
100: return(0);
101: }
103: static PetscErrorCode PetscSFDestroy_Neighbor(PetscSF sf)
104: {
108: PetscSFReset_Neighbor(sf);
109: PetscFree(sf->data);
110: return(0);
111: }
113: static PetscErrorCode PetscSFBcastAndOpBegin_Neighbor(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
114: {
115: PetscErrorCode ierr;
116: PetscSFLink link;
117: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
118: MPI_Comm distcomm;
119: void *rootbuf = NULL,*leafbuf = NULL;
120: MPI_Request *req;
123: PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_BCAST,&link);
124: PetscSFLinkPackRootData(sf,link,PETSCSF_REMOTE,rootdata);
125: /* Do neighborhood alltoallv for remote ranks */
126: PetscSFGetDistComm_Neighbor(sf,PETSCSF_../../../../../..2LEAF,&distcomm);
127: PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_../../../../../..2LEAF,&rootbuf,&leafbuf,&req,NULL);
128: MPI_Start_ineighbor_alltoallv(dat->rootdegree,dat->leafdegree,rootbuf,dat->rootcounts,dat->rootdispls,unit,leafbuf,dat->leafcounts,dat->leafdispls,unit,distcomm,req);
129: PetscSFLinkBcastAndOpLocal(sf,link,rootdata,leafdata,op);
130: return(0);
131: }
133: PETSC_STATIC_INLINE PetscErrorCode PetscSFLeafToRootBegin_Neighbor(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op,PetscSFOperation sfop,PetscSFLink *out)
134: {
135: PetscErrorCode ierr;
136: PetscSFLink link;
137: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
138: MPI_Comm distcomm = MPI_COMM_NULL;
139: void *rootbuf = NULL,*leafbuf = NULL;
140: MPI_Request *req = NULL;
143: PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,sfop,&link);
144: PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);
145: /* Do neighborhood alltoallv for remote ranks */
146: PetscSFGetDistComm_Neighbor(sf,PETSCSF_LEAF2../../../../../..,&distcomm);
147: PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_LEAF2../../../../../..,&rootbuf,&leafbuf,&req,NULL);
148: MPI_Start_ineighbor_alltoallv(dat->leafdegree,dat->rootdegree,leafbuf,dat->leafcounts,dat->leafdispls,unit,rootbuf,dat->rootcounts,dat->rootdispls,unit,distcomm,req);
149: *out = link;
150: return(0);
151: }
153: static PetscErrorCode PetscSFReduceBegin_Neighbor(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
154: {
155: PetscErrorCode ierr;
156: PetscSFLink link = NULL;
159: PetscSFLeafToRootBegin_Neighbor(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op,PETSCSF_REDUCE,&link);
160: PetscSFLinkReduceLocal(sf,link,leafdata,rootdata,op);
161: return(0);
162: }
164: static PetscErrorCode PetscSFFetchAndOpBegin_Neighbor(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,void *rootdata,PetscMemType leafmtype,const void *leafdata,void *leafupdate,MPI_Op op)
165: {
166: PetscErrorCode ierr;
167: PetscSFLink link = NULL;
170: PetscSFLeafToRootBegin_Neighbor(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op,PETSCSF_FETCH,&link);
171: PetscSFLinkFetchAndOpLocal(sf,link,rootdata,leafdata,leafupdate,op);
172: return(0);
173: }
175: static PetscErrorCode PetscSFFetchAndOpEnd_Neighbor(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
176: {
177: PetscErrorCode ierr;
178: PetscSFLink link = NULL;
179: MPI_Comm comm = MPI_COMM_NULL;
180: PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
181: void *rootbuf = NULL,*leafbuf = NULL;
184: PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
185: PetscSFLinkMPIWaitall(sf,link,PETSCSF_LEAF2../../../../../..);
186: /* Process remote fetch-and-op */
187: PetscSFLinkFetchRootData(sf,link,PETSCSF_REMOTE,rootdata,op);
189: /* Bcast the updated rootbuf back to leaves */
190: PetscSFGetDistComm_Neighbor(sf,PETSCSF_../../../../../..2LEAF,&comm);
191: PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_../../../../../..2LEAF,&rootbuf,&leafbuf,NULL,NULL);
192: MPI_Start_neighbor_alltoallv(dat->rootdegree,dat->leafdegree,rootbuf,dat->rootcounts,dat->rootdispls,unit,leafbuf,dat->leafcounts,dat->leafdispls,unit,comm);
193: PetscSFLinkUnpackLeafData(sf,link,PETSCSF_REMOTE,leafupdate,MPIU_REPLACE);
194: PetscSFLinkReclaim(sf,&link);
195: return(0);
196: }
198: PETSC_INTERN PetscErrorCode PetscSFCreate_Neighbor(PetscSF sf)
199: {
200: PetscErrorCode ierr;
201: PetscSF_Neighbor *dat;
204: sf->ops->CreateEmbeddedSF = PetscSFCreateEmbeddedSF_Basic;
205: sf->ops->BcastAndOpEnd = PetscSFBcastAndOpEnd_Basic;
206: sf->ops->ReduceEnd = PetscSFReduceEnd_Basic;
207: sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Basic;
208: sf->ops->View = PetscSFView_Basic;
210: sf->ops->SetUp = PetscSFSetUp_Neighbor;
211: sf->ops->Reset = PetscSFReset_Neighbor;
212: sf->ops->Destroy = PetscSFDestroy_Neighbor;
213: sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Neighbor;
214: sf->ops->ReduceBegin = PetscSFReduceBegin_Neighbor;
215: sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Neighbor;
216: sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Neighbor;
218: PetscNewLog(sf,&dat);
219: sf->data = (void*)dat;
220: return(0);
221: }
222: #endif