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 PetscReal x1[3], x2[3], n[3], norm; 231 PetscReal 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] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 240 x2[d] = PetscRealPart(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 (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 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 270 R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 271 R[6] = 0.0; R[7] = 0.0; R[8] = -1.0; 272 } else { 273 coords[2] = x1[0]; 274 coords[3] = x1[1]; 275 coords[4] = x2[0]; 276 coords[5] = x2[1]; 277 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 278 R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 279 R[6] = 0.0; R[7] = 0.0; R[8] = 1.0; 280 } 281 PetscFunctionReturn(0); 282 } 283 alpha = 1.0/sqrtz; 284 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 285 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 286 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 287 for (d = 0; d < dim; ++d) { 288 x1p[d] = 0.0; 289 x2p[d] = 0.0; 290 for (e = 0; e < dim; ++e) { 291 x1p[d] += R[d*dim+e]*x1[e]; 292 x2p[d] += R[d*dim+e]*x2[e]; 293 } 294 } 295 if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 296 if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 297 /* 2) Project to (x, y) */ 298 coords[0] = 0.0; 299 coords[1] = 0.0; 300 coords[2] = x1p[0]; 301 coords[3] = x1p[1]; 302 coords[4] = x2p[0]; 303 coords[5] = x2p[1]; 304 /* Output R^T which rotates \hat z to the input normal */ 305 for (d = 0; d < dim; ++d) { 306 for (e = d+1; e < dim; ++e) { 307 PetscReal tmp; 308 309 tmp = R[d*dim+e]; 310 R[d*dim+e] = R[e*dim+d]; 311 R[e*dim+d] = tmp; 312 } 313 } 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "Invert2D_Internal" 319 PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 320 { 321 const PetscReal invDet = 1.0/detJ; 322 323 invJ[0] = invDet*J[3]; 324 invJ[1] = -invDet*J[1]; 325 invJ[2] = -invDet*J[2]; 326 invJ[3] = invDet*J[0]; 327 PetscLogFlops(5.0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "Invert3D_Internal" 332 PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 333 { 334 const PetscReal invDet = 1.0/detJ; 335 336 invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 337 invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 338 invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 339 invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 340 invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 341 invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 342 invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 343 invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 344 invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 345 PetscLogFlops(37.0); 346 } 347 348 #undef __FUNCT__ 349 #define __FUNCT__ "Det2D_Internal" 350 PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[]) 351 { 352 *detJ = J[0]*J[3] - J[1]*J[2]; 353 PetscLogFlops(3.0); 354 } 355 356 #undef __FUNCT__ 357 #define __FUNCT__ "Det3D_Internal" 358 PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[]) 359 { 360 *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 361 J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 362 J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 363 PetscLogFlops(12.0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "Volume_Triangle_Internal" 368 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 369 { 370 /* Signed volume is 1/2 the determinant 371 372 | 1 1 1 | 373 | x0 x1 x2 | 374 | y0 y1 y2 | 375 376 but if x0,y0 is the origin, we have 377 378 | x1 x2 | 379 | y1 y2 | 380 */ 381 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 382 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 383 PetscReal M[4], detM; 384 M[0] = x1; M[1] = x2; 385 M[2] = y1; M[3] = y2; 386 Det2D_Internal(&detM, M); 387 *vol = 0.5*detM; 388 PetscLogFlops(5.0); 389 } 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "Volume_Triangle_Origin_Internal" 393 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 394 { 395 Det2D_Internal(vol, coords); 396 *vol *= 0.5; 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "Volume_Tetrahedron_Internal" 401 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 402 { 403 /* Signed volume is 1/6th of the determinant 404 405 | 1 1 1 1 | 406 | x0 x1 x2 x3 | 407 | y0 y1 y2 y3 | 408 | z0 z1 z2 z3 | 409 410 but if x0,y0,z0 is the origin, we have 411 412 | x1 x2 x3 | 413 | y1 y2 y3 | 414 | z1 z2 z3 | 415 */ 416 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 417 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 418 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 419 PetscReal M[9], detM; 420 M[0] = x1; M[1] = x2; M[2] = x3; 421 M[3] = y1; M[4] = y2; M[5] = y3; 422 M[6] = z1; M[7] = z2; M[8] = z3; 423 Det3D_Internal(&detM, M); 424 *vol = -0.16666666666666666666666*detM; 425 PetscLogFlops(10.0); 426 } 427 428 #undef __FUNCT__ 429 #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 430 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 431 { 432 Det3D_Internal(vol, coords); 433 *vol *= -0.16666666666666666666666; 434 } 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 438 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 439 { 440 PetscSection coordSection; 441 Vec coordinates; 442 PetscScalar *coords; 443 PetscInt numCoords, d; 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 448 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 449 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 450 *detJ = 0.0; 451 if (numCoords == 4) { 452 const PetscInt dim = 2; 453 PetscReal R[4], J0; 454 455 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 456 ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 457 if (J) { 458 J0 = 0.5*PetscRealPart(coords[1]); 459 J[0] = R[0]*J0; J[1] = R[1]; 460 J[2] = R[2]*J0; J[3] = R[3]; 461 Det2D_Internal(detJ, J); 462 } 463 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 464 } else if (numCoords == 2) { 465 const PetscInt dim = 1; 466 467 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 468 if (J) { 469 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 470 *detJ = J[0]; 471 PetscLogFlops(2.0); 472 } 473 if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 474 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords); 475 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 481 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 482 { 483 PetscSection coordSection; 484 Vec coordinates; 485 PetscScalar *coords; 486 PetscInt numCoords, d, f, g; 487 PetscErrorCode ierr; 488 489 PetscFunctionBegin; 490 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 491 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 492 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 493 *detJ = 0.0; 494 if (numCoords == 9) { 495 const PetscInt dim = 3; 496 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 497 498 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 499 ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr); 500 if (J) { 501 const PetscInt pdim = 2; 502 503 for (d = 0; d < pdim; d++) { 504 for (f = 0; f < pdim; f++) { 505 J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 506 } 507 } 508 PetscLogFlops(8.0); 509 Det3D_Internal(detJ, J0); 510 for (d = 0; d < dim; d++) { 511 for (f = 0; f < dim; f++) { 512 J[d*dim+f] = 0.0; 513 for (g = 0; g < dim; g++) { 514 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 515 } 516 } 517 } 518 PetscLogFlops(18.0); 519 } 520 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 521 } else if (numCoords == 6) { 522 const PetscInt dim = 2; 523 524 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 525 if (J) { 526 for (d = 0; d < dim; d++) { 527 for (f = 0; f < dim; f++) { 528 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 529 } 530 } 531 PetscLogFlops(8.0); 532 Det2D_Internal(detJ, J); 533 } 534 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 535 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords); 536 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 542 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 543 { 544 PetscSection coordSection; 545 Vec coordinates; 546 PetscScalar *coords; 547 const PetscInt dim = 2; 548 PetscInt d, f; 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 553 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 554 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 555 *detJ = 0.0; 556 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 557 if (J) { 558 for (d = 0; d < dim; d++) { 559 for (f = 0; f < dim; f++) { 560 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 561 } 562 } 563 PetscLogFlops(8.0); 564 Det2D_Internal(detJ, J); 565 } 566 if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 567 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 568 PetscFunctionReturn(0); 569 } 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 573 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 574 { 575 PetscSection coordSection; 576 Vec coordinates; 577 PetscScalar *coords; 578 const PetscInt dim = 3; 579 PetscInt d, f; 580 PetscErrorCode ierr; 581 582 PetscFunctionBegin; 583 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 584 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 585 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 586 *detJ = 0.0; 587 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 588 if (J) { 589 for (d = 0; d < dim; d++) { 590 for (f = 0; f < dim; f++) { 591 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 592 } 593 } 594 PetscLogFlops(18.0); 595 Det3D_Internal(detJ, J); 596 *detJ = -*detJ; 597 } 598 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 599 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 605 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 606 { 607 PetscSection coordSection; 608 Vec coordinates; 609 PetscScalar *coords; 610 const PetscInt dim = 3; 611 PetscInt d; 612 PetscErrorCode ierr; 613 614 PetscFunctionBegin; 615 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 616 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 617 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 618 *detJ = 0.0; 619 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 620 if (J) { 621 for (d = 0; d < dim; d++) { 622 J[d*dim+0] = 0.5*(PetscRealPart(coords[(2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 623 J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 624 J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 625 } 626 PetscLogFlops(18.0); 627 Det3D_Internal(detJ, J); 628 *detJ = -*detJ; 629 } 630 if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 631 *detJ *= -8.0; 632 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 633 PetscFunctionReturn(0); 634 } 635 636 #undef __FUNCT__ 637 #define __FUNCT__ "DMPlexComputeCellGeometry" 638 /*@C 639 DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 640 641 Collective on DM 642 643 Input Arguments: 644 + dm - the DM 645 - cell - the cell 646 647 Output Arguments: 648 + v0 - the translation part of this affine transform 649 . J - the Jacobian of the transform from the reference element 650 . invJ - the inverse of the Jacobian 651 - detJ - the Jacobian determinant 652 653 Level: advanced 654 655 Fortran Notes: 656 Since it returns arrays, this routine is only available in Fortran 90, and you must 657 include petsc.h90 in your code. 658 659 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 660 @*/ 661 PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 662 { 663 PetscInt depth, dim, coneSize; 664 PetscErrorCode ierr; 665 666 PetscFunctionBegin; 667 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 668 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 669 if (depth == 1) { 670 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 671 switch (dim) { 672 case 1: 673 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 674 break; 675 case 2: 676 switch (coneSize) { 677 case 3: 678 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 679 break; 680 case 4: 681 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 682 break; 683 default: 684 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 685 } 686 break; 687 case 3: 688 switch (coneSize) { 689 case 4: 690 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 691 break; 692 case 8: 693 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 694 break; 695 default: 696 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 697 } 698 break; 699 default: 700 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 701 } 702 } else { 703 /* We need to keep a pointer to the depth label */ 704 ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr); 705 /* Cone size is now the number of faces */ 706 switch (dim) { 707 case 1: 708 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 709 break; 710 case 2: 711 switch (coneSize) { 712 case 3: 713 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 714 break; 715 case 4: 716 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 717 break; 718 default: 719 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 720 } 721 break; 722 case 3: 723 switch (coneSize) { 724 case 4: 725 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 726 break; 727 case 6: 728 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 729 break; 730 default: 731 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 732 } 733 break; 734 default: 735 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 736 } 737 } 738 PetscFunctionReturn(0); 739 } 740 741 #undef __FUNCT__ 742 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 743 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 744 { 745 PetscSection coordSection; 746 Vec coordinates; 747 PetscScalar *coords; 748 PetscInt coordSize; 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 753 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 754 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 755 if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 756 if (centroid) { 757 centroid[0] = 0.5*PetscRealPart(coords[0] + coords[dim+0]); 758 centroid[1] = 0.5*PetscRealPart(coords[1] + coords[dim+1]); 759 } 760 if (normal) { 761 normal[0] = PetscRealPart(coords[1] - coords[dim+1]); 762 normal[1] = -PetscRealPart(coords[0] - coords[dim+0]); 763 } 764 if (vol) { 765 *vol = sqrt(PetscSqr(PetscRealPart(coords[0] - coords[dim+0])) + PetscSqr(PetscRealPart(coords[1] - coords[dim+1]))); 766 } 767 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 773 /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 774 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 775 { 776 PetscSection coordSection; 777 Vec coordinates; 778 PetscScalar *coords = NULL; 779 PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 780 PetscInt tdim = 2, coordSize, numCorners, p, d, e; 781 PetscErrorCode ierr; 782 783 PetscFunctionBegin; 784 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 785 ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 786 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 787 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 788 dim = coordSize/numCorners; 789 if (normal) { 790 if (dim > 2) { 791 const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 792 const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 793 const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 794 PetscReal norm; 795 796 v0[0] = PetscRealPart(coords[0]); 797 v0[1] = PetscRealPart(coords[1]); 798 v0[2] = PetscRealPart(coords[2]); 799 normal[0] = y0*z1 - z0*y1; 800 normal[1] = z0*x1 - x0*z1; 801 normal[2] = x0*y1 - y0*x1; 802 norm = sqrt(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 803 normal[0] /= norm; 804 normal[1] /= norm; 805 normal[2] /= norm; 806 } else { 807 for (d = 0; d < dim; ++d) normal[d] = 0.0; 808 } 809 } 810 if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr);} 811 for (p = 0; p < numCorners; ++p) { 812 /* Need to do this copy to get types right */ 813 for (d = 0; d < tdim; ++d) { 814 ctmp[d] = PetscRealPart(coords[p*tdim+d]); 815 ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 816 } 817 Volume_Triangle_Origin_Internal(&vtmp, ctmp); 818 vsum += vtmp; 819 for (d = 0; d < tdim; ++d) { 820 csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 821 } 822 } 823 for (d = 0; d < tdim; ++d) { 824 csum[d] /= (tdim+1)*vsum; 825 } 826 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 827 if (vol) *vol = PetscAbsScalar(vsum); 828 if (centroid) { 829 if (dim > 2) { 830 for (d = 0; d < dim; ++d) { 831 centroid[d] = v0[d]; 832 for (e = 0; e < dim; ++e) { 833 centroid[d] += R[d*dim+e]*csum[e]; 834 } 835 } 836 } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 837 } 838 PetscFunctionReturn(0); 839 } 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 843 /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 844 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 845 { 846 PetscSection coordSection; 847 Vec coordinates; 848 PetscScalar *coords = NULL; 849 PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 850 const PetscInt *faces; 851 PetscInt numFaces, f, coordSize, numCorners, p, d; 852 PetscErrorCode ierr; 853 854 PetscFunctionBegin; 855 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 856 ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 857 858 ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 859 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 860 for (f = 0; f < numFaces; ++f) { 861 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 862 numCorners = coordSize/dim; 863 switch (numCorners) { 864 case 3: 865 for (d = 0; d < dim; ++d) { 866 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 867 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 868 coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 869 } 870 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 871 vsum += vtmp; 872 if (centroid) { 873 for (d = 0; d < dim; ++d) { 874 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 875 } 876 } 877 break; 878 case 4: 879 /* DO FOR PYRAMID */ 880 /* First tet */ 881 for (d = 0; d < dim; ++d) { 882 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 883 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 884 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 885 } 886 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 887 vsum += vtmp; 888 if (centroid) { 889 for (d = 0; d < dim; ++d) { 890 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 891 } 892 } 893 /* Second tet */ 894 for (d = 0; d < dim; ++d) { 895 coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 896 coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 897 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 898 } 899 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 900 vsum += vtmp; 901 if (centroid) { 902 for (d = 0; d < dim; ++d) { 903 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 904 } 905 } 906 break; 907 default: 908 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %d vertices", numCorners); 909 } 910 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 911 } 912 if (vol) *vol = PetscAbsScalar(vsum); 913 if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 914 PetscFunctionReturn(0); 915 } 916 917 #undef __FUNCT__ 918 #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 919 /*@C 920 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 921 922 Collective on DM 923 924 Input Arguments: 925 + dm - the DM 926 - cell - the cell 927 928 Output Arguments: 929 + volume - the cell volume 930 . centroid - the cell centroid 931 - normal - the cell normal, if appropriate 932 933 Level: advanced 934 935 Fortran Notes: 936 Since it returns arrays, this routine is only available in Fortran 90, and you must 937 include petsc.h90 in your code. 938 939 .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 940 @*/ 941 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 942 { 943 PetscInt depth, dim; 944 PetscErrorCode ierr; 945 946 PetscFunctionBegin; 947 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 948 ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 949 if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 950 /* We need to keep a pointer to the depth label */ 951 ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 952 /* Cone size is now the number of faces */ 953 switch (depth) { 954 case 1: 955 ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 956 break; 957 case 2: 958 ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 959 break; 960 case 3: 961 ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 962 break; 963 default: 964 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 965 } 966 PetscFunctionReturn(0); 967 } 968