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 static PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 10 { 11 const PetscInt *cone = NULL; 12 PetscInt *facesTmp; 13 PetscInt maxConeSize, maxSupportSize, coneSize; 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 18 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 19 ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr); 20 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 21 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 22 switch (dim) { 23 case 2: 24 switch (coneSize) { 25 case 3: 26 if (faces) { 27 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 28 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 29 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 30 *faces = facesTmp; 31 } 32 if (numFaces) *numFaces = 3; 33 if (faceSize) *faceSize = 2; 34 break; 35 case 4: 36 /* Vertices follow right hand rule */ 37 if (faces) { 38 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 39 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 40 facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 41 facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 42 *faces = facesTmp; 43 } 44 if (numFaces) *numFaces = 4; 45 if (faceSize) *faceSize = 2; 46 if (faces) *faces = facesTmp; 47 break; 48 default: 49 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 50 } 51 break; 52 case 3: 53 switch (coneSize) { 54 case 3: 55 if (faces) { 56 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 57 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 58 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 59 *faces = facesTmp; 60 } 61 if (numFaces) *numFaces = 3; 62 if (faceSize) *faceSize = 2; 63 if (faces) *faces = facesTmp; 64 break; 65 case 4: 66 /* Vertices of first face follow right hand rule and normal points away from last vertex */ 67 if (faces) { 68 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 69 facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 70 facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 71 facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 72 *faces = facesTmp; 73 } 74 if (numFaces) *numFaces = 4; 75 if (faceSize) *faceSize = 3; 76 if (faces) *faces = facesTmp; 77 break; 78 case 8: 79 if (faces) { 80 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 81 facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; 82 facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; 83 facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; 84 facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; 85 facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; 86 *faces = facesTmp; 87 } 88 if (numFaces) *numFaces = 6; 89 if (faceSize) *faceSize = 4; 90 break; 91 default: 92 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 93 } 94 break; 95 default: 96 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 97 } 98 ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, &facesTmp);CHKERRQ(ierr); 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "DMPlexInterpolateFaces_Internal" 104 /* This interpolates faces for cells at some stratum */ 105 static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 106 { 107 PetscHashIJKL faceTable; 108 PetscInt *pStart, *pEnd; 109 PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 114 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 115 ++depth; 116 ++cellDepth; 117 cellDim -= depth - cellDepth; 118 ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr); 119 for (d = depth-1; d >= faceDepth; --d) { 120 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 121 } 122 ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 123 pEnd[faceDepth] = pStart[faceDepth]; 124 for (d = faceDepth-1; d >= 0; --d) { 125 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 126 } 127 if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 128 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 129 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 130 ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 131 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 132 const PetscInt *cellFaces; 133 PetscInt numCellFaces, faceSize, cf, f; 134 135 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 136 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 137 for (cf = 0; cf < numCellFaces; ++cf) { 138 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 139 PetscHashIJKLKey key; 140 141 if (faceSize == 2) { 142 key.i = PetscMin(cellFace[0], cellFace[1]); 143 key.j = PetscMax(cellFace[0], cellFace[1]); 144 } else { 145 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 146 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 147 } 148 ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 149 if (f < 0) { 150 ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 151 f = face++; 152 } 153 } 154 } 155 pEnd[faceDepth] = face; 156 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 157 /* Count new points */ 158 for (d = 0; d <= depth; ++d) { 159 numPoints += pEnd[d]-pStart[d]; 160 } 161 ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 162 /* Set cone sizes */ 163 for (d = 0; d <= depth; ++d) { 164 PetscInt coneSize, p; 165 166 if (d == faceDepth) { 167 for (p = pStart[d]; p < pEnd[d]; ++p) { 168 /* I see no way to do this if we admit faces of different shapes */ 169 ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 170 } 171 } else if (d == cellDepth) { 172 for (p = pStart[d]; p < pEnd[d]; ++p) { 173 /* Number of cell faces may be different from number of cell vertices*/ 174 ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 175 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 176 } 177 } else { 178 for (p = pStart[d]; p < pEnd[d]; ++p) { 179 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 180 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 181 } 182 } 183 } 184 ierr = DMSetUp(idm);CHKERRQ(ierr); 185 /* Get face cones from subsets of cell vertices */ 186 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 187 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 188 ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 189 for (d = depth; d > cellDepth; --d) { 190 const PetscInt *cone; 191 PetscInt p; 192 193 for (p = pStart[d]; p < pEnd[d]; ++p) { 194 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 195 ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 196 ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 197 ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 198 } 199 } 200 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 201 const PetscInt *cellFaces; 202 PetscInt numCellFaces, faceSize, cf, f; 203 204 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 205 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 206 for (cf = 0; cf < numCellFaces; ++cf) { 207 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 208 PetscHashIJKLKey key; 209 210 if (faceSize == 2) { 211 key.i = PetscMin(cellFace[0], cellFace[1]); 212 key.j = PetscMax(cellFace[0], cellFace[1]); 213 } else { 214 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 215 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 216 } 217 ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 218 if (f < 0) { 219 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 220 ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 221 f = face++; 222 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 223 } else { 224 const PetscInt *cone; 225 PetscInt coneSize, ornt, i, j; 226 227 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 228 /* Orient face: Do not allow reverse orientation at the first vertex */ 229 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 230 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 231 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); 232 /* - First find the initial vertex */ 233 for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 234 /* - Try forward comparison */ 235 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 236 if (j == faceSize) { 237 if ((faceSize == 2) && (i == 1)) ornt = -2; 238 else ornt = i; 239 } else { 240 /* - Try backward comparison */ 241 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 242 if (j == faceSize) { 243 if (i == 0) ornt = -faceSize; 244 else ornt = -(i+1); 245 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 246 } 247 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 248 } 249 } 250 } 251 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]); 252 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 253 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 254 ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 255 ierr = DMPlexStratify(idm);CHKERRQ(ierr); 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "DMPlexInterpolate" 261 /*@ 262 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 263 264 Collective on DM 265 266 Input Parameter: 267 . dmA - The DMPlex object with only cells and vertices 268 269 Output Parameter: 270 . dmB - The complete DMPlex object 271 272 Level: intermediate 273 274 .keywords: mesh 275 .seealso: DMPlexCreateFromCellList() 276 @*/ 277 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 278 { 279 DM idm, odm = dm; 280 PetscInt depth, dim, d; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 285 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 286 for (d = 1; d < dim; ++d) { 287 /* Create interpolated mesh */ 288 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 289 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 290 ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr); 291 if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);} 292 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 293 odm = idm; 294 } 295 *dmInt = idm; 296 PetscFunctionReturn(0); 297 } 298 299 #undef __FUNCT__ 300 #define __FUNCT__ "DMPlexCopyCoordinates" 301 /*@ 302 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 303 304 Collective on DM 305 306 Input Parameter: 307 . dmA - The DMPlex object with initial coordinates 308 309 Output Parameter: 310 . dmB - The DMPlex object with copied coordinates 311 312 Level: intermediate 313 314 Note: This is typically used when adding pieces other than vertices to a mesh 315 316 .keywords: mesh 317 .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 318 @*/ 319 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 320 { 321 Vec coordinatesA, coordinatesB; 322 PetscSection coordSectionA, coordSectionB; 323 PetscScalar *coordsA, *coordsB; 324 PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 329 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 330 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); 331 ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 332 ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 333 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 334 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 335 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 336 ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 337 for (v = vStartB; v < vEndB; ++v) { 338 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 339 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 340 } 341 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 342 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 343 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 344 ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr); 345 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 346 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 347 ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr); 348 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 349 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 350 for (v = 0; v < vEndB-vStartB; ++v) { 351 for (d = 0; d < spaceDim; ++d) { 352 coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 353 } 354 } 355 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 356 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 357 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 358 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361