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