1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2ccd2543fSMatthew G Knepley 3ccd2543fSMatthew G Knepley #undef __FUNCT__ 4fea14342SMatthew G. Knepley #define __FUNCT__ "DMPlexGetLineIntersection_2D_Internal" 5fea14342SMatthew G. Knepley static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 6fea14342SMatthew G. Knepley { 7fea14342SMatthew G. Knepley const PetscReal p0_x = segmentA[0*2+0]; 8fea14342SMatthew G. Knepley const PetscReal p0_y = segmentA[0*2+1]; 9fea14342SMatthew G. Knepley const PetscReal p1_x = segmentA[1*2+0]; 10fea14342SMatthew G. Knepley const PetscReal p1_y = segmentA[1*2+1]; 11fea14342SMatthew G. Knepley const PetscReal p2_x = segmentB[0*2+0]; 12fea14342SMatthew G. Knepley const PetscReal p2_y = segmentB[0*2+1]; 13fea14342SMatthew G. Knepley const PetscReal p3_x = segmentB[1*2+0]; 14fea14342SMatthew G. Knepley const PetscReal p3_y = segmentB[1*2+1]; 15fea14342SMatthew G. Knepley const PetscReal s1_x = p1_x - p0_x; 16fea14342SMatthew G. Knepley const PetscReal s1_y = p1_y - p0_y; 17fea14342SMatthew G. Knepley const PetscReal s2_x = p3_x - p2_x; 18fea14342SMatthew G. Knepley const PetscReal s2_y = p3_y - p2_y; 19fea14342SMatthew G. Knepley const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 20fea14342SMatthew G. Knepley 21fea14342SMatthew G. Knepley PetscFunctionBegin; 22fea14342SMatthew G. Knepley *hasIntersection = PETSC_FALSE; 23fea14342SMatthew G. Knepley /* Non-parallel lines */ 24fea14342SMatthew G. Knepley if (denom != 0.0) { 25fea14342SMatthew G. Knepley const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 26fea14342SMatthew G. Knepley const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 27fea14342SMatthew G. Knepley 28fea14342SMatthew G. Knepley if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 29fea14342SMatthew G. Knepley *hasIntersection = PETSC_TRUE; 30fea14342SMatthew G. Knepley if (intersection) { 31fea14342SMatthew G. Knepley intersection[0] = p0_x + (t * s1_x); 32fea14342SMatthew G. Knepley intersection[1] = p0_y + (t * s1_y); 33fea14342SMatthew G. Knepley } 34fea14342SMatthew G. Knepley } 35fea14342SMatthew G. Knepley } 36fea14342SMatthew G. Knepley PetscFunctionReturn(0); 37fea14342SMatthew G. Knepley } 38fea14342SMatthew G. Knepley 39fea14342SMatthew G. Knepley #undef __FUNCT__ 40ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal" 41ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 42ccd2543fSMatthew G Knepley { 43ccd2543fSMatthew G Knepley const PetscInt embedDim = 2; 44f5ebc837SMatthew G. Knepley const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 45ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 46ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 47ccd2543fSMatthew G Knepley PetscReal v0[2], J[4], invJ[4], detJ; 48ccd2543fSMatthew G Knepley PetscReal xi, eta; 49ccd2543fSMatthew G Knepley PetscErrorCode ierr; 50ccd2543fSMatthew G Knepley 51ccd2543fSMatthew G Knepley PetscFunctionBegin; 528e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 53ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 54ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 55ccd2543fSMatthew G Knepley 56f5ebc837SMatthew G. Knepley if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c; 57ccd2543fSMatthew G Knepley else *cell = -1; 58ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 59ccd2543fSMatthew G Knepley } 60ccd2543fSMatthew G Knepley 61ccd2543fSMatthew G Knepley #undef __FUNCT__ 62ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal" 63ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 64ccd2543fSMatthew G Knepley { 65ccd2543fSMatthew G Knepley PetscSection coordSection; 66ccd2543fSMatthew G Knepley Vec coordsLocal; 67a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 68ccd2543fSMatthew G Knepley const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 69ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 70ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 71ccd2543fSMatthew G Knepley PetscInt crossings = 0, f; 72ccd2543fSMatthew G Knepley PetscErrorCode ierr; 73ccd2543fSMatthew G Knepley 74ccd2543fSMatthew G Knepley PetscFunctionBegin; 75ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 7669d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 77ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 78ccd2543fSMatthew G Knepley for (f = 0; f < 4; ++f) { 79ccd2543fSMatthew G Knepley PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 80ccd2543fSMatthew G Knepley PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 81ccd2543fSMatthew G Knepley PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 82ccd2543fSMatthew G Knepley PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 83ccd2543fSMatthew G Knepley PetscReal slope = (y_j - y_i) / (x_j - x_i); 84ccd2543fSMatthew G Knepley PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 85ccd2543fSMatthew G Knepley PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 86ccd2543fSMatthew G Knepley PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 87ccd2543fSMatthew G Knepley if ((cond1 || cond2) && above) ++crossings; 88ccd2543fSMatthew G Knepley } 89ccd2543fSMatthew G Knepley if (crossings % 2) *cell = c; 90ccd2543fSMatthew G Knepley else *cell = -1; 91ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 92ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 93ccd2543fSMatthew G Knepley } 94ccd2543fSMatthew G Knepley 95ccd2543fSMatthew G Knepley #undef __FUNCT__ 96ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal" 97ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 98ccd2543fSMatthew G Knepley { 99ccd2543fSMatthew G Knepley const PetscInt embedDim = 3; 100ccd2543fSMatthew G Knepley PetscReal v0[3], J[9], invJ[9], detJ; 101ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 102ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 103ccd2543fSMatthew G Knepley PetscReal z = PetscRealPart(point[2]); 104ccd2543fSMatthew G Knepley PetscReal xi, eta, zeta; 105ccd2543fSMatthew G Knepley PetscErrorCode ierr; 106ccd2543fSMatthew G Knepley 107ccd2543fSMatthew G Knepley PetscFunctionBegin; 1088e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 109ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 110ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 111ccd2543fSMatthew G Knepley zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 112ccd2543fSMatthew G Knepley 113ccd2543fSMatthew G Knepley if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 114ccd2543fSMatthew G Knepley else *cell = -1; 115ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 116ccd2543fSMatthew G Knepley } 117ccd2543fSMatthew G Knepley 118ccd2543fSMatthew G Knepley #undef __FUNCT__ 119ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal" 120ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 121ccd2543fSMatthew G Knepley { 122ccd2543fSMatthew G Knepley PetscSection coordSection; 123ccd2543fSMatthew G Knepley Vec coordsLocal; 1247c1f9639SMatthew G Knepley PetscScalar *coords; 125fb150da6SMatthew G. Knepley const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 126fb150da6SMatthew G. Knepley 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 127ccd2543fSMatthew G Knepley PetscBool found = PETSC_TRUE; 128ccd2543fSMatthew G Knepley PetscInt f; 129ccd2543fSMatthew G Knepley PetscErrorCode ierr; 130ccd2543fSMatthew G Knepley 131ccd2543fSMatthew G Knepley PetscFunctionBegin; 132ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 13369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 134ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 135ccd2543fSMatthew G Knepley for (f = 0; f < 6; ++f) { 136ccd2543fSMatthew G Knepley /* Check the point is under plane */ 137ccd2543fSMatthew G Knepley /* Get face normal */ 138ccd2543fSMatthew G Knepley PetscReal v_i[3]; 139ccd2543fSMatthew G Knepley PetscReal v_j[3]; 140ccd2543fSMatthew G Knepley PetscReal normal[3]; 141ccd2543fSMatthew G Knepley PetscReal pp[3]; 142ccd2543fSMatthew G Knepley PetscReal dot; 143ccd2543fSMatthew G Knepley 144ccd2543fSMatthew G Knepley v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 145ccd2543fSMatthew G Knepley v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 146ccd2543fSMatthew G Knepley v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 147ccd2543fSMatthew G Knepley v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 148ccd2543fSMatthew G Knepley v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 149ccd2543fSMatthew G Knepley v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 150ccd2543fSMatthew G Knepley normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 151ccd2543fSMatthew G Knepley normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 152ccd2543fSMatthew G Knepley normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 153ccd2543fSMatthew G Knepley pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 154ccd2543fSMatthew G Knepley pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 155ccd2543fSMatthew G Knepley pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 156ccd2543fSMatthew G Knepley dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 157ccd2543fSMatthew G Knepley 158ccd2543fSMatthew G Knepley /* Check that projected point is in face (2D location problem) */ 159ccd2543fSMatthew G Knepley if (dot < 0.0) { 160ccd2543fSMatthew G Knepley found = PETSC_FALSE; 161ccd2543fSMatthew G Knepley break; 162ccd2543fSMatthew G Knepley } 163ccd2543fSMatthew G Knepley } 164ccd2543fSMatthew G Knepley if (found) *cell = c; 165ccd2543fSMatthew G Knepley else *cell = -1; 166ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 167ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 168ccd2543fSMatthew G Knepley } 169ccd2543fSMatthew G Knepley 170ccd2543fSMatthew G Knepley #undef __FUNCT__ 171c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashInitialize_Internal" 172c4eade1cSMatthew G. Knepley static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 173c4eade1cSMatthew G. Knepley { 174c4eade1cSMatthew G. Knepley PetscInt d; 175c4eade1cSMatthew G. Knepley 176c4eade1cSMatthew G. Knepley PetscFunctionBegin; 177c4eade1cSMatthew G. Knepley box->dim = dim; 178c4eade1cSMatthew G. Knepley for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 179c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 180c4eade1cSMatthew G. Knepley } 181c4eade1cSMatthew G. Knepley 182c4eade1cSMatthew G. Knepley #undef __FUNCT__ 183c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashCreate" 184c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 185c4eade1cSMatthew G. Knepley { 186c4eade1cSMatthew G. Knepley PetscErrorCode ierr; 187c4eade1cSMatthew G. Knepley 188c4eade1cSMatthew G. Knepley PetscFunctionBegin; 189c4eade1cSMatthew G. Knepley ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 190c4eade1cSMatthew G. Knepley ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 191c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 192c4eade1cSMatthew G. Knepley } 193c4eade1cSMatthew G. Knepley 194c4eade1cSMatthew G. Knepley #undef __FUNCT__ 195c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashEnlarge" 196c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 197c4eade1cSMatthew G. Knepley { 198c4eade1cSMatthew G. Knepley PetscInt d; 199c4eade1cSMatthew G. Knepley 200c4eade1cSMatthew G. Knepley PetscFunctionBegin; 201c4eade1cSMatthew G. Knepley for (d = 0; d < box->dim; ++d) { 202c4eade1cSMatthew G. Knepley box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 203c4eade1cSMatthew G. Knepley box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 204c4eade1cSMatthew G. Knepley } 205c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 206c4eade1cSMatthew G. Knepley } 207c4eade1cSMatthew G. Knepley 208c4eade1cSMatthew G. Knepley #undef __FUNCT__ 209c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashSetGrid" 210c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 211c4eade1cSMatthew G. Knepley { 212c4eade1cSMatthew G. Knepley PetscInt d; 213c4eade1cSMatthew G. Knepley 214c4eade1cSMatthew G. Knepley PetscFunctionBegin; 215c4eade1cSMatthew G. Knepley for (d = 0; d < box->dim; ++d) { 216c4eade1cSMatthew G. Knepley box->extent[d] = box->upper[d] - box->lower[d]; 217c4eade1cSMatthew G. Knepley if (n[d] == PETSC_DETERMINE) { 218c4eade1cSMatthew G. Knepley box->h[d] = h[d]; 219c4eade1cSMatthew G. Knepley box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 220c4eade1cSMatthew G. Knepley } else { 221c4eade1cSMatthew G. Knepley box->n[d] = n[d]; 222c4eade1cSMatthew G. Knepley box->h[d] = box->extent[d]/n[d]; 223c4eade1cSMatthew G. Knepley } 224c4eade1cSMatthew G. Knepley } 225c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 226c4eade1cSMatthew G. Knepley } 227c4eade1cSMatthew G. Knepley 228c4eade1cSMatthew G. Knepley #undef __FUNCT__ 229c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashGetEnclosingBox" 2301c6dfc3eSMatthew G. Knepley PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 231c4eade1cSMatthew G. Knepley { 232c4eade1cSMatthew G. Knepley const PetscReal *lower = box->lower; 233c4eade1cSMatthew G. Knepley const PetscReal *upper = box->upper; 234c4eade1cSMatthew G. Knepley const PetscReal *h = box->h; 235c4eade1cSMatthew G. Knepley const PetscInt *n = box->n; 236c4eade1cSMatthew G. Knepley const PetscInt dim = box->dim; 237c4eade1cSMatthew G. Knepley PetscInt d, p; 238c4eade1cSMatthew G. Knepley 239c4eade1cSMatthew G. Knepley PetscFunctionBegin; 240c4eade1cSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 241c4eade1cSMatthew G. Knepley for (d = 0; d < dim; ++d) { 2421c6dfc3eSMatthew G. Knepley PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 243c4eade1cSMatthew G. Knepley 2441c6dfc3eSMatthew G. Knepley if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 245c4eade1cSMatthew G. Knepley 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", 2461c6dfc3eSMatthew G. Knepley p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0); 247c4eade1cSMatthew G. Knepley dboxes[p*dim+d] = dbox; 248c4eade1cSMatthew G. Knepley } 249c4eade1cSMatthew G. Knepley if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 250c4eade1cSMatthew G. Knepley } 251c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 252c4eade1cSMatthew G. Knepley } 253c4eade1cSMatthew G. Knepley 254c4eade1cSMatthew G. Knepley #undef __FUNCT__ 255c4eade1cSMatthew G. Knepley #define __FUNCT__ "PetscGridHashDestroy" 256c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 257c4eade1cSMatthew G. Knepley { 258c4eade1cSMatthew G. Knepley PetscErrorCode ierr; 259c4eade1cSMatthew G. Knepley 260c4eade1cSMatthew G. Knepley PetscFunctionBegin; 261c4eade1cSMatthew G. Knepley if (*box) { 262c4eade1cSMatthew G. Knepley ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 263c4eade1cSMatthew G. Knepley ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 264c4eade1cSMatthew G. Knepley ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 265c4eade1cSMatthew G. Knepley } 266c4eade1cSMatthew G. Knepley ierr = PetscFree(*box);CHKERRQ(ierr); 267c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 268c4eade1cSMatthew G. Knepley } 269c4eade1cSMatthew G. Knepley 270cafe43deSMatthew G. Knepley #undef __FUNCT__ 271cafe43deSMatthew G. Knepley #define __FUNCT__ "DMPlexLocatePoint_Internal" 272cafe43deSMatthew G. Knepley PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 273cafe43deSMatthew G. Knepley { 274cafe43deSMatthew G. Knepley PetscInt coneSize; 275cafe43deSMatthew G. Knepley PetscErrorCode ierr; 276cafe43deSMatthew G. Knepley 277cafe43deSMatthew G. Knepley PetscFunctionBegin; 278cafe43deSMatthew G. Knepley switch (dim) { 279cafe43deSMatthew G. Knepley case 2: 280cafe43deSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 281cafe43deSMatthew G. Knepley switch (coneSize) { 282cafe43deSMatthew G. Knepley case 3: 283cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 284cafe43deSMatthew G. Knepley break; 285cafe43deSMatthew G. Knepley case 4: 286cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 287cafe43deSMatthew G. Knepley break; 288cafe43deSMatthew G. Knepley default: 289cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 290cafe43deSMatthew G. Knepley } 291cafe43deSMatthew G. Knepley break; 292cafe43deSMatthew G. Knepley case 3: 293cafe43deSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 294cafe43deSMatthew G. Knepley switch (coneSize) { 295cafe43deSMatthew G. Knepley case 4: 296cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 297cafe43deSMatthew G. Knepley break; 298cafe43deSMatthew G. Knepley case 6: 299cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 300cafe43deSMatthew G. Knepley break; 301cafe43deSMatthew G. Knepley default: 302cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 303cafe43deSMatthew G. Knepley } 304cafe43deSMatthew G. Knepley break; 305cafe43deSMatthew G. Knepley default: 306cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 307cafe43deSMatthew G. Knepley } 308cafe43deSMatthew G. Knepley PetscFunctionReturn(0); 309cafe43deSMatthew G. Knepley } 310cafe43deSMatthew G. Knepley 311cafe43deSMatthew G. Knepley #undef __FUNCT__ 312cafe43deSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGridHash_Internal" 313cafe43deSMatthew G. Knepley PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 314cafe43deSMatthew G. Knepley { 315cafe43deSMatthew G. Knepley MPI_Comm comm; 316cafe43deSMatthew G. Knepley PetscGridHash lbox; 317cafe43deSMatthew G. Knepley Vec coordinates; 318cafe43deSMatthew G. Knepley PetscSection coordSection; 319cafe43deSMatthew G. Knepley Vec coordsLocal; 320cafe43deSMatthew G. Knepley const PetscScalar *coords; 321722d0f5cSMatthew G. Knepley PetscInt *dboxes, *boxes; 322cafe43deSMatthew G. Knepley PetscInt n[3] = {10, 10, 10}; 3231d0c6c94SMatthew G. Knepley PetscInt dim, N, cStart, cEnd, cMax, c, i; 324cafe43deSMatthew G. Knepley PetscErrorCode ierr; 325cafe43deSMatthew G. Knepley 326cafe43deSMatthew G. Knepley PetscFunctionBegin; 327cafe43deSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 328cafe43deSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 329cafe43deSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 330cafe43deSMatthew G. Knepley ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 331cafe43deSMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 332cafe43deSMatthew G. Knepley ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr); 333cafe43deSMatthew G. Knepley for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);} 334cafe43deSMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 335cafe43deSMatthew G. Knepley ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr); 336cafe43deSMatthew G. Knepley #if 0 337cafe43deSMatthew G. Knepley /* Could define a custom reduction to merge these */ 338b2566f29SBarry Smith ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 339b2566f29SBarry Smith ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr); 340cafe43deSMatthew G. Knepley #endif 341cafe43deSMatthew G. Knepley /* Is there a reason to snap the local bounding box to a division of the global box? */ 342cafe43deSMatthew G. Knepley /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */ 343cafe43deSMatthew G. Knepley /* Create label */ 344cafe43deSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 3451d0c6c94SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 3461d0c6c94SMatthew G. Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 347cafe43deSMatthew G. Knepley ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 348cafe43deSMatthew G. Knepley ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 349722d0f5cSMatthew G. Knepley /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 350cafe43deSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 351cafe43deSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 35238353de4SMatthew G. Knepley ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr); 353cafe43deSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 354cafe43deSMatthew G. Knepley const PetscReal *h = lbox->h; 355cafe43deSMatthew G. Knepley PetscScalar *ccoords = NULL; 35638353de4SMatthew G. Knepley PetscInt csize = 0; 357cafe43deSMatthew G. Knepley PetscScalar point[3]; 358cafe43deSMatthew G. Knepley PetscInt dlim[6], d, e, i, j, k; 359cafe43deSMatthew G. Knepley 360cafe43deSMatthew G. Knepley /* Find boxes enclosing each vertex */ 36138353de4SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr); 36238353de4SMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr); 363722d0f5cSMatthew G. Knepley /* Mark cells containing the vertices */ 36438353de4SMatthew G. Knepley for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 365cafe43deSMatthew G. Knepley /* Get grid of boxes containing these */ 366cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 3672291669eSMatthew G. Knepley for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 368cafe43deSMatthew G. Knepley for (e = 1; e < dim+1; ++e) { 369cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) { 370cafe43deSMatthew G. Knepley dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 371cafe43deSMatthew G. Knepley dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 372cafe43deSMatthew G. Knepley } 373cafe43deSMatthew G. Knepley } 374fea14342SMatthew G. Knepley /* Check for intersection of box with cell */ 375cafe43deSMatthew G. Knepley for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) { 376cafe43deSMatthew G. Knepley for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) { 377cafe43deSMatthew G. Knepley for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) { 378cafe43deSMatthew G. Knepley const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 379cafe43deSMatthew G. Knepley PetscScalar cpoint[3]; 380fea14342SMatthew G. Knepley PetscInt cell, edge, ii, jj, kk; 381cafe43deSMatthew G. Knepley 382fea14342SMatthew G. Knepley /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 383cafe43deSMatthew G. Knepley for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 384cafe43deSMatthew G. Knepley for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 385cafe43deSMatthew G. Knepley for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 386cafe43deSMatthew G. Knepley 387cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 388cafe43deSMatthew G. Knepley if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 389cafe43deSMatthew G. Knepley } 390cafe43deSMatthew G. Knepley } 391cafe43deSMatthew G. Knepley } 392fea14342SMatthew G. Knepley /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 393fea14342SMatthew G. Knepley for (edge = 0; edge < dim+1; ++edge) { 394fea14342SMatthew G. Knepley PetscReal segA[6], segB[6]; 395fea14342SMatthew G. Knepley 396fea14342SMatthew G. Knepley for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);} 397fea14342SMatthew G. Knepley for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 3989a128ed2SMatthew G. Knepley if (dim > 2) {segB[2] = PetscRealPart(point[2]); 3999a128ed2SMatthew G. Knepley segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 400fea14342SMatthew G. Knepley for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 4019a128ed2SMatthew G. Knepley if (dim > 1) {segB[1] = PetscRealPart(point[1]); 4029a128ed2SMatthew G. Knepley segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 403fea14342SMatthew G. Knepley for (ii = 0; ii < 2; ++ii) { 404fea14342SMatthew G. Knepley PetscBool intersects; 405fea14342SMatthew G. Knepley 4069a128ed2SMatthew G. Knepley segB[0] = PetscRealPart(point[0]); 4079a128ed2SMatthew G. Knepley segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 408fea14342SMatthew G. Knepley ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 409fea14342SMatthew G. Knepley if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 410cafe43deSMatthew G. Knepley } 411cafe43deSMatthew G. Knepley } 412cafe43deSMatthew G. Knepley } 413cafe43deSMatthew G. Knepley } 414fea14342SMatthew G. Knepley } 415fea14342SMatthew G. Knepley } 416fea14342SMatthew G. Knepley } 417fea14342SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 418fea14342SMatthew G. Knepley } 419722d0f5cSMatthew G. Knepley ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 420cafe43deSMatthew G. Knepley ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 421cafe43deSMatthew G. Knepley ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 422cafe43deSMatthew G. Knepley *localBox = lbox; 423cafe43deSMatthew G. Knepley PetscFunctionReturn(0); 424cafe43deSMatthew G. Knepley } 425cafe43deSMatthew G. Knepley 426cafe43deSMatthew G. Knepley #undef __FUNCT__ 427ccd2543fSMatthew G Knepley #define __FUNCT__ "DMLocatePoints_Plex" 4283a93e3b7SToby Isaac PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, PetscSF cellSF) 429ccd2543fSMatthew G Knepley { 430cafe43deSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 431953fc75cSMatthew G. Knepley PetscBool hash = mesh->useHashLocation; 4323a93e3b7SToby Isaac PetscInt bs, numPoints, p, numFound, *found = NULL; 4331318edbeSMatthew G. Knepley PetscInt dim, cStart, cEnd, cMax, numCells, c; 434cafe43deSMatthew G. Knepley const PetscInt *boxCells; 4353a93e3b7SToby Isaac PetscSFNode *cells; 436ccd2543fSMatthew G Knepley PetscScalar *a; 4373a93e3b7SToby Isaac PetscMPIInt result; 438ccd2543fSMatthew G Knepley PetscErrorCode ierr; 439ccd2543fSMatthew G Knepley 440ccd2543fSMatthew G Knepley PetscFunctionBegin; 441cafe43deSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 442cafe43deSMatthew G. Knepley ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 4433a93e3b7SToby Isaac ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr); 4443a93e3b7SToby Isaac if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 445cafe43deSMatthew G. Knepley 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); 446ccd2543fSMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 447ccd2543fSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 448ccd2543fSMatthew G Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 449ccd2543fSMatthew G Knepley ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 450ccd2543fSMatthew G Knepley ierr = VecGetArray(v, &a);CHKERRQ(ierr); 451ccd2543fSMatthew G Knepley numPoints /= bs; 452785e854fSJed Brown ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 453953fc75cSMatthew G. Knepley if (hash) { 454ac6ec2abSMatthew G. Knepley if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 455cafe43deSMatthew G. Knepley /* Designate the local box for each point */ 456cafe43deSMatthew G. Knepley /* Send points to correct process */ 457cafe43deSMatthew G. Knepley /* Search cells that lie in each subbox */ 458cafe43deSMatthew G. Knepley /* Should we bin points before doing search? */ 459cafe43deSMatthew G. Knepley ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 460953fc75cSMatthew G. Knepley } 4613a93e3b7SToby Isaac for (p = 0, numFound = 0; p < numPoints; ++p) { 462ccd2543fSMatthew G Knepley const PetscScalar *point = &a[p*bs]; 463953fc75cSMatthew G. Knepley PetscInt dbin[3], bin, cell = -1, cellOffset; 464ccd2543fSMatthew G Knepley 4653a93e3b7SToby Isaac cells[p].rank = -1; 4663a93e3b7SToby Isaac cells[p].index = -1; 467953fc75cSMatthew G. Knepley if (hash) { 468cafe43deSMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 469cafe43deSMatthew G. Knepley /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 470cafe43deSMatthew G. Knepley ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 471cafe43deSMatthew G. Knepley ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 472cafe43deSMatthew G. Knepley for (c = cellOffset; c < cellOffset + numCells; ++c) { 473cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 4743a93e3b7SToby Isaac if (cell >= 0) { 4753a93e3b7SToby Isaac cells[p].rank = 0; 4763a93e3b7SToby Isaac cells[p].index = cell; 4773a93e3b7SToby Isaac numFound++; 4783a93e3b7SToby Isaac break; 479ccd2543fSMatthew G Knepley } 4803a93e3b7SToby Isaac } 481953fc75cSMatthew G. Knepley } else { 482953fc75cSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 483953fc75cSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 4843a93e3b7SToby Isaac if (cell >= 0) { 4853a93e3b7SToby Isaac cells[p].rank = 0; 4863a93e3b7SToby Isaac cells[p].index = cell; 4873a93e3b7SToby Isaac numFound++; 4883a93e3b7SToby Isaac break; 489953fc75cSMatthew G. Knepley } 490953fc75cSMatthew G. Knepley } 4913a93e3b7SToby Isaac } 492ccd2543fSMatthew G Knepley } 493953fc75cSMatthew G. Knepley if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);} 494cafe43deSMatthew G. Knepley /* Check for highest numbered proc that claims a point (do we care?) */ 495ccd2543fSMatthew G Knepley ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 4963a93e3b7SToby Isaac if (numFound < numPoints) { 4973a93e3b7SToby Isaac ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr); 4983a93e3b7SToby Isaac for (p = 0, numFound = 0; p < numPoints; p++) { 4993a93e3b7SToby Isaac if (cells[p].rank >= 0 && cells[p].index >= 0) { 5003a93e3b7SToby Isaac if (numFound < p) { 5013a93e3b7SToby Isaac cells[numFound] = cells[p]; 5023a93e3b7SToby Isaac } 5033a93e3b7SToby Isaac found[numFound++] = p; 5043a93e3b7SToby Isaac } 5053a93e3b7SToby Isaac } 5063a93e3b7SToby Isaac } 5073a93e3b7SToby Isaac ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 508ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 509ccd2543fSMatthew G Knepley } 510ccd2543fSMatthew G Knepley 511ccd2543fSMatthew G Knepley #undef __FUNCT__ 51217fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 51317fe8556SMatthew G. Knepley /* 51417fe8556SMatthew G. Knepley DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 51517fe8556SMatthew G. Knepley */ 5163beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 51717fe8556SMatthew G. Knepley { 51817fe8556SMatthew G. Knepley const PetscReal x = PetscRealPart(coords[2] - coords[0]); 51917fe8556SMatthew G. Knepley const PetscReal y = PetscRealPart(coords[3] - coords[1]); 5208b49ba18SBarry Smith const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 52117fe8556SMatthew G. Knepley 52217fe8556SMatthew G. Knepley PetscFunctionBegin; 5231c99cf0cSGeoffrey Irving R[0] = c; R[1] = -s; 5241c99cf0cSGeoffrey Irving R[2] = s; R[3] = c; 52517fe8556SMatthew G. Knepley coords[0] = 0.0; 5267f07f362SMatthew G. Knepley coords[1] = r; 52717fe8556SMatthew G. Knepley PetscFunctionReturn(0); 52817fe8556SMatthew G. Knepley } 52917fe8556SMatthew G. Knepley 53017fe8556SMatthew G. Knepley #undef __FUNCT__ 53128dbe442SToby Isaac #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal" 53228dbe442SToby Isaac /* 53328dbe442SToby Isaac DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D 53428dbe442SToby Isaac 53528dbe442SToby Isaac This uses the basis completion described by Frisvad, 53628dbe442SToby Isaac 53728dbe442SToby Isaac http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html 53828dbe442SToby Isaac DOI:10.1080/2165347X.2012.689606 53928dbe442SToby Isaac */ 5403beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 54128dbe442SToby Isaac { 54228dbe442SToby Isaac PetscReal x = PetscRealPart(coords[3] - coords[0]); 54328dbe442SToby Isaac PetscReal y = PetscRealPart(coords[4] - coords[1]); 54428dbe442SToby Isaac PetscReal z = PetscRealPart(coords[5] - coords[2]); 54528dbe442SToby Isaac PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 54628dbe442SToby Isaac PetscReal rinv = 1. / r; 54728dbe442SToby Isaac PetscFunctionBegin; 54828dbe442SToby Isaac 54928dbe442SToby Isaac x *= rinv; y *= rinv; z *= rinv; 55028dbe442SToby Isaac if (x > 0.) { 55128dbe442SToby Isaac PetscReal inv1pX = 1./ (1. + x); 55228dbe442SToby Isaac 55328dbe442SToby Isaac R[0] = x; R[1] = -y; R[2] = -z; 55428dbe442SToby Isaac R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 55528dbe442SToby Isaac R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 55628dbe442SToby Isaac } 55728dbe442SToby Isaac else { 55828dbe442SToby Isaac PetscReal inv1mX = 1./ (1. - x); 55928dbe442SToby Isaac 56028dbe442SToby Isaac R[0] = x; R[1] = z; R[2] = y; 56128dbe442SToby Isaac R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 56228dbe442SToby Isaac R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 56328dbe442SToby Isaac } 56428dbe442SToby Isaac coords[0] = 0.0; 56528dbe442SToby Isaac coords[1] = r; 56628dbe442SToby Isaac PetscFunctionReturn(0); 56728dbe442SToby Isaac } 56828dbe442SToby Isaac 56928dbe442SToby Isaac #undef __FUNCT__ 570ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 571ccd2543fSMatthew G Knepley /* 572ccd2543fSMatthew G Knepley DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 573ccd2543fSMatthew G Knepley */ 5743beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 575ccd2543fSMatthew G Knepley { 5761ee9d5ecSMatthew G. Knepley PetscReal x1[3], x2[3], n[3], norm; 57799dec3a6SMatthew G. Knepley PetscReal x1p[3], x2p[3], xnp[3]; 5784a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 579ccd2543fSMatthew G Knepley const PetscInt dim = 3; 58099dec3a6SMatthew G. Knepley PetscInt d, e, p; 581ccd2543fSMatthew G Knepley 582ccd2543fSMatthew G Knepley PetscFunctionBegin; 583ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 584ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 5851ee9d5ecSMatthew G. Knepley x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 5861ee9d5ecSMatthew G. Knepley x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 587ccd2543fSMatthew G Knepley } 588ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 589ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 590ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 5918b49ba18SBarry Smith norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 592ccd2543fSMatthew G Knepley n[0] /= norm; 593ccd2543fSMatthew G Knepley n[1] /= norm; 594ccd2543fSMatthew G Knepley n[2] /= norm; 595ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 596ccd2543fSMatthew G Knepley 597ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 598ccd2543fSMatthew G Knepley 599ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 600ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 601ccd2543fSMatthew G Knepley \ nx ny nz / 602ccd2543fSMatthew G Knepley 603ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 604ccd2543fSMatthew G Knepley */ 6058b49ba18SBarry Smith sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 60673868372SMatthew G. Knepley /* Check for n = z */ 60773868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 6087df32b8bSSanderA const PetscInt s = PetscSign(n[2]); 6097df32b8bSSanderA /* If nz < 0, rotate 180 degrees around x-axis */ 61099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 61199dec3a6SMatthew G. Knepley coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 6127df32b8bSSanderA coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s; 61373868372SMatthew G. Knepley } 61499dec3a6SMatthew G. Knepley coords[0] = 0.0; 61599dec3a6SMatthew G. Knepley coords[1] = 0.0; 6167df32b8bSSanderA coords[2] = x1[0]; 6177df32b8bSSanderA coords[3] = x1[1] * s; 6187df32b8bSSanderA coords[4] = x2[0]; 6197df32b8bSSanderA coords[5] = x2[1] * s; 6207df32b8bSSanderA R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 6217df32b8bSSanderA R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0; 6227df32b8bSSanderA R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s; 62373868372SMatthew G. Knepley PetscFunctionReturn(0); 62473868372SMatthew G. Knepley } 625da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 626ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 627ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 628ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 629ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 630ccd2543fSMatthew G Knepley x1p[d] = 0.0; 631ccd2543fSMatthew G Knepley x2p[d] = 0.0; 632ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 633ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 634ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 635ccd2543fSMatthew G Knepley } 636ccd2543fSMatthew G Knepley } 6378763be8eSMatthew G. Knepley if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 6388763be8eSMatthew G. Knepley if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 639ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 64099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 64199dec3a6SMatthew G. Knepley for (d = 0; d < dim; ++d) { 64299dec3a6SMatthew G. Knepley xnp[d] = 0.0; 64399dec3a6SMatthew G. Knepley for (e = 0; e < dim; ++e) { 64499dec3a6SMatthew G. Knepley xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 64599dec3a6SMatthew G. Knepley } 64699dec3a6SMatthew G. Knepley if (d < dim-1) coords[p*2+d] = xnp[d]; 64799dec3a6SMatthew G. Knepley } 64899dec3a6SMatthew G. Knepley } 649ccd2543fSMatthew G Knepley coords[0] = 0.0; 650ccd2543fSMatthew G Knepley coords[1] = 0.0; 651ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 652ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 653ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 654ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 6557f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 6567f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 6577f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 6587f07f362SMatthew G. Knepley PetscReal tmp; 6597f07f362SMatthew G. Knepley 6607f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 6617f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 6627f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 6637f07f362SMatthew G. Knepley } 6647f07f362SMatthew G. Knepley } 665ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 666ccd2543fSMatthew G Knepley } 667ccd2543fSMatthew G Knepley 668ccd2543fSMatthew G Knepley #undef __FUNCT__ 669834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal" 6706322fe33SJed Brown PETSC_UNUSED 671834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 672834e62ceSMatthew G. Knepley { 673834e62ceSMatthew G. Knepley /* Signed volume is 1/2 the determinant 674834e62ceSMatthew G. Knepley 675834e62ceSMatthew G. Knepley | 1 1 1 | 676834e62ceSMatthew G. Knepley | x0 x1 x2 | 677834e62ceSMatthew G. Knepley | y0 y1 y2 | 678834e62ceSMatthew G. Knepley 679834e62ceSMatthew G. Knepley but if x0,y0 is the origin, we have 680834e62ceSMatthew G. Knepley 681834e62ceSMatthew G. Knepley | x1 x2 | 682834e62ceSMatthew G. Knepley | y1 y2 | 683834e62ceSMatthew G. Knepley */ 684834e62ceSMatthew G. Knepley const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 685834e62ceSMatthew G. Knepley const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 686834e62ceSMatthew G. Knepley PetscReal M[4], detM; 687834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; 68886623015SMatthew G. Knepley M[2] = y1; M[3] = y2; 689923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(&detM, M); 690834e62ceSMatthew G. Knepley *vol = 0.5*detM; 6913bc0b13bSBarry Smith (void)PetscLogFlops(5.0); 692834e62ceSMatthew G. Knepley } 693834e62ceSMatthew G. Knepley 694834e62ceSMatthew G. Knepley #undef __FUNCT__ 695834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal" 696834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 697834e62ceSMatthew G. Knepley { 698923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(vol, coords); 699834e62ceSMatthew G. Knepley *vol *= 0.5; 700834e62ceSMatthew G. Knepley } 701834e62ceSMatthew G. Knepley 702834e62ceSMatthew G. Knepley #undef __FUNCT__ 703834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal" 7046322fe33SJed Brown PETSC_UNUSED 705834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 706834e62ceSMatthew G. Knepley { 707834e62ceSMatthew G. Knepley /* Signed volume is 1/6th of the determinant 708834e62ceSMatthew G. Knepley 709834e62ceSMatthew G. Knepley | 1 1 1 1 | 710834e62ceSMatthew G. Knepley | x0 x1 x2 x3 | 711834e62ceSMatthew G. Knepley | y0 y1 y2 y3 | 712834e62ceSMatthew G. Knepley | z0 z1 z2 z3 | 713834e62ceSMatthew G. Knepley 714834e62ceSMatthew G. Knepley but if x0,y0,z0 is the origin, we have 715834e62ceSMatthew G. Knepley 716834e62ceSMatthew G. Knepley | x1 x2 x3 | 717834e62ceSMatthew G. Knepley | y1 y2 y3 | 718834e62ceSMatthew G. Knepley | z1 z2 z3 | 719834e62ceSMatthew G. Knepley */ 720834e62ceSMatthew G. Knepley const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 721834e62ceSMatthew G. Knepley const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 722834e62ceSMatthew G. Knepley const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 723834e62ceSMatthew G. Knepley PetscReal M[9], detM; 724834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; M[2] = x3; 725834e62ceSMatthew G. Knepley M[3] = y1; M[4] = y2; M[5] = y3; 726834e62ceSMatthew G. Knepley M[6] = z1; M[7] = z2; M[8] = z3; 727923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(&detM, M); 728b7ad821dSMatthew G. Knepley *vol = -0.16666666666666666666666*detM; 7293bc0b13bSBarry Smith (void)PetscLogFlops(10.0); 730834e62ceSMatthew G. Knepley } 731834e62ceSMatthew G. Knepley 732834e62ceSMatthew G. Knepley #undef __FUNCT__ 7330ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 7340ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 7350ec8681fSMatthew G. Knepley { 736923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(vol, coords); 737b7ad821dSMatthew G. Knepley *vol *= -0.16666666666666666666666; 7380ec8681fSMatthew G. Knepley } 7390ec8681fSMatthew G. Knepley 7400ec8681fSMatthew G. Knepley #undef __FUNCT__ 74117fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 74217fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 74317fe8556SMatthew G. Knepley { 74417fe8556SMatthew G. Knepley PetscSection coordSection; 74517fe8556SMatthew G. Knepley Vec coordinates; 746a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 7477f07f362SMatthew G. Knepley PetscInt numCoords, d; 74817fe8556SMatthew G. Knepley PetscErrorCode ierr; 74917fe8556SMatthew G. Knepley 75017fe8556SMatthew G. Knepley PetscFunctionBegin; 75117fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 75269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 75317fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 754adac9986SMatthew G. Knepley if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 7557f07f362SMatthew G. Knepley *detJ = 0.0; 75628dbe442SToby Isaac if (numCoords == 6) { 75728dbe442SToby Isaac const PetscInt dim = 3; 75828dbe442SToby Isaac PetscReal R[9], J0; 75928dbe442SToby Isaac 76028dbe442SToby Isaac if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 76128dbe442SToby Isaac ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr); 76228dbe442SToby Isaac if (J) { 76328dbe442SToby Isaac J0 = 0.5*PetscRealPart(coords[1]); 76428dbe442SToby Isaac J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 76528dbe442SToby Isaac J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 76628dbe442SToby Isaac J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 76728dbe442SToby Isaac DMPlex_Det3D_Internal(detJ, J); 76828dbe442SToby Isaac if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 769adac9986SMatthew G. Knepley } 77028dbe442SToby Isaac } else if (numCoords == 4) { 7717f07f362SMatthew G. Knepley const PetscInt dim = 2; 7727f07f362SMatthew G. Knepley PetscReal R[4], J0; 7737f07f362SMatthew G. Knepley 7747f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 7757f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 77617fe8556SMatthew G. Knepley if (J) { 7777f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 7787f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 7797f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 780923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 781923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 782adac9986SMatthew G. Knepley } 7837f07f362SMatthew G. Knepley } else if (numCoords == 2) { 7847f07f362SMatthew G. Knepley const PetscInt dim = 1; 7857f07f362SMatthew G. Knepley 7867f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 7877f07f362SMatthew G. Knepley if (J) { 7887f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 78917fe8556SMatthew G. Knepley *detJ = J[0]; 7903bc0b13bSBarry Smith ierr = PetscLogFlops(2.0);CHKERRQ(ierr); 7913bc0b13bSBarry Smith if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);} 792adac9986SMatthew G. Knepley } 793796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 79417fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 79517fe8556SMatthew G. Knepley PetscFunctionReturn(0); 79617fe8556SMatthew G. Knepley } 79717fe8556SMatthew G. Knepley 79817fe8556SMatthew G. Knepley #undef __FUNCT__ 799ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 800ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 801ccd2543fSMatthew G Knepley { 802ccd2543fSMatthew G Knepley PetscSection coordSection; 803ccd2543fSMatthew G Knepley Vec coordinates; 804a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 8057f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 806ccd2543fSMatthew G Knepley PetscErrorCode ierr; 807ccd2543fSMatthew G Knepley 808ccd2543fSMatthew G Knepley PetscFunctionBegin; 809ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 81069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 811ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 8127f07f362SMatthew G. Knepley *detJ = 0.0; 813ccd2543fSMatthew G Knepley if (numCoords == 9) { 8147f07f362SMatthew G. Knepley const PetscInt dim = 3; 8157f07f362SMatthew G. Knepley PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 8167f07f362SMatthew G. Knepley 8177f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 81899dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 8197f07f362SMatthew G. Knepley if (J) { 820b7ad821dSMatthew G. Knepley const PetscInt pdim = 2; 821b7ad821dSMatthew G. Knepley 822b7ad821dSMatthew G. Knepley for (d = 0; d < pdim; d++) { 823b7ad821dSMatthew G. Knepley for (f = 0; f < pdim; f++) { 824b7ad821dSMatthew G. Knepley J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 825ccd2543fSMatthew G Knepley } 8267f07f362SMatthew G. Knepley } 8273bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 828923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 8297f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 8307f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 8317f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 8327f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 8337f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 8347f07f362SMatthew G. Knepley } 8357f07f362SMatthew G. Knepley } 8367f07f362SMatthew G. Knepley } 8373bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 8387f07f362SMatthew G. Knepley } 839923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 8407f07f362SMatthew G. Knepley } else if (numCoords == 6) { 8417f07f362SMatthew G. Knepley const PetscInt dim = 2; 8427f07f362SMatthew G. Knepley 8437f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 844ccd2543fSMatthew G Knepley if (J) { 845ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 846ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 847ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 848ccd2543fSMatthew G Knepley } 849ccd2543fSMatthew G Knepley } 8503bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 851923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 852ccd2543fSMatthew G Knepley } 853923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 854796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 855ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 856ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 857ccd2543fSMatthew G Knepley } 858ccd2543fSMatthew G Knepley 859ccd2543fSMatthew G Knepley #undef __FUNCT__ 860ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 861ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 862ccd2543fSMatthew G Knepley { 863ccd2543fSMatthew G Knepley PetscSection coordSection; 864ccd2543fSMatthew G Knepley Vec coordinates; 865a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 866*71f58de1SToby Isaac PetscInt numCoords, numSelfCoords, d, f, g; 867ccd2543fSMatthew G Knepley PetscErrorCode ierr; 868ccd2543fSMatthew G Knepley 869ccd2543fSMatthew G Knepley PetscFunctionBegin; 870ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 87169d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 872*71f58de1SToby Isaac ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr); 87399dec3a6SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 874*71f58de1SToby Isaac numCoords = numSelfCoords ? numSelfCoords : numCoords; 8757f07f362SMatthew G. Knepley *detJ = 0.0; 87699dec3a6SMatthew G. Knepley if (numCoords == 12) { 87799dec3a6SMatthew G. Knepley const PetscInt dim = 3; 87899dec3a6SMatthew G. Knepley PetscReal R[9], J0[9] = {1.0,0.0,0.0,0.0,1.0,0.0,0.0,0.0,1.0}; 87999dec3a6SMatthew G. Knepley 88099dec3a6SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 88199dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 88299dec3a6SMatthew G. Knepley if (J) { 88399dec3a6SMatthew G. Knepley const PetscInt pdim = 2; 88499dec3a6SMatthew G. Knepley 88599dec3a6SMatthew G. Knepley for (d = 0; d < pdim; d++) { 88699dec3a6SMatthew G. Knepley J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 88799dec3a6SMatthew G. Knepley J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 88899dec3a6SMatthew G. Knepley } 8893bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 890923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 89199dec3a6SMatthew G. Knepley for (d = 0; d < dim; d++) { 89299dec3a6SMatthew G. Knepley for (f = 0; f < dim; f++) { 89399dec3a6SMatthew G. Knepley J[d*dim+f] = 0.0; 89499dec3a6SMatthew G. Knepley for (g = 0; g < dim; g++) { 89599dec3a6SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 89699dec3a6SMatthew G. Knepley } 89799dec3a6SMatthew G. Knepley } 89899dec3a6SMatthew G. Knepley } 8993bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 90099dec3a6SMatthew G. Knepley } 901923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 902*71f58de1SToby Isaac } else if (numCoords == 8) { 90399dec3a6SMatthew G. Knepley const PetscInt dim = 2; 90499dec3a6SMatthew G. Knepley 9057f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 906ccd2543fSMatthew G Knepley if (J) { 907ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 90899dec3a6SMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 90999dec3a6SMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 910ccd2543fSMatthew G Knepley } 9113bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 912923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 913ccd2543fSMatthew G Knepley } 914923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 915796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 91699dec3a6SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 917ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 918ccd2543fSMatthew G Knepley } 919ccd2543fSMatthew G Knepley 920ccd2543fSMatthew G Knepley #undef __FUNCT__ 921ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 922ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 923ccd2543fSMatthew G Knepley { 924ccd2543fSMatthew G Knepley PetscSection coordSection; 925ccd2543fSMatthew G Knepley Vec coordinates; 926a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 927ccd2543fSMatthew G Knepley const PetscInt dim = 3; 92899dec3a6SMatthew G. Knepley PetscInt d; 929ccd2543fSMatthew G Knepley PetscErrorCode ierr; 930ccd2543fSMatthew G Knepley 931ccd2543fSMatthew G Knepley PetscFunctionBegin; 932ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 93369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 934ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 9357f07f362SMatthew G. Knepley *detJ = 0.0; 9367f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 937ccd2543fSMatthew G Knepley if (J) { 938ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 939f0df753eSMatthew G. Knepley /* I orient with outward face normals */ 940f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 941f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 942f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 943ccd2543fSMatthew G Knepley } 9443bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 945923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 946ccd2543fSMatthew G Knepley } 947923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 948ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 949ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 950ccd2543fSMatthew G Knepley } 951ccd2543fSMatthew G Knepley 952ccd2543fSMatthew G Knepley #undef __FUNCT__ 953ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 954ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 955ccd2543fSMatthew G Knepley { 956ccd2543fSMatthew G Knepley PetscSection coordSection; 957ccd2543fSMatthew G Knepley Vec coordinates; 958a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 959ccd2543fSMatthew G Knepley const PetscInt dim = 3; 960ccd2543fSMatthew G Knepley PetscInt d; 961ccd2543fSMatthew G Knepley PetscErrorCode ierr; 962ccd2543fSMatthew G Knepley 963ccd2543fSMatthew G Knepley PetscFunctionBegin; 964ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 96569d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 966ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 9677f07f362SMatthew G. Knepley *detJ = 0.0; 9687f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 969ccd2543fSMatthew G Knepley if (J) { 970ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 971f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 972f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 973f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 974ccd2543fSMatthew G Knepley } 9753bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 976923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 977ccd2543fSMatthew G Knepley } 978923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 979ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 980ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 981ccd2543fSMatthew G Knepley } 982ccd2543fSMatthew G Knepley 983ccd2543fSMatthew G Knepley #undef __FUNCT__ 9848e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 985ccd2543fSMatthew G Knepley /*@C 9868e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 987ccd2543fSMatthew G Knepley 988ccd2543fSMatthew G Knepley Collective on DM 989ccd2543fSMatthew G Knepley 990ccd2543fSMatthew G Knepley Input Arguments: 991ccd2543fSMatthew G Knepley + dm - the DM 992ccd2543fSMatthew G Knepley - cell - the cell 993ccd2543fSMatthew G Knepley 994ccd2543fSMatthew G Knepley Output Arguments: 995ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 996ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 997ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 998ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 999ccd2543fSMatthew G Knepley 1000ccd2543fSMatthew G Knepley Level: advanced 1001ccd2543fSMatthew G Knepley 1002ccd2543fSMatthew G Knepley Fortran Notes: 1003ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1004ccd2543fSMatthew G Knepley include petsc.h90 in your code. 1005ccd2543fSMatthew G Knepley 10068e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 1007ccd2543fSMatthew G Knepley @*/ 10088e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1009ccd2543fSMatthew G Knepley { 101049dc4407SMatthew G. Knepley PetscInt depth, dim, coneSize; 1011ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1012ccd2543fSMatthew G Knepley 1013ccd2543fSMatthew G Knepley PetscFunctionBegin; 1014139a35ccSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1015ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 101649dc4407SMatthew G. Knepley if (depth == 1) { 10178e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 10188e0841e0SMatthew G. Knepley } else { 10198e0841e0SMatthew G. Knepley DMLabel depth; 10208e0841e0SMatthew G. Knepley 10218e0841e0SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 10228e0841e0SMatthew G. Knepley ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 10238e0841e0SMatthew G. Knepley } 1024ccd2543fSMatthew G Knepley switch (dim) { 102517fe8556SMatthew G. Knepley case 1: 102617fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 102717fe8556SMatthew G. Knepley break; 1028ccd2543fSMatthew G Knepley case 2: 1029ccd2543fSMatthew G Knepley switch (coneSize) { 1030ccd2543fSMatthew G Knepley case 3: 1031ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1032ccd2543fSMatthew G Knepley break; 1033ccd2543fSMatthew G Knepley case 4: 1034ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1035ccd2543fSMatthew G Knepley break; 1036ccd2543fSMatthew G Knepley default: 10378e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1038ccd2543fSMatthew G Knepley } 1039ccd2543fSMatthew G Knepley break; 1040ccd2543fSMatthew G Knepley case 3: 1041ccd2543fSMatthew G Knepley switch (coneSize) { 1042ccd2543fSMatthew G Knepley case 4: 1043ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1044ccd2543fSMatthew G Knepley break; 10458e0841e0SMatthew G. Knepley case 6: /* Faces */ 10468e0841e0SMatthew G. Knepley case 8: /* Vertices */ 1047ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1048ccd2543fSMatthew G Knepley break; 1049ccd2543fSMatthew G Knepley default: 10508e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1051ccd2543fSMatthew G Knepley } 1052ccd2543fSMatthew G Knepley break; 1053ccd2543fSMatthew G Knepley default: 1054ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1055ccd2543fSMatthew G Knepley } 10568e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 10578e0841e0SMatthew G. Knepley } 10588e0841e0SMatthew G. Knepley 10598e0841e0SMatthew G. Knepley #undef __FUNCT__ 10608e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 10618e0841e0SMatthew G. Knepley static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 10628e0841e0SMatthew G. Knepley { 10638e0841e0SMatthew G. Knepley PetscQuadrature quad; 10648e0841e0SMatthew G. Knepley PetscSection coordSection; 10658e0841e0SMatthew G. Knepley Vec coordinates; 10668e0841e0SMatthew G. Knepley PetscScalar *coords = NULL; 10678e0841e0SMatthew G. Knepley const PetscReal *quadPoints; 10688e0841e0SMatthew G. Knepley PetscReal *basisDer; 10698e0841e0SMatthew G. Knepley PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 10708e0841e0SMatthew G. Knepley PetscErrorCode ierr; 10718e0841e0SMatthew G. Knepley 10728e0841e0SMatthew G. Knepley PetscFunctionBegin; 10738e0841e0SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 10748e0841e0SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 10758e0841e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 10768e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 10778e0841e0SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 10788e0841e0SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1079954b1791SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 10808e0841e0SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 10818e0841e0SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 10828e0841e0SMatthew G. Knepley *detJ = 0.0; 10838e0841e0SMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 10848e0841e0SMatthew G. Knepley 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); 10858e0841e0SMatthew G. Knepley if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 10868e0841e0SMatthew G. Knepley if (J) { 10870790e268SMatthew G. Knepley ierr = PetscMemzero(J, Nq*cdim*dim*sizeof(PetscReal));CHKERRQ(ierr); 10888e0841e0SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 10898e0841e0SMatthew G. Knepley PetscInt i, j, k, c, r; 10908e0841e0SMatthew G. Knepley 10918e0841e0SMatthew G. Knepley /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 10928e0841e0SMatthew G. Knepley for (k = 0; k < pdim; ++k) 10938e0841e0SMatthew G. Knepley for (j = 0; j < dim; ++j) 10948e0841e0SMatthew G. Knepley for (i = 0; i < cdim; ++i) 109571d6e60fSMatthew G. Knepley J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 10963bc0b13bSBarry Smith ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr); 10978e0841e0SMatthew G. Knepley if (cdim > dim) { 10988e0841e0SMatthew G. Knepley for (c = dim; c < cdim; ++c) 10998e0841e0SMatthew G. Knepley for (r = 0; r < cdim; ++r) 11008e0841e0SMatthew G. Knepley J[r*cdim+c] = r == c ? 1.0 : 0.0; 11018e0841e0SMatthew G. Knepley } 11028e0841e0SMatthew G. Knepley switch (cdim) { 11038e0841e0SMatthew G. Knepley case 3: 11048e0841e0SMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 11058e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 110617fe8556SMatthew G. Knepley break; 110749dc4407SMatthew G. Knepley case 2: 11088e0841e0SMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 11098e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 111049dc4407SMatthew G. Knepley break; 11118e0841e0SMatthew G. Knepley case 1: 11128e0841e0SMatthew G. Knepley *detJ = J[0]; 11138e0841e0SMatthew G. Knepley if (invJ) invJ[0] = 1.0/J[0]; 111449dc4407SMatthew G. Knepley } 111549dc4407SMatthew G. Knepley } 11168e0841e0SMatthew G. Knepley } 11178e0841e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 11188e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 11198e0841e0SMatthew G. Knepley } 11208e0841e0SMatthew G. Knepley 11218e0841e0SMatthew G. Knepley #undef __FUNCT__ 11228e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 11238e0841e0SMatthew G. Knepley /*@C 11248e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 11258e0841e0SMatthew G. Knepley 11268e0841e0SMatthew G. Knepley Collective on DM 11278e0841e0SMatthew G. Knepley 11288e0841e0SMatthew G. Knepley Input Arguments: 11298e0841e0SMatthew G. Knepley + dm - the DM 11308e0841e0SMatthew G. Knepley . cell - the cell 11318e0841e0SMatthew G. Knepley - fe - the finite element containing the quadrature 11328e0841e0SMatthew G. Knepley 11338e0841e0SMatthew G. Knepley Output Arguments: 11348e0841e0SMatthew G. Knepley + v0 - the translation part of this transform 11358e0841e0SMatthew G. Knepley . J - the Jacobian of the transform from the reference element at each quadrature point 11368e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point 11378e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point 11388e0841e0SMatthew G. Knepley 11398e0841e0SMatthew G. Knepley Level: advanced 11408e0841e0SMatthew G. Knepley 11418e0841e0SMatthew G. Knepley Fortran Notes: 11428e0841e0SMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 11438e0841e0SMatthew G. Knepley include petsc.h90 in your code. 11448e0841e0SMatthew G. Knepley 11458e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 11468e0841e0SMatthew G. Knepley @*/ 11478e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 11488e0841e0SMatthew G. Knepley { 11498e0841e0SMatthew G. Knepley PetscErrorCode ierr; 11508e0841e0SMatthew G. Knepley 11518e0841e0SMatthew G. Knepley PetscFunctionBegin; 11528e0841e0SMatthew G. Knepley if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 11538e0841e0SMatthew G. Knepley else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1154ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1155ccd2543fSMatthew G Knepley } 1156834e62ceSMatthew G. Knepley 1157834e62ceSMatthew G. Knepley #undef __FUNCT__ 1158cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 1159011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1160cc08537eSMatthew G. Knepley { 1161cc08537eSMatthew G. Knepley PetscSection coordSection; 1162cc08537eSMatthew G. Knepley Vec coordinates; 1163a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 116406e2781eSMatthew G. Knepley PetscScalar tmp[2]; 1165cc08537eSMatthew G. Knepley PetscInt coordSize; 1166cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1167cc08537eSMatthew G. Knepley 1168cc08537eSMatthew G. Knepley PetscFunctionBegin; 1169cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 117069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1171cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1172011ea5d8SMatthew G. Knepley if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 11732e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1174cc08537eSMatthew G. Knepley if (centroid) { 117506e2781eSMatthew G. Knepley centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 117606e2781eSMatthew G. Knepley centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1177cc08537eSMatthew G. Knepley } 1178cc08537eSMatthew G. Knepley if (normal) { 1179a60a936bSMatthew G. Knepley PetscReal norm; 1180a60a936bSMatthew G. Knepley 118106e2781eSMatthew G. Knepley normal[0] = -PetscRealPart(coords[1] - tmp[1]); 118206e2781eSMatthew G. Knepley normal[1] = PetscRealPart(coords[0] - tmp[0]); 1183a60a936bSMatthew G. Knepley norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1184a60a936bSMatthew G. Knepley normal[0] /= norm; 1185a60a936bSMatthew G. Knepley normal[1] /= norm; 1186cc08537eSMatthew G. Knepley } 1187cc08537eSMatthew G. Knepley if (vol) { 118806e2781eSMatthew G. Knepley *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1189cc08537eSMatthew G. Knepley } 1190cc08537eSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1191cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1192cc08537eSMatthew G. Knepley } 1193cc08537eSMatthew G. Knepley 1194cc08537eSMatthew G. Knepley #undef __FUNCT__ 1195cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 1196cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1197011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1198cc08537eSMatthew G. Knepley { 1199cc08537eSMatthew G. Knepley PetscSection coordSection; 1200cc08537eSMatthew G. Knepley Vec coordinates; 1201cc08537eSMatthew G. Knepley PetscScalar *coords = NULL; 12020a1d6728SMatthew G. Knepley PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 12030a1d6728SMatthew G. Knepley PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1204cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1205cc08537eSMatthew G. Knepley 1206cc08537eSMatthew G. Knepley PetscFunctionBegin; 1207cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 12080a1d6728SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 120969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1210cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 12110bce18caSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1212ceee4971SMatthew G. Knepley if (dim > 2 && centroid) { 1213ceee4971SMatthew G. Knepley v0[0] = PetscRealPart(coords[0]); 1214ceee4971SMatthew G. Knepley v0[1] = PetscRealPart(coords[1]); 1215ceee4971SMatthew G. Knepley v0[2] = PetscRealPart(coords[2]); 1216ceee4971SMatthew G. Knepley } 1217011ea5d8SMatthew G. Knepley if (normal) { 1218011ea5d8SMatthew G. Knepley if (dim > 2) { 12191ee9d5ecSMatthew G. Knepley const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 12201ee9d5ecSMatthew G. Knepley const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 12211ee9d5ecSMatthew G. Knepley const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 12220a1d6728SMatthew G. Knepley PetscReal norm; 12230a1d6728SMatthew G. Knepley 12240a1d6728SMatthew G. Knepley normal[0] = y0*z1 - z0*y1; 12250a1d6728SMatthew G. Knepley normal[1] = z0*x1 - x0*z1; 12260a1d6728SMatthew G. Knepley normal[2] = x0*y1 - y0*x1; 12278b49ba18SBarry Smith norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 12280a1d6728SMatthew G. Knepley normal[0] /= norm; 12290a1d6728SMatthew G. Knepley normal[1] /= norm; 12300a1d6728SMatthew G. Knepley normal[2] /= norm; 1231011ea5d8SMatthew G. Knepley } else { 1232011ea5d8SMatthew G. Knepley for (d = 0; d < dim; ++d) normal[d] = 0.0; 1233011ea5d8SMatthew G. Knepley } 1234011ea5d8SMatthew G. Knepley } 123599dec3a6SMatthew G. Knepley if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 12360a1d6728SMatthew G. Knepley for (p = 0; p < numCorners; ++p) { 12370a1d6728SMatthew G. Knepley /* Need to do this copy to get types right */ 12380a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12391ee9d5ecSMatthew G. Knepley ctmp[d] = PetscRealPart(coords[p*tdim+d]); 12401ee9d5ecSMatthew G. Knepley ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 12410a1d6728SMatthew G. Knepley } 12420a1d6728SMatthew G. Knepley Volume_Triangle_Origin_Internal(&vtmp, ctmp); 12430a1d6728SMatthew G. Knepley vsum += vtmp; 12440a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12450a1d6728SMatthew G. Knepley csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 12460a1d6728SMatthew G. Knepley } 12470a1d6728SMatthew G. Knepley } 12480a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12490a1d6728SMatthew G. Knepley csum[d] /= (tdim+1)*vsum; 12500a1d6728SMatthew G. Knepley } 12510a1d6728SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1252ee6bbdb2SSatish Balay if (vol) *vol = PetscAbsReal(vsum); 12530a1d6728SMatthew G. Knepley if (centroid) { 12540a1d6728SMatthew G. Knepley if (dim > 2) { 12550a1d6728SMatthew G. Knepley for (d = 0; d < dim; ++d) { 12560a1d6728SMatthew G. Knepley centroid[d] = v0[d]; 12570a1d6728SMatthew G. Knepley for (e = 0; e < dim; ++e) { 12580a1d6728SMatthew G. Knepley centroid[d] += R[d*dim+e]*csum[e]; 12590a1d6728SMatthew G. Knepley } 12600a1d6728SMatthew G. Knepley } 12610a1d6728SMatthew G. Knepley } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 12620a1d6728SMatthew G. Knepley } 1263cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1264cc08537eSMatthew G. Knepley } 1265cc08537eSMatthew G. Knepley 1266cc08537eSMatthew G. Knepley #undef __FUNCT__ 12670ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 12680ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1269011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 12700ec8681fSMatthew G. Knepley { 12710ec8681fSMatthew G. Knepley PetscSection coordSection; 12720ec8681fSMatthew G. Knepley Vec coordinates; 12730ec8681fSMatthew G. Knepley PetscScalar *coords = NULL; 127486623015SMatthew G. Knepley PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1275a7df9edeSMatthew G. Knepley const PetscInt *faces, *facesO; 12760ec8681fSMatthew G. Knepley PetscInt numFaces, f, coordSize, numCorners, p, d; 12770ec8681fSMatthew G. Knepley PetscErrorCode ierr; 12780ec8681fSMatthew G. Knepley 12790ec8681fSMatthew G. Knepley PetscFunctionBegin; 1280f6dae198SJed Brown if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 12810ec8681fSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 128269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 12830ec8681fSMatthew G. Knepley 1284d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 12850ec8681fSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 12860ec8681fSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1287a7df9edeSMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 12880ec8681fSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 1289011ea5d8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 12900ec8681fSMatthew G. Knepley numCorners = coordSize/dim; 12910ec8681fSMatthew G. Knepley switch (numCorners) { 12920ec8681fSMatthew G. Knepley case 3: 12930ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12941ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 12951ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 12961ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 12970ec8681fSMatthew G. Knepley } 12980ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1299a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 13000ec8681fSMatthew G. Knepley vsum += vtmp; 13014f25033aSJed Brown if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 13020ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13031ee9d5ecSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 13040ec8681fSMatthew G. Knepley } 13050ec8681fSMatthew G. Knepley } 13060ec8681fSMatthew G. Knepley break; 13070ec8681fSMatthew G. Knepley case 4: 13080ec8681fSMatthew G. Knepley /* DO FOR PYRAMID */ 13090ec8681fSMatthew G. Knepley /* First tet */ 13100ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13111ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 13121ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 13131ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 13140ec8681fSMatthew G. Knepley } 13150ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1316a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 13170ec8681fSMatthew G. Knepley vsum += vtmp; 13180ec8681fSMatthew G. Knepley if (centroid) { 13190ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13200ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 13210ec8681fSMatthew G. Knepley } 13220ec8681fSMatthew G. Knepley } 13230ec8681fSMatthew G. Knepley /* Second tet */ 13240ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13251ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 13261ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 13271ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 13280ec8681fSMatthew G. Knepley } 13290ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1330a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 13310ec8681fSMatthew G. Knepley vsum += vtmp; 13320ec8681fSMatthew G. Knepley if (centroid) { 13330ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13340ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 13350ec8681fSMatthew G. Knepley } 13360ec8681fSMatthew G. Knepley } 13370ec8681fSMatthew G. Knepley break; 13380ec8681fSMatthew G. Knepley default: 1339796f034aSJed Brown SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 13400ec8681fSMatthew G. Knepley } 13414f25033aSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 13420ec8681fSMatthew G. Knepley } 13438763be8eSMatthew G. Knepley if (vol) *vol = PetscAbsReal(vsum); 13440ec8681fSMatthew G. Knepley if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1345d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 13460ec8681fSMatthew G. Knepley PetscFunctionReturn(0); 13470ec8681fSMatthew G. Knepley } 13480ec8681fSMatthew G. Knepley 13490ec8681fSMatthew G. Knepley #undef __FUNCT__ 1350834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1351834e62ceSMatthew G. Knepley /*@C 1352834e62ceSMatthew G. Knepley DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1353834e62ceSMatthew G. Knepley 1354834e62ceSMatthew G. Knepley Collective on DM 1355834e62ceSMatthew G. Knepley 1356834e62ceSMatthew G. Knepley Input Arguments: 1357834e62ceSMatthew G. Knepley + dm - the DM 1358834e62ceSMatthew G. Knepley - cell - the cell 1359834e62ceSMatthew G. Knepley 1360834e62ceSMatthew G. Knepley Output Arguments: 1361834e62ceSMatthew G. Knepley + volume - the cell volume 1362cc08537eSMatthew G. Knepley . centroid - the cell centroid 1363cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate 1364834e62ceSMatthew G. Knepley 1365834e62ceSMatthew G. Knepley Level: advanced 1366834e62ceSMatthew G. Knepley 1367834e62ceSMatthew G. Knepley Fortran Notes: 1368834e62ceSMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1369834e62ceSMatthew G. Knepley include petsc.h90 in your code. 1370834e62ceSMatthew G. Knepley 137169d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1372834e62ceSMatthew G. Knepley @*/ 1373cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1374834e62ceSMatthew G. Knepley { 13750ec8681fSMatthew G. Knepley PetscInt depth, dim; 1376834e62ceSMatthew G. Knepley PetscErrorCode ierr; 1377834e62ceSMatthew G. Knepley 1378834e62ceSMatthew G. Knepley PetscFunctionBegin; 1379834e62ceSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1380c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1381834e62ceSMatthew G. Knepley if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1382834e62ceSMatthew G. Knepley /* We need to keep a pointer to the depth label */ 1383c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1384834e62ceSMatthew G. Knepley /* Cone size is now the number of faces */ 1385011ea5d8SMatthew G. Knepley switch (depth) { 1386cc08537eSMatthew G. Knepley case 1: 1387011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1388cc08537eSMatthew G. Knepley break; 1389834e62ceSMatthew G. Knepley case 2: 1390011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1391834e62ceSMatthew G. Knepley break; 1392834e62ceSMatthew G. Knepley case 3: 1393011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1394834e62ceSMatthew G. Knepley break; 1395834e62ceSMatthew G. Knepley default: 1396834e62ceSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1397834e62ceSMatthew G. Knepley } 1398834e62ceSMatthew G. Knepley PetscFunctionReturn(0); 1399834e62ceSMatthew G. Knepley } 1400113c68e6SMatthew G. Knepley 1401113c68e6SMatthew G. Knepley #undef __FUNCT__ 1402c0d900a5SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFEM" 1403c0d900a5SMatthew G. Knepley /* This should also take a PetscFE argument I think */ 1404c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1405c0d900a5SMatthew G. Knepley { 1406c0d900a5SMatthew G. Knepley DM dmCell; 1407c0d900a5SMatthew G. Knepley Vec coordinates; 1408c0d900a5SMatthew G. Knepley PetscSection coordSection, sectionCell; 1409c0d900a5SMatthew G. Knepley PetscScalar *cgeom; 1410c0d900a5SMatthew G. Knepley PetscInt cStart, cEnd, cMax, c; 1411c0d900a5SMatthew G. Knepley PetscErrorCode ierr; 1412c0d900a5SMatthew G. Knepley 1413c0d900a5SMatthew G. Knepley PetscFunctionBegin; 1414c0d900a5SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1415c0d900a5SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1416c0d900a5SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1417c0d900a5SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1418c0d900a5SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1419c0d900a5SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1420c0d900a5SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1421c0d900a5SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1422c0d900a5SMatthew G. Knepley cEnd = cMax < 0 ? cEnd : cMax; 1423c0d900a5SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1424c0d900a5SMatthew G. Knepley /* TODO This needs to be multiplied by Nq for non-affine */ 14259e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1426c0d900a5SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1427c0d900a5SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1428c0d900a5SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1429c0d900a5SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1430c0d900a5SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1431c0d900a5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1432c0d900a5SMatthew G. Knepley PetscFECellGeom *cg; 1433c0d900a5SMatthew G. Knepley 1434c0d900a5SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1435c0d900a5SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1436c0d900a5SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1437c0d900a5SMatthew G. Knepley if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1438c0d900a5SMatthew G. Knepley } 1439c0d900a5SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1440c0d900a5SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1441c0d900a5SMatthew G. Knepley PetscFunctionReturn(0); 1442c0d900a5SMatthew G. Knepley } 1443c0d900a5SMatthew G. Knepley 1444c0d900a5SMatthew G. Knepley #undef __FUNCT__ 1445113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM" 1446891a9168SMatthew G. Knepley /*@ 1447891a9168SMatthew G. Knepley DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1448891a9168SMatthew G. Knepley 1449891a9168SMatthew G. Knepley Input Parameter: 1450891a9168SMatthew G. Knepley . dm - The DM 1451891a9168SMatthew G. Knepley 1452891a9168SMatthew G. Knepley Output Parameters: 1453891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data 1454891a9168SMatthew G. Knepley . facegeom - A Vec of PetscFVFaceGeom data 1455891a9168SMatthew G. Knepley 1456891a9168SMatthew G. Knepley Level: developer 1457891a9168SMatthew G. Knepley 1458891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1459891a9168SMatthew G. Knepley @*/ 1460113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1461113c68e6SMatthew G. Knepley { 1462113c68e6SMatthew G. Knepley DM dmFace, dmCell; 1463113c68e6SMatthew G. Knepley DMLabel ghostLabel; 1464113c68e6SMatthew G. Knepley PetscSection sectionFace, sectionCell; 1465113c68e6SMatthew G. Knepley PetscSection coordSection; 1466113c68e6SMatthew G. Knepley Vec coordinates; 1467113c68e6SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1468113c68e6SMatthew G. Knepley PetscReal minradius, gminradius; 1469113c68e6SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1470113c68e6SMatthew G. Knepley PetscErrorCode ierr; 1471113c68e6SMatthew G. Knepley 1472113c68e6SMatthew G. Knepley PetscFunctionBegin; 1473113c68e6SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1474113c68e6SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1475113c68e6SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1476113c68e6SMatthew G. Knepley /* Make cell centroids and volumes */ 1477113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1478113c68e6SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1479113c68e6SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1480113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1481113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1482113c68e6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1483113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 14849e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1485113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1486113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1487113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1488113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 148906348e87SToby Isaac if (cEndInterior < 0) { 149006348e87SToby Isaac cEndInterior = cEnd; 149106348e87SToby Isaac } 1492113c68e6SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1493113c68e6SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 1494113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1495113c68e6SMatthew G. Knepley 1496113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1497113c68e6SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1498113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1499113c68e6SMatthew G. Knepley } 1500113c68e6SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 1501113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1502113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1503113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1504113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 15059e5edeeeSMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1506113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1507113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1508113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1509113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1510113c68e6SMatthew G. Knepley ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1511c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1512113c68e6SMatthew G. Knepley minradius = PETSC_MAX_REAL; 1513113c68e6SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 1514113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1515113c68e6SMatthew G. Knepley PetscReal area; 151650d63984SToby Isaac PetscInt ghost = -1, d, numChildren; 1517113c68e6SMatthew G. Knepley 15189ac3fadcSMatthew G. Knepley if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 151950d63984SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 152050d63984SToby Isaac if (ghost >= 0 || numChildren) continue; 1521113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1522113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1523113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1524113c68e6SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1525113c68e6SMatthew G. Knepley { 1526113c68e6SMatthew G. Knepley PetscFVCellGeom *cL, *cR; 152706348e87SToby Isaac PetscInt ncells; 1528113c68e6SMatthew G. Knepley const PetscInt *cells; 1529113c68e6SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 15300453c0cdSMatthew G. Knepley PetscReal l[3], r[3], v[3]; 1531113c68e6SMatthew G. Knepley 1532113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 153306348e87SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr); 1534113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1535113c68e6SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 153606348e87SToby Isaac if (ncells > 1) { 153706348e87SToby Isaac ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1538113c68e6SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 153906348e87SToby Isaac } 154006348e87SToby Isaac else { 154106348e87SToby Isaac rcentroid = fg->centroid; 154206348e87SToby Isaac } 15432e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 15442e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 15450453c0cdSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1546113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1547113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1548113c68e6SMatthew G. Knepley } 1549113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1550113c68e6SMatthew G. Knepley 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]); 1551113c68e6SMatthew G. Knepley 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]); 1552113c68e6SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1553113c68e6SMatthew G. Knepley } 1554113c68e6SMatthew G. Knepley if (cells[0] < cEndInterior) { 1555113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1556113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1557113c68e6SMatthew G. Knepley } 155806348e87SToby Isaac if (ncells > 1 && cells[1] < cEndInterior) { 1559113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1560113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1561113c68e6SMatthew G. Knepley } 1562113c68e6SMatthew G. Knepley } 1563113c68e6SMatthew G. Knepley } 1564b2566f29SBarry Smith ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1565113c68e6SMatthew G. Knepley ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1566113c68e6SMatthew G. Knepley /* Compute centroids of ghost cells */ 1567113c68e6SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 1568113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1569113c68e6SMatthew G. Knepley const PetscInt *cone, *support; 1570113c68e6SMatthew G. Knepley PetscInt coneSize, supportSize, s; 1571113c68e6SMatthew G. Knepley 1572113c68e6SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1573113c68e6SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1574113c68e6SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1575113c68e6SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 157650d63984SToby Isaac if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize); 1577113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1578113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1579113c68e6SMatthew G. Knepley for (s = 0; s < 2; ++s) { 1580113c68e6SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 1581113c68e6SMatthew G. Knepley if (support[s] == c) { 1582640bce14SSatish Balay PetscFVCellGeom *ci; 1583113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1584113c68e6SMatthew G. Knepley PetscReal c2f[3], a; 1585113c68e6SMatthew G. Knepley 1586113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1587113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1588113c68e6SMatthew G. Knepley a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1589113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1590113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1591113c68e6SMatthew G. Knepley cg->volume = ci->volume; 1592113c68e6SMatthew G. Knepley } 1593113c68e6SMatthew G. Knepley } 1594113c68e6SMatthew G. Knepley } 1595113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1596113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1597113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1598113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1599113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1600113c68e6SMatthew G. Knepley } 1601113c68e6SMatthew G. Knepley 1602113c68e6SMatthew G. Knepley #undef __FUNCT__ 1603113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexGetMinRadius" 1604113c68e6SMatthew G. Knepley /*@C 1605113c68e6SMatthew G. Knepley DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1606113c68e6SMatthew G. Knepley 1607113c68e6SMatthew G. Knepley Not collective 1608113c68e6SMatthew G. Knepley 1609113c68e6SMatthew G. Knepley Input Argument: 1610113c68e6SMatthew G. Knepley . dm - the DM 1611113c68e6SMatthew G. Knepley 1612113c68e6SMatthew G. Knepley Output Argument: 1613113c68e6SMatthew G. Knepley . minradius - the minium cell radius 1614113c68e6SMatthew G. Knepley 1615113c68e6SMatthew G. Knepley Level: developer 1616113c68e6SMatthew G. Knepley 1617113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates() 1618113c68e6SMatthew G. Knepley @*/ 1619113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1620113c68e6SMatthew G. Knepley { 1621113c68e6SMatthew G. Knepley PetscFunctionBegin; 1622113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1623113c68e6SMatthew G. Knepley PetscValidPointer(minradius,2); 1624113c68e6SMatthew G. Knepley *minradius = ((DM_Plex*) dm->data)->minradius; 1625113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1626113c68e6SMatthew G. Knepley } 1627113c68e6SMatthew G. Knepley 1628113c68e6SMatthew G. Knepley #undef __FUNCT__ 1629113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexSetMinRadius" 1630113c68e6SMatthew G. Knepley /*@C 1631113c68e6SMatthew G. Knepley DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1632113c68e6SMatthew G. Knepley 1633113c68e6SMatthew G. Knepley Logically collective 1634113c68e6SMatthew G. Knepley 1635113c68e6SMatthew G. Knepley Input Arguments: 1636113c68e6SMatthew G. Knepley + dm - the DM 1637113c68e6SMatthew G. Knepley - minradius - the minium cell radius 1638113c68e6SMatthew G. Knepley 1639113c68e6SMatthew G. Knepley Level: developer 1640113c68e6SMatthew G. Knepley 1641113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates() 1642113c68e6SMatthew G. Knepley @*/ 1643113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1644113c68e6SMatthew G. Knepley { 1645113c68e6SMatthew G. Knepley PetscFunctionBegin; 1646113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1647113c68e6SMatthew G. Knepley ((DM_Plex*) dm->data)->minradius = minradius; 1648113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1649113c68e6SMatthew G. Knepley } 1650856ac710SMatthew G. Knepley 1651856ac710SMatthew G. Knepley #undef __FUNCT__ 1652856ac710SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction_Internal" 1653856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1654856ac710SMatthew G. Knepley { 1655856ac710SMatthew G. Knepley DMLabel ghostLabel; 1656856ac710SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 1657856ac710SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1658856ac710SMatthew G. Knepley PetscErrorCode ierr; 1659856ac710SMatthew G. Knepley 1660856ac710SMatthew G. Knepley PetscFunctionBegin; 1661856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1662856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1663856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1664856ac710SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1665856ac710SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1666c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1667856ac710SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1668856ac710SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 1669856ac710SMatthew G. Knepley const PetscInt *faces; 1670856ac710SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 1671640bce14SSatish Balay PetscFVCellGeom *cg; 1672856ac710SMatthew G. Knepley PetscBool boundary; 1673856ac710SMatthew G. Knepley PetscInt ghost; 1674856ac710SMatthew G. Knepley 1675856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1676856ac710SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1677856ac710SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1678856ac710SMatthew G. Knepley if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1679856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1680640bce14SSatish Balay PetscFVCellGeom *cg1; 1681856ac710SMatthew G. Knepley PetscFVFaceGeom *fg; 1682856ac710SMatthew G. Knepley const PetscInt *fcells; 1683856ac710SMatthew G. Knepley PetscInt ncell, side; 1684856ac710SMatthew G. Knepley 1685856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1686a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1687856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1688856ac710SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1689856ac710SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1690856ac710SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 1691856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1692856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1693856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1694856ac710SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1695856ac710SMatthew G. Knepley } 1696856ac710SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1697856ac710SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1698856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1699856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1700a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1701856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1702856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1703856ac710SMatthew G. Knepley ++usedFaces; 1704856ac710SMatthew G. Knepley } 1705856ac710SMatthew G. Knepley } 1706856ac710SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1707856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1708856ac710SMatthew G. Knepley } 1709856ac710SMatthew G. Knepley 1710856ac710SMatthew G. Knepley #undef __FUNCT__ 1711b81db932SToby Isaac #define __FUNCT__ "BuildGradientReconstruction_Internal_Tree" 1712b81db932SToby Isaac static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1713b81db932SToby Isaac { 1714b81db932SToby Isaac DMLabel ghostLabel; 1715b81db932SToby Isaac PetscScalar *dx, *grad, **gref; 1716b81db932SToby Isaac PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 1717b81db932SToby Isaac PetscSection neighSec; 1718b81db932SToby Isaac PetscInt (*neighbors)[2]; 1719b81db932SToby Isaac PetscInt *counter; 1720b81db932SToby Isaac PetscErrorCode ierr; 1721b81db932SToby Isaac 1722b81db932SToby Isaac PetscFunctionBegin; 1723b81db932SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1724b81db932SToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1725b81db932SToby Isaac ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 17265bc680faSToby Isaac if (cEndInterior < 0) { 17275bc680faSToby Isaac cEndInterior = cEnd; 17285bc680faSToby Isaac } 1729b81db932SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr); 1730b81db932SToby Isaac ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr); 1731b81db932SToby Isaac ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1732c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1733b81db932SToby Isaac for (f = fStart; f < fEnd; f++) { 1734b81db932SToby Isaac const PetscInt *fcells; 1735b81db932SToby Isaac PetscBool boundary; 17365bc680faSToby Isaac PetscInt ghost = -1; 1737b81db932SToby Isaac PetscInt numChildren, numCells, c; 1738b81db932SToby Isaac 173906348e87SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1740a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1741b81db932SToby Isaac ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1742b81db932SToby Isaac if ((ghost >= 0) || boundary || numChildren) continue; 1743b81db932SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 174406348e87SToby Isaac if (numCells == 2) { 1745b81db932SToby Isaac ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1746b81db932SToby Isaac for (c = 0; c < 2; c++) { 1747b81db932SToby Isaac PetscInt cell = fcells[c]; 1748b81db932SToby Isaac 1749e6885bbbSToby Isaac if (cell >= cStart && cell < cEndInterior) { 1750b81db932SToby Isaac ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr); 1751b81db932SToby Isaac } 1752b81db932SToby Isaac } 1753b81db932SToby Isaac } 175406348e87SToby Isaac } 1755b81db932SToby Isaac ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr); 1756b81db932SToby Isaac ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr); 1757b81db932SToby Isaac ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1758b81db932SToby Isaac nStart = 0; 1759b81db932SToby Isaac ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr); 1760b81db932SToby Isaac ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr); 1761b81db932SToby Isaac ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr); 1762b81db932SToby Isaac for (f = fStart; f < fEnd; f++) { 1763b81db932SToby Isaac const PetscInt *fcells; 1764b81db932SToby Isaac PetscBool boundary; 17655bc680faSToby Isaac PetscInt ghost = -1; 1766b81db932SToby Isaac PetscInt numChildren, numCells, c; 1767b81db932SToby Isaac 176806348e87SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1769a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 1770b81db932SToby Isaac ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 1771b81db932SToby Isaac if ((ghost >= 0) || boundary || numChildren) continue; 1772b81db932SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 177306348e87SToby Isaac if (numCells == 2) { 1774b81db932SToby Isaac ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 1775b81db932SToby Isaac for (c = 0; c < 2; c++) { 1776b81db932SToby Isaac PetscInt cell = fcells[c], off; 1777b81db932SToby Isaac 1778e6885bbbSToby Isaac if (cell >= cStart && cell < cEndInterior) { 1779b81db932SToby Isaac ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr); 1780b81db932SToby Isaac off += counter[cell - cStart]++; 1781b81db932SToby Isaac neighbors[off][0] = f; 1782b81db932SToby Isaac neighbors[off][1] = fcells[1 - c]; 1783b81db932SToby Isaac } 1784b81db932SToby Isaac } 1785b81db932SToby Isaac } 178606348e87SToby Isaac } 1787b81db932SToby Isaac ierr = PetscFree(counter);CHKERRQ(ierr); 1788b81db932SToby Isaac ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1789b81db932SToby Isaac for (c = cStart; c < cEndInterior; c++) { 1790317218b9SToby Isaac PetscInt numFaces, f, d, off, ghost = -1; 1791640bce14SSatish Balay PetscFVCellGeom *cg; 1792b81db932SToby Isaac 1793b81db932SToby Isaac ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1794b81db932SToby Isaac ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr); 1795b81db932SToby Isaac ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr); 1796317218b9SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);} 1797317218b9SToby Isaac 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); 1798b81db932SToby Isaac for (f = 0; f < numFaces; ++f) { 1799640bce14SSatish Balay PetscFVCellGeom *cg1; 1800b81db932SToby Isaac PetscFVFaceGeom *fg; 1801b81db932SToby Isaac const PetscInt *fcells; 1802b81db932SToby Isaac PetscInt ncell, side, nface; 1803b81db932SToby Isaac 1804b81db932SToby Isaac nface = neighbors[off + f][0]; 1805b81db932SToby Isaac ncell = neighbors[off + f][1]; 1806b81db932SToby Isaac ierr = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr); 1807b81db932SToby Isaac side = (c != fcells[0]); 1808b81db932SToby Isaac ierr = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr); 1809b81db932SToby Isaac ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1810b81db932SToby Isaac for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1811b81db932SToby Isaac gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1812b81db932SToby Isaac } 1813b81db932SToby Isaac ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr); 1814b81db932SToby Isaac for (f = 0; f < numFaces; ++f) { 1815b81db932SToby Isaac for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d]; 1816b81db932SToby Isaac } 1817b81db932SToby Isaac } 1818b81db932SToby Isaac ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 18195fe94518SToby Isaac ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr); 1820b81db932SToby Isaac ierr = PetscFree(neighbors);CHKERRQ(ierr); 1821b81db932SToby Isaac PetscFunctionReturn(0); 1822b81db932SToby Isaac } 1823b81db932SToby Isaac 1824b81db932SToby Isaac #undef __FUNCT__ 1825856ac710SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGradientFVM" 1826856ac710SMatthew G. Knepley /*@ 1827856ac710SMatthew G. Knepley DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 1828856ac710SMatthew G. Knepley 1829856ac710SMatthew G. Knepley Collective on DM 1830856ac710SMatthew G. Knepley 1831856ac710SMatthew G. Knepley Input Arguments: 1832856ac710SMatthew G. Knepley + dm - The DM 1833856ac710SMatthew G. Knepley . fvm - The PetscFV 1834856ac710SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 1835856ac710SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 1836856ac710SMatthew G. Knepley 1837856ac710SMatthew G. Knepley Output Parameters: 1838856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted 1839856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data 1840856ac710SMatthew G. Knepley 1841856ac710SMatthew G. Knepley Level: developer 1842856ac710SMatthew G. Knepley 1843856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 1844856ac710SMatthew G. Knepley @*/ 1845856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 1846856ac710SMatthew G. Knepley { 1847856ac710SMatthew G. Knepley DM dmFace, dmCell; 1848856ac710SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1849b81db932SToby Isaac PetscSection sectionGrad, parentSection; 1850856ac710SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 1851856ac710SMatthew G. Knepley PetscErrorCode ierr; 1852856ac710SMatthew G. Knepley 1853856ac710SMatthew G. Knepley PetscFunctionBegin; 1854856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1855856ac710SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1856856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1857856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1858856ac710SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 1859856ac710SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1860856ac710SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1861856ac710SMatthew G. Knepley ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1862856ac710SMatthew G. Knepley ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1863b81db932SToby Isaac ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1864b81db932SToby Isaac if (!parentSection) { 1865856ac710SMatthew G. Knepley ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1866b81db932SToby Isaac } 1867b81db932SToby Isaac else { 1868b81db932SToby Isaac ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1869b81db932SToby Isaac } 1870856ac710SMatthew G. Knepley ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1871856ac710SMatthew G. Knepley ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1872856ac710SMatthew G. Knepley /* Create storage for gradients */ 1873856ac710SMatthew G. Knepley ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 1874856ac710SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 1875856ac710SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 1876856ac710SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 1877856ac710SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 1878856ac710SMatthew G. Knepley ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 1879856ac710SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 1880856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1881856ac710SMatthew G. Knepley } 1882