1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 /*@ 4 DMPlexSetMigrationSF - Sets the `PetscSF` for migrating from a parent `DM` into this `DM` 5 6 Logically Collective 7 8 Input Parameters: 9 + dm - The `DM` 10 - migrationSF - The `PetscSF` 11 12 Level: intermediate 13 14 Note: 15 It is necessary to call this in order to have `DMCreateSubDM()` or `DMCreateSuperDM()` build the Global-To-Natural map 16 17 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()` 18 @*/ 19 PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) 20 { 21 PetscFunctionBegin; 22 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 23 if (migrationSF) PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2); 24 PetscCall(PetscObjectReference((PetscObject)migrationSF)); 25 PetscCall(PetscSFDestroy(&dm->sfMigration)); 26 dm->sfMigration = migrationSF; 27 PetscFunctionReturn(PETSC_SUCCESS); 28 } 29 30 /*@ 31 DMPlexGetMigrationSF - Gets the `PetscSF` for migrating from a parent `DM` into this `DM` 32 33 Note Collective 34 35 Input Parameter: 36 . dm - The `DM` 37 38 Output Parameter: 39 . migrationSF - The `PetscSF` 40 41 Level: intermediate 42 43 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF` 44 @*/ 45 PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) 46 { 47 PetscFunctionBegin; 48 *migrationSF = dm->sfMigration; 49 PetscFunctionReturn(PETSC_SUCCESS); 50 } 51 52 /*@ 53 DMPlexSetGlobalToNaturalSF - Sets the `PetscSF` for mapping Global `Vec` to the Natural `Vec` 54 55 Input Parameters: 56 + dm - The `DM` 57 - naturalSF - The `PetscSF` 58 59 Level: intermediate 60 61 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()` 62 @*/ 63 PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) 64 { 65 PetscFunctionBegin; 66 dm->sfNatural = naturalSF; 67 PetscCall(PetscObjectReference((PetscObject)naturalSF)); 68 dm->useNatural = naturalSF ? PETSC_TRUE : PETSC_FALSE; 69 PetscFunctionReturn(PETSC_SUCCESS); 70 } 71 72 /*@ 73 DMPlexGetGlobalToNaturalSF - Gets the `PetscSF` for mapping Global `Vec` to the Natural `Vec` 74 75 Input Parameter: 76 . dm - The `DM` 77 78 Output Parameter: 79 . naturalSF - The `PetscSF` 80 81 Level: intermediate 82 83 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF` 84 @*/ 85 PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) 86 { 87 PetscFunctionBegin; 88 *naturalSF = dm->sfNatural; 89 PetscFunctionReturn(PETSC_SUCCESS); 90 } 91 92 /*@ 93 DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec` 94 95 Input Parameters: 96 + dm - The redistributed `DM` 97 . section - The local `PetscSection` describing the `Vec` before the mesh was distributed, or `NULL` if not available 98 - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed 99 100 Output Parameter: 101 . sfNatural - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering 102 103 Level: intermediate 104 105 Note: 106 This is not typically called by the user. 107 108 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()` 109 @*/ 110 PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 111 { 112 MPI_Comm comm; 113 PetscSF sf, sfEmbed, sfField; 114 PetscSection gSection, sectionDist, gLocSection; 115 PetscInt *spoints, *remoteOffsets; 116 PetscInt ssize, pStart, pEnd, p, localSize, maxStorageSize; 117 PetscBool destroyFlag = PETSC_FALSE, debug = PETSC_FALSE; 118 119 PetscFunctionBegin; 120 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 121 if (!sfMigration) { 122 /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */ 123 *sfNatural = NULL; 124 PetscFunctionReturn(PETSC_SUCCESS); 125 } else if (!section) { 126 /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */ 127 PetscSF sfMigrationInv; 128 PetscSection localSection; 129 130 PetscCall(DMGetLocalSection(dm, &localSection)); 131 PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv)); 132 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 133 PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section)); 134 PetscCall(PetscSFDestroy(&sfMigrationInv)); 135 destroyFlag = PETSC_TRUE; 136 } 137 if (debug) PetscCall(PetscSFView(sfMigration, NULL)); 138 /* Create a new section from distributing the original section */ 139 PetscCall(PetscSectionCreate(comm, §ionDist)); 140 PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist)); 141 PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section")); 142 if (debug) PetscCall(PetscSectionView(sectionDist, NULL)); 143 PetscCall(DMSetLocalSection(dm, sectionDist)); 144 /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */ 145 PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize)); 146 PetscCallMPI(MPIU_Allreduce(&localSize, &maxStorageSize, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 147 if (maxStorageSize) { 148 const PetscInt *leaves; 149 PetscInt *sortleaves, *indices; 150 PetscInt Nl; 151 152 /* Get a pruned version of migration SF */ 153 PetscCall(DMGetGlobalSection(dm, &gSection)); 154 if (debug) PetscCall(PetscSectionView(gSection, NULL)); 155 PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL)); 156 PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 157 for (p = pStart, ssize = 0; p < pEnd; ++p) { 158 PetscInt dof, off; 159 160 PetscCall(PetscSectionGetDof(gSection, p, &dof)); 161 PetscCall(PetscSectionGetOffset(gSection, p, &off)); 162 if ((dof > 0) && (off >= 0)) ++ssize; 163 } 164 PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices)); 165 for (p = 0; p < Nl; ++p) { 166 sortleaves[p] = leaves ? leaves[p] : p; 167 indices[p] = p; 168 } 169 PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices)); 170 for (p = pStart, ssize = 0; p < pEnd; ++p) { 171 PetscInt dof, off, loc; 172 173 PetscCall(PetscSectionGetDof(gSection, p, &dof)); 174 PetscCall(PetscSectionGetOffset(gSection, p, &off)); 175 if ((dof > 0) && (off >= 0)) { 176 PetscCall(PetscFindInt(p, Nl, sortleaves, &loc)); 177 PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " with nonzero dof is not a leaf of the migration SF", p); 178 spoints[ssize++] = indices[loc]; 179 } 180 } 181 PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed)); 182 PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF")); 183 PetscCall(PetscFree3(spoints, sortleaves, indices)); 184 if (debug) PetscCall(PetscSFView(sfEmbed, NULL)); 185 /* Create the SF associated with this section 186 Roots are natural dofs, leaves are global dofs */ 187 PetscCall(DMGetPointSF(dm, &sf)); 188 PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &gLocSection)); 189 PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField)); 190 PetscCall(PetscSFDestroy(&sfEmbed)); 191 PetscCall(PetscSectionDestroy(&gLocSection)); 192 PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF")); 193 if (debug) PetscCall(PetscSFView(sfField, NULL)); 194 /* Invert the field SF 195 Roots are global dofs, leaves are natural dofs */ 196 PetscCall(PetscSFCreateInverseSF(sfField, sfNatural)); 197 PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF")); 198 PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view")); 199 /* Clean up */ 200 PetscCall(PetscSFDestroy(&sfField)); 201 } else { 202 *sfNatural = NULL; 203 } 204 PetscCall(PetscSectionDestroy(§ionDist)); 205 PetscCall(PetscFree(remoteOffsets)); 206 if (destroyFlag) PetscCall(PetscSectionDestroy(§ion)); 207 PetscFunctionReturn(PETSC_SUCCESS); 208 } 209 210 /*@ 211 DMPlexMigrateGlobalToNaturalSF - Migrates the input `sfNatural` based on sfMigration 212 213 Input Parameters: 214 + dmOld - The original `DM` 215 . dmNew - The `DM` to be migrated to 216 . sfNaturalOld - The sfNatural for the `dmOld` 217 - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed 218 219 Output Parameter: 220 . sfNaturalNew - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering 221 222 Level: intermediate 223 224 Notes: 225 `sfNaturalOld` maps from the old Global section (roots) to the natural Vec layout (leaves, may or may not be described by a PetscSection). 226 `DMPlexMigrateGlobalToNaturalSF` creates an SF to map from the old global section to the new global section (generated from `sfMigration`). 227 That SF is then composed with the `sfNaturalOld` to generate `sfNaturalNew`. 228 This also distributes and sets the local section for `dmNew`. 229 230 This is not typically called by the user. 231 232 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()` 233 @*/ 234 PetscErrorCode DMPlexMigrateGlobalToNaturalSF(DM dmOld, DM dmNew, PetscSF sfNaturalOld, PetscSF sfMigration, PetscSF *sfNaturalNew) 235 { 236 MPI_Comm comm; 237 PetscSection oldGlobalSection, newGlobalSection; 238 PetscInt *remoteOffsets; 239 PetscBool debug = PETSC_FALSE; 240 241 PetscFunctionBegin; 242 PetscCall(PetscObjectGetComm((PetscObject)dmNew, &comm)); 243 if (!sfMigration) { 244 /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */ 245 *sfNaturalNew = NULL; 246 PetscFunctionReturn(PETSC_SUCCESS); 247 } 248 if (debug) PetscCall(PetscSFView(sfMigration, NULL)); 249 250 { // Create oldGlobalSection and newGlobalSection *with* localOffsets 251 PetscSection oldLocalSection, newLocalSection; 252 PetscSF pointSF; 253 254 PetscCall(DMGetLocalSection(dmOld, &oldLocalSection)); 255 PetscCall(DMGetPointSF(dmOld, &pointSF)); 256 PetscCall(PetscSectionCreateGlobalSection(oldLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &oldGlobalSection)); 257 258 PetscCall(PetscSectionCreate(comm, &newLocalSection)); 259 PetscCall(PetscSFDistributeSection(sfMigration, oldLocalSection, NULL, newLocalSection)); 260 PetscCall(DMSetLocalSection(dmNew, newLocalSection)); 261 262 PetscCall(PetscSectionCreate(comm, &newGlobalSection)); 263 PetscCall(DMGetPointSF(dmNew, &pointSF)); 264 PetscCall(PetscSectionCreateGlobalSection(newLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &newGlobalSection)); 265 266 PetscCall(PetscObjectSetName((PetscObject)oldLocalSection, "Old Local Section")); 267 if (debug) PetscCall(PetscSectionView(oldLocalSection, NULL)); 268 PetscCall(PetscObjectSetName((PetscObject)oldGlobalSection, "Old Global Section")); 269 if (debug) PetscCall(PetscSectionView(oldGlobalSection, NULL)); 270 PetscCall(PetscObjectSetName((PetscObject)newLocalSection, "New Local Section")); 271 if (debug) PetscCall(PetscSectionView(newLocalSection, NULL)); 272 PetscCall(PetscObjectSetName((PetscObject)newGlobalSection, "New Global Section")); 273 if (debug) PetscCall(PetscSectionView(newGlobalSection, NULL)); 274 PetscCall(PetscSectionDestroy(&newLocalSection)); 275 } 276 277 { // Create remoteOffsets array, mapping the oldGlobalSection offsets to the local points (according to sfMigration) 278 PetscInt lpStart, lpEnd, rpStart, rpEnd; 279 280 PetscCall(PetscSectionGetChart(oldGlobalSection, &rpStart, &rpEnd)); 281 PetscCall(PetscSectionGetChart(newGlobalSection, &lpStart, &lpEnd)); 282 283 // in `PetscSFDistributeSection` (where this is taken from), it possibly makes a new embedded SF. Should possibly do that here? 284 PetscCall(PetscMalloc1(lpEnd - lpStart, &remoteOffsets)); 285 PetscCall(PetscSFBcastBegin(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE)); 286 PetscCall(PetscSFBcastEnd(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE)); 287 if (debug) { 288 PetscViewer viewer; 289 290 PetscCall(PetscPrintf(comm, "Remote Offsets:\n")); 291 PetscCall(PetscViewerASCIIGetStdout(comm, &viewer)); 292 PetscCall(PetscIntView(lpEnd - lpStart, remoteOffsets, viewer)); 293 } 294 } 295 296 { // Create SF from oldGlobalSection to newGlobalSection and compose with sfNaturalOld 297 PetscSF oldglob_to_newglob_sf, newglob_to_oldglob_sf; 298 299 PetscCall(PetscSFCreateSectionSF(sfMigration, oldGlobalSection, remoteOffsets, newGlobalSection, &oldglob_to_newglob_sf)); 300 PetscCall(PetscObjectSetName((PetscObject)oldglob_to_newglob_sf, "OldGlobal-to-NewGlobal SF")); 301 if (debug) PetscCall(PetscSFView(oldglob_to_newglob_sf, NULL)); 302 303 PetscCall(PetscSFCreateInverseSF(oldglob_to_newglob_sf, &newglob_to_oldglob_sf)); 304 PetscCall(PetscObjectSetName((PetscObject)newglob_to_oldglob_sf, "NewGlobal-to-OldGlobal SF")); 305 PetscCall(PetscObjectViewFromOptions((PetscObject)newglob_to_oldglob_sf, (PetscObject)dmOld, "-natural_migrate_sf_view")); 306 PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNaturalOld, sfNaturalNew)); 307 308 PetscCall(PetscSFDestroy(&oldglob_to_newglob_sf)); 309 PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf)); 310 } 311 312 PetscCall(PetscSectionDestroy(&oldGlobalSection)); 313 PetscCall(PetscSectionDestroy(&newGlobalSection)); 314 PetscCall(PetscFree(remoteOffsets)); 315 PetscFunctionReturn(PETSC_SUCCESS); 316 } 317 318 /*@ 319 DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order. 320 321 Collective 322 323 Input Parameters: 324 + dm - The distributed `DMPLEX` 325 - gv - The global `Vec` 326 327 Output Parameter: 328 . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv` 329 330 Level: intermediate 331 332 Note: 333 The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`. 334 335 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 336 @*/ 337 PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 338 { 339 const PetscScalar *inarray; 340 PetscScalar *outarray; 341 MPI_Comm comm; 342 PetscMPIInt size; 343 344 PetscFunctionBegin; 345 PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 346 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 347 PetscCallMPI(MPI_Comm_size(comm, &size)); 348 if (dm->sfNatural) { 349 if (PetscDefined(USE_DEBUG)) { 350 PetscSection gs; 351 PetscInt Nl, n; 352 353 PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL)); 354 PetscCall(VecGetLocalSize(nv, &n)); 355 PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Natural vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of natural section", n, Nl); 356 357 PetscCall(DMGetGlobalSection(dm, &gs)); 358 PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl)); 359 PetscCall(VecGetLocalSize(gv, &n)); 360 PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Global vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of global section", n, Nl); 361 } 362 PetscCall(VecGetArray(nv, &outarray)); 363 PetscCall(VecGetArrayRead(gv, &inarray)); 364 PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 365 PetscCall(VecRestoreArrayRead(gv, &inarray)); 366 PetscCall(VecRestoreArray(nv, &outarray)); 367 } else if (size == 1) { 368 PetscCall(VecCopy(gv, nv)); 369 } else { 370 PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 371 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute()."); 372 } 373 PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 374 PetscFunctionReturn(PETSC_SUCCESS); 375 } 376 377 /*@ 378 DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order. 379 380 Collective 381 382 Input Parameters: 383 + dm - The distributed `DMPLEX` 384 - gv - The global `Vec` 385 386 Output Parameter: 387 . nv - The natural `Vec` 388 389 Level: intermediate 390 391 Note: 392 The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`. 393 394 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 395 @*/ 396 PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 397 { 398 const PetscScalar *inarray; 399 PetscScalar *outarray; 400 PetscMPIInt size; 401 402 PetscFunctionBegin; 403 PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 404 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 405 if (dm->sfNatural) { 406 PetscCall(VecGetArrayRead(gv, &inarray)); 407 PetscCall(VecGetArray(nv, &outarray)); 408 PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 409 PetscCall(VecRestoreArrayRead(gv, &inarray)); 410 PetscCall(VecRestoreArray(nv, &outarray)); 411 } else if (size == 1) { 412 } else { 413 PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 414 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute()."); 415 } 416 PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 417 PetscFunctionReturn(PETSC_SUCCESS); 418 } 419 420 /*@ 421 DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order. 422 423 Collective 424 425 Input Parameters: 426 + dm - The distributed `DMPLEX` 427 - nv - The natural `Vec` 428 429 Output Parameter: 430 . gv - The global `Vec` 431 432 Level: intermediate 433 434 Note: 435 The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`. 436 437 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()` 438 @*/ 439 PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 440 { 441 const PetscScalar *inarray; 442 PetscScalar *outarray; 443 PetscMPIInt size; 444 445 PetscFunctionBegin; 446 PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 447 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 448 if (dm->sfNatural) { 449 /* We only have access to the SF that goes from Global to Natural. 450 Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 451 Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 452 PetscCall(VecZeroEntries(gv)); 453 PetscCall(VecGetArray(gv, &outarray)); 454 PetscCall(VecGetArrayRead(nv, &inarray)); 455 PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 456 PetscCall(VecRestoreArrayRead(nv, &inarray)); 457 PetscCall(VecRestoreArray(gv, &outarray)); 458 } else if (size == 1) { 459 PetscCall(VecCopy(nv, gv)); 460 } else { 461 PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 462 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute()."); 463 } 464 PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 465 PetscFunctionReturn(PETSC_SUCCESS); 466 } 467 468 /*@ 469 DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order. 470 471 Collective 472 473 Input Parameters: 474 + dm - The distributed `DMPLEX` 475 - nv - The natural `Vec` 476 477 Output Parameter: 478 . gv - The global `Vec` 479 480 Level: intermediate 481 482 Note: 483 The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`. 484 485 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 486 @*/ 487 PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 488 { 489 const PetscScalar *inarray; 490 PetscScalar *outarray; 491 PetscMPIInt size; 492 493 PetscFunctionBegin; 494 PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 495 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 496 if (dm->sfNatural) { 497 PetscCall(VecGetArrayRead(nv, &inarray)); 498 PetscCall(VecGetArray(gv, &outarray)); 499 PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 500 PetscCall(VecRestoreArrayRead(nv, &inarray)); 501 PetscCall(VecRestoreArray(gv, &outarray)); 502 } else if (size == 1) { 503 } else { 504 PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 505 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute()."); 506 } 507 PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 508 PetscFunctionReturn(PETSC_SUCCESS); 509 } 510 511 /*@ 512 DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution. 513 514 Collective 515 516 Input Parameter: 517 . dm - The distributed `DMPLEX` 518 519 Output Parameter: 520 . nv - The natural `Vec` 521 522 Level: intermediate 523 524 Note: 525 The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`. 526 527 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 528 @*/ 529 PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv) 530 { 531 PetscMPIInt size; 532 533 PetscFunctionBegin; 534 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 535 if (dm->sfNatural) { 536 PetscInt nleaves, bs, maxbs; 537 Vec v; 538 539 /* 540 Setting the natural vector block size. 541 We can't get it from a global vector because of constraints, and the block size in the local vector 542 may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1 543 */ 544 PetscCall(DMGetLocalVector(dm, &v)); 545 PetscCall(VecGetBlockSize(v, &bs)); 546 PetscCallMPI(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 547 if (bs == 1 && maxbs > 1) bs = maxbs; 548 PetscCall(DMRestoreLocalVector(dm, &v)); 549 550 PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL)); 551 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv)); 552 PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE)); 553 PetscCall(VecSetBlockSize(*nv, bs)); 554 PetscCall(VecSetType(*nv, dm->vectype)); 555 PetscCall(VecSetDM(*nv, dm)); 556 } else if (size == 1) { 557 PetscCall(DMCreateLocalVector(dm, nv)); 558 } else { 559 PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 560 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute()."); 561 } 562 PetscFunctionReturn(PETSC_SUCCESS); 563 } 564