1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMPlexGetLineIntersection_2D_Internal" 5 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 6 { 7 const PetscReal p0_x = segmentA[0*2+0]; 8 const PetscReal p0_y = segmentA[0*2+1]; 9 const PetscReal p1_x = segmentA[1*2+0]; 10 const PetscReal p1_y = segmentA[1*2+1]; 11 const PetscReal p2_x = segmentB[0*2+0]; 12 const PetscReal p2_y = segmentB[0*2+1]; 13 const PetscReal p3_x = segmentB[1*2+0]; 14 const PetscReal p3_y = segmentB[1*2+1]; 15 const PetscReal s1_x = p1_x - p0_x; 16 const PetscReal s1_y = p1_y - p0_y; 17 const PetscReal s2_x = p3_x - p2_x; 18 const PetscReal s2_y = p3_y - p2_y; 19 const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 20 21 PetscFunctionBegin; 22 *hasIntersection = PETSC_FALSE; 23 /* Non-parallel lines */ 24 if (denom != 0.0) { 25 const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 26 const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 27 28 if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 29 *hasIntersection = PETSC_TRUE; 30 if (intersection) { 31 intersection[0] = p0_x + (t * s1_x); 32 intersection[1] = p0_y + (t * s1_y); 33 } 34 } 35 } 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal" 41 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 42 { 43 const PetscInt embedDim = 2; 44 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 45 PetscReal x = PetscRealPart(point[0]); 46 PetscReal y = PetscRealPart(point[1]); 47 PetscReal v0[2], J[4], invJ[4], detJ; 48 PetscReal xi, eta; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 53 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 54 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 55 56 if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c; 57 else *cell = -1; 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "DMPlexClosestPoint_Simplex_2D_Internal" 63 static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[]) 64 { 65 const PetscInt embedDim = 2; 66 PetscReal x = PetscRealPart(point[0]); 67 PetscReal y = PetscRealPart(point[1]); 68 PetscReal v0[2], J[4], invJ[4], detJ; 69 PetscReal xi, eta, r; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 74 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 75 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 76 77 xi = PetscMax(xi, 0.0); 78 eta = PetscMax(eta, 0.0); 79 r = (xi + eta)/2.0; 80 if (xi + eta > 2.0) { 81 r = (xi + eta)/2.0; 82 xi /= r; 83 eta /= r; 84 } 85 cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0]; 86 cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1]; 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal" 92 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 93 { 94 PetscSection coordSection; 95 Vec coordsLocal; 96 PetscScalar *coords = NULL; 97 const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 98 PetscReal x = PetscRealPart(point[0]); 99 PetscReal y = PetscRealPart(point[1]); 100 PetscInt crossings = 0, f; 101 PetscErrorCode ierr; 102 103 PetscFunctionBegin; 104 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 105 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 106 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 107 for (f = 0; f < 4; ++f) { 108 PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 109 PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 110 PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 111 PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 112 PetscReal slope = (y_j - y_i) / (x_j - x_i); 113 PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 114 PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 115 PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 116 if ((cond1 || cond2) && above) ++crossings; 117 } 118 if (crossings % 2) *cell = c; 119 else *cell = -1; 120 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 #undef __FUNCT__ 125 #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal" 126 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 127 { 128 const PetscInt embedDim = 3; 129 PetscReal v0[3], J[9], invJ[9], detJ; 130 PetscReal x = PetscRealPart(point[0]); 131 PetscReal y = PetscRealPart(point[1]); 132 PetscReal z = PetscRealPart(point[2]); 133 PetscReal xi, eta, zeta; 134 PetscErrorCode ierr; 135 136 PetscFunctionBegin; 137 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 138 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 139 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 140 zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 141 142 if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 143 else *cell = -1; 144 PetscFunctionReturn(0); 145 } 146 147 #undef __FUNCT__ 148 #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal" 149 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 150 { 151 PetscSection coordSection; 152 Vec coordsLocal; 153 PetscScalar *coords; 154 const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 155 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 156 PetscBool found = PETSC_TRUE; 157 PetscInt f; 158 PetscErrorCode ierr; 159 160 PetscFunctionBegin; 161 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 162 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 163 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 164 for (f = 0; f < 6; ++f) { 165 /* Check the point is under plane */ 166 /* Get face normal */ 167 PetscReal v_i[3]; 168 PetscReal v_j[3]; 169 PetscReal normal[3]; 170 PetscReal pp[3]; 171 PetscReal dot; 172 173 v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 174 v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 175 v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 176 v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 177 v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 178 v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 179 normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 180 normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 181 normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 182 pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 183 pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 184 pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 185 dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 186 187 /* Check that projected point is in face (2D location problem) */ 188 if (dot < 0.0) { 189 found = PETSC_FALSE; 190 break; 191 } 192 } 193 if (found) *cell = c; 194 else *cell = -1; 195 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "PetscGridHashInitialize_Internal" 201 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 202 { 203 PetscInt d; 204 205 PetscFunctionBegin; 206 box->dim = dim; 207 for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNCT__ 212 #define __FUNCT__ "PetscGridHashCreate" 213 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 219 ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "PetscGridHashEnlarge" 225 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 226 { 227 PetscInt d; 228 229 PetscFunctionBegin; 230 for (d = 0; d < box->dim; ++d) { 231 box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 232 box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 233 } 234 PetscFunctionReturn(0); 235 } 236 237 #undef __FUNCT__ 238 #define __FUNCT__ "PetscGridHashSetGrid" 239 /* 240 PetscGridHashSetGrid - Divide the grid into boxes 241 242 Not collective 243 244 Input Parameters: 245 + box - The grid hash object 246 . n - The number of boxes in each dimension, or PETSC_DETERMINE 247 - h - The box size in each dimension, only used if n[d] == PETSC_DETERMINE 248 249 Level: developer 250 251 .seealso: PetscGridHashCreate() 252 */ 253 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 254 { 255 PetscInt d; 256 257 PetscFunctionBegin; 258 for (d = 0; d < box->dim; ++d) { 259 box->extent[d] = box->upper[d] - box->lower[d]; 260 if (n[d] == PETSC_DETERMINE) { 261 box->h[d] = h[d]; 262 box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 263 } else { 264 box->n[d] = n[d]; 265 box->h[d] = box->extent[d]/n[d]; 266 } 267 } 268 PetscFunctionReturn(0); 269 } 270 271 #undef __FUNCT__ 272 #define __FUNCT__ "PetscGridHashGetEnclosingBox" 273 /* 274 PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point 275 276 Not collective 277 278 Input Parameters: 279 + box - The grid hash object 280 . numPoints - The number of input points 281 - points - The input point coordinates 282 283 Output Parameters: 284 + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 285 - boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL 286 287 Level: developer 288 289 .seealso: PetscGridHashCreate() 290 */ 291 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 292 { 293 const PetscReal *lower = box->lower; 294 const PetscReal *upper = box->upper; 295 const PetscReal *h = box->h; 296 const PetscInt *n = box->n; 297 const PetscInt dim = box->dim; 298 PetscInt d, p; 299 300 PetscFunctionBegin; 301 for (p = 0; p < numPoints; ++p) { 302 for (d = 0; d < dim; ++d) { 303 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 304 305 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 306 if (dbox < 0 || dbox >= n[d]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %d (%g, %g, %g) is outside of our bounding box", 307 p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0); 308 dboxes[p*dim+d] = dbox; 309 } 310 if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 311 } 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "PetscGridHashDestroy" 317 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 318 { 319 PetscErrorCode ierr; 320 321 PetscFunctionBegin; 322 if (*box) { 323 ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 324 ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 325 ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 326 } 327 ierr = PetscFree(*box);CHKERRQ(ierr); 328 PetscFunctionReturn(0); 329 } 330 331 #undef __FUNCT__ 332 #define __FUNCT__ "DMPlexLocatePoint_Internal" 333 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 334 { 335 PetscInt coneSize; 336 PetscErrorCode ierr; 337 338 PetscFunctionBegin; 339 switch (dim) { 340 case 2: 341 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 342 switch (coneSize) { 343 case 3: 344 ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 345 break; 346 case 4: 347 ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 348 break; 349 default: 350 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 351 } 352 break; 353 case 3: 354 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 355 switch (coneSize) { 356 case 4: 357 ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 358 break; 359 case 6: 360 ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 361 break; 362 default: 363 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 364 } 365 break; 366 default: 367 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 368 } 369 PetscFunctionReturn(0); 370 } 371 372 #undef __FUNCT__ 373 #define __FUNCT__ "DMPlexClosestPoint_Internal" 374 /* 375 DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point 376 */ 377 PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[]) 378 { 379 PetscInt coneSize; 380 PetscErrorCode ierr; 381 382 PetscFunctionBegin; 383 switch (dim) { 384 case 2: 385 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 386 switch (coneSize) { 387 case 3: 388 ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 389 break; 390 #if 0 391 case 4: 392 ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 393 break; 394 #endif 395 default: 396 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 397 } 398 break; 399 #if 0 400 case 3: 401 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 402 switch (coneSize) { 403 case 4: 404 ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 405 break; 406 case 6: 407 ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 408 break; 409 default: 410 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 411 } 412 break; 413 #endif 414 default: 415 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim); 416 } 417 PetscFunctionReturn(0); 418 } 419 420 #undef __FUNCT__ 421 #define __FUNCT__ "DMPlexComputeGridHash_Internal" 422 /* 423 DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex 424 425 Collective on DM 426 427 Input Parameter: 428 . dm - The Plex 429 430 Output Parameter: 431 . localBox - The grid hash object 432 433 Level: developer 434 435 .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox() 436 */ 437 PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 438 { 439 MPI_Comm comm; 440 PetscGridHash lbox; 441 Vec coordinates; 442 PetscSection coordSection; 443 Vec coordsLocal; 444 const PetscScalar *coords; 445 PetscInt *dboxes, *boxes; 446 PetscInt n[3] = {10, 10, 10}; 447 PetscInt dim, N, cStart, cEnd, cMax, c, i; 448 PetscErrorCode ierr; 449 450 PetscFunctionBegin; 451 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 452 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 453 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 454 if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D"); 455 ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 456 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 457 ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr); 458 for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);} 459 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 460 ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr); 461 #if 0 462 /* Could define a custom reduction to merge these */ 463 ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 464 ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr); 465 #endif 466 /* Is there a reason to snap the local bounding box to a division of the global box? */ 467 /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */ 468 /* Create label */ 469 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 470 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 471 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 472 ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 473 ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 474 /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 475 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 476 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 477 ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr); 478 for (c = cStart; c < cEnd; ++c) { 479 const PetscReal *h = lbox->h; 480 PetscScalar *ccoords = NULL; 481 PetscInt csize = 0; 482 PetscScalar point[3]; 483 PetscInt dlim[6], d, e, i, j, k; 484 485 /* Find boxes enclosing each vertex */ 486 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr); 487 ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr); 488 /* Mark cells containing the vertices */ 489 for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 490 /* Get grid of boxes containing these */ 491 for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 492 for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 493 for (e = 1; e < dim+1; ++e) { 494 for (d = 0; d < dim; ++d) { 495 dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 496 dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 497 } 498 } 499 /* Check for intersection of box with cell */ 500 for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) { 501 for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) { 502 for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) { 503 const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 504 PetscScalar cpoint[3]; 505 PetscInt cell, edge, ii, jj, kk; 506 507 /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 508 for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 509 for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 510 for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 511 512 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 513 if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 514 } 515 } 516 } 517 /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 518 for (edge = 0; edge < dim+1; ++edge) { 519 PetscReal segA[6], segB[6]; 520 521 for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);} 522 for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 523 if (dim > 2) {segB[2] = PetscRealPart(point[2]); 524 segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 525 for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 526 if (dim > 1) {segB[1] = PetscRealPart(point[1]); 527 segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 528 for (ii = 0; ii < 2; ++ii) { 529 PetscBool intersects; 530 531 segB[0] = PetscRealPart(point[0]); 532 segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 533 ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 534 if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 535 } 536 } 537 } 538 } 539 } 540 } 541 } 542 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 543 } 544 ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 545 ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 546 ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 547 *localBox = lbox; 548 PetscFunctionReturn(0); 549 } 550 551 #undef __FUNCT__ 552 #define __FUNCT__ "DMLocatePoints_Plex" 553 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF) 554 { 555 DM_Plex *mesh = (DM_Plex *) dm->data; 556 PetscBool hash = mesh->useHashLocation; 557 PetscInt bs, numPoints, p, numFound, *found = NULL; 558 PetscInt dim, cStart, cEnd, cMax, numCells, c; 559 const PetscInt *boxCells; 560 PetscSFNode *cells; 561 PetscScalar *a; 562 PetscMPIInt result; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing."); 567 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 568 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 569 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr); 570 if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 571 if (bs != dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Block size for point vector %D must be the mesh coordinate dimension %D", bs, dim); 572 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 573 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 574 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 575 ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 576 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 577 numPoints /= bs; 578 ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 579 if (hash) { 580 if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 581 /* Designate the local box for each point */ 582 /* Send points to correct process */ 583 /* Search cells that lie in each subbox */ 584 /* Should we bin points before doing search? */ 585 ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 586 } 587 for (p = 0, numFound = 0; p < numPoints; ++p) { 588 const PetscScalar *point = &a[p*bs]; 589 PetscInt dbin[3], bin, cell = -1, cellOffset; 590 591 cells[p].rank = -1; 592 cells[p].index = -1; 593 if (hash) { 594 ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 595 /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 596 ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 597 ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 598 for (c = cellOffset; c < cellOffset + numCells; ++c) { 599 ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 600 if (cell >= 0) { 601 cells[p].rank = 0; 602 cells[p].index = cell; 603 numFound++; 604 break; 605 } 606 } 607 } else { 608 for (c = cStart; c < cEnd; ++c) { 609 ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 610 if (cell >= 0) { 611 cells[p].rank = 0; 612 cells[p].index = cell; 613 numFound++; 614 break; 615 } 616 } 617 } 618 } 619 if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);} 620 if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) { 621 for (p = 0; p < numPoints; p++) { 622 const PetscScalar *point = &a[p*bs]; 623 PetscReal cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL; 624 PetscInt dbin[3], bin, cellOffset, d; 625 626 if (cells[p].rank < 0) { 627 ++numFound; 628 ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 629 ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 630 ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 631 for (c = cellOffset; c < cellOffset + numCells; ++c) { 632 ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr); 633 for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]); 634 dist = DMPlex_NormD_Internal(dim, diff); 635 if (dist < distMax) { 636 for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d]; 637 cells[p].rank = 0; 638 cells[p].index = boxCells[c]; 639 distMax = dist; 640 break; 641 } 642 } 643 } 644 } 645 } 646 /* This code is only be relevant when interfaced to parallel point location */ 647 /* Check for highest numbered proc that claims a point (do we care?) */ 648 if (numFound < numPoints) { 649 if (ltype == DM_POINTLOCATION_NEAREST) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location does not support parallel point location."); 650 ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr); 651 for (p = 0, numFound = 0; p < numPoints; p++) { 652 if (cells[p].rank >= 0 && cells[p].index >= 0) { 653 if (numFound < p) { 654 cells[numFound] = cells[p]; 655 } 656 found[numFound++] = p; 657 } 658 } 659 } 660 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 661 ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 662 PetscFunctionReturn(0); 663 } 664 665 #undef __FUNCT__ 666 #define __FUNCT__ "DMPlexComputeProjection2Dto1D" 667 /*@C 668 DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates 669 670 Not collective 671 672 Input Parameter: 673 . coords - The coordinates of a segment 674 675 Output Parameters: 676 + coords - The new y-coordinate, and 0 for x 677 - R - The rotation which accomplishes the projection 678 679 Level: developer 680 681 .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D() 682 @*/ 683 PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[]) 684 { 685 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 686 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 687 const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 688 689 PetscFunctionBegin; 690 R[0] = c; R[1] = -s; 691 R[2] = s; R[3] = c; 692 coords[0] = 0.0; 693 coords[1] = r; 694 PetscFunctionReturn(0); 695 } 696 697 #undef __FUNCT__ 698 #define __FUNCT__ "DMPlexComputeProjection3Dto1D" 699 /*@C 700 DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates 701 702 Not collective 703 704 Input Parameter: 705 . coords - The coordinates of a segment 706 707 Output Parameters: 708 + coords - The new y-coordinate, and 0 for x and z 709 - R - The rotation which accomplishes the projection 710 711 Note: This uses the basis completion described by Frisvad in http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html, DOI:10.1080/2165347X.2012.689606 712 713 Level: developer 714 715 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D() 716 @*/ 717 PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[]) 718 { 719 PetscReal x = PetscRealPart(coords[3] - coords[0]); 720 PetscReal y = PetscRealPart(coords[4] - coords[1]); 721 PetscReal z = PetscRealPart(coords[5] - coords[2]); 722 PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 723 PetscReal rinv = 1. / r; 724 PetscFunctionBegin; 725 726 x *= rinv; y *= rinv; z *= rinv; 727 if (x > 0.) { 728 PetscReal inv1pX = 1./ (1. + x); 729 730 R[0] = x; R[1] = -y; R[2] = -z; 731 R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 732 R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 733 } 734 else { 735 PetscReal inv1mX = 1./ (1. - x); 736 737 R[0] = x; R[1] = z; R[2] = y; 738 R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 739 R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 740 } 741 coords[0] = 0.0; 742 coords[1] = r; 743 PetscFunctionReturn(0); 744 } 745 746 #undef __FUNCT__ 747 #define __FUNCT__ "DMPlexComputeProjection3Dto2D" 748 /*@ 749 DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates 750 751 Not collective 752 753 Input Parameter: 754 . coords - The coordinates of a segment 755 756 Output Parameters: 757 + coords - The new y- and z-coordinates, and 0 for x 758 - R - The rotation which accomplishes the projection 759 760 Level: developer 761 762 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D() 763 @*/ 764 PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 765 { 766 PetscReal x1[3], x2[3], n[3], norm; 767 PetscReal x1p[3], x2p[3], xnp[3]; 768 PetscReal sqrtz, alpha; 769 const PetscInt dim = 3; 770 PetscInt d, e, p; 771 772 PetscFunctionBegin; 773 /* 0) Calculate normal vector */ 774 for (d = 0; d < dim; ++d) { 775 x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 776 x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 777 } 778 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 779 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 780 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 781 norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 782 n[0] /= norm; 783 n[1] /= norm; 784 n[2] /= norm; 785 /* 1) Take the normal vector and rotate until it is \hat z 786 787 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 788 789 R = / alpha nx nz alpha ny nz -1/alpha \ 790 | -alpha ny alpha nx 0 | 791 \ nx ny nz / 792 793 will rotate the normal vector to \hat z 794 */ 795 sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 796 /* Check for n = z */ 797 if (sqrtz < 1.0e-10) { 798 const PetscInt s = PetscSign(n[2]); 799 /* If nz < 0, rotate 180 degrees around x-axis */ 800 for (p = 3; p < coordSize/3; ++p) { 801 coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 802 coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s; 803 } 804 coords[0] = 0.0; 805 coords[1] = 0.0; 806 coords[2] = x1[0]; 807 coords[3] = x1[1] * s; 808 coords[4] = x2[0]; 809 coords[5] = x2[1] * s; 810 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 811 R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0; 812 R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s; 813 PetscFunctionReturn(0); 814 } 815 alpha = 1.0/sqrtz; 816 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 817 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 818 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 819 for (d = 0; d < dim; ++d) { 820 x1p[d] = 0.0; 821 x2p[d] = 0.0; 822 for (e = 0; e < dim; ++e) { 823 x1p[d] += R[d*dim+e]*x1[e]; 824 x2p[d] += R[d*dim+e]*x2[e]; 825 } 826 } 827 if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 828 if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 829 /* 2) Project to (x, y) */ 830 for (p = 3; p < coordSize/3; ++p) { 831 for (d = 0; d < dim; ++d) { 832 xnp[d] = 0.0; 833 for (e = 0; e < dim; ++e) { 834 xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 835 } 836 if (d < dim-1) coords[p*2+d] = xnp[d]; 837 } 838 } 839 coords[0] = 0.0; 840 coords[1] = 0.0; 841 coords[2] = x1p[0]; 842 coords[3] = x1p[1]; 843 coords[4] = x2p[0]; 844 coords[5] = x2p[1]; 845 /* Output R^T which rotates \hat z to the input normal */ 846 for (d = 0; d < dim; ++d) { 847 for (e = d+1; e < dim; ++e) { 848 PetscReal tmp; 849 850 tmp = R[d*dim+e]; 851 R[d*dim+e] = R[e*dim+d]; 852 R[e*dim+d] = tmp; 853 } 854 } 855 PetscFunctionReturn(0); 856 } 857 858 #undef __FUNCT__ 859 #define __FUNCT__ "Volume_Triangle_Internal" 860 PETSC_UNUSED 861 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 862 { 863 /* Signed volume is 1/2 the determinant 864 865 | 1 1 1 | 866 | x0 x1 x2 | 867 | y0 y1 y2 | 868 869 but if x0,y0 is the origin, we have 870 871 | x1 x2 | 872 | y1 y2 | 873 */ 874 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 875 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 876 PetscReal M[4], detM; 877 M[0] = x1; M[1] = x2; 878 M[2] = y1; M[3] = y2; 879 DMPlex_Det2D_Internal(&detM, M); 880 *vol = 0.5*detM; 881 (void)PetscLogFlops(5.0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "Volume_Triangle_Origin_Internal" 886 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 887 { 888 DMPlex_Det2D_Internal(vol, coords); 889 *vol *= 0.5; 890 } 891 892 #undef __FUNCT__ 893 #define __FUNCT__ "Volume_Tetrahedron_Internal" 894 PETSC_UNUSED 895 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 896 { 897 /* Signed volume is 1/6th of the determinant 898 899 | 1 1 1 1 | 900 | x0 x1 x2 x3 | 901 | y0 y1 y2 y3 | 902 | z0 z1 z2 z3 | 903 904 but if x0,y0,z0 is the origin, we have 905 906 | x1 x2 x3 | 907 | y1 y2 y3 | 908 | z1 z2 z3 | 909 */ 910 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 911 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 912 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 913 PetscReal M[9], detM; 914 M[0] = x1; M[1] = x2; M[2] = x3; 915 M[3] = y1; M[4] = y2; M[5] = y3; 916 M[6] = z1; M[7] = z2; M[8] = z3; 917 DMPlex_Det3D_Internal(&detM, M); 918 *vol = -0.16666666666666666666666*detM; 919 (void)PetscLogFlops(10.0); 920 } 921 922 #undef __FUNCT__ 923 #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 924 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 925 { 926 DMPlex_Det3D_Internal(vol, coords); 927 *vol *= -0.16666666666666666666666; 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 932 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 933 { 934 PetscSection coordSection; 935 Vec coordinates; 936 PetscScalar *coords = NULL; 937 PetscInt numCoords, d, pStart, pEnd, numSelfCoords = 0; 938 PetscErrorCode ierr; 939 940 PetscFunctionBegin; 941 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 942 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 943 ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 944 if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 945 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 946 numCoords = numSelfCoords ? numSelfCoords : numCoords; 947 if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 948 *detJ = 0.0; 949 if (numCoords == 6) { 950 const PetscInt dim = 3; 951 PetscReal R[9], J0; 952 953 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 954 ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr); 955 if (J) { 956 J0 = 0.5*PetscRealPart(coords[1]); 957 J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 958 J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 959 J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 960 DMPlex_Det3D_Internal(detJ, J); 961 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 962 } 963 } else if (numCoords == 4) { 964 const PetscInt dim = 2; 965 PetscReal R[4], J0; 966 967 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 968 ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr); 969 if (J) { 970 J0 = 0.5*PetscRealPart(coords[1]); 971 J[0] = R[0]*J0; J[1] = R[1]; 972 J[2] = R[2]*J0; J[3] = R[3]; 973 DMPlex_Det2D_Internal(detJ, J); 974 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 975 } 976 } else if (numCoords == 2) { 977 const PetscInt dim = 1; 978 979 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 980 if (J) { 981 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 982 *detJ = J[0]; 983 ierr = PetscLogFlops(2.0);CHKERRQ(ierr); 984 if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);} 985 } 986 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 987 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 988 PetscFunctionReturn(0); 989 } 990 991 #undef __FUNCT__ 992 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 993 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 994 { 995 PetscSection coordSection; 996 Vec coordinates; 997 PetscScalar *coords = NULL; 998 PetscInt numCoords, d, f, g; 999 PetscErrorCode ierr; 1000 1001 PetscFunctionBegin; 1002 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1003 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1004 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1005 *detJ = 0.0; 1006 if (numCoords == 9) { 1007 const PetscInt dim = 3; 1008 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 1009 1010 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1011 ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 1012 if (J) { 1013 const PetscInt pdim = 2; 1014 1015 for (d = 0; d < pdim; d++) { 1016 for (f = 0; f < pdim; f++) { 1017 J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1018 } 1019 } 1020 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1021 DMPlex_Det3D_Internal(detJ, J0); 1022 for (d = 0; d < dim; d++) { 1023 for (f = 0; f < dim; f++) { 1024 J[d*dim+f] = 0.0; 1025 for (g = 0; g < dim; g++) { 1026 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 1027 } 1028 } 1029 } 1030 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1031 } 1032 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1033 } else if (numCoords == 6) { 1034 const PetscInt dim = 2; 1035 1036 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1037 if (J) { 1038 for (d = 0; d < dim; d++) { 1039 for (f = 0; f < dim; f++) { 1040 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 1041 } 1042 } 1043 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1044 DMPlex_Det2D_Internal(detJ, J); 1045 } 1046 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1047 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 1048 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1049 PetscFunctionReturn(0); 1050 } 1051 1052 #undef __FUNCT__ 1053 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 1054 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1055 { 1056 PetscSection coordSection; 1057 Vec coordinates; 1058 PetscScalar *coords = NULL; 1059 PetscInt numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd; 1060 PetscErrorCode ierr; 1061 1062 PetscFunctionBegin; 1063 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1064 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1065 ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 1066 if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 1067 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1068 numCoords = numSelfCoords ? numSelfCoords : numCoords; 1069 *detJ = 0.0; 1070 if (numCoords == 12) { 1071 const PetscInt dim = 3; 1072 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 1073 1074 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1075 ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 1076 if (J) { 1077 const PetscInt pdim = 2; 1078 1079 for (d = 0; d < pdim; d++) { 1080 J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1081 J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1082 } 1083 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1084 DMPlex_Det3D_Internal(detJ, J0); 1085 for (d = 0; d < dim; d++) { 1086 for (f = 0; f < dim; f++) { 1087 J[d*dim+f] = 0.0; 1088 for (g = 0; g < dim; g++) { 1089 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 1090 } 1091 } 1092 } 1093 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1094 } 1095 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1096 } else if (numCoords == 8) { 1097 const PetscInt dim = 2; 1098 1099 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1100 if (J) { 1101 for (d = 0; d < dim; d++) { 1102 J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1103 J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1104 } 1105 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1106 DMPlex_Det2D_Internal(detJ, J); 1107 } 1108 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1109 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 1110 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 1116 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1117 { 1118 PetscSection coordSection; 1119 Vec coordinates; 1120 PetscScalar *coords = NULL; 1121 const PetscInt dim = 3; 1122 PetscInt d; 1123 PetscErrorCode ierr; 1124 1125 PetscFunctionBegin; 1126 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1127 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1128 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1129 *detJ = 0.0; 1130 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1131 if (J) { 1132 for (d = 0; d < dim; d++) { 1133 /* I orient with outward face normals */ 1134 J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 1135 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1136 J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1137 } 1138 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1139 DMPlex_Det3D_Internal(detJ, J); 1140 } 1141 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1142 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1143 PetscFunctionReturn(0); 1144 } 1145 1146 #undef __FUNCT__ 1147 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 1148 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1149 { 1150 PetscSection coordSection; 1151 Vec coordinates; 1152 PetscScalar *coords = NULL; 1153 const PetscInt dim = 3; 1154 PetscInt d; 1155 PetscErrorCode ierr; 1156 1157 PetscFunctionBegin; 1158 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1159 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1160 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1161 *detJ = 0.0; 1162 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1163 if (J) { 1164 for (d = 0; d < dim; d++) { 1165 J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1166 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1167 J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 1168 } 1169 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1170 DMPlex_Det3D_Internal(detJ, J); 1171 } 1172 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1173 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1174 PetscFunctionReturn(0); 1175 } 1176 1177 #undef __FUNCT__ 1178 #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 1179 /*@C 1180 DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 1181 1182 Collective on DM 1183 1184 Input Arguments: 1185 + dm - the DM 1186 - cell - the cell 1187 1188 Output Arguments: 1189 + v0 - the translation part of this affine transform 1190 . J - the Jacobian of the transform from the reference element 1191 . invJ - the inverse of the Jacobian 1192 - detJ - the Jacobian determinant 1193 1194 Level: advanced 1195 1196 Fortran Notes: 1197 Since it returns arrays, this routine is only available in Fortran 90, and you must 1198 include petsc.h90 in your code. 1199 1200 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 1201 @*/ 1202 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1203 { 1204 PetscInt depth, dim, coneSize; 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1209 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1210 if (depth == 1) { 1211 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1212 } else { 1213 DMLabel depth; 1214 1215 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1216 ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 1217 } 1218 switch (dim) { 1219 case 1: 1220 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1221 break; 1222 case 2: 1223 switch (coneSize) { 1224 case 3: 1225 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1226 break; 1227 case 4: 1228 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1229 break; 1230 default: 1231 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1232 } 1233 break; 1234 case 3: 1235 switch (coneSize) { 1236 case 4: 1237 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1238 break; 1239 case 6: /* Faces */ 1240 case 8: /* Vertices */ 1241 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1242 break; 1243 default: 1244 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1245 } 1246 break; 1247 default: 1248 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1249 } 1250 PetscFunctionReturn(0); 1251 } 1252 1253 #undef __FUNCT__ 1254 #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 1255 static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1256 { 1257 PetscQuadrature quad; 1258 PetscSection coordSection; 1259 Vec coordinates; 1260 PetscScalar *coords = NULL; 1261 const PetscReal *quadPoints; 1262 PetscReal *basisDer; 1263 PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 1264 PetscErrorCode ierr; 1265 1266 PetscFunctionBegin; 1267 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1268 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1269 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1270 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1271 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 1272 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1273 ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 1274 ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1275 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 1276 *detJ = 0.0; 1277 if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 1278 if (numCoords != pdim*cdim) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "There are %d coordinates for point %d != %d*%d", numCoords, point, pdim, cdim); 1279 if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 1280 if (J) { 1281 ierr = PetscMemzero(J, Nq*cdim*dim*sizeof(PetscReal));CHKERRQ(ierr); 1282 for (q = 0; q < Nq; ++q) { 1283 PetscInt i, j, k, c, r; 1284 1285 /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 1286 for (k = 0; k < pdim; ++k) 1287 for (j = 0; j < dim; ++j) 1288 for (i = 0; i < cdim; ++i) 1289 J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 1290 ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr); 1291 if (cdim > dim) { 1292 for (c = dim; c < cdim; ++c) 1293 for (r = 0; r < cdim; ++r) 1294 J[r*cdim+c] = r == c ? 1.0 : 0.0; 1295 } 1296 switch (cdim) { 1297 case 3: 1298 DMPlex_Det3D_Internal(detJ, J); 1299 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1300 break; 1301 case 2: 1302 DMPlex_Det2D_Internal(detJ, J); 1303 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1304 break; 1305 case 1: 1306 *detJ = J[0]; 1307 if (invJ) invJ[0] = 1.0/J[0]; 1308 } 1309 } 1310 } 1311 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1312 PetscFunctionReturn(0); 1313 } 1314 1315 #undef __FUNCT__ 1316 #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 1317 /*@C 1318 DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 1319 1320 Collective on DM 1321 1322 Input Arguments: 1323 + dm - the DM 1324 . cell - the cell 1325 - fe - the finite element containing the quadrature 1326 1327 Output Arguments: 1328 + v0 - the translation part of this transform 1329 . J - the Jacobian of the transform from the reference element at each quadrature point 1330 . invJ - the inverse of the Jacobian at each quadrature point 1331 - detJ - the Jacobian determinant at each quadrature point 1332 1333 Level: advanced 1334 1335 Fortran Notes: 1336 Since it returns arrays, this routine is only available in Fortran 90, and you must 1337 include petsc.h90 in your code. 1338 1339 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1340 @*/ 1341 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1342 { 1343 PetscErrorCode ierr; 1344 1345 PetscFunctionBegin; 1346 if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1347 else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1348 PetscFunctionReturn(0); 1349 } 1350 1351 #undef __FUNCT__ 1352 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 1353 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1354 { 1355 PetscSection coordSection; 1356 Vec coordinates; 1357 PetscScalar *coords = NULL; 1358 PetscScalar tmp[2]; 1359 PetscInt coordSize; 1360 PetscErrorCode ierr; 1361 1362 PetscFunctionBegin; 1363 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1364 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1365 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1366 if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 1367 ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1368 if (centroid) { 1369 centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 1370 centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1371 } 1372 if (normal) { 1373 PetscReal norm; 1374 1375 normal[0] = -PetscRealPart(coords[1] - tmp[1]); 1376 normal[1] = PetscRealPart(coords[0] - tmp[0]); 1377 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1378 normal[0] /= norm; 1379 normal[1] /= norm; 1380 } 1381 if (vol) { 1382 *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1383 } 1384 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1385 PetscFunctionReturn(0); 1386 } 1387 1388 #undef __FUNCT__ 1389 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 1390 /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1391 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1392 { 1393 PetscSection coordSection; 1394 Vec coordinates; 1395 PetscScalar *coords = NULL; 1396 PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 1397 PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1398 PetscErrorCode ierr; 1399 1400 PetscFunctionBegin; 1401 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1402 ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 1403 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1404 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1405 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1406 if (dim > 2 && centroid) { 1407 v0[0] = PetscRealPart(coords[0]); 1408 v0[1] = PetscRealPart(coords[1]); 1409 v0[2] = PetscRealPart(coords[2]); 1410 } 1411 if (normal) { 1412 if (dim > 2) { 1413 const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 1414 const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 1415 const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 1416 PetscReal norm; 1417 1418 normal[0] = y0*z1 - z0*y1; 1419 normal[1] = z0*x1 - x0*z1; 1420 normal[2] = x0*y1 - y0*x1; 1421 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 1422 normal[0] /= norm; 1423 normal[1] /= norm; 1424 normal[2] /= norm; 1425 } else { 1426 for (d = 0; d < dim; ++d) normal[d] = 0.0; 1427 } 1428 } 1429 if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);} 1430 for (p = 0; p < numCorners; ++p) { 1431 /* Need to do this copy to get types right */ 1432 for (d = 0; d < tdim; ++d) { 1433 ctmp[d] = PetscRealPart(coords[p*tdim+d]); 1434 ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 1435 } 1436 Volume_Triangle_Origin_Internal(&vtmp, ctmp); 1437 vsum += vtmp; 1438 for (d = 0; d < tdim; ++d) { 1439 csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 1440 } 1441 } 1442 for (d = 0; d < tdim; ++d) { 1443 csum[d] /= (tdim+1)*vsum; 1444 } 1445 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1446 if (vol) *vol = PetscAbsReal(vsum); 1447 if (centroid) { 1448 if (dim > 2) { 1449 for (d = 0; d < dim; ++d) { 1450 centroid[d] = v0[d]; 1451 for (e = 0; e < dim; ++e) { 1452 centroid[d] += R[d*dim+e]*csum[e]; 1453 } 1454 } 1455 } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 1456 } 1457 PetscFunctionReturn(0); 1458 } 1459 1460 #undef __FUNCT__ 1461 #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 1462 /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1463 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1464 { 1465 PetscSection coordSection; 1466 Vec coordinates; 1467 PetscScalar *coords = NULL; 1468 PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1469 const PetscInt *faces, *facesO; 1470 PetscInt numFaces, f, coordSize, numCorners, p, d; 1471 PetscErrorCode ierr; 1472 1473 PetscFunctionBegin; 1474 if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 1475 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1476 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1477 1478 if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 1479 ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 1480 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1481 ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 1482 for (f = 0; f < numFaces; ++f) { 1483 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1484 numCorners = coordSize/dim; 1485 switch (numCorners) { 1486 case 3: 1487 for (d = 0; d < dim; ++d) { 1488 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1489 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1490 coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 1491 } 1492 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1493 if (facesO[f] < 0) vtmp = -vtmp; 1494 vsum += vtmp; 1495 if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 1496 for (d = 0; d < dim; ++d) { 1497 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1498 } 1499 } 1500 break; 1501 case 4: 1502 /* DO FOR PYRAMID */ 1503 /* First tet */ 1504 for (d = 0; d < dim; ++d) { 1505 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1506 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1507 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1508 } 1509 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1510 if (facesO[f] < 0) vtmp = -vtmp; 1511 vsum += vtmp; 1512 if (centroid) { 1513 for (d = 0; d < dim; ++d) { 1514 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1515 } 1516 } 1517 /* Second tet */ 1518 for (d = 0; d < dim; ++d) { 1519 coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 1520 coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 1521 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1522 } 1523 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1524 if (facesO[f] < 0) vtmp = -vtmp; 1525 vsum += vtmp; 1526 if (centroid) { 1527 for (d = 0; d < dim; ++d) { 1528 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1529 } 1530 } 1531 break; 1532 default: 1533 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 1534 } 1535 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1536 } 1537 if (vol) *vol = PetscAbsReal(vsum); 1538 if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1539 if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 1540 PetscFunctionReturn(0); 1541 } 1542 1543 #undef __FUNCT__ 1544 #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1545 /*@C 1546 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1547 1548 Collective on DM 1549 1550 Input Arguments: 1551 + dm - the DM 1552 - cell - the cell 1553 1554 Output Arguments: 1555 + volume - the cell volume 1556 . centroid - the cell centroid 1557 - normal - the cell normal, if appropriate 1558 1559 Level: advanced 1560 1561 Fortran Notes: 1562 Since it returns arrays, this routine is only available in Fortran 90, and you must 1563 include petsc.h90 in your code. 1564 1565 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1566 @*/ 1567 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1568 { 1569 PetscInt depth, dim; 1570 PetscErrorCode ierr; 1571 1572 PetscFunctionBegin; 1573 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1574 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1575 if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1576 /* We need to keep a pointer to the depth label */ 1577 ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1578 /* Cone size is now the number of faces */ 1579 switch (depth) { 1580 case 1: 1581 ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1582 break; 1583 case 2: 1584 ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1585 break; 1586 case 3: 1587 ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1588 break; 1589 default: 1590 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1591 } 1592 PetscFunctionReturn(0); 1593 } 1594 1595 #undef __FUNCT__ 1596 #define __FUNCT__ "DMPlexComputeGeometryFEM" 1597 /* This should also take a PetscFE argument I think */ 1598 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1599 { 1600 DM dmCell; 1601 Vec coordinates; 1602 PetscSection coordSection, sectionCell; 1603 PetscScalar *cgeom; 1604 PetscInt cStart, cEnd, cMax, c; 1605 PetscErrorCode ierr; 1606 1607 PetscFunctionBegin; 1608 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1609 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1610 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1611 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1612 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1613 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1614 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1615 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1616 cEnd = cMax < 0 ? cEnd : cMax; 1617 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1618 /* TODO This needs to be multiplied by Nq for non-affine */ 1619 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1620 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1621 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1622 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1623 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1624 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1625 for (c = cStart; c < cEnd; ++c) { 1626 PetscFECellGeom *cg; 1627 1628 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1629 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1630 ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1631 if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1632 } 1633 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1634 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1635 PetscFunctionReturn(0); 1636 } 1637 1638 #undef __FUNCT__ 1639 #define __FUNCT__ "DMPlexComputeGeometryFVM" 1640 /*@ 1641 DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1642 1643 Input Parameter: 1644 . dm - The DM 1645 1646 Output Parameters: 1647 + cellgeom - A Vec of PetscFVCellGeom data 1648 . facegeom - A Vec of PetscFVFaceGeom data 1649 1650 Level: developer 1651 1652 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1653 @*/ 1654 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1655 { 1656 DM dmFace, dmCell; 1657 DMLabel ghostLabel; 1658 PetscSection sectionFace, sectionCell; 1659 PetscSection coordSection; 1660 Vec coordinates; 1661 PetscScalar *fgeom, *cgeom; 1662 PetscReal minradius, gminradius; 1663 PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1664 PetscErrorCode ierr; 1665 1666 PetscFunctionBegin; 1667 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1668 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1669 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1670 /* Make cell centroids and volumes */ 1671 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1672 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1673 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1674 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1675 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1676 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1677 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1678 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1679 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1680 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1681 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1682 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1683 if (cEndInterior < 0) { 1684 cEndInterior = cEnd; 1685 } 1686 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1687 for (c = cStart; c < cEndInterior; ++c) { 1688 PetscFVCellGeom *cg; 1689 1690 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1691 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1692 ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1693 } 1694 /* Compute face normals and minimum cell radius */ 1695 ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1696 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1697 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1698 ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 1699 for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1700 ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1701 ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1702 ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1703 ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1704 ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1705 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1706 minradius = PETSC_MAX_REAL; 1707 for (f = fStart; f < fEnd; ++f) { 1708 PetscFVFaceGeom *fg; 1709 PetscReal area; 1710 PetscInt ghost = -1, d, numChildren; 1711 1712 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1713 ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 1714 if (ghost >= 0 || numChildren) continue; 1715 ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1716 ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1717 for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1718 /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1719 { 1720 PetscFVCellGeom *cL, *cR; 1721 PetscInt ncells; 1722 const PetscInt *cells; 1723 PetscReal *lcentroid, *rcentroid; 1724 PetscReal l[3], r[3], v[3]; 1725 1726 ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 1727 ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr); 1728 ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1729 lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 1730 if (ncells > 1) { 1731 ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1732 rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 1733 } 1734 else { 1735 rcentroid = fg->centroid; 1736 } 1737 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 1738 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 1739 DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1740 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1741 for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1742 } 1743 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1744 if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]); 1745 if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]); 1746 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1747 } 1748 if (cells[0] < cEndInterior) { 1749 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1750 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1751 } 1752 if (ncells > 1 && cells[1] < cEndInterior) { 1753 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1754 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1755 } 1756 } 1757 } 1758 ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1759 ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1760 /* Compute centroids of ghost cells */ 1761 for (c = cEndInterior; c < cEnd; ++c) { 1762 PetscFVFaceGeom *fg; 1763 const PetscInt *cone, *support; 1764 PetscInt coneSize, supportSize, s; 1765 1766 ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1767 if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1768 ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1769 ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 1770 if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize); 1771 ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1772 ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1773 for (s = 0; s < 2; ++s) { 1774 /* Reflect ghost centroid across plane of face */ 1775 if (support[s] == c) { 1776 PetscFVCellGeom *ci; 1777 PetscFVCellGeom *cg; 1778 PetscReal c2f[3], a; 1779 1780 ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1781 DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1782 a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1783 ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1784 DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1785 cg->volume = ci->volume; 1786 } 1787 } 1788 } 1789 ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1790 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1791 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1792 ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1793 PetscFunctionReturn(0); 1794 } 1795 1796 #undef __FUNCT__ 1797 #define __FUNCT__ "DMPlexGetMinRadius" 1798 /*@C 1799 DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1800 1801 Not collective 1802 1803 Input Argument: 1804 . dm - the DM 1805 1806 Output Argument: 1807 . minradius - the minium cell radius 1808 1809 Level: developer 1810 1811 .seealso: DMGetCoordinates() 1812 @*/ 1813 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1814 { 1815 PetscFunctionBegin; 1816 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1817 PetscValidPointer(minradius,2); 1818 *minradius = ((DM_Plex*) dm->data)->minradius; 1819 PetscFunctionReturn(0); 1820 } 1821 1822 #undef __FUNCT__ 1823 #define __FUNCT__ "DMPlexSetMinRadius" 1824 /*@C 1825 DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1826 1827 Logically collective 1828 1829 Input Arguments: 1830 + dm - the DM 1831 - minradius - the minium cell radius 1832 1833 Level: developer 1834 1835 .seealso: DMSetCoordinates() 1836 @*/ 1837 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1838 { 1839 PetscFunctionBegin; 1840 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1841 ((DM_Plex*) dm->data)->minradius = minradius; 1842 PetscFunctionReturn(0); 1843 } 1844 1845 #undef __FUNCT__ 1846 #define __FUNCT__ "BuildGradientReconstruction_Internal" 1847 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1848 { 1849 DMLabel ghostLabel; 1850 PetscScalar *dx, *grad, **gref; 1851 PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1852 PetscErrorCode ierr; 1853 1854 PetscFunctionBegin; 1855 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1856 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1857 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1858 ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1859 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1860 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1861 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1862 for (c = cStart; c < cEndInterior; c++) { 1863 const PetscInt *faces; 1864 PetscInt numFaces, usedFaces, f, d; 1865 PetscFVCellGeom *cg; 1866 PetscBool boundary; 1867 PetscInt ghost; 1868 1869 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1870 ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1871 ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1872 if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1873 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1874 PetscFVCellGeom *cg1; 1875 PetscFVFaceGeom *fg; 1876 const PetscInt *fcells; 1877 PetscInt ncell, side; 1878 1879 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1880 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1881 if ((ghost >= 0) || boundary) continue; 1882 ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1883 side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1884 ncell = fcells[!side]; /* the neighbor */ 1885 ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1886 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1887 for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1888 gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1889 } 1890 if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1891 ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1892 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1893 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1894 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1895 if ((ghost >= 0) || boundary) continue; 1896 for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1897 ++usedFaces; 1898 } 1899 } 1900 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1901 PetscFunctionReturn(0); 1902 } 1903 1904 #undef __FUNCT__ 1905 #define __FUNCT__ "BuildGradientReconstruction_Internal_Tree" 1906 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1907 { 1908 DMLabel ghostLabel; 1909 PetscScalar *dx, *grad, **gref; 1910 PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 1911 PetscSection neighSec; 1912 PetscInt (*neighbors)[2]; 1913 PetscInt *counter; 1914 PetscErrorCode ierr; 1915 1916 PetscFunctionBegin; 1917 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1918 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1919 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1920 if (cEndInterior < 0) { 1921 cEndInterior = cEnd; 1922 } 1923 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr); 1924 ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr); 1925 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1926 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1927 for (f = fStart; f < fEnd; f++) { 1928 const PetscInt *fcells; 1929 PetscBool boundary; 1930 PetscInt ghost = -1; 1931 PetscInt numChildren, numCells, c; 1932 1933 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1934 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1935 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1936 if ((ghost >= 0) || boundary || numChildren) continue; 1937 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 1938 if (numCells == 2) { 1939 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1940 for (c = 0; c < 2; c++) { 1941 PetscInt cell = fcells[c]; 1942 1943 if (cell >= cStart && cell < cEndInterior) { 1944 ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr); 1945 } 1946 } 1947 } 1948 } 1949 ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr); 1950 ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr); 1951 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1952 nStart = 0; 1953 ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr); 1954 ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr); 1955 ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr); 1956 for (f = fStart; f < fEnd; f++) { 1957 const PetscInt *fcells; 1958 PetscBool boundary; 1959 PetscInt ghost = -1; 1960 PetscInt numChildren, numCells, c; 1961 1962 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1963 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1964 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1965 if ((ghost >= 0) || boundary || numChildren) continue; 1966 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 1967 if (numCells == 2) { 1968 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1969 for (c = 0; c < 2; c++) { 1970 PetscInt cell = fcells[c], off; 1971 1972 if (cell >= cStart && cell < cEndInterior) { 1973 ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr); 1974 off += counter[cell - cStart]++; 1975 neighbors[off][0] = f; 1976 neighbors[off][1] = fcells[1 - c]; 1977 } 1978 } 1979 } 1980 } 1981 ierr = PetscFree(counter);CHKERRQ(ierr); 1982 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1983 for (c = cStart; c < cEndInterior; c++) { 1984 PetscInt numFaces, f, d, off, ghost = -1; 1985 PetscFVCellGeom *cg; 1986 1987 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1988 ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr); 1989 ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr); 1990 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);} 1991 if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1992 for (f = 0; f < numFaces; ++f) { 1993 PetscFVCellGeom *cg1; 1994 PetscFVFaceGeom *fg; 1995 const PetscInt *fcells; 1996 PetscInt ncell, side, nface; 1997 1998 nface = neighbors[off + f][0]; 1999 ncell = neighbors[off + f][1]; 2000 ierr = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr); 2001 side = (c != fcells[0]); 2002 ierr = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr); 2003 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 2004 for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d]; 2005 gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 2006 } 2007 ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr); 2008 for (f = 0; f < numFaces; ++f) { 2009 for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d]; 2010 } 2011 } 2012 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 2013 ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr); 2014 ierr = PetscFree(neighbors);CHKERRQ(ierr); 2015 PetscFunctionReturn(0); 2016 } 2017 2018 #undef __FUNCT__ 2019 #define __FUNCT__ "DMPlexComputeGradientFVM" 2020 /*@ 2021 DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 2022 2023 Collective on DM 2024 2025 Input Arguments: 2026 + dm - The DM 2027 . fvm - The PetscFV 2028 . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 2029 - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 2030 2031 Output Parameters: 2032 + faceGeometry - The geometric factors for gradient calculation are inserted 2033 - dmGrad - The DM describing the layout of gradient data 2034 2035 Level: developer 2036 2037 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 2038 @*/ 2039 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 2040 { 2041 DM dmFace, dmCell; 2042 PetscScalar *fgeom, *cgeom; 2043 PetscSection sectionGrad, parentSection; 2044 PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 2045 PetscErrorCode ierr; 2046 2047 PetscFunctionBegin; 2048 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2049 ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 2050 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2051 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2052 /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 2053 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 2054 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 2055 ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2056 ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2057 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2058 if (!parentSection) { 2059 ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2060 } else { 2061 ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2062 } 2063 ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2064 ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2065 /* Create storage for gradients */ 2066 ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 2067 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 2068 ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 2069 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 2070 ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 2071 ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 2072 ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 2073 PetscFunctionReturn(0); 2074 } 2075 2076 #undef __FUNCT__ 2077 #define __FUNCT__ "DMPlexGetDataFVM" 2078 PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM) 2079 { 2080 PetscObject cellgeomobj, facegeomobj; 2081 PetscErrorCode ierr; 2082 2083 PetscFunctionBegin; 2084 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2085 if (!cellgeomobj) { 2086 Vec cellgeomInt, facegeomInt; 2087 2088 ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr); 2089 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr); 2090 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr); 2091 ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr); 2092 ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr); 2093 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2094 } 2095 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr); 2096 if (cellgeom) *cellgeom = (Vec) cellgeomobj; 2097 if (facegeom) *facegeom = (Vec) facegeomobj; 2098 if (gradDM) { 2099 PetscObject gradobj; 2100 PetscBool computeGradients; 2101 2102 ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr); 2103 if (!computeGradients) { 2104 *gradDM = NULL; 2105 PetscFunctionReturn(0); 2106 } 2107 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2108 if (!gradobj) { 2109 DM dmGradInt; 2110 2111 ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr); 2112 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr); 2113 ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr); 2114 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2115 } 2116 *gradDM = (DM) gradobj; 2117 } 2118 PetscFunctionReturn(0); 2119 } 2120