1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMPlexGetLineIntersection_2D_Internal" 5 static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 6 { 7 const PetscReal p0_x = segmentA[0*2+0]; 8 const PetscReal p0_y = segmentA[0*2+1]; 9 const PetscReal p1_x = segmentA[1*2+0]; 10 const PetscReal p1_y = segmentA[1*2+1]; 11 const PetscReal p2_x = segmentB[0*2+0]; 12 const PetscReal p2_y = segmentB[0*2+1]; 13 const PetscReal p3_x = segmentB[1*2+0]; 14 const PetscReal p3_y = segmentB[1*2+1]; 15 const PetscReal s1_x = p1_x - p0_x; 16 const PetscReal s1_y = p1_y - p0_y; 17 const PetscReal s2_x = p3_x - p2_x; 18 const PetscReal s2_y = p3_y - p2_y; 19 const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 20 21 PetscFunctionBegin; 22 *hasIntersection = PETSC_FALSE; 23 /* Non-parallel lines */ 24 if (denom != 0.0) { 25 const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 26 const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 27 28 if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 29 *hasIntersection = PETSC_TRUE; 30 if (intersection) { 31 intersection[0] = p0_x + (t * s1_x); 32 intersection[1] = p0_y + (t * s1_y); 33 } 34 } 35 } 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal" 41 static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 42 { 43 const PetscInt embedDim = 2; 44 const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 45 PetscReal x = PetscRealPart(point[0]); 46 PetscReal y = PetscRealPart(point[1]); 47 PetscReal v0[2], J[4], invJ[4], detJ; 48 PetscReal xi, eta; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 53 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 54 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 55 56 if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c; 57 else *cell = -1; 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal" 63 static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 64 { 65 PetscSection coordSection; 66 Vec coordsLocal; 67 PetscScalar *coords = NULL; 68 const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 69 PetscReal x = PetscRealPart(point[0]); 70 PetscReal y = PetscRealPart(point[1]); 71 PetscInt crossings = 0, f; 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 76 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 77 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 78 for (f = 0; f < 4; ++f) { 79 PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 80 PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 81 PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 82 PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 83 PetscReal slope = (y_j - y_i) / (x_j - x_i); 84 PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 85 PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 86 PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 87 if ((cond1 || cond2) && above) ++crossings; 88 } 89 if (crossings % 2) *cell = c; 90 else *cell = -1; 91 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal" 97 static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 98 { 99 const PetscInt embedDim = 3; 100 PetscReal v0[3], J[9], invJ[9], detJ; 101 PetscReal x = PetscRealPart(point[0]); 102 PetscReal y = PetscRealPart(point[1]); 103 PetscReal z = PetscRealPart(point[2]); 104 PetscReal xi, eta, zeta; 105 PetscErrorCode ierr; 106 107 PetscFunctionBegin; 108 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 109 xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 110 eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 111 zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 112 113 if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 114 else *cell = -1; 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal" 120 static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 121 { 122 PetscSection coordSection; 123 Vec coordsLocal; 124 PetscScalar *coords; 125 const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 126 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 127 PetscBool found = PETSC_TRUE; 128 PetscInt f; 129 PetscErrorCode ierr; 130 131 PetscFunctionBegin; 132 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 133 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 134 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 135 for (f = 0; f < 6; ++f) { 136 /* Check the point is under plane */ 137 /* Get face normal */ 138 PetscReal v_i[3]; 139 PetscReal v_j[3]; 140 PetscReal normal[3]; 141 PetscReal pp[3]; 142 PetscReal dot; 143 144 v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 145 v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 146 v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 147 v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 148 v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 149 v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 150 normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 151 normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 152 normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 153 pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 154 pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 155 pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 156 dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 157 158 /* Check that projected point is in face (2D location problem) */ 159 if (dot < 0.0) { 160 found = PETSC_FALSE; 161 break; 162 } 163 } 164 if (found) *cell = c; 165 else *cell = -1; 166 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 167 PetscFunctionReturn(0); 168 } 169 170 #undef __FUNCT__ 171 #define __FUNCT__ "PetscGridHashInitialize_Internal" 172 static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 173 { 174 PetscInt d; 175 176 PetscFunctionBegin; 177 box->dim = dim; 178 for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 179 PetscFunctionReturn(0); 180 } 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "PetscGridHashCreate" 184 PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 185 { 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 190 ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 194 #undef __FUNCT__ 195 #define __FUNCT__ "PetscGridHashEnlarge" 196 PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 197 { 198 PetscInt d; 199 200 PetscFunctionBegin; 201 for (d = 0; d < box->dim; ++d) { 202 box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 203 box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 204 } 205 PetscFunctionReturn(0); 206 } 207 208 #undef __FUNCT__ 209 #define __FUNCT__ "PetscGridHashSetGrid" 210 PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 211 { 212 PetscInt d; 213 214 PetscFunctionBegin; 215 for (d = 0; d < box->dim; ++d) { 216 box->extent[d] = box->upper[d] - box->lower[d]; 217 if (n[d] == PETSC_DETERMINE) { 218 box->h[d] = h[d]; 219 box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 220 } else { 221 box->n[d] = n[d]; 222 box->h[d] = box->extent[d]/n[d]; 223 } 224 } 225 PetscFunctionReturn(0); 226 } 227 228 #undef __FUNCT__ 229 #define __FUNCT__ "PetscGridHashGetEnclosingBox" 230 PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 231 { 232 const PetscReal *lower = box->lower; 233 const PetscReal *upper = box->upper; 234 const PetscReal *h = box->h; 235 const PetscInt *n = box->n; 236 const PetscInt dim = box->dim; 237 PetscInt d, p; 238 239 PetscFunctionBegin; 240 for (p = 0; p < numPoints; ++p) { 241 for (d = 0; d < dim; ++d) { 242 PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 243 244 if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 245 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", 246 p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0); 247 dboxes[p*dim+d] = dbox; 248 } 249 if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 250 } 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "PetscGridHashDestroy" 256 PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 257 { 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 if (*box) { 262 ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 263 ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 264 ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 265 } 266 ierr = PetscFree(*box);CHKERRQ(ierr); 267 PetscFunctionReturn(0); 268 } 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "DMPlexLocatePoint_Internal" 272 PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 273 { 274 PetscInt coneSize; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 switch (dim) { 279 case 2: 280 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 281 switch (coneSize) { 282 case 3: 283 ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 284 break; 285 case 4: 286 ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 287 break; 288 default: 289 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 290 } 291 break; 292 case 3: 293 ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 294 switch (coneSize) { 295 case 4: 296 ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 297 break; 298 case 6: 299 ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 300 break; 301 default: 302 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 303 } 304 break; 305 default: 306 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 307 } 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "DMPlexComputeGridHash_Internal" 313 PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 314 { 315 MPI_Comm comm; 316 PetscGridHash lbox; 317 Vec coordinates; 318 PetscSection coordSection; 319 Vec coordsLocal; 320 const PetscScalar *coords; 321 PetscInt *dboxes, *boxes; 322 PetscInt n[3] = {10, 10, 10}; 323 PetscInt dim, N, cStart, cEnd, cMax, c, i; 324 PetscErrorCode ierr; 325 326 PetscFunctionBegin; 327 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 328 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 329 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 330 ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 331 ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 332 ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr); 333 for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);} 334 ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 335 ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr); 336 #if 0 337 /* Could define a custom reduction to merge these */ 338 ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 339 ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr); 340 #endif 341 /* Is there a reason to snap the local bounding box to a division of the global box? */ 342 /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */ 343 /* Create label */ 344 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 345 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 346 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 347 ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 348 ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 349 /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 350 ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 351 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 352 ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr); 353 for (c = cStart; c < cEnd; ++c) { 354 const PetscReal *h = lbox->h; 355 PetscScalar *ccoords = NULL; 356 PetscInt csize = 0; 357 PetscScalar point[3]; 358 PetscInt dlim[6], d, e, i, j, k; 359 360 /* Find boxes enclosing each vertex */ 361 ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr); 362 ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr); 363 /* Mark cells containing the vertices */ 364 for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 365 /* Get grid of boxes containing these */ 366 for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 367 for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 368 for (e = 1; e < dim+1; ++e) { 369 for (d = 0; d < dim; ++d) { 370 dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 371 dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 372 } 373 } 374 /* Check for intersection of box with cell */ 375 for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) { 376 for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) { 377 for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) { 378 const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 379 PetscScalar cpoint[3]; 380 PetscInt cell, edge, ii, jj, kk; 381 382 /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 383 for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 384 for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 385 for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 386 387 ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 388 if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 389 } 390 } 391 } 392 /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 393 for (edge = 0; edge < dim+1; ++edge) { 394 PetscReal segA[6], segB[6]; 395 396 for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);} 397 for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 398 if (dim > 2) {segB[2] = PetscRealPart(point[2]); 399 segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 400 for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 401 if (dim > 1) {segB[1] = PetscRealPart(point[1]); 402 segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 403 for (ii = 0; ii < 2; ++ii) { 404 PetscBool intersects; 405 406 segB[0] = PetscRealPart(point[0]); 407 segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 408 ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 409 if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 410 } 411 } 412 } 413 } 414 } 415 } 416 } 417 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 418 } 419 ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 420 ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 421 ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 422 *localBox = lbox; 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "DMLocatePoints_Plex" 428 PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, PetscSF cellSF) 429 { 430 DM_Plex *mesh = (DM_Plex *) dm->data; 431 PetscBool hash = mesh->useHashLocation; 432 PetscInt bs, numPoints, p, numFound, *found = NULL; 433 PetscInt dim, cStart, cEnd, cMax, numCells, c; 434 const PetscInt *boxCells; 435 PetscSFNode *cells; 436 PetscScalar *a; 437 PetscMPIInt result; 438 PetscErrorCode ierr; 439 440 PetscFunctionBegin; 441 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 442 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 443 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr); 444 if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 445 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); 446 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 447 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 448 if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 449 ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 450 ierr = VecGetArray(v, &a);CHKERRQ(ierr); 451 numPoints /= bs; 452 ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 453 if (hash) { 454 if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 455 /* Designate the local box for each point */ 456 /* Send points to correct process */ 457 /* Search cells that lie in each subbox */ 458 /* Should we bin points before doing search? */ 459 ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 460 } 461 for (p = 0, numFound = 0; p < numPoints; ++p) { 462 const PetscScalar *point = &a[p*bs]; 463 PetscInt dbin[3], bin, cell = -1, cellOffset; 464 465 cells[p].rank = -1; 466 cells[p].index = -1; 467 if (hash) { 468 ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 469 /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 470 ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 471 ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 472 for (c = cellOffset; c < cellOffset + numCells; ++c) { 473 ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 474 if (cell >= 0) { 475 cells[p].rank = 0; 476 cells[p].index = cell; 477 numFound++; 478 break; 479 } 480 } 481 } else { 482 for (c = cStart; c < cEnd; ++c) { 483 ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 484 if (cell >= 0) { 485 cells[p].rank = 0; 486 cells[p].index = cell; 487 numFound++; 488 break; 489 } 490 } 491 } 492 } 493 if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);} 494 /* Check for highest numbered proc that claims a point (do we care?) */ 495 ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 496 if (numFound < numPoints) { 497 ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr); 498 for (p = 0, numFound = 0; p < numPoints; p++) { 499 if (cells[p].rank >= 0 && cells[p].index >= 0) { 500 if (numFound < p) { 501 cells[numFound] = cells[p]; 502 } 503 found[numFound++] = p; 504 } 505 } 506 } 507 ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 508 PetscFunctionReturn(0); 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 513 /* 514 DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 515 */ 516 PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 517 { 518 const PetscReal x = PetscRealPart(coords[2] - coords[0]); 519 const PetscReal y = PetscRealPart(coords[3] - coords[1]); 520 const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 521 522 PetscFunctionBegin; 523 R[0] = c; R[1] = -s; 524 R[2] = s; R[3] = c; 525 coords[0] = 0.0; 526 coords[1] = r; 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal" 532 /* 533 DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D 534 535 This uses the basis completion described by Frisvad, 536 537 http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html 538 DOI:10.1080/2165347X.2012.689606 539 */ 540 PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 541 { 542 PetscReal x = PetscRealPart(coords[3] - coords[0]); 543 PetscReal y = PetscRealPart(coords[4] - coords[1]); 544 PetscReal z = PetscRealPart(coords[5] - coords[2]); 545 PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 546 PetscReal rinv = 1. / r; 547 PetscFunctionBegin; 548 549 x *= rinv; y *= rinv; z *= rinv; 550 if (x > 0.) { 551 PetscReal inv1pX = 1./ (1. + x); 552 553 R[0] = x; R[1] = -y; R[2] = -z; 554 R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 555 R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 556 } 557 else { 558 PetscReal inv1mX = 1./ (1. - x); 559 560 R[0] = x; R[1] = z; R[2] = y; 561 R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 562 R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 563 } 564 coords[0] = 0.0; 565 coords[1] = r; 566 PetscFunctionReturn(0); 567 } 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 571 /* 572 DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 573 */ 574 PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 575 { 576 PetscReal x1[3], x2[3], n[3], norm; 577 PetscReal x1p[3], x2p[3], xnp[3]; 578 PetscReal sqrtz, alpha; 579 const PetscInt dim = 3; 580 PetscInt d, e, p; 581 582 PetscFunctionBegin; 583 /* 0) Calculate normal vector */ 584 for (d = 0; d < dim; ++d) { 585 x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 586 x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 587 } 588 n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 589 n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 590 n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 591 norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 592 n[0] /= norm; 593 n[1] /= norm; 594 n[2] /= norm; 595 /* 1) Take the normal vector and rotate until it is \hat z 596 597 Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 598 599 R = / alpha nx nz alpha ny nz -1/alpha \ 600 | -alpha ny alpha nx 0 | 601 \ nx ny nz / 602 603 will rotate the normal vector to \hat z 604 */ 605 sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 606 /* Check for n = z */ 607 if (sqrtz < 1.0e-10) { 608 const PetscInt s = PetscSign(n[2]); 609 /* If nz < 0, rotate 180 degrees around x-axis */ 610 for (p = 3; p < coordSize/3; ++p) { 611 coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 612 coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s; 613 } 614 coords[0] = 0.0; 615 coords[1] = 0.0; 616 coords[2] = x1[0]; 617 coords[3] = x1[1] * s; 618 coords[4] = x2[0]; 619 coords[5] = x2[1] * s; 620 R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 621 R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0; 622 R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s; 623 PetscFunctionReturn(0); 624 } 625 alpha = 1.0/sqrtz; 626 R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 627 R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 628 R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 629 for (d = 0; d < dim; ++d) { 630 x1p[d] = 0.0; 631 x2p[d] = 0.0; 632 for (e = 0; e < dim; ++e) { 633 x1p[d] += R[d*dim+e]*x1[e]; 634 x2p[d] += R[d*dim+e]*x2[e]; 635 } 636 } 637 if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 638 if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 639 /* 2) Project to (x, y) */ 640 for (p = 3; p < coordSize/3; ++p) { 641 for (d = 0; d < dim; ++d) { 642 xnp[d] = 0.0; 643 for (e = 0; e < dim; ++e) { 644 xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 645 } 646 if (d < dim-1) coords[p*2+d] = xnp[d]; 647 } 648 } 649 coords[0] = 0.0; 650 coords[1] = 0.0; 651 coords[2] = x1p[0]; 652 coords[3] = x1p[1]; 653 coords[4] = x2p[0]; 654 coords[5] = x2p[1]; 655 /* Output R^T which rotates \hat z to the input normal */ 656 for (d = 0; d < dim; ++d) { 657 for (e = d+1; e < dim; ++e) { 658 PetscReal tmp; 659 660 tmp = R[d*dim+e]; 661 R[d*dim+e] = R[e*dim+d]; 662 R[e*dim+d] = tmp; 663 } 664 } 665 PetscFunctionReturn(0); 666 } 667 668 #undef __FUNCT__ 669 #define __FUNCT__ "Volume_Triangle_Internal" 670 PETSC_UNUSED 671 PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 672 { 673 /* Signed volume is 1/2 the determinant 674 675 | 1 1 1 | 676 | x0 x1 x2 | 677 | y0 y1 y2 | 678 679 but if x0,y0 is the origin, we have 680 681 | x1 x2 | 682 | y1 y2 | 683 */ 684 const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 685 const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 686 PetscReal M[4], detM; 687 M[0] = x1; M[1] = x2; 688 M[2] = y1; M[3] = y2; 689 DMPlex_Det2D_Internal(&detM, M); 690 *vol = 0.5*detM; 691 (void)PetscLogFlops(5.0); 692 } 693 694 #undef __FUNCT__ 695 #define __FUNCT__ "Volume_Triangle_Origin_Internal" 696 PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 697 { 698 DMPlex_Det2D_Internal(vol, coords); 699 *vol *= 0.5; 700 } 701 702 #undef __FUNCT__ 703 #define __FUNCT__ "Volume_Tetrahedron_Internal" 704 PETSC_UNUSED 705 PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 706 { 707 /* Signed volume is 1/6th of the determinant 708 709 | 1 1 1 1 | 710 | x0 x1 x2 x3 | 711 | y0 y1 y2 y3 | 712 | z0 z1 z2 z3 | 713 714 but if x0,y0,z0 is the origin, we have 715 716 | x1 x2 x3 | 717 | y1 y2 y3 | 718 | z1 z2 z3 | 719 */ 720 const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 721 const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 722 const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 723 PetscReal M[9], detM; 724 M[0] = x1; M[1] = x2; M[2] = x3; 725 M[3] = y1; M[4] = y2; M[5] = y3; 726 M[6] = z1; M[7] = z2; M[8] = z3; 727 DMPlex_Det3D_Internal(&detM, M); 728 *vol = -0.16666666666666666666666*detM; 729 (void)PetscLogFlops(10.0); 730 } 731 732 #undef __FUNCT__ 733 #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 734 PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 735 { 736 DMPlex_Det3D_Internal(vol, coords); 737 *vol *= -0.16666666666666666666666; 738 } 739 740 #undef __FUNCT__ 741 #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 742 static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 743 { 744 PetscSection coordSection; 745 Vec coordinates; 746 PetscScalar *coords = NULL; 747 PetscInt numCoords, d; 748 PetscErrorCode ierr; 749 750 PetscFunctionBegin; 751 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 752 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 753 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 754 if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 755 *detJ = 0.0; 756 if (numCoords == 6) { 757 const PetscInt dim = 3; 758 PetscReal R[9], J0; 759 760 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 761 ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr); 762 if (J) { 763 J0 = 0.5*PetscRealPart(coords[1]); 764 J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 765 J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 766 J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 767 DMPlex_Det3D_Internal(detJ, J); 768 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 769 } 770 } else if (numCoords == 4) { 771 const PetscInt dim = 2; 772 PetscReal R[4], J0; 773 774 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 775 ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 776 if (J) { 777 J0 = 0.5*PetscRealPart(coords[1]); 778 J[0] = R[0]*J0; J[1] = R[1]; 779 J[2] = R[2]*J0; J[3] = R[3]; 780 DMPlex_Det2D_Internal(detJ, J); 781 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 782 } 783 } else if (numCoords == 2) { 784 const PetscInt dim = 1; 785 786 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 787 if (J) { 788 J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 789 *detJ = J[0]; 790 ierr = PetscLogFlops(2.0);CHKERRQ(ierr); 791 if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);} 792 } 793 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 794 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 795 PetscFunctionReturn(0); 796 } 797 798 #undef __FUNCT__ 799 #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 800 static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 801 { 802 PetscSection coordSection; 803 Vec coordinates; 804 PetscScalar *coords = NULL; 805 PetscInt numCoords, d, f, g; 806 PetscErrorCode ierr; 807 808 PetscFunctionBegin; 809 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 810 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 811 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 812 *detJ = 0.0; 813 if (numCoords == 9) { 814 const PetscInt dim = 3; 815 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 816 817 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 818 ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 819 if (J) { 820 const PetscInt pdim = 2; 821 822 for (d = 0; d < pdim; d++) { 823 for (f = 0; f < pdim; f++) { 824 J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 825 } 826 } 827 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 828 DMPlex_Det3D_Internal(detJ, J0); 829 for (d = 0; d < dim; d++) { 830 for (f = 0; f < dim; f++) { 831 J[d*dim+f] = 0.0; 832 for (g = 0; g < dim; g++) { 833 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 834 } 835 } 836 } 837 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 838 } 839 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 840 } else if (numCoords == 6) { 841 const PetscInt dim = 2; 842 843 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 844 if (J) { 845 for (d = 0; d < dim; d++) { 846 for (f = 0; f < dim; f++) { 847 J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 848 } 849 } 850 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 851 DMPlex_Det2D_Internal(detJ, J); 852 } 853 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 854 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 855 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 856 PetscFunctionReturn(0); 857 } 858 859 #undef __FUNCT__ 860 #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 861 static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 862 { 863 PetscSection coordSection; 864 Vec coordinates; 865 PetscScalar *coords = NULL; 866 PetscInt numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd; 867 PetscErrorCode ierr; 868 869 PetscFunctionBegin; 870 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 871 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 872 ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 873 if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 874 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 875 numCoords = numSelfCoords ? numSelfCoords : numCoords; 876 *detJ = 0.0; 877 if (numCoords == 12) { 878 const PetscInt dim = 3; 879 PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 880 881 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 882 ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 883 if (J) { 884 const PetscInt pdim = 2; 885 886 for (d = 0; d < pdim; d++) { 887 J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 888 J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 889 } 890 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 891 DMPlex_Det3D_Internal(detJ, J0); 892 for (d = 0; d < dim; d++) { 893 for (f = 0; f < dim; f++) { 894 J[d*dim+f] = 0.0; 895 for (g = 0; g < dim; g++) { 896 J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 897 } 898 } 899 } 900 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 901 } 902 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 903 } else if (numCoords == 8) { 904 const PetscInt dim = 2; 905 906 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 907 if (J) { 908 for (d = 0; d < dim; d++) { 909 J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 910 J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 911 } 912 ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 913 DMPlex_Det2D_Internal(detJ, J); 914 } 915 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 916 } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 917 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 918 PetscFunctionReturn(0); 919 } 920 921 #undef __FUNCT__ 922 #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 923 static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 924 { 925 PetscSection coordSection; 926 Vec coordinates; 927 PetscScalar *coords = NULL; 928 const PetscInt dim = 3; 929 PetscInt d; 930 PetscErrorCode ierr; 931 932 PetscFunctionBegin; 933 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 934 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 935 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 936 *detJ = 0.0; 937 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 938 if (J) { 939 for (d = 0; d < dim; d++) { 940 /* I orient with outward face normals */ 941 J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 942 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 943 J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 944 } 945 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 946 DMPlex_Det3D_Internal(detJ, J); 947 } 948 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 949 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 950 PetscFunctionReturn(0); 951 } 952 953 #undef __FUNCT__ 954 #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 955 static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 956 { 957 PetscSection coordSection; 958 Vec coordinates; 959 PetscScalar *coords = NULL; 960 const PetscInt dim = 3; 961 PetscInt d; 962 PetscErrorCode ierr; 963 964 PetscFunctionBegin; 965 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 966 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 967 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 968 *detJ = 0.0; 969 if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 970 if (J) { 971 for (d = 0; d < dim; d++) { 972 J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 973 J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 974 J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 975 } 976 ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 977 DMPlex_Det3D_Internal(detJ, J); 978 } 979 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 980 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 981 PetscFunctionReturn(0); 982 } 983 984 #undef __FUNCT__ 985 #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 986 /*@C 987 DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 988 989 Collective on DM 990 991 Input Arguments: 992 + dm - the DM 993 - cell - the cell 994 995 Output Arguments: 996 + v0 - the translation part of this affine transform 997 . J - the Jacobian of the transform from the reference element 998 . invJ - the inverse of the Jacobian 999 - detJ - the Jacobian determinant 1000 1001 Level: advanced 1002 1003 Fortran Notes: 1004 Since it returns arrays, this routine is only available in Fortran 90, and you must 1005 include petsc.h90 in your code. 1006 1007 .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 1008 @*/ 1009 PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1010 { 1011 PetscInt depth, dim, coneSize; 1012 PetscErrorCode ierr; 1013 1014 PetscFunctionBegin; 1015 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1016 ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1017 if (depth == 1) { 1018 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1019 } else { 1020 DMLabel depth; 1021 1022 ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 1023 ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 1024 } 1025 switch (dim) { 1026 case 1: 1027 ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1028 break; 1029 case 2: 1030 switch (coneSize) { 1031 case 3: 1032 ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1033 break; 1034 case 4: 1035 ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1036 break; 1037 default: 1038 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1039 } 1040 break; 1041 case 3: 1042 switch (coneSize) { 1043 case 4: 1044 ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1045 break; 1046 case 6: /* Faces */ 1047 case 8: /* Vertices */ 1048 ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1049 break; 1050 default: 1051 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1052 } 1053 break; 1054 default: 1055 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1056 } 1057 PetscFunctionReturn(0); 1058 } 1059 1060 #undef __FUNCT__ 1061 #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 1062 static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1063 { 1064 PetscQuadrature quad; 1065 PetscSection coordSection; 1066 Vec coordinates; 1067 PetscScalar *coords = NULL; 1068 const PetscReal *quadPoints; 1069 PetscReal *basisDer; 1070 PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 1071 PetscErrorCode ierr; 1072 1073 PetscFunctionBegin; 1074 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1075 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1076 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1077 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1078 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 1079 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1080 ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 1081 ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1082 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 1083 *detJ = 0.0; 1084 if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 1085 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); 1086 if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 1087 if (J) { 1088 ierr = PetscMemzero(J, Nq*cdim*dim*sizeof(PetscReal));CHKERRQ(ierr); 1089 for (q = 0; q < Nq; ++q) { 1090 PetscInt i, j, k, c, r; 1091 1092 /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 1093 for (k = 0; k < pdim; ++k) 1094 for (j = 0; j < dim; ++j) 1095 for (i = 0; i < cdim; ++i) 1096 J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 1097 ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr); 1098 if (cdim > dim) { 1099 for (c = dim; c < cdim; ++c) 1100 for (r = 0; r < cdim; ++r) 1101 J[r*cdim+c] = r == c ? 1.0 : 0.0; 1102 } 1103 switch (cdim) { 1104 case 3: 1105 DMPlex_Det3D_Internal(detJ, J); 1106 if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1107 break; 1108 case 2: 1109 DMPlex_Det2D_Internal(detJ, J); 1110 if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1111 break; 1112 case 1: 1113 *detJ = J[0]; 1114 if (invJ) invJ[0] = 1.0/J[0]; 1115 } 1116 } 1117 } 1118 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 1119 PetscFunctionReturn(0); 1120 } 1121 1122 #undef __FUNCT__ 1123 #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 1124 /*@C 1125 DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 1126 1127 Collective on DM 1128 1129 Input Arguments: 1130 + dm - the DM 1131 . cell - the cell 1132 - fe - the finite element containing the quadrature 1133 1134 Output Arguments: 1135 + v0 - the translation part of this transform 1136 . J - the Jacobian of the transform from the reference element at each quadrature point 1137 . invJ - the inverse of the Jacobian at each quadrature point 1138 - detJ - the Jacobian determinant at each quadrature point 1139 1140 Level: advanced 1141 1142 Fortran Notes: 1143 Since it returns arrays, this routine is only available in Fortran 90, and you must 1144 include petsc.h90 in your code. 1145 1146 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1147 @*/ 1148 PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1149 { 1150 PetscErrorCode ierr; 1151 1152 PetscFunctionBegin; 1153 if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1154 else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1155 PetscFunctionReturn(0); 1156 } 1157 1158 #undef __FUNCT__ 1159 #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 1160 static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1161 { 1162 PetscSection coordSection; 1163 Vec coordinates; 1164 PetscScalar *coords = NULL; 1165 PetscScalar tmp[2]; 1166 PetscInt coordSize; 1167 PetscErrorCode ierr; 1168 1169 PetscFunctionBegin; 1170 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1171 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1172 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1173 if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 1174 ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1175 if (centroid) { 1176 centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 1177 centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1178 } 1179 if (normal) { 1180 PetscReal norm; 1181 1182 normal[0] = -PetscRealPart(coords[1] - tmp[1]); 1183 normal[1] = PetscRealPart(coords[0] - tmp[0]); 1184 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1185 normal[0] /= norm; 1186 normal[1] /= norm; 1187 } 1188 if (vol) { 1189 *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1190 } 1191 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1192 PetscFunctionReturn(0); 1193 } 1194 1195 #undef __FUNCT__ 1196 #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 1197 /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1198 static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1199 { 1200 PetscSection coordSection; 1201 Vec coordinates; 1202 PetscScalar *coords = NULL; 1203 PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 1204 PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1209 ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 1210 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1211 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1212 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1213 if (dim > 2 && centroid) { 1214 v0[0] = PetscRealPart(coords[0]); 1215 v0[1] = PetscRealPart(coords[1]); 1216 v0[2] = PetscRealPart(coords[2]); 1217 } 1218 if (normal) { 1219 if (dim > 2) { 1220 const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 1221 const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 1222 const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 1223 PetscReal norm; 1224 1225 normal[0] = y0*z1 - z0*y1; 1226 normal[1] = z0*x1 - x0*z1; 1227 normal[2] = x0*y1 - y0*x1; 1228 norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 1229 normal[0] /= norm; 1230 normal[1] /= norm; 1231 normal[2] /= norm; 1232 } else { 1233 for (d = 0; d < dim; ++d) normal[d] = 0.0; 1234 } 1235 } 1236 if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 1237 for (p = 0; p < numCorners; ++p) { 1238 /* Need to do this copy to get types right */ 1239 for (d = 0; d < tdim; ++d) { 1240 ctmp[d] = PetscRealPart(coords[p*tdim+d]); 1241 ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 1242 } 1243 Volume_Triangle_Origin_Internal(&vtmp, ctmp); 1244 vsum += vtmp; 1245 for (d = 0; d < tdim; ++d) { 1246 csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 1247 } 1248 } 1249 for (d = 0; d < tdim; ++d) { 1250 csum[d] /= (tdim+1)*vsum; 1251 } 1252 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1253 if (vol) *vol = PetscAbsReal(vsum); 1254 if (centroid) { 1255 if (dim > 2) { 1256 for (d = 0; d < dim; ++d) { 1257 centroid[d] = v0[d]; 1258 for (e = 0; e < dim; ++e) { 1259 centroid[d] += R[d*dim+e]*csum[e]; 1260 } 1261 } 1262 } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 1263 } 1264 PetscFunctionReturn(0); 1265 } 1266 1267 #undef __FUNCT__ 1268 #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 1269 /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1270 static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1271 { 1272 PetscSection coordSection; 1273 Vec coordinates; 1274 PetscScalar *coords = NULL; 1275 PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1276 const PetscInt *faces, *facesO; 1277 PetscInt numFaces, f, coordSize, numCorners, p, d; 1278 PetscErrorCode ierr; 1279 1280 PetscFunctionBegin; 1281 if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 1282 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1283 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1284 1285 if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 1286 ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 1287 ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1288 ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 1289 for (f = 0; f < numFaces; ++f) { 1290 ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1291 numCorners = coordSize/dim; 1292 switch (numCorners) { 1293 case 3: 1294 for (d = 0; d < dim; ++d) { 1295 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1296 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1297 coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 1298 } 1299 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1300 if (facesO[f] < 0) vtmp = -vtmp; 1301 vsum += vtmp; 1302 if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 1303 for (d = 0; d < dim; ++d) { 1304 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1305 } 1306 } 1307 break; 1308 case 4: 1309 /* DO FOR PYRAMID */ 1310 /* First tet */ 1311 for (d = 0; d < dim; ++d) { 1312 coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 1313 coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 1314 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1315 } 1316 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1317 if (facesO[f] < 0) vtmp = -vtmp; 1318 vsum += vtmp; 1319 if (centroid) { 1320 for (d = 0; d < dim; ++d) { 1321 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1322 } 1323 } 1324 /* Second tet */ 1325 for (d = 0; d < dim; ++d) { 1326 coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 1327 coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 1328 coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 1329 } 1330 Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1331 if (facesO[f] < 0) vtmp = -vtmp; 1332 vsum += vtmp; 1333 if (centroid) { 1334 for (d = 0; d < dim; ++d) { 1335 for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 1336 } 1337 } 1338 break; 1339 default: 1340 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 1341 } 1342 ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 1343 } 1344 if (vol) *vol = PetscAbsReal(vsum); 1345 if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1346 if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 1347 PetscFunctionReturn(0); 1348 } 1349 1350 #undef __FUNCT__ 1351 #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1352 /*@C 1353 DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1354 1355 Collective on DM 1356 1357 Input Arguments: 1358 + dm - the DM 1359 - cell - the cell 1360 1361 Output Arguments: 1362 + volume - the cell volume 1363 . centroid - the cell centroid 1364 - normal - the cell normal, if appropriate 1365 1366 Level: advanced 1367 1368 Fortran Notes: 1369 Since it returns arrays, this routine is only available in Fortran 90, and you must 1370 include petsc.h90 in your code. 1371 1372 .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1373 @*/ 1374 PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1375 { 1376 PetscInt depth, dim; 1377 PetscErrorCode ierr; 1378 1379 PetscFunctionBegin; 1380 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1381 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1382 if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1383 /* We need to keep a pointer to the depth label */ 1384 ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1385 /* Cone size is now the number of faces */ 1386 switch (depth) { 1387 case 1: 1388 ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1389 break; 1390 case 2: 1391 ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1392 break; 1393 case 3: 1394 ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1395 break; 1396 default: 1397 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1398 } 1399 PetscFunctionReturn(0); 1400 } 1401 1402 #undef __FUNCT__ 1403 #define __FUNCT__ "DMPlexComputeGeometryFEM" 1404 /* This should also take a PetscFE argument I think */ 1405 PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1406 { 1407 DM dmCell; 1408 Vec coordinates; 1409 PetscSection coordSection, sectionCell; 1410 PetscScalar *cgeom; 1411 PetscInt cStart, cEnd, cMax, c; 1412 PetscErrorCode ierr; 1413 1414 PetscFunctionBegin; 1415 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1416 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1417 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1418 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1419 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1420 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1421 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1422 ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1423 cEnd = cMax < 0 ? cEnd : cMax; 1424 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1425 /* TODO This needs to be multiplied by Nq for non-affine */ 1426 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1427 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1428 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1429 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1430 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1431 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1432 for (c = cStart; c < cEnd; ++c) { 1433 PetscFECellGeom *cg; 1434 1435 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1436 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1437 ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1438 if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1439 } 1440 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1441 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1442 PetscFunctionReturn(0); 1443 } 1444 1445 #undef __FUNCT__ 1446 #define __FUNCT__ "DMPlexComputeGeometryFVM" 1447 /*@ 1448 DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1449 1450 Input Parameter: 1451 . dm - The DM 1452 1453 Output Parameters: 1454 + cellgeom - A Vec of PetscFVCellGeom data 1455 . facegeom - A Vec of PetscFVFaceGeom data 1456 1457 Level: developer 1458 1459 .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1460 @*/ 1461 PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1462 { 1463 DM dmFace, dmCell; 1464 DMLabel ghostLabel; 1465 PetscSection sectionFace, sectionCell; 1466 PetscSection coordSection; 1467 Vec coordinates; 1468 PetscScalar *fgeom, *cgeom; 1469 PetscReal minradius, gminradius; 1470 PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1471 PetscErrorCode ierr; 1472 1473 PetscFunctionBegin; 1474 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1475 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1476 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1477 /* Make cell centroids and volumes */ 1478 ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1479 ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1480 ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1481 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1482 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1483 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1484 ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1485 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1486 ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1487 ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1488 ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1489 ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1490 if (cEndInterior < 0) { 1491 cEndInterior = cEnd; 1492 } 1493 ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1494 for (c = cStart; c < cEndInterior; ++c) { 1495 PetscFVCellGeom *cg; 1496 1497 ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1498 ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1499 ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1500 } 1501 /* Compute face normals and minimum cell radius */ 1502 ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1503 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1504 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1505 ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 1506 for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1507 ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1508 ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1509 ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1510 ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1511 ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1512 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1513 minradius = PETSC_MAX_REAL; 1514 for (f = fStart; f < fEnd; ++f) { 1515 PetscFVFaceGeom *fg; 1516 PetscReal area; 1517 PetscInt ghost = -1, d, numChildren; 1518 1519 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1520 ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 1521 if (ghost >= 0 || numChildren) continue; 1522 ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1523 ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1524 for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1525 /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1526 { 1527 PetscFVCellGeom *cL, *cR; 1528 PetscInt ncells; 1529 const PetscInt *cells; 1530 PetscReal *lcentroid, *rcentroid; 1531 PetscReal l[3], r[3], v[3]; 1532 1533 ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 1534 ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr); 1535 ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1536 lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 1537 if (ncells > 1) { 1538 ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1539 rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 1540 } 1541 else { 1542 rcentroid = fg->centroid; 1543 } 1544 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 1545 ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 1546 DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1547 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1548 for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1549 } 1550 if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1551 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]); 1552 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]); 1553 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1554 } 1555 if (cells[0] < cEndInterior) { 1556 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1557 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1558 } 1559 if (ncells > 1 && cells[1] < cEndInterior) { 1560 DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1561 minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1562 } 1563 } 1564 } 1565 ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1566 ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1567 /* Compute centroids of ghost cells */ 1568 for (c = cEndInterior; c < cEnd; ++c) { 1569 PetscFVFaceGeom *fg; 1570 const PetscInt *cone, *support; 1571 PetscInt coneSize, supportSize, s; 1572 1573 ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1574 if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1575 ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1576 ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 1577 if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize); 1578 ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1579 ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1580 for (s = 0; s < 2; ++s) { 1581 /* Reflect ghost centroid across plane of face */ 1582 if (support[s] == c) { 1583 PetscFVCellGeom *ci; 1584 PetscFVCellGeom *cg; 1585 PetscReal c2f[3], a; 1586 1587 ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1588 DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1589 a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1590 ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1591 DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1592 cg->volume = ci->volume; 1593 } 1594 } 1595 } 1596 ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1597 ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1598 ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1599 ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1600 PetscFunctionReturn(0); 1601 } 1602 1603 #undef __FUNCT__ 1604 #define __FUNCT__ "DMPlexGetMinRadius" 1605 /*@C 1606 DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1607 1608 Not collective 1609 1610 Input Argument: 1611 . dm - the DM 1612 1613 Output Argument: 1614 . minradius - the minium cell radius 1615 1616 Level: developer 1617 1618 .seealso: DMGetCoordinates() 1619 @*/ 1620 PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1621 { 1622 PetscFunctionBegin; 1623 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1624 PetscValidPointer(minradius,2); 1625 *minradius = ((DM_Plex*) dm->data)->minradius; 1626 PetscFunctionReturn(0); 1627 } 1628 1629 #undef __FUNCT__ 1630 #define __FUNCT__ "DMPlexSetMinRadius" 1631 /*@C 1632 DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1633 1634 Logically collective 1635 1636 Input Arguments: 1637 + dm - the DM 1638 - minradius - the minium cell radius 1639 1640 Level: developer 1641 1642 .seealso: DMSetCoordinates() 1643 @*/ 1644 PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1645 { 1646 PetscFunctionBegin; 1647 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1648 ((DM_Plex*) dm->data)->minradius = minradius; 1649 PetscFunctionReturn(0); 1650 } 1651 1652 #undef __FUNCT__ 1653 #define __FUNCT__ "BuildGradientReconstruction_Internal" 1654 static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1655 { 1656 DMLabel ghostLabel; 1657 PetscScalar *dx, *grad, **gref; 1658 PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1659 PetscErrorCode ierr; 1660 1661 PetscFunctionBegin; 1662 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1663 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1664 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1665 ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1666 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1667 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1668 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1669 for (c = cStart; c < cEndInterior; c++) { 1670 const PetscInt *faces; 1671 PetscInt numFaces, usedFaces, f, d; 1672 PetscFVCellGeom *cg; 1673 PetscBool boundary; 1674 PetscInt ghost; 1675 1676 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1677 ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1678 ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1679 if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1680 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1681 PetscFVCellGeom *cg1; 1682 PetscFVFaceGeom *fg; 1683 const PetscInt *fcells; 1684 PetscInt ncell, side; 1685 1686 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1687 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1688 if ((ghost >= 0) || boundary) continue; 1689 ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1690 side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1691 ncell = fcells[!side]; /* the neighbor */ 1692 ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1693 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1694 for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1695 gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1696 } 1697 if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1698 ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1699 for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1700 ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1701 ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1702 if ((ghost >= 0) || boundary) continue; 1703 for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1704 ++usedFaces; 1705 } 1706 } 1707 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1708 PetscFunctionReturn(0); 1709 } 1710 1711 #undef __FUNCT__ 1712 #define __FUNCT__ "BuildGradientReconstruction_Internal_Tree" 1713 static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1714 { 1715 DMLabel ghostLabel; 1716 PetscScalar *dx, *grad, **gref; 1717 PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 1718 PetscSection neighSec; 1719 PetscInt (*neighbors)[2]; 1720 PetscInt *counter; 1721 PetscErrorCode ierr; 1722 1723 PetscFunctionBegin; 1724 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1725 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1726 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1727 if (cEndInterior < 0) { 1728 cEndInterior = cEnd; 1729 } 1730 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr); 1731 ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr); 1732 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1733 ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1734 for (f = fStart; f < fEnd; f++) { 1735 const PetscInt *fcells; 1736 PetscBool boundary; 1737 PetscInt ghost = -1; 1738 PetscInt numChildren, numCells, c; 1739 1740 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1741 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1742 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1743 if ((ghost >= 0) || boundary || numChildren) continue; 1744 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 1745 if (numCells == 2) { 1746 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1747 for (c = 0; c < 2; c++) { 1748 PetscInt cell = fcells[c]; 1749 1750 if (cell >= cStart && cell < cEndInterior) { 1751 ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr); 1752 } 1753 } 1754 } 1755 } 1756 ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr); 1757 ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr); 1758 ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1759 nStart = 0; 1760 ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr); 1761 ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr); 1762 ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr); 1763 for (f = fStart; f < fEnd; f++) { 1764 const PetscInt *fcells; 1765 PetscBool boundary; 1766 PetscInt ghost = -1; 1767 PetscInt numChildren, numCells, c; 1768 1769 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1770 ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1771 ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1772 if ((ghost >= 0) || boundary || numChildren) continue; 1773 ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 1774 if (numCells == 2) { 1775 ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1776 for (c = 0; c < 2; c++) { 1777 PetscInt cell = fcells[c], off; 1778 1779 if (cell >= cStart && cell < cEndInterior) { 1780 ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr); 1781 off += counter[cell - cStart]++; 1782 neighbors[off][0] = f; 1783 neighbors[off][1] = fcells[1 - c]; 1784 } 1785 } 1786 } 1787 } 1788 ierr = PetscFree(counter);CHKERRQ(ierr); 1789 ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1790 for (c = cStart; c < cEndInterior; c++) { 1791 PetscInt numFaces, f, d, off, ghost = -1; 1792 PetscFVCellGeom *cg; 1793 1794 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1795 ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr); 1796 ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr); 1797 if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);} 1798 if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1799 for (f = 0; f < numFaces; ++f) { 1800 PetscFVCellGeom *cg1; 1801 PetscFVFaceGeom *fg; 1802 const PetscInt *fcells; 1803 PetscInt ncell, side, nface; 1804 1805 nface = neighbors[off + f][0]; 1806 ncell = neighbors[off + f][1]; 1807 ierr = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr); 1808 side = (c != fcells[0]); 1809 ierr = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr); 1810 ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1811 for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1812 gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1813 } 1814 ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr); 1815 for (f = 0; f < numFaces; ++f) { 1816 for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d]; 1817 } 1818 } 1819 ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1820 ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr); 1821 ierr = PetscFree(neighbors);CHKERRQ(ierr); 1822 PetscFunctionReturn(0); 1823 } 1824 1825 #undef __FUNCT__ 1826 #define __FUNCT__ "DMPlexComputeGradientFVM" 1827 /*@ 1828 DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 1829 1830 Collective on DM 1831 1832 Input Arguments: 1833 + dm - The DM 1834 . fvm - The PetscFV 1835 . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 1836 - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 1837 1838 Output Parameters: 1839 + faceGeometry - The geometric factors for gradient calculation are inserted 1840 - dmGrad - The DM describing the layout of gradient data 1841 1842 Level: developer 1843 1844 .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 1845 @*/ 1846 PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 1847 { 1848 DM dmFace, dmCell; 1849 PetscScalar *fgeom, *cgeom; 1850 PetscSection sectionGrad, parentSection; 1851 PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 1852 PetscErrorCode ierr; 1853 1854 PetscFunctionBegin; 1855 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1856 ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1857 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1858 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1859 /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 1860 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1861 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1862 ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1863 ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1864 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1865 if (!parentSection) { 1866 ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1867 } 1868 else { 1869 ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1870 } 1871 ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1872 ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1873 /* Create storage for gradients */ 1874 ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 1875 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 1876 ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 1877 for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 1878 ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 1879 ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 1880 ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 1881 PetscFunctionReturn(0); 1882 } 1883