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