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