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