1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/hashmapi.h> 3 #include <petsc/private/hashmapij.h> 4 5 /* HashIJKL */ 6 7 #include <petsc/private/hashmap.h> 8 9 typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey; 10 11 #define PetscHashIJKLKeyHash(key) \ 12 PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \ 13 PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l))) 14 15 #define PetscHashIJKLKeyEqual(k1,k2) \ 16 (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0) 17 18 PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1) 19 20 21 /* 22 DMPlexGetFaces_Internal - Gets groups of vertices that correspond to faces for the given cell 23 This assumes that the mesh is not interpolated from the depth of point p to the vertices 24 */ 25 PetscErrorCode DMPlexGetFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 26 { 27 const PetscInt *cone = NULL; 28 PetscInt coneSize; 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 33 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 34 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 35 ierr = DMPlexGetRawFaces_Internal(dm, dim, coneSize, cone, numFaces, faceSize, faces);CHKERRQ(ierr); 36 PetscFunctionReturn(0); 37 } 38 39 /* 40 DMPlexRestoreFaces_Internal - Restores the array 41 */ 42 PetscErrorCode DMPlexRestoreFaces_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 43 { 44 PetscErrorCode ierr; 45 46 PetscFunctionBegin; 47 if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); } 48 PetscFunctionReturn(0); 49 } 50 51 /* 52 DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone 53 */ 54 PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *faceSize, const PetscInt *faces[]) 55 { 56 PetscInt *facesTmp; 57 PetscInt maxConeSize, maxSupportSize; 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 62 if (faces && coneSize) PetscValidIntPointer(cone,4); 63 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 64 if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 65 switch (dim) { 66 case 1: 67 switch (coneSize) { 68 case 2: 69 if (faces) { 70 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 71 *faces = facesTmp; 72 } 73 if (numFaces) *numFaces = 2; 74 if (faceSize) *faceSize = 1; 75 break; 76 default: 77 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 78 } 79 break; 80 case 2: 81 switch (coneSize) { 82 case 3: 83 if (faces) { 84 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 85 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 86 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 87 *faces = facesTmp; 88 } 89 if (numFaces) *numFaces = 3; 90 if (faceSize) *faceSize = 2; 91 break; 92 case 4: 93 /* Vertices follow right hand rule */ 94 if (faces) { 95 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 96 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 97 facesTmp[4] = cone[2]; facesTmp[5] = cone[3]; 98 facesTmp[6] = cone[3]; facesTmp[7] = cone[0]; 99 *faces = facesTmp; 100 } 101 if (numFaces) *numFaces = 4; 102 if (faceSize) *faceSize = 2; 103 break; 104 default: 105 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 106 } 107 break; 108 case 3: 109 switch (coneSize) { 110 case 3: 111 if (faces) { 112 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 113 facesTmp[2] = cone[1]; facesTmp[3] = cone[2]; 114 facesTmp[4] = cone[2]; facesTmp[5] = cone[0]; 115 *faces = facesTmp; 116 } 117 if (numFaces) *numFaces = 3; 118 if (faceSize) *faceSize = 2; 119 break; 120 case 4: 121 /* Vertices of first face follow right hand rule and normal points away from last vertex */ 122 if (faces) { 123 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; 124 facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1]; 125 facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3]; 126 facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3]; 127 *faces = facesTmp; 128 } 129 if (numFaces) *numFaces = 4; 130 if (faceSize) *faceSize = 3; 131 break; 132 case 8: 133 /* 7--------6 134 /| /| 135 / | / | 136 4--------5 | 137 | | | | 138 | | | | 139 | 1--------2 140 | / | / 141 |/ |/ 142 0--------3 143 */ 144 if (faces) { 145 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */ 146 facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */ 147 facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */ 148 facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */ 149 facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */ 150 facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */ 151 *faces = facesTmp; 152 } 153 if (numFaces) *numFaces = 6; 154 if (faceSize) *faceSize = 4; 155 break; 156 default: 157 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 158 } 159 break; 160 default: 161 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 162 } 163 PetscFunctionReturn(0); 164 } 165 166 /* 167 DMPlexGetRawFacesHybrid_Internal - Gets groups of vertices that correspond to faces for the given cone using hybrid ordering (prisms) 168 */ 169 static PetscErrorCode DMPlexGetRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 170 { 171 PetscInt *facesTmp; 172 PetscInt maxConeSize, maxSupportSize; 173 PetscErrorCode ierr; 174 175 PetscFunctionBegin; 176 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 177 if (faces && coneSize) PetscValidIntPointer(cone,4); 178 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 179 if (faces) {ierr = DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);CHKERRQ(ierr);} 180 switch (dim) { 181 case 1: 182 switch (coneSize) { 183 case 2: 184 if (faces) { 185 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 186 *faces = facesTmp; 187 } 188 if (numFaces) *numFaces = 2; 189 if (numFacesNotH) *numFacesNotH = 2; 190 if (faceSize) *faceSize = 1; 191 break; 192 default: 193 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 194 } 195 break; 196 case 2: 197 switch (coneSize) { 198 case 4: 199 if (faces) { 200 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; 201 facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; 202 facesTmp[4] = cone[0]; facesTmp[5] = cone[2]; 203 facesTmp[6] = cone[1]; facesTmp[7] = cone[3]; 204 *faces = facesTmp; 205 } 206 if (numFaces) *numFaces = 4; 207 if (numFacesNotH) *numFacesNotH = 2; 208 if (faceSize) *faceSize = 2; 209 break; 210 default: 211 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 212 } 213 break; 214 case 3: 215 switch (coneSize) { 216 case 6: /* triangular prism */ 217 if (faces) { 218 facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = -1; /* Bottom */ 219 facesTmp[4] = cone[3]; facesTmp[5] = cone[4]; facesTmp[6] = cone[5]; facesTmp[7] = -1; /* Top */ 220 facesTmp[8] = cone[0]; facesTmp[9] = cone[1]; facesTmp[10] = cone[3]; facesTmp[11] = cone[4]; /* Back left */ 221 facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[4]; facesTmp[15] = cone[5]; /* Back right */ 222 facesTmp[16] = cone[2]; facesTmp[17] = cone[0]; facesTmp[18] = cone[5]; facesTmp[19] = cone[3]; /* Front */ 223 *faces = facesTmp; 224 } 225 if (numFaces) *numFaces = 5; 226 if (numFacesNotH) *numFacesNotH = 2; 227 if (faceSize) *faceSize = -4; 228 break; 229 default: 230 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cone size %D not supported for dimension %D", coneSize, dim); 231 } 232 break; 233 default: 234 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim); 235 } 236 PetscFunctionReturn(0); 237 } 238 239 static PetscErrorCode DMPlexRestoreRawFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt coneSize, const PetscInt cone[], PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 240 { 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 if (faces) { ierr = DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);CHKERRQ(ierr); } 245 PetscFunctionReturn(0); 246 } 247 248 static PetscErrorCode DMPlexGetFacesHybrid_Internal(DM dm, PetscInt dim, PetscInt p, PetscInt *numFaces, PetscInt *numFacesNotH, PetscInt *faceSize, const PetscInt *faces[]) 249 { 250 const PetscInt *cone = NULL; 251 PetscInt coneSize; 252 PetscErrorCode ierr; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 256 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 257 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 258 ierr = DMPlexGetRawFacesHybrid_Internal(dm, dim, coneSize, cone, numFaces, numFacesNotH, faceSize, faces);CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 /* This interpolates faces for cells at some stratum */ 263 static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm) 264 { 265 DMLabel subpointMap; 266 PetscHashIJKL faceTable; 267 PetscInt *pStart, *pEnd; 268 PetscInt cellDim, depth, faceDepth = cellDepth, numPoints = 0, faceSizeAll = 0, face, c, d; 269 PetscInt coneSizeH = 0, faceSizeAllH = 0, numCellFacesH = 0, faceH, pMax = -1, dim, outerloop; 270 PetscInt cMax, fMax, eMax, vMax; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = DMGetDimension(dm, &cellDim);CHKERRQ(ierr); 275 /* HACK: I need a better way to determine face dimension, or an alternative to GetFaces() */ 276 ierr = DMPlexGetSubpointMap(dm, &subpointMap);CHKERRQ(ierr); 277 if (subpointMap) ++cellDim; 278 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 279 ++depth; 280 ++cellDepth; 281 cellDim -= depth - cellDepth; 282 ierr = PetscMalloc2(depth+1,&pStart,depth+1,&pEnd);CHKERRQ(ierr); 283 for (d = depth-1; d >= faceDepth; --d) { 284 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d+1], &pEnd[d+1]);CHKERRQ(ierr); 285 } 286 ierr = DMPlexGetDepthStratum(dm, -1, NULL, &pStart[faceDepth]);CHKERRQ(ierr); 287 pEnd[faceDepth] = pStart[faceDepth]; 288 for (d = faceDepth-1; d >= 0; --d) { 289 ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr); 290 } 291 cMax = fMax = eMax = vMax = PETSC_DETERMINE; 292 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 293 if (cellDim == dim) { 294 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 295 pMax = cMax; 296 } else if (cellDim == dim -1) { 297 ierr = DMPlexGetHybridBounds(dm, &cMax, &fMax, NULL, NULL);CHKERRQ(ierr); 298 pMax = fMax; 299 } 300 pMax = pMax < 0 ? pEnd[cellDepth] : pMax; 301 if (pMax < pEnd[cellDepth]) { 302 const PetscInt *cellFaces, *cone; 303 PetscInt numCellFacesT, faceSize, cf; 304 305 ierr = DMPlexGetConeSize(dm, pMax, &coneSizeH);CHKERRQ(ierr); 306 ierr = DMPlexGetCone(dm, pMax, &cone);CHKERRQ(ierr); 307 ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr); 308 if (faceSize < 0) { 309 PetscInt *sizes, minv, maxv; 310 311 /* count vertices of hybrid and non-hybrid faces */ 312 ierr = PetscCalloc1(numCellFacesH, &sizes);CHKERRQ(ierr); 313 for (cf = 0; cf < numCellFacesT; ++cf) { /* These are the non-hybrid faces */ 314 const PetscInt *cellFace = &cellFaces[-cf*faceSize]; 315 PetscInt f; 316 317 for (f = 0; f < -faceSize; ++f) sizes[cf] += (cellFace[f] >= 0 ? 1 : 0); 318 } 319 ierr = PetscSortInt(numCellFacesT, sizes);CHKERRQ(ierr); 320 minv = sizes[0]; 321 maxv = sizes[PetscMax(numCellFacesT-1, 0)]; 322 if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for non-hybrid face %D != %D", minv, maxv); 323 faceSizeAll = minv; 324 ierr = PetscMemzero(sizes, numCellFacesH*sizeof(PetscInt));CHKERRQ(ierr); 325 for (cf = numCellFacesT; cf < numCellFacesH; ++cf) { /* These are the hybrid faces */ 326 const PetscInt *cellFace = &cellFaces[-cf*faceSize]; 327 PetscInt f; 328 329 for (f = 0; f < -faceSize; ++f) sizes[cf-numCellFacesT] += (cellFace[f] >= 0 ? 1 : 0); 330 } 331 ierr = PetscSortInt(numCellFacesH - numCellFacesT, sizes);CHKERRQ(ierr); 332 minv = sizes[0]; 333 maxv = sizes[PetscMax(numCellFacesH - numCellFacesT-1, 0)]; 334 ierr = PetscFree(sizes);CHKERRQ(ierr); 335 if (minv != maxv) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Different number of vertices for hybrid face %D != %D", minv, maxv); 336 faceSizeAllH = minv; 337 } else { /* the size of the faces in hybrid cells is the same */ 338 faceSizeAll = faceSizeAllH = faceSize; 339 } 340 ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, cone, &numCellFacesH, &numCellFacesT, &faceSize, &cellFaces);CHKERRQ(ierr); 341 } else if (pEnd[cellDepth] > pStart[cellDepth]) { 342 ierr = DMPlexGetFaces_Internal(dm, cellDim, pStart[cellDepth], NULL, &faceSizeAll, NULL);CHKERRQ(ierr); 343 } 344 if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 345 346 /* With hybrid grids, we first iterate on hybrid cells and start numbering the non-hybrid faces 347 Then, faces for non-hybrid cells are numbered. 348 This is to guarantee consistent orientations (all 0) of all the points in the cone of the hybrid cells */ 349 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 350 for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) { 351 PetscInt start, end; 352 353 start = outerloop == 0 ? pMax : pStart[cellDepth]; 354 end = outerloop == 0 ? pEnd[cellDepth] : pMax; 355 for (c = start; c < end; ++c) { 356 const PetscInt *cellFaces; 357 PetscInt numCellFaces, faceSize, faceSizeInc, cf; 358 359 if (c < pMax) { 360 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 361 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 362 } else { /* Hybrid cell */ 363 const PetscInt *cone; 364 PetscInt numCellFacesN, coneSize; 365 366 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 367 if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH); 368 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 369 ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 370 if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c); 371 faceSize = PetscMax(faceSize, -faceSize); 372 if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize); 373 numCellFaces = numCellFacesN; /* process only non-hybrid faces */ 374 } 375 faceSizeInc = faceSize; 376 for (cf = 0; cf < numCellFaces; ++cf) { 377 const PetscInt *cellFace = &cellFaces[cf*faceSizeInc]; 378 PetscInt faceSizeH = faceSize; 379 PetscHashIJKLKey key; 380 PetscHashIter iter; 381 PetscBool missing; 382 383 if (faceSizeInc == 2) { 384 key.i = PetscMin(cellFace[0], cellFace[1]); 385 key.j = PetscMax(cellFace[0], cellFace[1]); 386 key.k = PETSC_MAX_INT; 387 key.l = PETSC_MAX_INT; 388 } else { 389 key.i = cellFace[0]; 390 key.j = cellFace[1]; 391 key.k = cellFace[2]; 392 key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 393 ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 394 } 395 /* this check is redundant for non-hybrid meshes */ 396 if (faceSizeH != faceSizeAll) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAll); 397 ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 398 if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);} 399 } 400 if (c < pMax) { 401 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 402 } else { 403 ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr); 404 } 405 } 406 } 407 pEnd[faceDepth] = face; 408 409 /* Second pass for hybrid meshes: number hybrid faces */ 410 for (c = pMax; c < pEnd[cellDepth]; ++c) { 411 const PetscInt *cellFaces, *cone; 412 PetscInt numCellFaces, numCellFacesN, faceSize, cf, coneSize; 413 414 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 415 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 416 ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 417 if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH); 418 faceSize = PetscMax(faceSize, -faceSize); 419 for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */ 420 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 421 PetscHashIJKLKey key; 422 PetscHashIter iter; 423 PetscBool missing; 424 PetscInt faceSizeH = faceSize; 425 426 if (faceSize == 2) { 427 key.i = PetscMin(cellFace[0], cellFace[1]); 428 key.j = PetscMax(cellFace[0], cellFace[1]); 429 key.k = PETSC_MAX_INT; 430 key.l = PETSC_MAX_INT; 431 } else { 432 key.i = cellFace[0]; 433 key.j = cellFace[1]; 434 key.k = cellFace[2]; 435 key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 436 ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 437 } 438 if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH); 439 ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 440 if (missing) {ierr = PetscHashIJKLIterSet(faceTable, iter, face++);CHKERRQ(ierr);} 441 } 442 ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 443 } 444 faceH = face - pEnd[faceDepth]; 445 if (faceH) { 446 if (fMax == PETSC_DETERMINE) fMax = pEnd[faceDepth]; 447 else if (eMax == PETSC_DETERMINE) eMax = pEnd[faceDepth]; 448 else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of unassigned hybrid facets %D for cellDim %D and dimension %D", faceH, cellDim, dim); 449 } 450 pEnd[faceDepth] = face; 451 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 452 /* Count new points */ 453 for (d = 0; d <= depth; ++d) { 454 numPoints += pEnd[d]-pStart[d]; 455 } 456 ierr = DMPlexSetChart(idm, 0, numPoints);CHKERRQ(ierr); 457 /* Set cone sizes */ 458 for (d = 0; d <= depth; ++d) { 459 PetscInt coneSize, p; 460 461 if (d == faceDepth) { 462 /* I see no way to do this if we admit faces of different shapes */ 463 for (p = pStart[d]; p < pEnd[d]-faceH; ++p) { 464 ierr = DMPlexSetConeSize(idm, p, faceSizeAll);CHKERRQ(ierr); 465 } 466 for (p = pEnd[d]-faceH; p < pEnd[d]; ++p) { 467 ierr = DMPlexSetConeSize(idm, p, faceSizeAllH);CHKERRQ(ierr); 468 } 469 } else if (d == cellDepth) { 470 for (p = pStart[d]; p < pEnd[d]; ++p) { 471 /* Number of cell faces may be different from number of cell vertices*/ 472 if (p < pMax) { 473 ierr = DMPlexGetFaces_Internal(dm, cellDim, p, &coneSize, NULL, NULL);CHKERRQ(ierr); 474 } else { 475 ierr = DMPlexGetFacesHybrid_Internal(dm, cellDim, p, &coneSize, NULL, NULL, NULL);CHKERRQ(ierr); 476 } 477 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 478 } 479 } else { 480 for (p = pStart[d]; p < pEnd[d]; ++p) { 481 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 482 ierr = DMPlexSetConeSize(idm, p, coneSize);CHKERRQ(ierr); 483 } 484 } 485 } 486 ierr = DMSetUp(idm);CHKERRQ(ierr); 487 /* Get face cones from subsets of cell vertices */ 488 if (faceSizeAll > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSizeAll); 489 ierr = PetscHashIJKLCreate(&faceTable);CHKERRQ(ierr); 490 for (d = depth; d > cellDepth; --d) { 491 const PetscInt *cone; 492 PetscInt p; 493 494 for (p = pStart[d]; p < pEnd[d]; ++p) { 495 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 496 ierr = DMPlexSetCone(idm, p, cone);CHKERRQ(ierr); 497 ierr = DMPlexGetConeOrientation(dm, p, &cone);CHKERRQ(ierr); 498 ierr = DMPlexSetConeOrientation(idm, p, cone);CHKERRQ(ierr); 499 } 500 } 501 for (outerloop = 0, face = pStart[faceDepth]; outerloop < 2; outerloop++) { 502 PetscInt start, end; 503 504 start = outerloop == 0 ? pMax : pStart[cellDepth]; 505 end = outerloop == 0 ? pEnd[cellDepth] : pMax; 506 for (c = start; c < end; ++c) { 507 const PetscInt *cellFaces; 508 PetscInt numCellFaces, faceSize, faceSizeInc, cf; 509 510 if (c < pMax) { 511 ierr = DMPlexGetFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 512 if (faceSize != faceSizeAll) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent face for cell %D of size %D != %D", c, faceSize, faceSizeAll); 513 } else { 514 const PetscInt *cone; 515 PetscInt numCellFacesN, coneSize; 516 517 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 518 if (coneSize != coneSizeH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid coneSize %D != %D", coneSize, coneSizeH); 519 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 520 ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 521 if (numCellFaces != numCellFacesH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected numCellFaces %D != %D for hybrid cell %D", numCellFaces, numCellFacesH, c); 522 faceSize = PetscMax(faceSize, -faceSize); 523 if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not support interpolation of meshes with faces of %D vertices", faceSize); 524 numCellFaces = numCellFacesN; /* process only non-hybrid faces */ 525 } 526 faceSizeInc = faceSize; 527 for (cf = 0; cf < numCellFaces; ++cf) { 528 const PetscInt *cellFace = &cellFaces[cf*faceSizeInc]; 529 PetscHashIJKLKey key; 530 PetscHashIter iter; 531 PetscBool missing; 532 533 if (faceSizeInc == 2) { 534 key.i = PetscMin(cellFace[0], cellFace[1]); 535 key.j = PetscMax(cellFace[0], cellFace[1]); 536 key.k = PETSC_MAX_INT; 537 key.l = PETSC_MAX_INT; 538 } else { 539 key.i = cellFace[0]; 540 key.j = cellFace[1]; 541 key.k = cellFace[2]; 542 key.l = faceSizeInc > 3 ? (cellFace[3] < 0 ? faceSize = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 543 ierr = PetscSortInt(faceSizeInc, (PetscInt *) &key);CHKERRQ(ierr); 544 } 545 ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 546 if (missing) { 547 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 548 ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr); 549 ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 550 } else { 551 const PetscInt *cone; 552 PetscInt coneSize, ornt, i, j, f; 553 554 ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 555 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 556 /* Orient face: Do not allow reverse orientation at the first vertex */ 557 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 558 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 559 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); 560 /* - First find the initial vertex */ 561 for (i = 0; i < faceSize; ++i) if (cellFace[0] == cone[i]) break; 562 /* - Try forward comparison */ 563 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+j)%faceSize]) break; 564 if (j == faceSize) { 565 if ((faceSize == 2) && (i == 1)) ornt = -2; 566 else ornt = i; 567 } else { 568 /* - Try backward comparison */ 569 for (j = 0; j < faceSize; ++j) if (cellFace[j] != cone[(i+faceSize-j)%faceSize]) break; 570 if (j == faceSize) { 571 if (i == 0) ornt = -faceSize; 572 else ornt = -i; 573 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 574 } 575 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 576 } 577 } 578 if (c < pMax) { 579 ierr = DMPlexRestoreFaces_Internal(dm, cellDim, c, &numCellFaces, &faceSize, &cellFaces);CHKERRQ(ierr); 580 } else { 581 ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSizeH, NULL, NULL, NULL, NULL, &cellFaces);CHKERRQ(ierr); 582 } 583 } 584 } 585 /* Second pass for hybrid meshes: orient hybrid faces */ 586 for (c = pMax; c < pEnd[cellDepth]; ++c) { 587 const PetscInt *cellFaces, *cone; 588 PetscInt numCellFaces, numCellFacesN, faceSize, cf, coneSize; 589 590 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 591 ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr); 592 ierr = DMPlexGetRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 593 if (numCellFaces != numCellFacesH) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected hybrid numCellFaces %D != %D", numCellFaces, numCellFacesH); 594 faceSize = PetscMax(faceSize, -faceSize); 595 for (cf = numCellFacesN; cf < numCellFaces; ++cf) { /* These are the hybrid faces */ 596 const PetscInt *cellFace = &cellFaces[cf*faceSize]; 597 PetscHashIJKLKey key; 598 PetscHashIter iter; 599 PetscBool missing; 600 PetscInt faceSizeH = faceSize; 601 602 if (faceSize == 2) { 603 key.i = PetscMin(cellFace[0], cellFace[1]); 604 key.j = PetscMax(cellFace[0], cellFace[1]); 605 key.k = PETSC_MAX_INT; 606 key.l = PETSC_MAX_INT; 607 } else { 608 key.i = cellFace[0]; 609 key.j = cellFace[1]; 610 key.k = cellFace[2]; 611 key.l = faceSize > 3 ? (cellFace[3] < 0 ? faceSizeH = 3, PETSC_MAX_INT : cellFace[3]) : PETSC_MAX_INT; 612 ierr = PetscSortInt(faceSize, (PetscInt *) &key);CHKERRQ(ierr); 613 } 614 if (faceSizeH != faceSizeAllH) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unexpected number of vertices for hybrid face %D of point %D -> %D != %D", cf, c, faceSizeH, faceSizeAllH); 615 ierr = PetscHashIJKLPut(faceTable, key, &iter, &missing);CHKERRQ(ierr); 616 if (missing) { 617 ierr = DMPlexSetCone(idm, face, cellFace);CHKERRQ(ierr); 618 ierr = PetscHashIJKLIterSet(faceTable, iter, face);CHKERRQ(ierr); 619 ierr = DMPlexInsertCone(idm, c, cf, face++);CHKERRQ(ierr); 620 } else { 621 const PetscInt *cone; 622 PetscInt coneSize, ornt, i, j, f; 623 624 ierr = PetscHashIJKLIterGet(faceTable, iter, &f);CHKERRQ(ierr); 625 ierr = DMPlexInsertCone(idm, c, cf, f);CHKERRQ(ierr); 626 /* Orient face: Do not allow reverse orientation at the first vertex */ 627 ierr = DMPlexGetConeSize(idm, f, &coneSize);CHKERRQ(ierr); 628 ierr = DMPlexGetCone(idm, f, &cone);CHKERRQ(ierr); 629 if (coneSize != faceSizeH) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSizeH); 630 /* - First find the initial vertex */ 631 for (i = 0; i < faceSizeH; ++i) if (cellFace[0] == cone[i]) break; 632 /* - Try forward comparison */ 633 for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+j)%faceSizeH]) break; 634 if (j == faceSizeH) { 635 if ((faceSizeH == 2) && (i == 1)) ornt = -2; 636 else ornt = i; 637 } else { 638 /* - Try backward comparison */ 639 for (j = 0; j < faceSizeH; ++j) if (cellFace[j] != cone[(i+faceSizeH-j)%faceSizeH]) break; 640 if (j == faceSizeH) { 641 if (i == 0) ornt = -faceSizeH; 642 else ornt = -i; 643 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine face orientation"); 644 } 645 ierr = DMPlexInsertConeOrientation(idm, c, cf, ornt);CHKERRQ(ierr); 646 } 647 } 648 ierr = DMPlexRestoreRawFacesHybrid_Internal(dm, cellDim, coneSize, cone, &numCellFaces, &numCellFacesN, &faceSize, &cellFaces);CHKERRQ(ierr); 649 } 650 if (face != pEnd[faceDepth]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of faces %D should be %D", face-pStart[faceDepth], pEnd[faceDepth]-pStart[faceDepth]); 651 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 652 ierr = PetscHashIJKLDestroy(&faceTable);CHKERRQ(ierr); 653 ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr); 654 ierr = DMPlexSetHybridBounds(idm, cMax, fMax, eMax, vMax);CHKERRQ(ierr); 655 ierr = DMPlexSymmetrize(idm);CHKERRQ(ierr); 656 ierr = DMPlexStratify(idm);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 PetscErrorCode DMPlexOrientCell(DM dm, PetscInt p, PetscInt masterConeSize, const PetscInt masterCone[]) 661 { 662 PetscInt coneSize; 663 PetscInt start1=0; 664 PetscBool reverse1=PETSC_FALSE; 665 const PetscInt *cone=NULL; 666 PetscErrorCode ierr; 667 668 PetscFunctionBegin; 669 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 670 if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 671 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 672 ierr = DMPlexFixFaceOrientations_Orient_Private(coneSize, masterConeSize, masterCone, cone, &start1, &reverse1);CHKERRQ(ierr); 673 #if defined(PETSC_USE_DEBUG) 674 if (PetscUnlikely(cone[start1] != masterCone[0])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[0]", start1, cone[start1], masterCone[0]); 675 #endif 676 ierr = DMPlexOrientCell_Internal(dm, p, start1, reverse1);CHKERRQ(ierr); 677 #if defined(PETSC_USE_DEBUG) 678 { 679 PetscInt c; 680 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 681 for (c = 0; c < 2; c++) { 682 if (PetscUnlikely(cone[c] != masterCone[c])) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "The algorithm above is wrong as cone[%d] = %d != %d = masterCone[%d]", c, cone[c], masterCone[c], c); 683 } 684 } 685 #endif 686 PetscFunctionReturn(0); 687 } 688 689 PetscErrorCode DMPlexOrientCell_Internal(DM dm, PetscInt p, PetscInt start1, PetscBool reverse1) 690 { 691 PetscInt i, j, k, maxConeSize, coneSize, coneConeSize, supportSize, supportConeSize; 692 PetscInt start0, start; 693 PetscBool reverse0, reverse; 694 PetscInt newornt; 695 const PetscInt *cone=NULL, *support=NULL, *supportCone=NULL, *ornts=NULL; 696 PetscInt *newcone=NULL, *newornts=NULL; 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 if (!start1 && !reverse1) PetscFunctionReturn(0); 701 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 702 if (!coneSize) PetscFunctionReturn(0); /* do nothing for points with no cone */ 703 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 704 ierr = DMPlexGetMaxSizes(dm, &maxConeSize, NULL);CHKERRQ(ierr); 705 /* permute p's cone and orientations */ 706 ierr = DMPlexGetConeOrientation(dm, p, &ornts);CHKERRQ(ierr); 707 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 708 ierr = DMGetWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 709 ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, cone, start1, reverse1, newcone);CHKERRQ(ierr); 710 ierr = DMPlexFixFaceOrientations_Permute_Private(coneSize, ornts, start1, reverse1, newornts);CHKERRQ(ierr); 711 /* if direction of p (face) is flipped, flip also p's cone points (edges) */ 712 if (reverse1) { 713 for (i=0; i<coneSize; i++) { 714 ierr = DMPlexGetConeSize(dm, cone[i], &coneConeSize);CHKERRQ(ierr); 715 ierr = DMPlexFixFaceOrientations_Translate_Private(newornts[i], &start0, &reverse0);CHKERRQ(ierr); 716 ierr = DMPlexFixFaceOrientations_Combine_Private(coneConeSize, start0, reverse0, 1, PETSC_FALSE, &start, &reverse);CHKERRQ(ierr); 717 ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneConeSize, start, reverse, &newornts[i]);CHKERRQ(ierr); 718 } 719 } 720 ierr = DMPlexSetConeOrientation(dm, p, newornts);CHKERRQ(ierr); 721 /* fix oriention of p within cones of p's support points */ 722 ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 723 ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 724 for (j=0; j<supportSize; j++) { 725 ierr = DMPlexGetCone(dm, support[j], &supportCone);CHKERRQ(ierr); 726 ierr = DMPlexGetConeSize(dm, support[j], &supportConeSize);CHKERRQ(ierr); 727 ierr = DMPlexGetConeOrientation(dm, support[j], &ornts);CHKERRQ(ierr); 728 for (k=0; k<supportConeSize; k++) { 729 if (supportCone[k] != p) continue; 730 ierr = DMPlexFixFaceOrientations_Translate_Private(ornts[k], &start0, &reverse0);CHKERRQ(ierr); 731 ierr = DMPlexFixFaceOrientations_Combine_Private(coneSize, start0, reverse0, start1, reverse1, &start, &reverse);CHKERRQ(ierr); 732 ierr = DMPlexFixFaceOrientations_TranslateBack_Private(coneSize, start, reverse, &newornt);CHKERRQ(ierr); 733 ierr = DMPlexInsertConeOrientation(dm, support[j], k, newornt);CHKERRQ(ierr); 734 } 735 } 736 /* rewrite cone */ 737 ierr = DMPlexSetCone(dm, p, newcone);CHKERRQ(ierr); 738 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newcone);CHKERRQ(ierr); 739 ierr = DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &newornts);CHKERRQ(ierr); 740 PetscFunctionReturn(0); 741 } 742 743 static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 744 { 745 PetscInt nleaves; 746 PetscInt nranks; 747 const PetscMPIInt *ranks=NULL; 748 const PetscInt *roffset=NULL, *rmine=NULL, *rremote=NULL; 749 PetscInt n, o, r; 750 PetscErrorCode ierr; 751 752 PetscFunctionBegin; 753 ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr); 754 nleaves = roffset[nranks]; 755 ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr); 756 for (r=0; r<nranks; r++) { 757 /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 758 - to unify order with the other side */ 759 o = roffset[r]; 760 n = roffset[r+1] - o; 761 ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr); 762 ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr); 763 ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr); 764 } 765 PetscFunctionReturn(0); 766 } 767 768 PetscErrorCode DMPlexOrientInterface(DM dm) 769 { 770 PetscSF sf=NULL; 771 PetscInt (*roots)[2], (*leaves)[2]; 772 PetscMPIInt (*rootsRanks)[2], (*leavesRanks)[2]; 773 const PetscInt *locals=NULL; 774 const PetscSFNode *remotes=NULL; 775 PetscInt nroots, nleaves, p, c; 776 PetscInt nranks, n, o, r; 777 const PetscMPIInt *ranks=NULL; 778 const PetscInt *roffset=NULL; 779 PetscInt *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 780 const PetscInt *cone=NULL; 781 PetscInt coneSize, ind0; 782 MPI_Comm comm; 783 PetscMPIInt rank; 784 PetscInt debug = 0; 785 PetscErrorCode ierr; 786 787 PetscFunctionBegin; 788 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 789 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr); 790 if (nroots < 0) PetscFunctionReturn(0); 791 ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 792 ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr); 793 #if defined(PETSC_USE_DEBUG) 794 ierr = DMViewFromOptions(dm, NULL, "-before_fix_dm_view");CHKERRQ(ierr); 795 ierr = DMPlexCheckPointSF(dm);CHKERRQ(ierr); 796 #endif 797 ierr = SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr); 798 ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr); 799 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 800 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 801 if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);} 802 for (p = 0; p < nroots; ++p) { 803 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 804 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 805 /* Translate all points to root numbering */ 806 for (c = 0; c < 2; c++) { 807 if (coneSize > 1) { 808 ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 809 if (ind0 < 0) { 810 roots[p][c] = cone[c]; 811 rootsRanks[p][c] = rank; 812 } else { 813 roots[p][c] = remotes[ind0].index; 814 rootsRanks[p][c] = remotes[ind0].rank; 815 } 816 } else { 817 roots[p][c] = -1; 818 rootsRanks[p][c] = -1; 819 } 820 } 821 } 822 if (debug) { 823 for (p = 0; p < nroots; ++p) { 824 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 825 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 826 if (coneSize > 1) { 827 ierr = PetscSynchronizedPrintf(comm, "[%d] %D: cone=[%D %D] roots=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], roots[p][0], roots[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr); 828 } 829 } 830 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 831 } 832 for (p = 0; p < nroots; ++p) { 833 for (c = 0; c < 2; c++) { 834 leaves[p][c] = -2; 835 leavesRanks[p][c] = -2; 836 } 837 } 838 ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 839 ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 840 ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 841 ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 842 if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 843 if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);} 844 for (p = 0; p < nroots; ++p) { 845 if (leaves[p][0] < 0) continue; 846 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 847 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 848 if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %D: cone=[%D %D] leaves=[%D %D] roots=[%D %D] leavesRanks=[%D %D] rootsRanks=[%D %D]\n", rank, p, cone[0], cone[1], leaves[p][0], leaves[p][1], roots[p][0], roots[p][1], leavesRanks[p][0], leavesRanks[p][1], rootsRanks[p][0], rootsRanks[p][1]);CHKERRQ(ierr);} 849 if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][0] != rootsRanks[p][0])) { 850 PetscInt masterCone[2]; 851 /* Translate these two cone points back to leave numbering */ 852 for (c = 0; c < 2; c++) { 853 if (leavesRanks[p][c] == rank) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - remote rank of point %D is the same rank",leavesRanks[p][c]); 854 /* Find index of rank leavesRanks[p][c] among remote ranks */ 855 /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */ 856 ierr = PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);CHKERRQ(ierr); 857 if (PetscUnlikely(r < 0)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "this should never happen - rank %D not found among remote ranks",leavesRanks[p][c]); 858 /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */ 859 o = roffset[r]; 860 n = roffset[r+1] - o; 861 ierr = PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);CHKERRQ(ierr); 862 if (PetscUnlikely(ind0 < 0)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "No cone point of %D is connected to (%D, %D) - it seems there is missing connection in point SF!",p,ranks[r],leaves[p][c]); 863 /* Get the corresponding local point */ 864 masterCone[c] = rmine1[o+ind0];CHKERRQ(ierr); 865 } 866 if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %D: masterCone=[%D %D]\n", rank, p, masterCone[0], masterCone[1]);CHKERRQ(ierr);} 867 /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */ 868 ierr = DMPlexOrientCell(dm, p, 2, masterCone);CHKERRQ(ierr); 869 } 870 } 871 #if defined(PETSC_USE_DEBUG) 872 ierr = DMViewFromOptions(dm, NULL, "-after_fix_dm_view");CHKERRQ(ierr); 873 for (r = 0; r < nleaves; ++r) { 874 p = locals[r]; 875 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 876 if (!coneSize) continue; 877 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 878 for (c = 0; c < 2; c++) { 879 ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 880 if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing its cone point cone[%D] = %D!", p, c, cone[c]); 881 if (leaves[p][c] != remotes[ind0].index || leavesRanks[p][c] != remotes[ind0].rank) { 882 if (leavesRanks[p][c] == rank) { 883 PetscInt ind1; 884 ierr = PetscFindInt(leaves[p][c], nleaves, locals, &ind1);CHKERRQ(ierr); 885 if (ind1 < 0) { 886 SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). The latter was not even found among the local SF points - it is probably broken!", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]); 887 } else { 888 SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced %D --> (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leaves[p][c], remotes[ind1].rank, remotes[ind1].index); 889 } 890 } else { 891 SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D = locals[%d]: cone[%D]=%D --> (%D, %D) differs from the enforced (%D, %D). Is the algorithm above or the point SF broken?", p, r, c, cone[c], remotes[ind0].rank, remotes[ind0].index, leavesRanks[p][c], leaves[p][c]); 892 } 893 } 894 } 895 } 896 #endif 897 if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 898 ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr); 899 ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr); 900 PetscFunctionReturn(0); 901 } 902 903 static PetscErrorCode IntArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], const char valname[], PetscInt n, const PetscInt a[]) 904 { 905 PetscInt idx; 906 PetscMPIInt rank; 907 PetscBool flg; 908 PetscErrorCode ierr; 909 910 PetscFunctionBegin; 911 ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 912 if (!flg) PetscFunctionReturn(0); 913 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 914 ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 915 for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D %s %D\n", rank, idxname, idx, valname, a[idx]);CHKERRQ(ierr);} 916 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 917 PetscFunctionReturn(0); 918 } 919 920 static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[]) 921 { 922 PetscInt idx; 923 PetscMPIInt rank; 924 PetscBool flg; 925 PetscErrorCode ierr; 926 927 PetscFunctionBegin; 928 ierr = PetscOptionsHasName(NULL, NULL, opt, &flg);CHKERRQ(ierr); 929 if (!flg) PetscFunctionReturn(0); 930 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 931 ierr = PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);CHKERRQ(ierr); 932 if (idxname) { 933 for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]%s %D rank %D index %D\n", rank, idxname, idx, a[idx].rank, a[idx].index);CHKERRQ(ierr);} 934 } else { 935 for (idx = 0; idx < n; ++idx) {ierr = PetscSynchronizedPrintf(comm, "[%d]rank %D index %D\n", rank, a[idx].rank, a[idx].index);CHKERRQ(ierr);} 936 } 937 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 938 PetscFunctionReturn(0); 939 } 940 941 static PetscErrorCode DMPlexMapToLocalPoint(PetscHMapIJ roothash, const PetscInt localPoints[], PetscMPIInt rank, PetscSFNode remotePoint, PetscInt *localPoint) 942 { 943 PetscErrorCode ierr; 944 945 PetscFunctionBegin; 946 if (remotePoint.rank == rank) { 947 *localPoint = remotePoint.index; 948 } else { 949 PetscHashIJKey key; 950 PetscInt root; 951 952 key.i = remotePoint.index; 953 key.j = remotePoint.rank; 954 ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 955 if (root >= 0) { 956 *localPoint = localPoints[root]; 957 } else PetscFunctionReturn(1); 958 } 959 PetscFunctionReturn(0); 960 } 961 962 /*@ 963 DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation 964 965 Collective on DM 966 967 Input Parameters: 968 + dm - The interpolated DM 969 - pointSF - The initial SF without interpolated points 970 971 Output Parameter: 972 . pointSF - The SF including interpolated points 973 974 Level: intermediate 975 976 Note: All debugging for this process can be turned on with the options: -dm_interp_pre_view -petscsf_interp_pre_view -petscsection_interp_candidate_view -petscsection_interp_candidate_remote_view -petscsection_interp_claim_view -petscsf_interp_pre_view -dmplex_interp_debug 977 978 .keywords: mesh 979 .seealso: DMPlexInterpolate(), DMPlexUninterpolate() 980 @*/ 981 PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF) 982 { 983 /* 984 Okay, the algorithm is: 985 - Take each point in the overlap (root) 986 - Look at the neighboring points in the overlap (candidates) 987 - Send these candidate points to neighbors 988 - Neighbor checks for edge between root and candidate 989 - If edge is found, it replaces candidate point with edge point 990 - Send back the overwritten candidates (claims) 991 - Original guy checks for edges, different from original candidate, and gets its own edge 992 - This pair is put into SF 993 994 We need a new algorithm that tolerates groups larger than 2. 995 - Take each point in the overlap (root) 996 - Find all collections of points in the overlap which make faces (do early join) 997 - Send collections as candidates (add size as first number) 998 - Make sure to send collection to all owners of all overlap points in collection 999 - Neighbor check for face in collections 1000 - If face is found, it replaces candidate point with face point 1001 - Send back the overwritten candidates (claims) 1002 - Original guy checks for faces, different from original candidate, and gets its own face 1003 - This pair is put into SF 1004 */ 1005 PetscHMapI leafhash; 1006 PetscHMapIJ roothash; 1007 const PetscInt *localPoints, *rootdegree; 1008 const PetscSFNode *remotePoints; 1009 PetscSFNode *candidates, *candidatesRemote, *claims; 1010 PetscSection candidateSection, candidateSectionRemote, claimSection; 1011 PetscInt numLeaves, l, numRoots, r, candidatesSize, candidatesRemoteSize; 1012 PetscMPIInt size, rank; 1013 PetscHashIJKey key; 1014 PetscBool debug = PETSC_FALSE; 1015 PetscErrorCode ierr; 1016 1017 PetscFunctionBegin; 1018 ierr = PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);CHKERRQ(ierr); 1019 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 1020 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1021 ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 1022 if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 1023 ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");CHKERRQ(ierr); 1024 ierr = PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");CHKERRQ(ierr); 1025 ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 1026 /* Build hashes of points in the SF for efficient lookup */ 1027 ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr); 1028 ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr); 1029 for (l = 0; l < numLeaves; ++l) { 1030 key.i = remotePoints[l].index; 1031 key.j = remotePoints[l].rank; 1032 ierr = PetscHMapISet(leafhash, localPoints[l], l);CHKERRQ(ierr); 1033 ierr = PetscHMapIJSet(roothash, key, l);CHKERRQ(ierr); 1034 } 1035 /* Compute root degree to identify shared points */ 1036 ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 1037 ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 1038 ierr = IntArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-interp_root_degree_view", "Root degree", "point", "degree", numRoots, rootdegree);CHKERRQ(ierr); 1039 /* Build a section / SFNode array of candidate points (face bd points) in the cone(support(leaf)), 1040 where each candidate is defined by a set of remote points (roots) for the other points that define the face. */ 1041 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 1042 ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 1043 { 1044 PetscHMapIJ facehash; 1045 1046 ierr = PetscHMapIJCreate(&facehash);CHKERRQ(ierr); 1047 for (l = 0; l < numLeaves; ++l) { 1048 const PetscInt localPoint = localPoints[l]; 1049 const PetscInt *support; 1050 PetscInt supportSize, s; 1051 1052 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);} 1053 ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr); 1054 ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr); 1055 for (s = 0; s < supportSize; ++s) { 1056 const PetscInt face = support[s]; 1057 const PetscInt *cone; 1058 PetscInt coneSize, c, f, root; 1059 PetscBool isFace = PETSC_TRUE; 1060 1061 /* Only add face once */ 1062 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Support point %D\n", rank, face);CHKERRQ(ierr);} 1063 key.i = localPoint; 1064 key.j = face; 1065 ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr); 1066 if (f >= 0) continue; 1067 ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 1068 ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr); 1069 /* If a cone point does not map to leaves on any proc, then do not put face in SF */ 1070 for (c = 0; c < coneSize; ++c) { 1071 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);} 1072 ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 1073 if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;} 1074 } 1075 if (isFace) { 1076 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found shared face %D\n", rank, face);CHKERRQ(ierr);} 1077 ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr); 1078 ierr = PetscSectionAddDof(candidateSection, localPoint, coneSize);CHKERRQ(ierr); 1079 } 1080 } 1081 } 1082 if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 1083 ierr = PetscHMapIJClear(facehash);CHKERRQ(ierr); 1084 ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 1085 ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 1086 ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 1087 for (l = 0; l < numLeaves; ++l) { 1088 const PetscInt localPoint = localPoints[l]; 1089 const PetscInt *support; 1090 PetscInt supportSize, s, offset, idx = 0; 1091 1092 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local point %D\n", rank, localPoint);CHKERRQ(ierr);} 1093 ierr = PetscSectionGetOffset(candidateSection, localPoint, &offset);CHKERRQ(ierr); 1094 ierr = DMPlexGetSupportSize(dm, localPoint, &supportSize);CHKERRQ(ierr); 1095 ierr = DMPlexGetSupport(dm, localPoint, &support);CHKERRQ(ierr); 1096 for (s = 0; s < supportSize; ++s) { 1097 const PetscInt face = support[s]; 1098 const PetscInt *cone; 1099 PetscInt coneSize, c, f, root; 1100 PetscBool isFace = PETSC_TRUE; 1101 1102 /* Only add face once */ 1103 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Support point %D\n", rank, face);CHKERRQ(ierr);} 1104 key.i = localPoint; 1105 key.j = face; 1106 ierr = PetscHMapIJGet(facehash, key, &f);CHKERRQ(ierr); 1107 if (f >= 0) continue; 1108 ierr = DMPlexGetConeSize(dm, face, &coneSize);CHKERRQ(ierr); 1109 ierr = DMPlexGetCone(dm, face, &cone);CHKERRQ(ierr); 1110 /* If a cone point does not map to leaves on any proc, then do not put face in SF */ 1111 for (c = 0; c < coneSize; ++c) { 1112 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Cone point %D\n", rank, cone[c]);CHKERRQ(ierr);} 1113 ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 1114 if (!rootdegree[cone[c]] && (root < 0)) {isFace = PETSC_FALSE; break;} 1115 } 1116 if (isFace) { 1117 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Adding shared face %D at idx %D\n", rank, face, idx);CHKERRQ(ierr);} 1118 ierr = PetscHMapIJSet(facehash, key, l);CHKERRQ(ierr); 1119 candidates[offset+idx].rank = -1; 1120 candidates[offset+idx++].index = coneSize-1; 1121 for (c = 0; c < coneSize; ++c) { 1122 if (cone[c] == localPoint) continue; 1123 if (rootdegree[cone[c]]) { 1124 candidates[offset+idx].rank = rank; 1125 candidates[offset+idx++].index = cone[c]; 1126 } else { 1127 ierr = PetscHMapIGet(leafhash, cone[c], &root);CHKERRQ(ierr); 1128 if (root < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot locate local point %D in SF", cone[c]); 1129 candidates[offset+idx++] = remotePoints[root]; 1130 } 1131 } 1132 } 1133 } 1134 } 1135 if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 1136 ierr = PetscHMapIJDestroy(&facehash);CHKERRQ(ierr); 1137 ierr = PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");CHKERRQ(ierr); 1138 ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);CHKERRQ(ierr); 1139 } 1140 /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 1141 /* Note that this section is indexed by offsets into leaves, not by point number */ 1142 { 1143 PetscSF sfMulti, sfInverse, sfCandidates; 1144 PetscInt *remoteOffsets; 1145 1146 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 1147 ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 1148 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 1149 ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 1150 ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 1151 ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 1152 ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 1153 ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 1154 ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 1155 ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 1156 ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 1157 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1158 1159 ierr = PetscObjectViewFromOptions((PetscObject) candidateSectionRemote, NULL, "-petscsection_interp_candidate_remote_view");CHKERRQ(ierr); 1160 ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);CHKERRQ(ierr); 1161 } 1162 /* */ 1163 { 1164 PetscInt idx; 1165 /* There is a section point for every leaf attached to a given root point */ 1166 for (r = 0, idx = 0; r < numRoots; ++r) { 1167 PetscInt deg; 1168 for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) { 1169 PetscInt offset, dof, d; 1170 1171 ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 1172 ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 1173 for (d = 0; d < dof; ++d) { 1174 const PetscInt sizeInd = offset+d; 1175 const PetscInt numPoints = candidatesRemote[sizeInd].index; 1176 const PetscInt *join = NULL; 1177 PetscInt points[1024], p, joinSize; 1178 1179 points[0] = r; 1180 for (p = 0; p < numPoints; ++p) { 1181 ierr = DMPlexMapToLocalPoint(roothash, localPoints, rank, candidatesRemote[offset+(++d)], &points[p+1]); 1182 if (ierr) {d += numPoints-1 - p; break;} /* We got a point not in our overlap */ 1183 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local candidate %D\n", rank, points[p+1]);CHKERRQ(ierr);} 1184 } 1185 if (ierr) continue; 1186 ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 1187 if (joinSize == 1) { 1188 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Adding face %D at idx %D\n", rank, join[0], sizeInd);CHKERRQ(ierr);} 1189 candidatesRemote[sizeInd].rank = rank; 1190 candidatesRemote[sizeInd].index = join[0]; 1191 } 1192 ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 1193 } 1194 } 1195 } 1196 if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 1197 } 1198 /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 1199 { 1200 PetscSF sfMulti, sfClaims, sfPointNew; 1201 PetscSFNode *remotePointsNew; 1202 PetscHMapI claimshash; 1203 PetscInt *remoteOffsets, *localPointsNew; 1204 PetscInt claimsSize, pStart, pEnd, root, numLocalNew, p, d; 1205 1206 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 1207 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 1208 ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 1209 ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 1210 ierr = PetscSectionGetStorageSize(claimSection, &claimsSize);CHKERRQ(ierr); 1211 ierr = PetscMalloc1(claimsSize, &claims);CHKERRQ(ierr); 1212 ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 1213 ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 1214 ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 1215 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 1216 ierr = PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");CHKERRQ(ierr); 1217 ierr = SFNodeArrayViewFromOptions(PetscObjectComm((PetscObject) dm), "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);CHKERRQ(ierr); 1218 /* Walk the original section of local supports and add an SF entry for each updated item */ 1219 ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr); 1220 for (p = 0; p < numRoots; ++p) { 1221 PetscInt dof, offset; 1222 1223 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking root for claims %D\n", rank, p);CHKERRQ(ierr);} 1224 ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 1225 ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 1226 for (d = 0; d < dof;) { 1227 if (claims[offset+d].rank >= 0) { 1228 const PetscInt faceInd = offset+d; 1229 const PetscInt numPoints = candidates[faceInd].index; 1230 const PetscInt *join = NULL; 1231 PetscInt joinSize, points[1024], c; 1232 1233 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);CHKERRQ(ierr);} 1234 points[0] = p; 1235 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] point %D\n", rank, points[0]);CHKERRQ(ierr);} 1236 for (c = 0, ++d; c < numPoints; ++c, ++d) { 1237 key.i = candidates[offset+d].index; 1238 key.j = candidates[offset+d].rank; 1239 ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 1240 points[c+1] = localPoints[root]; 1241 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] point %D\n", rank, points[c+1]);CHKERRQ(ierr);} 1242 } 1243 ierr = DMPlexGetJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 1244 if (joinSize == 1) { 1245 if (debug) {ierr = PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Found local face %D\n", rank, join[0]);CHKERRQ(ierr);} 1246 ierr = PetscHMapISet(claimshash, join[0], faceInd);CHKERRQ(ierr); 1247 } 1248 ierr = DMPlexRestoreJoin(dm, numPoints+1, points, &joinSize, &join);CHKERRQ(ierr); 1249 } else d += claims[offset+d].index+1; 1250 } 1251 } 1252 if (debug) {ierr = PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);CHKERRQ(ierr);} 1253 /* Create new pointSF from hashed claims */ 1254 ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr); 1255 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1256 ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 1257 ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 1258 for (p = 0; p < numLeaves; ++p) { 1259 localPointsNew[p] = localPoints[p]; 1260 remotePointsNew[p].index = remotePoints[p].index; 1261 remotePointsNew[p].rank = remotePoints[p].rank; 1262 } 1263 p = numLeaves; 1264 ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 1265 ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr); 1266 for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 1267 PetscInt offset; 1268 ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr); 1269 remotePointsNew[p] = claims[offset]; 1270 } 1271 ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 1272 ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 1273 ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 1274 ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 1275 ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr); 1276 } 1277 ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr); 1278 ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr); 1279 ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 1280 ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 1281 ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 1282 ierr = PetscFree(candidates);CHKERRQ(ierr); 1283 ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 1284 ierr = PetscFree(claims);CHKERRQ(ierr); 1285 ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 1286 PetscFunctionReturn(0); 1287 } 1288 1289 /*@C 1290 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 1291 1292 Collective on DM 1293 1294 Input Parameters: 1295 + dm - The DMPlex object with only cells and vertices 1296 - dmInt - The interpolated DM 1297 1298 Output Parameter: 1299 . dmInt - The complete DMPlex object 1300 1301 Level: intermediate 1302 1303 Notes: 1304 It does not copy over the coordinates. 1305 1306 .keywords: mesh 1307 .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 1308 @*/ 1309 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 1310 { 1311 DM idm, odm = dm; 1312 PetscSF sfPoint; 1313 PetscInt depth, dim, d; 1314 const char *name; 1315 PetscBool flg=PETSC_TRUE; 1316 PetscErrorCode ierr; 1317 1318 PetscFunctionBegin; 1319 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1320 PetscValidPointer(dmInt, 2); 1321 ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 1322 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1323 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1324 if ((depth == dim) || (dim <= 1)) { 1325 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 1326 idm = dm; 1327 } else { 1328 for (d = 1; d < dim; ++d) { 1329 /* Create interpolated mesh */ 1330 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 1331 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 1332 ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 1333 if (depth > 0) { 1334 ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 1335 ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 1336 ierr = DMPlexInterpolatePointSF(idm, sfPoint);CHKERRQ(ierr); 1337 } 1338 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 1339 odm = idm; 1340 } 1341 ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 1342 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 1343 ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr); 1344 ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr); 1345 ierr = PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);CHKERRQ(ierr); 1346 if (flg) {ierr = DMPlexOrientInterface(idm);CHKERRQ(ierr);} 1347 } 1348 { 1349 PetscBool isper; 1350 const PetscReal *maxCell, *L; 1351 const DMBoundaryType *bd; 1352 1353 ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 1354 ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr); 1355 } 1356 *dmInt = idm; 1357 ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 1358 PetscFunctionReturn(0); 1359 } 1360 1361 /*@ 1362 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 1363 1364 Collective on DM 1365 1366 Input Parameter: 1367 . dmA - The DMPlex object with initial coordinates 1368 1369 Output Parameter: 1370 . dmB - The DMPlex object with copied coordinates 1371 1372 Level: intermediate 1373 1374 Note: This is typically used when adding pieces other than vertices to a mesh 1375 1376 .keywords: mesh 1377 .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 1378 @*/ 1379 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 1380 { 1381 Vec coordinatesA, coordinatesB; 1382 VecType vtype; 1383 PetscSection coordSectionA, coordSectionB; 1384 PetscScalar *coordsA, *coordsB; 1385 PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 1386 PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE; 1387 PetscBool lc = PETSC_FALSE; 1388 PetscErrorCode ierr; 1389 1390 PetscFunctionBegin; 1391 PetscValidHeaderSpecific(dmA, DM_CLASSID, 1); 1392 PetscValidHeaderSpecific(dmB, DM_CLASSID, 2); 1393 if (dmA == dmB) PetscFunctionReturn(0); 1394 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 1395 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 1396 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); 1397 ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr); 1398 ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr); 1399 ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 1400 ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 1401 if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 1402 ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr); 1403 if (!Nf) PetscFunctionReturn(0); 1404 if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf); 1405 if (!coordSectionB) { 1406 PetscInt dim; 1407 1408 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 1409 ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 1410 ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 1411 ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 1412 } 1413 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 1414 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 1415 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 1416 ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr); 1417 if (cStartA <= cS && cS < cEndA) { /* localized coordinates */ 1418 if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cells in first DM %D != %D in the second DM", cEndA-cStartA, cEndB-cStartB); 1419 cS = cS - cStartA + cStartB; 1420 cE = vEndB; 1421 lc = PETSC_TRUE; 1422 } else { 1423 cS = vStartB; 1424 cE = vEndB; 1425 } 1426 ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr); 1427 for (v = vStartB; v < vEndB; ++v) { 1428 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 1429 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 1430 } 1431 if (lc) { /* localized coordinates */ 1432 PetscInt c; 1433 1434 for (c = cS-cStartB; c < cEndB-cStartB; c++) { 1435 PetscInt dof; 1436 1437 ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 1438 ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr); 1439 ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr); 1440 } 1441 } 1442 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 1443 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 1444 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 1445 ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 1446 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 1447 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 1448 ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr); 1449 ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr); 1450 ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr); 1451 ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr); 1452 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 1453 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 1454 for (v = 0; v < vEndB-vStartB; ++v) { 1455 PetscInt offA, offB; 1456 1457 ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr); 1458 ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr); 1459 for (d = 0; d < spaceDim; ++d) { 1460 coordsB[offB+d] = coordsA[offA+d]; 1461 } 1462 } 1463 if (lc) { /* localized coordinates */ 1464 PetscInt c; 1465 1466 for (c = cS-cStartB; c < cEndB-cStartB; c++) { 1467 PetscInt dof, offA, offB; 1468 1469 ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr); 1470 ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr); 1471 ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 1472 ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr); 1473 } 1474 } 1475 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 1476 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 1477 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 1478 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 1479 PetscFunctionReturn(0); 1480 } 1481 1482 /*@ 1483 DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 1484 1485 Collective on DM 1486 1487 Input Parameter: 1488 . dm - The complete DMPlex object 1489 1490 Output Parameter: 1491 . dmUnint - The DMPlex object with only cells and vertices 1492 1493 Level: intermediate 1494 1495 Notes: 1496 It does not copy over the coordinates. 1497 1498 .keywords: mesh 1499 .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 1500 @*/ 1501 PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 1502 { 1503 DM udm; 1504 PetscInt dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone; 1505 PetscErrorCode ierr; 1506 1507 PetscFunctionBegin; 1508 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1509 PetscValidPointer(dmUnint, 2); 1510 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1511 if (dim <= 1) { 1512 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 1513 *dmUnint = dm; 1514 PetscFunctionReturn(0); 1515 } 1516 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1517 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1518 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1519 ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 1520 ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 1521 ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 1522 ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 1523 for (c = cStart; c < cEnd; ++c) { 1524 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 1525 1526 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1527 for (cl = 0; cl < closureSize*2; cl += 2) { 1528 const PetscInt p = closure[cl]; 1529 1530 if ((p >= vStart) && (p < vEnd)) ++coneSize; 1531 } 1532 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1533 ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 1534 maxConeSize = PetscMax(maxConeSize, coneSize); 1535 } 1536 ierr = DMSetUp(udm);CHKERRQ(ierr); 1537 ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 1538 for (c = cStart; c < cEnd; ++c) { 1539 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 1540 1541 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1542 for (cl = 0; cl < closureSize*2; cl += 2) { 1543 const PetscInt p = closure[cl]; 1544 1545 if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 1546 } 1547 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1548 ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 1549 } 1550 ierr = PetscFree(cone);CHKERRQ(ierr); 1551 ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1552 ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 1553 ierr = DMPlexStratify(udm);CHKERRQ(ierr); 1554 /* Reduce SF */ 1555 { 1556 PetscSF sfPoint, sfPointUn; 1557 const PetscSFNode *remotePoints; 1558 const PetscInt *localPoints; 1559 PetscSFNode *remotePointsUn; 1560 PetscInt *localPointsUn; 1561 PetscInt vEnd, numRoots, numLeaves, l; 1562 PetscInt numLeavesUn = 0, n = 0; 1563 PetscErrorCode ierr; 1564 1565 /* Get original SF information */ 1566 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1567 ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 1568 ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 1569 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 1570 /* Allocate space for cells and vertices */ 1571 for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 1572 /* Fill in leaves */ 1573 if (vEnd >= 0) { 1574 ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 1575 ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 1576 for (l = 0; l < numLeaves; l++) { 1577 if (localPoints[l] < vEnd) { 1578 localPointsUn[n] = localPoints[l]; 1579 remotePointsUn[n].rank = remotePoints[l].rank; 1580 remotePointsUn[n].index = remotePoints[l].index; 1581 ++n; 1582 } 1583 } 1584 if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 1585 ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 1586 } 1587 } 1588 { 1589 PetscBool isper; 1590 const PetscReal *maxCell, *L; 1591 const DMBoundaryType *bd; 1592 1593 ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 1594 ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr); 1595 } 1596 1597 *dmUnint = udm; 1598 PetscFunctionReturn(0); 1599 } 1600