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 DMLabel subpointMap; 108 PetscHashIJKL faceTable; 109 PetscInt *pStart, *pEnd; 110 PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 111 PetscErrorCode ierr; 112 113 PetscFunctionBegin; 114 ierr = DMPlexGetDimension(dm, &cellDim);CHKERRQ(ierr); 115 /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 116 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 117 if (subpointMap) ++cellDim; 118 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 119 ++depth; 120 ++cellDepth; 121 cellDim -= depth - cellDepth; 122 ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr); 123 for (d = depth-1; d >= faceDepth; --d) { 124 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 125 } 126 ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 127 pEnd[faceDepth] = pStart[faceDepth]; 128 for (d = faceDepth-1; d >= 0; --d) { 129 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 130 } 131 if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 132 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 133 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 134 ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 135 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 136 const PetscInt *cellFaces; 137 PetscInt numCellFaces, faceSize, cf, f; 138 139 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 140 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 141 for (cf = 0; cf < numCellFaces; ++cf) { 142 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 143 PetscHashIJKLKey key; 144 145 if (faceSize == 2) { 146 key.i = PetscMin(cellFace[0], cellFace[1]); 147 key.j = PetscMax(cellFace[0], cellFace[1]); 148 } else { 149 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 150 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 151 } 152 ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 153 if (f < 0) { 154 ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 155 f = face++; 156 } 157 } 158 } 159 pEnd[faceDepth] = face; 160 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 161 /* Count new points */ 162 for (d = 0; d <= depth; ++d) { 163 numPoints += pEnd[d]-pStart[d]; 164 } 165 ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 166 /* Set cone sizes */ 167 for (d = 0; d <= depth; ++d) { 168 PetscInt coneSize, p; 169 170 if (d == faceDepth) { 171 for (p = pStart[d]; p < pEnd[d]; ++p) { 172 /* I see no way to do this if we admit faces of different shapes */ 173 ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 174 } 175 } else if (d == cellDepth) { 176 for (p = pStart[d]; p < pEnd[d]; ++p) { 177 /* Number of cell faces may be different from number of cell vertices*/ 178 ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 179 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 180 } 181 } else { 182 for (p = pStart[d]; p < pEnd[d]; ++p) { 183 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 184 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 185 } 186 } 187 } 188 ierr = DMSetUp(idm);CHKERRQ(ierr); 189 /* Get face cones from subsets of cell vertices */ 190 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 191 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 192 ierr = PetscHashIJKLSetMultivalued(faceTable, PETSC_FALSE);CHKERRQ(ierr); 193 for (d = depth; d > cellDepth; --d) { 194 const PetscInt *cone; 195 PetscInt p; 196 197 for (p = pStart[d]; p < pEnd[d]; ++p) { 198 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 199 ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 200 ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 201 ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 202 } 203 } 204 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 205 const PetscInt *cellFaces; 206 PetscInt numCellFaces, faceSize, cf, f; 207 208 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 209 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 210 for (cf = 0; cf < numCellFaces; ++cf) { 211 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 212 PetscHashIJKLKey key; 213 214 if (faceSize == 2) { 215 key.i = PetscMin(cellFace[0], cellFace[1]); 216 key.j = PetscMax(cellFace[0], cellFace[1]); 217 } else { 218 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 219 ierr = PetscSortInt(faceSize, (PetscInt *) &key); 220 } 221 ierr = PetscHashIJKLGet(faceTable, key, &f);CHKERRQ(ierr); 222 if (f < 0) { 223 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 224 ierr = PetscHashIJKLAdd(faceTable, key, face);CHKERRQ(ierr); 225 f = face++; 226 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 227 } else { 228 const PetscInt *cone; 229 PetscInt coneSize, ornt, i, j; 230 231 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 232 /* Orient face: Do not allow reverse orientation at the first vertex */ 233 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 234 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 235 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); 236 /* - First find the initial vertex */ 237 for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 238 /* - Try forward comparison */ 239 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 240 if (j == faceSize) { 241 if ((faceSize == 2) && (i == 1)) ornt = -2; 242 else ornt = i; 243 } else { 244 /* - Try backward comparison */ 245 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 246 if (j == faceSize) { 247 if (i == 0) ornt = -faceSize; 248 else ornt = -(i+1); 249 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 250 } 251 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 252 } 253 } 254 } 255 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]); 256 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 257 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 258 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 259 ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 260 ierr = DMPlexStratify(idm);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "DMPlexInterpolate" 266 /*@ 267 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 268 269 Collective on DM 270 271 Input Parameter: 272 . dmA - The DMPlex object with only cells and vertices 273 274 Output Parameter: 275 . dmB - The complete DMPlex object 276 277 Level: intermediate 278 279 .keywords: mesh 280 .seealso: DMPlexCreateFromCellList() 281 @*/ 282 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 283 { 284 DM idm, odm = dm; 285 PetscInt depth, dim, d; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 290 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 291 if (dim <= 1) { 292 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 293 idm = dm; 294 } 295 for (d = 1; d < dim; ++d) { 296 /* Create interpolated mesh */ 297 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 298 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 299 ierr = DMPlexSetDimension(idm, dim);CHKERRQ(ierr); 300 if (depth > 0) {ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr);} 301 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 302 odm = idm; 303 } 304 *dmInt = idm; 305 PetscFunctionReturn(0); 306 } 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "DMPlexCopyCoordinates" 310 /*@ 311 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 312 313 Collective on DM 314 315 Input Parameter: 316 . dmA - The DMPlex object with initial coordinates 317 318 Output Parameter: 319 . dmB - The DMPlex object with copied coordinates 320 321 Level: intermediate 322 323 Note: This is typically used when adding pieces other than vertices to a mesh 324 325 .keywords: mesh 326 .seealso: DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMPlexGetCoordinateSection() 327 @*/ 328 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 329 { 330 Vec coordinatesA, coordinatesB; 331 PetscSection coordSectionA, coordSectionB; 332 PetscScalar *coordsA, *coordsB; 333 PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 334 PetscErrorCode ierr; 335 336 PetscFunctionBegin; 337 if (dmA == dmB) PetscFunctionReturn(0); 338 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 339 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 340 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); 341 ierr = DMPlexGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 342 ierr = DMPlexGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 343 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 344 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 345 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 346 ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 347 for (v = vStartB; v < vEndB; ++v) { 348 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 349 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 350 } 351 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 352 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 353 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 354 ierr = VecCreate(PetscObjectComm((PetscObject) dmB), &coordinatesB);CHKERRQ(ierr); 355 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 356 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 357 ierr = VecSetFromOptions(coordinatesB);CHKERRQ(ierr); 358 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 359 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 360 for (v = 0; v < vEndB-vStartB; ++v) { 361 for (d = 0; d < spaceDim; ++d) { 362 coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 363 } 364 } 365 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 366 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 367 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 368 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371