1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <../src/sys/utils/hash.h> 3 4 /* 5 DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 6 */ 7 PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 8 { 9 const PetscInt *cone = NULL; 10 PetscInt maxConeSize, maxSupportSize, coneSize; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 16 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 17 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 18 ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 /* 23 DMPlexRestoreFaces_Internal - Restores the array 24 */ 25 PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 26 { 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, (void *) faces);CHKERRQ(ierr); 31 PetscFunctionReturn(0); 32 } 33 34 /* 35 DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 36 */ 37 PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 38 { 39 PetscInt *facesTmp; 40 PetscInt maxConeSize, maxSupportSize; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 45 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 46 if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), PETSC_INT, &facesTmp);CHKERRQ(ierr);} 47 switch (dim) { 48 case 1: 49 switch (coneSize) { 50 case 2: 51 if (faces) { 52 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 53 *faces = facesTmp; 54 } 55 if (numFaces) *numFaces = 2; 56 if (faceSize) *faceSize = 1; 57 break; 58 default: 59 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 60 } 61 break; 62 case 2: 63 switch (coneSize) { 64 case 3: 65 if (faces) { 66 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 67 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 68 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 69 *faces = facesTmp; 70 } 71 if (numFaces) *numFaces = 3; 72 if (faceSize) *faceSize = 2; 73 break; 74 case 4: 75 /* Vertices follow right hand rule */ 76 if (faces) { 77 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 78 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 79 facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 80 facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 81 *faces = facesTmp; 82 } 83 if (numFaces) *numFaces = 4; 84 if (faceSize) *faceSize = 2; 85 break; 86 default: 87 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 88 } 89 break; 90 case 3: 91 switch (coneSize) { 92 case 3: 93 if (faces) { 94 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 95 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 96 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 97 *faces = facesTmp; 98 } 99 if (numFaces) *numFaces = 3; 100 if (faceSize) *faceSize = 2; 101 break; 102 case 4: 103 /* Vertices of first face follow right hand rule and normal points away from last vertex */ 104 if (faces) { 105 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 106 facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 107 facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 108 facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 109 *faces = facesTmp; 110 } 111 if (numFaces) *numFaces = 4; 112 if (faceSize) *faceSize = 3; 113 break; 114 case 8: 115 if (faces) { 116 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 117 facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 118 facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 119 facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 120 facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 121 facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 122 *faces = facesTmp; 123 } 124 if (numFaces) *numFaces = 6; 125 if (faceSize) *faceSize = 4; 126 break; 127 default: 128 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 129 } 130 break; 131 default: 132 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 133 } 134 PetscFunctionReturn(0); 135 } 136 137 /* This interpolates faces for cells at some stratum */ 138 static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 139 { 140 DMLabel subpointMap; 141 PetscHashIJKL faceTable; 142 PetscInt *pStart, *pEnd; 143 PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr); 148 /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 149 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 150 if (subpointMap) ++cellDim; 151 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 152 ++depth; 153 ++cellDepth; 154 cellDim -= depth - cellDepth; 155 ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 156 for (d = depth-1; d >= faceDepth; --d) { 157 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 158 } 159 ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 160 pEnd[faceDepth] = pStart[faceDepth]; 161 for (d = faceDepth-1; d >= 0; --d) { 162 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 163 } 164 if (pEnd[cellDepth] > pStart[cellDepth]) {ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr);} 165 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 166 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 167 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 168 const PetscInt *cellFaces; 169 PetscInt numCellFaces, faceSize, cf; 170 171 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 172 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 173 for (cf = 0; cf < numCellFaces; ++cf) { 174 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 175 PetscHashIJKLKey key; 176 PetscHashIJKLIter missing, iter; 177 178 if (faceSize == 2) { 179 key.i = PetscMin(cellFace[0], cellFace[1]); 180 key.j = PetscMax(cellFace[0], cellFace[1]); 181 key.k = 0; 182 key.l = 0; 183 } else { 184 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 185 ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 186 } 187 ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 188 if (missing) {ierr = PetscHashIJKLSet(faceTable, iter, face++);CHKERRQ(ierr);} 189 } 190 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 191 } 192 pEnd[faceDepth] = face; 193 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 194 /* Count new points */ 195 for (d = 0; d <= depth; ++d) { 196 numPoints += pEnd[d]-pStart[d]; 197 } 198 ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 199 /* Set cone sizes */ 200 for (d = 0; d <= depth; ++d) { 201 PetscInt coneSize, p; 202 203 if (d == faceDepth) { 204 for (p = pStart[d]; p < pEnd[d]; ++p) { 205 /* I see no way to do this if we admit faces of different shapes */ 206 ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 207 } 208 } else if (d == cellDepth) { 209 for (p = pStart[d]; p < pEnd[d]; ++p) { 210 /* Number of cell faces may be different from number of cell vertices*/ 211 ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 212 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 213 } 214 } else { 215 for (p = pStart[d]; p < pEnd[d]; ++p) { 216 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 217 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 218 } 219 } 220 } 221 ierr = DMSetUp(idm);CHKERRQ(ierr); 222 /* Get face cones from subsets of cell vertices */ 223 if (faceSizeAll > 4) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 224 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 225 for (d = depth; d > cellDepth; --d) { 226 const PetscInt *cone; 227 PetscInt p; 228 229 for (p = pStart[d]; p < pEnd[d]; ++p) { 230 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 231 ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 232 ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 233 ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 234 } 235 } 236 for (c = pStart[cellDepth], face = pStart[faceDepth]; c < pEnd[cellDepth]; ++c) { 237 const PetscInt *cellFaces; 238 PetscInt numCellFaces, faceSize, cf; 239 240 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 241 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 242 for (cf = 0; cf < numCellFaces; ++cf) { 243 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 244 PetscHashIJKLKey key; 245 PetscHashIJKLIter missing, iter; 246 247 if (faceSize == 2) { 248 key.i = PetscMin(cellFace[0], cellFace[1]); 249 key.j = PetscMax(cellFace[0], cellFace[1]); 250 key.k = 0; 251 key.l = 0; 252 } else { 253 key.i = cellFace[0]; key.j = cellFace[1]; key.k = cellFace[2]; key.l = faceSize > 3 ? cellFace[3] : 0; 254 ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 255 } 256 ierr = PetscHashIJKLPut(faceTable, key, &missing, &iter);CHKERRQ(ierr); 257 if (missing) { 258 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 259 ierr = PetscHashIJKLSet(faceTable, iter, face);CHKERRQ(ierr); 260 ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 261 } else { 262 const PetscInt *cone; 263 PetscInt coneSize, ornt, i, j, f; 264 265 ierr = PetscHashIJKLGet(faceTable, iter, &f);CHKERRQ(ierr); 266 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 267 /* Orient face: Do not allow reverse orientation at the first vertex */ 268 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 269 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 270 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); 271 /* - First find the initial vertex */ 272 for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 273 /* - Try forward comparison */ 274 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 275 if (j == faceSize) { 276 if ((faceSize == 2) && (i == 1)) ornt = -2; 277 else ornt = i; 278 } else { 279 /* - Try backward comparison */ 280 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 281 if (j == faceSize) { 282 if (i == 0) ornt = -faceSize; 283 else ornt = -i; 284 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 285 } 286 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 287 } 288 } 289 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 290 } 291 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]); 292 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 293 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 294 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 295 ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 296 ierr = DMPlexStratify(idm);CHKERRQ(ierr); 297 PetscFunctionReturn(0); 298 } 299 300 /* This interpolates the PointSF in parallel following local interpolation */ 301 static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth) 302 { 303 PetscMPIInt size, rank; 304 PetscInt p, c, d, dof, offset; 305 PetscInt numLeaves, numRoots, candidatesSize, candidatesRemoteSize; 306 const PetscInt *localPoints; 307 const PetscSFNode *remotePoints; 308 PetscSFNode *candidates, *candidatesRemote, *claims; 309 PetscSection candidateSection, candidateSectionRemote, claimSection; 310 PetscHashI leafhash; 311 PetscHashIJ roothash; 312 PetscHashIJKey key; 313 PetscErrorCode ierr; 314 315 PetscFunctionBegin; 316 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 317 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 318 ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 319 if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 320 ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 321 /* Build hashes of points in the SF for efficient lookup */ 322 PetscHashICreate(leafhash); 323 PetscHashIJCreate(&roothash); 324 ierr = PetscHashIJSetMultivalued(roothash, PETSC_FALSE);CHKERRQ(ierr); 325 for (p = 0; p < numLeaves; ++p) { 326 PetscHashIAdd(leafhash, localPoints[p], p); 327 key.i = remotePoints[p].index; key.j = remotePoints[p].rank; 328 PetscHashIJAdd(roothash, key, p); 329 } 330 /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves, 331 where each candidate is defined by the root entry for the other vertex that defines the edge. */ 332 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 333 ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 334 { 335 PetscInt leaf, root, idx, a, *adj = NULL; 336 for (p = 0; p < numLeaves; ++p) { 337 PetscInt adjSize = PETSC_DETERMINE; 338 ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 339 for (a = 0; a < adjSize; ++a) { 340 PetscHashIMap(leafhash, adj[a], leaf); 341 if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);} 342 } 343 } 344 ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 345 ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 346 ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 347 for (p = 0; p < numLeaves; ++p) { 348 PetscInt adjSize = PETSC_DETERMINE; 349 ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr); 350 ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 351 for (idx = 0, a = 0; a < adjSize; ++a) { 352 PetscHashIMap(leafhash, adj[a], root); 353 if (root >= 0) candidates[offset+idx++] = remotePoints[root]; 354 } 355 } 356 ierr = PetscFree(adj);CHKERRQ(ierr); 357 } 358 /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 359 { 360 PetscSF sfMulti, sfInverse, sfCandidates; 361 PetscInt *remoteOffsets; 362 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 363 ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 364 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 365 ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 366 ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 367 ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 368 ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 369 ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 370 ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 371 ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 372 ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 373 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 374 } 375 /* Walk local roots and check for each remote candidate whether we know all required points, 376 either from owning it or having a root entry in the point SF. If we do we place a claim 377 by replacing the vertex number with our edge ID. */ 378 { 379 PetscInt idx, root, joinSize, vertices[2]; 380 const PetscInt *rootdegree, *join = NULL; 381 ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 382 ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 383 /* Loop remote edge connections and put in a claim if both vertices are known */ 384 for (idx = 0, p = 0; p < numRoots; ++p) { 385 for (d = 0; d < rootdegree[p]; ++d) { 386 ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 387 ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 388 for (c = 0; c < dof; ++c) { 389 /* We own both vertices, so we claim the edge by replacing vertex with edge */ 390 if (candidatesRemote[offset+c].rank == rank) { 391 vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index; 392 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 393 if (joinSize == 1) candidatesRemote[offset+c].index = join[0]; 394 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 395 continue; 396 } 397 /* If we own one vertex and share a root with the other, we claim it */ 398 key.i = candidatesRemote[offset+c].index; key.j = candidatesRemote[offset+c].rank; 399 PetscHashIJGet(roothash, key, &root); 400 if (root >= 0) { 401 vertices[0] = p; vertices[1] = localPoints[root]; 402 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 403 if (joinSize == 1) { 404 candidatesRemote[offset+c].index = join[0]; 405 candidatesRemote[offset+c].rank = rank; 406 } 407 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 408 } 409 } 410 idx++; 411 } 412 } 413 } 414 /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 415 { 416 PetscSF sfMulti, sfClaims, sfPointNew; 417 PetscHashI claimshash; 418 PetscInt size, pStart, pEnd, root, joinSize, numLocalNew; 419 PetscInt *remoteOffsets, *localPointsNew, vertices[2]; 420 const PetscInt *join = NULL; 421 PetscSFNode *remotePointsNew; 422 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 423 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 424 ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 425 ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 426 ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr); 427 ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr); 428 ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 429 ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 430 ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 431 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 432 /* Walk the original section of local supports and add an SF entry for each updated item */ 433 PetscHashICreate(claimshash); 434 for (p = 0; p < numRoots; ++p) { 435 ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 436 ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 437 for (d = 0; d < dof; ++d) { 438 if (candidates[offset+d].index != claims[offset+d].index) { 439 key.i = candidates[offset+d].index; key.j = candidates[offset+d].rank; 440 PetscHashIJGet(roothash, key, &root); 441 if (root >= 0) { 442 vertices[0] = p; vertices[1] = localPoints[root]; 443 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 444 if (joinSize == 1) PetscHashIAdd(claimshash, join[0], offset+d); 445 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 446 } 447 } 448 } 449 } 450 /* Create new pointSF from hashed claims */ 451 PetscHashISize(claimshash, numLocalNew); 452 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 453 ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 454 ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 455 for (p = 0; p < numLeaves; ++p) { 456 localPointsNew[p] = localPoints[p]; 457 remotePointsNew[p].index = remotePoints[p].index; 458 remotePointsNew[p].rank = remotePoints[p].rank; 459 } 460 p = numLeaves; ierr = PetscHashIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 461 for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 462 PetscHashIMap(claimshash, localPointsNew[p], offset); 463 remotePointsNew[p] = claims[offset]; 464 } 465 ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 466 ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 467 ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 468 ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 469 PetscHashIDestroy(claimshash); 470 } 471 PetscHashIDestroy(leafhash); 472 ierr = PetscHashIJDestroy(&roothash);CHKERRQ(ierr); 473 ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 474 ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 475 ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 476 ierr = PetscFree(candidates);CHKERRQ(ierr); 477 ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 478 ierr = PetscFree(claims);CHKERRQ(ierr); 479 ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 480 PetscFunctionReturn(0); 481 } 482 483 /*@C 484 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 485 486 Collective on DM 487 488 Input Parameters: 489 + dm - The DMPlex object with only cells and vertices 490 - dmInt - If NULL a new DM is created, otherwise the interpolated DM is put into the given DM 491 492 Output Parameter: 493 . dmInt - The complete DMPlex object 494 495 Level: intermediate 496 497 .keywords: mesh 498 .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList() 499 @*/ 500 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 501 { 502 DM idm, odm = dm; 503 PetscSF sfPoint; 504 PetscInt depth, dim, d; 505 PetscErrorCode ierr; 506 507 PetscFunctionBegin; 508 ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 509 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 510 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 511 if (dim <= 1) { 512 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 513 idm = dm; 514 } 515 for (d = 1; d < dim; ++d) { 516 /* Create interpolated mesh */ 517 if ((d == dim-1) && *dmInt) {idm = *dmInt;} 518 else {ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr);} 519 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 520 ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 521 if (depth > 0) { 522 ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 523 ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 524 ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr); 525 } 526 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 527 odm = idm; 528 } 529 *dmInt = idm; 530 ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 531 PetscFunctionReturn(0); 532 } 533 534 /*@ 535 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 536 537 Collective on DM 538 539 Input Parameter: 540 . dmA - The DMPlex object with initial coordinates 541 542 Output Parameter: 543 . dmB - The DMPlex object with copied coordinates 544 545 Level: intermediate 546 547 Note: This is typically used when adding pieces other than vertices to a mesh 548 549 .keywords: mesh 550 .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 551 @*/ 552 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 553 { 554 Vec coordinatesA, coordinatesB; 555 PetscSection coordSectionA, coordSectionB; 556 PetscScalar *coordsA, *coordsB; 557 PetscInt spaceDim, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 if (dmA == dmB) PetscFunctionReturn(0); 562 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 563 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 564 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); 565 ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 566 ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 567 if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 568 if (!coordSectionB) { 569 PetscInt dim; 570 571 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 572 ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 573 ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 574 ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 575 } 576 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 577 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 578 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 579 ierr = PetscSectionSetChart(coordSectionB, vStartB, vEndB);CHKERRQ(ierr); 580 for (v = vStartB; v < vEndB; ++v) { 581 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 582 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 583 } 584 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 585 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 586 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 587 ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 588 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 589 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 590 ierr = VecSetType(coordinatesB,VECSTANDARD);CHKERRQ(ierr); 591 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 592 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 593 for (v = 0; v < vEndB-vStartB; ++v) { 594 for (d = 0; d < spaceDim; ++d) { 595 coordsB[v*spaceDim+d] = coordsA[v*spaceDim+d]; 596 } 597 } 598 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 599 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 600 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 601 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 602 PetscFunctionReturn(0); 603 } 604 605 /*@ 606 DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 607 608 Collective on DM 609 610 Input Parameter: 611 . dm - The complete DMPlex object 612 613 Output Parameter: 614 . dmUnint - The DMPlex object with only cells and vertices 615 616 Level: intermediate 617 618 .keywords: mesh 619 .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList() 620 @*/ 621 PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 622 { 623 DM udm; 624 PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 629 if (dim <= 1) { 630 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 631 *dmUnint = dm; 632 PetscFunctionReturn(0); 633 } 634 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 635 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 636 ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 637 ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 638 ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 639 ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 640 for (c = cStart; c < cEnd; ++c) { 641 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 642 643 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 644 for (cl = 0; cl < closureSize*2; cl += 2) { 645 const PetscInt p = closure[cl]; 646 647 if ((p >= vStart) && (p < vEnd)) ++coneSize; 648 } 649 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 650 ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 651 maxConeSize = PetscMax(maxConeSize, coneSize); 652 } 653 ierr = DMSetUp(udm);CHKERRQ(ierr); 654 ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 655 for (c = cStart; c < cEnd; ++c) { 656 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 657 658 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 659 for (cl = 0; cl < closureSize*2; cl += 2) { 660 const PetscInt p = closure[cl]; 661 662 if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 663 } 664 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 665 ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 666 } 667 ierr = PetscFree(cone);CHKERRQ(ierr); 668 ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 669 ierr = DMPlexStratify(udm);CHKERRQ(ierr); 670 /* Reduce SF */ 671 { 672 PetscSF sfPoint, sfPointUn; 673 const PetscSFNode *remotePoints; 674 const PetscInt *localPoints; 675 PetscSFNode *remotePointsUn; 676 PetscInt *localPointsUn; 677 PetscInt vEnd, numRoots, numLeaves, l; 678 PetscInt numLeavesUn = 0, n = 0; 679 PetscErrorCode ierr; 680 681 /* Get original SF information */ 682 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 683 ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 684 ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 685 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 686 /* Allocate space for cells and vertices */ 687 for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 688 /* Fill in leaves */ 689 if (vEnd >= 0) { 690 ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 691 ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 692 for (l = 0; l < numLeaves; l++) { 693 if (localPoints[l] < vEnd) { 694 localPointsUn[n] = localPoints[l]; 695 remotePointsUn[n].rank = remotePoints[l].rank; 696 remotePointsUn[n].index = remotePoints[l].index; 697 ++n; 698 } 699 } 700 if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 701 ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 702 } 703 } 704 *dmUnint = udm; 705 PetscFunctionReturn(0); 706 } 707