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