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__ "DMPlexComputeProjection3Dto2D_Internal" 206 /* 207 DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 208 */ 209 static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[]) 210 { 211 PetscScalar x1[3], x2[3], n[3], norm; 212 PetscScalar R[9], x1p[3], x2p[3]; 213 PetscScalar sqrtz, alpha; 214 const PetscInt dim = 3; 215 PetscInt d, e; 216 217 PetscFunctionBegin; 218 /* 0) Calculate normal vector */ 219 for (d = 0; d < dim; ++d) { 220 x1[d] = coords[1*dim+d] - coords[0*dim+d]; 221 x2[d] = coords[2*dim+d] - coords[0*dim+d]; 222 } 223 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 224 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 225 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 226 norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 227 n[0] /= norm; 228 n[1] /= norm; 229 n[2] /= norm; 230 /* 1) Take the normal vector and rotate until it is \hat z 231 232 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 233 234 R = / alpha nx nz alpha ny nz -1/alpha \ 235 | -alpha ny alpha nx 0 | 236 \ nx ny nz / 237 238 will rotate the normal vector to \hat z 239 */ 240 sqrtz = sqrt(1.0 - n[2]*n[2]); 241 alpha = 1.0/sqrtz; 242 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 243 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 244 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 245 for (d = 0; d < dim; ++d) { 246 x1p[d] = 0.0; 247 x2p[d] = 0.0; 248 for (e = 0; e < dim; ++e) { 249 x1p[d] += R[d*dim+e]*x1[e]; 250 x2p[d] += R[d*dim+e]*x2[e]; 251 } 252 } 253 if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 254 if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 255 /* 2) Project to (x, y) */ 256 coords[0] = 0.0; 257 coords[1] = 0.0; 258 coords[2] = x1p[0]; 259 coords[3] = x1p[1]; 260 coords[4] = x2p[0]; 261 coords[5] = x2p[1]; 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 267 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 268 { 269 PetscSection coordSection; 270 Vec coordinates; 271 PetscScalar *coords; 272 const PetscInt dim = 2; 273 PetscInt numCoords, d, f; 274 PetscErrorCode ierr; 275 276 PetscFunctionBegin; 277 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 278 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 279 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 280 if (numCoords == 9) { 281 ierr = DMPlexComputeProjection3Dto2D_Internal(coords);CHKERRQ(ierr); 282 } else if (numCoords != 6) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords); 283 if (v0) { 284 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 285 } 286 if (J) { 287 for (d = 0; d < dim; d++) { 288 for (f = 0; f < dim; f++) { 289 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 290 } 291 } 292 *detJ = J[0]*J[3] - J[1]*J[2]; 293 #if 0 294 if (detJ < 0.0) { 295 const PetscReal xLength = mesh->periodicity[0]; 296 297 if (xLength != 0.0) { 298 PetscReal v0x = coords[0*dim+0]; 299 300 if (v0x == 0.0) v0x = v0[0] = xLength; 301 for (f = 0; f < dim; f++) { 302 const PetscReal px = coords[(f+1)*dim+0] == 0.0 ? xLength : coords[(f+1)*dim+0]; 303 304 J[0*dim+f] = 0.5*(px - v0x); 305 } 306 } 307 detJ = J[0]*J[3] - J[1]*J[2]; 308 } 309 #endif 310 PetscLogFlops(8.0 + 3.0); 311 } else { 312 *detJ = 0.0; 313 } 314 if (invJ) { 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 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 324 PetscFunctionReturn(0); 325 } 326 327 #undef __FUNCT__ 328 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 329 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 330 { 331 PetscSection coordSection; 332 Vec coordinates; 333 PetscScalar *coords; 334 const PetscInt dim = 2; 335 PetscInt d, f; 336 PetscErrorCode ierr; 337 338 PetscFunctionBegin; 339 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 340 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 341 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 342 if (v0) { 343 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 344 } 345 if (J) { 346 for (d = 0; d < dim; d++) { 347 for (f = 0; f < dim; f++) { 348 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 349 } 350 } 351 *detJ = J[0]*J[3] - J[1]*J[2]; 352 PetscLogFlops(8.0 + 3.0); 353 } else { 354 *detJ = 0.0; 355 } 356 if (invJ) { 357 const PetscReal invDet = 1.0/(*detJ); 358 359 invJ[0] = invDet*J[3]; 360 invJ[1] = -invDet*J[1]; 361 invJ[2] = -invDet*J[2]; 362 invJ[3] = invDet*J[0]; 363 PetscLogFlops(5.0); 364 } 365 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 366 PetscFunctionReturn(0); 367 } 368 369 #undef __FUNCT__ 370 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 371 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 372 { 373 PetscSection coordSection; 374 Vec coordinates; 375 PetscScalar *coords; 376 const PetscInt dim = 3; 377 PetscInt d, f; 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 382 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 383 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 384 if (v0) { 385 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 386 } 387 if (J) { 388 for (d = 0; d < dim; d++) { 389 for (f = 0; f < dim; f++) { 390 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 391 } 392 } 393 *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) + 394 J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) + 395 J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1])); 396 PetscLogFlops(18.0 + 12.0); 397 } 398 if (invJ) { 399 const PetscReal invDet = 1.0/(*detJ); 400 401 invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 402 invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 403 invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 404 invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 405 invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 406 invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 407 invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 408 invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 409 invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 410 PetscLogFlops(37.0); 411 } 412 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 413 PetscFunctionReturn(0); 414 } 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 418 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 419 { 420 PetscSection coordSection; 421 Vec coordinates; 422 PetscScalar *coords; 423 const PetscInt dim = 3; 424 PetscInt d; 425 PetscErrorCode ierr; 426 427 PetscFunctionBegin; 428 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 429 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 430 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 431 if (v0) { 432 for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 433 } 434 if (J) { 435 for (d = 0; d < dim; d++) { 436 J[d*dim+0] = 0.5*(PetscRealPart(coords[(0+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 437 J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 438 J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 439 } 440 *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) + 441 J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) + 442 J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1])); 443 PetscLogFlops(18.0 + 12.0); 444 } else { 445 *detJ = 0.0; 446 } 447 if (invJ) { 448 const PetscReal invDet = -1.0/(*detJ); 449 450 invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 451 invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 452 invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 453 invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 454 invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 455 invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 456 invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 457 invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 458 invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 459 PetscLogFlops(37.0); 460 } 461 *detJ *= 8.0; 462 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "DMPlexComputeCellGeometry" 468 /*@C 469 DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 470 471 Collective on DM 472 473 Input Arguments: 474 + dm - the DM 475 - cell - the cell 476 477 Output Arguments: 478 + v0 - the translation part of this affine transform 479 . J - the Jacobian of the transform from the reference element 480 . invJ - the inverse of the Jacobian 481 - detJ - the Jacobian determinant 482 483 Level: advanced 484 485 Fortran Notes: 486 Since it returns arrays, this routine is only available in Fortran 90, and you must 487 include petsc.h90 in your code. 488 489 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 490 @*/ 491 PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 492 { 493 PetscInt depth, dim, coneSize; 494 PetscErrorCode ierr; 495 496 PetscFunctionBegin; 497 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 498 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 499 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 500 if (depth == 1) { 501 switch (dim) { 502 case 2: 503 switch (coneSize) { 504 case 3: 505 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 506 break; 507 case 4: 508 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 509 break; 510 default: 511 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 512 } 513 break; 514 case 3: 515 switch (coneSize) { 516 case 4: 517 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 518 break; 519 case 8: 520 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 521 break; 522 default: 523 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 524 } 525 break; 526 default: 527 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 528 } 529 } else { 530 /* Cone size is now the number of faces */ 531 switch (dim) { 532 case 2: 533 switch (coneSize) { 534 case 3: 535 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 536 break; 537 case 4: 538 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 539 break; 540 default: 541 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 542 } 543 break; 544 case 3: 545 switch (coneSize) { 546 case 4: 547 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 548 break; 549 case 6: 550 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 551 break; 552 default: 553 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 554 } 555 break; 556 default: 557 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 558 } 559 } 560 PetscFunctionReturn(0); 561 } 562