1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 3 #include <petscblaslapack.h> 4 #include <petsctime.h> 5 6 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 7 { 8 const PetscReal p0_x = segmentA[0*2+0]; 9 const PetscReal p0_y = segmentA[0*2+1]; 10 const PetscReal p1_x = segmentA[1*2+0]; 11 const PetscReal p1_y = segmentA[1*2+1]; 12 const PetscReal p2_x = segmentB[0*2+0]; 13 const PetscReal p2_y = segmentB[0*2+1]; 14 const PetscReal p3_x = segmentB[1*2+0]; 15 const PetscReal p3_y = segmentB[1*2+1]; 16 const PetscReal s1_x = p1_x - p0_x; 17 const PetscReal s1_y = p1_y - p0_y; 18 const PetscReal s2_x = p3_x - p2_x; 19 const PetscReal s2_y = p3_y - p2_y; 20 const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 21 22 PetscFunctionBegin; 23 *hasIntersection = PETSC_FALSE; 24 /* Non-parallel lines */ 25 if (denom != 0.0) { 26 const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 27 const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 28 29 if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 30 *hasIntersection = PETSC_TRUE; 31 if (intersection) { 32 intersection[0] = p0_x + (t * s1_x); 33 intersection[1] = p0_y + (t * s1_y); 34 } 35 } 36 } 37 PetscFunctionReturn(0); 38 } 39 40 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 41 { 42 const PetscInt embedDim = 2; 43 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 44 PetscReal x = PetscRealPart(point[0]); 45 PetscReal y = PetscRealPart(point[1]); 46 PetscReal v0[2], J[4], invJ[4], detJ; 47 PetscReal xi, eta; 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 52 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 53 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 54 55 if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c; 56 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 57 PetscFunctionReturn(0); 58 } 59 60 static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[]) 61 { 62 const PetscInt embedDim = 2; 63 PetscReal x = PetscRealPart(point[0]); 64 PetscReal y = PetscRealPart(point[1]); 65 PetscReal v0[2], J[4], invJ[4], detJ; 66 PetscReal xi, eta, r; 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 71 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 72 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 73 74 xi = PetscMax(xi, 0.0); 75 eta = PetscMax(eta, 0.0); 76 r = (xi + eta)/2.0; 77 if (xi + eta > 2.0) { 78 r = (xi + eta)/2.0; 79 xi /= r; 80 eta /= r; 81 } 82 cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0]; 83 cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1]; 84 PetscFunctionReturn(0); 85 } 86 87 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 88 { 89 PetscSection coordSection; 90 Vec coordsLocal; 91 PetscScalar *coords = NULL; 92 const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 93 PetscReal x = PetscRealPart(point[0]); 94 PetscReal y = PetscRealPart(point[1]); 95 PetscInt crossings = 0, f; 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 100 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 101 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 102 for (f = 0; f < 4; ++f) { 103 PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 104 PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 105 PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 106 PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 107 PetscReal slope = (y_j - y_i) / (x_j - x_i); 108 PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 109 PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 110 PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 111 if ((cond1 || cond2) && above) ++crossings; 112 } 113 if (crossings % 2) *cell = c; 114 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 115 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 116 PetscFunctionReturn(0); 117 } 118 119 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 120 { 121 const PetscInt embedDim = 3; 122 PetscReal v0[3], J[9], invJ[9], detJ; 123 PetscReal x = PetscRealPart(point[0]); 124 PetscReal y = PetscRealPart(point[1]); 125 PetscReal z = PetscRealPart(point[2]); 126 PetscReal xi, eta, zeta; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 131 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 132 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 133 zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 134 135 if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 136 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 137 PetscFunctionReturn(0); 138 } 139 140 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 141 { 142 PetscSection coordSection; 143 Vec coordsLocal; 144 PetscScalar *coords; 145 const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 146 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 147 PetscBool found = PETSC_TRUE; 148 PetscInt f; 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 153 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 154 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 155 for (f = 0; f < 6; ++f) { 156 /* Check the point is under plane */ 157 /* Get face normal */ 158 PetscReal v_i[3]; 159 PetscReal v_j[3]; 160 PetscReal normal[3]; 161 PetscReal pp[3]; 162 PetscReal dot; 163 164 v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 165 v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 166 v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 167 v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 168 v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 169 v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 170 normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 171 normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 172 normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 173 pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 174 pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 175 pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 176 dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 177 178 /* Check that projected point is in face (2D location problem) */ 179 if (dot < 0.0) { 180 found = PETSC_FALSE; 181 break; 182 } 183 } 184 if (found) *cell = c; 185 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 186 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 191 { 192 PetscInt d; 193 194 PetscFunctionBegin; 195 box->dim = dim; 196 for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 197 PetscFunctionReturn(0); 198 } 199 200 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 206 ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 207 PetscFunctionReturn(0); 208 } 209 210 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 211 { 212 PetscInt d; 213 214 PetscFunctionBegin; 215 for (d = 0; d < box->dim; ++d) { 216 box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 217 box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 218 } 219 PetscFunctionReturn(0); 220 } 221 222 /* 223 PetscGridHashSetGrid - Divide the grid into boxes 224 225 Not collective 226 227 Input Parameters: 228 + box - The grid hash object 229 . n - The number of boxes in each dimension, or PETSC_DETERMINE 230 - h - The box size in each dimension, only used if n[d] == PETSC_DETERMINE 231 232 Level: developer 233 234 .seealso: PetscGridHashCreate() 235 */ 236 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 237 { 238 PetscInt d; 239 240 PetscFunctionBegin; 241 for (d = 0; d < box->dim; ++d) { 242 box->extent[d] = box->upper[d] - box->lower[d]; 243 if (n[d] == PETSC_DETERMINE) { 244 box->h[d] = h[d]; 245 box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 246 } else { 247 box->n[d] = n[d]; 248 box->h[d] = box->extent[d]/n[d]; 249 } 250 } 251 PetscFunctionReturn(0); 252 } 253 254 /* 255 PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point 256 257 Not collective 258 259 Input Parameters: 260 + box - The grid hash object 261 . numPoints - The number of input points 262 - points - The input point coordinates 263 264 Output Parameters: 265 + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 266 - boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL 267 268 Level: developer 269 270 .seealso: PetscGridHashCreate() 271 */ 272 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 273 { 274 const PetscReal *lower = box->lower; 275 const PetscReal *upper = box->upper; 276 const PetscReal *h = box->h; 277 const PetscInt *n = box->n; 278 const PetscInt dim = box->dim; 279 PetscInt d, p; 280 281 PetscFunctionBegin; 282 for (p = 0; p < numPoints; ++p) { 283 for (d = 0; d < dim; ++d) { 284 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 285 286 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 287 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", 288 p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0); 289 dboxes[p*dim+d] = dbox; 290 } 291 if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 292 } 293 PetscFunctionReturn(0); 294 } 295 296 /* 297 PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point 298 299 Not collective 300 301 Input Parameters: 302 + box - The grid hash object 303 . numPoints - The number of input points 304 - points - The input point coordinates 305 306 Output Parameters: 307 + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 308 . boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL 309 - found - Flag indicating if point was located within a box 310 311 Level: developer 312 313 .seealso: PetscGridHashGetEnclosingBox() 314 */ 315 PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found) 316 { 317 const PetscReal *lower = box->lower; 318 const PetscReal *upper = box->upper; 319 const PetscReal *h = box->h; 320 const PetscInt *n = box->n; 321 const PetscInt dim = box->dim; 322 PetscInt d, p; 323 324 PetscFunctionBegin; 325 *found = PETSC_FALSE; 326 for (p = 0; p < numPoints; ++p) { 327 for (d = 0; d < dim; ++d) { 328 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 329 330 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 331 if (dbox < 0 || dbox >= n[d]) { 332 PetscFunctionReturn(0); 333 } 334 dboxes[p*dim+d] = dbox; 335 } 336 if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 337 } 338 *found = PETSC_TRUE; 339 PetscFunctionReturn(0); 340 } 341 342 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 343 { 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 if (*box) { 348 ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 349 ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 350 ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 351 } 352 ierr = PetscFree(*box);CHKERRQ(ierr); 353 PetscFunctionReturn(0); 354 } 355 356 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 357 { 358 PetscInt coneSize; 359 PetscErrorCode ierr; 360 361 PetscFunctionBegin; 362 switch (dim) { 363 case 2: 364 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 365 switch (coneSize) { 366 case 3: 367 ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 368 break; 369 case 4: 370 ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 371 break; 372 default: 373 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 374 } 375 break; 376 case 3: 377 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 378 switch (coneSize) { 379 case 4: 380 ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 381 break; 382 case 6: 383 ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 384 break; 385 default: 386 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 387 } 388 break; 389 default: 390 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 391 } 392 PetscFunctionReturn(0); 393 } 394 395 /* 396 DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point 397 */ 398 PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[]) 399 { 400 PetscInt coneSize; 401 PetscErrorCode ierr; 402 403 PetscFunctionBegin; 404 switch (dim) { 405 case 2: 406 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 407 switch (coneSize) { 408 case 3: 409 ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 410 break; 411 #if 0 412 case 4: 413 ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 414 break; 415 #endif 416 default: 417 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 418 } 419 break; 420 #if 0 421 case 3: 422 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 423 switch (coneSize) { 424 case 4: 425 ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 426 break; 427 case 6: 428 ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 429 break; 430 default: 431 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 432 } 433 break; 434 #endif 435 default: 436 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim); 437 } 438 PetscFunctionReturn(0); 439 } 440 441 /* 442 DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex 443 444 Collective on DM 445 446 Input Parameter: 447 . dm - The Plex 448 449 Output Parameter: 450 . localBox - The grid hash object 451 452 Level: developer 453 454 .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox() 455 */ 456 PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 457 { 458 MPI_Comm comm; 459 PetscGridHash lbox; 460 Vec coordinates; 461 PetscSection coordSection; 462 Vec coordsLocal; 463 const PetscScalar *coords; 464 PetscInt *dboxes, *boxes; 465 PetscInt n[3] = {10, 10, 10}; 466 PetscInt dim, N, cStart, cEnd, cMax, c, i; 467 PetscErrorCode ierr; 468 469 PetscFunctionBegin; 470 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 471 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 472 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 473 if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D"); 474 ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 475 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 476 ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr); 477 for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);} 478 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 479 ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr); 480 #if 0 481 /* Could define a custom reduction to merge these */ 482 ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 483 ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr); 484 #endif 485 /* Is there a reason to snap the local bounding box to a division of the global box? */ 486 /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */ 487 /* Create label */ 488 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 489 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 490 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 491 ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 492 ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 493 /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 494 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 495 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 496 ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr); 497 for (c = cStart; c < cEnd; ++c) { 498 const PetscReal *h = lbox->h; 499 PetscScalar *ccoords = NULL; 500 PetscInt csize = 0; 501 PetscScalar point[3]; 502 PetscInt dlim[6], d, e, i, j, k; 503 504 /* Find boxes enclosing each vertex */ 505 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr); 506 ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr); 507 /* Mark cells containing the vertices */ 508 for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 509 /* Get grid of boxes containing these */ 510 for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 511 for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 512 for (e = 1; e < dim+1; ++e) { 513 for (d = 0; d < dim; ++d) { 514 dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 515 dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 516 } 517 } 518 /* Check for intersection of box with cell */ 519 for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) { 520 for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) { 521 for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) { 522 const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 523 PetscScalar cpoint[3]; 524 PetscInt cell, edge, ii, jj, kk; 525 526 /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 527 for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 528 for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 529 for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 530 531 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 532 if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 533 } 534 } 535 } 536 /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 537 for (edge = 0; edge < dim+1; ++edge) { 538 PetscReal segA[6], segB[6]; 539 540 for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);} 541 for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 542 if (dim > 2) {segB[2] = PetscRealPart(point[2]); 543 segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 544 for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 545 if (dim > 1) {segB[1] = PetscRealPart(point[1]); 546 segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 547 for (ii = 0; ii < 2; ++ii) { 548 PetscBool intersects; 549 550 segB[0] = PetscRealPart(point[0]); 551 segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 552 ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 553 if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 554 } 555 } 556 } 557 } 558 } 559 } 560 } 561 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 562 } 563 ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 564 ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 565 ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 566 *localBox = lbox; 567 PetscFunctionReturn(0); 568 } 569 570 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF) 571 { 572 DM_Plex *mesh = (DM_Plex *) dm->data; 573 PetscBool hash = mesh->useHashLocation, reuse = PETSC_FALSE; 574 PetscInt bs, numPoints, p, numFound, *found = NULL; 575 PetscInt dim, cStart, cEnd, cMax, numCells, c; 576 const PetscInt *boxCells; 577 PetscSFNode *cells; 578 PetscScalar *a; 579 PetscMPIInt result; 580 PetscLogDouble t0,t1; 581 PetscErrorCode ierr; 582 583 PetscFunctionBegin; 584 ierr = PetscTime(&t0);CHKERRQ(ierr); 585 if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it."); 586 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 587 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 588 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr); 589 if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 590 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); 591 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 592 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 593 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 594 ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 595 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 596 numPoints /= bs; 597 { 598 const PetscSFNode *sf_cells; 599 600 ierr = PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);CHKERRQ(ierr); 601 if (sf_cells) { 602 ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");CHKERRQ(ierr); 603 cells = (PetscSFNode*)sf_cells; 604 reuse = PETSC_TRUE; 605 } else { 606 ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");CHKERRQ(ierr); 607 ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 608 /* initialize cells if created */ 609 for (p=0; p<numPoints; p++) { 610 cells[p].rank = 0; 611 cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 612 } 613 } 614 } 615 if (hash) { 616 if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 617 /* Designate the local box for each point */ 618 /* Send points to correct process */ 619 /* Search cells that lie in each subbox */ 620 /* Should we bin points before doing search? */ 621 ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 622 } 623 for (p = 0, numFound = 0; p < numPoints; ++p) { 624 const PetscScalar *point = &a[p*bs]; 625 PetscInt dbin[3], bin, cell = -1, cellOffset; 626 627 /* check initial values in cells[].index - abort early if found */ 628 if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 629 c = cells[p].index; 630 cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 631 ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 632 if (cell >= 0) { 633 cells[p].rank = 0; 634 cells[p].index = cell; 635 numFound++; 636 } 637 } 638 if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { continue; } 639 640 if (hash) { 641 PetscBool found_box; 642 643 //ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 644 /* allow for case that point is outside box - abort early */ 645 ierr = PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);CHKERRQ(ierr); 646 if (found_box) { 647 /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 648 ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 649 ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 650 for (c = cellOffset; c < cellOffset + numCells; ++c) { 651 ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 652 if (cell >= 0) { 653 cells[p].rank = 0; 654 cells[p].index = cell; 655 numFound++; 656 break; 657 } 658 } 659 } 660 } else { 661 for (c = cStart; c < cEnd; ++c) { 662 ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 663 if (cell >= 0) { 664 cells[p].rank = 0; 665 cells[p].index = cell; 666 numFound++; 667 break; 668 } 669 } 670 } 671 } 672 if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);} 673 if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) { 674 for (p = 0; p < numPoints; p++) { 675 const PetscScalar *point = &a[p*bs]; 676 PetscReal cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL; 677 PetscInt dbin[3], bin, cellOffset, d; 678 679 if (cells[p].index < 0) { 680 ++numFound; 681 ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 682 ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 683 ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 684 for (c = cellOffset; c < cellOffset + numCells; ++c) { 685 ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr); 686 for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]); 687 dist = DMPlex_NormD_Internal(dim, diff); 688 if (dist < distMax) { 689 for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d]; 690 cells[p].rank = 0; 691 cells[p].index = boxCells[c]; 692 distMax = dist; 693 break; 694 } 695 } 696 } 697 } 698 } 699 /* This code is only be relevant when interfaced to parallel point location */ 700 /* Check for highest numbered proc that claims a point (do we care?) */ 701 if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) { 702 ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr); 703 for (p = 0, numFound = 0; p < numPoints; p++) { 704 if (cells[p].rank >= 0 && cells[p].index >= 0) { 705 if (numFound < p) { 706 cells[numFound] = cells[p]; 707 } 708 found[numFound++] = p; 709 } 710 } 711 } 712 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 713 if (!reuse) { 714 ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 715 } 716 ierr = PetscTime(&t1);CHKERRQ(ierr); 717 ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] npoints %D : time(rank0) %1.2e (sec): points/sec %1.4e\n",numPoints,t1-t0,(double)((double)numPoints/(t1-t0)));CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720 721 /*@C 722 DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates 723 724 Not collective 725 726 Input Parameter: 727 . coords - The coordinates of a segment 728 729 Output Parameters: 730 + coords - The new y-coordinate, and 0 for x 731 - R - The rotation which accomplishes the projection 732 733 Level: developer 734 735 .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D() 736 @*/ 737 PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[]) 738 { 739 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 740 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 741 const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 742 743 PetscFunctionBegin; 744 R[0] = c; R[1] = -s; 745 R[2] = s; R[3] = c; 746 coords[0] = 0.0; 747 coords[1] = r; 748 PetscFunctionReturn(0); 749 } 750 751 /*@C 752 DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates 753 754 Not collective 755 756 Input Parameter: 757 . coords - The coordinates of a segment 758 759 Output Parameters: 760 + coords - The new y-coordinate, and 0 for x and z 761 - R - The rotation which accomplishes the projection 762 763 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 764 765 Level: developer 766 767 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D() 768 @*/ 769 PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[]) 770 { 771 PetscReal x = PetscRealPart(coords[3] - coords[0]); 772 PetscReal y = PetscRealPart(coords[4] - coords[1]); 773 PetscReal z = PetscRealPart(coords[5] - coords[2]); 774 PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 775 PetscReal rinv = 1. / r; 776 PetscFunctionBegin; 777 778 x *= rinv; y *= rinv; z *= rinv; 779 if (x > 0.) { 780 PetscReal inv1pX = 1./ (1. + x); 781 782 R[0] = x; R[1] = -y; R[2] = -z; 783 R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 784 R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 785 } 786 else { 787 PetscReal inv1mX = 1./ (1. - x); 788 789 R[0] = x; R[1] = z; R[2] = y; 790 R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 791 R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 792 } 793 coords[0] = 0.0; 794 coords[1] = r; 795 PetscFunctionReturn(0); 796 } 797 798 /*@ 799 DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates 800 801 Not collective 802 803 Input Parameter: 804 . coords - The coordinates of a segment 805 806 Output Parameters: 807 + coords - The new y- and z-coordinates, and 0 for x 808 - R - The rotation which accomplishes the projection 809 810 Level: developer 811 812 .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D() 813 @*/ 814 PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 815 { 816 PetscReal x1[3], x2[3], n[3], norm; 817 PetscReal x1p[3], x2p[3], xnp[3]; 818 PetscReal sqrtz, alpha; 819 const PetscInt dim = 3; 820 PetscInt d, e, p; 821 822 PetscFunctionBegin; 823 /* 0) Calculate normal vector */ 824 for (d = 0; d < dim; ++d) { 825 x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 826 x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 827 } 828 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 829 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 830 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 831 norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 832 n[0] /= norm; 833 n[1] /= norm; 834 n[2] /= norm; 835 /* 1) Take the normal vector and rotate until it is \hat z 836 837 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 838 839 R = / alpha nx nz alpha ny nz -1/alpha \ 840 | -alpha ny alpha nx 0 | 841 \ nx ny nz / 842 843 will rotate the normal vector to \hat z 844 */ 845 sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 846 /* Check for n = z */ 847 if (sqrtz < 1.0e-10) { 848 const PetscInt s = PetscSign(n[2]); 849 /* If nz < 0, rotate 180 degrees around x-axis */ 850 for (p = 3; p < coordSize/3; ++p) { 851 coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 852 coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s; 853 } 854 coords[0] = 0.0; 855 coords[1] = 0.0; 856 coords[2] = x1[0]; 857 coords[3] = x1[1] * s; 858 coords[4] = x2[0]; 859 coords[5] = x2[1] * s; 860 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 861 R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0; 862 R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s; 863 PetscFunctionReturn(0); 864 } 865 alpha = 1.0/sqrtz; 866 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 867 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 868 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 869 for (d = 0; d < dim; ++d) { 870 x1p[d] = 0.0; 871 x2p[d] = 0.0; 872 for (e = 0; e < dim; ++e) { 873 x1p[d] += R[d*dim+e]*x1[e]; 874 x2p[d] += R[d*dim+e]*x2[e]; 875 } 876 } 877 if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 878 if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 879 /* 2) Project to (x, y) */ 880 for (p = 3; p < coordSize/3; ++p) { 881 for (d = 0; d < dim; ++d) { 882 xnp[d] = 0.0; 883 for (e = 0; e < dim; ++e) { 884 xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 885 } 886 if (d < dim-1) coords[p*2+d] = xnp[d]; 887 } 888 } 889 coords[0] = 0.0; 890 coords[1] = 0.0; 891 coords[2] = x1p[0]; 892 coords[3] = x1p[1]; 893 coords[4] = x2p[0]; 894 coords[5] = x2p[1]; 895 /* Output R^T which rotates \hat z to the input normal */ 896 for (d = 0; d < dim; ++d) { 897 for (e = d+1; e < dim; ++e) { 898 PetscReal tmp; 899 900 tmp = R[d*dim+e]; 901 R[d*dim+e] = R[e*dim+d]; 902 R[e*dim+d] = tmp; 903 } 904 } 905 PetscFunctionReturn(0); 906 } 907 908 PETSC_UNUSED 909 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 910 { 911 /* Signed volume is 1/2 the determinant 912 913 | 1 1 1 | 914 | x0 x1 x2 | 915 | y0 y1 y2 | 916 917 but if x0,y0 is the origin, we have 918 919 | x1 x2 | 920 | y1 y2 | 921 */ 922 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 923 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 924 PetscReal M[4], detM; 925 M[0] = x1; M[1] = x2; 926 M[2] = y1; M[3] = y2; 927 DMPlex_Det2D_Internal(&detM, M); 928 *vol = 0.5*detM; 929 (void)PetscLogFlops(5.0); 930 } 931 932 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 933 { 934 DMPlex_Det2D_Internal(vol, coords); 935 *vol *= 0.5; 936 } 937 938 PETSC_UNUSED 939 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 940 { 941 /* Signed volume is 1/6th of the determinant 942 943 | 1 1 1 1 | 944 | x0 x1 x2 x3 | 945 | y0 y1 y2 y3 | 946 | z0 z1 z2 z3 | 947 948 but if x0,y0,z0 is the origin, we have 949 950 | x1 x2 x3 | 951 | y1 y2 y3 | 952 | z1 z2 z3 | 953 */ 954 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 955 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 956 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 957 PetscReal M[9], detM; 958 M[0] = x1; M[1] = x2; M[2] = x3; 959 M[3] = y1; M[4] = y2; M[5] = y3; 960 M[6] = z1; M[7] = z2; M[8] = z3; 961 DMPlex_Det3D_Internal(&detM, M); 962 *vol = -0.16666666666666666666666*detM; 963 (void)PetscLogFlops(10.0); 964 } 965 966 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 967 { 968 DMPlex_Det3D_Internal(vol, coords); 969 *vol *= -0.16666666666666666666666; 970 } 971 972 static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 973 { 974 PetscSection coordSection; 975 Vec coordinates; 976 const PetscScalar *coords; 977 PetscInt dim, d, off; 978 PetscErrorCode ierr; 979 980 PetscFunctionBegin; 981 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 982 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 983 ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr); 984 if (!dim) PetscFunctionReturn(0); 985 ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr); 986 ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr); 987 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);} 988 ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr); 989 *detJ = 1.; 990 if (J) { 991 for (d = 0; d < dim * dim; d++) J[d] = 0.; 992 for (d = 0; d < dim; d++) J[d * dim + d] = 1.; 993 if (invJ) { 994 for (d = 0; d < dim * dim; d++) invJ[d] = 0.; 995 for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.; 996 } 997 } 998 PetscFunctionReturn(0); 999 } 1000 1001 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1002 { 1003 PetscSection coordSection; 1004 Vec coordinates; 1005 PetscScalar *coords = NULL; 1006 PetscInt numCoords, d, pStart, pEnd, numSelfCoords = 0; 1007 PetscErrorCode ierr; 1008 1009 PetscFunctionBegin; 1010 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1011 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1012 ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 1013 if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 1014 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1015 numCoords = numSelfCoords ? numSelfCoords : numCoords; 1016 if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 1017 *detJ = 0.0; 1018 if (numCoords == 6) { 1019 const PetscInt dim = 3; 1020 PetscReal R[9], J0; 1021 1022 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1023 ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr); 1024 if (J) { 1025 J0 = 0.5*PetscRealPart(coords[1]); 1026 J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 1027 J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 1028 J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 1029 DMPlex_Det3D_Internal(detJ, J); 1030 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1031 } 1032 } else if (numCoords == 4) { 1033 const PetscInt dim = 2; 1034 PetscReal R[4], J0; 1035 1036 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1037 ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr); 1038 if (J) { 1039 J0 = 0.5*PetscRealPart(coords[1]); 1040 J[0] = R[0]*J0; J[1] = R[1]; 1041 J[2] = R[2]*J0; J[3] = R[3]; 1042 DMPlex_Det2D_Internal(detJ, J); 1043 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1044 } 1045 } else if (numCoords == 2) { 1046 const PetscInt dim = 1; 1047 1048 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1049 if (J) { 1050 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 1051 *detJ = J[0]; 1052 ierr = PetscLogFlops(2.0);CHKERRQ(ierr); 1053 if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);} 1054 } 1055 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 1056 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 1060 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1061 { 1062 PetscSection coordSection; 1063 Vec coordinates; 1064 PetscScalar *coords = NULL; 1065 PetscInt numCoords, d, f, g; 1066 PetscErrorCode ierr; 1067 1068 PetscFunctionBegin; 1069 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1070 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1071 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1072 *detJ = 0.0; 1073 if (numCoords == 9) { 1074 const PetscInt dim = 3; 1075 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 1076 1077 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1078 ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 1079 if (J) { 1080 const PetscInt pdim = 2; 1081 1082 for (d = 0; d < pdim; d++) { 1083 for (f = 0; f < pdim; f++) { 1084 J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1085 } 1086 } 1087 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1088 DMPlex_Det3D_Internal(detJ, J0); 1089 for (d = 0; d < dim; d++) { 1090 for (f = 0; f < dim; f++) { 1091 J[d*dim+f] = 0.0; 1092 for (g = 0; g < dim; g++) { 1093 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 1094 } 1095 } 1096 } 1097 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1098 } 1099 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1100 } else if (numCoords == 6) { 1101 const PetscInt dim = 2; 1102 1103 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1104 if (J) { 1105 for (d = 0; d < dim; d++) { 1106 for (f = 0; f < dim; f++) { 1107 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 1108 } 1109 } 1110 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1111 DMPlex_Det2D_Internal(detJ, J); 1112 } 1113 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1114 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 1115 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1116 PetscFunctionReturn(0); 1117 } 1118 1119 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1120 { 1121 PetscSection coordSection; 1122 Vec coordinates; 1123 PetscScalar *coords = NULL; 1124 PetscInt numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd; 1125 PetscErrorCode ierr; 1126 1127 PetscFunctionBegin; 1128 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1129 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1130 ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 1131 if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 1132 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1133 numCoords = numSelfCoords ? numSelfCoords : numCoords; 1134 if (!Nq) { 1135 *detJ = 0.0; 1136 if (numCoords == 12) { 1137 const PetscInt dim = 3; 1138 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 1139 1140 if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1141 ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 1142 if (J) { 1143 const PetscInt pdim = 2; 1144 1145 for (d = 0; d < pdim; d++) { 1146 J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1147 J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1148 } 1149 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1150 DMPlex_Det3D_Internal(detJ, J0); 1151 for (d = 0; d < dim; d++) { 1152 for (f = 0; f < dim; f++) { 1153 J[d*dim+f] = 0.0; 1154 for (g = 0; g < dim; g++) { 1155 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 1156 } 1157 } 1158 } 1159 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1160 } 1161 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1162 } else if (numCoords == 8) { 1163 const PetscInt dim = 2; 1164 1165 if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1166 if (J) { 1167 for (d = 0; d < dim; d++) { 1168 J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1169 J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1170 } 1171 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1172 DMPlex_Det2D_Internal(detJ, J); 1173 } 1174 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1175 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 1176 } else { 1177 const PetscInt Nv = 4; 1178 const PetscInt dimR = 2; 1179 const PetscInt zToPlex[4] = {0, 1, 3, 2}; 1180 PetscReal zOrder[12]; 1181 PetscReal zCoeff[12]; 1182 PetscInt i, j, k, l, dim; 1183 1184 if (numCoords == 12) { 1185 dim = 3; 1186 } else if (numCoords == 8) { 1187 dim = 2; 1188 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 1189 for (i = 0; i < Nv; i++) { 1190 PetscInt zi = zToPlex[i]; 1191 1192 for (j = 0; j < dim; j++) { 1193 zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 1194 } 1195 } 1196 for (j = 0; j < dim; j++) { 1197 zCoeff[dim * 0 + j] = 0.25 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1198 zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1199 zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1200 zCoeff[dim * 3 + j] = 0.25 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1201 } 1202 for (i = 0; i < Nq; i++) { 1203 PetscReal xi = points[dimR * i], eta = points[dimR * i + 1]; 1204 1205 if (v) { 1206 PetscReal extPoint[4]; 1207 1208 extPoint[0] = 1.; 1209 extPoint[1] = xi; 1210 extPoint[2] = eta; 1211 extPoint[3] = xi * eta; 1212 for (j = 0; j < dim; j++) { 1213 PetscReal val = 0.; 1214 1215 for (k = 0; k < Nv; k++) { 1216 val += extPoint[k] * zCoeff[dim * k + j]; 1217 } 1218 v[i * dim + j] = val; 1219 } 1220 } 1221 if (J) { 1222 PetscReal extJ[8]; 1223 1224 extJ[0] = 0.; 1225 extJ[1] = 0.; 1226 extJ[2] = 1.; 1227 extJ[3] = 0.; 1228 extJ[4] = 0.; 1229 extJ[5] = 1.; 1230 extJ[6] = eta; 1231 extJ[7] = xi; 1232 for (j = 0; j < dim; j++) { 1233 for (k = 0; k < dimR; k++) { 1234 PetscReal val = 0.; 1235 1236 for (l = 0; l < Nv; l++) { 1237 val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 1238 } 1239 J[i * dim * dim + dim * j + k] = val; 1240 } 1241 } 1242 if (dim == 3) { /* put the cross product in the third component of the Jacobian */ 1243 PetscReal x, y, z; 1244 PetscReal *iJ = &J[i * dim * dim]; 1245 PetscReal norm; 1246 1247 x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0]; 1248 y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1]; 1249 z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0]; 1250 norm = PetscSqrtReal(x * x + y * y + z * z); 1251 iJ[2] = x / norm; 1252 iJ[5] = y / norm; 1253 iJ[8] = z / norm; 1254 DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 1255 if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1256 } else { 1257 DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]); 1258 if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1259 } 1260 } 1261 } 1262 } 1263 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1264 PetscFunctionReturn(0); 1265 } 1266 1267 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1268 { 1269 PetscSection coordSection; 1270 Vec coordinates; 1271 PetscScalar *coords = NULL; 1272 const PetscInt dim = 3; 1273 PetscInt d; 1274 PetscErrorCode ierr; 1275 1276 PetscFunctionBegin; 1277 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1278 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1279 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1280 *detJ = 0.0; 1281 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1282 if (J) { 1283 for (d = 0; d < dim; d++) { 1284 /* I orient with outward face normals */ 1285 J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 1286 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1287 J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1288 } 1289 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1290 DMPlex_Det3D_Internal(detJ, J); 1291 } 1292 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1293 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1294 PetscFunctionReturn(0); 1295 } 1296 1297 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1298 { 1299 PetscSection coordSection; 1300 Vec coordinates; 1301 PetscScalar *coords = NULL; 1302 const PetscInt dim = 3; 1303 PetscInt d; 1304 PetscErrorCode ierr; 1305 1306 PetscFunctionBegin; 1307 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1308 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1309 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1310 if (!Nq) { 1311 *detJ = 0.0; 1312 if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1313 if (J) { 1314 for (d = 0; d < dim; d++) { 1315 J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1316 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1317 J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 1318 } 1319 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1320 DMPlex_Det3D_Internal(detJ, J); 1321 } 1322 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1323 } else { 1324 const PetscInt Nv = 8; 1325 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 1326 const PetscInt dim = 3; 1327 const PetscInt dimR = 3; 1328 PetscReal zOrder[24]; 1329 PetscReal zCoeff[24]; 1330 PetscInt i, j, k, l; 1331 1332 for (i = 0; i < Nv; i++) { 1333 PetscInt zi = zToPlex[i]; 1334 1335 for (j = 0; j < dim; j++) { 1336 zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 1337 } 1338 } 1339 for (j = 0; j < dim; j++) { 1340 zCoeff[dim * 0 + j] = 0.125 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1341 zCoeff[dim * 1 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1342 zCoeff[dim * 2 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1343 zCoeff[dim * 3 + j] = 0.125 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1344 zCoeff[dim * 4 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1345 zCoeff[dim * 5 + j] = 0.125 * (+ zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1346 zCoeff[dim * 6 + j] = 0.125 * (+ zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1347 zCoeff[dim * 7 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1348 } 1349 for (i = 0; i < Nq; i++) { 1350 PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2]; 1351 1352 if (v) { 1353 PetscReal extPoint[8]; 1354 1355 extPoint[0] = 1.; 1356 extPoint[1] = xi; 1357 extPoint[2] = eta; 1358 extPoint[3] = xi * eta; 1359 extPoint[4] = theta; 1360 extPoint[5] = theta * xi; 1361 extPoint[6] = theta * eta; 1362 extPoint[7] = theta * eta * xi; 1363 for (j = 0; j < dim; j++) { 1364 PetscReal val = 0.; 1365 1366 for (k = 0; k < Nv; k++) { 1367 val += extPoint[k] * zCoeff[dim * k + j]; 1368 } 1369 v[i * dim + j] = val; 1370 } 1371 } 1372 if (J) { 1373 PetscReal extJ[24]; 1374 1375 extJ[0] = 0. ; extJ[1] = 0. ; extJ[2] = 0. ; 1376 extJ[3] = 1. ; extJ[4] = 0. ; extJ[5] = 0. ; 1377 extJ[6] = 0. ; extJ[7] = 1. ; extJ[8] = 0. ; 1378 extJ[9] = eta ; extJ[10] = xi ; extJ[11] = 0. ; 1379 extJ[12] = 0. ; extJ[13] = 0. ; extJ[14] = 1. ; 1380 extJ[15] = theta ; extJ[16] = 0. ; extJ[17] = xi ; 1381 extJ[18] = 0. ; extJ[19] = theta ; extJ[20] = eta ; 1382 extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi; 1383 1384 for (j = 0; j < dim; j++) { 1385 for (k = 0; k < dimR; k++) { 1386 PetscReal val = 0.; 1387 1388 for (l = 0; l < Nv; l++) { 1389 val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 1390 } 1391 J[i * dim * dim + dim * j + k] = val; 1392 } 1393 } 1394 DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 1395 if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1396 } 1397 } 1398 } 1399 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1400 PetscFunctionReturn(0); 1401 } 1402 1403 static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1404 { 1405 PetscInt depth, dim, coordDim, coneSize, i; 1406 PetscInt Nq = 0; 1407 const PetscReal *points = NULL; 1408 DMLabel depthLabel; 1409 PetscReal v0[3] = {-1.}, J0[9] = {-1.}, detJ0 = -1.; 1410 PetscBool isAffine = PETSC_TRUE; 1411 PetscErrorCode ierr; 1412 1413 PetscFunctionBegin; 1414 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1415 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1416 ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1417 ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr); 1418 if (depth == 1 && dim == 1) { 1419 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1420 } 1421 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 1422 if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim); 1423 if (quad) {ierr = PetscQuadratureGetData(quad, NULL, &Nq, &points, NULL);CHKERRQ(ierr);} 1424 switch (dim) { 1425 case 0: 1426 ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 1427 isAffine = PETSC_FALSE; 1428 break; 1429 case 1: 1430 if (Nq) { 1431 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 1432 } else { 1433 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 1434 } 1435 break; 1436 case 2: 1437 switch (coneSize) { 1438 case 3: 1439 if (Nq) { 1440 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 1441 } else { 1442 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 1443 } 1444 break; 1445 case 4: 1446 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr); 1447 isAffine = PETSC_FALSE; 1448 break; 1449 default: 1450 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1451 } 1452 break; 1453 case 3: 1454 switch (coneSize) { 1455 case 4: 1456 if (Nq) { 1457 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 1458 } else { 1459 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 1460 } 1461 break; 1462 case 6: /* Faces */ 1463 case 8: /* Vertices */ 1464 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr); 1465 isAffine = PETSC_FALSE; 1466 break; 1467 default: 1468 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1469 } 1470 break; 1471 default: 1472 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1473 } 1474 if (isAffine && Nq) { 1475 if (v) { 1476 for (i = 0; i < Nq; i++) { 1477 CoordinatesRefToReal(coordDim,dim,v0, J0, &points[dim * i], &v[coordDim * i]); 1478 } 1479 } 1480 if (detJ) { 1481 for (i = 0; i < Nq; i++) { 1482 detJ[i] = detJ0; 1483 } 1484 } 1485 if (J) { 1486 PetscInt k; 1487 1488 for (i = 0, k = 0; i < Nq; i++) { 1489 PetscInt j; 1490 1491 for (j = 0; j < coordDim * coordDim; j++, k++) { 1492 J[k] = J0[j]; 1493 } 1494 } 1495 } 1496 if (invJ) { 1497 PetscInt k; 1498 switch (coordDim) { 1499 case 0: 1500 break; 1501 case 1: 1502 invJ[0] = 1./J0[0]; 1503 break; 1504 case 2: 1505 DMPlex_Invert2D_Internal(invJ, J0, detJ0); 1506 break; 1507 case 3: 1508 DMPlex_Invert3D_Internal(invJ, J0, detJ0); 1509 break; 1510 } 1511 for (i = 1, k = coordDim * coordDim; i < Nq; i++) { 1512 PetscInt j; 1513 1514 for (j = 0; j < coordDim * coordDim; j++, k++) { 1515 invJ[k] = invJ[j]; 1516 } 1517 } 1518 } 1519 } 1520 PetscFunctionReturn(0); 1521 } 1522 1523 /*@C 1524 DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 1525 1526 Collective on DM 1527 1528 Input Arguments: 1529 + dm - the DM 1530 - cell - the cell 1531 1532 Output Arguments: 1533 + v0 - the translation part of this affine transform 1534 . J - the Jacobian of the transform from the reference element 1535 . invJ - the inverse of the Jacobian 1536 - detJ - the Jacobian determinant 1537 1538 Level: advanced 1539 1540 Fortran Notes: 1541 Since it returns arrays, this routine is only available in Fortran 90, and you must 1542 include petsc.h90 in your code. 1543 1544 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 1545 @*/ 1546 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1547 { 1548 PetscErrorCode ierr; 1549 1550 PetscFunctionBegin; 1551 ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr); 1552 PetscFunctionReturn(0); 1553 } 1554 1555 static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1556 { 1557 PetscQuadrature feQuad; 1558 PetscSection coordSection; 1559 Vec coordinates; 1560 PetscScalar *coords = NULL; 1561 const PetscReal *quadPoints; 1562 PetscReal *basisDer, *basis, detJt; 1563 PetscInt dim, cdim, pdim, qdim, Nq, numCoords, q; 1564 PetscErrorCode ierr; 1565 1566 PetscFunctionBegin; 1567 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1568 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1569 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1570 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1571 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 1572 if (!quad) { /* use the first point of the first functional of the dual space */ 1573 PetscDualSpace dsp; 1574 1575 ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr); 1576 ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr); 1577 ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1578 Nq = 1; 1579 } else { 1580 ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1581 } 1582 ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 1583 ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr); 1584 if (feQuad == quad) { 1585 ierr = PetscFEGetDefaultTabulation(fe, &basis, &basisDer, NULL);CHKERRQ(ierr); 1586 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); 1587 } else { 1588 ierr = PetscFEGetTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr); 1589 } 1590 if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 1591 if (v) { 1592 ierr = PetscMemzero(v, Nq*cdim*sizeof(PetscReal));CHKERRQ(ierr); 1593 for (q = 0; q < Nq; ++q) { 1594 PetscInt i, k; 1595 1596 for (k = 0; k < pdim; ++k) 1597 for (i = 0; i < cdim; ++i) 1598 v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]); 1599 ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr); 1600 } 1601 } 1602 if (J) { 1603 ierr = PetscMemzero(J, Nq*cdim*cdim*sizeof(PetscReal));CHKERRQ(ierr); 1604 for (q = 0; q < Nq; ++q) { 1605 PetscInt i, j, k, c, r; 1606 1607 /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 1608 for (k = 0; k < pdim; ++k) 1609 for (j = 0; j < dim; ++j) 1610 for (i = 0; i < cdim; ++i) 1611 J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 1612 ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr); 1613 if (cdim > dim) { 1614 for (c = dim; c < cdim; ++c) 1615 for (r = 0; r < cdim; ++r) 1616 J[r*cdim+c] = r == c ? 1.0 : 0.0; 1617 } 1618 if (!detJ && !invJ) continue; 1619 detJt = 0.; 1620 switch (cdim) { 1621 case 3: 1622 DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]); 1623 if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);} 1624 break; 1625 case 2: 1626 DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]); 1627 if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);} 1628 break; 1629 case 1: 1630 detJt = J[q*cdim*dim]; 1631 if (invJ) invJ[q*cdim*dim] = 1.0/detJt; 1632 } 1633 if (detJ) detJ[q] = detJt; 1634 } 1635 } 1636 else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ"); 1637 if (feQuad != quad) { 1638 ierr = PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr); 1639 } 1640 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1641 PetscFunctionReturn(0); 1642 } 1643 1644 /*@C 1645 DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 1646 1647 Collective on DM 1648 1649 Input Arguments: 1650 + dm - the DM 1651 . cell - the cell 1652 - quad - the quadrature containing the points in the reference element where the geometry will be evaluated. If quad == NULL, geometry will be 1653 evaluated at the first vertex of the reference element 1654 1655 Output Arguments: 1656 + v - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element 1657 . J - the Jacobian of the transform from the reference element at each quadrature point 1658 . invJ - the inverse of the Jacobian at each quadrature point 1659 - detJ - the Jacobian determinant at each quadrature point 1660 1661 Level: advanced 1662 1663 Fortran Notes: 1664 Since it returns arrays, this routine is only available in Fortran 90, and you must 1665 include petsc.h90 in your code. 1666 1667 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1668 @*/ 1669 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1670 { 1671 PetscFE fe = NULL; 1672 PetscErrorCode ierr; 1673 1674 PetscFunctionBegin; 1675 if (dm->coordinateDM) { 1676 PetscClassId id; 1677 PetscInt numFields; 1678 PetscDS prob = dm->coordinateDM->prob; 1679 PetscObject disc; 1680 1681 ierr = PetscDSGetNumFields(prob, &numFields);CHKERRQ(ierr); 1682 if (numFields) { 1683 ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr); 1684 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 1685 if (id == PETSCFE_CLASSID) { 1686 fe = (PetscFE) disc; 1687 } 1688 } 1689 } 1690 if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);} 1691 else {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);} 1692 PetscFunctionReturn(0); 1693 } 1694 1695 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1696 { 1697 PetscSection coordSection; 1698 Vec coordinates; 1699 PetscScalar *coords = NULL; 1700 PetscScalar tmp[2]; 1701 PetscInt coordSize; 1702 PetscErrorCode ierr; 1703 1704 PetscFunctionBegin; 1705 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1706 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1707 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1708 if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 1709 ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1710 if (centroid) { 1711 centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 1712 centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1713 } 1714 if (normal) { 1715 PetscReal norm; 1716 1717 normal[0] = -PetscRealPart(coords[1] - tmp[1]); 1718 normal[1] = PetscRealPart(coords[0] - tmp[0]); 1719 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1720 normal[0] /= norm; 1721 normal[1] /= norm; 1722 } 1723 if (vol) { 1724 *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1725 } 1726 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1727 PetscFunctionReturn(0); 1728 } 1729 1730 /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1731 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1732 { 1733 PetscSection coordSection; 1734 Vec coordinates; 1735 PetscScalar *coords = NULL; 1736 PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 1737 PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1738 PetscErrorCode ierr; 1739 1740 PetscFunctionBegin; 1741 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1742 ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 1743 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1744 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1745 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1746 if (dim > 2 && centroid) { 1747 v0[0] = PetscRealPart(coords[0]); 1748 v0[1] = PetscRealPart(coords[1]); 1749 v0[2] = PetscRealPart(coords[2]); 1750 } 1751 if (normal) { 1752 if (dim > 2) { 1753 const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 1754 const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 1755 const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 1756 PetscReal norm; 1757 1758 normal[0] = y0*z1 - z0*y1; 1759 normal[1] = z0*x1 - x0*z1; 1760 normal[2] = x0*y1 - y0*x1; 1761 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 1762 normal[0] /= norm; 1763 normal[1] /= norm; 1764 normal[2] /= norm; 1765 } else { 1766 for (d = 0; d < dim; ++d) normal[d] = 0.0; 1767 } 1768 } 1769 if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);} 1770 for (p = 0; p < numCorners; ++p) { 1771 /* Need to do this copy to get types right */ 1772 for (d = 0; d < tdim; ++d) { 1773 ctmp[d] = PetscRealPart(coords[p*tdim+d]); 1774 ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 1775 } 1776 Volume_Triangle_Origin_Internal(&vtmp, ctmp); 1777 vsum += vtmp; 1778 for (d = 0; d < tdim; ++d) { 1779 csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 1780 } 1781 } 1782 for (d = 0; d < tdim; ++d) { 1783 csum[d] /= (tdim+1)*vsum; 1784 } 1785 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1786 if (vol) *vol = PetscAbsReal(vsum); 1787 if (centroid) { 1788 if (dim > 2) { 1789 for (d = 0; d < dim; ++d) { 1790 centroid[d] = v0[d]; 1791 for (e = 0; e < dim; ++e) { 1792 centroid[d] += R[d*dim+e]*csum[e]; 1793 } 1794 } 1795 } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 1796 } 1797 PetscFunctionReturn(0); 1798 } 1799 1800 /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1801 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1802 { 1803 PetscSection coordSection; 1804 Vec coordinates; 1805 PetscScalar *coords = NULL; 1806 PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1807 const PetscInt *faces, *facesO; 1808 PetscInt numFaces, f, coordSize, numCorners, p, d; 1809 PetscErrorCode ierr; 1810 1811 PetscFunctionBegin; 1812 if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 1813 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1814 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1815 1816 if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 1817 ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 1818 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1819 ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 1820 for (f = 0; f < numFaces; ++f) { 1821 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1822 numCorners = coordSize/dim; 1823 switch (numCorners) { 1824 case 3: 1825 for (d = 0; d < dim; ++d) { 1826 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1827 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1828 coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 1829 } 1830 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1831 if (facesO[f] < 0) vtmp = -vtmp; 1832 vsum += vtmp; 1833 if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 1834 for (d = 0; d < dim; ++d) { 1835 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1836 } 1837 } 1838 break; 1839 case 4: 1840 /* DO FOR PYRAMID */ 1841 /* First tet */ 1842 for (d = 0; d < dim; ++d) { 1843 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1844 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1845 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1846 } 1847 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1848 if (facesO[f] < 0) vtmp = -vtmp; 1849 vsum += vtmp; 1850 if (centroid) { 1851 for (d = 0; d < dim; ++d) { 1852 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1853 } 1854 } 1855 /* Second tet */ 1856 for (d = 0; d < dim; ++d) { 1857 coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 1858 coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 1859 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1860 } 1861 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1862 if (facesO[f] < 0) vtmp = -vtmp; 1863 vsum += vtmp; 1864 if (centroid) { 1865 for (d = 0; d < dim; ++d) { 1866 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1867 } 1868 } 1869 break; 1870 default: 1871 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 1872 } 1873 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1874 } 1875 if (vol) *vol = PetscAbsReal(vsum); 1876 if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1877 if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 1878 PetscFunctionReturn(0); 1879 } 1880 1881 /*@C 1882 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1883 1884 Collective on DM 1885 1886 Input Arguments: 1887 + dm - the DM 1888 - cell - the cell 1889 1890 Output Arguments: 1891 + volume - the cell volume 1892 . centroid - the cell centroid 1893 - normal - the cell normal, if appropriate 1894 1895 Level: advanced 1896 1897 Fortran Notes: 1898 Since it returns arrays, this routine is only available in Fortran 90, and you must 1899 include petsc.h90 in your code. 1900 1901 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1902 @*/ 1903 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1904 { 1905 PetscInt depth, dim; 1906 PetscErrorCode ierr; 1907 1908 PetscFunctionBegin; 1909 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1910 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1911 if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1912 /* We need to keep a pointer to the depth label */ 1913 ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1914 /* Cone size is now the number of faces */ 1915 switch (depth) { 1916 case 1: 1917 ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1918 break; 1919 case 2: 1920 ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1921 break; 1922 case 3: 1923 ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1924 break; 1925 default: 1926 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1927 } 1928 PetscFunctionReturn(0); 1929 } 1930 1931 /*@ 1932 DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh 1933 1934 Collective on dm 1935 1936 Input Parameter: 1937 . dm - The DMPlex 1938 1939 Output Parameter: 1940 . cellgeom - A vector with the cell geometry data for each cell 1941 1942 Level: beginner 1943 1944 .keywords: DMPlexComputeCellGeometryFEM() 1945 @*/ 1946 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1947 { 1948 DM dmCell; 1949 Vec coordinates; 1950 PetscSection coordSection, sectionCell; 1951 PetscScalar *cgeom; 1952 PetscInt cStart, cEnd, cMax, c; 1953 PetscErrorCode ierr; 1954 1955 PetscFunctionBegin; 1956 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1957 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1958 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1959 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1960 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1961 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1962 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1963 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1964 cEnd = cMax < 0 ? cEnd : cMax; 1965 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1966 /* TODO This needs to be multiplied by Nq for non-affine */ 1967 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1968 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1969 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1970 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1971 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1972 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1973 for (c = cStart; c < cEnd; ++c) { 1974 PetscFECellGeom *cg; 1975 1976 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1977 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1978 ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1979 if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1980 } 1981 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1982 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1983 PetscFunctionReturn(0); 1984 } 1985 1986 /*@ 1987 DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1988 1989 Input Parameter: 1990 . dm - The DM 1991 1992 Output Parameters: 1993 + cellgeom - A Vec of PetscFVCellGeom data 1994 . facegeom - A Vec of PetscFVFaceGeom data 1995 1996 Level: developer 1997 1998 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1999 @*/ 2000 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 2001 { 2002 DM dmFace, dmCell; 2003 DMLabel ghostLabel; 2004 PetscSection sectionFace, sectionCell; 2005 PetscSection coordSection; 2006 Vec coordinates; 2007 PetscScalar *fgeom, *cgeom; 2008 PetscReal minradius, gminradius; 2009 PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 2010 PetscErrorCode ierr; 2011 2012 PetscFunctionBegin; 2013 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2014 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2015 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2016 /* Make cell centroids and volumes */ 2017 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 2018 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 2019 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 2020 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 2021 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2022 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2023 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 2024 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 2025 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 2026 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 2027 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 2028 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 2029 if (cEndInterior < 0) { 2030 cEndInterior = cEnd; 2031 } 2032 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 2033 for (c = cStart; c < cEndInterior; ++c) { 2034 PetscFVCellGeom *cg; 2035 2036 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2037 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 2038 ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 2039 } 2040 /* Compute face normals and minimum cell radius */ 2041 ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 2042 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 2043 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2044 ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 2045 for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 2046 ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 2047 ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 2048 ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 2049 ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 2050 ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 2051 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2052 minradius = PETSC_MAX_REAL; 2053 for (f = fStart; f < fEnd; ++f) { 2054 PetscFVFaceGeom *fg; 2055 PetscReal area; 2056 PetscInt ghost = -1, d, numChildren; 2057 2058 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 2059 ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 2060 if (ghost >= 0 || numChildren) continue; 2061 ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 2062 ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 2063 for (d = 0; d < dim; ++d) fg->normal[d] *= area; 2064 /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 2065 { 2066 PetscFVCellGeom *cL, *cR; 2067 PetscInt ncells; 2068 const PetscInt *cells; 2069 PetscReal *lcentroid, *rcentroid; 2070 PetscReal l[3], r[3], v[3]; 2071 2072 ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 2073 ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr); 2074 ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 2075 lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 2076 if (ncells > 1) { 2077 ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 2078 rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 2079 } 2080 else { 2081 rcentroid = fg->centroid; 2082 } 2083 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 2084 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 2085 DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 2086 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 2087 for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 2088 } 2089 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 2090 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]); 2091 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]); 2092 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 2093 } 2094 if (cells[0] < cEndInterior) { 2095 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 2096 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 2097 } 2098 if (ncells > 1 && cells[1] < cEndInterior) { 2099 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 2100 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 2101 } 2102 } 2103 } 2104 ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 2105 ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 2106 /* Compute centroids of ghost cells */ 2107 for (c = cEndInterior; c < cEnd; ++c) { 2108 PetscFVFaceGeom *fg; 2109 const PetscInt *cone, *support; 2110 PetscInt coneSize, supportSize, s; 2111 2112 ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 2113 if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 2114 ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 2115 ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 2116 if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize); 2117 ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 2118 ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 2119 for (s = 0; s < 2; ++s) { 2120 /* Reflect ghost centroid across plane of face */ 2121 if (support[s] == c) { 2122 PetscFVCellGeom *ci; 2123 PetscFVCellGeom *cg; 2124 PetscReal c2f[3], a; 2125 2126 ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 2127 DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 2128 a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 2129 ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 2130 DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 2131 cg->volume = ci->volume; 2132 } 2133 } 2134 } 2135 ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 2136 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 2137 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 2138 ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 2139 PetscFunctionReturn(0); 2140 } 2141 2142 /*@C 2143 DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 2144 2145 Not collective 2146 2147 Input Argument: 2148 . dm - the DM 2149 2150 Output Argument: 2151 . minradius - the minium cell radius 2152 2153 Level: developer 2154 2155 .seealso: DMGetCoordinates() 2156 @*/ 2157 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 2158 { 2159 PetscFunctionBegin; 2160 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2161 PetscValidPointer(minradius,2); 2162 *minradius = ((DM_Plex*) dm->data)->minradius; 2163 PetscFunctionReturn(0); 2164 } 2165 2166 /*@C 2167 DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 2168 2169 Logically collective 2170 2171 Input Arguments: 2172 + dm - the DM 2173 - minradius - the minium cell radius 2174 2175 Level: developer 2176 2177 .seealso: DMSetCoordinates() 2178 @*/ 2179 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 2180 { 2181 PetscFunctionBegin; 2182 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2183 ((DM_Plex*) dm->data)->minradius = minradius; 2184 PetscFunctionReturn(0); 2185 } 2186 2187 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2188 { 2189 DMLabel ghostLabel; 2190 PetscScalar *dx, *grad, **gref; 2191 PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 2192 PetscErrorCode ierr; 2193 2194 PetscFunctionBegin; 2195 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2196 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2197 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2198 ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 2199 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 2200 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2201 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 2202 for (c = cStart; c < cEndInterior; c++) { 2203 const PetscInt *faces; 2204 PetscInt numFaces, usedFaces, f, d; 2205 PetscFVCellGeom *cg; 2206 PetscBool boundary; 2207 PetscInt ghost; 2208 2209 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2210 ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 2211 ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 2212 if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 2213 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 2214 PetscFVCellGeom *cg1; 2215 PetscFVFaceGeom *fg; 2216 const PetscInt *fcells; 2217 PetscInt ncell, side; 2218 2219 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2220 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2221 if ((ghost >= 0) || boundary) continue; 2222 ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 2223 side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 2224 ncell = fcells[!side]; /* the neighbor */ 2225 ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 2226 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 2227 for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 2228 gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 2229 } 2230 if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 2231 ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 2232 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 2233 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2234 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2235 if ((ghost >= 0) || boundary) continue; 2236 for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 2237 ++usedFaces; 2238 } 2239 } 2240 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 2241 PetscFunctionReturn(0); 2242 } 2243 2244 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2245 { 2246 DMLabel ghostLabel; 2247 PetscScalar *dx, *grad, **gref; 2248 PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 2249 PetscSection neighSec; 2250 PetscInt (*neighbors)[2]; 2251 PetscInt *counter; 2252 PetscErrorCode ierr; 2253 2254 PetscFunctionBegin; 2255 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2256 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2257 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2258 if (cEndInterior < 0) { 2259 cEndInterior = cEnd; 2260 } 2261 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr); 2262 ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr); 2263 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2264 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2265 for (f = fStart; f < fEnd; f++) { 2266 const PetscInt *fcells; 2267 PetscBool boundary; 2268 PetscInt ghost = -1; 2269 PetscInt numChildren, numCells, c; 2270 2271 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 2272 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 2273 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 2274 if ((ghost >= 0) || boundary || numChildren) continue; 2275 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 2276 if (numCells == 2) { 2277 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 2278 for (c = 0; c < 2; c++) { 2279 PetscInt cell = fcells[c]; 2280 2281 if (cell >= cStart && cell < cEndInterior) { 2282 ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr); 2283 } 2284 } 2285 } 2286 } 2287 ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr); 2288 ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr); 2289 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 2290 nStart = 0; 2291 ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr); 2292 ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr); 2293 ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr); 2294 for (f = fStart; f < fEnd; f++) { 2295 const PetscInt *fcells; 2296 PetscBool boundary; 2297 PetscInt ghost = -1; 2298 PetscInt numChildren, numCells, c; 2299 2300 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 2301 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 2302 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 2303 if ((ghost >= 0) || boundary || numChildren) continue; 2304 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 2305 if (numCells == 2) { 2306 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 2307 for (c = 0; c < 2; c++) { 2308 PetscInt cell = fcells[c], off; 2309 2310 if (cell >= cStart && cell < cEndInterior) { 2311 ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr); 2312 off += counter[cell - cStart]++; 2313 neighbors[off][0] = f; 2314 neighbors[off][1] = fcells[1 - c]; 2315 } 2316 } 2317 } 2318 } 2319 ierr = PetscFree(counter);CHKERRQ(ierr); 2320 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 2321 for (c = cStart; c < cEndInterior; c++) { 2322 PetscInt numFaces, f, d, off, ghost = -1; 2323 PetscFVCellGeom *cg; 2324 2325 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2326 ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr); 2327 ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr); 2328 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);} 2329 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); 2330 for (f = 0; f < numFaces; ++f) { 2331 PetscFVCellGeom *cg1; 2332 PetscFVFaceGeom *fg; 2333 const PetscInt *fcells; 2334 PetscInt ncell, side, nface; 2335 2336 nface = neighbors[off + f][0]; 2337 ncell = neighbors[off + f][1]; 2338 ierr = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr); 2339 side = (c != fcells[0]); 2340 ierr = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr); 2341 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 2342 for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d]; 2343 gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 2344 } 2345 ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr); 2346 for (f = 0; f < numFaces; ++f) { 2347 for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d]; 2348 } 2349 } 2350 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 2351 ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr); 2352 ierr = PetscFree(neighbors);CHKERRQ(ierr); 2353 PetscFunctionReturn(0); 2354 } 2355 2356 /*@ 2357 DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 2358 2359 Collective on DM 2360 2361 Input Arguments: 2362 + dm - The DM 2363 . fvm - The PetscFV 2364 . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM() 2365 - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM() 2366 2367 Output Parameters: 2368 + faceGeometry - The geometric factors for gradient calculation are inserted 2369 - dmGrad - The DM describing the layout of gradient data 2370 2371 Level: developer 2372 2373 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 2374 @*/ 2375 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 2376 { 2377 DM dmFace, dmCell; 2378 PetscScalar *fgeom, *cgeom; 2379 PetscSection sectionGrad, parentSection; 2380 PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 2381 PetscErrorCode ierr; 2382 2383 PetscFunctionBegin; 2384 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2385 ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 2386 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2387 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2388 /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 2389 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 2390 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 2391 ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2392 ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2393 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2394 if (!parentSection) { 2395 ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2396 } else { 2397 ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2398 } 2399 ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2400 ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2401 /* Create storage for gradients */ 2402 ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 2403 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 2404 ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 2405 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 2406 ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 2407 ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 2408 ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 2409 PetscFunctionReturn(0); 2410 } 2411 2412 /*@ 2413 DMPlexGetDataFVM - Retrieve precomputed cell geometry 2414 2415 Collective on DM 2416 2417 Input Arguments: 2418 + dm - The DM 2419 - fvm - The PetscFV 2420 2421 Output Parameters: 2422 + cellGeometry - The cell geometry 2423 . faceGeometry - The face geometry 2424 - dmGrad - The gradient matrices 2425 2426 Level: developer 2427 2428 .seealso: DMPlexComputeGeometryFVM() 2429 @*/ 2430 PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM) 2431 { 2432 PetscObject cellgeomobj, facegeomobj; 2433 PetscErrorCode ierr; 2434 2435 PetscFunctionBegin; 2436 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2437 if (!cellgeomobj) { 2438 Vec cellgeomInt, facegeomInt; 2439 2440 ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr); 2441 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr); 2442 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr); 2443 ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr); 2444 ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr); 2445 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2446 } 2447 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr); 2448 if (cellgeom) *cellgeom = (Vec) cellgeomobj; 2449 if (facegeom) *facegeom = (Vec) facegeomobj; 2450 if (gradDM) { 2451 PetscObject gradobj; 2452 PetscBool computeGradients; 2453 2454 ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr); 2455 if (!computeGradients) { 2456 *gradDM = NULL; 2457 PetscFunctionReturn(0); 2458 } 2459 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2460 if (!gradobj) { 2461 DM dmGradInt; 2462 2463 ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr); 2464 ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr); 2465 ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr); 2466 ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2467 } 2468 *gradDM = (DM) gradobj; 2469 } 2470 PetscFunctionReturn(0); 2471 } 2472 2473 static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work, PetscReal *resNeg, PetscReal *guess) 2474 { 2475 PetscInt l, m; 2476 2477 PetscFunctionBeginHot; 2478 if (dimC == dimR && dimR <= 3) { 2479 /* invert Jacobian, multiply */ 2480 PetscScalar det, idet; 2481 2482 switch (dimR) { 2483 case 1: 2484 invJ[0] = 1./ J[0]; 2485 break; 2486 case 2: 2487 det = J[0] * J[3] - J[1] * J[2]; 2488 idet = 1./det; 2489 invJ[0] = J[3] * idet; 2490 invJ[1] = -J[1] * idet; 2491 invJ[2] = -J[2] * idet; 2492 invJ[3] = J[0] * idet; 2493 break; 2494 case 3: 2495 { 2496 invJ[0] = J[4] * J[8] - J[5] * J[7]; 2497 invJ[1] = J[2] * J[7] - J[1] * J[8]; 2498 invJ[2] = J[1] * J[5] - J[2] * J[4]; 2499 det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6]; 2500 idet = 1./det; 2501 invJ[0] *= idet; 2502 invJ[1] *= idet; 2503 invJ[2] *= idet; 2504 invJ[3] = idet * (J[5] * J[6] - J[3] * J[8]); 2505 invJ[4] = idet * (J[0] * J[8] - J[2] * J[6]); 2506 invJ[5] = idet * (J[2] * J[3] - J[0] * J[5]); 2507 invJ[6] = idet * (J[3] * J[7] - J[4] * J[6]); 2508 invJ[7] = idet * (J[1] * J[6] - J[0] * J[7]); 2509 invJ[8] = idet * (J[0] * J[4] - J[1] * J[3]); 2510 } 2511 break; 2512 } 2513 for (l = 0; l < dimR; l++) { 2514 for (m = 0; m < dimC; m++) { 2515 guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m]; 2516 } 2517 } 2518 } else { 2519 #if defined(PETSC_USE_COMPLEX) 2520 char transpose = 'C'; 2521 #else 2522 char transpose = 'T'; 2523 #endif 2524 PetscBLASInt m = dimR; 2525 PetscBLASInt n = dimC; 2526 PetscBLASInt one = 1; 2527 PetscBLASInt worksize = dimR * dimC, info; 2528 2529 for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];} 2530 2531 PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info)); 2532 if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS"); 2533 2534 for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);} 2535 } 2536 PetscFunctionReturn(0); 2537 } 2538 2539 static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 2540 { 2541 PetscInt coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR); 2542 PetscScalar *coordsScalar = NULL; 2543 PetscReal *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg; 2544 PetscScalar *J, *invJ, *work; 2545 PetscErrorCode ierr; 2546 2547 PetscFunctionBegin; 2548 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2549 ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 2550 if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr); 2551 ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr); 2552 ierr = DMGetWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr); 2553 cellCoords = &cellData[0]; 2554 cellCoeffs = &cellData[coordSize]; 2555 extJ = &cellData[2 * coordSize]; 2556 resNeg = &cellData[2 * coordSize + dimR]; 2557 invJ = &J[dimR * dimC]; 2558 work = &J[2 * dimR * dimC]; 2559 if (dimR == 2) { 2560 const PetscInt zToPlex[4] = {0, 1, 3, 2}; 2561 2562 for (i = 0; i < 4; i++) { 2563 PetscInt plexI = zToPlex[i]; 2564 2565 for (j = 0; j < dimC; j++) { 2566 cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 2567 } 2568 } 2569 } else if (dimR == 3) { 2570 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 2571 2572 for (i = 0; i < 8; i++) { 2573 PetscInt plexI = zToPlex[i]; 2574 2575 for (j = 0; j < dimC; j++) { 2576 cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 2577 } 2578 } 2579 } else { 2580 for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);} 2581 } 2582 /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 2583 for (i = 0; i < dimR; i++) { 2584 PetscReal *swap; 2585 2586 for (j = 0; j < (numV / 2); j++) { 2587 for (k = 0; k < dimC; k++) { 2588 cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 2589 cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 2590 } 2591 } 2592 2593 if (i < dimR - 1) { 2594 swap = cellCoeffs; 2595 cellCoeffs = cellCoords; 2596 cellCoords = swap; 2597 } 2598 } 2599 ierr = PetscMemzero(refCoords,numPoints * dimR * sizeof (PetscReal));CHKERRQ(ierr); 2600 for (j = 0; j < numPoints; j++) { 2601 for (i = 0; i < maxIts; i++) { 2602 PetscReal *guess = &refCoords[dimR * j]; 2603 2604 /* compute -residual and Jacobian */ 2605 for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];} 2606 for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;} 2607 for (k = 0; k < numV; k++) { 2608 PetscReal extCoord = 1.; 2609 for (l = 0; l < dimR; l++) { 2610 PetscReal coord = guess[l]; 2611 PetscInt dep = (k & (1 << l)) >> l; 2612 2613 extCoord *= dep * coord + !dep; 2614 extJ[l] = dep; 2615 2616 for (m = 0; m < dimR; m++) { 2617 PetscReal coord = guess[m]; 2618 PetscInt dep = ((k & (1 << m)) >> m) && (m != l); 2619 PetscReal mult = dep * coord + !dep; 2620 2621 extJ[l] *= mult; 2622 } 2623 } 2624 for (l = 0; l < dimC; l++) { 2625 PetscReal coeff = cellCoeffs[dimC * k + l]; 2626 2627 resNeg[l] -= coeff * extCoord; 2628 for (m = 0; m < dimR; m++) { 2629 J[dimR * l + m] += coeff * extJ[m]; 2630 } 2631 } 2632 } 2633 #if 0 && defined(PETSC_USE_DEBUG) 2634 { 2635 PetscReal maxAbs = 0.; 2636 2637 for (l = 0; l < dimC; l++) { 2638 maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l])); 2639 } 2640 ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr); 2641 } 2642 #endif 2643 2644 ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr); 2645 } 2646 } 2647 ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr); 2648 ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr); 2649 ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 2650 PetscFunctionReturn(0); 2651 } 2652 2653 static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 2654 { 2655 PetscInt coordSize, i, j, k, l, numV = (1 << dimR); 2656 PetscScalar *coordsScalar = NULL; 2657 PetscReal *cellData, *cellCoords, *cellCoeffs; 2658 PetscErrorCode ierr; 2659 2660 PetscFunctionBegin; 2661 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2662 ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 2663 if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr); 2664 ierr = DMGetWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr); 2665 cellCoords = &cellData[0]; 2666 cellCoeffs = &cellData[coordSize]; 2667 if (dimR == 2) { 2668 const PetscInt zToPlex[4] = {0, 1, 3, 2}; 2669 2670 for (i = 0; i < 4; i++) { 2671 PetscInt plexI = zToPlex[i]; 2672 2673 for (j = 0; j < dimC; j++) { 2674 cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 2675 } 2676 } 2677 } else if (dimR == 3) { 2678 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 2679 2680 for (i = 0; i < 8; i++) { 2681 PetscInt plexI = zToPlex[i]; 2682 2683 for (j = 0; j < dimC; j++) { 2684 cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 2685 } 2686 } 2687 } else { 2688 for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);} 2689 } 2690 /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 2691 for (i = 0; i < dimR; i++) { 2692 PetscReal *swap; 2693 2694 for (j = 0; j < (numV / 2); j++) { 2695 for (k = 0; k < dimC; k++) { 2696 cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 2697 cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 2698 } 2699 } 2700 2701 if (i < dimR - 1) { 2702 swap = cellCoeffs; 2703 cellCoeffs = cellCoords; 2704 cellCoords = swap; 2705 } 2706 } 2707 ierr = PetscMemzero(realCoords,numPoints * dimC * sizeof (PetscReal));CHKERRQ(ierr); 2708 for (j = 0; j < numPoints; j++) { 2709 const PetscReal *guess = &refCoords[dimR * j]; 2710 PetscReal *mapped = &realCoords[dimC * j]; 2711 2712 for (k = 0; k < numV; k++) { 2713 PetscReal extCoord = 1.; 2714 for (l = 0; l < dimR; l++) { 2715 PetscReal coord = guess[l]; 2716 PetscInt dep = (k & (1 << l)) >> l; 2717 2718 extCoord *= dep * coord + !dep; 2719 } 2720 for (l = 0; l < dimC; l++) { 2721 PetscReal coeff = cellCoeffs[dimC * k + l]; 2722 2723 mapped[l] += coeff * extCoord; 2724 } 2725 } 2726 } 2727 ierr = DMRestoreWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr); 2728 ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 2729 PetscFunctionReturn(0); 2730 } 2731 2732 static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 2733 { 2734 PetscInt numComp, numDof, i, j, k, l, m, maxIter = 7, coordSize; 2735 PetscScalar *nodes = NULL; 2736 PetscReal *invV, *modes; 2737 PetscReal *B, *D, *resNeg; 2738 PetscScalar *J, *invJ, *work; 2739 PetscErrorCode ierr; 2740 2741 PetscFunctionBegin; 2742 ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr); 2743 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 2744 if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC); 2745 ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 2746 /* convert nodes to values in the stable evaluation basis */ 2747 ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 2748 invV = fe->invV; 2749 for (i = 0; i < numDof; i++) { 2750 for (j = 0; j < dimC; j++) { 2751 modes[i * dimC + j] = 0.; 2752 for (k = 0; k < numDof; k++) { 2753 modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]); 2754 } 2755 } 2756 } 2757 ierr = DMGetWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr); 2758 D = &B[numDof]; 2759 resNeg = &D[numDof * dimR]; 2760 ierr = DMGetWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr); 2761 invJ = &J[dimC * dimR]; 2762 work = &invJ[dimC * dimR]; 2763 for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;} 2764 for (j = 0; j < numPoints; j++) { 2765 for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */ 2766 PetscReal *guess = &refCoords[j * dimR]; 2767 ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr); 2768 for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[j * dimC + k];} 2769 for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;} 2770 for (k = 0; k < numDof; k++) { 2771 for (l = 0; l < dimC; l++) { 2772 resNeg[l] -= modes[k * dimC + l] * B[k]; 2773 for (m = 0; m < dimR; m++) { 2774 J[l * dimR + m] += modes[k * dimC + l] * D[k * dimR + m]; 2775 } 2776 } 2777 } 2778 #if 0 && defined(PETSC_USE_DEBUG) 2779 { 2780 PetscReal maxAbs = 0.; 2781 2782 for (l = 0; l < dimC; l++) { 2783 maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l])); 2784 } 2785 ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr); 2786 } 2787 #endif 2788 ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr); 2789 } 2790 } 2791 ierr = DMRestoreWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr); 2792 ierr = DMRestoreWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr); 2793 ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 2794 ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 2795 PetscFunctionReturn(0); 2796 } 2797 2798 static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 2799 { 2800 PetscInt numComp, numDof, i, j, k, l, coordSize; 2801 PetscScalar *nodes = NULL; 2802 PetscReal *invV, *modes; 2803 PetscReal *B; 2804 PetscErrorCode ierr; 2805 2806 PetscFunctionBegin; 2807 ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr); 2808 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 2809 if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC); 2810 ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 2811 /* convert nodes to values in the stable evaluation basis */ 2812 ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 2813 invV = fe->invV; 2814 for (i = 0; i < numDof; i++) { 2815 for (j = 0; j < dimC; j++) { 2816 modes[i * dimC + j] = 0.; 2817 for (k = 0; k < numDof; k++) { 2818 modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]); 2819 } 2820 } 2821 } 2822 ierr = DMGetWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr); 2823 for (i = 0; i < numPoints * dimC; i++) {realCoords[i] = 0.;} 2824 for (j = 0; j < numPoints; j++) { 2825 const PetscReal *guess = &refCoords[j * dimR]; 2826 PetscReal *mapped = &realCoords[j * dimC]; 2827 2828 ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, NULL, NULL);CHKERRQ(ierr); 2829 for (k = 0; k < numDof; k++) { 2830 for (l = 0; l < dimC; l++) { 2831 mapped[l] += modes[k * dimC + l] * B[k]; 2832 } 2833 } 2834 } 2835 ierr = DMRestoreWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr); 2836 ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 2837 ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 2838 PetscFunctionReturn(0); 2839 } 2840 2841 /*@ 2842 DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element 2843 map. This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not 2844 extend uniquely outside the reference cell (e.g, most non-affine maps) 2845 2846 Not collective 2847 2848 Input Parameters: 2849 + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or 2850 implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 2851 as a multilinear map for tensor-product elements 2852 . cell - the cell whose map is used. 2853 . numPoints - the number of points to locate 2854 - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim()) 2855 2856 Output Parameters: 2857 . refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension()) 2858 2859 Level: intermediate 2860 @*/ 2861 PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[]) 2862 { 2863 PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i; 2864 DM coordDM = NULL; 2865 Vec coords; 2866 PetscFE fe = NULL; 2867 PetscErrorCode ierr; 2868 2869 PetscFunctionBegin; 2870 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2871 ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr); 2872 ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr); 2873 if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0); 2874 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 2875 ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 2876 ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 2877 if (coordDM) { 2878 PetscInt coordFields; 2879 2880 ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr); 2881 if (coordFields) { 2882 PetscClassId id; 2883 PetscObject disc; 2884 2885 ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr); 2886 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 2887 if (id == PETSCFE_CLASSID) { 2888 fe = (PetscFE) disc; 2889 } 2890 } 2891 } 2892 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 2893 ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr); 2894 cEnd = cEndInterior > 0 ? cEndInterior : cEnd; 2895 if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr); 2896 if (!fe) { /* implicit discretization: affine or multilinear */ 2897 PetscInt coneSize; 2898 PetscBool isSimplex, isTensor; 2899 2900 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 2901 isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 2902 isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 2903 if (isSimplex) { 2904 PetscReal detJ, *v0, *J, *invJ; 2905 2906 ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 2907 J = &v0[dimC]; 2908 invJ = &J[dimC * dimC]; 2909 ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr); 2910 for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */ 2911 CoordinatesRealToRef(dimC, dimR, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]); 2912 } 2913 ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 2914 } else if (isTensor) { 2915 ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr); 2916 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize); 2917 } else { 2918 ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr); 2919 } 2920 PetscFunctionReturn(0); 2921 } 2922 2923 /*@ 2924 DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map. 2925 2926 Not collective 2927 2928 Input Parameters: 2929 + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or 2930 implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 2931 as a multilinear map for tensor-product elements 2932 . cell - the cell whose map is used. 2933 . numPoints - the number of points to locate 2934 + refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension()) 2935 2936 Output Parameters: 2937 . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim()) 2938 2939 Level: intermediate 2940 @*/ 2941 PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[]) 2942 { 2943 PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i; 2944 DM coordDM = NULL; 2945 Vec coords; 2946 PetscFE fe = NULL; 2947 PetscErrorCode ierr; 2948 2949 PetscFunctionBegin; 2950 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2951 ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr); 2952 ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr); 2953 if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0); 2954 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 2955 ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 2956 ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 2957 if (coordDM) { 2958 PetscInt coordFields; 2959 2960 ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr); 2961 if (coordFields) { 2962 PetscClassId id; 2963 PetscObject disc; 2964 2965 ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr); 2966 ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 2967 if (id == PETSCFE_CLASSID) { 2968 fe = (PetscFE) disc; 2969 } 2970 } 2971 } 2972 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 2973 ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr); 2974 cEnd = cEndInterior > 0 ? cEndInterior : cEnd; 2975 if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr); 2976 if (!fe) { /* implicit discretization: affine or multilinear */ 2977 PetscInt coneSize; 2978 PetscBool isSimplex, isTensor; 2979 2980 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 2981 isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 2982 isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 2983 if (isSimplex) { 2984 PetscReal detJ, *v0, *J; 2985 2986 ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 2987 J = &v0[dimC]; 2988 ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr); 2989 for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */ 2990 CoordinatesRefToReal(dimC, dimR, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]); 2991 } 2992 ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 2993 } else if (isTensor) { 2994 ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr); 2995 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize); 2996 } else { 2997 ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr); 2998 } 2999 PetscFunctionReturn(0); 3000 } 3001