1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/petscfeimpl.h> /* For PetscFEInterpolate_Static() */ 3 4 #include <petscdmplextransform.h> /*I "petscdmplextransform.h" I*/ 5 #include <petscsf.h> 6 7 /*@ 8 DMPlexCreateProcessSF - Create an `PetscSF` which just has process connectivity 9 10 Collective 11 12 Input Parameters: 13 + dm - The `DM` 14 - sfPoint - The `PetscSF` which encodes point connectivity 15 16 Output Parameters: 17 + processRanks - A list of process neighbors, or `NULL` 18 - sfProcess - An `PetscSF` encoding the process connectivity, or `NULL` 19 20 Level: developer 21 22 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSFCreate()`, `DMPlexCreateTwoSidedProcessSF()` 23 @*/ 24 PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess) 25 { 26 PetscInt numRoots, numLeaves, l; 27 const PetscInt *localPoints; 28 const PetscSFNode *remotePoints; 29 PetscInt *localPointsNew; 30 PetscSFNode *remotePointsNew; 31 PetscMPIInt *ranks; 32 PetscInt *ranksNew; 33 PetscMPIInt size; 34 35 PetscFunctionBegin; 36 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 37 PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 38 if (processRanks) PetscAssertPointer(processRanks, 3); 39 if (sfProcess) PetscAssertPointer(sfProcess, 4); 40 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 41 PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints)); 42 PetscCall(PetscMalloc1(numLeaves, &ranks)); 43 for (l = 0; l < numLeaves; ++l) ranks[l] = (PetscMPIInt)remotePoints[l].rank; 44 PetscCall(PetscSortRemoveDupsMPIInt(&numLeaves, ranks)); 45 PetscCall(PetscMalloc1(numLeaves, &ranksNew)); 46 PetscCall(PetscMalloc1(numLeaves, &localPointsNew)); 47 PetscCall(PetscMalloc1(numLeaves, &remotePointsNew)); 48 for (l = 0; l < numLeaves; ++l) { 49 ranksNew[l] = ranks[l]; 50 localPointsNew[l] = l; 51 remotePointsNew[l].index = 0; 52 remotePointsNew[l].rank = ranksNew[l]; 53 } 54 PetscCall(PetscFree(ranks)); 55 if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks)); 56 else PetscCall(PetscFree(ranksNew)); 57 if (sfProcess) { 58 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess)); 59 PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF")); 60 PetscCall(PetscSFSetFromOptions(*sfProcess)); 61 PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER)); 62 } 63 PetscFunctionReturn(PETSC_SUCCESS); 64 } 65 66 /*@ 67 DMPlexCreateCoarsePointIS - Creates an `IS` covering the coarse `DM` chart with the fine points as data 68 69 Collective 70 71 Input Parameter: 72 . dm - The coarse `DM` 73 74 Output Parameter: 75 . fpointIS - The `IS` of all the fine points which exist in the original coarse mesh 76 77 Level: developer 78 79 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `IS`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()` 80 @*/ 81 PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS) 82 { 83 DMPlexTransform tr; 84 PetscInt *fpoints; 85 PetscInt pStart, pEnd, p, vStart, vEnd, v; 86 87 PetscFunctionBegin; 88 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 89 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 90 PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 91 PetscCall(DMPlexTransformSetUp(tr)); 92 PetscCall(PetscMalloc1(pEnd - pStart, &fpoints)); 93 for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1; 94 for (v = vStart; v < vEnd; ++v) { 95 PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */ 96 97 PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew)); 98 fpoints[v - pStart] = vNew; 99 } 100 PetscCall(DMPlexTransformDestroy(&tr)); 101 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS)); 102 PetscFunctionReturn(PETSC_SUCCESS); 103 } 104 105 /*@ 106 DMPlexSetTransformType - Set the transform type for uniform refinement 107 108 Input Parameters: 109 + dm - The `DM` 110 - type - The transform type for uniform refinement 111 112 Level: developer 113 114 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()` 115 @*/ 116 PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type) 117 { 118 DM_Plex *mesh = (DM_Plex *)dm->data; 119 120 PetscFunctionBegin; 121 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 122 if (type) PetscAssertPointer(type, 2); 123 PetscCall(PetscFree(mesh->transformType)); 124 PetscCall(PetscStrallocpy(type, &mesh->transformType)); 125 PetscFunctionReturn(PETSC_SUCCESS); 126 } 127 128 /*@C 129 DMPlexGetTransformType - Retrieve the transform type for uniform refinement 130 131 Input Parameter: 132 . dm - The `DM` 133 134 Output Parameter: 135 . type - The transform type for uniform refinement 136 137 Level: developer 138 139 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()` 140 @*/ 141 PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type) 142 { 143 DM_Plex *mesh = (DM_Plex *)dm->data; 144 145 PetscFunctionBegin; 146 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 147 PetscAssertPointer(type, 2); 148 *type = mesh->transformType; 149 PetscFunctionReturn(PETSC_SUCCESS); 150 } 151 152 /*@ 153 DMPlexSetRefinementUniform - Set the flag for uniform refinement 154 155 Input Parameters: 156 + dm - The `DM` 157 - refinementUniform - The flag for uniform refinement 158 159 Level: developer 160 161 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 162 @*/ 163 PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform) 164 { 165 DM_Plex *mesh = (DM_Plex *)dm->data; 166 167 PetscFunctionBegin; 168 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 169 mesh->refinementUniform = refinementUniform; 170 PetscFunctionReturn(PETSC_SUCCESS); 171 } 172 173 /*@ 174 DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement 175 176 Input Parameter: 177 . dm - The `DM` 178 179 Output Parameter: 180 . refinementUniform - The flag for uniform refinement 181 182 Level: developer 183 184 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 185 @*/ 186 PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform) 187 { 188 DM_Plex *mesh = (DM_Plex *)dm->data; 189 190 PetscFunctionBegin; 191 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 192 PetscAssertPointer(refinementUniform, 2); 193 *refinementUniform = mesh->refinementUniform; 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 /*@ 198 DMPlexSetRefinementLimit - Set the maximum cell volume for refinement 199 200 Input Parameters: 201 + dm - The `DM` 202 - refinementLimit - The maximum cell volume in the refined mesh 203 204 Level: developer 205 206 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()` 207 @*/ 208 PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit) 209 { 210 DM_Plex *mesh = (DM_Plex *)dm->data; 211 212 PetscFunctionBegin; 213 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 214 mesh->refinementLimit = refinementLimit; 215 PetscFunctionReturn(PETSC_SUCCESS); 216 } 217 218 /*@ 219 DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement 220 221 Input Parameter: 222 . dm - The `DM` 223 224 Output Parameter: 225 . refinementLimit - The maximum cell volume in the refined mesh 226 227 Level: developer 228 229 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()` 230 @*/ 231 PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit) 232 { 233 DM_Plex *mesh = (DM_Plex *)dm->data; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 237 PetscAssertPointer(refinementLimit, 2); 238 /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */ 239 *refinementLimit = mesh->refinementLimit; 240 PetscFunctionReturn(PETSC_SUCCESS); 241 } 242 243 /*@C 244 DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement 245 246 Input Parameters: 247 + dm - The `DM` 248 - refinementFunc - Function giving the maximum cell volume in the refined mesh 249 250 Calling Sequence of `refinementFunc`: 251 + coords - Coordinates of the current point, usually a cell centroid 252 - limit - The maximum cell volume for a cell containing this point 253 254 Level: developer 255 256 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 257 @*/ 258 PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal coords[], PetscReal *limit)) 259 { 260 DM_Plex *mesh = (DM_Plex *)dm->data; 261 262 PetscFunctionBegin; 263 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 264 mesh->refinementFunc = refinementFunc; 265 PetscFunctionReturn(PETSC_SUCCESS); 266 } 267 268 /*@C 269 DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement 270 271 Input Parameter: 272 . dm - The `DM` 273 274 Output Parameter: 275 . refinementFunc - Function giving the maximum cell volume in the refined mesh 276 277 Calling Sequence of `refinementFunc`: 278 + coords - Coordinates of the current point, usually a cell centroid 279 - limit - The maximum cell volume for a cell containing this point 280 281 Level: developer 282 283 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 284 @*/ 285 PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal coords[], PetscReal *limit)) 286 { 287 DM_Plex *mesh = (DM_Plex *)dm->data; 288 289 PetscFunctionBegin; 290 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 291 PetscAssertPointer(refinementFunc, 2); 292 *refinementFunc = mesh->refinementFunc; 293 PetscFunctionReturn(PETSC_SUCCESS); 294 } 295 296 PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm) 297 { 298 PetscBool isUniform; 299 300 PetscFunctionBegin; 301 PetscCall(DMPlexGetRefinementUniform(dm, &isUniform)); 302 PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view")); 303 if (isUniform) { 304 DMPlexTransform tr; 305 DM cdm, rcdm; 306 DMPlexTransformType trType; 307 const char *prefix; 308 PetscOptions options; 309 PetscInt cDegree; 310 PetscBool useCeed; 311 312 PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 313 PetscCall(DMPlexTransformSetDM(tr, dm)); 314 PetscCall(DMPlexGetTransformType(dm, &trType)); 315 if (trType) PetscCall(DMPlexTransformSetType(tr, trType)); 316 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 317 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 318 PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 319 PetscCall(PetscObjectSetOptions((PetscObject)tr, options)); 320 PetscCall(DMPlexTransformSetFromOptions(tr)); 321 PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL)); 322 PetscCall(DMPlexTransformSetUp(tr)); 323 PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view")); 324 PetscCall(DMPlexTransformApply(tr, dm, rdm)); 325 PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE)); 326 PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 327 PetscCall(DMPlexSetUseCeed(*rdm, useCeed)); 328 PetscCall(DMSetMatType(*rdm, dm->mattype)); 329 PetscCall(DMCopyDisc(dm, *rdm)); 330 PetscCall(DMGetCoordinateDM(dm, &cdm)); 331 PetscCall(DMGetCoordinateDM(*rdm, &rcdm)); 332 PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree)); 333 { 334 PetscDS cds, rcds; 335 336 PetscCall(DMPlexCreateCoordinateSpace(*rdm, cDegree, PETSC_TRUE, NULL)); 337 PetscCall(DMGetCoordinateDM(*rdm, &rcdm)); 338 PetscCall(DMGetDS(cdm, &cds)); 339 PetscCall(DMGetDS(rcdm, &rcds)); 340 PetscCall(PetscDSCopyConstants(cds, rcds)); 341 } 342 PetscCall(DMPlexGetUseCeed(cdm, &useCeed)); 343 PetscCall(DMPlexSetUseCeed(rcdm, useCeed)); 344 if (useCeed) { 345 PetscCall(DMPlexSetUseMatClosurePermutation(rcdm, PETSC_FALSE)); 346 PetscCall(DMUseTensorOrder(rcdm, PETSC_TRUE)); 347 } 348 PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm)); 349 PetscCall(DMPlexTransformDestroy(&tr)); 350 } else { 351 PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm)); 352 } 353 if (*rdm) { 354 ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 355 ((DM_Plex *)(*rdm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 356 } 357 PetscCall(DMViewFromOptions(*rdm, NULL, "-postref_dm_view")); 358 PetscFunctionReturn(PETSC_SUCCESS); 359 } 360 361 PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[]) 362 { 363 DM cdm = dm; 364 PetscInt r; 365 PetscBool isUniform, localized, useCeed; 366 367 PetscFunctionBegin; 368 PetscCall(DMPlexGetRefinementUniform(dm, &isUniform)); 369 PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 370 if (isUniform) { 371 for (r = 0; r < nlevels; ++r) { 372 DMPlexTransform tr; 373 DM codm, rcodm; 374 const char *prefix; 375 376 PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr)); 377 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix)); 378 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 379 PetscCall(DMPlexTransformSetDM(tr, cdm)); 380 PetscCall(DMPlexTransformSetFromOptions(tr)); 381 PetscCall(DMPlexTransformSetUp(tr)); 382 PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r])); 383 PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown)); 384 PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1)); 385 PetscCall(DMSetMatType(rdm[r], dm->mattype)); 386 PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 387 PetscCall(DMPlexSetUseCeed(rdm[r], useCeed)); 388 PetscCall(DMCopyDisc(cdm, rdm[r])); 389 PetscCall(DMGetCoordinateDM(dm, &codm)); 390 PetscCall(DMGetCoordinateDM(rdm[r], &rcodm)); 391 PetscCall(DMCopyDisc(codm, rcodm)); 392 PetscCall(DMPlexGetUseCeed(codm, &useCeed)); 393 PetscCall(DMPlexSetUseCeed(rcodm, useCeed)); 394 if (useCeed) { 395 PetscCall(DMPlexSetUseMatClosurePermutation(rcodm, PETSC_FALSE)); 396 PetscCall(DMUseTensorOrder(rcodm, PETSC_TRUE)); 397 } 398 PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r])); 399 PetscCall(DMSetCoarseDM(rdm[r], cdm)); 400 PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE)); 401 if (rdm[r]) { 402 ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 403 ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 404 } 405 cdm = rdm[r]; 406 PetscCall(DMPlexTransformDestroy(&tr)); 407 } 408 } else { 409 for (r = 0; r < nlevels; ++r) { 410 PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r])); 411 PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 412 PetscCall(DMPlexSetUseCeed(rdm[r], useCeed)); 413 PetscCall(DMCopyDisc(cdm, rdm[r])); 414 if (localized) PetscCall(DMLocalizeCoordinates(rdm[r])); 415 PetscCall(DMSetCoarseDM(rdm[r], cdm)); 416 if (rdm[r]) { 417 ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 418 ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 419 } 420 cdm = rdm[r]; 421 } 422 } 423 PetscFunctionReturn(PETSC_SUCCESS); 424 } 425