1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <../src/sys/utils/hash.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMPlexGetFaces_Internal" 6 /* 7 DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 8 */ 9 PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 10 { 11 const PetscInt *cone = NULL; 12 PetscInt maxConeSize, maxSupportSize, coneSize; 13 PetscErrorCode ierr; 14 15 PetscFunctionBegin; 16 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 17 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 18 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 19 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 20 ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 21 PetscFunctionReturn(0); 22 } 23 24 #undef __FUNCT__ 25 #define __FUNCT__ "DMPlexRestoreFaces_Internal" 26 /* 27 DMPlexRestoreFaces_Internal - Restores the array 28 */ 29 PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 30 { 31 PetscErrorCode ierr; 32 33 PetscFunctionBegin; 34 ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, (void *) faces);CHKERRQ(ierr); 35 PetscFunctionReturn(0); 36 } 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "DMPlexGetRawFaces_Internal" 40 /* 41 DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 42 */ 43 PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 44 { 45 PetscInt *facesTmp; 46 PetscInt maxConeSize, maxSupportSize; 47 PetscErrorCode ierr; 48 49 PetscFunctionBegin; 50 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 51 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 52 if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);} 53 switch (dim) { 54 case 2: 55 switch (coneSize) { 56 case 3: 57 if (faces) { 58 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 59 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 60 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 61 *faces = facesTmp; 62 } 63 if (numFaces) *numFaces = 3; 64 if (faceSize) *faceSize = 2; 65 break; 66 case 4: 67 /* Vertices follow right hand rule */ 68 if (faces) { 69 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 70 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 71 facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 72 facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 73 *faces = facesTmp; 74 } 75 if (numFaces) *numFaces = 4; 76 if (faceSize) *faceSize = 2; 77 break; 78 default: 79 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 80 } 81 break; 82 case 3: 83 switch (coneSize) { 84 case 3: 85 if (faces) { 86 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 87 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 88 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 89 *faces = facesTmp; 90 } 91 if (numFaces) *numFaces = 3; 92 if (faceSize) *faceSize = 2; 93 break; 94 case 4: 95 /* Vertices of first face follow right hand rule and normal points away from last vertex */ 96 if (faces) { 97 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 98 facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 99 facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 100 facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 101 *faces = facesTmp; 102 } 103 if (numFaces) *numFaces = 4; 104 if (faceSize) *faceSize = 3; 105 break; 106 case 8: 107 if (faces) { 108 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 109 facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; 110 facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; 111 facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; 112 facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; 113 facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; 114 *faces = facesTmp; 115 } 116 if (numFaces) *numFaces = 6; 117 if (faceSize) *faceSize = 4; 118 break; 119 default: 120 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 121 } 122 break; 123 default: 124 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 125 } 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "DMPlexInterpolateFaces_Internal" 131 /* This interpolates faces for cells at some stratum */ 132 static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 133 { 134 DMLabel subpointMap; 135 PetscHashIJKL faceTable; 136 PetscInt *pStart, *pEnd; 137 PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 142 /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 143 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 144 if (subpointMap) ++cellDim; 145 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 146 ++depth; 147 ++cellDepth; 148 cellDim -= depth - cellDepth; 149 ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 150 for (d = depth-1; d >= faceDepth; --d) { 151 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 152 } 153 ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 154 pEnd[faceDepth] = pStart[faceDepth]; 155 for (d = faceDepth-1; d >= 0; --d) { 156 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 157 } 158 if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 159 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 160 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 161 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 162 const PetscInt *cellFaces; 163 PetscInt numCellFaces, faceSize, cf; 164 165 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 166 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 167 for (cf = 0; cf < numCellFaces; ++cf) { 168 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 169 PetscHashIJKLKey key; 170 PetscHashIJKLIter missing, iter; 171 172 if (faceSize == 2) { 173 key.i = PetscMin(cellFace[0], cellFace[1]); 174 key.j = PetscMax(cellFace[0], cellFace[1]); 175 key.k = 0; 176 key.l = 0; 177 } else { 178 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 179 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 180 } 181 ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 182 if (missing) {ierr = PetscHashIJKLSet(faceTable, iter, face++);CHKERRQ(ierr);} 183 } 184 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 185 } 186 pEnd[faceDepth] = face; 187 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 188 /* Count new points */ 189 for (d = 0; d <= depth; ++d) { 190 numPoints += pEnd[d]-pStart[d]; 191 } 192 ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 193 /* Set cone sizes */ 194 for (d = 0; d <= depth; ++d) { 195 PetscInt coneSize, p; 196 197 if (d == faceDepth) { 198 for (p = pStart[d]; p < pEnd[d]; ++p) { 199 /* I see no way to do this if we admit faces of different shapes */ 200 ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 201 } 202 } else if (d == cellDepth) { 203 for (p = pStart[d]; p < pEnd[d]; ++p) { 204 /* Number of cell faces may be different from number of cell vertices*/ 205 ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 206 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 207 } 208 } else { 209 for (p = pStart[d]; p < pEnd[d]; ++p) { 210 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 211 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 212 } 213 } 214 } 215 ierr = DMSetUp(idm);CHKERRQ(ierr); 216 /* Get face cones from subsets of cell vertices */ 217 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 218 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 219 for (d = depth; d > cellDepth; --d) { 220 const PetscInt *cone; 221 PetscInt p; 222 223 for (p = pStart[d]; p < pEnd[d]; ++p) { 224 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 225 ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 226 ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 227 ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 228 } 229 } 230 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 231 const PetscInt *cellFaces; 232 PetscInt numCellFaces, faceSize, cf; 233 234 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 235 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 236 for (cf = 0; cf < numCellFaces; ++cf) { 237 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 238 PetscHashIJKLKey key; 239 PetscHashIJKLIter missing, iter; 240 241 if (faceSize == 2) { 242 key.i = PetscMin(cellFace[0], cellFace[1]); 243 key.j = PetscMax(cellFace[0], cellFace[1]); 244 key.k = 0; 245 key.l = 0; 246 } else { 247 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 248 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 249 } 250 ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 251 if (missing) { 252 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 253 ierr = PetscHashIJKLSet(faceTable, iter, face);CHKERRQ(ierr); 254 ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 255 } else { 256 const PetscInt *cone; 257 PetscInt coneSize, ornt, i, j, f; 258 259 ierr = PetscHashIJKLGet(faceTable, iter, &f);CHKERRQ(ierr); 260 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 261 /* Orient face: Do not allow reverse orientation at the first vertex */ 262 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 263 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 264 if (coneSize != faceSize) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSize); 265 /* - First find the initial vertex */ 266 for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 267 /* - Try forward comparison */ 268 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 269 if (j == faceSize) { 270 if ((faceSize == 2) && (i == 1)) ornt = -2; 271 else ornt = i; 272 } else { 273 /* - Try backward comparison */ 274 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 275 if (j == faceSize) { 276 if (i == 0) ornt = -faceSize; 277 else ornt = -i; 278 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 279 } 280 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 281 } 282 } 283 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 284 } 285 if (face != pEnd[faceDepth]) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]); 286 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 287 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 288 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 289 ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 290 ierr = DMPlexStratify(idm);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "DMPlexInterpolate" 296 /*@ 297 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 298 299 Collective on DM 300 301 Input Parameter: 302 . dm - The DMPlex object with only cells and vertices 303 304 Output Parameter: 305 . dmInt - The complete DMPlex object 306 307 Level: intermediate 308 309 .keywords: mesh 310 .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList() 311 @*/ 312 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 313 { 314 DM idm, odm = dm; 315 PetscInt depth, dim, d; 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 320 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 321 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 322 if (dim <= 1) { 323 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 324 idm = dm; 325 } 326 for (d = 1; d < dim; ++d) { 327 /* Create interpolated mesh */ 328 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 329 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 330 ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr); 331 if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);} 332 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 333 odm = idm; 334 } 335 *dmInt = idm; 336 ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNCT__ 341 #define __FUNCT__ "DMPlexCopyCoordinates" 342 /*@ 343 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 344 345 Collective on DM 346 347 Input Parameter: 348 . dmA - The DMPlex object with initial coordinates 349 350 Output Parameter: 351 . dmB - The DMPlex object with copied coordinates 352 353 Level: intermediate 354 355 Note: This is typically used when adding pieces other than vertices to a mesh 356 357 .keywords: mesh 358 .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 359 @*/ 360 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 361 { 362 Vec coordinatesA, coordinatesB; 363 PetscSection coordSectionA, coordSectionB; 364 PetscScalar *coordsA, *coordsB; 365 PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 if (dmA == dmB) PetscFunctionReturn(0); 370 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 371 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 372 if ((vEndA-vStartA) != (vEndB-vStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of vertices in first DM %d != %d in the second DM", vEndA-vStartA, vEndB-vStartB); 373 ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 374 ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 375 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 376 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 377 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 378 ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 379 for (v = vStartB; v < vEndB; ++v) { 380 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 381 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 382 } 383 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 384 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 385 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 386 ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr); 387 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 388 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 389 ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr); 390 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 391 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 392 for (v = 0; v < vEndB-vStartB; ++v) { 393 for (d = 0; d < spaceDim; ++d) { 394 coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 395 } 396 } 397 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 398 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 399 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 400 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 401 PetscFunctionReturn(0); 402 } 403 404 #undef __FUNCT__ 405 #define __FUNCT__ "DMPlexCopyLabels" 406 /*@ 407 DMPlexCopyLabels - Copy labels from one mesh to another with a superset of the points 408 409 Collective on DM 410 411 Input Parameter: 412 . dmA - The DMPlex object with initial labels 413 414 Output Parameter: 415 . dmB - The DMPlex object with copied labels 416 417 Level: intermediate 418 419 Note: This is typically used when interpolating or otherwise adding to a mesh 420 421 .keywords: mesh 422 .seealso: DMCopyCoordinates(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 423 @*/ 424 PetscErrorCode DMPlexCopyLabels(DM dmA, DM dmB) 425 { 426 PetscInt numLabels, l; 427 PetscErrorCode ierr; 428 429 PetscFunctionBegin; 430 if (dmA == dmB) PetscFunctionReturn(0); 431 ierr = DMPlexGetNumLabels(dmA, &numLabels);CHKERRQ(ierr); 432 for (l = 0; l < numLabels; ++l) { 433 DMLabel label, labelNew; 434 const char *name; 435 PetscBool flg; 436 437 ierr = DMPlexGetLabelName(dmA, l, &name);CHKERRQ(ierr); 438 ierr = PetscStrcmp(name, "depth", &flg);CHKERRQ(ierr); 439 if (flg) continue; 440 ierr = DMPlexGetLabel(dmA, name, &label);CHKERRQ(ierr); 441 ierr = DMLabelDuplicate(label, &labelNew);CHKERRQ(ierr); 442 ierr = DMPlexAddLabel(dmB, labelNew);CHKERRQ(ierr); 443 } 444 PetscFunctionReturn(0); 445 } 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "DMPlexUninterpolate" 449 /*@ 450 DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 451 452 Collective on DM 453 454 Input Parameter: 455 . dm - The complete DMPlex object 456 457 Output Parameter: 458 . dmUnint - The DMPlex object with only cells and vertices 459 460 Level: intermediate 461 462 .keywords: mesh 463 .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList() 464 @*/ 465 PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 466 { 467 DM udm; 468 PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 469 PetscErrorCode ierr; 470 471 PetscFunctionBegin; 472 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 473 if (dim <= 1) { 474 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 475 *dmUnint = dm; 476 PetscFunctionReturn(0); 477 } 478 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 479 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 480 ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 481 ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 482 ierr = DMPlexSetDimension(udm, dim);CHKERRQ(ierr); 483 ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 484 for (c = cStart; c < cEnd; ++c) { 485 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 486 487 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 488 for (cl = 0; cl < closureSize*2; cl += 2) { 489 const PetscInt p = closure[cl]; 490 491 if ((p >= vStart) && (p < vEnd)) ++coneSize; 492 } 493 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 494 ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 495 maxConeSize = PetscMax(maxConeSize, coneSize); 496 } 497 ierr = DMSetUp(udm);CHKERRQ(ierr); 498 ierr = PetscMalloc(maxConeSize * sizeof(PetscInt), &cone);CHKERRQ(ierr); 499 for (c = cStart; c < cEnd; ++c) { 500 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 501 502 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 503 for (cl = 0; cl < closureSize*2; cl += 2) { 504 const PetscInt p = closure[cl]; 505 506 if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 507 } 508 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 509 ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 510 } 511 ierr = PetscFree(cone);CHKERRQ(ierr); 512 ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 513 ierr = DMPlexStratify(udm);CHKERRQ(ierr); 514 *dmUnint = udm; 515 PetscFunctionReturn(0); 516 } 517