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 grids 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 static PetscErrorCode DMPlexOrientPointSF_Internal(DM dm, PetscSF sf) 661 { 662 PetscInt (*roots)[2], (*leaves)[2]; 663 PetscMPIInt (*rootsRanks)[2], (*leavesRanks)[2]; 664 const PetscInt *locals; 665 const PetscSFNode *remotes; 666 PetscInt nroots, nleaves, p, c; 667 const PetscInt *cone; 668 PetscInt coneSize, ind0, ind1; 669 MPI_Comm comm; 670 PetscMPIInt rank; 671 PetscInt debug = 0; 672 PetscErrorCode ierr; 673 674 PetscFunctionBegin; 675 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);CHKERRQ(ierr); 676 if (nroots < 0) PetscFunctionReturn(0); 677 ierr = PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);CHKERRQ(ierr); 678 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 679 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 680 if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Roots\n");CHKERRQ(ierr);} 681 for (p = 0; p < nroots; ++p) { 682 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 683 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 684 /* Translate all points to root numbering */ 685 for (c = 0; c < 2; c++) { 686 if (coneSize > 1) { 687 ierr = PetscFindInt(cone[c], nleaves, locals, &ind0);CHKERRQ(ierr); 688 if (ind0 < 0) { 689 roots[p][c] = cone[c]; 690 rootsRanks[p][c] = rank; 691 } else { 692 roots[p][c] = remotes[ind0].index; 693 rootsRanks[p][c] = remotes[ind0].rank; 694 } 695 } else { 696 roots[p][c] = -1; 697 rootsRanks[p][c] = -1; 698 } 699 } 700 } 701 if (debug) { 702 for (p = 0; p < nroots; ++p) { 703 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 704 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 705 if (coneSize > 1) { 706 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); 707 } 708 } 709 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 710 } 711 for (p = 0; p < nroots; ++p) { 712 leaves[p][0] = -2; 713 leaves[p][1] = -2; 714 } 715 ierr = PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 716 ierr = PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 717 ierr = PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);CHKERRQ(ierr); 718 ierr = PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);CHKERRQ(ierr); 719 if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 720 if (debug && rank == 0) {ierr = PetscSynchronizedPrintf(comm, "Referred leaves\n");CHKERRQ(ierr);} 721 for (p = 0; p < nroots; ++p) { 722 if (leaves[p][0] < 0) continue; 723 if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] %D: %D %D\n", rank, p, leaves[p][0], leaves[p][1]);CHKERRQ(ierr);} 724 725 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 726 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 727 if (coneSize > 1) { 728 ierr = PetscFindInt(cone[0], nleaves, locals, &ind0);CHKERRQ(ierr); 729 ierr = PetscFindInt(cone[1], nleaves, locals, &ind1);CHKERRQ(ierr); 730 if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 0, cone[0]); 731 if (ind1 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 1, cone[1]); 732 if ((leaves[p][0] != remotes[ind0].index) || (leaves[p][1] != remotes[ind1].index)) { 733 const PetscInt *newcone; 734 735 /* TODO Generalize this to correct an arbitrary orientation */ 736 ierr = DMPlexReverseCell(dm, p);CHKERRQ(ierr); 737 if (debug) {ierr = PetscSynchronizedPrintf(comm, "[%d] reversed point %D\n", rank, p);CHKERRQ(ierr);} 738 ierr = DMPlexGetCone(dm, p, &newcone);CHKERRQ(ierr); 739 ierr = PetscFindInt(newcone[0], nleaves, locals, &ind0);CHKERRQ(ierr); 740 ierr = PetscFindInt(newcone[1], nleaves, locals, &ind1);CHKERRQ(ierr); 741 if (ind0 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 0, cone[0]); 742 if (ind1 < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but is missing cone point %D = %D", p, 1, cone[1]); 743 if (leaves[p][0] != remotes[ind0].index) 744 SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D Cone[%D] %D --> (%D, %D) != %D Cone[%D] SF Point (%D, %D)", p, 0, cone[0], remotes[ind0].rank, remotes[ind0].index, leaves[p][0], 0, remotes[p].rank, remotes[p].index); 745 if (leaves[p][1] != remotes[ind1].index) 746 SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D Cone[%D] %D --> (%D, %D) != %D Cone[%D] SF Point (%D, %D)", p, 1, cone[1], remotes[ind1].rank, remotes[ind1].index, leaves[p][1], 1, remotes[p].rank, remotes[p].index); 747 } 748 } 749 } 750 if (debug) {ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr);} 751 ierr = PetscFree4(roots, leaves, rootsRanks, leavesRanks);CHKERRQ(ierr); 752 PetscFunctionReturn(0); 753 } 754 755 /* This interpolates the PointSF in parallel following local interpolation */ 756 static PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF, PetscInt depth) 757 { 758 PetscMPIInt size, rank; 759 PetscInt p, c, d, dof, offset; 760 PetscInt numLeaves, numRoots, candidatesSize, candidatesRemoteSize; 761 const PetscInt *localPoints; 762 const PetscSFNode *remotePoints; 763 PetscSFNode *candidates, *candidatesRemote, *claims; 764 PetscSection candidateSection, candidateSectionRemote, claimSection; 765 PetscHMapI leafhash; 766 PetscHMapIJ roothash; 767 PetscHashIJKey key; 768 PetscErrorCode ierr; 769 770 PetscFunctionBegin; 771 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 772 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 773 ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 774 if (size < 2 || numRoots < 0) PetscFunctionReturn(0); 775 ierr = PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 776 /* Build hashes of points in the SF for efficient lookup */ 777 ierr = PetscHMapICreate(&leafhash);CHKERRQ(ierr); 778 ierr = PetscHMapIJCreate(&roothash);CHKERRQ(ierr); 779 for (p = 0; p < numLeaves; ++p) { 780 ierr = PetscHMapISet(leafhash, localPoints[p], p);CHKERRQ(ierr); 781 key.i = remotePoints[p].index; 782 key.j = remotePoints[p].rank; 783 ierr = PetscHMapIJSet(roothash, key, p);CHKERRQ(ierr); 784 } 785 /* Build a section / SFNode array of candidate points in the single-level adjacency of leaves, 786 where each candidate is defined by the root entry for the other vertex that defines the edge. */ 787 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSection);CHKERRQ(ierr); 788 ierr = PetscSectionSetChart(candidateSection, 0, numRoots);CHKERRQ(ierr); 789 { 790 PetscInt leaf, root, idx, a, *adj = NULL; 791 for (p = 0; p < numLeaves; ++p) { 792 PetscInt adjSize = PETSC_DETERMINE; 793 ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 794 for (a = 0; a < adjSize; ++a) { 795 ierr = PetscHMapIGet(leafhash, adj[a], &leaf);CHKERRQ(ierr); 796 if (leaf >= 0) {ierr = PetscSectionAddDof(candidateSection, localPoints[p], 1);CHKERRQ(ierr);} 797 } 798 } 799 ierr = PetscSectionSetUp(candidateSection);CHKERRQ(ierr); 800 ierr = PetscSectionGetStorageSize(candidateSection, &candidatesSize);CHKERRQ(ierr); 801 ierr = PetscMalloc1(candidatesSize, &candidates);CHKERRQ(ierr); 802 for (p = 0; p < numLeaves; ++p) { 803 PetscInt adjSize = PETSC_DETERMINE; 804 ierr = PetscSectionGetOffset(candidateSection, localPoints[p], &offset);CHKERRQ(ierr); 805 ierr = DMPlexGetAdjacency_Internal(dm, localPoints[p], PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);CHKERRQ(ierr); 806 for (idx = 0, a = 0; a < adjSize; ++a) { 807 ierr = PetscHMapIGet(leafhash, adj[a], &root);CHKERRQ(ierr); 808 if (root >= 0) candidates[offset+idx++] = remotePoints[root]; 809 } 810 } 811 ierr = PetscFree(adj);CHKERRQ(ierr); 812 } 813 /* Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */ 814 { 815 PetscSF sfMulti, sfInverse, sfCandidates; 816 PetscInt *remoteOffsets; 817 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 818 ierr = PetscSFCreateInverseSF(sfMulti, &sfInverse);CHKERRQ(ierr); 819 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &candidateSectionRemote);CHKERRQ(ierr); 820 ierr = PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateSectionRemote);CHKERRQ(ierr); 821 ierr = PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateSectionRemote, &sfCandidates);CHKERRQ(ierr); 822 ierr = PetscSectionGetStorageSize(candidateSectionRemote, &candidatesRemoteSize);CHKERRQ(ierr); 823 ierr = PetscMalloc1(candidatesRemoteSize, &candidatesRemote);CHKERRQ(ierr); 824 ierr = PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 825 ierr = PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);CHKERRQ(ierr); 826 ierr = PetscSFDestroy(&sfInverse);CHKERRQ(ierr); 827 ierr = PetscSFDestroy(&sfCandidates);CHKERRQ(ierr); 828 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 829 } 830 /* Walk local roots and check for each remote candidate whether we know all required points, 831 either from owning it or having a root entry in the point SF. If we do we place a claim 832 by replacing the vertex number with our edge ID. */ 833 { 834 PetscInt idx, root, joinSize, vertices[2]; 835 const PetscInt *rootdegree, *join = NULL; 836 ierr = PetscSFComputeDegreeBegin(pointSF, &rootdegree);CHKERRQ(ierr); 837 ierr = PetscSFComputeDegreeEnd(pointSF, &rootdegree);CHKERRQ(ierr); 838 /* Loop remote edge connections and put in a claim if both vertices are known */ 839 for (idx = 0, p = 0; p < numRoots; ++p) { 840 for (d = 0; d < rootdegree[p]; ++d) { 841 ierr = PetscSectionGetDof(candidateSectionRemote, idx, &dof);CHKERRQ(ierr); 842 ierr = PetscSectionGetOffset(candidateSectionRemote, idx, &offset);CHKERRQ(ierr); 843 for (c = 0; c < dof; ++c) { 844 /* We own both vertices, so we claim the edge by replacing vertex with edge */ 845 if (candidatesRemote[offset+c].rank == rank) { 846 vertices[0] = p; vertices[1] = candidatesRemote[offset+c].index; 847 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 848 if (joinSize == 1) candidatesRemote[offset+c].index = join[0]; 849 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 850 continue; 851 } 852 /* If we own one vertex and share a root with the other, we claim it */ 853 key.i = candidatesRemote[offset+c].index; 854 key.j = candidatesRemote[offset+c].rank; 855 ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 856 if (root >= 0) { 857 vertices[0] = p; vertices[1] = localPoints[root]; 858 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 859 if (joinSize == 1) { 860 candidatesRemote[offset+c].index = join[0]; 861 candidatesRemote[offset+c].rank = rank; 862 } 863 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 864 } 865 } 866 idx++; 867 } 868 } 869 } 870 /* Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */ 871 { 872 PetscSF sfMulti, sfClaims, sfPointNew; 873 PetscHMapI claimshash; 874 PetscInt size, pStart, pEnd, root, joinSize, numLocalNew; 875 PetscInt *remoteOffsets, *localPointsNew, vertices[2]; 876 const PetscInt *join = NULL; 877 PetscSFNode *remotePointsNew; 878 ierr = PetscSFGetMultiSF(pointSF, &sfMulti);CHKERRQ(ierr); 879 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &claimSection);CHKERRQ(ierr); 880 ierr = PetscSFDistributeSection(sfMulti, candidateSectionRemote, &remoteOffsets, claimSection);CHKERRQ(ierr); 881 ierr = PetscSFCreateSectionSF(sfMulti, candidateSectionRemote, remoteOffsets, claimSection, &sfClaims);CHKERRQ(ierr); 882 ierr = PetscSectionGetStorageSize(claimSection, &size);CHKERRQ(ierr); 883 ierr = PetscMalloc1(size, &claims);CHKERRQ(ierr); 884 ierr = PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 885 ierr = PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);CHKERRQ(ierr); 886 ierr = PetscSFDestroy(&sfClaims);CHKERRQ(ierr); 887 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 888 /* Walk the original section of local supports and add an SF entry for each updated item */ 889 ierr = PetscHMapICreate(&claimshash);CHKERRQ(ierr); 890 for (p = 0; p < numRoots; ++p) { 891 ierr = PetscSectionGetDof(candidateSection, p, &dof);CHKERRQ(ierr); 892 ierr = PetscSectionGetOffset(candidateSection, p, &offset);CHKERRQ(ierr); 893 for (d = 0; d < dof; ++d) { 894 if (candidates[offset+d].index != claims[offset+d].index) { 895 key.i = candidates[offset+d].index; 896 key.j = candidates[offset+d].rank; 897 ierr = PetscHMapIJGet(roothash, key, &root);CHKERRQ(ierr); 898 if (root >= 0) { 899 vertices[0] = p; vertices[1] = localPoints[root]; 900 ierr = DMPlexGetJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 901 if (joinSize == 1) {ierr = PetscHMapISet(claimshash, join[0], offset+d);CHKERRQ(ierr);} 902 ierr = DMPlexRestoreJoin(dm, 2, vertices, &joinSize, &join);CHKERRQ(ierr); 903 } 904 } 905 } 906 } 907 /* Create new pointSF from hashed claims */ 908 ierr = PetscHMapIGetSize(claimshash, &numLocalNew);CHKERRQ(ierr); 909 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 910 ierr = PetscMalloc1(numLeaves + numLocalNew, &localPointsNew);CHKERRQ(ierr); 911 ierr = PetscMalloc1(numLeaves + numLocalNew, &remotePointsNew);CHKERRQ(ierr); 912 for (p = 0; p < numLeaves; ++p) { 913 localPointsNew[p] = localPoints[p]; 914 remotePointsNew[p].index = remotePoints[p].index; 915 remotePointsNew[p].rank = remotePoints[p].rank; 916 } 917 p = numLeaves; 918 ierr = PetscHMapIGetKeys(claimshash, &p, localPointsNew);CHKERRQ(ierr); 919 ierr = PetscSortInt(numLocalNew, &localPointsNew[numLeaves]);CHKERRQ(ierr); 920 for (p = numLeaves; p < numLeaves + numLocalNew; ++p) { 921 ierr = PetscHMapIGet(claimshash, localPointsNew[p], &offset);CHKERRQ(ierr); 922 remotePointsNew[p] = claims[offset]; 923 } 924 ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), &sfPointNew);CHKERRQ(ierr); 925 ierr = PetscSFSetGraph(sfPointNew, pEnd-pStart, numLeaves+numLocalNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 926 ierr = DMSetPointSF(dm, sfPointNew);CHKERRQ(ierr); 927 ierr = PetscSFDestroy(&sfPointNew);CHKERRQ(ierr); 928 ierr = PetscHMapIDestroy(&claimshash);CHKERRQ(ierr); 929 } 930 ierr = PetscHMapIDestroy(&leafhash);CHKERRQ(ierr); 931 ierr = PetscHMapIJDestroy(&roothash);CHKERRQ(ierr); 932 ierr = PetscSectionDestroy(&candidateSection);CHKERRQ(ierr); 933 ierr = PetscSectionDestroy(&candidateSectionRemote);CHKERRQ(ierr); 934 ierr = PetscSectionDestroy(&claimSection);CHKERRQ(ierr); 935 ierr = PetscFree(candidates);CHKERRQ(ierr); 936 ierr = PetscFree(candidatesRemote);CHKERRQ(ierr); 937 ierr = PetscFree(claims);CHKERRQ(ierr); 938 ierr = PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);CHKERRQ(ierr); 939 PetscFunctionReturn(0); 940 } 941 942 /*@C 943 DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc. 944 945 Collective on DM 946 947 Input Parameters: 948 + dm - The DMPlex object with only cells and vertices 949 - dmInt - The interpolated DM 950 951 Output Parameter: 952 . dmInt - The complete DMPlex object 953 954 Level: intermediate 955 956 Notes: 957 It does not copy over the coordinates. 958 959 .keywords: mesh 960 .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 961 @*/ 962 PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt) 963 { 964 DM idm, odm = dm; 965 PetscSF sfPoint; 966 PetscInt depth, dim, d; 967 const char *name; 968 PetscErrorCode ierr; 969 970 PetscFunctionBegin; 971 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 972 PetscValidPointer(dmInt, 2); 973 ierr = PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 974 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 975 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 976 if ((depth == dim) || (dim <= 1)) { 977 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 978 idm = dm; 979 } else { 980 for (d = 1; d < dim; ++d) { 981 /* Create interpolated mesh */ 982 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &idm);CHKERRQ(ierr); 983 ierr = DMSetType(idm, DMPLEX);CHKERRQ(ierr); 984 ierr = DMSetDimension(idm, dim);CHKERRQ(ierr); 985 if (depth > 0) { 986 ierr = DMPlexInterpolateFaces_Internal(odm, 1, idm);CHKERRQ(ierr); 987 ierr = DMGetPointSF(odm, &sfPoint);CHKERRQ(ierr); 988 ierr = DMPlexInterpolatePointSF(idm, sfPoint, depth);CHKERRQ(ierr); 989 } 990 if (odm != dm) {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 991 odm = idm; 992 } 993 ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 994 ierr = PetscObjectSetName((PetscObject) idm, name);CHKERRQ(ierr); 995 ierr = DMPlexCopyCoordinates(dm, idm);CHKERRQ(ierr); 996 ierr = DMCopyLabels(dm, idm);CHKERRQ(ierr); 997 { 998 /* TODO temporary */ 999 PetscBool flg=PETSC_FALSE; 1000 ierr = PetscOptionsGetBool(NULL,NULL, "-dm_plex_check_point_sf", &flg, NULL);CHKERRQ(ierr); 1001 if (flg) {ierr = DMPlexCheckPointSF(idm);CHKERRQ(ierr);} 1002 ierr = DMViewFromOptions(idm, NULL, "-before_fix_dm_view");CHKERRQ(ierr); 1003 } 1004 ierr = DMGetPointSF(idm, &sfPoint);CHKERRQ(ierr); 1005 ierr = DMPlexOrientPointSF_Internal(idm, sfPoint);CHKERRQ(ierr); 1006 } 1007 { 1008 PetscBool isper; 1009 const PetscReal *maxCell, *L; 1010 const DMBoundaryType *bd; 1011 1012 ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 1013 ierr = DMSetPeriodicity(idm,isper,maxCell,L,bd);CHKERRQ(ierr); 1014 } 1015 *dmInt = idm; 1016 ierr = PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);CHKERRQ(ierr); 1017 PetscFunctionReturn(0); 1018 } 1019 1020 /*@ 1021 DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices 1022 1023 Collective on DM 1024 1025 Input Parameter: 1026 . dmA - The DMPlex object with initial coordinates 1027 1028 Output Parameter: 1029 . dmB - The DMPlex object with copied coordinates 1030 1031 Level: intermediate 1032 1033 Note: This is typically used when adding pieces other than vertices to a mesh 1034 1035 .keywords: mesh 1036 .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 1037 @*/ 1038 PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB) 1039 { 1040 Vec coordinatesA, coordinatesB; 1041 VecType vtype; 1042 PetscSection coordSectionA, coordSectionB; 1043 PetscScalar *coordsA, *coordsB; 1044 PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d; 1045 PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE; 1046 PetscBool lc = PETSC_FALSE; 1047 PetscErrorCode ierr; 1048 1049 PetscFunctionBegin; 1050 PetscValidHeaderSpecific(dmA, DM_CLASSID, 1); 1051 PetscValidHeaderSpecific(dmB, DM_CLASSID, 2); 1052 if (dmA == dmB) PetscFunctionReturn(0); 1053 ierr = DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);CHKERRQ(ierr); 1054 ierr = DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);CHKERRQ(ierr); 1055 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); 1056 ierr = DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);CHKERRQ(ierr); 1057 ierr = DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);CHKERRQ(ierr); 1058 ierr = DMGetCoordinateSection(dmA, &coordSectionA);CHKERRQ(ierr); 1059 ierr = DMGetCoordinateSection(dmB, &coordSectionB);CHKERRQ(ierr); 1060 if (coordSectionA == coordSectionB) PetscFunctionReturn(0); 1061 ierr = PetscSectionGetNumFields(coordSectionA, &Nf);CHKERRQ(ierr); 1062 if (!Nf) PetscFunctionReturn(0); 1063 if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf); 1064 if (!coordSectionB) { 1065 PetscInt dim; 1066 1067 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);CHKERRQ(ierr); 1068 ierr = DMGetCoordinateDim(dmA, &dim);CHKERRQ(ierr); 1069 ierr = DMSetCoordinateSection(dmB, dim, coordSectionB);CHKERRQ(ierr); 1070 ierr = PetscObjectDereference((PetscObject) coordSectionB);CHKERRQ(ierr); 1071 } 1072 ierr = PetscSectionSetNumFields(coordSectionB, 1);CHKERRQ(ierr); 1073 ierr = PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);CHKERRQ(ierr); 1074 ierr = PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);CHKERRQ(ierr); 1075 ierr = PetscSectionGetChart(coordSectionA, &cS, &cE);CHKERRQ(ierr); 1076 if (cStartA <= cS && cS < cEndA) { /* localized coordinates */ 1077 if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cellls in first DM %d != %d in the second DM", cEndA-cStartA, cEndB-cStartB); 1078 cS = cS - cStartA + cStartB; 1079 cE = vEndB; 1080 lc = PETSC_TRUE; 1081 } else { 1082 cS = vStartB; 1083 cE = vEndB; 1084 } 1085 ierr = PetscSectionSetChart(coordSectionB, cS, cE);CHKERRQ(ierr); 1086 for (v = vStartB; v < vEndB; ++v) { 1087 ierr = PetscSectionSetDof(coordSectionB, v, spaceDim);CHKERRQ(ierr); 1088 ierr = PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);CHKERRQ(ierr); 1089 } 1090 if (lc) { /* localized coordinates */ 1091 PetscInt c; 1092 1093 for (c = cS-cStartB; c < cEndB-cStartB; c++) { 1094 PetscInt dof; 1095 1096 ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 1097 ierr = PetscSectionSetDof(coordSectionB, c + cStartB, dof);CHKERRQ(ierr); 1098 ierr = PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);CHKERRQ(ierr); 1099 } 1100 } 1101 ierr = PetscSectionSetUp(coordSectionB);CHKERRQ(ierr); 1102 ierr = PetscSectionGetStorageSize(coordSectionB, &coordSizeB);CHKERRQ(ierr); 1103 ierr = DMGetCoordinatesLocal(dmA, &coordinatesA);CHKERRQ(ierr); 1104 ierr = VecCreate(PETSC_COMM_SELF, &coordinatesB);CHKERRQ(ierr); 1105 ierr = PetscObjectSetName((PetscObject) coordinatesB, "coordinates");CHKERRQ(ierr); 1106 ierr = VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);CHKERRQ(ierr); 1107 ierr = VecGetBlockSize(coordinatesA, &d);CHKERRQ(ierr); 1108 ierr = VecSetBlockSize(coordinatesB, d);CHKERRQ(ierr); 1109 ierr = VecGetType(coordinatesA, &vtype);CHKERRQ(ierr); 1110 ierr = VecSetType(coordinatesB, vtype);CHKERRQ(ierr); 1111 ierr = VecGetArray(coordinatesA, &coordsA);CHKERRQ(ierr); 1112 ierr = VecGetArray(coordinatesB, &coordsB);CHKERRQ(ierr); 1113 for (v = 0; v < vEndB-vStartB; ++v) { 1114 PetscInt offA, offB; 1115 1116 ierr = PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);CHKERRQ(ierr); 1117 ierr = PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);CHKERRQ(ierr); 1118 for (d = 0; d < spaceDim; ++d) { 1119 coordsB[offB+d] = coordsA[offA+d]; 1120 } 1121 } 1122 if (lc) { /* localized coordinates */ 1123 PetscInt c; 1124 1125 for (c = cS-cStartB; c < cEndB-cStartB; c++) { 1126 PetscInt dof, offA, offB; 1127 1128 ierr = PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);CHKERRQ(ierr); 1129 ierr = PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);CHKERRQ(ierr); 1130 ierr = PetscSectionGetDof(coordSectionA, c + cStartA, &dof);CHKERRQ(ierr); 1131 ierr = PetscMemcpy(coordsB + offB,coordsA + offA,dof*sizeof(*coordsB));CHKERRQ(ierr); 1132 } 1133 } 1134 ierr = VecRestoreArray(coordinatesA, &coordsA);CHKERRQ(ierr); 1135 ierr = VecRestoreArray(coordinatesB, &coordsB);CHKERRQ(ierr); 1136 ierr = DMSetCoordinatesLocal(dmB, coordinatesB);CHKERRQ(ierr); 1137 ierr = VecDestroy(&coordinatesB);CHKERRQ(ierr); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*@ 1142 DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh 1143 1144 Collective on DM 1145 1146 Input Parameter: 1147 . dm - The complete DMPlex object 1148 1149 Output Parameter: 1150 . dmUnint - The DMPlex object with only cells and vertices 1151 1152 Level: intermediate 1153 1154 Notes: 1155 It does not copy over the coordinates. 1156 1157 .keywords: mesh 1158 .seealso: DMPlexInterpolate(), DMPlexCreateFromCellList(), DMPlexCopyCoordinates() 1159 @*/ 1160 PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint) 1161 { 1162 DM udm; 1163 PetscInt dim, vStart, vEnd, cStart, cEnd, cMax, c, maxConeSize = 0, *cone; 1164 PetscErrorCode ierr; 1165 1166 PetscFunctionBegin; 1167 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1168 PetscValidPointer(dmUnint, 2); 1169 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1170 if (dim <= 1) { 1171 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 1172 *dmUnint = dm; 1173 PetscFunctionReturn(0); 1174 } 1175 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1176 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1177 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1178 ierr = DMCreate(PetscObjectComm((PetscObject) dm), &udm);CHKERRQ(ierr); 1179 ierr = DMSetType(udm, DMPLEX);CHKERRQ(ierr); 1180 ierr = DMSetDimension(udm, dim);CHKERRQ(ierr); 1181 ierr = DMPlexSetChart(udm, cStart, vEnd);CHKERRQ(ierr); 1182 for (c = cStart; c < cEnd; ++c) { 1183 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 1184 1185 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1186 for (cl = 0; cl < closureSize*2; cl += 2) { 1187 const PetscInt p = closure[cl]; 1188 1189 if ((p >= vStart) && (p < vEnd)) ++coneSize; 1190 } 1191 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1192 ierr = DMPlexSetConeSize(udm, c, coneSize);CHKERRQ(ierr); 1193 maxConeSize = PetscMax(maxConeSize, coneSize); 1194 } 1195 ierr = DMSetUp(udm);CHKERRQ(ierr); 1196 ierr = PetscMalloc1(maxConeSize, &cone);CHKERRQ(ierr); 1197 for (c = cStart; c < cEnd; ++c) { 1198 PetscInt *closure = NULL, closureSize, cl, coneSize = 0; 1199 1200 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1201 for (cl = 0; cl < closureSize*2; cl += 2) { 1202 const PetscInt p = closure[cl]; 1203 1204 if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p; 1205 } 1206 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1207 ierr = DMPlexSetCone(udm, c, cone);CHKERRQ(ierr); 1208 } 1209 ierr = PetscFree(cone);CHKERRQ(ierr); 1210 ierr = DMPlexSetHybridBounds(udm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1211 ierr = DMPlexSymmetrize(udm);CHKERRQ(ierr); 1212 ierr = DMPlexStratify(udm);CHKERRQ(ierr); 1213 /* Reduce SF */ 1214 { 1215 PetscSF sfPoint, sfPointUn; 1216 const PetscSFNode *remotePoints; 1217 const PetscInt *localPoints; 1218 PetscSFNode *remotePointsUn; 1219 PetscInt *localPointsUn; 1220 PetscInt vEnd, numRoots, numLeaves, l; 1221 PetscInt numLeavesUn = 0, n = 0; 1222 PetscErrorCode ierr; 1223 1224 /* Get original SF information */ 1225 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1226 ierr = DMGetPointSF(udm, &sfPointUn);CHKERRQ(ierr); 1227 ierr = DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);CHKERRQ(ierr); 1228 ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 1229 /* Allocate space for cells and vertices */ 1230 for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++; 1231 /* Fill in leaves */ 1232 if (vEnd >= 0) { 1233 ierr = PetscMalloc1(numLeavesUn, &remotePointsUn);CHKERRQ(ierr); 1234 ierr = PetscMalloc1(numLeavesUn, &localPointsUn);CHKERRQ(ierr); 1235 for (l = 0; l < numLeaves; l++) { 1236 if (localPoints[l] < vEnd) { 1237 localPointsUn[n] = localPoints[l]; 1238 remotePointsUn[n].rank = remotePoints[l].rank; 1239 remotePointsUn[n].index = remotePoints[l].index; 1240 ++n; 1241 } 1242 } 1243 if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn); 1244 ierr = PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);CHKERRQ(ierr); 1245 } 1246 } 1247 { 1248 PetscBool isper; 1249 const PetscReal *maxCell, *L; 1250 const DMBoundaryType *bd; 1251 1252 ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 1253 ierr = DMSetPeriodicity(udm,isper,maxCell,L,bd);CHKERRQ(ierr); 1254 } 1255 1256 *dmUnint = udm; 1257 PetscFunctionReturn(0); 1258 } 1259