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