1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 /*@ 4 DMPlexSetMigrationSF - Sets the SF for migrating from a parent DM into this DM 5 6 Input Parameters: 7 + dm - The DM 8 - naturalSF - The PetscSF 9 10 Note: It is necessary to call this in order to have DMCreateSubDM() or DMCreateSuperDM() build the Global-To-Natural map 11 12 Level: intermediate 13 14 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexGetMigrationSF() 15 @*/ 16 PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) 17 { 18 PetscErrorCode ierr; 19 PetscFunctionBegin; 20 dm->sfMigration = migrationSF; 21 ierr = PetscObjectReference((PetscObject) migrationSF);CHKERRQ(ierr); 22 PetscFunctionReturn(0); 23 } 24 25 /*@ 26 DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM 27 28 Input Parameter: 29 . dm - The DM 30 31 Output Parameter: 32 . migrationSF - The PetscSF 33 34 Level: intermediate 35 36 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexSetMigrationSF 37 @*/ 38 PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) 39 { 40 PetscFunctionBegin; 41 *migrationSF = dm->sfMigration; 42 PetscFunctionReturn(0); 43 } 44 45 /*@ 46 DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec 47 48 Input Parameters: 49 + dm - The DM 50 - naturalSF - The PetscSF 51 52 Level: intermediate 53 54 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexGetGlobaltoNaturalSF() 55 @*/ 56 PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) 57 { 58 PetscErrorCode ierr; 59 PetscFunctionBegin; 60 dm->sfNatural = naturalSF; 61 ierr = PetscObjectReference((PetscObject) naturalSF);CHKERRQ(ierr); 62 dm->useNatural = PETSC_TRUE; 63 PetscFunctionReturn(0); 64 } 65 66 /*@ 67 DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec 68 69 Input Parameter: 70 . dm - The DM 71 72 Output Parameter: 73 . naturalSF - The PetscSF 74 75 Level: intermediate 76 77 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexSetGlobaltoNaturalSF 78 @*/ 79 PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) 80 { 81 PetscFunctionBegin; 82 *naturalSF = dm->sfNatural; 83 PetscFunctionReturn(0); 84 } 85 86 /*@ 87 DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec 88 89 Input Parameters: 90 + dm - The DM 91 . section - The PetscSection describing the Vec before the mesh was distributed, 92 or NULL if not available 93 - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed 94 95 Output Parameter: 96 . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 97 98 Note: This is not typically called by the user. 99 100 Level: intermediate 101 102 .seealso: DMPlexDistribute(), DMPlexDistributeField() 103 @*/ 104 PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 105 { 106 MPI_Comm comm; 107 Vec gv, tmpVec; 108 PetscSF sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv; 109 PetscSection gSection, sectionDist, gLocSection; 110 PetscInt *spoints, *remoteOffsets; 111 PetscInt ssize, pStart, pEnd, p, globalSize; 112 PetscLayout map; 113 PetscBool destroyFlag = PETSC_FALSE; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 118 if (!sfMigration) { 119 /* If sfMigration is missing, 120 sfNatural cannot be computed and is set to NULL */ 121 *sfNatural = NULL; 122 PetscFunctionReturn(0); 123 } else if (!section) { 124 /* If the sequential section is not provided (NULL), 125 it is reconstructed from the parallel section */ 126 PetscSF sfMigrationInv; 127 PetscSection localSection; 128 129 ierr = DMGetLocalSection(dm, &localSection);CHKERRQ(ierr); 130 ierr = PetscSFCreateInverseSF(sfMigration, &sfMigrationInv);CHKERRQ(ierr); 131 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 132 ierr = PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section);CHKERRQ(ierr); 133 ierr = PetscSFDestroy(&sfMigrationInv);CHKERRQ(ierr); 134 destroyFlag = PETSC_TRUE; 135 } 136 /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr); 137 ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */ 138 /* Create a new section from distributing the original section */ 139 ierr = PetscSectionCreate(comm, §ionDist);CHKERRQ(ierr); 140 ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr); 141 /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr); 142 ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 143 ierr = DMSetLocalSection(dm, sectionDist);CHKERRQ(ierr); 144 /* If a sequential section is provided but no dof is affected, 145 sfNatural cannot be computed and is set to NULL */ 146 ierr = DMCreateGlobalVector(dm, &tmpVec);CHKERRQ(ierr); 147 ierr = VecGetSize(tmpVec, &globalSize);CHKERRQ(ierr); 148 ierr = DMRestoreGlobalVector(dm, &tmpVec);CHKERRQ(ierr); 149 if (!globalSize) { 150 *sfNatural = NULL; 151 if (destroyFlag) {ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr);} 152 PetscFunctionReturn(0); 153 } 154 /* Get a pruned version of migration SF */ 155 ierr = DMGetGlobalSection(dm, &gSection);CHKERRQ(ierr); 156 ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr); 157 for (p = pStart, ssize = 0; p < pEnd; ++p) { 158 PetscInt dof, off; 159 160 ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 161 ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 162 if ((dof > 0) && (off >= 0)) ++ssize; 163 } 164 ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr); 165 for (p = pStart, ssize = 0; p < pEnd; ++p) { 166 PetscInt dof, off; 167 168 ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 169 ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 170 if ((dof > 0) && (off >= 0)) spoints[ssize++] = p; 171 } 172 ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr); 173 ierr = PetscFree(spoints);CHKERRQ(ierr); 174 /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr); 175 ierr = PetscSFView(sfEmbed, 0);CHKERRQ(ierr); */ 176 /* Create the SF for seq to natural */ 177 ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr); 178 ierr = VecGetLayout(gv,&map);CHKERRQ(ierr); 179 /* Note that entries of gv are leaves in sfSeqToNatural, entries of the seq vec are roots */ 180 ierr = PetscSFCreate(comm, &sfSeqToNatural);CHKERRQ(ierr); 181 ierr = PetscSFSetGraphWithPattern(sfSeqToNatural, map, PETSCSF_PATTERN_GATHER);CHKERRQ(ierr); 182 ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr); 183 /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr); 184 ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */ 185 /* Create the SF associated with this section */ 186 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 187 ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr); 188 ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr); 189 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 190 ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); 191 ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr); 192 ierr = PetscSectionDestroy(§ionDist);CHKERRQ(ierr); 193 /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr); 194 ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */ 195 /* Invert the field SF so it's now from distributed to sequential */ 196 ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr); 197 ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr); 198 /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr); 199 ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */ 200 /* Multiply the sfFieldInv with the */ 201 ierr = PetscSFComposeInverse(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr); 202 ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr); 203 /* Clean up */ 204 ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr); 205 ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr); 206 if (destroyFlag) {ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr);} 207 PetscFunctionReturn(0); 208 } 209 210 /*@ 211 DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 212 213 Collective on dm 214 215 Input Parameters: 216 + dm - The distributed DMPlex 217 - gv - The global Vec 218 219 Output Parameters: 220 . nv - Vec in the canonical ordering distributed over all processors associated with gv 221 222 Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 223 224 Level: intermediate 225 226 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd() 227 @*/ 228 PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 229 { 230 const PetscScalar *inarray; 231 PetscScalar *outarray; 232 PetscMPIInt size; 233 PetscErrorCode ierr; 234 235 PetscFunctionBegin; 236 ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 237 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 238 if (dm->sfNatural) { 239 ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 240 ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 241 ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 242 ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 243 ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 244 } else if (size == 1) { 245 ierr = VecCopy(gv, nv);CHKERRQ(ierr); 246 } else if (dm->useNatural) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.\n"); 247 else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n"); 248 ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 249 PetscFunctionReturn(0); 250 } 251 252 /*@ 253 DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 254 255 Collective on dm 256 257 Input Parameters: 258 + dm - The distributed DMPlex 259 - gv - The global Vec 260 261 Output Parameters: 262 . nv - The natural Vec 263 264 Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 265 266 Level: intermediate 267 268 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 269 @*/ 270 PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 271 { 272 const PetscScalar *inarray; 273 PetscScalar *outarray; 274 PetscMPIInt size; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 279 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 280 if (dm->sfNatural) { 281 ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 282 ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 283 ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 284 ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 285 ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 286 } else if (size == 1) { 287 } else if (dm->useNatural) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.\n"); 288 else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n"); 289 ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 290 PetscFunctionReturn(0); 291 } 292 293 /*@ 294 DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 295 296 Collective on dm 297 298 Input Parameters: 299 + dm - The distributed DMPlex 300 - nv - The natural Vec 301 302 Output Parameters: 303 . gv - The global Vec 304 305 Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 306 307 Level: intermediate 308 309 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd() 310 @*/ 311 PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 312 { 313 const PetscScalar *inarray; 314 PetscScalar *outarray; 315 PetscMPIInt size; 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 320 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 321 if (dm->sfNatural) { 322 /* We only have acces to the SF that goes from Global to Natural. 323 Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 324 Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 325 ierr = VecZeroEntries(gv);CHKERRQ(ierr); 326 ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 327 ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 328 ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 329 ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 330 ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 331 } else if (size == 1) { 332 ierr = VecCopy(nv, gv);CHKERRQ(ierr); 333 } else if (dm->useNatural) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.\n"); 334 else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n"); 335 ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 336 PetscFunctionReturn(0); 337 } 338 339 /*@ 340 DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 341 342 Collective on dm 343 344 Input Parameters: 345 + dm - The distributed DMPlex 346 - nv - The natural Vec 347 348 Output Parameters: 349 . gv - The global Vec 350 351 Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 352 353 Level: intermediate 354 355 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 356 @*/ 357 PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 358 { 359 const PetscScalar *inarray; 360 PetscScalar *outarray; 361 PetscMPIInt size; 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 366 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 367 if (dm->sfNatural) { 368 ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 369 ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 370 ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 371 ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 372 ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 373 } else if (size == 1) { 374 } else if (dm->useNatural) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.\n"); 375 else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n"); 376 ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379