1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal" 5 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 6 { 7 const PetscInt embedDim = 2; 8 PetscReal x = PetscRealPart(point[0]); 9 PetscReal y = PetscRealPart(point[1]); 10 PetscReal v0[2], J[4], invJ[4], detJ; 11 PetscReal xi, eta; 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 16 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 17 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 18 19 if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) *cell = c; 20 else *cell = -1; 21 PetscFunctionReturn(0); 22 } 23 24 #undef __FUNCT__ 25 #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal" 26 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 27 { 28 PetscSection coordSection; 29 Vec coordsLocal; 30 PetscScalar *coords; 31 const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 32 PetscReal x = PetscRealPart(point[0]); 33 PetscReal y = PetscRealPart(point[1]); 34 PetscInt crossings = 0, f; 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 39 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 40 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 41 for (f = 0; f < 4; ++f) { 42 PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 43 PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 44 PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 45 PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 46 PetscReal slope = (y_j - y_i) / (x_j - x_i); 47 PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 48 PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 49 PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 50 if ((cond1 || cond2) && above) ++crossings; 51 } 52 if (crossings % 2) *cell = c; 53 else *cell = -1; 54 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal" 60 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 61 { 62 const PetscInt embedDim = 3; 63 PetscReal v0[3], J[9], invJ[9], detJ; 64 PetscReal x = PetscRealPart(point[0]); 65 PetscReal y = PetscRealPart(point[1]); 66 PetscReal z = PetscRealPart(point[2]); 67 PetscReal xi, eta, zeta; 68 PetscErrorCode ierr; 69 70 PetscFunctionBegin; 71 ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 72 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 73 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 74 zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 75 76 if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 77 else *cell = -1; 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal" 83 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 84 { 85 PetscSection coordSection; 86 Vec coordsLocal; 87 PetscScalar *coords; 88 const PetscInt faces[24] = {0, 1, 2, 3, 5, 4, 7, 6, 1, 0, 4, 5, 89 3, 2, 6, 7, 1, 5, 6, 2, 0, 3, 7, 4}; 90 PetscBool found = PETSC_TRUE; 91 PetscInt f; 92 PetscErrorCode ierr; 93 94 PetscFunctionBegin; 95 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 96 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 97 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 98 for (f = 0; f < 6; ++f) { 99 /* Check the point is under plane */ 100 /* Get face normal */ 101 PetscReal v_i[3]; 102 PetscReal v_j[3]; 103 PetscReal normal[3]; 104 PetscReal pp[3]; 105 PetscReal dot; 106 107 v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 108 v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 109 v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 110 v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 111 v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 112 v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 113 normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 114 normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 115 normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 116 pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 117 pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 118 pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 119 dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 120 121 /* Check that projected point is in face (2D location problem) */ 122 if (dot < 0.0) { 123 found = PETSC_FALSE; 124 break; 125 } 126 } 127 if (found) *cell = c; 128 else *cell = -1; 129 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "DMLocatePoints_Plex" 135 /* 136 Need to implement using the guess 137 */ 138 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS) 139 { 140 PetscInt cell = -1 /*, guess = -1*/; 141 PetscInt bs, numPoints, p; 142 PetscInt dim, cStart, cEnd, cMax, c, coneSize; 143 PetscInt *cells; 144 PetscScalar *a; 145 PetscErrorCode ierr; 146 147 PetscFunctionBegin; 148 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 149 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 150 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 151 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 152 ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 153 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 154 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 155 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); 156 numPoints /= bs; 157 ierr = PetscMalloc(numPoints * sizeof(PetscInt), &cells);CHKERRQ(ierr); 158 for (p = 0; p < numPoints; ++p) { 159 const PetscScalar *point = &a[p*bs]; 160 161 switch (dim) { 162 case 2: 163 for (c = cStart; c < cEnd; ++c) { 164 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 165 switch (coneSize) { 166 case 3: 167 ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 168 break; 169 case 4: 170 ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 171 break; 172 default: 173 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize); 174 } 175 if (cell >= 0) break; 176 } 177 break; 178 case 3: 179 for (c = cStart; c < cEnd; ++c) { 180 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 181 switch (coneSize) { 182 case 4: 183 ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 184 break; 185 case 8: 186 ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 187 break; 188 default: 189 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %d", coneSize); 190 } 191 if (cell >= 0) break; 192 } 193 break; 194 default: 195 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %d", dim); 196 } 197 cells[p] = cell; 198 } 199 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 200 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 206 /* 207 DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 208 */ 209 static PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 210 { 211 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 212 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 213 const PetscReal r = sqrt(x*x + y*y), c = x/r, s = y/r; 214 215 PetscFunctionBegin; 216 R[0] = c; R[1] = s; 217 R[2] = -s; R[3] = c; 218 coords[0] = 0.0; 219 coords[1] = r; 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 225 /* 226 DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 227 */ 228 static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[], PetscReal R[]) 229 { 230 PetscScalar x1[3], x2[3], n[3], norm; 231 PetscScalar x1p[3], x2p[3]; 232 PetscReal sqrtz, alpha; 233 const PetscInt dim = 3; 234 PetscInt d, e; 235 236 PetscFunctionBegin; 237 /* 0) Calculate normal vector */ 238 for (d = 0; d < dim; ++d) { 239 x1[d] = coords[1*dim+d] - coords[0*dim+d]; 240 x2[d] = coords[2*dim+d] - coords[0*dim+d]; 241 } 242 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 243 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 244 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 245 norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 246 n[0] /= norm; 247 n[1] /= norm; 248 n[2] /= norm; 249 /* 1) Take the normal vector and rotate until it is \hat z 250 251 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 252 253 R = / alpha nx nz alpha ny nz -1/alpha \ 254 | -alpha ny alpha nx 0 | 255 \ nx ny nz / 256 257 will rotate the normal vector to \hat z 258 */ 259 sqrtz = sqrt(1.0 - PetscAbsScalar(n[2]*n[2])); 260 /* Check for n = z */ 261 if (sqrtz < 1.0e-10) { 262 coords[0] = 0.0; 263 coords[1] = 0.0; 264 if (PetscRealPart(n[2]) < 0.0) { 265 coords[2] = x2[0]; 266 coords[3] = x2[1]; 267 coords[4] = x1[0]; 268 coords[5] = x1[1]; 269 } else { 270 coords[2] = x1[0]; 271 coords[3] = x1[1]; 272 coords[4] = x2[0]; 273 coords[5] = x2[1]; 274 } 275 PetscFunctionReturn(0); 276 } 277 alpha = 1.0/sqrtz; 278 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 279 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 280 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 281 for (d = 0; d < dim; ++d) { 282 x1p[d] = 0.0; 283 x2p[d] = 0.0; 284 for (e = 0; e < dim; ++e) { 285 x1p[d] += R[d*dim+e]*x1[e]; 286 x2p[d] += R[d*dim+e]*x2[e]; 287 } 288 } 289 if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 290 if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 291 /* 2) Project to (x, y) */ 292 coords[0] = 0.0; 293 coords[1] = 0.0; 294 coords[2] = x1p[0]; 295 coords[3] = x1p[1]; 296 coords[4] = x2p[0]; 297 coords[5] = x2p[1]; 298 /* Output R^T which rotates \hat z to the input normal */ 299 for (d = 0; d < dim; ++d) { 300 for (e = d+1; e < dim; ++e) { 301 PetscReal tmp; 302 303 tmp = R[d*dim+e]; 304 R[d*dim+e] = R[e*dim+d]; 305 R[e*dim+d] = tmp; 306 } 307 } 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "Invert2D_Internal" 313 PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 314 { 315 const PetscReal invDet = 1.0/detJ; 316 317 invJ[0] = invDet*J[3]; 318 invJ[1] = -invDet*J[1]; 319 invJ[2] = -invDet*J[2]; 320 invJ[3] = invDet*J[0]; 321 PetscLogFlops(5.0); 322 } 323 324 #undef __FUNCT__ 325 #define __FUNCT__ "Invert3D_Internal" 326 PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 327 { 328 const PetscReal invDet = 1.0/detJ; 329 330 invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 331 invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 332 invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 333 invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 334 invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 335 invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 336 invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 337 invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 338 invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 339 PetscLogFlops(37.0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "Det2D_Internal" 344 PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[]) 345 { 346 *detJ = J[0]*J[3] - J[1]*J[2]; 347 PetscLogFlops(3.0); 348 } 349 350 #undef __FUNCT__ 351 #define __FUNCT__ "Det3D_Internal" 352 PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[]) 353 { 354 *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) + 355 J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) + 356 J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1])); 357 PetscLogFlops(12.0); 358 } 359 360 #undef __FUNCT__ 361 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 362 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 363 { 364 PetscSection coordSection; 365 Vec coordinates; 366 PetscScalar *coords; 367 PetscInt numCoords, d; 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 372 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 373 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 374 *detJ = 0.0; 375 if (numCoords == 4) { 376 const PetscInt dim = 2; 377 PetscReal R[4], J0; 378 379 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 380 ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 381 if (J) { 382 J0 = 0.5*PetscRealPart(coords[1]); 383 J[0] = R[0]*J0; J[1] = R[1]; 384 J[2] = R[2]*J0; J[3] = R[3]; 385 Det2D_Internal(detJ, J); 386 } 387 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 388 } else if (numCoords == 2) { 389 const PetscInt dim = 1; 390 391 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 392 if (J) { 393 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 394 *detJ = J[0]; 395 PetscLogFlops(2.0); 396 } 397 if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 398 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords); 399 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 400 PetscFunctionReturn(0); 401 } 402 403 #undef __FUNCT__ 404 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 405 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 406 { 407 PetscSection coordSection; 408 Vec coordinates; 409 PetscScalar *coords; 410 PetscInt numCoords, d, f, g; 411 PetscErrorCode ierr; 412 413 PetscFunctionBegin; 414 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 415 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 416 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 417 *detJ = 0.0; 418 if (numCoords == 9) { 419 const PetscInt dim = 3; 420 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 421 422 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 423 ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr); 424 if (J) { 425 for (d = 0; d < dim-1; d++) { 426 for (f = 0; f < dim-1; f++) { 427 J0[d*(dim+1)+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 428 } 429 } 430 PetscLogFlops(8.0); 431 for (d = 0; d < dim; d++) { 432 for (f = 0; f < dim; f++) { 433 J[d*dim+f] = 0.0; 434 for (g = 0; g < dim; g++) { 435 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 436 } 437 } 438 } 439 PetscLogFlops(18.0); 440 Det3D_Internal(detJ, J); 441 } 442 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 443 } else if (numCoords == 6) { 444 const PetscInt dim = 2; 445 446 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 447 if (J) { 448 for (d = 0; d < dim; d++) { 449 for (f = 0; f < dim; f++) { 450 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 451 } 452 } 453 PetscLogFlops(8.0); 454 Det2D_Internal(detJ, J); 455 } 456 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 457 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords); 458 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 459 PetscFunctionReturn(0); 460 } 461 462 #undef __FUNCT__ 463 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 464 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 465 { 466 PetscSection coordSection; 467 Vec coordinates; 468 PetscScalar *coords; 469 const PetscInt dim = 2; 470 PetscInt d, f; 471 PetscErrorCode ierr; 472 473 PetscFunctionBegin; 474 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 475 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 476 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 477 *detJ = 0.0; 478 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 479 if (J) { 480 for (d = 0; d < dim; d++) { 481 for (f = 0; f < dim; f++) { 482 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 483 } 484 } 485 PetscLogFlops(8.0); 486 Det2D_Internal(detJ, J); 487 } 488 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 489 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 490 PetscFunctionReturn(0); 491 } 492 493 #undef __FUNCT__ 494 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 495 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 496 { 497 PetscSection coordSection; 498 Vec coordinates; 499 PetscScalar *coords; 500 const PetscInt dim = 3; 501 PetscInt d, f; 502 PetscErrorCode ierr; 503 504 PetscFunctionBegin; 505 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 506 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 507 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 508 *detJ = 0.0; 509 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 510 if (J) { 511 for (d = 0; d < dim; d++) { 512 for (f = 0; f < dim; f++) { 513 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 514 } 515 } 516 PetscLogFlops(18.0); 517 Det3D_Internal(detJ, J); 518 } 519 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 520 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 521 PetscFunctionReturn(0); 522 } 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 526 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 527 { 528 PetscSection coordSection; 529 Vec coordinates; 530 PetscScalar *coords; 531 const PetscInt dim = 3; 532 PetscInt d; 533 PetscErrorCode ierr; 534 535 PetscFunctionBegin; 536 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 537 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 538 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 539 *detJ = 0.0; 540 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 541 if (J) { 542 for (d = 0; d < dim; d++) { 543 J[d*dim+0] = 0.5*(PetscRealPart(coords[(0+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 544 J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 545 J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 546 } 547 PetscLogFlops(18.0); 548 Det3D_Internal(detJ, J); 549 } 550 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 551 *detJ *= -8.0; 552 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 553 PetscFunctionReturn(0); 554 } 555 556 #undef __FUNCT__ 557 #define __FUNCT__ "DMPlexComputeCellGeometry" 558 /*@C 559 DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 560 561 Collective on DM 562 563 Input Arguments: 564 + dm - the DM 565 - cell - the cell 566 567 Output Arguments: 568 + v0 - the translation part of this affine transform 569 . J - the Jacobian of the transform from the reference element 570 . invJ - the inverse of the Jacobian 571 - detJ - the Jacobian determinant 572 573 Level: advanced 574 575 Fortran Notes: 576 Since it returns arrays, this routine is only available in Fortran 90, and you must 577 include petsc.h90 in your code. 578 579 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 580 @*/ 581 PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 582 { 583 PetscInt depth, dim, coneSize; 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 588 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 589 if (depth == 1) { 590 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 591 switch (dim) { 592 case 1: 593 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 594 break; 595 case 2: 596 switch (coneSize) { 597 case 3: 598 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 599 break; 600 case 4: 601 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 602 break; 603 default: 604 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 605 } 606 break; 607 case 3: 608 switch (coneSize) { 609 case 4: 610 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 611 break; 612 case 8: 613 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 614 break; 615 default: 616 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 617 } 618 break; 619 default: 620 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 621 } 622 } else { 623 /* We need to keep a pointer to the depth label */ 624 ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr); 625 /* Cone size is now the number of faces */ 626 switch (dim) { 627 case 1: 628 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 629 break; 630 case 2: 631 switch (coneSize) { 632 case 3: 633 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 634 break; 635 case 4: 636 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 637 break; 638 default: 639 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 640 } 641 break; 642 case 3: 643 switch (coneSize) { 644 case 4: 645 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 646 break; 647 case 6: 648 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 649 break; 650 default: 651 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 652 } 653 break; 654 default: 655 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 656 } 657 } 658 PetscFunctionReturn(0); 659 } 660