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