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