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