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