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 const char *const DMPlexCoordMaps[] = {"none", "shear", "flare", "annulus", "shell", "unknown", "DMPlexCoordMap", "DM_COORD_MAP_", NULL}; 7 8 /*@ 9 DMPlexFindVertices - Try to find DAG points based on their coordinates. 10 11 Not Collective (provided `DMGetCoordinatesLocalSetUp()` has been already called) 12 13 Input Parameters: 14 + dm - The `DMPLEX` object 15 . coordinates - The `Vec` of coordinates of the sought points 16 - eps - The tolerance or `PETSC_DEFAULT` 17 18 Output Parameter: 19 . points - The `IS` of found DAG points or -1 20 21 Level: intermediate 22 23 Notes: 24 The length of `Vec` coordinates must be npoints * dim where dim is the spatial dimension returned by `DMGetCoordinateDim()` and npoints is the number of sought points. 25 26 The output `IS` is living on `PETSC_COMM_SELF` and its length is npoints. 27 Each rank does the search independently. 28 If this rank's local `DMPLEX` portion contains the DAG point corresponding to the i-th tuple of coordinates, the i-th entry of the output `IS` is set to that DAG point, otherwise to -1. 29 30 The output `IS` must be destroyed by user. 31 32 The tolerance is interpreted as the maximum Euclidean (L2) distance of the sought point from the specified coordinates. 33 34 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 if needed. 35 36 .seealso: `DMPLEX`, `DMPlexCreate()`, `DMGetCoordinatesLocal()` 37 @*/ 38 PetscErrorCode DMPlexFindVertices(DM dm, Vec coordinates, PetscReal eps, IS *points) 39 { 40 PetscInt c, cdim, i, j, o, p, vStart, vEnd; 41 PetscInt npoints; 42 const PetscScalar *coord; 43 Vec allCoordsVec; 44 const PetscScalar *allCoords; 45 PetscInt *dagPoints; 46 47 PetscFunctionBegin; 48 if (eps < 0) eps = PETSC_SQRT_MACHINE_EPSILON; 49 PetscCall(DMGetCoordinateDim(dm, &cdim)); 50 { 51 PetscInt n; 52 53 PetscCall(VecGetLocalSize(coordinates, &n)); 54 PetscCheck(n % cdim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Given coordinates Vec has local length %" PetscInt_FMT " not divisible by coordinate dimension %" PetscInt_FMT " of given DM", n, cdim); 55 npoints = n / cdim; 56 } 57 PetscCall(DMGetCoordinatesLocal(dm, &allCoordsVec)); 58 PetscCall(VecGetArrayRead(allCoordsVec, &allCoords)); 59 PetscCall(VecGetArrayRead(coordinates, &coord)); 60 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 61 if (PetscDefined(USE_DEBUG)) { 62 /* check coordinate section is consistent with DM dimension */ 63 PetscSection cs; 64 PetscInt ndof; 65 66 PetscCall(DMGetCoordinateSection(dm, &cs)); 67 for (p = vStart; p < vEnd; p++) { 68 PetscCall(PetscSectionGetDof(cs, p, &ndof)); 69 PetscCheck(ndof == cdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "point %" PetscInt_FMT ": ndof = %" PetscInt_FMT " != %" PetscInt_FMT " = cdim", p, ndof, cdim); 70 } 71 } 72 PetscCall(PetscMalloc1(npoints, &dagPoints)); 73 if (eps == 0.0) { 74 for (i = 0, j = 0; i < npoints; i++, j += cdim) { 75 dagPoints[i] = -1; 76 for (p = vStart, o = 0; p < vEnd; p++, o += cdim) { 77 for (c = 0; c < cdim; c++) { 78 if (coord[j + c] != allCoords[o + c]) break; 79 } 80 if (c == cdim) { 81 dagPoints[i] = p; 82 break; 83 } 84 } 85 } 86 } else { 87 for (i = 0, j = 0; i < npoints; i++, j += cdim) { 88 PetscReal norm; 89 90 dagPoints[i] = -1; 91 for (p = vStart, o = 0; p < vEnd; p++, o += cdim) { 92 norm = 0.0; 93 for (c = 0; c < cdim; c++) norm += PetscRealPart(PetscSqr(coord[j + c] - allCoords[o + c])); 94 norm = PetscSqrtReal(norm); 95 if (norm <= eps) { 96 dagPoints[i] = p; 97 break; 98 } 99 } 100 } 101 } 102 PetscCall(VecRestoreArrayRead(allCoordsVec, &allCoords)); 103 PetscCall(VecRestoreArrayRead(coordinates, &coord)); 104 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, dagPoints, PETSC_OWN_POINTER, points)); 105 PetscFunctionReturn(PETSC_SUCCESS); 106 } 107 108 #if 0 109 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 110 { 111 const PetscReal p0_x = segmentA[0 * 2 + 0]; 112 const PetscReal p0_y = segmentA[0 * 2 + 1]; 113 const PetscReal p1_x = segmentA[1 * 2 + 0]; 114 const PetscReal p1_y = segmentA[1 * 2 + 1]; 115 const PetscReal p2_x = segmentB[0 * 2 + 0]; 116 const PetscReal p2_y = segmentB[0 * 2 + 1]; 117 const PetscReal p3_x = segmentB[1 * 2 + 0]; 118 const PetscReal p3_y = segmentB[1 * 2 + 1]; 119 const PetscReal s1_x = p1_x - p0_x; 120 const PetscReal s1_y = p1_y - p0_y; 121 const PetscReal s2_x = p3_x - p2_x; 122 const PetscReal s2_y = p3_y - p2_y; 123 const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 124 125 PetscFunctionBegin; 126 *hasIntersection = PETSC_FALSE; 127 /* Non-parallel lines */ 128 if (denom != 0.0) { 129 const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 130 const PetscReal t = (s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 131 132 if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 133 *hasIntersection = PETSC_TRUE; 134 if (intersection) { 135 intersection[0] = p0_x + (t * s1_x); 136 intersection[1] = p0_y + (t * s1_y); 137 } 138 } 139 } 140 PetscFunctionReturn(PETSC_SUCCESS); 141 } 142 143 /* The plane is segmentB x segmentC: https://en.wikipedia.org/wiki/Line%E2%80%93plane_intersection */ 144 static PetscErrorCode DMPlexGetLinePlaneIntersection_3D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], const PetscReal segmentC[], PetscReal intersection[], PetscBool *hasIntersection) 145 { 146 const PetscReal p0_x = segmentA[0 * 3 + 0]; 147 const PetscReal p0_y = segmentA[0 * 3 + 1]; 148 const PetscReal p0_z = segmentA[0 * 3 + 2]; 149 const PetscReal p1_x = segmentA[1 * 3 + 0]; 150 const PetscReal p1_y = segmentA[1 * 3 + 1]; 151 const PetscReal p1_z = segmentA[1 * 3 + 2]; 152 const PetscReal q0_x = segmentB[0 * 3 + 0]; 153 const PetscReal q0_y = segmentB[0 * 3 + 1]; 154 const PetscReal q0_z = segmentB[0 * 3 + 2]; 155 const PetscReal q1_x = segmentB[1 * 3 + 0]; 156 const PetscReal q1_y = segmentB[1 * 3 + 1]; 157 const PetscReal q1_z = segmentB[1 * 3 + 2]; 158 const PetscReal r0_x = segmentC[0 * 3 + 0]; 159 const PetscReal r0_y = segmentC[0 * 3 + 1]; 160 const PetscReal r0_z = segmentC[0 * 3 + 2]; 161 const PetscReal r1_x = segmentC[1 * 3 + 0]; 162 const PetscReal r1_y = segmentC[1 * 3 + 1]; 163 const PetscReal r1_z = segmentC[1 * 3 + 2]; 164 const PetscReal s0_x = p1_x - p0_x; 165 const PetscReal s0_y = p1_y - p0_y; 166 const PetscReal s0_z = p1_z - p0_z; 167 const PetscReal s1_x = q1_x - q0_x; 168 const PetscReal s1_y = q1_y - q0_y; 169 const PetscReal s1_z = q1_z - q0_z; 170 const PetscReal s2_x = r1_x - r0_x; 171 const PetscReal s2_y = r1_y - r0_y; 172 const PetscReal s2_z = r1_z - r0_z; 173 const PetscReal s3_x = s1_y * s2_z - s1_z * s2_y; /* s1 x s2 */ 174 const PetscReal s3_y = s1_z * s2_x - s1_x * s2_z; 175 const PetscReal s3_z = s1_x * s2_y - s1_y * s2_x; 176 const PetscReal s4_x = s0_y * s2_z - s0_z * s2_y; /* s0 x s2 */ 177 const PetscReal s4_y = s0_z * s2_x - s0_x * s2_z; 178 const PetscReal s4_z = s0_x * s2_y - s0_y * s2_x; 179 const PetscReal s5_x = s1_y * s0_z - s1_z * s0_y; /* s1 x s0 */ 180 const PetscReal s5_y = s1_z * s0_x - s1_x * s0_z; 181 const PetscReal s5_z = s1_x * s0_y - s1_y * s0_x; 182 const PetscReal denom = -(s0_x * s3_x + s0_y * s3_y + s0_z * s3_z); /* -s0 . (s1 x s2) */ 183 184 PetscFunctionBegin; 185 *hasIntersection = PETSC_FALSE; 186 /* Line not parallel to plane */ 187 if (denom != 0.0) { 188 const PetscReal t = (s3_x * (p0_x - q0_x) + s3_y * (p0_y - q0_y) + s3_z * (p0_z - q0_z)) / denom; 189 const PetscReal u = (s4_x * (p0_x - q0_x) + s4_y * (p0_y - q0_y) + s4_z * (p0_z - q0_z)) / denom; 190 const PetscReal v = (s5_x * (p0_x - q0_x) + s5_y * (p0_y - q0_y) + s5_z * (p0_z - q0_z)) / denom; 191 192 if (t >= 0 && t <= 1 && u >= 0 && u <= 1 && v >= 0 && v <= 1) { 193 *hasIntersection = PETSC_TRUE; 194 if (intersection) { 195 intersection[0] = p0_x + (t * s0_x); 196 intersection[1] = p0_y + (t * s0_y); 197 intersection[2] = p0_z + (t * s0_z); 198 } 199 } 200 } 201 PetscFunctionReturn(PETSC_SUCCESS); 202 } 203 #endif 204 205 static PetscErrorCode DMPlexGetPlaneSimplexIntersection_Coords_Internal(DM dm, PetscInt dim, PetscInt cdim, const PetscScalar coords[], const PetscReal p[], const PetscReal normal[], PetscBool *pos, PetscInt *Nint, PetscReal intPoints[]) 206 { 207 PetscReal d[4]; // distance of vertices to the plane 208 PetscReal dp; // distance from origin to the plane 209 PetscInt n = 0; 210 211 PetscFunctionBegin; 212 if (pos) *pos = PETSC_FALSE; 213 if (Nint) *Nint = 0; 214 if (PetscDefined(USE_DEBUG)) { 215 PetscReal mag = DMPlex_NormD_Internal(cdim, normal); 216 PetscCheck(PetscAbsReal(mag - (PetscReal)1.0) < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Normal vector is not normalized: %g", (double)mag); 217 } 218 219 dp = DMPlex_DotRealD_Internal(cdim, normal, p); 220 for (PetscInt v = 0; v < dim + 1; ++v) { 221 // d[v] is positive, zero, or negative if vertex i is above, on, or below the plane 222 #if defined(PETSC_USE_COMPLEX) 223 PetscReal c[4]; 224 for (PetscInt i = 0; i < cdim; ++i) c[i] = PetscRealPart(coords[v * cdim + i]); 225 d[v] = DMPlex_DotRealD_Internal(cdim, normal, c); 226 #else 227 d[v] = DMPlex_DotRealD_Internal(cdim, normal, &coords[v * cdim]); 228 #endif 229 d[v] -= dp; 230 } 231 232 // If all d are positive or negative, no intersection 233 { 234 PetscInt v; 235 for (v = 0; v < dim + 1; ++v) 236 if (d[v] >= 0.) break; 237 if (v == dim + 1) PetscFunctionReturn(PETSC_SUCCESS); 238 for (v = 0; v < dim + 1; ++v) 239 if (d[v] <= 0.) break; 240 if (v == dim + 1) { 241 if (pos) *pos = PETSC_TRUE; 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 } 245 246 for (PetscInt v = 0; v < dim + 1; ++v) { 247 // Points with zero distance are automatically added to the list. 248 if (PetscAbsReal(d[v]) < PETSC_MACHINE_EPSILON) { 249 for (PetscInt i = 0; i < cdim; ++i) intPoints[n * cdim + i] = PetscRealPart(coords[v * cdim + i]); 250 ++n; 251 } else { 252 // For each point with nonzero distance, seek another point with opposite sign 253 // and higher index, and compute the intersection of the line between those 254 // points and the plane. 255 for (PetscInt w = v + 1; w < dim + 1; ++w) { 256 if (d[v] * d[w] < 0.) { 257 PetscReal inv_dist = 1. / (d[v] - d[w]); 258 for (PetscInt i = 0; i < cdim; ++i) intPoints[n * cdim + i] = (d[v] * PetscRealPart(coords[w * cdim + i]) - d[w] * PetscRealPart(coords[v * cdim + i])) * inv_dist; 259 ++n; 260 } 261 } 262 } 263 } 264 // TODO order output points if there are 4 265 *Nint = n; 266 PetscFunctionReturn(PETSC_SUCCESS); 267 } 268 269 static PetscErrorCode DMPlexGetPlaneSimplexIntersection_Internal(DM dm, PetscInt dim, PetscInt c, const PetscReal p[], const PetscReal normal[], PetscBool *pos, PetscInt *Nint, PetscReal intPoints[]) 270 { 271 const PetscScalar *array; 272 PetscScalar *coords = NULL; 273 PetscInt numCoords; 274 PetscBool isDG; 275 PetscInt cdim; 276 277 PetscFunctionBegin; 278 PetscCall(DMGetCoordinateDim(dm, &cdim)); 279 PetscCheck(cdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM has coordinates in %" PetscInt_FMT "D instead of %" PetscInt_FMT "D", cdim, dim); 280 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 281 PetscCheck(numCoords == dim * (dim + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Tetrahedron should have %" PetscInt_FMT " coordinates, not %" PetscInt_FMT, dim * (dim + 1), numCoords); 282 PetscCall(PetscArrayzero(intPoints, dim * (dim + 1))); 283 284 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, coords, p, normal, pos, Nint, intPoints)); 285 286 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 287 PetscFunctionReturn(PETSC_SUCCESS); 288 } 289 290 static PetscErrorCode DMPlexGetPlaneQuadIntersection_Internal(DM dm, PetscInt dim, PetscInt c, const PetscReal p[], const PetscReal normal[], PetscBool *pos, PetscInt *Nint, PetscReal intPoints[]) 291 { 292 const PetscScalar *array; 293 PetscScalar *coords = NULL; 294 PetscInt numCoords; 295 PetscBool isDG; 296 PetscInt cdim; 297 PetscScalar tcoords[6] = {0., 0., 0., 0., 0., 0.}; 298 const PetscInt vertsA[3] = {0, 1, 3}; 299 const PetscInt vertsB[3] = {1, 2, 3}; 300 PetscInt NintA, NintB; 301 302 PetscFunctionBegin; 303 PetscCall(DMGetCoordinateDim(dm, &cdim)); 304 PetscCheck(cdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM has coordinates in %" PetscInt_FMT "D instead of %" PetscInt_FMT "D", cdim, dim); 305 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 306 PetscCheck(numCoords == dim * 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Quadrilateral should have %" PetscInt_FMT " coordinates, not %" PetscInt_FMT, dim * 4, numCoords); 307 PetscCall(PetscArrayzero(intPoints, dim * 4)); 308 309 for (PetscInt v = 0; v < 3; ++v) 310 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsA[v] * cdim + d]; 311 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintA, intPoints)); 312 for (PetscInt v = 0; v < 3; ++v) 313 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsB[v] * cdim + d]; 314 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintB, &intPoints[NintA * cdim])); 315 *Nint = NintA + NintB; 316 317 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 static PetscErrorCode DMPlexGetPlaneHexIntersection_Internal(DM dm, PetscInt dim, PetscInt c, const PetscReal p[], const PetscReal normal[], PetscBool *pos, PetscInt *Nint, PetscReal intPoints[]) 322 { 323 const PetscScalar *array; 324 PetscScalar *coords = NULL; 325 PetscInt numCoords; 326 PetscBool isDG; 327 PetscInt cdim; 328 PetscScalar tcoords[12] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}; 329 // We split using the (2, 4) main diagonal, so all tets contain those vertices 330 const PetscInt vertsA[4] = {0, 1, 2, 4}; 331 const PetscInt vertsB[4] = {0, 2, 3, 4}; 332 const PetscInt vertsC[4] = {1, 7, 2, 4}; 333 const PetscInt vertsD[4] = {2, 7, 6, 4}; 334 const PetscInt vertsE[4] = {3, 5, 4, 2}; 335 const PetscInt vertsF[4] = {4, 5, 6, 2}; 336 PetscInt NintA, NintB, NintC, NintD, NintE, NintF, Nsum = 0; 337 338 PetscFunctionBegin; 339 PetscCall(DMGetCoordinateDim(dm, &cdim)); 340 PetscCheck(cdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DM has coordinates in %" PetscInt_FMT "D instead of %" PetscInt_FMT "D", cdim, dim); 341 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 342 PetscCheck(numCoords == dim * 8, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Hexahedron should have %" PetscInt_FMT " coordinates, not %" PetscInt_FMT, dim * 8, numCoords); 343 PetscCall(PetscArrayzero(intPoints, dim * 18)); 344 345 for (PetscInt v = 0; v < 4; ++v) 346 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsA[v] * cdim + d]; 347 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintA, &intPoints[Nsum * cdim])); 348 Nsum += NintA; 349 for (PetscInt v = 0; v < 4; ++v) 350 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsB[v] * cdim + d]; 351 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintB, &intPoints[Nsum * cdim])); 352 Nsum += NintB; 353 for (PetscInt v = 0; v < 4; ++v) 354 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsC[v] * cdim + d]; 355 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintC, &intPoints[Nsum * cdim])); 356 Nsum += NintC; 357 for (PetscInt v = 0; v < 4; ++v) 358 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsD[v] * cdim + d]; 359 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintD, &intPoints[Nsum * cdim])); 360 Nsum += NintD; 361 for (PetscInt v = 0; v < 4; ++v) 362 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsE[v] * cdim + d]; 363 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintE, &intPoints[Nsum * cdim])); 364 Nsum += NintE; 365 for (PetscInt v = 0; v < 4; ++v) 366 for (PetscInt d = 0; d < cdim; ++d) tcoords[v * cdim + d] = coords[vertsF[v] * cdim + d]; 367 PetscCall(DMPlexGetPlaneSimplexIntersection_Coords_Internal(dm, dim, cdim, tcoords, p, normal, pos, &NintF, &intPoints[Nsum * cdim])); 368 Nsum += NintF; 369 *Nint = Nsum; 370 371 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 372 PetscFunctionReturn(PETSC_SUCCESS); 373 } 374 375 /* 376 DMPlexGetPlaneCellIntersection_Internal - Finds the intersection of a plane with a cell 377 378 Not collective 379 380 Input Parameters: 381 + dm - the DM 382 . c - the mesh point 383 . p - a point on the plane. 384 - normal - a normal vector to the plane, must be normalized 385 386 Output Parameters: 387 . pos - `PETSC_TRUE` is the cell is on the positive side of the plane, `PETSC_FALSE` on the negative side 388 + Nint - the number of intersection points, in [0, 4] 389 - intPoints - the coordinates of the intersection points, should be length at least 12 390 391 Note: The `pos` argument is only meaningful if the number of intersections is 0. The algorithmic idea comes from https://github.com/chrisk314/tet-plane-intersection. 392 393 Level: developer 394 395 .seealso: 396 @*/ 397 static PetscErrorCode DMPlexGetPlaneCellIntersection_Internal(DM dm, PetscInt c, const PetscReal p[], const PetscReal normal[], PetscBool *pos, PetscInt *Nint, PetscReal intPoints[]) 398 { 399 DMPolytopeType ct; 400 401 PetscFunctionBegin; 402 PetscCall(DMPlexGetCellType(dm, c, &ct)); 403 switch (ct) { 404 case DM_POLYTOPE_SEGMENT: 405 case DM_POLYTOPE_TRIANGLE: 406 case DM_POLYTOPE_TETRAHEDRON: 407 PetscCall(DMPlexGetPlaneSimplexIntersection_Internal(dm, DMPolytopeTypeGetDim(ct), c, p, normal, pos, Nint, intPoints)); 408 break; 409 case DM_POLYTOPE_QUADRILATERAL: 410 PetscCall(DMPlexGetPlaneQuadIntersection_Internal(dm, DMPolytopeTypeGetDim(ct), c, p, normal, pos, Nint, intPoints)); 411 break; 412 case DM_POLYTOPE_HEXAHEDRON: 413 PetscCall(DMPlexGetPlaneHexIntersection_Internal(dm, DMPolytopeTypeGetDim(ct), c, p, normal, pos, Nint, intPoints)); 414 break; 415 default: 416 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No plane intersection for cell %" PetscInt_FMT " with type %s", c, DMPolytopeTypes[ct]); 417 } 418 PetscFunctionReturn(PETSC_SUCCESS); 419 } 420 421 static PetscErrorCode DMPlexLocatePoint_Simplex_1D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 422 { 423 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 424 const PetscReal x = PetscRealPart(point[0]); 425 PetscReal v0, J, invJ, detJ; 426 PetscReal xi; 427 428 PetscFunctionBegin; 429 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ)); 430 xi = invJ * (x - v0); 431 432 if ((xi >= -eps) && (xi <= 2. + eps)) *cell = c; 433 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 434 PetscFunctionReturn(PETSC_SUCCESS); 435 } 436 437 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 438 { 439 const PetscInt embedDim = 2; 440 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 441 PetscReal x = PetscRealPart(point[0]); 442 PetscReal y = PetscRealPart(point[1]); 443 PetscReal v0[2], J[4], invJ[4], detJ; 444 PetscReal xi, eta; 445 446 PetscFunctionBegin; 447 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 448 xi = invJ[0 * embedDim + 0] * (x - v0[0]) + invJ[0 * embedDim + 1] * (y - v0[1]); 449 eta = invJ[1 * embedDim + 0] * (x - v0[0]) + invJ[1 * embedDim + 1] * (y - v0[1]); 450 451 if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0 + eps)) *cell = c; 452 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 453 PetscFunctionReturn(PETSC_SUCCESS); 454 } 455 456 static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[]) 457 { 458 const PetscInt embedDim = 2; 459 PetscReal x = PetscRealPart(point[0]); 460 PetscReal y = PetscRealPart(point[1]); 461 PetscReal v0[2], J[4], invJ[4], detJ; 462 PetscReal xi, eta, r; 463 464 PetscFunctionBegin; 465 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 466 xi = invJ[0 * embedDim + 0] * (x - v0[0]) + invJ[0 * embedDim + 1] * (y - v0[1]); 467 eta = invJ[1 * embedDim + 0] * (x - v0[0]) + invJ[1 * embedDim + 1] * (y - v0[1]); 468 469 xi = PetscMax(xi, 0.0); 470 eta = PetscMax(eta, 0.0); 471 if (xi + eta > 2.0) { 472 r = (xi + eta) / 2.0; 473 xi /= r; 474 eta /= r; 475 } 476 cpoint[0] = J[0 * embedDim + 0] * xi + J[0 * embedDim + 1] * eta + v0[0]; 477 cpoint[1] = J[1 * embedDim + 0] * xi + J[1 * embedDim + 1] * eta + v0[1]; 478 PetscFunctionReturn(PETSC_SUCCESS); 479 } 480 481 // This is the ray-casting, or even-odd algorithm: https://en.wikipedia.org/wiki/Even%E2%80%93odd_rule 482 static PetscErrorCode DMPlexLocatePoint_Quad_2D_Linear_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 483 { 484 const PetscScalar *array; 485 PetscScalar *coords = NULL; 486 const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 487 PetscReal x = PetscRealPart(point[0]); 488 PetscReal y = PetscRealPart(point[1]); 489 PetscInt crossings = 0, numCoords, f; 490 PetscBool isDG; 491 492 PetscFunctionBegin; 493 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 494 PetscCheck(numCoords == 8, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Quadrilateral should have 8 coordinates, not %" PetscInt_FMT, numCoords); 495 for (f = 0; f < 4; ++f) { 496 PetscReal x_i = PetscRealPart(coords[faces[2 * f + 0] * 2 + 0]); 497 PetscReal y_i = PetscRealPart(coords[faces[2 * f + 0] * 2 + 1]); 498 PetscReal x_j = PetscRealPart(coords[faces[2 * f + 1] * 2 + 0]); 499 PetscReal y_j = PetscRealPart(coords[faces[2 * f + 1] * 2 + 1]); 500 501 if ((x == x_j) && (y == y_j)) { 502 // point is a corner 503 crossings = 1; 504 break; 505 } 506 if ((y_j > y) != (y_i > y)) { 507 PetscReal slope = (x - x_j) * (y_i - y_j) - (x_i - x_j) * (y - y_j); 508 if (slope == 0) { 509 // point is a corner 510 crossings = 1; 511 break; 512 } 513 if ((slope < 0) != (y_i < y_j)) ++crossings; 514 } 515 } 516 if (crossings % 2) *cell = c; 517 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 518 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 519 PetscFunctionReturn(PETSC_SUCCESS); 520 } 521 522 static PetscErrorCode DMPlexLocatePoint_Quad_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 523 { 524 DM cdm; 525 PetscInt degree, dimR, dimC; 526 PetscFE fe; 527 PetscClassId id; 528 PetscSpace sp; 529 PetscReal ref[2], error; 530 Vec coords; 531 PetscBool found = PETSC_FALSE; 532 533 PetscFunctionBegin; 534 PetscCall(DMGetDimension(dm, &dimR)); 535 PetscCall(DMGetCoordinateDM(dm, &cdm)); 536 PetscCall(DMGetDimension(cdm, &dimC)); 537 PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 538 PetscCall(PetscObjectGetClassId((PetscObject)fe, &id)); 539 if (id != PETSCFE_CLASSID) degree = 1; 540 else { 541 PetscCall(PetscFEGetBasisSpace(fe, &sp)); 542 PetscCall(PetscSpaceGetDegree(sp, °ree, NULL)); 543 } 544 if (degree == 1) { 545 /* Use simple location method for linear elements*/ 546 PetscCall(DMPlexLocatePoint_Quad_2D_Linear_Internal(dm, point, c, cell)); 547 PetscFunctionReturn(PETSC_SUCCESS); 548 } 549 /* Otherwise, we have to solve for the real to reference coordinates */ 550 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 551 error = PETSC_SQRT_MACHINE_EPSILON; 552 PetscCall(DMPlexCoordinatesToReference_FE(cdm, fe, c, 1, point, ref, coords, dimC, dimR, 10, &error)); 553 if (error < PETSC_SQRT_MACHINE_EPSILON) found = PETSC_TRUE; 554 if ((ref[0] > 1.0 + PETSC_SMALL) || (ref[0] < -1.0 - PETSC_SMALL) || (ref[1] > 1.0 + PETSC_SMALL) || (ref[1] < -1.0 - PETSC_SMALL)) found = PETSC_FALSE; 555 if (PetscDefined(USE_DEBUG) && found) { 556 PetscReal real[2], inverseError = 0; 557 558 PetscCall(DMPlexReferenceToCoordinates_FE(cdm, fe, c, 1, ref, real, coords, dimC, dimR)); 559 for (PetscInt l = 0; l < dimC; l++) inverseError += (real[l] - PetscRealPart(point[l])) * (real[l] - PetscRealPart(point[l])); 560 inverseError = PetscSqrtReal(inverseError); 561 if (inverseError > PETSC_SQRT_MACHINE_EPSILON) found = PETSC_FALSE; 562 if (!found) PetscCall(PetscInfo(dm, "Point (%g, %g, %g) != Mapped Ref Coords (%g, %g, %g) with error %g\n", (double)point[0], (double)point[1], (double)point[2], (double)real[0], (double)real[1], (double)real[2], (double)inverseError)); 563 } 564 if (found) *cell = c; 565 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 566 PetscFunctionReturn(PETSC_SUCCESS); 567 } 568 569 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 570 { 571 const PetscInt embedDim = 3; 572 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 573 PetscReal v0[3], J[9], invJ[9], detJ; 574 PetscReal x = PetscRealPart(point[0]); 575 PetscReal y = PetscRealPart(point[1]); 576 PetscReal z = PetscRealPart(point[2]); 577 PetscReal xi, eta, zeta; 578 579 PetscFunctionBegin; 580 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 581 xi = invJ[0 * embedDim + 0] * (x - v0[0]) + invJ[0 * embedDim + 1] * (y - v0[1]) + invJ[0 * embedDim + 2] * (z - v0[2]); 582 eta = invJ[1 * embedDim + 0] * (x - v0[0]) + invJ[1 * embedDim + 1] * (y - v0[1]) + invJ[1 * embedDim + 2] * (z - v0[2]); 583 zeta = invJ[2 * embedDim + 0] * (x - v0[0]) + invJ[2 * embedDim + 1] * (y - v0[1]) + invJ[2 * embedDim + 2] * (z - v0[2]); 584 585 if ((xi >= -eps) && (eta >= -eps) && (zeta >= -eps) && (xi + eta + zeta <= 2.0 + eps)) *cell = c; 586 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 587 PetscFunctionReturn(PETSC_SUCCESS); 588 } 589 590 static PetscErrorCode DMPlexLocatePoint_Hex_3D_Linear_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 591 { 592 const PetscScalar *array; 593 PetscScalar *coords = NULL; 594 const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 595 PetscBool found = PETSC_TRUE; 596 PetscInt numCoords, f; 597 PetscBool isDG; 598 599 PetscFunctionBegin; 600 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 601 PetscCheck(numCoords == 24, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Quadrilateral should have 8 coordinates, not %" PetscInt_FMT, numCoords); 602 for (f = 0; f < 6; ++f) { 603 /* Check the point is under plane */ 604 /* Get face normal */ 605 PetscReal v_i[3]; 606 PetscReal v_j[3]; 607 PetscReal normal[3]; 608 PetscReal pp[3]; 609 PetscReal dot; 610 611 v_i[0] = PetscRealPart(coords[faces[f * 4 + 3] * 3 + 0] - coords[faces[f * 4 + 0] * 3 + 0]); 612 v_i[1] = PetscRealPart(coords[faces[f * 4 + 3] * 3 + 1] - coords[faces[f * 4 + 0] * 3 + 1]); 613 v_i[2] = PetscRealPart(coords[faces[f * 4 + 3] * 3 + 2] - coords[faces[f * 4 + 0] * 3 + 2]); 614 v_j[0] = PetscRealPart(coords[faces[f * 4 + 1] * 3 + 0] - coords[faces[f * 4 + 0] * 3 + 0]); 615 v_j[1] = PetscRealPart(coords[faces[f * 4 + 1] * 3 + 1] - coords[faces[f * 4 + 0] * 3 + 1]); 616 v_j[2] = PetscRealPart(coords[faces[f * 4 + 1] * 3 + 2] - coords[faces[f * 4 + 0] * 3 + 2]); 617 normal[0] = v_i[1] * v_j[2] - v_i[2] * v_j[1]; 618 normal[1] = v_i[2] * v_j[0] - v_i[0] * v_j[2]; 619 normal[2] = v_i[0] * v_j[1] - v_i[1] * v_j[0]; 620 pp[0] = PetscRealPart(coords[faces[f * 4 + 0] * 3 + 0] - point[0]); 621 pp[1] = PetscRealPart(coords[faces[f * 4 + 0] * 3 + 1] - point[1]); 622 pp[2] = PetscRealPart(coords[faces[f * 4 + 0] * 3 + 2] - point[2]); 623 dot = normal[0] * pp[0] + normal[1] * pp[1] + normal[2] * pp[2]; 624 625 /* Check that projected point is in face (2D location problem) */ 626 if (dot < 0.0) { 627 found = PETSC_FALSE; 628 break; 629 } 630 } 631 if (found) *cell = c; 632 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 633 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 634 PetscFunctionReturn(PETSC_SUCCESS); 635 } 636 637 static PetscErrorCode DMPlexLocatePoint_Hex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 638 { 639 DM cdm; 640 PetscInt degree, dimR, dimC; 641 PetscFE fe; 642 PetscClassId id; 643 PetscSpace sp; 644 PetscReal ref[3], error; 645 Vec coords; 646 PetscBool found = PETSC_FALSE; 647 648 PetscFunctionBegin; 649 PetscCall(DMGetDimension(dm, &dimR)); 650 PetscCall(DMGetCoordinateDM(dm, &cdm)); 651 PetscCall(DMGetDimension(cdm, &dimC)); 652 PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 653 PetscCall(PetscObjectGetClassId((PetscObject)fe, &id)); 654 if (id != PETSCFE_CLASSID) degree = 1; 655 else { 656 PetscCall(PetscFEGetBasisSpace(fe, &sp)); 657 PetscCall(PetscSpaceGetDegree(sp, °ree, NULL)); 658 } 659 if (degree == 1) { 660 /* Use simple location method for linear elements*/ 661 PetscCall(DMPlexLocatePoint_Hex_3D_Linear_Internal(dm, point, c, cell)); 662 PetscFunctionReturn(PETSC_SUCCESS); 663 } 664 /* Otherwise, we have to solve for the real to reference coordinates */ 665 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 666 error = PETSC_SQRT_MACHINE_EPSILON; 667 PetscCall(DMPlexCoordinatesToReference_FE(cdm, fe, c, 1, point, ref, coords, dimC, dimR, 10, &error)); 668 if (error < PETSC_SQRT_MACHINE_EPSILON) found = PETSC_TRUE; 669 if ((ref[0] > 1.0 + PETSC_SMALL) || (ref[0] < -1.0 - PETSC_SMALL) || (ref[1] > 1.0 + PETSC_SMALL) || (ref[1] < -1.0 - PETSC_SMALL) || (ref[2] > 1.0 + PETSC_SMALL) || (ref[2] < -1.0 - PETSC_SMALL)) found = PETSC_FALSE; 670 if (PetscDefined(USE_DEBUG) && found) { 671 PetscReal real[3], inverseError = 0; 672 673 PetscCall(DMPlexReferenceToCoordinates_FE(cdm, fe, c, 1, ref, real, coords, dimC, dimR)); 674 for (PetscInt l = 0; l < dimC; l++) inverseError += (real[l] - PetscRealPart(point[l])) * (real[l] - PetscRealPart(point[l])); 675 inverseError = PetscSqrtReal(inverseError); 676 if (inverseError > PETSC_SQRT_MACHINE_EPSILON) found = PETSC_FALSE; 677 if (!found) PetscCall(PetscInfo(dm, "Point (%g, %g, %g) != Mapped Ref Coords (%g, %g, %g) with error %g\n", (double)point[0], (double)point[1], (double)point[2], (double)real[0], (double)real[1], (double)real[2], (double)inverseError)); 678 } 679 if (found) *cell = c; 680 else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 681 PetscFunctionReturn(PETSC_SUCCESS); 682 } 683 684 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 685 { 686 PetscInt d; 687 688 PetscFunctionBegin; 689 box->dim = dim; 690 for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = point ? PetscRealPart(point[d]) : 0.; 691 PetscFunctionReturn(PETSC_SUCCESS); 692 } 693 694 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 695 { 696 PetscFunctionBegin; 697 PetscCall(PetscCalloc1(1, box)); 698 PetscCall(PetscGridHashInitialize_Internal(*box, dim, point)); 699 PetscFunctionReturn(PETSC_SUCCESS); 700 } 701 702 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 703 { 704 PetscInt d; 705 706 PetscFunctionBegin; 707 for (d = 0; d < box->dim; ++d) { 708 box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 709 box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 710 } 711 PetscFunctionReturn(PETSC_SUCCESS); 712 } 713 714 static PetscErrorCode DMPlexCreateGridHash(DM dm, PetscGridHash *box) 715 { 716 Vec coordinates; 717 const PetscScalar *a; 718 PetscInt cdim, cStart, cEnd; 719 720 PetscFunctionBegin; 721 PetscCall(DMGetCoordinateDim(dm, &cdim)); 722 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 723 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 724 725 PetscCall(VecGetArrayRead(coordinates, &a)); 726 PetscCall(PetscGridHashCreate(PetscObjectComm((PetscObject)dm), cdim, a, box)); 727 PetscCall(VecRestoreArrayRead(coordinates, &a)); 728 for (PetscInt c = cStart; c < cEnd; ++c) { 729 const PetscScalar *array; 730 PetscScalar *coords = NULL; 731 PetscInt numCoords; 732 PetscBool isDG; 733 734 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 735 for (PetscInt i = 0; i < numCoords / cdim; ++i) PetscCall(PetscGridHashEnlarge(*box, &coords[i * cdim])); 736 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 737 } 738 PetscFunctionReturn(PETSC_SUCCESS); 739 } 740 741 /*@C 742 PetscGridHashSetGrid - Divide the grid into boxes 743 744 Not Collective 745 746 Input Parameters: 747 + box - The grid hash object 748 . n - The number of boxes in each dimension, may use `PETSC_DETERMINE` for the entries 749 - h - The box size in each dimension, only used if n[d] == `PETSC_DETERMINE`, if not needed you can pass in `NULL` 750 751 Level: developer 752 753 .seealso: `DMPLEX`, `PetscGridHashCreate()` 754 @*/ 755 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 756 { 757 PetscInt d; 758 759 PetscFunctionBegin; 760 PetscAssertPointer(n, 2); 761 if (h) PetscAssertPointer(h, 3); 762 for (d = 0; d < box->dim; ++d) { 763 box->extent[d] = box->upper[d] - box->lower[d]; 764 if (n[d] == PETSC_DETERMINE) { 765 PetscCheck(h, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Missing h"); 766 box->h[d] = h[d]; 767 box->n[d] = PetscCeilReal(box->extent[d] / h[d]); 768 } else { 769 box->n[d] = n[d]; 770 box->h[d] = box->extent[d] / n[d]; 771 } 772 } 773 PetscFunctionReturn(PETSC_SUCCESS); 774 } 775 776 /*@C 777 PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point 778 779 Not Collective 780 781 Input Parameters: 782 + box - The grid hash object 783 . numPoints - The number of input points 784 - points - The input point coordinates 785 786 Output Parameters: 787 + dboxes - An array of `numPoints` x `dim` integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 788 - boxes - An array of `numPoints` integers expressing the enclosing box as single number, or `NULL` 789 790 Level: developer 791 792 Note: 793 This only guarantees that a box contains a point, not that a cell does. 794 795 .seealso: `DMPLEX`, `PetscGridHashCreate()` 796 @*/ 797 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 798 { 799 const PetscReal *lower = box->lower; 800 const PetscReal *upper = box->upper; 801 const PetscReal *h = box->h; 802 const PetscInt *n = box->n; 803 const PetscInt dim = box->dim; 804 PetscInt d, p; 805 806 PetscFunctionBegin; 807 for (p = 0; p < numPoints; ++p) { 808 for (d = 0; d < dim; ++d) { 809 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p * dim + d]) - lower[d]) / h[d]); 810 811 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p * dim + d]) - upper[d]) < 1.0e-9) dbox = n[d] - 1; 812 if (dbox == -1 && PetscAbsReal(PetscRealPart(points[p * dim + d]) - lower[d]) < 1.0e-9) dbox = 0; 813 PetscCheck(dbox >= 0 && dbox < n[d], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %" PetscInt_FMT " (%g, %g, %g) is outside of our bounding box (%g, %g, %g) - (%g, %g, %g)", 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, (double)lower[0], (double)lower[1], (double)lower[2], (double)upper[0], (double)upper[1], (double)upper[2]); 814 dboxes[p * dim + d] = dbox; 815 } 816 if (boxes) 817 for (d = dim - 2, boxes[p] = dboxes[p * dim + dim - 1]; d >= 0; --d) boxes[p] = boxes[p] * n[d] + dboxes[p * dim + d]; 818 } 819 PetscFunctionReturn(PETSC_SUCCESS); 820 } 821 822 /* 823 PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point 824 825 Not Collective 826 827 Input Parameters: 828 + box - The grid hash object 829 . cellSection - The PetscSection mapping cells to boxes 830 . numPoints - The number of input points 831 - points - The input point coordinates 832 833 Output Parameters: 834 + dboxes - An array of `numPoints`*`dim` integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 835 . boxes - An array of `numPoints` integers expressing the enclosing box as single number, or `NULL` 836 - found - Flag indicating if point was located within a box 837 838 Level: developer 839 840 Note: 841 This does an additional check that a cell actually contains the point, and found is `PETSC_FALSE` if no cell does. Thus, this function requires that `cellSection` is already constructed. 842 843 .seealso: `DMPLEX`, `PetscGridHashGetEnclosingBox()` 844 */ 845 static PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscSection cellSection, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[], PetscBool *found) 846 { 847 const PetscReal *lower = box->lower; 848 const PetscReal *upper = box->upper; 849 const PetscReal *h = box->h; 850 const PetscInt *n = box->n; 851 const PetscInt dim = box->dim; 852 PetscInt bStart, bEnd, d, p; 853 854 PetscFunctionBegin; 855 PetscValidHeaderSpecific(cellSection, PETSC_SECTION_CLASSID, 2); 856 *found = PETSC_FALSE; 857 PetscCall(PetscSectionGetChart(box->cellSection, &bStart, &bEnd)); 858 for (p = 0; p < numPoints; ++p) { 859 for (d = 0; d < dim; ++d) { 860 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p * dim + d]) - lower[d]) / h[d]); 861 862 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p * dim + d]) - upper[d]) < 1.0e-9) dbox = n[d] - 1; 863 if (dbox < 0 || dbox >= n[d]) PetscFunctionReturn(PETSC_SUCCESS); 864 dboxes[p * dim + d] = dbox; 865 } 866 if (boxes) 867 for (d = dim - 2, boxes[p] = dboxes[p * dim + dim - 1]; d >= 0; --d) boxes[p] = boxes[p] * n[d] + dboxes[p * dim + d]; 868 // It is possible for a box to overlap no grid cells 869 if (boxes[p] < bStart || boxes[p] >= bEnd) PetscFunctionReturn(PETSC_SUCCESS); 870 } 871 *found = PETSC_TRUE; 872 PetscFunctionReturn(PETSC_SUCCESS); 873 } 874 875 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 876 { 877 PetscFunctionBegin; 878 if (*box) { 879 PetscCall(PetscSectionDestroy(&(*box)->cellSection)); 880 PetscCall(ISDestroy(&(*box)->cells)); 881 PetscCall(DMLabelDestroy(&(*box)->cellsSparse)); 882 } 883 PetscCall(PetscFree(*box)); 884 PetscFunctionReturn(PETSC_SUCCESS); 885 } 886 887 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 888 { 889 DMPolytopeType ct; 890 891 PetscFunctionBegin; 892 PetscCall(DMPlexGetCellType(dm, cellStart, &ct)); 893 switch (ct) { 894 case DM_POLYTOPE_SEGMENT: 895 PetscCall(DMPlexLocatePoint_Simplex_1D_Internal(dm, point, cellStart, cell)); 896 break; 897 case DM_POLYTOPE_TRIANGLE: 898 PetscCall(DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell)); 899 break; 900 case DM_POLYTOPE_QUADRILATERAL: 901 PetscCall(DMPlexLocatePoint_Quad_2D_Internal(dm, point, cellStart, cell)); 902 break; 903 case DM_POLYTOPE_TETRAHEDRON: 904 PetscCall(DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell)); 905 break; 906 case DM_POLYTOPE_HEXAHEDRON: 907 PetscCall(DMPlexLocatePoint_Hex_3D_Internal(dm, point, cellStart, cell)); 908 break; 909 default: 910 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell %" PetscInt_FMT " with type %s", cellStart, DMPolytopeTypes[ct]); 911 } 912 PetscFunctionReturn(PETSC_SUCCESS); 913 } 914 915 /* 916 DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point 917 */ 918 static PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[]) 919 { 920 DMPolytopeType ct; 921 922 PetscFunctionBegin; 923 PetscCall(DMPlexGetCellType(dm, cell, &ct)); 924 switch (ct) { 925 case DM_POLYTOPE_TRIANGLE: 926 PetscCall(DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint)); 927 break; 928 #if 0 929 case DM_POLYTOPE_QUADRILATERAL: 930 PetscCall(DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint));break; 931 case DM_POLYTOPE_TETRAHEDRON: 932 PetscCall(DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint));break; 933 case DM_POLYTOPE_HEXAHEDRON: 934 PetscCall(DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint));break; 935 #endif 936 default: 937 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell %" PetscInt_FMT " with type %s", cell, DMPolytopeTypes[ct]); 938 } 939 PetscFunctionReturn(PETSC_SUCCESS); 940 } 941 942 /* 943 DMPlexComputeGridHash_Internal - Create a grid hash structure covering the `DMPLEX` 944 945 Collective 946 947 Input Parameter: 948 . dm - The `DMPLEX` 949 950 Output Parameter: 951 . localBox - The grid hash object 952 953 Level: developer 954 955 Notes: 956 How do we determine all boxes intersecting a given cell? 957 958 1) Get convex body enclosing cell. We will use a box called the box-hull. 959 960 2) Get smallest brick of boxes enclosing the box-hull 961 962 3) Each box is composed of 6 planes, 3 lower and 3 upper. We loop over dimensions, and 963 for each new plane determine whether the cell is on the negative side, positive side, or intersects it. 964 965 a) If the cell is on the negative side of the lower planes, it is not in the box 966 967 b) If the cell is on the positive side of the upper planes, it is not in the box 968 969 c) If there is no intersection, it is in the box 970 971 d) If any intersection point is within the box limits, it is in the box 972 973 .seealso: `DMPLEX`, `PetscGridHashCreate()`, `PetscGridHashGetEnclosingBox()` 974 */ 975 static PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 976 { 977 PetscInt debug = ((DM_Plex *)dm->data)->printLocate; 978 PetscGridHash lbox; 979 PetscSF sf; 980 const PetscInt *leaves; 981 PetscInt *dboxes, *boxes; 982 PetscInt cdim, cStart, cEnd, Nl = -1; 983 PetscBool flg; 984 985 PetscFunctionBegin; 986 PetscCall(DMGetCoordinateDim(dm, &cdim)); 987 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 988 PetscCall(DMPlexCreateGridHash(dm, &lbox)); 989 { 990 PetscInt n[3], d; 991 992 PetscCall(PetscOptionsGetIntArray(NULL, ((PetscObject)dm)->prefix, "-dm_plex_hash_box_faces", n, &d, &flg)); 993 if (flg) { 994 for (PetscInt i = d; i < cdim; ++i) n[i] = n[d - 1]; 995 } else { 996 for (PetscInt i = 0; i < cdim; ++i) n[i] = PetscMax(2, PetscFloorReal(PetscPowReal((PetscReal)(cEnd - cStart), 1.0 / cdim) * 0.8)); 997 } 998 PetscCall(PetscGridHashSetGrid(lbox, n, NULL)); 999 if (debug) 1000 PetscCall(PetscPrintf(PETSC_COMM_SELF, "GridHash:\n (%g, %g, %g) -- (%g, %g, %g)\n n %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n h %g %g %g\n", (double)lbox->lower[0], (double)lbox->lower[1], cdim > 2 ? (double)lbox->lower[2] : 0., 1001 (double)lbox->upper[0], (double)lbox->upper[1], cdim > 2 ? (double)lbox->upper[2] : 0, n[0], n[1], cdim > 2 ? n[2] : 0, (double)lbox->h[0], (double)lbox->h[1], cdim > 2 ? (double)lbox->h[2] : 0.)); 1002 } 1003 1004 PetscCall(DMGetPointSF(dm, &sf)); 1005 if (sf) PetscCall(PetscSFGetGraph(sf, NULL, &Nl, &leaves, NULL)); 1006 Nl = PetscMax(Nl, 0); 1007 PetscCall(PetscCalloc2(16 * cdim, &dboxes, 16, &boxes)); 1008 1009 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "cells", &lbox->cellsSparse)); 1010 PetscCall(DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd)); 1011 for (PetscInt c = cStart; c < cEnd; ++c) { 1012 PetscReal intPoints[6 * 6 * 6 * 3]; 1013 const PetscScalar *array; 1014 PetscScalar *coords = NULL; 1015 const PetscReal *h = lbox->h; 1016 PetscReal normal[9] = {1., 0., 0., 0., 1., 0., 0., 0., 1.}; 1017 PetscReal *lowerIntPoints[3] = {&intPoints[0 * 6 * 6 * 3], &intPoints[1 * 6 * 6 * 3], &intPoints[2 * 6 * 6 * 3]}; 1018 PetscReal *upperIntPoints[3] = {&intPoints[3 * 6 * 6 * 3], &intPoints[4 * 6 * 6 * 3], &intPoints[5 * 6 * 6 * 3]}; 1019 PetscReal lp[3], up[3], *tmp; 1020 PetscInt numCoords, idx, dlim[6], lowerInt[3], upperInt[3]; 1021 PetscBool isDG, lower[3], upper[3]; 1022 1023 PetscCall(PetscFindInt(c, Nl, leaves, &idx)); 1024 if (idx >= 0) continue; 1025 // Get grid of boxes containing the cell 1026 PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 1027 PetscCall(PetscGridHashGetEnclosingBox(lbox, numCoords / cdim, coords, dboxes, boxes)); 1028 PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &coords)); 1029 for (PetscInt d = 0; d < cdim; ++d) dlim[d * 2 + 0] = dlim[d * 2 + 1] = dboxes[d]; 1030 for (PetscInt d = cdim; d < 3; ++d) dlim[d * 2 + 0] = dlim[d * 2 + 1] = 0; 1031 for (PetscInt e = 1; e < numCoords / cdim; ++e) { 1032 for (PetscInt d = 0; d < cdim; ++d) { 1033 dlim[d * 2 + 0] = PetscMin(dlim[d * 2 + 0], dboxes[e * cdim + d]); 1034 dlim[d * 2 + 1] = PetscMax(dlim[d * 2 + 1], dboxes[e * cdim + d]); 1035 } 1036 } 1037 if (debug > 4) { 1038 for (PetscInt d = 0; d < cdim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " direction %" PetscInt_FMT " box limits %" PetscInt_FMT "--%" PetscInt_FMT "\n", c, d, dlim[d * 2 + 0], dlim[d * 2 + 1])); 1039 } 1040 // Initialize with lower planes for first box 1041 for (PetscInt d = 0; d < cdim; ++d) { 1042 lp[d] = lbox->lower[d] + dlim[d * 2 + 0] * h[d]; 1043 up[d] = lp[d] + h[d]; 1044 } 1045 for (PetscInt d = 0; d < cdim; ++d) { 1046 PetscCall(DMPlexGetPlaneCellIntersection_Internal(dm, c, lp, &normal[d * 3], &lower[d], &lowerInt[d], lowerIntPoints[d])); 1047 if (debug > 4) { 1048 if (!lowerInt[d]) 1049 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " lower direction %" PetscInt_FMT " (%g, %g, %g) does not intersect %s\n", c, d, (double)lp[0], (double)lp[1], cdim > 2 ? (double)lp[2] : 0., lower[d] ? "positive" : "negative")); 1050 else PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " lower direction %" PetscInt_FMT " (%g, %g, %g) intersects %" PetscInt_FMT " times\n", c, d, (double)lp[0], (double)lp[1], cdim > 2 ? (double)lp[2] : 0., lowerInt[d])); 1051 } 1052 } 1053 // Loop over grid 1054 for (PetscInt k = dlim[2 * 2 + 0]; k <= dlim[2 * 2 + 1]; ++k, lp[2] = up[2], up[2] += h[2]) { 1055 if (cdim > 2) PetscCall(DMPlexGetPlaneCellIntersection_Internal(dm, c, up, &normal[3 * 2], &upper[2], &upperInt[2], upperIntPoints[2])); 1056 if (cdim > 2 && debug > 4) { 1057 if (!upperInt[2]) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 2 (%g, %g, %g) does not intersect %s\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upper[2] ? "positive" : "negative")); 1058 else PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 2 (%g, %g, %g) intersects %" PetscInt_FMT " times\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upperInt[2])); 1059 } 1060 for (PetscInt j = dlim[1 * 2 + 0]; j <= dlim[1 * 2 + 1]; ++j, lp[1] = up[1], up[1] += h[1]) { 1061 if (cdim > 1) PetscCall(DMPlexGetPlaneCellIntersection_Internal(dm, c, up, &normal[3 * 1], &upper[1], &upperInt[1], upperIntPoints[1])); 1062 if (cdim > 1 && debug > 4) { 1063 if (!upperInt[1]) 1064 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 1 (%g, %g, %g) does not intersect %s\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upper[1] ? "positive" : "negative")); 1065 else PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 1 (%g, %g, %g) intersects %" PetscInt_FMT " times\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upperInt[1])); 1066 } 1067 for (PetscInt i = dlim[0 * 2 + 0]; i <= dlim[0 * 2 + 1]; ++i, lp[0] = up[0], up[0] += h[0]) { 1068 const PetscInt box = (k * lbox->n[1] + j) * lbox->n[0] + i; 1069 PetscBool excNeg = PETSC_TRUE; 1070 PetscBool excPos = PETSC_TRUE; 1071 PetscInt NlInt = 0; 1072 PetscInt NuInt = 0; 1073 1074 PetscCall(DMPlexGetPlaneCellIntersection_Internal(dm, c, up, &normal[3 * 0], &upper[0], &upperInt[0], upperIntPoints[0])); 1075 if (debug > 4) { 1076 if (!upperInt[0]) 1077 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 0 (%g, %g, %g) does not intersect %s\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upper[0] ? "positive" : "negative")); 1078 else PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " upper direction 0 (%g, %g, %g) intersects %" PetscInt_FMT " times\n", c, (double)up[0], (double)up[1], cdim > 2 ? (double)up[2] : 0., upperInt[0])); 1079 } 1080 for (PetscInt d = 0; d < cdim; ++d) { 1081 NlInt += lowerInt[d]; 1082 NuInt += upperInt[d]; 1083 } 1084 // If there is no intersection... 1085 if (!NlInt && !NuInt) { 1086 // If the cell is on the negative side of the lower planes, it is not in the box 1087 for (PetscInt d = 0; d < cdim; ++d) 1088 if (lower[d]) { 1089 excNeg = PETSC_FALSE; 1090 break; 1091 } 1092 // If the cell is on the positive side of the upper planes, it is not in the box 1093 for (PetscInt d = 0; d < cdim; ++d) 1094 if (!upper[d]) { 1095 excPos = PETSC_FALSE; 1096 break; 1097 } 1098 if (excNeg || excPos) { 1099 if (debug && excNeg) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " is on the negative side of the lower plane\n", c)); 1100 if (debug && excPos) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " is on the positive side of the upper plane\n", c)); 1101 continue; 1102 } 1103 // Otherwise it is in the box 1104 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " is contained in box %" PetscInt_FMT "\n", c, box)); 1105 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1106 continue; 1107 } 1108 /* 1109 If any intersection point is within the box limits, it is in the box 1110 We need to have tolerances here since intersection point calculations can introduce errors 1111 Initialize a count to track which planes have intersection outside the box. 1112 if two adjacent planes have intersection points upper and lower all outside the box, look 1113 first at if another plane has intersection points outside the box, if so, it is inside the cell 1114 look next if no intersection points exist on the other planes, and check if the planes are on the 1115 outside of the intersection points but on opposite ends. If so, the box cuts through the cell. 1116 */ 1117 PetscInt outsideCount[6] = {0, 0, 0, 0, 0, 0}; 1118 for (PetscInt plane = 0; plane < cdim; ++plane) { 1119 for (PetscInt ip = 0; ip < lowerInt[plane]; ++ip) { 1120 PetscInt d; 1121 1122 for (d = 0; d < cdim; ++d) { 1123 if ((lowerIntPoints[plane][ip * cdim + d] < (lp[d] - PETSC_SMALL)) || (lowerIntPoints[plane][ip * cdim + d] > (up[d] + PETSC_SMALL))) { 1124 if (lowerIntPoints[plane][ip * cdim + d] < (lp[d] - PETSC_SMALL)) outsideCount[d]++; // The lower point is to the left of this box, and we count it 1125 break; 1126 } 1127 } 1128 if (d == cdim) { 1129 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " intersected lower plane %" PetscInt_FMT " of box %" PetscInt_FMT "\n", c, plane, box)); 1130 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1131 goto end; 1132 } 1133 } 1134 for (PetscInt ip = 0; ip < upperInt[plane]; ++ip) { 1135 PetscInt d; 1136 1137 for (d = 0; d < cdim; ++d) { 1138 if ((upperIntPoints[plane][ip * cdim + d] < (lp[d] - PETSC_SMALL)) || (upperIntPoints[plane][ip * cdim + d] > (up[d] + PETSC_SMALL))) { 1139 if (upperIntPoints[plane][ip * cdim + d] > (up[d] + PETSC_SMALL)) outsideCount[cdim + d]++; // The upper point is to the right of this box, and we count it 1140 break; 1141 } 1142 } 1143 if (d == cdim) { 1144 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Cell %" PetscInt_FMT " intersected upper plane %" PetscInt_FMT " of box %" PetscInt_FMT "\n", c, plane, box)); 1145 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1146 goto end; 1147 } 1148 } 1149 } 1150 /* 1151 Check the planes with intersections 1152 in 2D, check if the square falls in the middle of a cell 1153 ie all four planes have intersection points outside of the box 1154 You do not want to be doing this, because it means your grid hashing is finer than your grid, 1155 but we should still support it I guess 1156 */ 1157 if (cdim == 2) { 1158 PetscInt nIntersects = 0; 1159 for (PetscInt d = 0; d < cdim; ++d) nIntersects += (outsideCount[d] + outsideCount[d + cdim]); 1160 // if the count adds up to 8, that means each plane has 2 external intersections and thus it is in the cell 1161 if (nIntersects == 8) { 1162 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1163 goto end; 1164 } 1165 } 1166 /* 1167 In 3 dimensions, if two adjacent planes have at least 3 intersections outside the cell in the appropriate direction, 1168 we then check the 3rd planar dimension. If a plane falls between intersection points, the cell belongs to that box. 1169 If the planes are on opposite sides of the intersection points, the cell belongs to that box and it passes through the cell. 1170 */ 1171 if (cdim == 3) { 1172 PetscInt faces[3] = {0, 0, 0}, checkInternalFace = 0; 1173 // Find two adjacent planes with at least 3 intersection points in the upper and lower 1174 // if the third plane has 3 intersection points or more, a pyramid base is formed on that plane and it is in the cell 1175 for (PetscInt d = 0; d < cdim; ++d) 1176 if (outsideCount[d] >= 3 && outsideCount[cdim + d] >= 3) { 1177 faces[d]++; 1178 checkInternalFace++; 1179 } 1180 if (checkInternalFace == 3) { 1181 // All planes have 3 intersection points, add it. 1182 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1183 goto end; 1184 } 1185 // Gross, figure out which adjacent faces have at least 3 points 1186 PetscInt nonIntersectingFace = -1; 1187 if (faces[0] == faces[1]) nonIntersectingFace = 2; 1188 if (faces[0] == faces[2]) nonIntersectingFace = 1; 1189 if (faces[1] == faces[2]) nonIntersectingFace = 0; 1190 if (nonIntersectingFace >= 0) { 1191 for (PetscInt plane = 0; plane < cdim; ++plane) { 1192 if (!lowerInt[nonIntersectingFace] && !upperInt[nonIntersectingFace]) continue; 1193 // If we have 2 adjacent sides with pyramids of intersection outside of them, and there is a point between the end caps at all, it must be between the two non intersecting ends, and the box is inside the cell. 1194 for (PetscInt ip = 0; ip < lowerInt[nonIntersectingFace]; ++ip) { 1195 if (lowerIntPoints[plane][ip * cdim + nonIntersectingFace] > lp[nonIntersectingFace] - PETSC_SMALL || lowerIntPoints[plane][ip * cdim + nonIntersectingFace] < up[nonIntersectingFace] + PETSC_SMALL) goto setpoint; 1196 } 1197 for (PetscInt ip = 0; ip < upperInt[nonIntersectingFace]; ++ip) { 1198 if (upperIntPoints[plane][ip * cdim + nonIntersectingFace] > lp[nonIntersectingFace] - PETSC_SMALL || upperIntPoints[plane][ip * cdim + nonIntersectingFace] < up[nonIntersectingFace] + PETSC_SMALL) goto setpoint; 1199 } 1200 goto end; 1201 } 1202 // The points are within the bonds of the non intersecting planes, add it. 1203 setpoint: 1204 PetscCall(DMLabelSetValue(lbox->cellsSparse, c, box)); 1205 goto end; 1206 } 1207 } 1208 end: 1209 lower[0] = upper[0]; 1210 lowerInt[0] = upperInt[0]; 1211 tmp = lowerIntPoints[0]; 1212 lowerIntPoints[0] = upperIntPoints[0]; 1213 upperIntPoints[0] = tmp; 1214 } 1215 lp[0] = lbox->lower[0] + dlim[0 * 2 + 0] * h[0]; 1216 up[0] = lp[0] + h[0]; 1217 lower[1] = upper[1]; 1218 lowerInt[1] = upperInt[1]; 1219 tmp = lowerIntPoints[1]; 1220 lowerIntPoints[1] = upperIntPoints[1]; 1221 upperIntPoints[1] = tmp; 1222 } 1223 lp[1] = lbox->lower[1] + dlim[1 * 2 + 0] * h[1]; 1224 up[1] = lp[1] + h[1]; 1225 lower[2] = upper[2]; 1226 lowerInt[2] = upperInt[2]; 1227 tmp = lowerIntPoints[2]; 1228 lowerIntPoints[2] = upperIntPoints[2]; 1229 upperIntPoints[2] = tmp; 1230 } 1231 } 1232 PetscCall(PetscFree2(dboxes, boxes)); 1233 1234 if (debug) PetscCall(DMLabelView(lbox->cellsSparse, PETSC_VIEWER_STDOUT_SELF)); 1235 PetscCall(DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells)); 1236 PetscCall(DMLabelDestroy(&lbox->cellsSparse)); 1237 *localBox = lbox; 1238 PetscFunctionReturn(PETSC_SUCCESS); 1239 } 1240 1241 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF) 1242 { 1243 PetscInt debug = ((DM_Plex *)dm->data)->printLocate; 1244 DM_Plex *mesh = (DM_Plex *)dm->data; 1245 PetscBool hash = mesh->useHashLocation, reuse = PETSC_FALSE; 1246 PetscInt bs, numPoints, p, numFound, *found = NULL; 1247 PetscInt dim, Nl = 0, cStart, cEnd, numCells, c, d; 1248 PetscSF sf; 1249 const PetscInt *leaves; 1250 const PetscInt *boxCells; 1251 PetscSFNode *cells; 1252 PetscScalar *a; 1253 PetscMPIInt result; 1254 PetscLogDouble t0, t1; 1255 PetscReal gmin[3], gmax[3]; 1256 PetscInt terminating_query_type[] = {0, 0, 0}; 1257 PetscMPIInt rank; 1258 1259 PetscFunctionBegin; 1260 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1261 PetscCall(PetscLogEventBegin(DMPLEX_LocatePoints, 0, 0, 0, 0)); 1262 PetscCall(PetscTime(&t0)); 1263 PetscCheck(ltype != DM_POINTLOCATION_NEAREST || hash, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it."); 1264 PetscCall(DMGetCoordinateDim(dm, &dim)); 1265 PetscCall(VecGetBlockSize(v, &bs)); 1266 PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF), PETSC_COMM_SELF, &result)); 1267 PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PetscObjectComm((PetscObject)cellSF), PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 1268 // We ignore extra coordinates 1269 PetscCheck(bs >= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Block size for point vector %" PetscInt_FMT " must be the mesh coordinate dimension %" PetscInt_FMT, bs, dim); 1270 PetscCall(DMGetCoordinatesLocalSetUp(dm)); 1271 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 1272 PetscCall(DMGetPointSF(dm, &sf)); 1273 if (sf) PetscCall(PetscSFGetGraph(sf, NULL, &Nl, &leaves, NULL)); 1274 Nl = PetscMax(Nl, 0); 1275 PetscCall(VecGetLocalSize(v, &numPoints)); 1276 PetscCall(VecGetArray(v, &a)); 1277 numPoints /= bs; 1278 { 1279 const PetscSFNode *sf_cells; 1280 1281 PetscCall(PetscSFGetGraph(cellSF, NULL, NULL, NULL, &sf_cells)); 1282 if (sf_cells) { 1283 PetscCall(PetscInfo(dm, "[DMLocatePoints_Plex] Re-using existing StarForest node list\n")); 1284 cells = (PetscSFNode *)sf_cells; 1285 reuse = PETSC_TRUE; 1286 } else { 1287 PetscCall(PetscInfo(dm, "[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n")); 1288 PetscCall(PetscMalloc1(numPoints, &cells)); 1289 /* initialize cells if created */ 1290 for (p = 0; p < numPoints; p++) { 1291 cells[p].rank = 0; 1292 cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 1293 } 1294 } 1295 } 1296 PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 1297 if (hash) { 1298 if (!mesh->lbox) { 1299 PetscCall(PetscInfo(dm, "Initializing grid hashing\n")); 1300 PetscCall(DMPlexComputeGridHash_Internal(dm, &mesh->lbox)); 1301 } 1302 /* Designate the local box for each point */ 1303 /* Send points to correct process */ 1304 /* Search cells that lie in each subbox */ 1305 /* Should we bin points before doing search? */ 1306 PetscCall(ISGetIndices(mesh->lbox->cells, &boxCells)); 1307 } 1308 for (p = 0, numFound = 0; p < numPoints; ++p) { 1309 const PetscScalar *point = &a[p * bs]; 1310 PetscInt dbin[3] = {-1, -1, -1}, bin, cell = -1, cellOffset; 1311 PetscBool point_outside_domain = PETSC_FALSE; 1312 1313 /* check bounding box of domain */ 1314 for (d = 0; d < dim; d++) { 1315 if (PetscRealPart(point[d]) < gmin[d]) { 1316 point_outside_domain = PETSC_TRUE; 1317 break; 1318 } 1319 if (PetscRealPart(point[d]) > gmax[d]) { 1320 point_outside_domain = PETSC_TRUE; 1321 break; 1322 } 1323 } 1324 if (point_outside_domain) { 1325 cells[p].rank = 0; 1326 cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 1327 terminating_query_type[0]++; 1328 continue; 1329 } 1330 1331 /* check initial values in cells[].index - abort early if found */ 1332 if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1333 c = cells[p].index; 1334 cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 1335 PetscCall(DMPlexLocatePoint_Internal(dm, dim, point, c, &cell)); 1336 if (cell >= 0) { 1337 cells[p].rank = 0; 1338 cells[p].index = cell; 1339 numFound++; 1340 } 1341 } 1342 if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1343 terminating_query_type[1]++; 1344 continue; 1345 } 1346 1347 if (hash) { 1348 PetscBool found_box; 1349 1350 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Checking point %" PetscInt_FMT " (%.2g, %.2g, %.2g)\n", rank, p, (double)PetscRealPart(point[0]), (double)PetscRealPart(point[1]), dim > 2 ? (double)PetscRealPart(point[2]) : 0.)); 1351 /* allow for case that point is outside box - abort early */ 1352 PetscCall(PetscGridHashGetEnclosingBoxQuery(mesh->lbox, mesh->lbox->cellSection, 1, point, dbin, &bin, &found_box)); 1353 if (found_box) { 1354 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Found point in box %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", rank, bin, dbin[0], dbin[1], dim > 2 ? dbin[2] : 0)); 1355 /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 1356 PetscCall(PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells)); 1357 PetscCall(PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset)); 1358 for (c = cellOffset; c < cellOffset + numCells; ++c) { 1359 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Checking for point in cell %" PetscInt_FMT "\n", rank, boxCells[c])); 1360 PetscCall(DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell)); 1361 if (cell >= 0) { 1362 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] FOUND in cell %" PetscInt_FMT "\n", rank, cell)); 1363 cells[p].rank = 0; 1364 cells[p].index = cell; 1365 numFound++; 1366 terminating_query_type[2]++; 1367 break; 1368 } 1369 } 1370 } 1371 } else { 1372 PetscBool found = PETSC_FALSE; 1373 for (c = cStart; c < cEnd; ++c) { 1374 PetscInt idx; 1375 1376 PetscCall(PetscFindInt(c, Nl, leaves, &idx)); 1377 if (idx >= 0) continue; 1378 PetscCall(DMPlexLocatePoint_Internal(dm, dim, point, c, &cell)); 1379 if (cell >= 0) { 1380 cells[p].rank = 0; 1381 cells[p].index = cell; 1382 numFound++; 1383 terminating_query_type[2]++; 1384 found = PETSC_TRUE; 1385 break; 1386 } 1387 } 1388 if (!found) terminating_query_type[0]++; 1389 } 1390 } 1391 if (hash) PetscCall(ISRestoreIndices(mesh->lbox->cells, &boxCells)); 1392 if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) { 1393 for (p = 0; p < numPoints; p++) { 1394 const PetscScalar *point = &a[p * bs]; 1395 PetscReal cpoint[3] = {0, 0, 0}, diff[3], best[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}, dist, distMax = PETSC_MAX_REAL; 1396 PetscInt dbin[3] = {-1, -1, -1}, bin, cellOffset, d, bestc = -1; 1397 1398 if (cells[p].index < 0) { 1399 PetscCall(PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin)); 1400 PetscCall(PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells)); 1401 PetscCall(PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset)); 1402 for (c = cellOffset; c < cellOffset + numCells; ++c) { 1403 PetscCall(DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint)); 1404 for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]); 1405 dist = DMPlex_NormD_Internal(dim, diff); 1406 if (dist < distMax) { 1407 for (d = 0; d < dim; ++d) best[d] = cpoint[d]; 1408 bestc = boxCells[c]; 1409 distMax = dist; 1410 } 1411 } 1412 if (distMax < PETSC_MAX_REAL) { 1413 ++numFound; 1414 cells[p].rank = 0; 1415 cells[p].index = bestc; 1416 for (d = 0; d < dim; ++d) a[p * bs + d] = best[d]; 1417 } 1418 } 1419 } 1420 } 1421 /* This code is only be relevant when interfaced to parallel point location */ 1422 /* Check for highest numbered proc that claims a point (do we care?) */ 1423 if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) { 1424 PetscCall(PetscMalloc1(numFound, &found)); 1425 for (p = 0, numFound = 0; p < numPoints; p++) { 1426 if (cells[p].rank >= 0 && cells[p].index >= 0) { 1427 if (numFound < p) cells[numFound] = cells[p]; 1428 found[numFound++] = p; 1429 } 1430 } 1431 } 1432 PetscCall(VecRestoreArray(v, &a)); 1433 if (!reuse) PetscCall(PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER)); 1434 PetscCall(PetscTime(&t1)); 1435 if (hash) { 1436 PetscCall(PetscInfo(dm, "[DMLocatePoints_Plex] terminating_query_type : %" PetscInt_FMT " [outside domain] : %" PetscInt_FMT " [inside initial cell] : %" PetscInt_FMT " [hash]\n", terminating_query_type[0], terminating_query_type[1], terminating_query_type[2])); 1437 } else { 1438 PetscCall(PetscInfo(dm, "[DMLocatePoints_Plex] terminating_query_type : %" PetscInt_FMT " [outside domain] : %" PetscInt_FMT " [inside initial cell] : %" PetscInt_FMT " [brute-force]\n", terminating_query_type[0], terminating_query_type[1], terminating_query_type[2])); 1439 } 1440 PetscCall(PetscInfo(dm, "[DMLocatePoints_Plex] npoints %" PetscInt_FMT " : time(rank0) %1.2e (sec): points/sec %1.4e\n", numPoints, t1 - t0, numPoints / (t1 - t0))); 1441 PetscCall(PetscLogEventEnd(DMPLEX_LocatePoints, 0, 0, 0, 0)); 1442 PetscFunctionReturn(PETSC_SUCCESS); 1443 } 1444 1445 /*@ 1446 DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates 1447 1448 Not Collective 1449 1450 Input/Output Parameter: 1451 . coords - The coordinates of a segment, on output the new y-coordinate, and 0 for x, an array of size 4, last two entries are unchanged 1452 1453 Output Parameter: 1454 . R - The rotation which accomplishes the projection, array of size 4 1455 1456 Level: developer 1457 1458 .seealso: `DMPLEX`, `DMPlexComputeProjection3Dto1D()`, `DMPlexComputeProjection3Dto2D()` 1459 @*/ 1460 PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[]) 1461 { 1462 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 1463 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 1464 const PetscReal r = PetscSqrtReal(x * x + y * y), c = x / r, s = y / r; 1465 1466 PetscFunctionBegin; 1467 R[0] = c; 1468 R[1] = -s; 1469 R[2] = s; 1470 R[3] = c; 1471 coords[0] = 0.0; 1472 coords[1] = r; 1473 PetscFunctionReturn(PETSC_SUCCESS); 1474 } 1475 1476 /*@ 1477 DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates 1478 1479 Not Collective 1480 1481 Input/Output Parameter: 1482 . coords - The coordinates of a segment; on output, the new y-coordinate, and 0 for x and z, an array of size 6, the other entries are unchanged 1483 1484 Output Parameter: 1485 . R - The rotation which accomplishes the projection, an array of size 9 1486 1487 Level: developer 1488 1489 Note: 1490 This uses the basis completion described by Frisvad {cite}`frisvad2012building` 1491 1492 .seealso: `DMPLEX`, `DMPlexComputeProjection2Dto1D()`, `DMPlexComputeProjection3Dto2D()` 1493 @*/ 1494 PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[]) 1495 { 1496 PetscReal x = PetscRealPart(coords[3] - coords[0]); 1497 PetscReal y = PetscRealPart(coords[4] - coords[1]); 1498 PetscReal z = PetscRealPart(coords[5] - coords[2]); 1499 PetscReal r = PetscSqrtReal(x * x + y * y + z * z); 1500 PetscReal rinv = 1. / r; 1501 1502 PetscFunctionBegin; 1503 x *= rinv; 1504 y *= rinv; 1505 z *= rinv; 1506 if (x > 0.) { 1507 PetscReal inv1pX = 1. / (1. + x); 1508 1509 R[0] = x; 1510 R[1] = -y; 1511 R[2] = -z; 1512 R[3] = y; 1513 R[4] = 1. - y * y * inv1pX; 1514 R[5] = -y * z * inv1pX; 1515 R[6] = z; 1516 R[7] = -y * z * inv1pX; 1517 R[8] = 1. - z * z * inv1pX; 1518 } else { 1519 PetscReal inv1mX = 1. / (1. - x); 1520 1521 R[0] = x; 1522 R[1] = z; 1523 R[2] = y; 1524 R[3] = y; 1525 R[4] = -y * z * inv1mX; 1526 R[5] = 1. - y * y * inv1mX; 1527 R[6] = z; 1528 R[7] = 1. - z * z * inv1mX; 1529 R[8] = -y * z * inv1mX; 1530 } 1531 coords[0] = 0.0; 1532 coords[1] = r; 1533 coords[2] = 0.0; 1534 PetscFunctionReturn(PETSC_SUCCESS); 1535 } 1536 1537 /*@ 1538 DMPlexComputeProjection3Dto2D - Rewrite coordinates of 3 or more coplanar 3D points to a common 2D basis for the 1539 plane. The normal is defined by positive orientation of the first 3 points. 1540 1541 Not Collective 1542 1543 Input Parameter: 1544 . coordSize - Length of coordinate array (3x number of points); must be at least 9 (3 points) 1545 1546 Input/Output Parameter: 1547 . coords - The interlaced coordinates of each coplanar 3D point; on output the first 1548 2*coordSize/3 entries contain interlaced 2D points, with the rest undefined 1549 1550 Output Parameter: 1551 . 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. 1552 1553 Level: developer 1554 1555 .seealso: `DMPLEX`, `DMPlexComputeProjection2Dto1D()`, `DMPlexComputeProjection3Dto1D()` 1556 @*/ 1557 PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 1558 { 1559 PetscReal x1[3], x2[3], n[3], c[3], norm; 1560 const PetscInt dim = 3; 1561 PetscInt d, p; 1562 1563 PetscFunctionBegin; 1564 /* 0) Calculate normal vector */ 1565 for (d = 0; d < dim; ++d) { 1566 x1[d] = PetscRealPart(coords[1 * dim + d] - coords[0 * dim + d]); 1567 x2[d] = PetscRealPart(coords[2 * dim + d] - coords[0 * dim + d]); 1568 } 1569 // n = x1 \otimes x2 1570 n[0] = x1[1] * x2[2] - x1[2] * x2[1]; 1571 n[1] = x1[2] * x2[0] - x1[0] * x2[2]; 1572 n[2] = x1[0] * x2[1] - x1[1] * x2[0]; 1573 norm = PetscSqrtReal(n[0] * n[0] + n[1] * n[1] + n[2] * n[2]); 1574 for (d = 0; d < dim; d++) n[d] /= norm; 1575 norm = PetscSqrtReal(x1[0] * x1[0] + x1[1] * x1[1] + x1[2] * x1[2]); 1576 for (d = 0; d < dim; d++) x1[d] /= norm; 1577 // x2 = n \otimes x1 1578 x2[0] = n[1] * x1[2] - n[2] * x1[1]; 1579 x2[1] = n[2] * x1[0] - n[0] * x1[2]; 1580 x2[2] = n[0] * x1[1] - n[1] * x1[0]; 1581 for (d = 0; d < dim; d++) { 1582 R[d * dim + 0] = x1[d]; 1583 R[d * dim + 1] = x2[d]; 1584 R[d * dim + 2] = n[d]; 1585 c[d] = PetscRealPart(coords[0 * dim + d]); 1586 } 1587 for (p = 0; p < coordSize / dim; p++) { 1588 PetscReal y[3]; 1589 for (d = 0; d < dim; d++) y[d] = PetscRealPart(coords[p * dim + d]) - c[d]; 1590 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]; 1591 } 1592 PetscFunctionReturn(PETSC_SUCCESS); 1593 } 1594 1595 PETSC_UNUSED static inline void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 1596 { 1597 /* Signed volume is 1/2 the determinant 1598 1599 | 1 1 1 | 1600 | x0 x1 x2 | 1601 | y0 y1 y2 | 1602 1603 but if x0,y0 is the origin, we have 1604 1605 | x1 x2 | 1606 | y1 y2 | 1607 */ 1608 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 1609 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 1610 PetscReal M[4], detM; 1611 M[0] = x1; 1612 M[1] = x2; 1613 M[2] = y1; 1614 M[3] = y2; 1615 DMPlex_Det2D_Internal(&detM, M); 1616 *vol = 0.5 * detM; 1617 (void)PetscLogFlops(5.0); 1618 } 1619 1620 PETSC_UNUSED static inline void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 1621 { 1622 /* Signed volume is 1/6th of the determinant 1623 1624 | 1 1 1 1 | 1625 | x0 x1 x2 x3 | 1626 | y0 y1 y2 y3 | 1627 | z0 z1 z2 z3 | 1628 1629 but if x0,y0,z0 is the origin, we have 1630 1631 | x1 x2 x3 | 1632 | y1 y2 y3 | 1633 | z1 z2 z3 | 1634 */ 1635 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 1636 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 1637 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 1638 const PetscReal onesixth = ((PetscReal)1. / (PetscReal)6.); 1639 PetscReal M[9], detM; 1640 M[0] = x1; 1641 M[1] = x2; 1642 M[2] = x3; 1643 M[3] = y1; 1644 M[4] = y2; 1645 M[5] = y3; 1646 M[6] = z1; 1647 M[7] = z2; 1648 M[8] = z3; 1649 DMPlex_Det3D_Internal(&detM, M); 1650 *vol = -onesixth * detM; 1651 (void)PetscLogFlops(10.0); 1652 } 1653 1654 static inline void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 1655 { 1656 const PetscReal onesixth = ((PetscReal)1. / (PetscReal)6.); 1657 DMPlex_Det3D_Internal(vol, coords); 1658 *vol *= -onesixth; 1659 } 1660 1661 static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1662 { 1663 PetscSection coordSection; 1664 Vec coordinates; 1665 const PetscScalar *coords; 1666 PetscInt dim, d, off; 1667 1668 PetscFunctionBegin; 1669 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1670 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1671 PetscCall(PetscSectionGetDof(coordSection, e, &dim)); 1672 if (!dim) PetscFunctionReturn(PETSC_SUCCESS); 1673 PetscCall(PetscSectionGetOffset(coordSection, e, &off)); 1674 PetscCall(VecGetArrayRead(coordinates, &coords)); 1675 if (v0) { 1676 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]); 1677 } 1678 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 1679 *detJ = 1.; 1680 if (J) { 1681 for (d = 0; d < dim * dim; d++) J[d] = 0.; 1682 for (d = 0; d < dim; d++) J[d * dim + d] = 1.; 1683 if (invJ) { 1684 for (d = 0; d < dim * dim; d++) invJ[d] = 0.; 1685 for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.; 1686 } 1687 } 1688 PetscFunctionReturn(PETSC_SUCCESS); 1689 } 1690 1691 /*@C 1692 DMPlexGetCellCoordinates - Get coordinates for a cell, taking into account periodicity 1693 1694 Not Collective 1695 1696 Input Parameters: 1697 + dm - The `DMPLEX` 1698 - cell - The cell number 1699 1700 Output Parameters: 1701 + isDG - Using cellwise coordinates 1702 . Nc - The number of coordinates 1703 . array - The coordinate array 1704 - coords - The cell coordinates 1705 1706 Level: developer 1707 1708 .seealso: `DMPLEX`, `DMPlexRestoreCellCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCellCoordinatesLocal()` 1709 @*/ 1710 PetscErrorCode DMPlexGetCellCoordinates(DM dm, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *coords[]) 1711 { 1712 DM cdm; 1713 Vec coordinates; 1714 PetscSection cs; 1715 const PetscScalar *ccoords; 1716 PetscInt pStart, pEnd; 1717 1718 PetscFunctionBeginHot; 1719 *isDG = PETSC_FALSE; 1720 *Nc = 0; 1721 *array = NULL; 1722 *coords = NULL; 1723 /* Check for cellwise coordinates */ 1724 PetscCall(DMGetCellCoordinateSection(dm, &cs)); 1725 if (!cs) goto cg; 1726 /* Check that the cell exists in the cellwise section */ 1727 PetscCall(PetscSectionGetChart(cs, &pStart, &pEnd)); 1728 if (cell < pStart || cell >= pEnd) goto cg; 1729 /* Check for cellwise coordinates for this cell */ 1730 PetscCall(PetscSectionGetDof(cs, cell, Nc)); 1731 if (!*Nc) goto cg; 1732 /* Check for cellwise coordinates */ 1733 PetscCall(DMGetCellCoordinatesLocalNoncollective(dm, &coordinates)); 1734 if (!coordinates) goto cg; 1735 /* Get cellwise coordinates */ 1736 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 1737 PetscCall(VecGetArrayRead(coordinates, array)); 1738 PetscCall(DMPlexPointLocalRead(cdm, cell, *array, &ccoords)); 1739 PetscCall(DMGetWorkArray(cdm, *Nc, MPIU_SCALAR, coords)); 1740 PetscCall(PetscArraycpy(*coords, ccoords, *Nc)); 1741 PetscCall(VecRestoreArrayRead(coordinates, array)); 1742 *isDG = PETSC_TRUE; 1743 PetscFunctionReturn(PETSC_SUCCESS); 1744 cg: 1745 /* Use continuous coordinates */ 1746 PetscCall(DMGetCoordinateDM(dm, &cdm)); 1747 PetscCall(DMGetCoordinateSection(dm, &cs)); 1748 PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coordinates)); 1749 PetscCall(DMPlexVecGetOrientedClosure_Internal(cdm, cs, PETSC_FALSE, coordinates, cell, 0, Nc, coords)); 1750 PetscFunctionReturn(PETSC_SUCCESS); 1751 } 1752 1753 /*@C 1754 DMPlexRestoreCellCoordinates - Get coordinates for a cell, taking into account periodicity 1755 1756 Not Collective 1757 1758 Input Parameters: 1759 + dm - The `DMPLEX` 1760 - cell - The cell number 1761 1762 Output Parameters: 1763 + isDG - Using cellwise coordinates 1764 . Nc - The number of coordinates 1765 . array - The coordinate array 1766 - coords - The cell coordinates 1767 1768 Level: developer 1769 1770 .seealso: `DMPLEX`, `DMPlexGetCellCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCellCoordinatesLocal()` 1771 @*/ 1772 PetscErrorCode DMPlexRestoreCellCoordinates(DM dm, PetscInt cell, PetscBool *isDG, PetscInt *Nc, const PetscScalar *array[], PetscScalar *coords[]) 1773 { 1774 DM cdm; 1775 PetscSection cs; 1776 Vec coordinates; 1777 1778 PetscFunctionBeginHot; 1779 if (*isDG) { 1780 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 1781 PetscCall(DMRestoreWorkArray(cdm, *Nc, MPIU_SCALAR, coords)); 1782 } else { 1783 PetscCall(DMGetCoordinateDM(dm, &cdm)); 1784 PetscCall(DMGetCoordinateSection(dm, &cs)); 1785 PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coordinates)); 1786 PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, cell, Nc, coords)); 1787 } 1788 PetscFunctionReturn(PETSC_SUCCESS); 1789 } 1790 1791 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1792 { 1793 const PetscScalar *array; 1794 PetscScalar *coords = NULL; 1795 PetscInt numCoords, d; 1796 PetscBool isDG; 1797 1798 PetscFunctionBegin; 1799 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 1800 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 1801 *detJ = 0.0; 1802 if (numCoords == 6) { 1803 const PetscInt dim = 3; 1804 PetscReal R[9], J0; 1805 1806 if (v0) { 1807 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 1808 } 1809 PetscCall(DMPlexComputeProjection3Dto1D(coords, R)); 1810 if (J) { 1811 J0 = 0.5 * PetscRealPart(coords[1]); 1812 J[0] = R[0] * J0; 1813 J[1] = R[1]; 1814 J[2] = R[2]; 1815 J[3] = R[3] * J0; 1816 J[4] = R[4]; 1817 J[5] = R[5]; 1818 J[6] = R[6] * J0; 1819 J[7] = R[7]; 1820 J[8] = R[8]; 1821 DMPlex_Det3D_Internal(detJ, J); 1822 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 1823 } 1824 } else if (numCoords == 4) { 1825 const PetscInt dim = 2; 1826 PetscReal R[4], J0; 1827 1828 if (v0) { 1829 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 1830 } 1831 PetscCall(DMPlexComputeProjection2Dto1D(coords, R)); 1832 if (J) { 1833 J0 = 0.5 * PetscRealPart(coords[1]); 1834 J[0] = R[0] * J0; 1835 J[1] = R[1]; 1836 J[2] = R[2] * J0; 1837 J[3] = R[3]; 1838 DMPlex_Det2D_Internal(detJ, J); 1839 if (invJ) DMPlex_Invert2D_Internal(invJ, J, *detJ); 1840 } 1841 } else if (numCoords == 2) { 1842 const PetscInt dim = 1; 1843 1844 if (v0) { 1845 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 1846 } 1847 if (J) { 1848 J[0] = 0.5 * (PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 1849 *detJ = J[0]; 1850 PetscCall(PetscLogFlops(2.0)); 1851 if (invJ) { 1852 invJ[0] = 1.0 / J[0]; 1853 PetscCall(PetscLogFlops(1.0)); 1854 } 1855 } 1856 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for segment %" PetscInt_FMT " is %" PetscInt_FMT " != 2 or 4 or 6", e, numCoords); 1857 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 1858 PetscFunctionReturn(PETSC_SUCCESS); 1859 } 1860 1861 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1862 { 1863 const PetscScalar *array; 1864 PetscScalar *coords = NULL; 1865 PetscInt numCoords, d; 1866 PetscBool isDG; 1867 1868 PetscFunctionBegin; 1869 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 1870 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 1871 *detJ = 0.0; 1872 if (numCoords == 9) { 1873 const PetscInt dim = 3; 1874 PetscReal R[9], J0[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 1875 1876 if (v0) { 1877 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 1878 } 1879 PetscCall(DMPlexComputeProjection3Dto2D(numCoords, coords, R)); 1880 if (J) { 1881 const PetscInt pdim = 2; 1882 1883 for (d = 0; d < pdim; d++) { 1884 for (PetscInt f = 0; f < pdim; f++) J0[d * dim + f] = 0.5 * (PetscRealPart(coords[(f + 1) * pdim + d]) - PetscRealPart(coords[0 * pdim + d])); 1885 } 1886 PetscCall(PetscLogFlops(8.0)); 1887 DMPlex_Det3D_Internal(detJ, J0); 1888 for (d = 0; d < dim; d++) { 1889 for (PetscInt f = 0; f < dim; f++) { 1890 J[d * dim + f] = 0.0; 1891 for (PetscInt g = 0; g < dim; g++) J[d * dim + f] += R[d * dim + g] * J0[g * dim + f]; 1892 } 1893 } 1894 PetscCall(PetscLogFlops(18.0)); 1895 } 1896 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 1897 } else if (numCoords == 6) { 1898 const PetscInt dim = 2; 1899 1900 if (v0) { 1901 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 1902 } 1903 if (J) { 1904 for (d = 0; d < dim; d++) { 1905 for (PetscInt f = 0; f < dim; f++) J[d * dim + f] = 0.5 * (PetscRealPart(coords[(f + 1) * dim + d]) - PetscRealPart(coords[0 * dim + d])); 1906 } 1907 PetscCall(PetscLogFlops(8.0)); 1908 DMPlex_Det2D_Internal(detJ, J); 1909 } 1910 if (invJ) DMPlex_Invert2D_Internal(invJ, J, *detJ); 1911 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %" PetscInt_FMT " != 6 or 9", numCoords); 1912 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 1913 PetscFunctionReturn(PETSC_SUCCESS); 1914 } 1915 1916 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscBool isTensor, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1917 { 1918 const PetscScalar *array; 1919 PetscScalar *coords = NULL; 1920 PetscInt numCoords, d; 1921 PetscBool isDG; 1922 1923 PetscFunctionBegin; 1924 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 1925 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 1926 if (!Nq) { 1927 PetscInt vorder[4] = {0, 1, 2, 3}; 1928 1929 if (isTensor) { 1930 vorder[2] = 3; 1931 vorder[3] = 2; 1932 } 1933 *detJ = 0.0; 1934 if (numCoords == 12) { 1935 const PetscInt dim = 3; 1936 PetscReal R[9], J0[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 1937 1938 if (v) { 1939 for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]); 1940 } 1941 PetscCall(DMPlexComputeProjection3Dto2D(numCoords, coords, R)); 1942 if (J) { 1943 const PetscInt pdim = 2; 1944 1945 for (d = 0; d < pdim; d++) { 1946 J0[d * dim + 0] = 0.5 * (PetscRealPart(coords[vorder[1] * pdim + d]) - PetscRealPart(coords[vorder[0] * pdim + d])); 1947 J0[d * dim + 1] = 0.5 * (PetscRealPart(coords[vorder[2] * pdim + d]) - PetscRealPart(coords[vorder[1] * pdim + d])); 1948 } 1949 PetscCall(PetscLogFlops(8.0)); 1950 DMPlex_Det3D_Internal(detJ, J0); 1951 for (d = 0; d < dim; d++) { 1952 for (PetscInt f = 0; f < dim; f++) { 1953 J[d * dim + f] = 0.0; 1954 for (PetscInt g = 0; g < dim; g++) J[d * dim + f] += R[d * dim + g] * J0[g * dim + f]; 1955 } 1956 } 1957 PetscCall(PetscLogFlops(18.0)); 1958 } 1959 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 1960 } else if (numCoords == 8) { 1961 const PetscInt dim = 2; 1962 1963 if (v) { 1964 for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]); 1965 } 1966 if (J) { 1967 for (d = 0; d < dim; d++) { 1968 J[d * dim + 0] = 0.5 * (PetscRealPart(coords[vorder[1] * dim + d]) - PetscRealPart(coords[vorder[0] * dim + d])); 1969 J[d * dim + 1] = 0.5 * (PetscRealPart(coords[vorder[3] * dim + d]) - PetscRealPart(coords[vorder[0] * dim + d])); 1970 } 1971 PetscCall(PetscLogFlops(8.0)); 1972 DMPlex_Det2D_Internal(detJ, J); 1973 } 1974 if (invJ) DMPlex_Invert2D_Internal(invJ, J, *detJ); 1975 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %" PetscInt_FMT " != 8 or 12", numCoords); 1976 } else { 1977 const PetscInt Nv = 4; 1978 const PetscInt dimR = 2; 1979 PetscInt zToPlex[4] = {0, 1, 3, 2}; 1980 PetscReal zOrder[12]; 1981 PetscReal zCoeff[12]; 1982 PetscInt i, j, k, l, dim; 1983 1984 if (isTensor) { 1985 zToPlex[2] = 2; 1986 zToPlex[3] = 3; 1987 } 1988 if (numCoords == 12) { 1989 dim = 3; 1990 } else if (numCoords == 8) { 1991 dim = 2; 1992 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %" PetscInt_FMT " != 8 or 12", numCoords); 1993 for (i = 0; i < Nv; i++) { 1994 PetscInt zi = zToPlex[i]; 1995 1996 for (j = 0; j < dim; j++) zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 1997 } 1998 for (j = 0; j < dim; j++) { 1999 /* Nodal basis for evaluation at the vertices: (1 \mp xi) (1 \mp eta): 2000 \phi^0 = (1 - xi - eta + xi eta) --> 1 = 1/4 ( \phi^0 + \phi^1 + \phi^2 + \phi^3) 2001 \phi^1 = (1 + xi - eta - xi eta) --> xi = 1/4 (-\phi^0 + \phi^1 - \phi^2 + \phi^3) 2002 \phi^2 = (1 - xi + eta - xi eta) --> eta = 1/4 (-\phi^0 - \phi^1 + \phi^2 + \phi^3) 2003 \phi^3 = (1 + xi + eta + xi eta) --> xi eta = 1/4 ( \phi^0 - \phi^1 - \phi^2 + \phi^3) 2004 */ 2005 zCoeff[dim * 0 + j] = 0.25 * (zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 2006 zCoeff[dim * 1 + j] = 0.25 * (-zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 2007 zCoeff[dim * 2 + j] = 0.25 * (-zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 2008 zCoeff[dim * 3 + j] = 0.25 * (zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 2009 } 2010 for (i = 0; i < Nq; i++) { 2011 PetscReal xi = points[dimR * i], eta = points[dimR * i + 1]; 2012 2013 if (v) { 2014 PetscReal extPoint[4]; 2015 2016 extPoint[0] = 1.; 2017 extPoint[1] = xi; 2018 extPoint[2] = eta; 2019 extPoint[3] = xi * eta; 2020 for (j = 0; j < dim; j++) { 2021 PetscReal val = 0.; 2022 2023 for (k = 0; k < Nv; k++) val += extPoint[k] * zCoeff[dim * k + j]; 2024 v[i * dim + j] = val; 2025 } 2026 } 2027 if (J) { 2028 PetscReal extJ[8]; 2029 2030 extJ[0] = 0.; 2031 extJ[1] = 0.; 2032 extJ[2] = 1.; 2033 extJ[3] = 0.; 2034 extJ[4] = 0.; 2035 extJ[5] = 1.; 2036 extJ[6] = eta; 2037 extJ[7] = xi; 2038 for (j = 0; j < dim; j++) { 2039 for (k = 0; k < dimR; k++) { 2040 PetscReal val = 0.; 2041 2042 for (l = 0; l < Nv; l++) val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 2043 J[i * dim * dim + dim * j + k] = val; 2044 } 2045 } 2046 if (dim == 3) { /* put the cross product in the third component of the Jacobian */ 2047 PetscReal x, y, z; 2048 PetscReal *iJ = &J[i * dim * dim]; 2049 PetscReal norm; 2050 2051 x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0]; 2052 y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1]; 2053 z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0]; 2054 norm = PetscSqrtReal(x * x + y * y + z * z); 2055 iJ[2] = x / norm; 2056 iJ[5] = y / norm; 2057 iJ[8] = z / norm; 2058 DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 2059 if (invJ) DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]); 2060 } else { 2061 DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]); 2062 if (invJ) DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]); 2063 } 2064 } 2065 } 2066 } 2067 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2068 PetscFunctionReturn(PETSC_SUCCESS); 2069 } 2070 2071 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 2072 { 2073 const PetscScalar *array; 2074 PetscScalar *coords = NULL; 2075 const PetscInt dim = 3; 2076 PetscInt numCoords, d; 2077 PetscBool isDG; 2078 2079 PetscFunctionBegin; 2080 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2081 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 2082 *detJ = 0.0; 2083 if (v0) { 2084 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 2085 } 2086 if (J) { 2087 for (d = 0; d < dim; d++) { 2088 /* I orient with outward face normals */ 2089 J[d * dim + 0] = 0.5 * (PetscRealPart(coords[2 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2090 J[d * dim + 1] = 0.5 * (PetscRealPart(coords[1 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2091 J[d * dim + 2] = 0.5 * (PetscRealPart(coords[3 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2092 } 2093 PetscCall(PetscLogFlops(18.0)); 2094 DMPlex_Det3D_Internal(detJ, J); 2095 } 2096 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 2097 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2098 PetscFunctionReturn(PETSC_SUCCESS); 2099 } 2100 2101 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 2102 { 2103 const PetscScalar *array; 2104 PetscScalar *coords = NULL; 2105 const PetscInt dim = 3; 2106 PetscInt numCoords, d; 2107 PetscBool isDG; 2108 2109 PetscFunctionBegin; 2110 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2111 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 2112 if (!Nq) { 2113 *detJ = 0.0; 2114 if (v) { 2115 for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]); 2116 } 2117 if (J) { 2118 for (d = 0; d < dim; d++) { 2119 J[d * dim + 0] = 0.5 * (PetscRealPart(coords[3 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2120 J[d * dim + 1] = 0.5 * (PetscRealPart(coords[1 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2121 J[d * dim + 2] = 0.5 * (PetscRealPart(coords[4 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2122 } 2123 PetscCall(PetscLogFlops(18.0)); 2124 DMPlex_Det3D_Internal(detJ, J); 2125 } 2126 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 2127 } else { 2128 const PetscInt Nv = 8; 2129 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 2130 const PetscInt dim = 3; 2131 const PetscInt dimR = 3; 2132 PetscReal zOrder[24]; 2133 PetscReal zCoeff[24]; 2134 PetscInt i, j, k, l; 2135 2136 for (i = 0; i < Nv; i++) { 2137 PetscInt zi = zToPlex[i]; 2138 2139 for (j = 0; j < dim; j++) zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 2140 } 2141 for (j = 0; j < dim; j++) { 2142 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]); 2143 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]); 2144 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]); 2145 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]); 2146 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]); 2147 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]); 2148 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]); 2149 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]); 2150 } 2151 for (i = 0; i < Nq; i++) { 2152 PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2]; 2153 2154 if (v) { 2155 PetscReal extPoint[8]; 2156 2157 extPoint[0] = 1.; 2158 extPoint[1] = xi; 2159 extPoint[2] = eta; 2160 extPoint[3] = xi * eta; 2161 extPoint[4] = theta; 2162 extPoint[5] = theta * xi; 2163 extPoint[6] = theta * eta; 2164 extPoint[7] = theta * eta * xi; 2165 for (j = 0; j < dim; j++) { 2166 PetscReal val = 0.; 2167 2168 for (k = 0; k < Nv; k++) val += extPoint[k] * zCoeff[dim * k + j]; 2169 v[i * dim + j] = val; 2170 } 2171 } 2172 if (J) { 2173 PetscReal extJ[24]; 2174 2175 extJ[0] = 0.; 2176 extJ[1] = 0.; 2177 extJ[2] = 0.; 2178 extJ[3] = 1.; 2179 extJ[4] = 0.; 2180 extJ[5] = 0.; 2181 extJ[6] = 0.; 2182 extJ[7] = 1.; 2183 extJ[8] = 0.; 2184 extJ[9] = eta; 2185 extJ[10] = xi; 2186 extJ[11] = 0.; 2187 extJ[12] = 0.; 2188 extJ[13] = 0.; 2189 extJ[14] = 1.; 2190 extJ[15] = theta; 2191 extJ[16] = 0.; 2192 extJ[17] = xi; 2193 extJ[18] = 0.; 2194 extJ[19] = theta; 2195 extJ[20] = eta; 2196 extJ[21] = theta * eta; 2197 extJ[22] = theta * xi; 2198 extJ[23] = eta * xi; 2199 2200 for (j = 0; j < dim; j++) { 2201 for (k = 0; k < dimR; k++) { 2202 PetscReal val = 0.; 2203 2204 for (l = 0; l < Nv; l++) val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 2205 J[i * dim * dim + dim * j + k] = val; 2206 } 2207 } 2208 DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 2209 if (invJ) DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]); 2210 } 2211 } 2212 } 2213 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2214 PetscFunctionReturn(PETSC_SUCCESS); 2215 } 2216 2217 static PetscErrorCode DMPlexComputeTriangularPrismGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 2218 { 2219 const PetscScalar *array; 2220 PetscScalar *coords = NULL; 2221 const PetscInt dim = 3; 2222 PetscInt numCoords, d; 2223 PetscBool isDG; 2224 2225 PetscFunctionBegin; 2226 PetscCall(DMPlexGetCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2227 PetscCheck(!invJ || J, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 2228 if (!Nq) { 2229 /* Assume that the map to the reference is affine */ 2230 *detJ = 0.0; 2231 if (v) { 2232 for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]); 2233 } 2234 if (J) { 2235 for (d = 0; d < dim; d++) { 2236 J[d * dim + 0] = 0.5 * (PetscRealPart(coords[2 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2237 J[d * dim + 1] = 0.5 * (PetscRealPart(coords[1 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2238 J[d * dim + 2] = 0.5 * (PetscRealPart(coords[4 * dim + d]) - PetscRealPart(coords[0 * dim + d])); 2239 } 2240 PetscCall(PetscLogFlops(18.0)); 2241 DMPlex_Det3D_Internal(detJ, J); 2242 } 2243 if (invJ) DMPlex_Invert3D_Internal(invJ, J, *detJ); 2244 } else { 2245 const PetscInt dim = 3; 2246 const PetscInt dimR = 3; 2247 const PetscInt Nv = 6; 2248 PetscReal verts[18]; 2249 PetscReal coeff[18]; 2250 PetscInt i, j, k, l; 2251 2252 for (i = 0; i < Nv; ++i) 2253 for (j = 0; j < dim; ++j) verts[dim * i + j] = PetscRealPart(coords[dim * i + j]); 2254 for (j = 0; j < dim; ++j) { 2255 /* Check for triangle, 2256 phi^0 = -1/2 (xi + eta) chi^0 = delta(-1, -1) x(xi) = \sum_k x_k phi^k(xi) = \sum_k chi^k(x) phi^k(xi) 2257 phi^1 = 1/2 (1 + xi) chi^1 = delta( 1, -1) y(xi) = \sum_k y_k phi^k(xi) = \sum_k chi^k(y) phi^k(xi) 2258 phi^2 = 1/2 (1 + eta) chi^2 = delta(-1, 1) 2259 2260 phi^0 + phi^1 + phi^2 = 1 coef_1 = 1/2 ( chi^1 + chi^2) 2261 -phi^0 + phi^1 - phi^2 = xi coef_xi = 1/2 (-chi^0 + chi^1) 2262 -phi^0 - phi^1 + phi^2 = eta coef_eta = 1/2 (-chi^0 + chi^2) 2263 2264 < chi_0 chi_1 chi_2> A / 1 1 1 \ / phi_0 \ <chi> I <phi>^T so we need the inverse transpose 2265 | -1 1 -1 | | phi_1 | = 2266 \ -1 -1 1 / \ phi_2 / 2267 2268 Check phi^0: 1/2 (phi^0 chi^1 + phi^0 chi^2 + phi^0 chi^0 - phi^0 chi^1 + phi^0 chi^0 - phi^0 chi^2) = phi^0 chi^0 2269 */ 2270 /* Nodal basis for evaluation at the vertices: {-xi - eta, 1 + xi, 1 + eta} (1 \mp zeta): 2271 \phi^0 = 1/4 ( -xi - eta + xi zeta + eta zeta) --> / 1 1 1 1 1 1 \ 1 2272 \phi^1 = 1/4 (1 + eta - zeta - eta zeta) --> | -1 1 -1 -1 -1 1 | eta 2273 \phi^2 = 1/4 (1 + xi - zeta - xi zeta) --> | -1 -1 1 -1 1 -1 | xi 2274 \phi^3 = 1/4 ( -xi - eta - xi zeta - eta zeta) --> | -1 -1 -1 1 1 1 | zeta 2275 \phi^4 = 1/4 (1 + xi + zeta + xi zeta) --> | 1 1 -1 -1 1 -1 | xi zeta 2276 \phi^5 = 1/4 (1 + eta + zeta + eta zeta) --> \ 1 -1 1 -1 -1 1 / eta zeta 2277 1/4 / 0 1 1 0 1 1 \ 2278 | -1 1 0 -1 0 1 | 2279 | -1 0 1 -1 1 0 | 2280 | 0 -1 -1 0 1 1 | 2281 | 1 0 -1 -1 1 0 | 2282 \ 1 -1 0 -1 0 1 / 2283 */ 2284 coeff[dim * 0 + j] = (1. / 4.) * (verts[dim * 1 + j] + verts[dim * 2 + j] + verts[dim * 4 + j] + verts[dim * 5 + j]); 2285 coeff[dim * 1 + j] = (1. / 4.) * (-verts[dim * 0 + j] + verts[dim * 1 + j] - verts[dim * 3 + j] + verts[dim * 5 + j]); 2286 coeff[dim * 2 + j] = (1. / 4.) * (-verts[dim * 0 + j] + verts[dim * 2 + j] - verts[dim * 3 + j] + verts[dim * 4 + j]); 2287 coeff[dim * 3 + j] = (1. / 4.) * (-verts[dim * 1 + j] - verts[dim * 2 + j] + verts[dim * 4 + j] + verts[dim * 5 + j]); 2288 coeff[dim * 4 + j] = (1. / 4.) * (verts[dim * 0 + j] - verts[dim * 2 + j] - verts[dim * 3 + j] + verts[dim * 4 + j]); 2289 coeff[dim * 5 + j] = (1. / 4.) * (verts[dim * 0 + j] - verts[dim * 1 + j] - verts[dim * 3 + j] + verts[dim * 5 + j]); 2290 /* For reference prism: 2291 {0, 0, 0} 2292 {0, 1, 0} 2293 {1, 0, 0} 2294 {0, 0, 1} 2295 {0, 0, 0} 2296 {0, 0, 0} 2297 */ 2298 } 2299 for (i = 0; i < Nq; ++i) { 2300 const PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], zeta = points[dimR * i + 2]; 2301 2302 if (v) { 2303 PetscReal extPoint[6]; 2304 PetscInt c; 2305 2306 extPoint[0] = 1.; 2307 extPoint[1] = eta; 2308 extPoint[2] = xi; 2309 extPoint[3] = zeta; 2310 extPoint[4] = xi * zeta; 2311 extPoint[5] = eta * zeta; 2312 for (c = 0; c < dim; ++c) { 2313 PetscReal val = 0.; 2314 2315 for (k = 0; k < Nv; ++k) val += extPoint[k] * coeff[k * dim + c]; 2316 v[i * dim + c] = val; 2317 } 2318 } 2319 if (J) { 2320 PetscReal extJ[18]; 2321 2322 extJ[0] = 0.; 2323 extJ[1] = 0.; 2324 extJ[2] = 0.; 2325 extJ[3] = 0.; 2326 extJ[4] = 1.; 2327 extJ[5] = 0.; 2328 extJ[6] = 1.; 2329 extJ[7] = 0.; 2330 extJ[8] = 0.; 2331 extJ[9] = 0.; 2332 extJ[10] = 0.; 2333 extJ[11] = 1.; 2334 extJ[12] = zeta; 2335 extJ[13] = 0.; 2336 extJ[14] = xi; 2337 extJ[15] = 0.; 2338 extJ[16] = zeta; 2339 extJ[17] = eta; 2340 2341 for (j = 0; j < dim; j++) { 2342 for (k = 0; k < dimR; k++) { 2343 PetscReal val = 0.; 2344 2345 for (l = 0; l < Nv; l++) val += coeff[dim * l + j] * extJ[dimR * l + k]; 2346 J[i * dim * dim + dim * j + k] = val; 2347 } 2348 } 2349 DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 2350 if (invJ) DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]); 2351 } 2352 } 2353 } 2354 PetscCall(DMPlexRestoreCellCoordinates(dm, e, &isDG, &numCoords, &array, &coords)); 2355 PetscFunctionReturn(PETSC_SUCCESS); 2356 } 2357 2358 static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 2359 { 2360 DMPolytopeType ct; 2361 PetscInt depth, dim, coordDim, coneSize, i; 2362 PetscInt Nq = 0; 2363 const PetscReal *points = NULL; 2364 DMLabel depthLabel; 2365 PetscReal xi0[3] = {-1., -1., -1.}, v0[3], J0[9], detJ0; 2366 PetscBool isAffine = PETSC_TRUE; 2367 2368 PetscFunctionBegin; 2369 PetscCall(DMPlexGetDepth(dm, &depth)); 2370 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2371 PetscCall(DMPlexGetDepthLabel(dm, &depthLabel)); 2372 PetscCall(DMLabelGetValue(depthLabel, cell, &dim)); 2373 if (depth == 1 && dim == 1) PetscCall(DMGetDimension(dm, &dim)); 2374 PetscCall(DMGetCoordinateDim(dm, &coordDim)); 2375 PetscCheck(coordDim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %" PetscInt_FMT " > 3", coordDim); 2376 if (quad) PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &points, NULL)); 2377 PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2378 switch (ct) { 2379 case DM_POLYTOPE_POINT: 2380 PetscCall(DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ)); 2381 isAffine = PETSC_FALSE; 2382 break; 2383 case DM_POLYTOPE_SEGMENT: 2384 case DM_POLYTOPE_POINT_PRISM_TENSOR: 2385 if (Nq) PetscCall(DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0)); 2386 else PetscCall(DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ)); 2387 break; 2388 case DM_POLYTOPE_TRIANGLE: 2389 if (Nq) PetscCall(DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0)); 2390 else PetscCall(DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ)); 2391 break; 2392 case DM_POLYTOPE_QUADRILATERAL: 2393 PetscCall(DMPlexComputeRectangleGeometry_Internal(dm, cell, PETSC_FALSE, Nq, points, v, J, invJ, detJ)); 2394 isAffine = PETSC_FALSE; 2395 break; 2396 case DM_POLYTOPE_SEG_PRISM_TENSOR: 2397 PetscCall(DMPlexComputeRectangleGeometry_Internal(dm, cell, PETSC_TRUE, Nq, points, v, J, invJ, detJ)); 2398 isAffine = PETSC_FALSE; 2399 break; 2400 case DM_POLYTOPE_TETRAHEDRON: 2401 if (Nq) PetscCall(DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0)); 2402 else PetscCall(DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ)); 2403 break; 2404 case DM_POLYTOPE_HEXAHEDRON: 2405 PetscCall(DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ)); 2406 isAffine = PETSC_FALSE; 2407 break; 2408 case DM_POLYTOPE_TRI_PRISM: 2409 PetscCall(DMPlexComputeTriangularPrismGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ)); 2410 isAffine = PETSC_FALSE; 2411 break; 2412 default: 2413 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No element geometry for cell %" PetscInt_FMT " with type %s", cell, DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]); 2414 } 2415 if (isAffine && Nq) { 2416 if (v) { 2417 for (i = 0; i < Nq; i++) CoordinatesRefToReal(coordDim, dim, xi0, v0, J0, &points[dim * i], &v[coordDim * i]); 2418 } 2419 if (detJ) { 2420 for (i = 0; i < Nq; i++) detJ[i] = detJ0; 2421 } 2422 if (J) { 2423 PetscInt k; 2424 2425 for (i = 0, k = 0; i < Nq; i++) { 2426 PetscInt j; 2427 2428 for (j = 0; j < coordDim * coordDim; j++, k++) J[k] = J0[j]; 2429 } 2430 } 2431 if (invJ) { 2432 PetscInt k; 2433 switch (coordDim) { 2434 case 0: 2435 break; 2436 case 1: 2437 invJ[0] = 1. / J0[0]; 2438 break; 2439 case 2: 2440 DMPlex_Invert2D_Internal(invJ, J0, detJ0); 2441 break; 2442 case 3: 2443 DMPlex_Invert3D_Internal(invJ, J0, detJ0); 2444 break; 2445 } 2446 for (i = 1, k = coordDim * coordDim; i < Nq; i++) { 2447 PetscInt j; 2448 2449 for (j = 0; j < coordDim * coordDim; j++, k++) invJ[k] = invJ[j]; 2450 } 2451 } 2452 } 2453 PetscFunctionReturn(PETSC_SUCCESS); 2454 } 2455 2456 /*@C 2457 DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 2458 2459 Collective 2460 2461 Input Parameters: 2462 + dm - the `DMPLEX` 2463 - cell - the cell 2464 2465 Output Parameters: 2466 + v0 - the translation part of this affine transform, meaning the translation to the origin (not the first vertex of the reference cell) 2467 . J - the Jacobian of the transform from the reference element 2468 . invJ - the inverse of the Jacobian 2469 - detJ - the Jacobian determinant 2470 2471 Level: advanced 2472 2473 .seealso: `DMPLEX`, `DMPlexComputeCellGeometryFEM()`, `DMGetCoordinateSection()`, `DMGetCoordinates()` 2474 @*/ 2475 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 2476 { 2477 PetscFunctionBegin; 2478 PetscCall(DMPlexComputeCellGeometryFEM_Implicit(dm, cell, NULL, v0, J, invJ, detJ)); 2479 PetscFunctionReturn(PETSC_SUCCESS); 2480 } 2481 2482 static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 2483 { 2484 const PetscScalar *array; 2485 PetscScalar *coords = NULL; 2486 PetscInt numCoords; 2487 PetscBool isDG; 2488 PetscQuadrature feQuad; 2489 const PetscReal *quadPoints; 2490 PetscTabulation T; 2491 PetscInt dim, cdim, pdim, qdim, Nq, q; 2492 2493 PetscFunctionBegin; 2494 PetscCall(DMGetDimension(dm, &dim)); 2495 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2496 PetscCall(DMPlexGetCellCoordinates(dm, point, &isDG, &numCoords, &array, &coords)); 2497 if (!quad) { /* use the first point of the first functional of the dual space */ 2498 PetscDualSpace dsp; 2499 2500 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 2501 PetscCall(PetscDualSpaceGetFunctional(dsp, 0, &quad)); 2502 PetscCall(PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL)); 2503 Nq = 1; 2504 } else { 2505 PetscCall(PetscQuadratureGetData(quad, &qdim, NULL, &Nq, &quadPoints, NULL)); 2506 } 2507 PetscCall(PetscFEGetDimension(fe, &pdim)); 2508 PetscCall(PetscFEGetQuadrature(fe, &feQuad)); 2509 if (feQuad == quad) { 2510 PetscCall(PetscFEGetCellTabulation(fe, J ? 1 : 0, &T)); 2511 PetscCheck(numCoords == pdim * cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "There are %" PetscInt_FMT " coordinates for point %" PetscInt_FMT " != %" PetscInt_FMT "*%" PetscInt_FMT, numCoords, point, pdim, cdim); 2512 } else { 2513 PetscCall(PetscFECreateTabulation(fe, 1, Nq, quadPoints, J ? 1 : 0, &T)); 2514 } 2515 PetscCheck(qdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %" PetscInt_FMT " != quadrature dimension %" PetscInt_FMT, dim, qdim); 2516 { 2517 const PetscReal *basis = T->T[0]; 2518 const PetscReal *basisDer = T->T[1]; 2519 PetscReal detJt; 2520 2521 PetscAssert(Nq == T->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Np %" PetscInt_FMT " != %" PetscInt_FMT, Nq, T->Np); 2522 PetscAssert(pdim == T->Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nb %" PetscInt_FMT " != %" PetscInt_FMT, pdim, T->Nb); 2523 PetscAssert(dim == T->Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nc %" PetscInt_FMT " != %" PetscInt_FMT, dim, T->Nc); 2524 PetscAssert(cdim == T->cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "cdim %" PetscInt_FMT " != %" PetscInt_FMT, cdim, T->cdim); 2525 if (v) { 2526 PetscCall(PetscArrayzero(v, Nq * cdim)); 2527 for (q = 0; q < Nq; ++q) { 2528 PetscInt i, k; 2529 2530 for (k = 0; k < pdim; ++k) { 2531 const PetscInt vertex = k / cdim; 2532 for (i = 0; i < cdim; ++i) v[q * cdim + i] += basis[(q * pdim + k) * cdim + i] * PetscRealPart(coords[vertex * cdim + i]); 2533 } 2534 PetscCall(PetscLogFlops(2.0 * pdim * cdim)); 2535 } 2536 } 2537 if (J) { 2538 PetscCall(PetscArrayzero(J, Nq * cdim * cdim)); 2539 for (q = 0; q < Nq; ++q) { 2540 PetscInt i, j, k, c, r; 2541 2542 /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 2543 for (k = 0; k < pdim; ++k) { 2544 const PetscInt vertex = k / cdim; 2545 for (j = 0; j < dim; ++j) { 2546 for (i = 0; i < cdim; ++i) J[(q * cdim + i) * cdim + j] += basisDer[((q * pdim + k) * cdim + i) * dim + j] * PetscRealPart(coords[vertex * cdim + i]); 2547 } 2548 } 2549 PetscCall(PetscLogFlops(2.0 * pdim * dim * cdim)); 2550 if (cdim > dim) { 2551 for (c = dim; c < cdim; ++c) 2552 for (r = 0; r < cdim; ++r) J[r * cdim + c] = r == c ? 1.0 : 0.0; 2553 } 2554 if (!detJ && !invJ) continue; 2555 detJt = 0.; 2556 switch (cdim) { 2557 case 3: 2558 DMPlex_Det3D_Internal(&detJt, &J[q * cdim * dim]); 2559 if (invJ) DMPlex_Invert3D_Internal(&invJ[q * cdim * dim], &J[q * cdim * dim], detJt); 2560 break; 2561 case 2: 2562 DMPlex_Det2D_Internal(&detJt, &J[q * cdim * dim]); 2563 if (invJ) DMPlex_Invert2D_Internal(&invJ[q * cdim * dim], &J[q * cdim * dim], detJt); 2564 break; 2565 case 1: 2566 detJt = J[q * cdim * dim]; 2567 if (invJ) invJ[q * cdim * dim] = 1.0 / detJt; 2568 } 2569 if (detJ) detJ[q] = detJt; 2570 } 2571 } else PetscCheck(!detJ && !invJ, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ"); 2572 } 2573 if (feQuad != quad) PetscCall(PetscTabulationDestroy(&T)); 2574 PetscCall(DMPlexRestoreCellCoordinates(dm, point, &isDG, &numCoords, &array, &coords)); 2575 PetscFunctionReturn(PETSC_SUCCESS); 2576 } 2577 2578 /*@C 2579 DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 2580 2581 Collective 2582 2583 Input Parameters: 2584 + dm - the `DMPLEX` 2585 . cell - the cell 2586 - quad - the quadrature containing the points in the reference element where the geometry will be evaluated. If `quad` is `NULL`, geometry will be 2587 evaluated at the first vertex of the reference element 2588 2589 Output Parameters: 2590 + v - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element 2591 . J - the Jacobian of the transform from the reference element at each quadrature point 2592 . invJ - the inverse of the Jacobian at each quadrature point 2593 - detJ - the Jacobian determinant at each quadrature point 2594 2595 Level: advanced 2596 2597 .seealso: `DMPLEX`, `DMGetCoordinateSection()`, `DMGetCoordinates()` 2598 @*/ 2599 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 2600 { 2601 DM cdm; 2602 PetscFE fe = NULL; 2603 2604 PetscFunctionBegin; 2605 PetscAssertPointer(detJ, 7); 2606 PetscCall(DMGetCoordinateDM(dm, &cdm)); 2607 if (cdm) { 2608 PetscClassId id; 2609 PetscInt numFields; 2610 PetscDS prob; 2611 PetscObject disc; 2612 2613 PetscCall(DMGetNumFields(cdm, &numFields)); 2614 if (numFields) { 2615 PetscCall(DMGetDS(cdm, &prob)); 2616 PetscCall(PetscDSGetDiscretization(prob, 0, &disc)); 2617 PetscCall(PetscObjectGetClassId(disc, &id)); 2618 if (id == PETSCFE_CLASSID) fe = (PetscFE)disc; 2619 } 2620 } 2621 if (!fe) PetscCall(DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ)); 2622 else PetscCall(DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ)); 2623 PetscFunctionReturn(PETSC_SUCCESS); 2624 } 2625 2626 static PetscErrorCode DMPlexComputeGeometryFVM_0D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 2627 { 2628 PetscSection coordSection; 2629 Vec coordinates; 2630 const PetscScalar *coords = NULL; 2631 PetscInt d, dof, off; 2632 2633 PetscFunctionBegin; 2634 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2635 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2636 PetscCall(VecGetArrayRead(coordinates, &coords)); 2637 2638 /* for a point the centroid is just the coord */ 2639 if (centroid) { 2640 PetscCall(PetscSectionGetDof(coordSection, cell, &dof)); 2641 PetscCall(PetscSectionGetOffset(coordSection, cell, &off)); 2642 for (d = 0; d < dof; d++) centroid[d] = PetscRealPart(coords[off + d]); 2643 } 2644 if (normal) { 2645 const PetscInt *support, *cones; 2646 PetscInt supportSize; 2647 PetscReal norm, sign; 2648 2649 /* compute the norm based upon the support centroids */ 2650 PetscCall(DMPlexGetSupportSize(dm, cell, &supportSize)); 2651 PetscCall(DMPlexGetSupport(dm, cell, &support)); 2652 PetscCall(DMPlexComputeCellGeometryFVM(dm, support[0], NULL, normal, NULL)); 2653 2654 /* Take the normal from the centroid of the support to the vertex*/ 2655 PetscCall(PetscSectionGetDof(coordSection, cell, &dof)); 2656 PetscCall(PetscSectionGetOffset(coordSection, cell, &off)); 2657 for (d = 0; d < dof; d++) normal[d] -= PetscRealPart(coords[off + d]); 2658 2659 /* Determine the sign of the normal based upon its location in the support */ 2660 PetscCall(DMPlexGetCone(dm, support[0], &cones)); 2661 sign = cones[0] == cell ? 1.0 : -1.0; 2662 2663 norm = DMPlex_NormD_Internal(dim, normal); 2664 for (d = 0; d < dim; ++d) normal[d] /= (norm * sign); 2665 } 2666 if (vol) *vol = 1.0; 2667 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 2668 PetscFunctionReturn(PETSC_SUCCESS); 2669 } 2670 2671 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 2672 { 2673 const PetscScalar *array; 2674 PetscScalar *coords = NULL; 2675 PetscInt cdim, coordSize, d; 2676 PetscBool isDG; 2677 2678 PetscFunctionBegin; 2679 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2680 PetscCall(DMPlexGetCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2681 PetscCheck(coordSize == cdim * 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge has %" PetscInt_FMT " coordinates != %" PetscInt_FMT, coordSize, cdim * 2); 2682 if (centroid) { 2683 for (d = 0; d < cdim; ++d) centroid[d] = 0.5 * PetscRealPart(coords[d] + coords[cdim + d]); 2684 } 2685 if (normal) { 2686 PetscReal norm; 2687 2688 switch (cdim) { 2689 case 3: 2690 normal[2] = 0.; /* fall through */ 2691 case 2: 2692 normal[0] = -PetscRealPart(coords[1] - coords[cdim + 1]); 2693 normal[1] = PetscRealPart(coords[0] - coords[cdim + 0]); 2694 break; 2695 case 1: 2696 normal[0] = 1.0; 2697 break; 2698 default: 2699 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension %" PetscInt_FMT " not supported", cdim); 2700 } 2701 norm = DMPlex_NormD_Internal(cdim, normal); 2702 for (d = 0; d < cdim; ++d) normal[d] /= norm; 2703 } 2704 if (vol) { 2705 *vol = 0.0; 2706 for (d = 0; d < cdim; ++d) *vol += PetscSqr(PetscRealPart(coords[d] - coords[cdim + d])); 2707 *vol = PetscSqrtReal(*vol); 2708 } 2709 PetscCall(DMPlexRestoreCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2710 PetscFunctionReturn(PETSC_SUCCESS); 2711 } 2712 2713 /* Centroid_i = (\sum_n A_n Cn_i) / A */ 2714 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 2715 { 2716 DMPolytopeType ct; 2717 const PetscScalar *array; 2718 PetscScalar *coords = NULL; 2719 PetscInt coordSize; 2720 PetscBool isDG; 2721 PetscInt fv[4] = {0, 1, 2, 3}; 2722 PetscInt cdim, numCorners, p, d; 2723 2724 PetscFunctionBegin; 2725 /* Must check for hybrid cells because prisms have a different orientation scheme */ 2726 PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2727 switch (ct) { 2728 case DM_POLYTOPE_SEG_PRISM_TENSOR: 2729 fv[2] = 3; 2730 fv[3] = 2; 2731 break; 2732 default: 2733 break; 2734 } 2735 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2736 PetscCall(DMPlexGetConeSize(dm, cell, &numCorners)); 2737 PetscCall(DMPlexGetCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2738 { 2739 PetscReal c[3] = {0., 0., 0.}, n[3] = {0., 0., 0.}, origin[3] = {0., 0., 0.}, norm; 2740 2741 for (d = 0; d < cdim; d++) origin[d] = PetscRealPart(coords[d]); 2742 for (p = 0; p < numCorners - 2; ++p) { 2743 PetscReal e0[3] = {0., 0., 0.}, e1[3] = {0., 0., 0.}; 2744 for (d = 0; d < cdim; d++) { 2745 e0[d] = PetscRealPart(coords[cdim * fv[p + 1] + d]) - origin[d]; 2746 e1[d] = PetscRealPart(coords[cdim * fv[p + 2] + d]) - origin[d]; 2747 } 2748 const PetscReal dx = e0[1] * e1[2] - e0[2] * e1[1]; 2749 const PetscReal dy = e0[2] * e1[0] - e0[0] * e1[2]; 2750 const PetscReal dz = e0[0] * e1[1] - e0[1] * e1[0]; 2751 const PetscReal a = PetscSqrtReal(dx * dx + dy * dy + dz * dz); 2752 2753 n[0] += dx; 2754 n[1] += dy; 2755 n[2] += dz; 2756 for (d = 0; d < cdim; d++) c[d] += a * PetscRealPart(origin[d] + coords[cdim * fv[p + 1] + d] + coords[cdim * fv[p + 2] + d]) / 3.; 2757 } 2758 norm = PetscSqrtReal(n[0] * n[0] + n[1] * n[1] + n[2] * n[2]); 2759 // Allow zero volume cells 2760 if (norm != 0) { 2761 n[0] /= norm; 2762 n[1] /= norm; 2763 n[2] /= norm; 2764 c[0] /= norm; 2765 c[1] /= norm; 2766 c[2] /= norm; 2767 } 2768 if (vol) *vol = 0.5 * norm; 2769 if (centroid) 2770 for (d = 0; d < cdim; ++d) centroid[d] = c[d]; 2771 if (normal) 2772 for (d = 0; d < cdim; ++d) normal[d] = n[d]; 2773 } 2774 PetscCall(DMPlexRestoreCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2775 PetscFunctionReturn(PETSC_SUCCESS); 2776 } 2777 2778 /* Centroid_i = (\sum_n V_n Cn_i) / V */ 2779 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 2780 { 2781 DMPolytopeType ct; 2782 const PetscScalar *array; 2783 PetscScalar *coords = NULL; 2784 PetscInt coordSize; 2785 PetscBool isDG; 2786 PetscReal vsum = 0.0, vtmp, coordsTmp[3 * 3], origin[3]; 2787 const PetscInt order[16] = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15}; 2788 const PetscInt *cone, *faceSizes, *faces; 2789 const DMPolytopeType *faceTypes; 2790 PetscBool isHybrid = PETSC_FALSE; 2791 PetscInt numFaces, f, fOff = 0, p, d; 2792 2793 PetscFunctionBegin; 2794 PetscCheck(dim <= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No support for dim %" PetscInt_FMT " > 3", dim); 2795 /* Must check for hybrid cells because prisms have a different orientation scheme */ 2796 PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2797 switch (ct) { 2798 case DM_POLYTOPE_POINT_PRISM_TENSOR: 2799 case DM_POLYTOPE_SEG_PRISM_TENSOR: 2800 case DM_POLYTOPE_TRI_PRISM_TENSOR: 2801 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 2802 isHybrid = PETSC_TRUE; 2803 default: 2804 break; 2805 } 2806 2807 if (centroid) 2808 for (d = 0; d < dim; ++d) centroid[d] = 0.0; 2809 PetscCall(DMPlexGetCone(dm, cell, &cone)); 2810 2811 // Using the closure of faces for coordinates does not work in periodic geometries, so we index into the cell coordinates 2812 PetscCall(DMPlexGetRawFaces_Internal(dm, ct, order, &numFaces, &faceTypes, &faceSizes, &faces)); 2813 PetscCall(DMPlexGetCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2814 for (f = 0; f < numFaces; ++f) { 2815 PetscBool flip = isHybrid && f == 0 ? PETSC_TRUE : PETSC_FALSE; /* The first hybrid face is reversed */ 2816 2817 // If using zero as the origin vertex for each tetrahedron, an element far from the origin will have positive and 2818 // negative volumes that nearly cancel, thus incurring rounding error. Here we define origin[] as the first vertex 2819 // so that all tetrahedra have positive volume. 2820 if (f == 0) 2821 for (d = 0; d < dim; d++) origin[d] = PetscRealPart(coords[d]); 2822 switch (faceTypes[f]) { 2823 case DM_POLYTOPE_TRIANGLE: 2824 for (d = 0; d < dim; ++d) { 2825 coordsTmp[0 * dim + d] = PetscRealPart(coords[faces[fOff + 0] * dim + d]) - origin[d]; 2826 coordsTmp[1 * dim + d] = PetscRealPart(coords[faces[fOff + 1] * dim + d]) - origin[d]; 2827 coordsTmp[2 * dim + d] = PetscRealPart(coords[faces[fOff + 2] * dim + d]) - origin[d]; 2828 } 2829 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 2830 if (flip) vtmp = -vtmp; 2831 vsum += vtmp; 2832 if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 2833 for (d = 0; d < dim; ++d) { 2834 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p * dim + d] * vtmp; 2835 } 2836 } 2837 break; 2838 case DM_POLYTOPE_QUADRILATERAL: 2839 case DM_POLYTOPE_SEG_PRISM_TENSOR: { 2840 PetscInt fv[4] = {0, 1, 2, 3}; 2841 2842 /* Side faces for hybrid cells are stored as tensor products */ 2843 if (isHybrid && f > 1) { 2844 fv[2] = 3; 2845 fv[3] = 2; 2846 } 2847 /* DO FOR PYRAMID */ 2848 /* First tet */ 2849 for (d = 0; d < dim; ++d) { 2850 coordsTmp[0 * dim + d] = PetscRealPart(coords[faces[fOff + fv[0]] * dim + d]) - origin[d]; 2851 coordsTmp[1 * dim + d] = PetscRealPart(coords[faces[fOff + fv[1]] * dim + d]) - origin[d]; 2852 coordsTmp[2 * dim + d] = PetscRealPart(coords[faces[fOff + fv[3]] * dim + d]) - origin[d]; 2853 } 2854 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 2855 if (flip) vtmp = -vtmp; 2856 vsum += vtmp; 2857 if (centroid) { 2858 for (d = 0; d < dim; ++d) { 2859 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p * dim + d] * vtmp; 2860 } 2861 } 2862 /* Second tet */ 2863 for (d = 0; d < dim; ++d) { 2864 coordsTmp[0 * dim + d] = PetscRealPart(coords[faces[fOff + fv[1]] * dim + d]) - origin[d]; 2865 coordsTmp[1 * dim + d] = PetscRealPart(coords[faces[fOff + fv[2]] * dim + d]) - origin[d]; 2866 coordsTmp[2 * dim + d] = PetscRealPart(coords[faces[fOff + fv[3]] * dim + d]) - origin[d]; 2867 } 2868 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 2869 if (flip) vtmp = -vtmp; 2870 vsum += vtmp; 2871 if (centroid) { 2872 for (d = 0; d < dim; ++d) { 2873 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p * dim + d] * vtmp; 2874 } 2875 } 2876 break; 2877 } 2878 default: 2879 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle face %" PetscInt_FMT " of type %s", cone[f], DMPolytopeTypes[ct]); 2880 } 2881 fOff += faceSizes[f]; 2882 } 2883 PetscCall(DMPlexRestoreRawFaces_Internal(dm, ct, order, &numFaces, &faceTypes, &faceSizes, &faces)); 2884 PetscCall(DMPlexRestoreCellCoordinates(dm, cell, &isDG, &coordSize, &array, &coords)); 2885 if (vol) *vol = PetscAbsReal(vsum); 2886 if (normal) 2887 for (d = 0; d < dim; ++d) normal[d] = 0.0; 2888 if (centroid) 2889 for (d = 0; d < dim; ++d) centroid[d] = centroid[d] / (vsum * 4) + origin[d]; 2890 PetscFunctionReturn(PETSC_SUCCESS); 2891 } 2892 2893 /*@C 2894 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 2895 2896 Collective 2897 2898 Input Parameters: 2899 + dm - the `DMPLEX` 2900 - cell - the cell 2901 2902 Output Parameters: 2903 + vol - the cell volume 2904 . centroid - the cell centroid 2905 - normal - the cell normal, if appropriate 2906 2907 Level: advanced 2908 2909 .seealso: `DMPLEX`, `DMGetCoordinateSection()`, `DMGetCoordinates()` 2910 @*/ 2911 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 2912 { 2913 PetscInt depth, dim; 2914 2915 PetscFunctionBegin; 2916 PetscCall(DMPlexGetDepth(dm, &depth)); 2917 PetscCall(DMGetDimension(dm, &dim)); 2918 PetscCheck(depth == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 2919 PetscCall(DMPlexGetPointDepth(dm, cell, &depth)); 2920 switch (depth) { 2921 case 0: 2922 PetscCall(DMPlexComputeGeometryFVM_0D_Internal(dm, dim, cell, vol, centroid, normal)); 2923 break; 2924 case 1: 2925 PetscCall(DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal)); 2926 break; 2927 case 2: 2928 PetscCall(DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal)); 2929 break; 2930 case 3: 2931 PetscCall(DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal)); 2932 break; 2933 default: 2934 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT " (depth %" PetscInt_FMT ") for element geometry computation", dim, depth); 2935 } 2936 PetscFunctionReturn(PETSC_SUCCESS); 2937 } 2938 2939 /*@ 2940 DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 2941 2942 Input Parameter: 2943 . dm - The `DMPLEX` 2944 2945 Output Parameters: 2946 + cellgeom - A `Vec` of `PetscFVCellGeom` data 2947 - facegeom - A `Vec` of `PetscFVFaceGeom` data 2948 2949 Level: developer 2950 2951 .seealso: `DMPLEX`, `PetscFVFaceGeom`, `PetscFVCellGeom` 2952 @*/ 2953 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 2954 { 2955 DM dmFace, dmCell; 2956 DMLabel ghostLabel; 2957 PetscSection sectionFace, sectionCell; 2958 PetscSection coordSection; 2959 Vec coordinates; 2960 PetscScalar *fgeom, *cgeom; 2961 PetscReal minradius, gminradius; 2962 PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 2963 2964 PetscFunctionBegin; 2965 PetscCall(DMGetDimension(dm, &dim)); 2966 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2967 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2968 /* Make cell centroids and volumes */ 2969 PetscCall(DMClone(dm, &dmCell)); 2970 PetscCall(DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection)); 2971 PetscCall(DMSetCoordinatesLocal(dmCell, coordinates)); 2972 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionCell)); 2973 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2974 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL)); 2975 PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd)); 2976 for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, (PetscInt)PetscCeilReal(((PetscReal)sizeof(PetscFVCellGeom)) / sizeof(PetscScalar)))); 2977 PetscCall(PetscSectionSetUp(sectionCell)); 2978 PetscCall(DMSetLocalSection(dmCell, sectionCell)); 2979 PetscCall(PetscSectionDestroy(§ionCell)); 2980 PetscCall(DMCreateLocalVector(dmCell, cellgeom)); 2981 if (cEndInterior < 0) cEndInterior = cEnd; 2982 PetscCall(VecGetArray(*cellgeom, &cgeom)); 2983 for (c = cStart; c < cEndInterior; ++c) { 2984 PetscFVCellGeom *cg; 2985 2986 PetscCall(DMPlexPointLocalRef(dmCell, c, cgeom, &cg)); 2987 PetscCall(PetscArrayzero(cg, 1)); 2988 PetscCall(DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL)); 2989 } 2990 /* Compute face normals and minimum cell radius */ 2991 PetscCall(DMClone(dm, &dmFace)); 2992 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionFace)); 2993 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 2994 PetscCall(PetscSectionSetChart(sectionFace, fStart, fEnd)); 2995 for (f = fStart; f < fEnd; ++f) PetscCall(PetscSectionSetDof(sectionFace, f, (PetscInt)PetscCeilReal(((PetscReal)sizeof(PetscFVFaceGeom)) / sizeof(PetscScalar)))); 2996 PetscCall(PetscSectionSetUp(sectionFace)); 2997 PetscCall(DMSetLocalSection(dmFace, sectionFace)); 2998 PetscCall(PetscSectionDestroy(§ionFace)); 2999 PetscCall(DMCreateLocalVector(dmFace, facegeom)); 3000 PetscCall(VecGetArray(*facegeom, &fgeom)); 3001 PetscCall(DMGetLabel(dm, "ghost", &ghostLabel)); 3002 minradius = PETSC_MAX_REAL; 3003 for (f = fStart; f < fEnd; ++f) { 3004 PetscFVFaceGeom *fg; 3005 PetscReal area; 3006 const PetscInt *cells; 3007 PetscInt ncells, ghost = -1, d, numChildren; 3008 3009 if (ghostLabel) PetscCall(DMLabelGetValue(ghostLabel, f, &ghost)); 3010 PetscCall(DMPlexGetTreeChildren(dm, f, &numChildren, NULL)); 3011 PetscCall(DMPlexGetSupport(dm, f, &cells)); 3012 PetscCall(DMPlexGetSupportSize(dm, f, &ncells)); 3013 /* It is possible to get a face with no support when using partition overlap */ 3014 if (!ncells || ghost >= 0 || numChildren) continue; 3015 PetscCall(DMPlexPointLocalRef(dmFace, f, fgeom, &fg)); 3016 PetscCall(DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal)); 3017 for (d = 0; d < dim; ++d) fg->normal[d] *= area; 3018 /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 3019 { 3020 PetscFVCellGeom *cL, *cR; 3021 PetscReal *lcentroid, *rcentroid; 3022 PetscReal l[3], r[3], v[3]; 3023 3024 PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL)); 3025 lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 3026 if (ncells > 1) { 3027 PetscCall(DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR)); 3028 rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 3029 } else { 3030 rcentroid = fg->centroid; 3031 } 3032 PetscCall(DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l)); 3033 PetscCall(DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r)); 3034 DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 3035 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 3036 for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 3037 } 3038 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 3039 PetscCheck(dim != 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " 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]); 3040 PetscCheck(dim != 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " 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]); 3041 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Direction for face %" PetscInt_FMT " could not be fixed", f); 3042 } 3043 if (cells[0] < cEndInterior) { 3044 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 3045 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 3046 } 3047 if (ncells > 1 && cells[1] < cEndInterior) { 3048 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 3049 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 3050 } 3051 } 3052 } 3053 PetscCallMPI(MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm))); 3054 PetscCall(DMPlexSetMinRadius(dm, gminradius)); 3055 /* Compute centroids of ghost cells */ 3056 for (c = cEndInterior; c < cEnd; ++c) { 3057 PetscFVFaceGeom *fg; 3058 const PetscInt *cone, *support; 3059 PetscInt coneSize, supportSize, s; 3060 3061 PetscCall(DMPlexGetConeSize(dmCell, c, &coneSize)); 3062 PetscCheck(coneSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %" PetscInt_FMT " has cone size %" PetscInt_FMT " != 1", c, coneSize); 3063 PetscCall(DMPlexGetCone(dmCell, c, &cone)); 3064 PetscCall(DMPlexGetSupportSize(dmCell, cone[0], &supportSize)); 3065 PetscCheck(supportSize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[0], supportSize); 3066 PetscCall(DMPlexGetSupport(dmCell, cone[0], &support)); 3067 PetscCall(DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg)); 3068 for (s = 0; s < 2; ++s) { 3069 /* Reflect ghost centroid across plane of face */ 3070 if (support[s] == c) { 3071 PetscFVCellGeom *ci; 3072 PetscFVCellGeom *cg; 3073 PetscReal c2f[3], a; 3074 3075 PetscCall(DMPlexPointLocalRead(dmCell, support[(s + 1) % 2], cgeom, &ci)); 3076 DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 3077 a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal) / DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 3078 PetscCall(DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg)); 3079 DMPlex_WaxpyD_Internal(dim, 2 * a, fg->normal, ci->centroid, cg->centroid); 3080 cg->volume = ci->volume; 3081 } 3082 } 3083 } 3084 PetscCall(VecRestoreArray(*facegeom, &fgeom)); 3085 PetscCall(VecRestoreArray(*cellgeom, &cgeom)); 3086 PetscCall(DMDestroy(&dmCell)); 3087 PetscCall(DMDestroy(&dmFace)); 3088 PetscFunctionReturn(PETSC_SUCCESS); 3089 } 3090 3091 /*@ 3092 DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 3093 3094 Not Collective 3095 3096 Input Parameter: 3097 . dm - the `DMPLEX` 3098 3099 Output Parameter: 3100 . minradius - the minimum cell radius 3101 3102 Level: developer 3103 3104 .seealso: `DMPLEX`, `DMGetCoordinates()` 3105 @*/ 3106 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 3107 { 3108 PetscFunctionBegin; 3109 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3110 PetscAssertPointer(minradius, 2); 3111 *minradius = ((DM_Plex *)dm->data)->minradius; 3112 PetscFunctionReturn(PETSC_SUCCESS); 3113 } 3114 3115 /*@ 3116 DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 3117 3118 Logically Collective 3119 3120 Input Parameters: 3121 + dm - the `DMPLEX` 3122 - minradius - the minimum cell radius 3123 3124 Level: developer 3125 3126 .seealso: `DMPLEX`, `DMSetCoordinates()` 3127 @*/ 3128 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 3129 { 3130 PetscFunctionBegin; 3131 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3132 ((DM_Plex *)dm->data)->minradius = minradius; 3133 PetscFunctionReturn(PETSC_SUCCESS); 3134 } 3135 3136 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 3137 { 3138 DMLabel ghostLabel; 3139 PetscScalar *dx, *grad, **gref; 3140 PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 3141 3142 PetscFunctionBegin; 3143 PetscCall(DMGetDimension(dm, &dim)); 3144 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3145 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL)); 3146 cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior; 3147 PetscCall(DMPlexGetMaxSizes(dm, &maxNumFaces, NULL)); 3148 PetscCall(PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces)); 3149 PetscCall(DMGetLabel(dm, "ghost", &ghostLabel)); 3150 PetscCall(PetscMalloc3(maxNumFaces * dim, &dx, maxNumFaces * dim, &grad, maxNumFaces, &gref)); 3151 for (c = cStart; c < cEndInterior; c++) { 3152 const PetscInt *faces; 3153 PetscInt numFaces, usedFaces, f, d; 3154 PetscFVCellGeom *cg; 3155 PetscBool boundary; 3156 PetscInt ghost; 3157 3158 // do not attempt to compute a gradient reconstruction stencil in a ghost cell. It will never be used 3159 PetscCall(DMLabelGetValue(ghostLabel, c, &ghost)); 3160 if (ghost >= 0) continue; 3161 3162 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 3163 PetscCall(DMPlexGetConeSize(dm, c, &numFaces)); 3164 PetscCall(DMPlexGetCone(dm, c, &faces)); 3165 PetscCheck(numFaces >= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cell %" PetscInt_FMT " has only %" PetscInt_FMT " faces, not enough for gradient reconstruction", c, numFaces); 3166 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 3167 PetscFVCellGeom *cg1; 3168 PetscFVFaceGeom *fg; 3169 const PetscInt *fcells; 3170 PetscInt ncell, side; 3171 3172 PetscCall(DMLabelGetValue(ghostLabel, faces[f], &ghost)); 3173 PetscCall(DMIsBoundaryPoint(dm, faces[f], &boundary)); 3174 if ((ghost >= 0) || boundary) continue; 3175 PetscCall(DMPlexGetSupport(dm, faces[f], &fcells)); 3176 side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 3177 ncell = fcells[!side]; /* the neighbor */ 3178 PetscCall(DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg)); 3179 PetscCall(DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1)); 3180 for (d = 0; d < dim; ++d) dx[usedFaces * dim + d] = cg1->centroid[d] - cg->centroid[d]; 3181 gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 3182 } 3183 PetscCheck(usedFaces, PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 3184 PetscCall(PetscFVComputeGradient(fvm, usedFaces, dx, grad)); 3185 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 3186 PetscCall(DMLabelGetValue(ghostLabel, faces[f], &ghost)); 3187 PetscCall(DMIsBoundaryPoint(dm, faces[f], &boundary)); 3188 if ((ghost >= 0) || boundary) continue; 3189 for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces * dim + d]; 3190 ++usedFaces; 3191 } 3192 } 3193 PetscCall(PetscFree3(dx, grad, gref)); 3194 PetscFunctionReturn(PETSC_SUCCESS); 3195 } 3196 3197 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 3198 { 3199 DMLabel ghostLabel; 3200 PetscScalar *dx, *grad, **gref; 3201 PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 3202 PetscSection neighSec; 3203 PetscInt(*neighbors)[2]; 3204 PetscInt *counter; 3205 3206 PetscFunctionBegin; 3207 PetscCall(DMGetDimension(dm, &dim)); 3208 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3209 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL)); 3210 if (cEndInterior < 0) cEndInterior = cEnd; 3211 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &neighSec)); 3212 PetscCall(PetscSectionSetChart(neighSec, cStart, cEndInterior)); 3213 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 3214 PetscCall(DMGetLabel(dm, "ghost", &ghostLabel)); 3215 for (f = fStart; f < fEnd; f++) { 3216 const PetscInt *fcells; 3217 PetscBool boundary; 3218 PetscInt ghost = -1; 3219 PetscInt numChildren, numCells, c; 3220 3221 if (ghostLabel) PetscCall(DMLabelGetValue(ghostLabel, f, &ghost)); 3222 PetscCall(DMIsBoundaryPoint(dm, f, &boundary)); 3223 PetscCall(DMPlexGetTreeChildren(dm, f, &numChildren, NULL)); 3224 if ((ghost >= 0) || boundary || numChildren) continue; 3225 PetscCall(DMPlexGetSupportSize(dm, f, &numCells)); 3226 if (numCells == 2) { 3227 PetscCall(DMPlexGetSupport(dm, f, &fcells)); 3228 for (c = 0; c < 2; c++) { 3229 PetscInt cell = fcells[c]; 3230 3231 if (cell >= cStart && cell < cEndInterior) PetscCall(PetscSectionAddDof(neighSec, cell, 1)); 3232 } 3233 } 3234 } 3235 PetscCall(PetscSectionSetUp(neighSec)); 3236 PetscCall(PetscSectionGetMaxDof(neighSec, &maxNumFaces)); 3237 PetscCall(PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces)); 3238 nStart = 0; 3239 PetscCall(PetscSectionGetStorageSize(neighSec, &nEnd)); 3240 PetscCall(PetscMalloc1(nEnd - nStart, &neighbors)); 3241 PetscCall(PetscCalloc1(cEndInterior - cStart, &counter)); 3242 for (f = fStart; f < fEnd; f++) { 3243 const PetscInt *fcells; 3244 PetscBool boundary; 3245 PetscInt ghost = -1; 3246 PetscInt numChildren, numCells, c; 3247 3248 if (ghostLabel) PetscCall(DMLabelGetValue(ghostLabel, f, &ghost)); 3249 PetscCall(DMIsBoundaryPoint(dm, f, &boundary)); 3250 PetscCall(DMPlexGetTreeChildren(dm, f, &numChildren, NULL)); 3251 if ((ghost >= 0) || boundary || numChildren) continue; 3252 PetscCall(DMPlexGetSupportSize(dm, f, &numCells)); 3253 if (numCells == 2) { 3254 PetscCall(DMPlexGetSupport(dm, f, &fcells)); 3255 for (c = 0; c < 2; c++) { 3256 PetscInt cell = fcells[c], off; 3257 3258 if (cell >= cStart && cell < cEndInterior) { 3259 PetscCall(PetscSectionGetOffset(neighSec, cell, &off)); 3260 off += counter[cell - cStart]++; 3261 neighbors[off][0] = f; 3262 neighbors[off][1] = fcells[1 - c]; 3263 } 3264 } 3265 } 3266 } 3267 PetscCall(PetscFree(counter)); 3268 PetscCall(PetscMalloc3(maxNumFaces * dim, &dx, maxNumFaces * dim, &grad, maxNumFaces, &gref)); 3269 for (c = cStart; c < cEndInterior; c++) { 3270 PetscInt numFaces, f, d, off, ghost = -1; 3271 PetscFVCellGeom *cg; 3272 3273 PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 3274 PetscCall(PetscSectionGetDof(neighSec, c, &numFaces)); 3275 PetscCall(PetscSectionGetOffset(neighSec, c, &off)); 3276 3277 // do not attempt to compute a gradient reconstruction stencil in a ghost cell. It will never be used 3278 if (ghostLabel) PetscCall(DMLabelGetValue(ghostLabel, c, &ghost)); 3279 if (ghost >= 0) continue; 3280 3281 PetscCheck(numFaces >= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cell %" PetscInt_FMT " has only %" PetscInt_FMT " faces, not enough for gradient reconstruction", c, numFaces); 3282 for (f = 0; f < numFaces; ++f) { 3283 PetscFVCellGeom *cg1; 3284 PetscFVFaceGeom *fg; 3285 const PetscInt *fcells; 3286 PetscInt ncell, side, nface; 3287 3288 nface = neighbors[off + f][0]; 3289 ncell = neighbors[off + f][1]; 3290 PetscCall(DMPlexGetSupport(dm, nface, &fcells)); 3291 side = (c != fcells[0]); 3292 PetscCall(DMPlexPointLocalRef(dmFace, nface, fgeom, &fg)); 3293 PetscCall(DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1)); 3294 for (d = 0; d < dim; ++d) dx[f * dim + d] = cg1->centroid[d] - cg->centroid[d]; 3295 gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 3296 } 3297 PetscCall(PetscFVComputeGradient(fvm, numFaces, dx, grad)); 3298 for (f = 0; f < numFaces; ++f) { 3299 for (d = 0; d < dim; ++d) gref[f][d] = grad[f * dim + d]; 3300 } 3301 } 3302 PetscCall(PetscFree3(dx, grad, gref)); 3303 PetscCall(PetscSectionDestroy(&neighSec)); 3304 PetscCall(PetscFree(neighbors)); 3305 PetscFunctionReturn(PETSC_SUCCESS); 3306 } 3307 3308 /*@ 3309 DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 3310 3311 Collective 3312 3313 Input Parameters: 3314 + dm - The `DMPLEX` 3315 . fvm - The `PetscFV` 3316 - cellGeometry - The face geometry from `DMPlexComputeCellGeometryFVM()` 3317 3318 Input/Output Parameter: 3319 . faceGeometry - The face geometry from `DMPlexComputeFaceGeometryFVM()`; on output 3320 the geometric factors for gradient calculation are inserted 3321 3322 Output Parameter: 3323 . dmGrad - The `DM` describing the layout of gradient data 3324 3325 Level: developer 3326 3327 .seealso: `DMPLEX`, `DMPlexGetFaceGeometryFVM()`, `DMPlexGetCellGeometryFVM()` 3328 @*/ 3329 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 3330 { 3331 DM dmFace, dmCell; 3332 PetscScalar *fgeom, *cgeom; 3333 PetscSection sectionGrad, parentSection; 3334 PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 3335 3336 PetscFunctionBegin; 3337 PetscCall(DMGetDimension(dm, &dim)); 3338 PetscCall(PetscFVGetNumComponents(fvm, &pdim)); 3339 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3340 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL)); 3341 /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 3342 PetscCall(VecGetDM(faceGeometry, &dmFace)); 3343 PetscCall(VecGetDM(cellGeometry, &dmCell)); 3344 PetscCall(VecGetArray(faceGeometry, &fgeom)); 3345 PetscCall(VecGetArray(cellGeometry, &cgeom)); 3346 PetscCall(DMPlexGetTree(dm, &parentSection, NULL, NULL, NULL, NULL)); 3347 if (!parentSection) { 3348 PetscCall(BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom)); 3349 } else { 3350 PetscCall(BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom)); 3351 } 3352 PetscCall(VecRestoreArray(faceGeometry, &fgeom)); 3353 PetscCall(VecRestoreArray(cellGeometry, &cgeom)); 3354 /* Create storage for gradients */ 3355 PetscCall(DMClone(dm, dmGrad)); 3356 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionGrad)); 3357 PetscCall(PetscSectionSetChart(sectionGrad, cStart, cEnd)); 3358 for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionGrad, c, pdim * dim)); 3359 PetscCall(PetscSectionSetUp(sectionGrad)); 3360 PetscCall(DMSetLocalSection(*dmGrad, sectionGrad)); 3361 PetscCall(PetscSectionDestroy(§ionGrad)); 3362 PetscFunctionReturn(PETSC_SUCCESS); 3363 } 3364 3365 /*@ 3366 DMPlexGetDataFVM - Retrieve precomputed cell geometry 3367 3368 Collective 3369 3370 Input Parameters: 3371 + dm - The `DM` 3372 - fv - The `PetscFV` 3373 3374 Output Parameters: 3375 + cellgeom - The cell geometry 3376 . facegeom - The face geometry 3377 - gradDM - The gradient matrices 3378 3379 Level: developer 3380 3381 .seealso: `DMPLEX`, `DMPlexComputeGeometryFVM()` 3382 @*/ 3383 PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM) 3384 { 3385 PetscObject cellgeomobj, facegeomobj; 3386 3387 PetscFunctionBegin; 3388 PetscCall(PetscObjectQuery((PetscObject)dm, "DMPlex_cellgeom_fvm", &cellgeomobj)); 3389 if (!cellgeomobj) { 3390 Vec cellgeomInt, facegeomInt; 3391 3392 PetscCall(DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt)); 3393 PetscCall(PetscObjectCompose((PetscObject)dm, "DMPlex_cellgeom_fvm", (PetscObject)cellgeomInt)); 3394 PetscCall(PetscObjectCompose((PetscObject)dm, "DMPlex_facegeom_fvm", (PetscObject)facegeomInt)); 3395 PetscCall(VecDestroy(&cellgeomInt)); 3396 PetscCall(VecDestroy(&facegeomInt)); 3397 PetscCall(PetscObjectQuery((PetscObject)dm, "DMPlex_cellgeom_fvm", &cellgeomobj)); 3398 } 3399 PetscCall(PetscObjectQuery((PetscObject)dm, "DMPlex_facegeom_fvm", &facegeomobj)); 3400 if (cellgeom) *cellgeom = (Vec)cellgeomobj; 3401 if (facegeom) *facegeom = (Vec)facegeomobj; 3402 if (gradDM) { 3403 PetscObject gradobj; 3404 PetscBool computeGradients; 3405 3406 PetscCall(PetscFVGetComputeGradients(fv, &computeGradients)); 3407 if (!computeGradients) { 3408 *gradDM = NULL; 3409 PetscFunctionReturn(PETSC_SUCCESS); 3410 } 3411 PetscCall(PetscObjectQuery((PetscObject)dm, "DMPlex_dmgrad_fvm", &gradobj)); 3412 if (!gradobj) { 3413 DM dmGradInt; 3414 3415 PetscCall(DMPlexComputeGradientFVM(dm, fv, (Vec)facegeomobj, (Vec)cellgeomobj, &dmGradInt)); 3416 PetscCall(PetscObjectCompose((PetscObject)dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt)); 3417 PetscCall(DMDestroy(&dmGradInt)); 3418 PetscCall(PetscObjectQuery((PetscObject)dm, "DMPlex_dmgrad_fvm", &gradobj)); 3419 } 3420 *gradDM = (DM)gradobj; 3421 } 3422 PetscFunctionReturn(PETSC_SUCCESS); 3423 } 3424 3425 static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work, PetscReal *resNeg, PetscReal *guess) 3426 { 3427 PetscInt l, m; 3428 3429 PetscFunctionBeginHot; 3430 if (dimC == dimR && dimR <= 3) { 3431 /* invert Jacobian, multiply */ 3432 PetscScalar det, idet; 3433 3434 switch (dimR) { 3435 case 1: 3436 invJ[0] = 1. / J[0]; 3437 break; 3438 case 2: 3439 det = J[0] * J[3] - J[1] * J[2]; 3440 idet = 1. / det; 3441 invJ[0] = J[3] * idet; 3442 invJ[1] = -J[1] * idet; 3443 invJ[2] = -J[2] * idet; 3444 invJ[3] = J[0] * idet; 3445 break; 3446 case 3: { 3447 invJ[0] = J[4] * J[8] - J[5] * J[7]; 3448 invJ[1] = J[2] * J[7] - J[1] * J[8]; 3449 invJ[2] = J[1] * J[5] - J[2] * J[4]; 3450 det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6]; 3451 idet = 1. / det; 3452 invJ[0] *= idet; 3453 invJ[1] *= idet; 3454 invJ[2] *= idet; 3455 invJ[3] = idet * (J[5] * J[6] - J[3] * J[8]); 3456 invJ[4] = idet * (J[0] * J[8] - J[2] * J[6]); 3457 invJ[5] = idet * (J[2] * J[3] - J[0] * J[5]); 3458 invJ[6] = idet * (J[3] * J[7] - J[4] * J[6]); 3459 invJ[7] = idet * (J[1] * J[6] - J[0] * J[7]); 3460 invJ[8] = idet * (J[0] * J[4] - J[1] * J[3]); 3461 } break; 3462 } 3463 for (l = 0; l < dimR; l++) { 3464 for (m = 0; m < dimC; m++) guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m]; 3465 } 3466 } else { 3467 #if defined(PETSC_USE_COMPLEX) 3468 char transpose = 'C'; 3469 #else 3470 char transpose = 'T'; 3471 #endif 3472 PetscBLASInt m, n, one = 1, worksize, info; 3473 3474 PetscCall(PetscBLASIntCast(dimR, &m)); 3475 PetscCall(PetscBLASIntCast(dimC, &n)); 3476 PetscCall(PetscBLASIntCast(dimC * dimC, &worksize)); 3477 for (l = 0; l < dimC; l++) invJ[l] = resNeg[l]; 3478 3479 PetscCallBLAS("LAPACKgels", LAPACKgels_(&transpose, &m, &n, &one, J, &m, invJ, &n, work, &worksize, &info)); 3480 PetscCheck(info == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELS %" PetscBLASInt_FMT, info); 3481 3482 for (l = 0; l < dimR; l++) guess[l] += PetscRealPart(invJ[l]); 3483 } 3484 PetscFunctionReturn(PETSC_SUCCESS); 3485 } 3486 3487 static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 3488 { 3489 PetscInt coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR); 3490 PetscScalar *coordsScalar = NULL; 3491 PetscReal *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg; 3492 PetscScalar *J, *invJ, *work; 3493 3494 PetscFunctionBegin; 3495 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3496 PetscCall(DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar)); 3497 PetscCheck(coordSize >= dimC * numV, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expecting at least %" PetscInt_FMT " coordinates, got %" PetscInt_FMT, dimC * (1 << dimR), coordSize); 3498 PetscCall(DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData)); 3499 PetscCall(DMGetWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J)); 3500 cellCoords = &cellData[0]; 3501 cellCoeffs = &cellData[coordSize]; 3502 extJ = &cellData[2 * coordSize]; 3503 resNeg = &cellData[2 * coordSize + dimR]; 3504 invJ = &J[dimR * dimC]; 3505 work = &J[2 * dimR * dimC]; 3506 if (dimR == 2) { 3507 const PetscInt zToPlex[4] = {0, 1, 3, 2}; 3508 3509 for (i = 0; i < 4; i++) { 3510 PetscInt plexI = zToPlex[i]; 3511 3512 for (j = 0; j < dimC; j++) cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 3513 } 3514 } else if (dimR == 3) { 3515 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 3516 3517 for (i = 0; i < 8; i++) { 3518 PetscInt plexI = zToPlex[i]; 3519 3520 for (j = 0; j < dimC; j++) cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 3521 } 3522 } else { 3523 for (i = 0; i < coordSize; i++) cellCoords[i] = PetscRealPart(coordsScalar[i]); 3524 } 3525 /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 3526 for (i = 0; i < dimR; i++) { 3527 PetscReal *swap; 3528 3529 for (j = 0; j < (numV / 2); j++) { 3530 for (k = 0; k < dimC; k++) { 3531 cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 3532 cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 3533 } 3534 } 3535 3536 if (i < dimR - 1) { 3537 swap = cellCoeffs; 3538 cellCoeffs = cellCoords; 3539 cellCoords = swap; 3540 } 3541 } 3542 PetscCall(PetscArrayzero(refCoords, numPoints * dimR)); 3543 for (j = 0; j < numPoints; j++) { 3544 for (i = 0; i < maxIts; i++) { 3545 PetscReal *guess = &refCoords[dimR * j]; 3546 3547 /* compute -residual and Jacobian */ 3548 for (k = 0; k < dimC; k++) resNeg[k] = realCoords[dimC * j + k]; 3549 for (k = 0; k < dimC * dimR; k++) J[k] = 0.; 3550 for (k = 0; k < numV; k++) { 3551 PetscReal extCoord = 1.; 3552 for (l = 0; l < dimR; l++) { 3553 PetscReal coord = guess[l]; 3554 PetscInt dep = (k & (1 << l)) >> l; 3555 3556 extCoord *= dep * coord + !dep; 3557 extJ[l] = dep; 3558 3559 for (m = 0; m < dimR; m++) { 3560 PetscReal coord = guess[m]; 3561 PetscInt dep = ((k & (1 << m)) >> m) && (m != l); 3562 PetscReal mult = dep * coord + !dep; 3563 3564 extJ[l] *= mult; 3565 } 3566 } 3567 for (l = 0; l < dimC; l++) { 3568 PetscReal coeff = cellCoeffs[dimC * k + l]; 3569 3570 resNeg[l] -= coeff * extCoord; 3571 for (m = 0; m < dimR; m++) J[dimR * l + m] += coeff * extJ[m]; 3572 } 3573 } 3574 if (0 && PetscDefined(USE_DEBUG)) { 3575 PetscReal maxAbs = 0.; 3576 3577 for (l = 0; l < dimC; l++) maxAbs = PetscMax(maxAbs, PetscAbsReal(resNeg[l])); 3578 PetscCall(PetscInfo(dm, "cell %" PetscInt_FMT ", point %" PetscInt_FMT ", iter %" PetscInt_FMT ": res %g\n", cell, j, i, (double)maxAbs)); 3579 } 3580 3581 PetscCall(DMPlexCoordinatesToReference_NewtonUpdate(dimC, dimR, J, invJ, work, resNeg, guess)); 3582 } 3583 } 3584 PetscCall(DMRestoreWorkArray(dm, 3 * dimR * dimC, MPIU_SCALAR, &J)); 3585 PetscCall(DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, MPIU_REAL, &cellData)); 3586 PetscCall(DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar)); 3587 PetscFunctionReturn(PETSC_SUCCESS); 3588 } 3589 3590 static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 3591 { 3592 PetscInt coordSize, i, j, k, l, numV = (1 << dimR); 3593 PetscScalar *coordsScalar = NULL; 3594 PetscReal *cellData, *cellCoords, *cellCoeffs; 3595 3596 PetscFunctionBegin; 3597 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3598 PetscCall(DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar)); 3599 PetscCheck(coordSize >= dimC * numV, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expecting at least %" PetscInt_FMT " coordinates, got %" PetscInt_FMT, dimC * (1 << dimR), coordSize); 3600 PetscCall(DMGetWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData)); 3601 cellCoords = &cellData[0]; 3602 cellCoeffs = &cellData[coordSize]; 3603 if (dimR == 2) { 3604 const PetscInt zToPlex[4] = {0, 1, 3, 2}; 3605 3606 for (i = 0; i < 4; i++) { 3607 PetscInt plexI = zToPlex[i]; 3608 3609 for (j = 0; j < dimC; j++) cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 3610 } 3611 } else if (dimR == 3) { 3612 const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 3613 3614 for (i = 0; i < 8; i++) { 3615 PetscInt plexI = zToPlex[i]; 3616 3617 for (j = 0; j < dimC; j++) cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 3618 } 3619 } else { 3620 for (i = 0; i < coordSize; i++) cellCoords[i] = PetscRealPart(coordsScalar[i]); 3621 } 3622 /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 3623 for (i = 0; i < dimR; i++) { 3624 PetscReal *swap; 3625 3626 for (j = 0; j < (numV / 2); j++) { 3627 for (k = 0; k < dimC; k++) { 3628 cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 3629 cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 3630 } 3631 } 3632 3633 if (i < dimR - 1) { 3634 swap = cellCoeffs; 3635 cellCoeffs = cellCoords; 3636 cellCoords = swap; 3637 } 3638 } 3639 PetscCall(PetscArrayzero(realCoords, numPoints * dimC)); 3640 for (j = 0; j < numPoints; j++) { 3641 const PetscReal *guess = &refCoords[dimR * j]; 3642 PetscReal *mapped = &realCoords[dimC * j]; 3643 3644 for (k = 0; k < numV; k++) { 3645 PetscReal extCoord = 1.; 3646 for (l = 0; l < dimR; l++) { 3647 PetscReal coord = guess[l]; 3648 PetscInt dep = (k & (1 << l)) >> l; 3649 3650 extCoord *= dep * coord + !dep; 3651 } 3652 for (l = 0; l < dimC; l++) { 3653 PetscReal coeff = cellCoeffs[dimC * k + l]; 3654 3655 mapped[l] += coeff * extCoord; 3656 } 3657 } 3658 } 3659 PetscCall(DMRestoreWorkArray(dm, 2 * coordSize, MPIU_REAL, &cellData)); 3660 PetscCall(DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar)); 3661 PetscFunctionReturn(PETSC_SUCCESS); 3662 } 3663 3664 PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt Nc, PetscInt dimR, PetscInt maxIter, PetscReal *tol) 3665 { 3666 PetscInt numComp, pdim, i, j, k, l, m, coordSize; 3667 PetscScalar *nodes = NULL; 3668 PetscReal *invV, *modes; 3669 PetscReal *B, *D, *resNeg; 3670 PetscScalar *J, *invJ, *work; 3671 PetscReal sq_tolerance = tol == NULL ? 0.0 : (*tol) * (*tol); 3672 3673 PetscFunctionBegin; 3674 PetscCall(PetscFEGetDimension(fe, &pdim)); 3675 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 3676 PetscCheck(numComp == Nc, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "coordinate discretization must have as many components (%" PetscInt_FMT ") as embedding dimension (!= %" PetscInt_FMT ")", numComp, Nc); 3677 /* we shouldn't apply inverse closure permutation, if one exists */ 3678 PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, coords, cell, 0, &coordSize, &nodes)); 3679 /* convert nodes to values in the stable evaluation basis */ 3680 PetscCall(DMGetWorkArray(dm, pdim, MPIU_REAL, &modes)); 3681 invV = fe->invV; 3682 for (i = 0; i < pdim; ++i) { 3683 modes[i] = 0.; 3684 for (j = 0; j < pdim; ++j) modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]); 3685 } 3686 PetscCall(DMGetWorkArray(dm, pdim * Nc + pdim * Nc * dimR + Nc, MPIU_REAL, &B)); 3687 D = &B[pdim * Nc]; 3688 resNeg = &D[pdim * Nc * dimR]; 3689 PetscCall(DMGetWorkArray(dm, 3 * Nc * dimR, MPIU_SCALAR, &J)); 3690 invJ = &J[Nc * dimR]; 3691 work = &invJ[Nc * dimR]; 3692 for (i = 0; i < numPoints * dimR; i++) refCoords[i] = 0.; 3693 for (j = 0; j < numPoints; j++) { 3694 for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */ 3695 PetscReal *guess = &refCoords[j * dimR], sq_error = 0; 3696 PetscCall(PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL)); 3697 for (k = 0; k < Nc; k++) resNeg[k] = realCoords[j * Nc + k]; 3698 for (k = 0; k < Nc * dimR; k++) J[k] = 0.; 3699 for (k = 0; k < pdim; k++) { 3700 for (l = 0; l < Nc; l++) { 3701 resNeg[l] -= modes[k] * B[k * Nc + l]; 3702 for (m = 0; m < dimR; m++) J[l * dimR + m] += modes[k] * D[(k * Nc + l) * dimR + m]; 3703 } 3704 } 3705 if (0 && PetscDefined(USE_DEBUG)) { 3706 PetscReal maxAbs = 0.; 3707 3708 for (l = 0; l < Nc; l++) maxAbs = PetscMax(maxAbs, PetscAbsReal(resNeg[l])); 3709 PetscCall(PetscInfo(dm, "cell %" PetscInt_FMT ", point %" PetscInt_FMT ", iter %" PetscInt_FMT ": res %g\n", cell, j, i, (double)maxAbs)); 3710 } 3711 for (l = 0; l < Nc; l++) sq_error += resNeg[l] * resNeg[l]; 3712 if (sq_error < sq_tolerance) { 3713 if (tol) *tol = PetscSqrtReal(sq_error); 3714 break; 3715 } 3716 PetscCall(DMPlexCoordinatesToReference_NewtonUpdate(Nc, dimR, J, invJ, work, resNeg, guess)); 3717 } 3718 } 3719 PetscCall(DMRestoreWorkArray(dm, 3 * Nc * dimR, MPIU_SCALAR, &J)); 3720 PetscCall(DMRestoreWorkArray(dm, pdim * Nc + pdim * Nc * dimR + Nc, MPIU_REAL, &B)); 3721 PetscCall(DMRestoreWorkArray(dm, pdim, MPIU_REAL, &modes)); 3722 PetscCall(DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes)); 3723 PetscFunctionReturn(PETSC_SUCCESS); 3724 } 3725 3726 /* TODO: TOBY please fix this for Nc > 1 */ 3727 PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt Nc, PetscInt dimR) 3728 { 3729 PetscInt numComp, pdim, i, j, k, l, coordSize; 3730 PetscScalar *nodes = NULL; 3731 PetscReal *invV, *modes; 3732 PetscReal *B; 3733 3734 PetscFunctionBegin; 3735 PetscCall(PetscFEGetDimension(fe, &pdim)); 3736 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 3737 PetscCheck(numComp == Nc, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "coordinate discretization must have as many components (%" PetscInt_FMT ") as embedding dimension (!= %" PetscInt_FMT ")", numComp, Nc); 3738 /* we shouldn't apply inverse closure permutation, if one exists */ 3739 PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, coords, cell, 0, &coordSize, &nodes)); 3740 /* convert nodes to values in the stable evaluation basis */ 3741 PetscCall(DMGetWorkArray(dm, pdim, MPIU_REAL, &modes)); 3742 invV = fe->invV; 3743 for (i = 0; i < pdim; ++i) { 3744 modes[i] = 0.; 3745 for (j = 0; j < pdim; ++j) modes[i] += invV[i * pdim + j] * PetscRealPart(nodes[j]); 3746 } 3747 PetscCall(DMGetWorkArray(dm, numPoints * pdim * Nc, MPIU_REAL, &B)); 3748 PetscCall(PetscSpaceEvaluate(fe->basisSpace, numPoints, refCoords, B, NULL, NULL)); 3749 for (i = 0; i < numPoints * Nc; i++) realCoords[i] = 0.; 3750 for (j = 0; j < numPoints; j++) { 3751 PetscReal *mapped = &realCoords[j * Nc]; 3752 3753 for (k = 0; k < pdim; k++) { 3754 for (l = 0; l < Nc; l++) mapped[l] += modes[k] * B[(j * pdim + k) * Nc + l]; 3755 } 3756 } 3757 PetscCall(DMRestoreWorkArray(dm, numPoints * pdim * Nc, MPIU_REAL, &B)); 3758 PetscCall(DMRestoreWorkArray(dm, pdim, MPIU_REAL, &modes)); 3759 PetscCall(DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes)); 3760 PetscFunctionReturn(PETSC_SUCCESS); 3761 } 3762 3763 /*@ 3764 DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element 3765 using a single element map. 3766 3767 Not Collective 3768 3769 Input Parameters: 3770 + dm - The mesh, with coordinate maps defined either by a `PetscDS` for the coordinate `DM` (see `DMGetCoordinateDM()`) or 3771 implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 3772 as a multilinear map for tensor-product elements 3773 . cell - the cell whose map is used. 3774 . numPoints - the number of points to locate 3775 - realCoords - (numPoints x coordinate dimension) array of coordinates (see `DMGetCoordinateDim()`) 3776 3777 Output Parameter: 3778 . refCoords - (`numPoints` x `dimension`) array of reference coordinates (see `DMGetDimension()`) 3779 3780 Level: intermediate 3781 3782 Notes: 3783 This inversion will be accurate inside the reference element, but may be inaccurate for 3784 mappings that do not extend uniquely outside the reference cell (e.g, most non-affine maps) 3785 3786 .seealso: `DMPLEX`, `DMPlexReferenceToCoordinates()` 3787 @*/ 3788 PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[]) 3789 { 3790 PetscInt dimC, dimR, depth, cStart, cEnd, i; 3791 DM coordDM = NULL; 3792 Vec coords; 3793 PetscFE fe = NULL; 3794 3795 PetscFunctionBegin; 3796 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3797 PetscCall(DMGetDimension(dm, &dimR)); 3798 PetscCall(DMGetCoordinateDim(dm, &dimC)); 3799 if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(PETSC_SUCCESS); 3800 PetscCall(DMPlexGetDepth(dm, &depth)); 3801 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 3802 PetscCall(DMGetCoordinateDM(dm, &coordDM)); 3803 if (coordDM) { 3804 PetscInt coordFields; 3805 3806 PetscCall(DMGetNumFields(coordDM, &coordFields)); 3807 if (coordFields) { 3808 PetscClassId id; 3809 PetscObject disc; 3810 3811 PetscCall(DMGetField(coordDM, 0, NULL, &disc)); 3812 PetscCall(PetscObjectGetClassId(disc, &id)); 3813 if (id == PETSCFE_CLASSID) fe = (PetscFE)disc; 3814 } 3815 } 3816 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 3817 PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " not in cell range [%" PetscInt_FMT ",%" PetscInt_FMT ")", cell, cStart, cEnd); 3818 if (!fe) { /* implicit discretization: affine or multilinear */ 3819 PetscInt coneSize; 3820 PetscBool isSimplex, isTensor; 3821 3822 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 3823 isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 3824 isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 3825 if (isSimplex) { 3826 PetscReal detJ, *v0, *J, *invJ; 3827 3828 PetscCall(DMGetWorkArray(dm, dimC + 2 * dimC * dimC, MPIU_REAL, &v0)); 3829 J = &v0[dimC]; 3830 invJ = &J[dimC * dimC]; 3831 PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ)); 3832 for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */ 3833 const PetscReal x0[3] = {-1., -1., -1.}; 3834 3835 CoordinatesRealToRef(dimC, dimR, x0, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]); 3836 } 3837 PetscCall(DMRestoreWorkArray(dm, dimC + 2 * dimC * dimC, MPIU_REAL, &v0)); 3838 } else if (isTensor) { 3839 PetscCall(DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR)); 3840 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unrecognized cone size %" PetscInt_FMT, coneSize); 3841 } else { 3842 PetscCall(DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR, 7, NULL)); 3843 } 3844 PetscFunctionReturn(PETSC_SUCCESS); 3845 } 3846 3847 /*@ 3848 DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the mesh for a single element map. 3849 3850 Not Collective 3851 3852 Input Parameters: 3853 + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate `DM` (see `DMGetCoordinateDM()`) or 3854 implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 3855 as a multilinear map for tensor-product elements 3856 . cell - the cell whose map is used. 3857 . numPoints - the number of points to locate 3858 - refCoords - (numPoints x dimension) array of reference coordinates (see `DMGetDimension()`) 3859 3860 Output Parameter: 3861 . realCoords - (numPoints x coordinate dimension) array of coordinates (see `DMGetCoordinateDim()`) 3862 3863 Level: intermediate 3864 3865 .seealso: `DMPLEX`, `DMPlexCoordinatesToReference()` 3866 @*/ 3867 PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[]) 3868 { 3869 PetscInt dimC, dimR, depth, cStart, cEnd, i; 3870 DM coordDM = NULL; 3871 Vec coords; 3872 PetscFE fe = NULL; 3873 3874 PetscFunctionBegin; 3875 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3876 PetscCall(DMGetDimension(dm, &dimR)); 3877 PetscCall(DMGetCoordinateDim(dm, &dimC)); 3878 if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(PETSC_SUCCESS); 3879 PetscCall(DMPlexGetDepth(dm, &depth)); 3880 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 3881 PetscCall(DMGetCoordinateDM(dm, &coordDM)); 3882 if (coordDM) { 3883 PetscInt coordFields; 3884 3885 PetscCall(DMGetNumFields(coordDM, &coordFields)); 3886 if (coordFields) { 3887 PetscClassId id; 3888 PetscObject disc; 3889 3890 PetscCall(DMGetField(coordDM, 0, NULL, &disc)); 3891 PetscCall(PetscObjectGetClassId(disc, &id)); 3892 if (id == PETSCFE_CLASSID) fe = (PetscFE)disc; 3893 } 3894 } 3895 PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 3896 PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " not in cell range [%" PetscInt_FMT ",%" PetscInt_FMT ")", cell, cStart, cEnd); 3897 if (!fe) { /* implicit discretization: affine or multilinear */ 3898 PetscInt coneSize; 3899 PetscBool isSimplex, isTensor; 3900 3901 PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 3902 isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 3903 isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 3904 if (isSimplex) { 3905 PetscReal detJ, *v0, *J; 3906 3907 PetscCall(DMGetWorkArray(dm, dimC + 2 * dimC * dimC, MPIU_REAL, &v0)); 3908 J = &v0[dimC]; 3909 PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ)); 3910 for (i = 0; i < numPoints; i++) { /* Apply the affine transformation for each point */ 3911 const PetscReal xi0[3] = {-1., -1., -1.}; 3912 3913 CoordinatesRefToReal(dimC, dimR, xi0, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]); 3914 } 3915 PetscCall(DMRestoreWorkArray(dm, dimC + 2 * dimC * dimC, MPIU_REAL, &v0)); 3916 } else if (isTensor) { 3917 PetscCall(DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR)); 3918 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unrecognized cone size %" PetscInt_FMT, coneSize); 3919 } else { 3920 PetscCall(DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR)); 3921 } 3922 PetscFunctionReturn(PETSC_SUCCESS); 3923 } 3924 3925 void coordMap_identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 3926 { 3927 const PetscInt Nc = uOff[1] - uOff[0]; 3928 PetscInt c; 3929 3930 for (c = 0; c < Nc; ++c) f0[c] = u[c]; 3931 } 3932 3933 /* Shear applies the transformation, assuming we fix z, 3934 / 1 0 m_0 \ 3935 | 0 1 m_1 | 3936 \ 0 0 1 / 3937 */ 3938 void coordMap_shear(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[]) 3939 { 3940 const PetscInt Nc = uOff[1] - uOff[0]; 3941 const PetscInt ax = (PetscInt)PetscRealPart(constants[0]); 3942 PetscInt c; 3943 3944 for (c = 0; c < Nc; ++c) coords[c] = u[c] + constants[c + 1] * u[ax]; 3945 } 3946 3947 /* Flare applies the transformation, assuming we fix x_f, 3948 3949 x_i = x_i * alpha_i x_f 3950 */ 3951 void coordMap_flare(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[]) 3952 { 3953 const PetscInt Nc = uOff[1] - uOff[0]; 3954 const PetscInt cf = (PetscInt)PetscRealPart(constants[0]); 3955 PetscInt c; 3956 3957 for (c = 0; c < Nc; ++c) coords[c] = u[c] * (c == cf ? 1.0 : constants[c + 1] * u[cf]); 3958 } 3959 3960 /* 3961 We would like to map the unit square to a quarter of the annulus between circles of radius 1 and 2. We start by mapping the straight sections, which 3962 will correspond to the top and bottom of our square. So 3963 3964 (0,0)--(1,0) ==> (1,0)--(2,0) Just a shift of (1,0) 3965 (0,1)--(1,1) ==> (0,1)--(0,2) Switch x and y 3966 3967 So it looks like we want to map each layer in y to a ray, so x is the radius and y is the angle: 3968 3969 (x, y) ==> (x+1, \pi/2 y) in (r', \theta') space 3970 ==> ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space 3971 */ 3972 void coordMap_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 3973 { 3974 const PetscReal ri = PetscRealPart(constants[0]); 3975 const PetscReal ro = PetscRealPart(constants[1]); 3976 3977 xp[0] = (x[0] * (ro - ri) + ri) * PetscCosReal(0.5 * PETSC_PI * x[1]); 3978 xp[1] = (x[0] * (ro - ri) + ri) * PetscSinReal(0.5 * PETSC_PI * x[1]); 3979 } 3980 3981 /* 3982 We would like to map the unit cube to a hemisphere of the spherical shell between balls of radius 1 and 2. We want to map the bottom surface onto the 3983 lower hemisphere and the upper surface onto the top, letting z be the radius. 3984 3985 (x, y) ==> ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x) in (r', \theta', \phi') space 3986 ==> ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space 3987 */ 3988 void coordMap_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) 3989 { 3990 const PetscReal pi4 = PETSC_PI / 4.0; 3991 const PetscReal ri = PetscRealPart(constants[0]); 3992 const PetscReal ro = PetscRealPart(constants[1]); 3993 const PetscReal rp = (x[2] + 1) * 0.5 * (ro - ri) + ri; 3994 const PetscReal phip = PetscAtan2Real(x[1], x[0]); 3995 const PetscReal thetap = 0.5 * PETSC_PI * (1.0 - ((((phip <= pi4) && (phip >= -pi4)) || ((phip >= 3.0 * pi4) || (phip <= -3.0 * pi4))) ? PetscAbsReal(x[0]) : PetscAbsReal(x[1]))); 3996 3997 xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip); 3998 xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip); 3999 xp[2] = rp * PetscSinReal(thetap); 4000 } 4001 4002 /*@C 4003 DMPlexRemapGeometry - This function maps the original `DM` coordinates to new coordinates. 4004 4005 Not Collective 4006 4007 Input Parameters: 4008 + dm - The `DM` 4009 . time - The time 4010 - func - The function transforming current coordinates to new coordinates 4011 4012 Calling sequence of `func`: 4013 + dim - The spatial dimension 4014 . Nf - The number of input fields (here 1) 4015 . NfAux - The number of input auxiliary fields 4016 . uOff - The offset of the coordinates in u[] (here 0) 4017 . uOff_x - The offset of the coordinates in u_x[] (here 0) 4018 . u - The coordinate values at this point in space 4019 . u_t - The coordinate time derivative at this point in space (here `NULL`) 4020 . u_x - The coordinate derivatives at this point in space 4021 . aOff - The offset of each auxiliary field in u[] 4022 . aOff_x - The offset of each auxiliary field in u_x[] 4023 . a - The auxiliary field values at this point in space 4024 . a_t - The auxiliary field time derivative at this point in space (or `NULL`) 4025 . a_x - The auxiliary field derivatives at this point in space 4026 . t - The current time 4027 . x - The coordinates of this point (here not used) 4028 . numConstants - The number of constants 4029 . constants - The value of each constant 4030 - f - The new coordinates at this point in space 4031 4032 Level: intermediate 4033 4034 .seealso: `DMPLEX`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMProjectFieldLocal()`, `DMProjectFieldLabelLocal()` 4035 @*/ 4036 PetscErrorCode DMPlexRemapGeometry(DM dm, PetscReal time, void (*func)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])) 4037 { 4038 DM cdm; 4039 PetscDS cds; 4040 DMField cf; 4041 PetscObject obj; 4042 PetscClassId id; 4043 Vec lCoords, tmpCoords; 4044 4045 PetscFunctionBegin; 4046 PetscCall(DMGetCoordinateDM(dm, &cdm)); 4047 PetscCall(DMGetCoordinatesLocal(dm, &lCoords)); 4048 PetscCall(DMGetDS(cdm, &cds)); 4049 PetscCall(PetscDSGetDiscretization(cds, 0, &obj)); 4050 PetscCall(PetscObjectGetClassId(obj, &id)); 4051 if (id != PETSCFE_CLASSID) { 4052 PetscSection cSection; 4053 const PetscScalar *constants; 4054 PetscScalar *coords, f[16]; 4055 PetscInt dim, cdim, Nc, vStart, vEnd; 4056 4057 PetscCall(DMGetDimension(dm, &dim)); 4058 PetscCall(DMGetCoordinateDim(dm, &cdim)); 4059 PetscCheck(cdim <= 16, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Affine version of DMPlexRemapGeometry is currently limited to dimensions <= 16, not %" PetscInt_FMT, cdim); 4060 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4061 PetscCall(DMGetCoordinateSection(dm, &cSection)); 4062 PetscCall(PetscDSGetConstants(cds, &Nc, &constants)); 4063 PetscCall(VecGetArrayWrite(lCoords, &coords)); 4064 for (PetscInt v = vStart; v < vEnd; ++v) { 4065 PetscInt uOff[2] = {0, cdim}; 4066 PetscInt off, c; 4067 4068 PetscCall(PetscSectionGetOffset(cSection, v, &off)); 4069 (*func)(dim, 1, 0, uOff, NULL, &coords[off], NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0.0, NULL, Nc, constants, f); 4070 for (c = 0; c < cdim; ++c) coords[off + c] = f[c]; 4071 } 4072 PetscCall(VecRestoreArrayWrite(lCoords, &coords)); 4073 } else { 4074 PetscCall(DMGetLocalVector(cdm, &tmpCoords)); 4075 PetscCall(VecCopy(lCoords, tmpCoords)); 4076 /* We have to do the coordinate field manually right now since the coordinate DM will not have its own */ 4077 PetscCall(DMGetCoordinateField(dm, &cf)); 4078 cdm->coordinates[0].field = cf; 4079 PetscCall(DMProjectFieldLocal(cdm, time, tmpCoords, &func, INSERT_VALUES, lCoords)); 4080 cdm->coordinates[0].field = NULL; 4081 PetscCall(DMRestoreLocalVector(cdm, &tmpCoords)); 4082 PetscCall(DMSetCoordinatesLocal(dm, lCoords)); 4083 } 4084 PetscFunctionReturn(PETSC_SUCCESS); 4085 } 4086 4087 /*@ 4088 DMPlexShearGeometry - This shears the domain, meaning adds a multiple of the shear coordinate to all other coordinates. 4089 4090 Not Collective 4091 4092 Input Parameters: 4093 + dm - The `DMPLEX` 4094 . direction - The shear coordinate direction, e.g. `DM_X` is the x-axis 4095 - multipliers - The multiplier m for each direction which is not the shear direction 4096 4097 Level: intermediate 4098 4099 .seealso: `DMPLEX`, `DMPlexRemapGeometry()`, `DMDirection`, `DM_X`, `DM_Y`, `DM_Z` 4100 @*/ 4101 PetscErrorCode DMPlexShearGeometry(DM dm, DMDirection direction, PetscReal multipliers[]) 4102 { 4103 DM cdm; 4104 PetscDS cds; 4105 PetscScalar *moduli; 4106 const PetscInt dir = (PetscInt)direction; 4107 PetscInt dE, d, e; 4108 4109 PetscFunctionBegin; 4110 PetscCall(DMGetCoordinateDM(dm, &cdm)); 4111 PetscCall(DMGetCoordinateDim(dm, &dE)); 4112 PetscCall(PetscMalloc1(dE + 1, &moduli)); 4113 moduli[0] = dir; 4114 for (d = 0, e = 0; d < dE; ++d) moduli[d + 1] = d == dir ? 0.0 : (multipliers ? multipliers[e++] : 1.0); 4115 PetscCall(DMGetDS(cdm, &cds)); 4116 PetscCall(PetscDSSetConstants(cds, dE + 1, moduli)); 4117 PetscCall(DMPlexRemapGeometry(dm, 0.0, coordMap_shear)); 4118 PetscCall(PetscFree(moduli)); 4119 PetscFunctionReturn(PETSC_SUCCESS); 4120 } 4121