1 #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2 #include <petsc/private/sectionimpl.h> 3 4 /*@C 5 PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout 6 7 Collective 8 9 Input Arguments: 10 + sf - star forest 11 . layout - PetscLayout defining the global space 12 . nleaves - number of leaf vertices on the current process, each of these references a root on any process 13 . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage 14 . localmode - copy mode for ilocal 15 - iremote - remote locations of root vertices for each leaf on the current process 16 17 Level: intermediate 18 19 Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they 20 encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not 21 needed 22 23 .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph() 24 @*/ 25 PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote) 26 { 27 PetscErrorCode ierr; 28 const PetscInt *range; 29 PetscInt i, nroots, ls = -1, ln = -1; 30 PetscMPIInt lr = -1; 31 PetscSFNode *remote; 32 33 PetscFunctionBegin; 34 ierr = PetscLayoutGetLocalSize(layout,&nroots);CHKERRQ(ierr); 35 ierr = PetscLayoutGetRanges(layout,&range);CHKERRQ(ierr); 36 ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr); 37 if (nleaves) { ls = iremote[0] + 1; } 38 for (i=0; i<nleaves; i++) { 39 const PetscInt idx = iremote[i] - ls; 40 if (idx < 0 || idx >= ln) { /* short-circuit the search */ 41 ierr = PetscLayoutFindOwnerIndex(layout,iremote[i],&lr,&remote[i].index);CHKERRQ(ierr); 42 remote[i].rank = lr; 43 ls = range[lr]; 44 ln = range[lr+1] - ls; 45 } else { 46 remote[i].rank = lr; 47 remote[i].index = idx; 48 } 49 } 50 ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 54 /*@ 55 PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections 56 describing the data layout. 57 58 Input Parameters: 59 + sf - The SF 60 . localSection - PetscSection describing the local data layout 61 - globalSection - PetscSection describing the global data layout 62 63 Level: developer 64 65 .seealso: PetscSFSetGraph(), PetscSFSetGraphLayout() 66 @*/ 67 PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) 68 { 69 MPI_Comm comm; 70 PetscLayout layout; 71 const PetscInt *ranges; 72 PetscInt *local; 73 PetscSFNode *remote; 74 PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 75 PetscMPIInt size, rank; 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 80 PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2); 81 PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3); 82 83 ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 84 ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 85 ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 86 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 87 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 88 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 89 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 90 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 91 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 92 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 93 for (p = pStart; p < pEnd; ++p) { 94 PetscInt gdof, gcdof; 95 96 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 97 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 98 if (gcdof > (gdof < 0 ? -(gdof+1) : gdof)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has %D constraints > %D dof", p, gcdof, (gdof < 0 ? -(gdof+1) : gdof)); 99 nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 100 } 101 ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr); 102 ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr); 103 for (p = pStart, l = 0; p < pEnd; ++p) { 104 const PetscInt *cind; 105 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 106 107 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 108 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 109 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 110 ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 111 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 112 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 113 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 114 if (!gdof) continue; /* Censored point */ 115 gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 116 if (gsize != dof-cdof) { 117 if (gsize != dof) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Global dof %D for point %D is neither the constrained size %D, nor the unconstrained %D", gsize, p, dof-cdof, dof); 118 cdof = 0; /* Ignore constraints */ 119 } 120 for (d = 0, c = 0; d < dof; ++d) { 121 if ((c < cdof) && (cind[c] == d)) {++c; continue;} 122 local[l+d-c] = off+d; 123 } 124 if (gdof < 0) { 125 for (d = 0; d < gsize; ++d, ++l) { 126 PetscInt offset = -(goff+1) + d, r; 127 128 ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr); 129 if (r < 0) r = -(r+2); 130 if ((r < 0) || (r >= size)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D mapped to invalid process %D (%D, %D)", p, r, gdof, goff); 131 remote[l].rank = r; 132 remote[l].index = offset - ranges[r]; 133 } 134 } else { 135 for (d = 0; d < gsize; ++d, ++l) { 136 remote[l].rank = rank; 137 remote[l].index = goff+d - ranges[rank]; 138 } 139 } 140 } 141 if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %D != nleaves %D", l, nleaves); 142 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 143 ierr = PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 144 PetscFunctionReturn(0); 145 } 146 147 static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s) 148 { 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 if (!s->bc) { 153 ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr); 154 ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr); 155 } 156 PetscFunctionReturn(0); 157 } 158 159 /*@C 160 PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF 161 162 Collective on sf 163 164 Input Parameters: 165 + sf - The SF 166 - rootSection - Section defined on root space 167 168 Output Parameters: 169 + remoteOffsets - root offsets in leaf storage, or NULL 170 - leafSection - Section defined on the leaf space 171 172 Level: advanced 173 174 .seealso: PetscSFCreate() 175 @*/ 176 PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) 177 { 178 PetscSF embedSF; 179 const PetscInt *indices; 180 IS selected; 181 PetscInt numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c; 182 PetscBool *sub, hasc; 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 ierr = PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr); 187 ierr = PetscSectionGetNumFields(rootSection, &numFields);CHKERRQ(ierr); 188 if (numFields) { 189 IS perm; 190 191 /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys 192 leafSection->perm. To keep this permutation set by the user, we grab 193 the reference before calling PetscSectionSetNumFields() and set it 194 back after. */ 195 ierr = PetscSectionGetPermutation(leafSection, &perm);CHKERRQ(ierr); 196 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 197 ierr = PetscSectionSetNumFields(leafSection, numFields);CHKERRQ(ierr); 198 ierr = PetscSectionSetPermutation(leafSection, perm);CHKERRQ(ierr); 199 ierr = ISDestroy(&perm);CHKERRQ(ierr); 200 } 201 ierr = PetscMalloc1(numFields+2, &sub);CHKERRQ(ierr); 202 sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE; 203 for (f = 0; f < numFields; ++f) { 204 PetscSectionSym sym; 205 const char *name = NULL; 206 PetscInt numComp = 0; 207 208 sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE; 209 ierr = PetscSectionGetFieldComponents(rootSection, f, &numComp);CHKERRQ(ierr); 210 ierr = PetscSectionGetFieldName(rootSection, f, &name);CHKERRQ(ierr); 211 ierr = PetscSectionGetFieldSym(rootSection, f, &sym);CHKERRQ(ierr); 212 ierr = PetscSectionSetFieldComponents(leafSection, f, numComp);CHKERRQ(ierr); 213 ierr = PetscSectionSetFieldName(leafSection, f, name);CHKERRQ(ierr); 214 ierr = PetscSectionSetFieldSym(leafSection, f, sym);CHKERRQ(ierr); 215 for (c = 0; c < rootSection->numFieldComponents[f]; ++c) { 216 ierr = PetscSectionGetComponentName(rootSection, f, c, &name);CHKERRQ(ierr); 217 ierr = PetscSectionSetComponentName(leafSection, f, c, name);CHKERRQ(ierr); 218 } 219 } 220 ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr); 221 ierr = PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);CHKERRQ(ierr); 222 rpEnd = PetscMin(rpEnd,nroots); 223 rpEnd = PetscMax(rpStart,rpEnd); 224 /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */ 225 sub[0] = (PetscBool)(nroots != rpEnd - rpStart); 226 ierr = MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf));CHKERRMPI(ierr); 227 if (sub[0]) { 228 ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr); 229 ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr); 230 ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr); 231 ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr); 232 ierr = ISDestroy(&selected);CHKERRQ(ierr); 233 } else { 234 ierr = PetscObjectReference((PetscObject)sf);CHKERRQ(ierr); 235 embedSF = sf; 236 } 237 ierr = PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd);CHKERRQ(ierr); 238 lpEnd++; 239 240 ierr = PetscSectionSetChart(leafSection, lpStart, lpEnd);CHKERRQ(ierr); 241 242 /* Constrained dof section */ 243 hasc = sub[1]; 244 for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]); 245 246 /* Could fuse these at the cost of copies and extra allocation */ 247 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 248 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 249 if (sub[1]) { 250 ierr = PetscSectionCheckConstraints_Static(rootSection);CHKERRQ(ierr); 251 ierr = PetscSectionCheckConstraints_Static(leafSection);CHKERRQ(ierr); 252 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 253 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 254 } 255 for (f = 0; f < numFields; ++f) { 256 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 257 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 258 if (sub[2+f]) { 259 ierr = PetscSectionCheckConstraints_Static(rootSection->field[f]);CHKERRQ(ierr); 260 ierr = PetscSectionCheckConstraints_Static(leafSection->field[f]);CHKERRQ(ierr); 261 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 262 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 263 } 264 } 265 if (remoteOffsets) { 266 ierr = PetscMalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr); 267 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 268 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 269 } 270 ierr = PetscSectionSetUp(leafSection);CHKERRQ(ierr); 271 if (hasc) { /* need to communicate bcIndices */ 272 PetscSF bcSF; 273 PetscInt *rOffBc; 274 275 ierr = PetscMalloc1(lpEnd - lpStart, &rOffBc);CHKERRQ(ierr); 276 if (sub[1]) { 277 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 278 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 279 ierr = PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF);CHKERRQ(ierr); 280 ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 281 ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 282 ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr); 283 } 284 for (f = 0; f < numFields; ++f) { 285 if (sub[2+f]) { 286 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 287 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 288 ierr = PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF);CHKERRQ(ierr); 289 ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 290 ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE);CHKERRQ(ierr); 291 ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr); 292 } 293 } 294 ierr = PetscFree(rOffBc);CHKERRQ(ierr); 295 } 296 ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr); 297 ierr = PetscFree(sub);CHKERRQ(ierr); 298 ierr = PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr); 299 PetscFunctionReturn(0); 300 } 301 302 /*@C 303 PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes 304 305 Collective on sf 306 307 Input Parameters: 308 + sf - The SF 309 . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots) 310 - leafSection - Data layout of local points for incoming data (this is layout for SF leaves) 311 312 Output Parameter: 313 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 314 315 Level: developer 316 317 .seealso: PetscSFCreate() 318 @*/ 319 PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) 320 { 321 PetscSF embedSF; 322 const PetscInt *indices; 323 IS selected; 324 PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 *remoteOffsets = NULL; 329 ierr = PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);CHKERRQ(ierr); 330 if (numRoots < 0) PetscFunctionReturn(0); 331 ierr = PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr); 332 ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr); 333 ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr); 334 ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr); 335 ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr); 336 ierr = PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr); 337 ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr); 338 ierr = ISDestroy(&selected);CHKERRQ(ierr); 339 ierr = PetscCalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr); 340 ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 341 ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE);CHKERRQ(ierr); 342 ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr); 343 ierr = PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr); 344 PetscFunctionReturn(0); 345 } 346 347 /*@C 348 PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points 349 350 Collective on sf 351 352 Input Parameters: 353 + sf - The SF 354 . rootSection - Data layout of remote points for outgoing data (this is usually the serial section) 355 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL 356 - leafSection - Data layout of local points for incoming data (this is the distributed section) 357 358 Output Parameters: 359 - sectionSF - The new SF 360 361 Note: Either rootSection or remoteOffsets can be specified 362 363 Level: advanced 364 365 .seealso: PetscSFCreate() 366 @*/ 367 PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) 368 { 369 MPI_Comm comm; 370 const PetscInt *localPoints; 371 const PetscSFNode *remotePoints; 372 PetscInt lpStart, lpEnd; 373 PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0; 374 PetscInt *localIndices; 375 PetscSFNode *remoteIndices; 376 PetscInt i, ind; 377 PetscErrorCode ierr; 378 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 381 PetscValidPointer(rootSection,2); 382 /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */ 383 PetscValidPointer(leafSection,4); 384 PetscValidPointer(sectionSF,5); 385 ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 386 ierr = PetscSFCreate(comm, sectionSF);CHKERRQ(ierr); 387 ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr); 388 ierr = PetscSectionGetStorageSize(rootSection, &numSectionRoots);CHKERRQ(ierr); 389 ierr = PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);CHKERRQ(ierr); 390 if (numRoots < 0) PetscFunctionReturn(0); 391 ierr = PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr); 392 for (i = 0; i < numPoints; ++i) { 393 PetscInt localPoint = localPoints ? localPoints[i] : i; 394 PetscInt dof; 395 396 if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 397 ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr); 398 numIndices += dof; 399 } 400 } 401 ierr = PetscMalloc1(numIndices, &localIndices);CHKERRQ(ierr); 402 ierr = PetscMalloc1(numIndices, &remoteIndices);CHKERRQ(ierr); 403 /* Create new index graph */ 404 for (i = 0, ind = 0; i < numPoints; ++i) { 405 PetscInt localPoint = localPoints ? localPoints[i] : i; 406 PetscInt rank = remotePoints[i].rank; 407 408 if ((localPoint >= lpStart) && (localPoint < lpEnd)) { 409 PetscInt remoteOffset = remoteOffsets[localPoint-lpStart]; 410 PetscInt loff, dof, d; 411 412 ierr = PetscSectionGetOffset(leafSection, localPoint, &loff);CHKERRQ(ierr); 413 ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr); 414 for (d = 0; d < dof; ++d, ++ind) { 415 localIndices[ind] = loff+d; 416 remoteIndices[ind].rank = rank; 417 remoteIndices[ind].index = remoteOffset+d; 418 } 419 } 420 } 421 if (numIndices != ind) SETERRQ2(comm, PETSC_ERR_PLIB, "Inconsistency in indices, %D should be %D", ind, numIndices); 422 ierr = PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);CHKERRQ(ierr); 423 ierr = PetscSFSetUp(*sectionSF);CHKERRQ(ierr); 424 ierr = PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr); 425 PetscFunctionReturn(0); 426 } 427 428 /*@C 429 PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects 430 431 Collective 432 433 Input Arguments: 434 + rmap - PetscLayout defining the global root space 435 - lmap - PetscLayout defining the global leaf space 436 437 Output Arguments: 438 . sf - The parallel star forest 439 440 Level: intermediate 441 442 .seealso: PetscSFCreate(), PetscLayoutCreate(), PetscSFSetGraphLayout() 443 @*/ 444 PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF* sf) 445 { 446 PetscErrorCode ierr; 447 PetscInt i,nroots,nleaves = 0; 448 PetscInt rN, lst, len; 449 PetscMPIInt owner = -1; 450 PetscSFNode *remote; 451 MPI_Comm rcomm = rmap->comm; 452 MPI_Comm lcomm = lmap->comm; 453 PetscMPIInt flg; 454 455 PetscFunctionBegin; 456 PetscValidPointer(sf,3); 457 if (!rmap->setupcalled) SETERRQ(rcomm,PETSC_ERR_ARG_WRONGSTATE,"Root layout not setup"); 458 if (!lmap->setupcalled) SETERRQ(lcomm,PETSC_ERR_ARG_WRONGSTATE,"Leaf layout not setup"); 459 ierr = MPI_Comm_compare(rcomm,lcomm,&flg);CHKERRMPI(ierr); 460 if (flg != MPI_CONGRUENT && flg != MPI_IDENT) SETERRQ(rcomm,PETSC_ERR_SUP,"cannot map two layouts with non-matching communicators"); 461 ierr = PetscSFCreate(rcomm,sf);CHKERRQ(ierr); 462 ierr = PetscLayoutGetLocalSize(rmap,&nroots);CHKERRQ(ierr); 463 ierr = PetscLayoutGetSize(rmap,&rN);CHKERRQ(ierr); 464 ierr = PetscLayoutGetRange(lmap,&lst,&len);CHKERRQ(ierr); 465 ierr = PetscMalloc1(len-lst,&remote);CHKERRQ(ierr); 466 for (i = lst; i < len && i < rN; i++) { 467 if (owner < -1 || i >= rmap->range[owner+1]) { 468 ierr = PetscLayoutFindOwner(rmap,i,&owner);CHKERRQ(ierr); 469 } 470 remote[nleaves].rank = owner; 471 remote[nleaves].index = i - rmap->range[owner]; 472 nleaves++; 473 } 474 ierr = PetscSFSetGraph(*sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,remote,PETSC_COPY_VALUES);CHKERRQ(ierr); 475 ierr = PetscFree(remote);CHKERRQ(ierr); 476 PetscFunctionReturn(0); 477 } 478 479 /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */ 480 PetscErrorCode PetscLayoutMapLocal(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 481 { 482 PetscInt *owners = map->range; 483 PetscInt n = map->n; 484 PetscSF sf; 485 PetscInt *lidxs,*work = NULL; 486 PetscSFNode *ridxs; 487 PetscMPIInt rank, p = 0; 488 PetscInt r, len = 0; 489 PetscErrorCode ierr; 490 491 PetscFunctionBegin; 492 if (on) *on = 0; /* squelch -Wmaybe-uninitialized */ 493 /* Create SF where leaves are input idxs and roots are owned idxs */ 494 ierr = MPI_Comm_rank(map->comm,&rank);CHKERRMPI(ierr); 495 ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 496 for (r = 0; r < n; ++r) lidxs[r] = -1; 497 ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 498 for (r = 0; r < N; ++r) { 499 const PetscInt idx = idxs[r]; 500 if (idx < 0 || map->N <= idx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range [0,%D)",idx,map->N); 501 if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 502 ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 503 } 504 ridxs[r].rank = p; 505 ridxs[r].index = idxs[r] - owners[p]; 506 } 507 ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 508 ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 509 ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 510 ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 511 if (ogidxs) { /* communicate global idxs */ 512 PetscInt cum = 0,start,*work2; 513 514 ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 515 ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 516 for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 517 ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRMPI(ierr); 518 start -= cum; 519 cum = 0; 520 for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 521 ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr); 522 ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPI_REPLACE);CHKERRQ(ierr); 523 ierr = PetscFree(work2);CHKERRQ(ierr); 524 } 525 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 526 /* Compress and put in indices */ 527 for (r = 0; r < n; ++r) 528 if (lidxs[r] >= 0) { 529 if (work) work[len] = work[r]; 530 lidxs[len++] = r; 531 } 532 if (on) *on = len; 533 if (oidxs) *oidxs = lidxs; 534 if (ogidxs) *ogidxs = work; 535 PetscFunctionReturn(0); 536 } 537 538 /*@ 539 PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices 540 541 Collective 542 543 Input Parameters: 544 + layout - PetscLayout defining the global index space and the rank that brokers each index 545 . numRootIndices - size of rootIndices 546 . rootIndices - PetscInt array of global indices of which this process requests ownership 547 . rootLocalIndices - root local index permutation (NULL if no permutation) 548 . rootLocalOffset - offset to be added to root local indices 549 . numLeafIndices - size of leafIndices 550 . leafIndices - PetscInt array of global indices with which this process requires data associated 551 . leafLocalIndices - leaf local index permutation (NULL if no permutation) 552 - leafLocalOffset - offset to be added to leaf local indices 553 554 Output Parameter: 555 + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed) 556 - sf - star forest representing the communication pattern from the root space to the leaf space 557 558 Example 1: 559 $ 560 $ rank : 0 1 2 561 $ rootIndices : [1 0 2] [3] [3] 562 $ rootLocalOffset : 100 200 300 563 $ layout : [0 1] [2] [3] 564 $ leafIndices : [0] [2] [0 3] 565 $ leafLocalOffset : 400 500 600 566 $ 567 would build the following SF: 568 $ 569 $ [0] 400 <- (0,101) 570 $ [1] 500 <- (0,102) 571 $ [2] 600 <- (0,101) 572 $ [2] 601 <- (2,300) 573 $ 574 Example 2: 575 $ 576 $ rank : 0 1 2 577 $ rootIndices : [1 0 2] [3] [3] 578 $ rootLocalOffset : 100 200 300 579 $ layout : [0 1] [2] [3] 580 $ leafIndices : rootIndices rootIndices rootIndices 581 $ leafLocalOffset : rootLocalOffset rootLocalOffset rootLocalOffset 582 $ 583 would build the following SF: 584 $ 585 $ [1] 200 <- (2,300) 586 $ 587 Example 3: 588 $ 589 $ No process requests ownership of global index 1, but no process needs it. 590 $ 591 $ rank : 0 1 2 592 $ numRootIndices : 2 1 1 593 $ rootIndices : [0 2] [3] [3] 594 $ rootLocalOffset : 100 200 300 595 $ layout : [0 1] [2] [3] 596 $ numLeafIndices : 1 1 2 597 $ leafIndices : [0] [2] [0 3] 598 $ leafLocalOffset : 400 500 600 599 $ 600 would build the following SF: 601 $ 602 $ [0] 400 <- (0,100) 603 $ [1] 500 <- (0,101) 604 $ [2] 600 <- (0,100) 605 $ [2] 601 <- (2,300) 606 $ 607 608 Notes: 609 The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its 610 local size can be set to PETSC_DECIDE. 611 If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests 612 ownership of x and sends its own rank and the local index of x to process i. 613 If multiple processes request ownership of x, the one with the highest rank is to own x. 614 Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the 615 ownership information of x. 616 The output sf is constructed by associating each leaf point to a root point in this way. 617 618 Suppose there is point data ordered according to the global indices and partitioned according to the given layout. 619 The optional output PetscSF object sfA can be used to push such data to leaf points. 620 621 All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices 622 must cover that of leafIndices, but need not cover the entire layout. 623 624 If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output 625 star forest is almost identity, so will only include non-trivial part of the map. 626 627 Developer Notes: 628 Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using 629 hash(rank, root_local_index) as the bid for the ownership determination. 630 631 Level: advanced 632 633 .seealso: PetscSFCreate() 634 @*/ 635 PetscErrorCode PetscSFCreateByMatchingIndices(PetscLayout layout, PetscInt numRootIndices, const PetscInt *rootIndices, const PetscInt *rootLocalIndices, PetscInt rootLocalOffset, PetscInt numLeafIndices, const PetscInt *leafIndices, const PetscInt *leafLocalIndices, PetscInt leafLocalOffset, PetscSF *sfA, PetscSF *sf) 636 { 637 MPI_Comm comm = layout->comm; 638 PetscMPIInt size, rank; 639 PetscSF sf1; 640 PetscSFNode *owners, *buffer, *iremote; 641 PetscInt *ilocal, nleaves, N, n, i; 642 #if defined(PETSC_USE_DEBUG) 643 PetscInt N1; 644 #endif 645 PetscBool flag; 646 PetscErrorCode ierr; 647 648 PetscFunctionBegin; 649 if (rootIndices) PetscValidIntPointer(rootIndices,3); 650 if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices,4); 651 if (leafIndices) PetscValidIntPointer(leafIndices,7); 652 if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices,8); 653 if (sfA) PetscValidPointer(sfA,10); 654 PetscValidPointer(sf,11); 655 if (numRootIndices < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%D) must be non-negative", numRootIndices); 656 if (numLeafIndices < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%D) must be non-negative", numLeafIndices); 657 ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 658 ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 659 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 660 ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 661 ierr = PetscLayoutGetLocalSize(layout, &n);CHKERRQ(ierr); 662 flag = (PetscBool)(leafIndices == rootIndices); 663 ierr = MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr); 664 if (flag && numLeafIndices != numRootIndices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%D) != numRootIndices(%D)", numLeafIndices, numRootIndices); 665 #if defined(PETSC_USE_DEBUG) 666 N1 = PETSC_MIN_INT; 667 for (i = 0; i < numRootIndices; i++) if (rootIndices[i] > N1) N1 = rootIndices[i]; 668 ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 669 if (N1 >= N) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%D) out of layout range [0,%D)", N1, N); 670 if (!flag) { 671 N1 = PETSC_MIN_INT; 672 for (i = 0; i < numLeafIndices; i++) if (leafIndices[i] > N1) N1 = leafIndices[i]; 673 ierr = MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 674 if (N1 >= N) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%D) out of layout range [0,%D)", N1, N); 675 } 676 #endif 677 /* Reduce: owners -> buffer */ 678 ierr = PetscMalloc1(n, &buffer);CHKERRQ(ierr); 679 ierr = PetscSFCreate(comm, &sf1);CHKERRQ(ierr); 680 ierr = PetscSFSetFromOptions(sf1);CHKERRQ(ierr); 681 ierr = PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices);CHKERRQ(ierr); 682 ierr = PetscMalloc1(numRootIndices, &owners);CHKERRQ(ierr); 683 for (i = 0; i < numRootIndices; ++i) { 684 owners[i].rank = rank; 685 owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i); 686 } 687 for (i = 0; i < n; ++i) { 688 buffer[i].index = -1; 689 buffer[i].rank = -1; 690 } 691 ierr = PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr); 692 ierr = PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC);CHKERRQ(ierr); 693 /* Bcast: buffer -> owners */ 694 if (!flag) { 695 /* leafIndices is different from rootIndices */ 696 ierr = PetscFree(owners);CHKERRQ(ierr); 697 ierr = PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices);CHKERRQ(ierr); 698 ierr = PetscMalloc1(numLeafIndices, &owners);CHKERRQ(ierr); 699 } 700 ierr = PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr); 701 ierr = PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE);CHKERRQ(ierr); 702 for (i = 0; i < numLeafIndices; ++i) if (owners[i].rank < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global point %D was unclaimed", leafIndices[i]); 703 ierr = PetscFree(buffer);CHKERRQ(ierr); 704 if (sfA) {*sfA = sf1;} 705 else {ierr = PetscSFDestroy(&sf1);CHKERRQ(ierr);} 706 /* Create sf */ 707 if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) { 708 /* leaf space == root space */ 709 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) if (owners[i].rank != rank) ++nleaves; 710 ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 711 ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); 712 for (i = 0, nleaves = 0; i < numLeafIndices; ++i) { 713 if (owners[i].rank != rank) { 714 ilocal[nleaves] = leafLocalOffset + i; 715 iremote[nleaves].rank = owners[i].rank; 716 iremote[nleaves].index = owners[i].index; 717 ++nleaves; 718 } 719 } 720 ierr = PetscFree(owners);CHKERRQ(ierr); 721 } else { 722 nleaves = numLeafIndices; 723 ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 724 for (i = 0; i < nleaves; ++i) {ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);} 725 iremote = owners; 726 } 727 ierr = PetscSFCreate(comm, sf);CHKERRQ(ierr); 728 ierr = PetscSFSetFromOptions(*sf);CHKERRQ(ierr); 729 ierr = PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 730 PetscFunctionReturn(0); 731 } 732