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 = DMPlexComputeCellGeometryFEM(dm, c, NULL, 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 = NULL; 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 = DMGetCoordinateSection(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 = DMPlexComputeCellGeometryFEM(dm, c, NULL, 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, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 89 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 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 = DMGetCoordinateSection(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__ "PetscGridHashInitialize_Internal" 135 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 136 { 137 PetscInt d; 138 139 PetscFunctionBegin; 140 box->dim = dim; 141 for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 142 PetscFunctionReturn(0); 143 } 144 145 #undef __FUNCT__ 146 #define __FUNCT__ "PetscGridHashCreate" 147 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 148 { 149 PetscErrorCode ierr; 150 151 PetscFunctionBegin; 152 ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 153 ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 154 PetscFunctionReturn(0); 155 } 156 157 #undef __FUNCT__ 158 #define __FUNCT__ "PetscGridHashEnlarge" 159 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 160 { 161 PetscInt d; 162 163 PetscFunctionBegin; 164 for (d = 0; d < box->dim; ++d) { 165 box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 166 box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 167 } 168 PetscFunctionReturn(0); 169 } 170 171 #undef __FUNCT__ 172 #define __FUNCT__ "PetscGridHashSetGrid" 173 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 174 { 175 PetscInt d; 176 177 PetscFunctionBegin; 178 for (d = 0; d < box->dim; ++d) { 179 box->extent[d] = box->upper[d] - box->lower[d]; 180 if (n[d] == PETSC_DETERMINE) { 181 box->h[d] = h[d]; 182 box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 183 } else { 184 box->n[d] = n[d]; 185 box->h[d] = box->extent[d]/n[d]; 186 } 187 } 188 PetscFunctionReturn(0); 189 } 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "PetscGridHashGetEnclosingBox" 193 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscReal points[], PetscInt dboxes[], PetscInt boxes[]) 194 { 195 const PetscReal *lower = box->lower; 196 const PetscReal *upper = box->upper; 197 const PetscReal *h = box->h; 198 const PetscInt *n = box->n; 199 const PetscInt dim = box->dim; 200 PetscInt d, p; 201 202 PetscFunctionBegin; 203 for (p = 0; p < numPoints; ++p) { 204 for (d = 0; d < dim; ++d) { 205 PetscInt dbox = PetscFloorReal((points[p*dim+d] - lower[d])/h[d]); 206 207 if (dbox == n[d] && PetscAbsReal(points[p*dim+d] - upper[d]) < 1.0e-9) dbox = n[d]-1; 208 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", 209 p, points[p*dim+0], dim > 1 ? points[p*dim+1] : 0.0, dim > 2 ? points[p*dim+2] : 0.0); 210 dboxes[p*dim+d] = dbox; 211 } 212 if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 213 } 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "PetscGridHashDestroy" 219 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 220 { 221 PetscErrorCode ierr; 222 223 PetscFunctionBegin; 224 if (*box) { 225 ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 226 ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 227 ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 228 } 229 ierr = PetscFree(*box);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 233 #define __FUNCT__ "DMLocatePoints_Plex" 234 /* 235 Need to implement using the guess 236 */ 237 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS) 238 { 239 PetscInt cell = -1 /*, guess = -1*/; 240 PetscInt bs, numPoints, p; 241 PetscInt dim, cStart, cEnd, cMax, c, coneSize; 242 PetscInt *cells; 243 PetscScalar *a; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 248 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 249 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 250 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 251 ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 252 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 253 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 254 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); 255 numPoints /= bs; 256 ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 257 for (p = 0; p < numPoints; ++p) { 258 const PetscScalar *point = &a[p*bs]; 259 260 switch (dim) { 261 case 2: 262 for (c = cStart; c < cEnd; ++c) { 263 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 264 switch (coneSize) { 265 case 3: 266 ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 267 break; 268 case 4: 269 ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 270 break; 271 default: 272 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 273 } 274 if (cell >= 0) break; 275 } 276 break; 277 case 3: 278 for (c = cStart; c < cEnd; ++c) { 279 ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 280 switch (coneSize) { 281 case 4: 282 ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 283 break; 284 case 6: 285 ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 286 break; 287 default: 288 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 289 } 290 if (cell >= 0) break; 291 } 292 break; 293 default: 294 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 295 } 296 cells[p] = cell; 297 } 298 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 299 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr); 300 PetscFunctionReturn(0); 301 } 302 303 #undef __FUNCT__ 304 #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 305 /* 306 DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 307 */ 308 PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 309 { 310 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 311 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 312 const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 313 314 PetscFunctionBegin; 315 R[0] = c; R[1] = -s; 316 R[2] = s; R[3] = c; 317 coords[0] = 0.0; 318 coords[1] = r; 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal" 324 /* 325 DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D 326 327 This uses the basis completion described by Frisvad, 328 329 http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html 330 DOI:10.1080/2165347X.2012.689606 331 */ 332 PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 333 { 334 PetscReal x = PetscRealPart(coords[3] - coords[0]); 335 PetscReal y = PetscRealPart(coords[4] - coords[1]); 336 PetscReal z = PetscRealPart(coords[5] - coords[2]); 337 PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 338 PetscReal rinv = 1. / r; 339 PetscFunctionBegin; 340 341 x *= rinv; y *= rinv; z *= rinv; 342 if (x > 0.) { 343 PetscReal inv1pX = 1./ (1. + x); 344 345 R[0] = x; R[1] = -y; R[2] = -z; 346 R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 347 R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 348 } 349 else { 350 PetscReal inv1mX = 1./ (1. - x); 351 352 R[0] = x; R[1] = z; R[2] = y; 353 R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 354 R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 355 } 356 coords[0] = 0.0; 357 coords[1] = r; 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 363 /* 364 DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 365 */ 366 PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 367 { 368 PetscReal x1[3], x2[3], n[3], norm; 369 PetscReal x1p[3], x2p[3], xnp[3]; 370 PetscReal sqrtz, alpha; 371 const PetscInt dim = 3; 372 PetscInt d, e, p; 373 374 PetscFunctionBegin; 375 /* 0) Calculate normal vector */ 376 for (d = 0; d < dim; ++d) { 377 x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 378 x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 379 } 380 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 381 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 382 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 383 norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 384 n[0] /= norm; 385 n[1] /= norm; 386 n[2] /= norm; 387 /* 1) Take the normal vector and rotate until it is \hat z 388 389 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 390 391 R = / alpha nx nz alpha ny nz -1/alpha \ 392 | -alpha ny alpha nx 0 | 393 \ nx ny nz / 394 395 will rotate the normal vector to \hat z 396 */ 397 sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 398 /* Check for n = z */ 399 if (sqrtz < 1.0e-10) { 400 if (n[2] < 0.0) { 401 if (coordSize > 9) { 402 coords[2] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]); 403 coords[3] = PetscRealPart(coords[3*dim+1] - coords[0*dim+1]); 404 coords[4] = x2[0]; 405 coords[5] = x2[1]; 406 coords[6] = x1[0]; 407 coords[7] = x1[1]; 408 } else { 409 coords[2] = x2[0]; 410 coords[3] = x2[1]; 411 coords[4] = x1[0]; 412 coords[5] = x1[1]; 413 } 414 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 415 R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 416 R[6] = 0.0; R[7] = 0.0; R[8] = -1.0; 417 } else { 418 for (p = 3; p < coordSize/3; ++p) { 419 coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 420 coords[p*2+1] = PetscRealPart(coords[p*dim+1] - coords[0*dim+1]); 421 } 422 coords[2] = x1[0]; 423 coords[3] = x1[1]; 424 coords[4] = x2[0]; 425 coords[5] = x2[1]; 426 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 427 R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 428 R[6] = 0.0; R[7] = 0.0; R[8] = 1.0; 429 } 430 coords[0] = 0.0; 431 coords[1] = 0.0; 432 PetscFunctionReturn(0); 433 } 434 alpha = 1.0/sqrtz; 435 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 436 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 437 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 438 for (d = 0; d < dim; ++d) { 439 x1p[d] = 0.0; 440 x2p[d] = 0.0; 441 for (e = 0; e < dim; ++e) { 442 x1p[d] += R[d*dim+e]*x1[e]; 443 x2p[d] += R[d*dim+e]*x2[e]; 444 } 445 } 446 if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 447 if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 448 /* 2) Project to (x, y) */ 449 for (p = 3; p < coordSize/3; ++p) { 450 for (d = 0; d < dim; ++d) { 451 xnp[d] = 0.0; 452 for (e = 0; e < dim; ++e) { 453 xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 454 } 455 if (d < dim-1) coords[p*2+d] = xnp[d]; 456 } 457 } 458 coords[0] = 0.0; 459 coords[1] = 0.0; 460 coords[2] = x1p[0]; 461 coords[3] = x1p[1]; 462 coords[4] = x2p[0]; 463 coords[5] = x2p[1]; 464 /* Output R^T which rotates \hat z to the input normal */ 465 for (d = 0; d < dim; ++d) { 466 for (e = d+1; e < dim; ++e) { 467 PetscReal tmp; 468 469 tmp = R[d*dim+e]; 470 R[d*dim+e] = R[e*dim+d]; 471 R[e*dim+d] = tmp; 472 } 473 } 474 PetscFunctionReturn(0); 475 } 476 477 #undef __FUNCT__ 478 #define __FUNCT__ "Volume_Triangle_Internal" 479 PETSC_UNUSED 480 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 481 { 482 /* Signed volume is 1/2 the determinant 483 484 | 1 1 1 | 485 | x0 x1 x2 | 486 | y0 y1 y2 | 487 488 but if x0,y0 is the origin, we have 489 490 | x1 x2 | 491 | y1 y2 | 492 */ 493 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 494 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 495 PetscReal M[4], detM; 496 M[0] = x1; M[1] = x2; 497 M[2] = y1; M[3] = y2; 498 DMPlex_Det2D_Internal(&detM, M); 499 *vol = 0.5*detM; 500 PetscLogFlops(5.0); 501 } 502 503 #undef __FUNCT__ 504 #define __FUNCT__ "Volume_Triangle_Origin_Internal" 505 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 506 { 507 DMPlex_Det2D_Internal(vol, coords); 508 *vol *= 0.5; 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "Volume_Tetrahedron_Internal" 513 PETSC_UNUSED 514 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 515 { 516 /* Signed volume is 1/6th of the determinant 517 518 | 1 1 1 1 | 519 | x0 x1 x2 x3 | 520 | y0 y1 y2 y3 | 521 | z0 z1 z2 z3 | 522 523 but if x0,y0,z0 is the origin, we have 524 525 | x1 x2 x3 | 526 | y1 y2 y3 | 527 | z1 z2 z3 | 528 */ 529 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 530 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 531 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 532 PetscReal M[9], detM; 533 M[0] = x1; M[1] = x2; M[2] = x3; 534 M[3] = y1; M[4] = y2; M[5] = y3; 535 M[6] = z1; M[7] = z2; M[8] = z3; 536 DMPlex_Det3D_Internal(&detM, M); 537 *vol = -0.16666666666666666666666*detM; 538 PetscLogFlops(10.0); 539 } 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 543 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 544 { 545 DMPlex_Det3D_Internal(vol, coords); 546 *vol *= -0.16666666666666666666666; 547 } 548 549 #undef __FUNCT__ 550 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 551 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 552 { 553 PetscSection coordSection; 554 Vec coordinates; 555 PetscScalar *coords = NULL; 556 PetscInt numCoords, d; 557 PetscErrorCode ierr; 558 559 PetscFunctionBegin; 560 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 561 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 562 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 563 *detJ = 0.0; 564 if (numCoords == 6) { 565 const PetscInt dim = 3; 566 PetscReal R[9], J0; 567 568 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 569 ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr); 570 if (J) { 571 J0 = 0.5*PetscRealPart(coords[1]); 572 J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 573 J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 574 J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 575 DMPlex_Det3D_Internal(detJ, J); 576 } 577 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 578 } else if (numCoords == 4) { 579 const PetscInt dim = 2; 580 PetscReal R[4], J0; 581 582 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 583 ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 584 if (J) { 585 J0 = 0.5*PetscRealPart(coords[1]); 586 J[0] = R[0]*J0; J[1] = R[1]; 587 J[2] = R[2]*J0; J[3] = R[3]; 588 DMPlex_Det2D_Internal(detJ, J); 589 } 590 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 591 } else if (numCoords == 2) { 592 const PetscInt dim = 1; 593 594 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 595 if (J) { 596 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 597 *detJ = J[0]; 598 PetscLogFlops(2.0); 599 } 600 if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 601 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 602 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 603 PetscFunctionReturn(0); 604 } 605 606 #undef __FUNCT__ 607 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 608 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 609 { 610 PetscSection coordSection; 611 Vec coordinates; 612 PetscScalar *coords = NULL; 613 PetscInt numCoords, d, f, g; 614 PetscErrorCode ierr; 615 616 PetscFunctionBegin; 617 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 618 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 619 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 620 *detJ = 0.0; 621 if (numCoords == 9) { 622 const PetscInt dim = 3; 623 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 624 625 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 626 ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 627 if (J) { 628 const PetscInt pdim = 2; 629 630 for (d = 0; d < pdim; d++) { 631 for (f = 0; f < pdim; f++) { 632 J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 633 } 634 } 635 PetscLogFlops(8.0); 636 DMPlex_Det3D_Internal(detJ, J0); 637 for (d = 0; d < dim; d++) { 638 for (f = 0; f < dim; f++) { 639 J[d*dim+f] = 0.0; 640 for (g = 0; g < dim; g++) { 641 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 642 } 643 } 644 } 645 PetscLogFlops(18.0); 646 } 647 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 648 } else if (numCoords == 6) { 649 const PetscInt dim = 2; 650 651 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 652 if (J) { 653 for (d = 0; d < dim; d++) { 654 for (f = 0; f < dim; f++) { 655 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 656 } 657 } 658 PetscLogFlops(8.0); 659 DMPlex_Det2D_Internal(detJ, J); 660 } 661 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 662 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 663 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 664 PetscFunctionReturn(0); 665 } 666 667 #undef __FUNCT__ 668 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 669 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 670 { 671 PetscSection coordSection; 672 Vec coordinates; 673 PetscScalar *coords = NULL; 674 PetscInt numCoords, d, f, g; 675 PetscErrorCode ierr; 676 677 PetscFunctionBegin; 678 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 679 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 680 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 681 *detJ = 0.0; 682 if (numCoords == 12) { 683 const PetscInt dim = 3; 684 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 685 686 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 687 ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 688 if (J) { 689 const PetscInt pdim = 2; 690 691 for (d = 0; d < pdim; d++) { 692 J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 693 J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 694 } 695 PetscLogFlops(8.0); 696 DMPlex_Det3D_Internal(detJ, J0); 697 for (d = 0; d < dim; d++) { 698 for (f = 0; f < dim; f++) { 699 J[d*dim+f] = 0.0; 700 for (g = 0; g < dim; g++) { 701 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 702 } 703 } 704 } 705 PetscLogFlops(18.0); 706 } 707 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 708 } else if ((numCoords == 8) || (numCoords == 16)) { 709 const PetscInt dim = 2; 710 711 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 712 if (J) { 713 for (d = 0; d < dim; d++) { 714 J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 715 J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 716 } 717 PetscLogFlops(8.0); 718 DMPlex_Det2D_Internal(detJ, J); 719 } 720 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 721 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 722 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 723 PetscFunctionReturn(0); 724 } 725 726 #undef __FUNCT__ 727 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 728 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 729 { 730 PetscSection coordSection; 731 Vec coordinates; 732 PetscScalar *coords = NULL; 733 const PetscInt dim = 3; 734 PetscInt d; 735 PetscErrorCode ierr; 736 737 PetscFunctionBegin; 738 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 739 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 740 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 741 *detJ = 0.0; 742 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 743 if (J) { 744 for (d = 0; d < dim; d++) { 745 /* I orient with outward face normals */ 746 J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 747 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 748 J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 749 } 750 PetscLogFlops(18.0); 751 DMPlex_Det3D_Internal(detJ, J); 752 } 753 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 754 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 755 PetscFunctionReturn(0); 756 } 757 758 #undef __FUNCT__ 759 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 760 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 761 { 762 PetscSection coordSection; 763 Vec coordinates; 764 PetscScalar *coords = NULL; 765 const PetscInt dim = 3; 766 PetscInt d; 767 PetscErrorCode ierr; 768 769 PetscFunctionBegin; 770 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 771 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 772 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 773 *detJ = 0.0; 774 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 775 if (J) { 776 for (d = 0; d < dim; d++) { 777 J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 778 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 779 J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 780 } 781 PetscLogFlops(18.0); 782 DMPlex_Det3D_Internal(detJ, J); 783 } 784 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 785 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 786 PetscFunctionReturn(0); 787 } 788 789 #undef __FUNCT__ 790 #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 791 /*@C 792 DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 793 794 Collective on DM 795 796 Input Arguments: 797 + dm - the DM 798 - cell - the cell 799 800 Output Arguments: 801 + v0 - the translation part of this affine transform 802 . J - the Jacobian of the transform from the reference element 803 . invJ - the inverse of the Jacobian 804 - detJ - the Jacobian determinant 805 806 Level: advanced 807 808 Fortran Notes: 809 Since it returns arrays, this routine is only available in Fortran 90, and you must 810 include petsc.h90 in your code. 811 812 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 813 @*/ 814 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 815 { 816 PetscInt depth, dim, coneSize; 817 PetscErrorCode ierr; 818 819 PetscFunctionBegin; 820 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 821 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 822 if (depth == 1) { 823 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 824 } else { 825 DMLabel depth; 826 827 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 828 ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 829 } 830 switch (dim) { 831 case 1: 832 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 833 break; 834 case 2: 835 switch (coneSize) { 836 case 3: 837 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 838 break; 839 case 4: 840 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 841 break; 842 default: 843 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 844 } 845 break; 846 case 3: 847 switch (coneSize) { 848 case 4: 849 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 850 break; 851 case 6: /* Faces */ 852 case 8: /* Vertices */ 853 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 854 break; 855 default: 856 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 857 } 858 break; 859 default: 860 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 861 } 862 PetscFunctionReturn(0); 863 } 864 865 #undef __FUNCT__ 866 #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 867 static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 868 { 869 PetscQuadrature quad; 870 PetscSection coordSection; 871 Vec coordinates; 872 PetscScalar *coords = NULL; 873 const PetscReal *quadPoints; 874 PetscReal *basisDer; 875 PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 876 PetscErrorCode ierr; 877 878 PetscFunctionBegin; 879 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 880 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 881 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 882 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 883 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 884 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 885 ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 886 ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 887 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 888 *detJ = 0.0; 889 if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 890 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); 891 if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 892 if (J) { 893 for (q = 0; q < Nq; ++q) { 894 PetscInt i, j, k, c, r; 895 896 /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 897 for (k = 0; k < pdim; ++k) 898 for (j = 0; j < dim; ++j) 899 for (i = 0; i < cdim; ++i) 900 J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 901 PetscLogFlops(2.0*pdim*dim*cdim); 902 if (cdim > dim) { 903 for (c = dim; c < cdim; ++c) 904 for (r = 0; r < cdim; ++r) 905 J[r*cdim+c] = r == c ? 1.0 : 0.0; 906 } 907 switch (cdim) { 908 case 3: 909 DMPlex_Det3D_Internal(detJ, J); 910 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 911 break; 912 case 2: 913 DMPlex_Det2D_Internal(detJ, J); 914 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 915 break; 916 case 1: 917 *detJ = J[0]; 918 if (invJ) invJ[0] = 1.0/J[0]; 919 } 920 } 921 } 922 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 923 PetscFunctionReturn(0); 924 } 925 926 #undef __FUNCT__ 927 #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 928 /*@C 929 DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 930 931 Collective on DM 932 933 Input Arguments: 934 + dm - the DM 935 . cell - the cell 936 - fe - the finite element containing the quadrature 937 938 Output Arguments: 939 + v0 - the translation part of this transform 940 . J - the Jacobian of the transform from the reference element at each quadrature point 941 . invJ - the inverse of the Jacobian at each quadrature point 942 - detJ - the Jacobian determinant at each quadrature point 943 944 Level: advanced 945 946 Fortran Notes: 947 Since it returns arrays, this routine is only available in Fortran 90, and you must 948 include petsc.h90 in your code. 949 950 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 951 @*/ 952 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 953 { 954 PetscErrorCode ierr; 955 956 PetscFunctionBegin; 957 if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 958 else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 959 PetscFunctionReturn(0); 960 } 961 962 #undef __FUNCT__ 963 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 964 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 965 { 966 PetscSection coordSection; 967 Vec coordinates; 968 PetscScalar *coords = NULL; 969 PetscScalar tmp[2]; 970 PetscInt coordSize; 971 PetscErrorCode ierr; 972 973 PetscFunctionBegin; 974 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 975 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 976 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 977 if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 978 ierr = DMPlexLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 979 if (centroid) { 980 centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 981 centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 982 } 983 if (normal) { 984 PetscReal norm; 985 986 normal[0] = -PetscRealPart(coords[1] - tmp[1]); 987 normal[1] = PetscRealPart(coords[0] - tmp[0]); 988 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 989 normal[0] /= norm; 990 normal[1] /= norm; 991 } 992 if (vol) { 993 *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 994 } 995 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 996 PetscFunctionReturn(0); 997 } 998 999 #undef __FUNCT__ 1000 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 1001 /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1002 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1003 { 1004 PetscSection coordSection; 1005 Vec coordinates; 1006 PetscScalar *coords = NULL; 1007 PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 1008 PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1013 ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 1014 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1015 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1016 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1017 if (normal) { 1018 if (dim > 2) { 1019 const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 1020 const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 1021 const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 1022 PetscReal norm; 1023 1024 v0[0] = PetscRealPart(coords[0]); 1025 v0[1] = PetscRealPart(coords[1]); 1026 v0[2] = PetscRealPart(coords[2]); 1027 normal[0] = y0*z1 - z0*y1; 1028 normal[1] = z0*x1 - x0*z1; 1029 normal[2] = x0*y1 - y0*x1; 1030 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 1031 normal[0] /= norm; 1032 normal[1] /= norm; 1033 normal[2] /= norm; 1034 } else { 1035 for (d = 0; d < dim; ++d) normal[d] = 0.0; 1036 } 1037 } 1038 if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 1039 for (p = 0; p < numCorners; ++p) { 1040 /* Need to do this copy to get types right */ 1041 for (d = 0; d < tdim; ++d) { 1042 ctmp[d] = PetscRealPart(coords[p*tdim+d]); 1043 ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 1044 } 1045 Volume_Triangle_Origin_Internal(&vtmp, ctmp); 1046 vsum += vtmp; 1047 for (d = 0; d < tdim; ++d) { 1048 csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 1049 } 1050 } 1051 for (d = 0; d < tdim; ++d) { 1052 csum[d] /= (tdim+1)*vsum; 1053 } 1054 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1055 if (vol) *vol = PetscAbsReal(vsum); 1056 if (centroid) { 1057 if (dim > 2) { 1058 for (d = 0; d < dim; ++d) { 1059 centroid[d] = v0[d]; 1060 for (e = 0; e < dim; ++e) { 1061 centroid[d] += R[d*dim+e]*csum[e]; 1062 } 1063 } 1064 } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 1065 } 1066 PetscFunctionReturn(0); 1067 } 1068 1069 #undef __FUNCT__ 1070 #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 1071 /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1072 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1073 { 1074 PetscSection coordSection; 1075 Vec coordinates; 1076 PetscScalar *coords = NULL; 1077 PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1078 const PetscInt *faces, *facesO; 1079 PetscInt numFaces, f, coordSize, numCorners, p, d; 1080 PetscErrorCode ierr; 1081 1082 PetscFunctionBegin; 1083 if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 1084 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1085 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1086 1087 if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 1088 ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 1089 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1090 ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 1091 for (f = 0; f < numFaces; ++f) { 1092 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1093 numCorners = coordSize/dim; 1094 switch (numCorners) { 1095 case 3: 1096 for (d = 0; d < dim; ++d) { 1097 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1098 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1099 coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 1100 } 1101 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1102 if (facesO[f] < 0) vtmp = -vtmp; 1103 vsum += vtmp; 1104 if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 1105 for (d = 0; d < dim; ++d) { 1106 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1107 } 1108 } 1109 break; 1110 case 4: 1111 /* DO FOR PYRAMID */ 1112 /* First tet */ 1113 for (d = 0; d < dim; ++d) { 1114 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1115 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1116 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1117 } 1118 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1119 if (facesO[f] < 0) vtmp = -vtmp; 1120 vsum += vtmp; 1121 if (centroid) { 1122 for (d = 0; d < dim; ++d) { 1123 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1124 } 1125 } 1126 /* Second tet */ 1127 for (d = 0; d < dim; ++d) { 1128 coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 1129 coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 1130 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1131 } 1132 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1133 if (facesO[f] < 0) vtmp = -vtmp; 1134 vsum += vtmp; 1135 if (centroid) { 1136 for (d = 0; d < dim; ++d) { 1137 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1138 } 1139 } 1140 break; 1141 default: 1142 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 1143 } 1144 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1145 } 1146 if (vol) *vol = PetscAbsReal(vsum); 1147 if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1148 if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 1149 PetscFunctionReturn(0); 1150 } 1151 1152 #undef __FUNCT__ 1153 #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1154 /*@C 1155 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1156 1157 Collective on DM 1158 1159 Input Arguments: 1160 + dm - the DM 1161 - cell - the cell 1162 1163 Output Arguments: 1164 + volume - the cell volume 1165 . centroid - the cell centroid 1166 - normal - the cell normal, if appropriate 1167 1168 Level: advanced 1169 1170 Fortran Notes: 1171 Since it returns arrays, this routine is only available in Fortran 90, and you must 1172 include petsc.h90 in your code. 1173 1174 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1175 @*/ 1176 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1177 { 1178 PetscInt depth, dim; 1179 PetscErrorCode ierr; 1180 1181 PetscFunctionBegin; 1182 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1183 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1184 if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1185 /* We need to keep a pointer to the depth label */ 1186 ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1187 /* Cone size is now the number of faces */ 1188 switch (depth) { 1189 case 1: 1190 ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1191 break; 1192 case 2: 1193 ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1194 break; 1195 case 3: 1196 ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1197 break; 1198 default: 1199 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1200 } 1201 PetscFunctionReturn(0); 1202 } 1203 1204 #undef __FUNCT__ 1205 #define __FUNCT__ "DMPlexComputeGeometryFEM" 1206 /* This should also take a PetscFE argument I think */ 1207 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1208 { 1209 DM dmCell; 1210 Vec coordinates; 1211 PetscSection coordSection, sectionCell; 1212 PetscScalar *cgeom; 1213 PetscInt cStart, cEnd, cMax, c; 1214 PetscErrorCode ierr; 1215 1216 PetscFunctionBegin; 1217 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1218 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1219 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1220 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1221 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1222 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1223 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1224 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1225 cEnd = cMax < 0 ? cEnd : cMax; 1226 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1227 /* TODO This needs to be multiplied by Nq for non-affine */ 1228 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1229 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1230 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1231 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1232 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1233 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1234 for (c = cStart; c < cEnd; ++c) { 1235 PetscFECellGeom *cg; 1236 1237 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1238 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1239 ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1240 if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1241 } 1242 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1243 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1244 PetscFunctionReturn(0); 1245 } 1246 1247 #undef __FUNCT__ 1248 #define __FUNCT__ "DMPlexComputeGeometryFVM" 1249 /*@ 1250 DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1251 1252 Input Parameter: 1253 . dm - The DM 1254 1255 Output Parameters: 1256 + cellgeom - A Vec of PetscFVCellGeom data 1257 . facegeom - A Vec of PetscFVFaceGeom data 1258 1259 Level: developer 1260 1261 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1262 @*/ 1263 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1264 { 1265 DM dmFace, dmCell; 1266 DMLabel ghostLabel; 1267 PetscSection sectionFace, sectionCell; 1268 PetscSection coordSection; 1269 Vec coordinates; 1270 PetscScalar *fgeom, *cgeom; 1271 PetscReal minradius, gminradius; 1272 PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1273 PetscErrorCode ierr; 1274 1275 PetscFunctionBegin; 1276 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1277 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1278 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1279 /* Make cell centroids and volumes */ 1280 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1281 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1282 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1283 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1284 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1285 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1286 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1287 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1288 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1289 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1290 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1291 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1292 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1293 for (c = cStart; c < cEndInterior; ++c) { 1294 PetscFVCellGeom *cg; 1295 1296 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1297 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1298 ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1299 } 1300 /* Compute face normals and minimum cell radius */ 1301 ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1302 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1303 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1304 ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 1305 for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1306 ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1307 ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1308 ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1309 ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1310 ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1311 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1312 minradius = PETSC_MAX_REAL; 1313 for (f = fStart; f < fEnd; ++f) { 1314 PetscFVFaceGeom *fg; 1315 PetscReal area; 1316 PetscInt ghost = -1, d; 1317 1318 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1319 if (ghost >= 0) continue; 1320 ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1321 ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1322 for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1323 /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1324 { 1325 PetscFVCellGeom *cL, *cR; 1326 const PetscInt *cells; 1327 PetscReal *lcentroid, *rcentroid; 1328 PetscReal l[3], r[3], v[3]; 1329 1330 ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 1331 ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1332 ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1333 lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 1334 rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 1335 ierr = DMPlexLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 1336 ierr = DMPlexLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 1337 DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1338 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1339 for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1340 } 1341 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1342 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]); 1343 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]); 1344 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1345 } 1346 if (cells[0] < cEndInterior) { 1347 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1348 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1349 } 1350 if (cells[1] < cEndInterior) { 1351 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1352 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1353 } 1354 } 1355 } 1356 ierr = MPI_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1357 ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1358 /* Compute centroids of ghost cells */ 1359 for (c = cEndInterior; c < cEnd; ++c) { 1360 PetscFVFaceGeom *fg; 1361 const PetscInt *cone, *support; 1362 PetscInt coneSize, supportSize, s; 1363 1364 ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1365 if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1366 ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1367 ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 1368 if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize); 1369 ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1370 ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1371 for (s = 0; s < 2; ++s) { 1372 /* Reflect ghost centroid across plane of face */ 1373 if (support[s] == c) { 1374 const PetscFVCellGeom *ci; 1375 PetscFVCellGeom *cg; 1376 PetscReal c2f[3], a; 1377 1378 ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1379 DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1380 a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1381 ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1382 DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1383 cg->volume = ci->volume; 1384 } 1385 } 1386 } 1387 ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1388 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1389 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1390 ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1391 PetscFunctionReturn(0); 1392 } 1393 1394 #undef __FUNCT__ 1395 #define __FUNCT__ "DMPlexGetMinRadius" 1396 /*@C 1397 DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1398 1399 Not collective 1400 1401 Input Argument: 1402 . dm - the DM 1403 1404 Output Argument: 1405 . minradius - the minium cell radius 1406 1407 Level: developer 1408 1409 .seealso: DMGetCoordinates() 1410 @*/ 1411 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1412 { 1413 PetscFunctionBegin; 1414 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1415 PetscValidPointer(minradius,2); 1416 *minradius = ((DM_Plex*) dm->data)->minradius; 1417 PetscFunctionReturn(0); 1418 } 1419 1420 #undef __FUNCT__ 1421 #define __FUNCT__ "DMPlexSetMinRadius" 1422 /*@C 1423 DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1424 1425 Logically collective 1426 1427 Input Arguments: 1428 + dm - the DM 1429 - minradius - the minium cell radius 1430 1431 Level: developer 1432 1433 .seealso: DMSetCoordinates() 1434 @*/ 1435 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1436 { 1437 PetscFunctionBegin; 1438 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1439 ((DM_Plex*) dm->data)->minradius = minradius; 1440 PetscFunctionReturn(0); 1441 } 1442 1443 #undef __FUNCT__ 1444 #define __FUNCT__ "BuildGradientReconstruction_Internal" 1445 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1446 { 1447 DMLabel ghostLabel; 1448 PetscScalar *dx, *grad, **gref; 1449 PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1450 PetscErrorCode ierr; 1451 1452 PetscFunctionBegin; 1453 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1454 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1455 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1456 ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1457 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1458 ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1459 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1460 for (c = cStart; c < cEndInterior; c++) { 1461 const PetscInt *faces; 1462 PetscInt numFaces, usedFaces, f, d; 1463 const PetscFVCellGeom *cg; 1464 PetscBool boundary; 1465 PetscInt ghost; 1466 1467 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1468 ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1469 ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1470 if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1471 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1472 const PetscFVCellGeom *cg1; 1473 PetscFVFaceGeom *fg; 1474 const PetscInt *fcells; 1475 PetscInt ncell, side; 1476 1477 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1478 ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1479 if ((ghost >= 0) || boundary) continue; 1480 ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1481 side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1482 ncell = fcells[!side]; /* the neighbor */ 1483 ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1484 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1485 for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1486 gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1487 } 1488 if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1489 ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1490 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1491 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1492 ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1493 if ((ghost >= 0) || boundary) continue; 1494 for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1495 ++usedFaces; 1496 } 1497 } 1498 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1499 PetscFunctionReturn(0); 1500 } 1501 1502 #undef __FUNCT__ 1503 #define __FUNCT__ "DMPlexComputeGradientFVM" 1504 /*@ 1505 DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 1506 1507 Collective on DM 1508 1509 Input Arguments: 1510 + dm - The DM 1511 . fvm - The PetscFV 1512 . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 1513 - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 1514 1515 Output Parameters: 1516 + faceGeometry - The geometric factors for gradient calculation are inserted 1517 - dmGrad - The DM describing the layout of gradient data 1518 1519 Level: developer 1520 1521 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 1522 @*/ 1523 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 1524 { 1525 DM dmFace, dmCell; 1526 PetscScalar *fgeom, *cgeom; 1527 PetscSection sectionGrad; 1528 PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 1529 PetscErrorCode ierr; 1530 1531 PetscFunctionBegin; 1532 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1533 ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1534 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1535 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1536 /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 1537 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1538 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1539 ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1540 ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1541 ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1542 ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1543 ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1544 /* Create storage for gradients */ 1545 ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 1546 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 1547 ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 1548 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 1549 ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 1550 ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 1551 ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 1552 PetscFunctionReturn(0); 1553 } 1554