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; 44*f5ebc837SMatthew 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 56*f5ebc837SMatthew 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}; 323cafe43deSMatthew G. Knepley PetscInt dim, N, cStart, cEnd, 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 */ 338cafe43deSMatthew G. Knepley ierr = MPI_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 339cafe43deSMatthew G. Knepley ierr = MPI_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); 345cafe43deSMatthew G. Knepley ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 346cafe43deSMatthew G. Knepley ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 347722d0f5cSMatthew G. Knepley /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 348cafe43deSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 349cafe43deSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 350722d0f5cSMatthew G. Knepley ierr = PetscCalloc2((dim+1) * dim, &dboxes, dim+1, &boxes);CHKERRQ(ierr); 351cafe43deSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 352cafe43deSMatthew G. Knepley const PetscReal *h = lbox->h; 353cafe43deSMatthew G. Knepley PetscScalar *ccoords = NULL; 354cafe43deSMatthew G. Knepley PetscScalar point[3]; 355cafe43deSMatthew G. Knepley PetscInt dlim[6], d, e, i, j, k; 356cafe43deSMatthew G. Knepley 357cafe43deSMatthew G. Knepley /* Find boxes enclosing each vertex */ 358cafe43deSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 359722d0f5cSMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(lbox, dim+1, ccoords, dboxes, boxes);CHKERRQ(ierr); 360722d0f5cSMatthew G. Knepley /* Mark cells containing the vertices */ 361722d0f5cSMatthew G. Knepley for (e = 0; e < dim+1; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 362cafe43deSMatthew G. Knepley /* Get grid of boxes containing these */ 363cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 3642291669eSMatthew G. Knepley for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 365cafe43deSMatthew G. Knepley for (e = 1; e < dim+1; ++e) { 366cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) { 367cafe43deSMatthew G. Knepley dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 368cafe43deSMatthew G. Knepley dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 369cafe43deSMatthew G. Knepley } 370cafe43deSMatthew G. Knepley } 371fea14342SMatthew G. Knepley /* Check for intersection of box with cell */ 372cafe43deSMatthew 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]) { 373cafe43deSMatthew 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]) { 374cafe43deSMatthew 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]) { 375cafe43deSMatthew G. Knepley const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 376cafe43deSMatthew G. Knepley PetscScalar cpoint[3]; 377fea14342SMatthew G. Knepley PetscInt cell, edge, ii, jj, kk; 378cafe43deSMatthew G. Knepley 379fea14342SMatthew G. Knepley /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 380cafe43deSMatthew G. Knepley for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 381cafe43deSMatthew G. Knepley for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 382cafe43deSMatthew G. Knepley for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 383cafe43deSMatthew G. Knepley 384cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 385cafe43deSMatthew G. Knepley if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 386cafe43deSMatthew G. Knepley } 387cafe43deSMatthew G. Knepley } 388cafe43deSMatthew G. Knepley } 389fea14342SMatthew G. Knepley /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 390fea14342SMatthew G. Knepley for (edge = 0; edge < dim+1; ++edge) { 391fea14342SMatthew G. Knepley PetscReal segA[6], segB[6]; 392fea14342SMatthew G. Knepley 393fea14342SMatthew 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]);} 394fea14342SMatthew G. Knepley for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 3959a128ed2SMatthew G. Knepley if (dim > 2) {segB[2] = PetscRealPart(point[2]); 3969a128ed2SMatthew G. Knepley segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 397fea14342SMatthew G. Knepley for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 3989a128ed2SMatthew G. Knepley if (dim > 1) {segB[1] = PetscRealPart(point[1]); 3999a128ed2SMatthew G. Knepley segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 400fea14342SMatthew G. Knepley for (ii = 0; ii < 2; ++ii) { 401fea14342SMatthew G. Knepley PetscBool intersects; 402fea14342SMatthew G. Knepley 4039a128ed2SMatthew G. Knepley segB[0] = PetscRealPart(point[0]); 4049a128ed2SMatthew G. Knepley segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 405fea14342SMatthew G. Knepley ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 406fea14342SMatthew G. Knepley if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 407cafe43deSMatthew G. Knepley } 408cafe43deSMatthew G. Knepley } 409cafe43deSMatthew G. Knepley } 410cafe43deSMatthew G. Knepley } 411fea14342SMatthew G. Knepley } 412fea14342SMatthew G. Knepley } 413fea14342SMatthew G. Knepley } 414fea14342SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 415fea14342SMatthew G. Knepley } 416722d0f5cSMatthew G. Knepley ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 417cafe43deSMatthew G. Knepley ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 418cafe43deSMatthew G. Knepley ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 419cafe43deSMatthew G. Knepley *localBox = lbox; 420cafe43deSMatthew G. Knepley PetscFunctionReturn(0); 421cafe43deSMatthew G. Knepley } 422cafe43deSMatthew G. Knepley 423cafe43deSMatthew G. Knepley #undef __FUNCT__ 424ccd2543fSMatthew G Knepley #define __FUNCT__ "DMLocatePoints_Plex" 425ccd2543fSMatthew G Knepley /* 426ccd2543fSMatthew G Knepley Need to implement using the guess 427ccd2543fSMatthew G Knepley */ 428ccd2543fSMatthew G Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS) 429ccd2543fSMatthew G Knepley { 430cafe43deSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 431ccd2543fSMatthew G Knepley PetscInt bs, numPoints, p; 4321318edbeSMatthew G. Knepley PetscInt dim, cStart, cEnd, cMax, numCells, c; 433cafe43deSMatthew G. Knepley const PetscInt *boxCells; 434ccd2543fSMatthew G Knepley PetscInt *cells; 435ccd2543fSMatthew G Knepley PetscScalar *a; 436ccd2543fSMatthew G Knepley PetscErrorCode ierr; 437ccd2543fSMatthew G Knepley 438ccd2543fSMatthew G Knepley PetscFunctionBegin; 439cafe43deSMatthew G. Knepley if (!mesh->lbox) {ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 440cafe43deSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 441cafe43deSMatthew G. Knepley ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 442cafe43deSMatthew 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); 443ccd2543fSMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 444ccd2543fSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 445ccd2543fSMatthew G Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 446ccd2543fSMatthew G Knepley ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 447ccd2543fSMatthew G Knepley ierr = VecGetArray(v, &a);CHKERRQ(ierr); 448ccd2543fSMatthew G Knepley numPoints /= bs; 449785e854fSJed Brown ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 450cafe43deSMatthew G. Knepley /* Designate the local box for each point */ 451cafe43deSMatthew G. Knepley /* Send points to correct process */ 452cafe43deSMatthew G. Knepley /* Search cells that lie in each subbox */ 453cafe43deSMatthew G. Knepley /* Should we bin points before doing search? */ 454cafe43deSMatthew G. Knepley ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 455ccd2543fSMatthew G Knepley for (p = 0; p < numPoints; ++p) { 456ccd2543fSMatthew G Knepley const PetscScalar *point = &a[p*bs]; 4571318edbeSMatthew G. Knepley PetscInt dbin[3], bin, cell, cellOffset; 458ccd2543fSMatthew G Knepley 459cafe43deSMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 460cafe43deSMatthew G. Knepley /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 461cafe43deSMatthew G. Knepley ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 462cafe43deSMatthew G. Knepley ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 463cafe43deSMatthew G. Knepley for (c = cellOffset; c < cellOffset + numCells; ++c) { 464cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 465ccd2543fSMatthew G Knepley if (cell >= 0) break; 466ccd2543fSMatthew G Knepley } 4674f6087d8SMatthew G. Knepley if (cell < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D not found in mesh", p); 468ccd2543fSMatthew G Knepley cells[p] = cell; 469ccd2543fSMatthew G Knepley } 470cafe43deSMatthew G. Knepley ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 471cafe43deSMatthew G. Knepley /* Check for highest numbered proc that claims a point (do we care?) */ 472ccd2543fSMatthew G Knepley ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 473ccd2543fSMatthew G Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr); 474ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 475ccd2543fSMatthew G Knepley } 476ccd2543fSMatthew G Knepley 477ccd2543fSMatthew G Knepley #undef __FUNCT__ 47817fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 47917fe8556SMatthew G. Knepley /* 48017fe8556SMatthew G. Knepley DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 48117fe8556SMatthew G. Knepley */ 4823beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 48317fe8556SMatthew G. Knepley { 48417fe8556SMatthew G. Knepley const PetscReal x = PetscRealPart(coords[2] - coords[0]); 48517fe8556SMatthew G. Knepley const PetscReal y = PetscRealPart(coords[3] - coords[1]); 4868b49ba18SBarry Smith const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 48717fe8556SMatthew G. Knepley 48817fe8556SMatthew G. Knepley PetscFunctionBegin; 4891c99cf0cSGeoffrey Irving R[0] = c; R[1] = -s; 4901c99cf0cSGeoffrey Irving R[2] = s; R[3] = c; 49117fe8556SMatthew G. Knepley coords[0] = 0.0; 4927f07f362SMatthew G. Knepley coords[1] = r; 49317fe8556SMatthew G. Knepley PetscFunctionReturn(0); 49417fe8556SMatthew G. Knepley } 49517fe8556SMatthew G. Knepley 49617fe8556SMatthew G. Knepley #undef __FUNCT__ 49728dbe442SToby Isaac #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal" 49828dbe442SToby Isaac /* 49928dbe442SToby Isaac DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D 50028dbe442SToby Isaac 50128dbe442SToby Isaac This uses the basis completion described by Frisvad, 50228dbe442SToby Isaac 50328dbe442SToby Isaac http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html 50428dbe442SToby Isaac DOI:10.1080/2165347X.2012.689606 50528dbe442SToby Isaac */ 5063beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 50728dbe442SToby Isaac { 50828dbe442SToby Isaac PetscReal x = PetscRealPart(coords[3] - coords[0]); 50928dbe442SToby Isaac PetscReal y = PetscRealPart(coords[4] - coords[1]); 51028dbe442SToby Isaac PetscReal z = PetscRealPart(coords[5] - coords[2]); 51128dbe442SToby Isaac PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 51228dbe442SToby Isaac PetscReal rinv = 1. / r; 51328dbe442SToby Isaac PetscFunctionBegin; 51428dbe442SToby Isaac 51528dbe442SToby Isaac x *= rinv; y *= rinv; z *= rinv; 51628dbe442SToby Isaac if (x > 0.) { 51728dbe442SToby Isaac PetscReal inv1pX = 1./ (1. + x); 51828dbe442SToby Isaac 51928dbe442SToby Isaac R[0] = x; R[1] = -y; R[2] = -z; 52028dbe442SToby Isaac R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 52128dbe442SToby Isaac R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 52228dbe442SToby Isaac } 52328dbe442SToby Isaac else { 52428dbe442SToby Isaac PetscReal inv1mX = 1./ (1. - x); 52528dbe442SToby Isaac 52628dbe442SToby Isaac R[0] = x; R[1] = z; R[2] = y; 52728dbe442SToby Isaac R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 52828dbe442SToby Isaac R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 52928dbe442SToby Isaac } 53028dbe442SToby Isaac coords[0] = 0.0; 53128dbe442SToby Isaac coords[1] = r; 53228dbe442SToby Isaac PetscFunctionReturn(0); 53328dbe442SToby Isaac } 53428dbe442SToby Isaac 53528dbe442SToby Isaac #undef __FUNCT__ 536ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 537ccd2543fSMatthew G Knepley /* 538ccd2543fSMatthew G Knepley DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 539ccd2543fSMatthew G Knepley */ 5403beb2758SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 541ccd2543fSMatthew G Knepley { 5421ee9d5ecSMatthew G. Knepley PetscReal x1[3], x2[3], n[3], norm; 54399dec3a6SMatthew G. Knepley PetscReal x1p[3], x2p[3], xnp[3]; 5444a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 545ccd2543fSMatthew G Knepley const PetscInt dim = 3; 54699dec3a6SMatthew G. Knepley PetscInt d, e, p; 547ccd2543fSMatthew G Knepley 548ccd2543fSMatthew G Knepley PetscFunctionBegin; 549ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 550ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 5511ee9d5ecSMatthew G. Knepley x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 5521ee9d5ecSMatthew G. Knepley x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 553ccd2543fSMatthew G Knepley } 554ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 555ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 556ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 5578b49ba18SBarry Smith norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 558ccd2543fSMatthew G Knepley n[0] /= norm; 559ccd2543fSMatthew G Knepley n[1] /= norm; 560ccd2543fSMatthew G Knepley n[2] /= norm; 561ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 562ccd2543fSMatthew G Knepley 563ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 564ccd2543fSMatthew G Knepley 565ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 566ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 567ccd2543fSMatthew G Knepley \ nx ny nz / 568ccd2543fSMatthew G Knepley 569ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 570ccd2543fSMatthew G Knepley */ 5718b49ba18SBarry Smith sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 57273868372SMatthew G. Knepley /* Check for n = z */ 57373868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 5741ee9d5ecSMatthew G. Knepley if (n[2] < 0.0) { 57599dec3a6SMatthew G. Knepley if (coordSize > 9) { 57699dec3a6SMatthew G. Knepley coords[2] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]); 57708b58242SMatthew G. Knepley coords[3] = PetscRealPart(coords[3*dim+1] - coords[0*dim+1]); 57899dec3a6SMatthew G. Knepley coords[4] = x2[0]; 57999dec3a6SMatthew G. Knepley coords[5] = x2[1]; 58099dec3a6SMatthew G. Knepley coords[6] = x1[0]; 58199dec3a6SMatthew G. Knepley coords[7] = x1[1]; 58299dec3a6SMatthew G. Knepley } else { 58373868372SMatthew G. Knepley coords[2] = x2[0]; 58473868372SMatthew G. Knepley coords[3] = x2[1]; 58573868372SMatthew G. Knepley coords[4] = x1[0]; 58673868372SMatthew G. Knepley coords[5] = x1[1]; 58799dec3a6SMatthew G. Knepley } 588b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 589b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 590b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = -1.0; 59173868372SMatthew G. Knepley } else { 59299dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 59399dec3a6SMatthew G. Knepley coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 59499dec3a6SMatthew G. Knepley coords[p*2+1] = PetscRealPart(coords[p*dim+1] - coords[0*dim+1]); 59599dec3a6SMatthew G. Knepley } 59673868372SMatthew G. Knepley coords[2] = x1[0]; 59773868372SMatthew G. Knepley coords[3] = x1[1]; 59873868372SMatthew G. Knepley coords[4] = x2[0]; 59973868372SMatthew G. Knepley coords[5] = x2[1]; 600b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 601b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 602b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = 1.0; 60373868372SMatthew G. Knepley } 60499dec3a6SMatthew G. Knepley coords[0] = 0.0; 60599dec3a6SMatthew G. Knepley coords[1] = 0.0; 60673868372SMatthew G. Knepley PetscFunctionReturn(0); 60773868372SMatthew G. Knepley } 608da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 609ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 610ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 611ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 612ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 613ccd2543fSMatthew G Knepley x1p[d] = 0.0; 614ccd2543fSMatthew G Knepley x2p[d] = 0.0; 615ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 616ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 617ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 618ccd2543fSMatthew G Knepley } 619ccd2543fSMatthew G Knepley } 6208763be8eSMatthew G. Knepley if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 6218763be8eSMatthew G. Knepley if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 622ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 62399dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 62499dec3a6SMatthew G. Knepley for (d = 0; d < dim; ++d) { 62599dec3a6SMatthew G. Knepley xnp[d] = 0.0; 62699dec3a6SMatthew G. Knepley for (e = 0; e < dim; ++e) { 62799dec3a6SMatthew G. Knepley xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 62899dec3a6SMatthew G. Knepley } 62999dec3a6SMatthew G. Knepley if (d < dim-1) coords[p*2+d] = xnp[d]; 63099dec3a6SMatthew G. Knepley } 63199dec3a6SMatthew G. Knepley } 632ccd2543fSMatthew G Knepley coords[0] = 0.0; 633ccd2543fSMatthew G Knepley coords[1] = 0.0; 634ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 635ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 636ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 637ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 6387f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 6397f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 6407f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 6417f07f362SMatthew G. Knepley PetscReal tmp; 6427f07f362SMatthew G. Knepley 6437f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 6447f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 6457f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 6467f07f362SMatthew G. Knepley } 6477f07f362SMatthew G. Knepley } 648ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 649ccd2543fSMatthew G Knepley } 650ccd2543fSMatthew G Knepley 651ccd2543fSMatthew G Knepley #undef __FUNCT__ 652834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal" 6536322fe33SJed Brown PETSC_UNUSED 654834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 655834e62ceSMatthew G. Knepley { 656834e62ceSMatthew G. Knepley /* Signed volume is 1/2 the determinant 657834e62ceSMatthew G. Knepley 658834e62ceSMatthew G. Knepley | 1 1 1 | 659834e62ceSMatthew G. Knepley | x0 x1 x2 | 660834e62ceSMatthew G. Knepley | y0 y1 y2 | 661834e62ceSMatthew G. Knepley 662834e62ceSMatthew G. Knepley but if x0,y0 is the origin, we have 663834e62ceSMatthew G. Knepley 664834e62ceSMatthew G. Knepley | x1 x2 | 665834e62ceSMatthew G. Knepley | y1 y2 | 666834e62ceSMatthew G. Knepley */ 667834e62ceSMatthew G. Knepley const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 668834e62ceSMatthew G. Knepley const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 669834e62ceSMatthew G. Knepley PetscReal M[4], detM; 670834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; 67186623015SMatthew G. Knepley M[2] = y1; M[3] = y2; 672923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(&detM, M); 673834e62ceSMatthew G. Knepley *vol = 0.5*detM; 674834e62ceSMatthew G. Knepley PetscLogFlops(5.0); 675834e62ceSMatthew G. Knepley } 676834e62ceSMatthew G. Knepley 677834e62ceSMatthew G. Knepley #undef __FUNCT__ 678834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal" 679834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 680834e62ceSMatthew G. Knepley { 681923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(vol, coords); 682834e62ceSMatthew G. Knepley *vol *= 0.5; 683834e62ceSMatthew G. Knepley } 684834e62ceSMatthew G. Knepley 685834e62ceSMatthew G. Knepley #undef __FUNCT__ 686834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal" 6876322fe33SJed Brown PETSC_UNUSED 688834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 689834e62ceSMatthew G. Knepley { 690834e62ceSMatthew G. Knepley /* Signed volume is 1/6th of the determinant 691834e62ceSMatthew G. Knepley 692834e62ceSMatthew G. Knepley | 1 1 1 1 | 693834e62ceSMatthew G. Knepley | x0 x1 x2 x3 | 694834e62ceSMatthew G. Knepley | y0 y1 y2 y3 | 695834e62ceSMatthew G. Knepley | z0 z1 z2 z3 | 696834e62ceSMatthew G. Knepley 697834e62ceSMatthew G. Knepley but if x0,y0,z0 is the origin, we have 698834e62ceSMatthew G. Knepley 699834e62ceSMatthew G. Knepley | x1 x2 x3 | 700834e62ceSMatthew G. Knepley | y1 y2 y3 | 701834e62ceSMatthew G. Knepley | z1 z2 z3 | 702834e62ceSMatthew G. Knepley */ 703834e62ceSMatthew G. Knepley const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 704834e62ceSMatthew G. Knepley const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 705834e62ceSMatthew G. Knepley const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 706834e62ceSMatthew G. Knepley PetscReal M[9], detM; 707834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; M[2] = x3; 708834e62ceSMatthew G. Knepley M[3] = y1; M[4] = y2; M[5] = y3; 709834e62ceSMatthew G. Knepley M[6] = z1; M[7] = z2; M[8] = z3; 710923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(&detM, M); 711b7ad821dSMatthew G. Knepley *vol = -0.16666666666666666666666*detM; 712834e62ceSMatthew G. Knepley PetscLogFlops(10.0); 713834e62ceSMatthew G. Knepley } 714834e62ceSMatthew G. Knepley 715834e62ceSMatthew G. Knepley #undef __FUNCT__ 7160ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 7170ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 7180ec8681fSMatthew G. Knepley { 719923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(vol, coords); 720b7ad821dSMatthew G. Knepley *vol *= -0.16666666666666666666666; 7210ec8681fSMatthew G. Knepley } 7220ec8681fSMatthew G. Knepley 7230ec8681fSMatthew G. Knepley #undef __FUNCT__ 72417fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 72517fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 72617fe8556SMatthew G. Knepley { 72717fe8556SMatthew G. Knepley PetscSection coordSection; 72817fe8556SMatthew G. Knepley Vec coordinates; 729a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 7307f07f362SMatthew G. Knepley PetscInt numCoords, d; 73117fe8556SMatthew G. Knepley PetscErrorCode ierr; 73217fe8556SMatthew G. Knepley 73317fe8556SMatthew G. Knepley PetscFunctionBegin; 73417fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 73569d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 73617fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 7377f07f362SMatthew G. Knepley *detJ = 0.0; 73828dbe442SToby Isaac if (numCoords == 6) { 73928dbe442SToby Isaac const PetscInt dim = 3; 74028dbe442SToby Isaac PetscReal R[9], J0; 74128dbe442SToby Isaac 74228dbe442SToby Isaac if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 74328dbe442SToby Isaac ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr); 74428dbe442SToby Isaac if (J) { 74528dbe442SToby Isaac J0 = 0.5*PetscRealPart(coords[1]); 74628dbe442SToby Isaac J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 74728dbe442SToby Isaac J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 74828dbe442SToby Isaac J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 74928dbe442SToby Isaac DMPlex_Det3D_Internal(detJ, J); 75028dbe442SToby Isaac } 75128dbe442SToby Isaac if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 75228dbe442SToby Isaac } else if (numCoords == 4) { 7537f07f362SMatthew G. Knepley const PetscInt dim = 2; 7547f07f362SMatthew G. Knepley PetscReal R[4], J0; 7557f07f362SMatthew G. Knepley 7567f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 7577f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 75817fe8556SMatthew G. Knepley if (J) { 7597f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 7607f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 7617f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 762923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 76317fe8556SMatthew G. Knepley } 764923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 7657f07f362SMatthew G. Knepley } else if (numCoords == 2) { 7667f07f362SMatthew G. Knepley const PetscInt dim = 1; 7677f07f362SMatthew G. Knepley 7687f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 7697f07f362SMatthew G. Knepley if (J) { 7707f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 77117fe8556SMatthew G. Knepley *detJ = J[0]; 77217fe8556SMatthew G. Knepley PetscLogFlops(2.0); 77317fe8556SMatthew G. Knepley } 7747f07f362SMatthew G. Knepley if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 775796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 77617fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 77717fe8556SMatthew G. Knepley PetscFunctionReturn(0); 77817fe8556SMatthew G. Knepley } 77917fe8556SMatthew G. Knepley 78017fe8556SMatthew G. Knepley #undef __FUNCT__ 781ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 782ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 783ccd2543fSMatthew G Knepley { 784ccd2543fSMatthew G Knepley PetscSection coordSection; 785ccd2543fSMatthew G Knepley Vec coordinates; 786a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 7877f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 788ccd2543fSMatthew G Knepley PetscErrorCode ierr; 789ccd2543fSMatthew G Knepley 790ccd2543fSMatthew G Knepley PetscFunctionBegin; 791ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 79269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 793ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 7947f07f362SMatthew G. Knepley *detJ = 0.0; 795ccd2543fSMatthew G Knepley if (numCoords == 9) { 7967f07f362SMatthew G. Knepley const PetscInt dim = 3; 7977f07f362SMatthew 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}; 7987f07f362SMatthew G. Knepley 7997f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 80099dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 8017f07f362SMatthew G. Knepley if (J) { 802b7ad821dSMatthew G. Knepley const PetscInt pdim = 2; 803b7ad821dSMatthew G. Knepley 804b7ad821dSMatthew G. Knepley for (d = 0; d < pdim; d++) { 805b7ad821dSMatthew G. Knepley for (f = 0; f < pdim; f++) { 806b7ad821dSMatthew G. Knepley J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 807ccd2543fSMatthew G Knepley } 8087f07f362SMatthew G. Knepley } 8097f07f362SMatthew G. Knepley PetscLogFlops(8.0); 810923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 8117f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 8127f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 8137f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 8147f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 8157f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 8167f07f362SMatthew G. Knepley } 8177f07f362SMatthew G. Knepley } 8187f07f362SMatthew G. Knepley } 8197f07f362SMatthew G. Knepley PetscLogFlops(18.0); 8207f07f362SMatthew G. Knepley } 821923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 8227f07f362SMatthew G. Knepley } else if (numCoords == 6) { 8237f07f362SMatthew G. Knepley const PetscInt dim = 2; 8247f07f362SMatthew G. Knepley 8257f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 826ccd2543fSMatthew G Knepley if (J) { 827ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 828ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 829ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 830ccd2543fSMatthew G Knepley } 831ccd2543fSMatthew G Knepley } 8327f07f362SMatthew G. Knepley PetscLogFlops(8.0); 833923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 834ccd2543fSMatthew G Knepley } 835923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 836796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 837ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 838ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 839ccd2543fSMatthew G Knepley } 840ccd2543fSMatthew G Knepley 841ccd2543fSMatthew G Knepley #undef __FUNCT__ 842ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 843ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 844ccd2543fSMatthew G Knepley { 845ccd2543fSMatthew G Knepley PetscSection coordSection; 846ccd2543fSMatthew G Knepley Vec coordinates; 847a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 84899dec3a6SMatthew G. Knepley PetscInt numCoords, d, f, g; 849ccd2543fSMatthew G Knepley PetscErrorCode ierr; 850ccd2543fSMatthew G Knepley 851ccd2543fSMatthew G Knepley PetscFunctionBegin; 852ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 85369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 85499dec3a6SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 8557f07f362SMatthew G. Knepley *detJ = 0.0; 85699dec3a6SMatthew G. Knepley if (numCoords == 12) { 85799dec3a6SMatthew G. Knepley const PetscInt dim = 3; 85899dec3a6SMatthew 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}; 85999dec3a6SMatthew G. Knepley 86099dec3a6SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 86199dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 86299dec3a6SMatthew G. Knepley if (J) { 86399dec3a6SMatthew G. Knepley const PetscInt pdim = 2; 86499dec3a6SMatthew G. Knepley 86599dec3a6SMatthew G. Knepley for (d = 0; d < pdim; d++) { 86699dec3a6SMatthew G. Knepley J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 86799dec3a6SMatthew G. Knepley J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 86899dec3a6SMatthew G. Knepley } 86999dec3a6SMatthew G. Knepley PetscLogFlops(8.0); 870923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 87199dec3a6SMatthew G. Knepley for (d = 0; d < dim; d++) { 87299dec3a6SMatthew G. Knepley for (f = 0; f < dim; f++) { 87399dec3a6SMatthew G. Knepley J[d*dim+f] = 0.0; 87499dec3a6SMatthew G. Knepley for (g = 0; g < dim; g++) { 87599dec3a6SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 87699dec3a6SMatthew G. Knepley } 87799dec3a6SMatthew G. Knepley } 87899dec3a6SMatthew G. Knepley } 87999dec3a6SMatthew G. Knepley PetscLogFlops(18.0); 88099dec3a6SMatthew G. Knepley } 881923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 88297f1a218SMatthew G. Knepley } else if ((numCoords == 8) || (numCoords == 16)) { 88399dec3a6SMatthew G. Knepley const PetscInt dim = 2; 88499dec3a6SMatthew G. Knepley 8857f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 886ccd2543fSMatthew G Knepley if (J) { 887ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 88899dec3a6SMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 88999dec3a6SMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 890ccd2543fSMatthew G Knepley } 8917f07f362SMatthew G. Knepley PetscLogFlops(8.0); 892923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 893ccd2543fSMatthew G Knepley } 894923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 895796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 89699dec3a6SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 897ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 898ccd2543fSMatthew G Knepley } 899ccd2543fSMatthew G Knepley 900ccd2543fSMatthew G Knepley #undef __FUNCT__ 901ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 902ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 903ccd2543fSMatthew G Knepley { 904ccd2543fSMatthew G Knepley PetscSection coordSection; 905ccd2543fSMatthew G Knepley Vec coordinates; 906a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 907ccd2543fSMatthew G Knepley const PetscInt dim = 3; 90899dec3a6SMatthew G. Knepley PetscInt d; 909ccd2543fSMatthew G Knepley PetscErrorCode ierr; 910ccd2543fSMatthew G Knepley 911ccd2543fSMatthew G Knepley PetscFunctionBegin; 912ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 91369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 914ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 9157f07f362SMatthew G. Knepley *detJ = 0.0; 9167f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 917ccd2543fSMatthew G Knepley if (J) { 918ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 919f0df753eSMatthew G. Knepley /* I orient with outward face normals */ 920f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 921f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 922f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 923ccd2543fSMatthew G Knepley } 9247f07f362SMatthew G. Knepley PetscLogFlops(18.0); 925923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 926ccd2543fSMatthew G Knepley } 927923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 928ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 929ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 930ccd2543fSMatthew G Knepley } 931ccd2543fSMatthew G Knepley 932ccd2543fSMatthew G Knepley #undef __FUNCT__ 933ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 934ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 935ccd2543fSMatthew G Knepley { 936ccd2543fSMatthew G Knepley PetscSection coordSection; 937ccd2543fSMatthew G Knepley Vec coordinates; 938a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 939ccd2543fSMatthew G Knepley const PetscInt dim = 3; 940ccd2543fSMatthew G Knepley PetscInt d; 941ccd2543fSMatthew G Knepley PetscErrorCode ierr; 942ccd2543fSMatthew G Knepley 943ccd2543fSMatthew G Knepley PetscFunctionBegin; 944ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 94569d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 946ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 9477f07f362SMatthew G. Knepley *detJ = 0.0; 9487f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 949ccd2543fSMatthew G Knepley if (J) { 950ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 951f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 952f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 953f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 954ccd2543fSMatthew G Knepley } 9557f07f362SMatthew G. Knepley PetscLogFlops(18.0); 956923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 957ccd2543fSMatthew G Knepley } 958923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 959ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 960ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 961ccd2543fSMatthew G Knepley } 962ccd2543fSMatthew G Knepley 963ccd2543fSMatthew G Knepley #undef __FUNCT__ 9648e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 965ccd2543fSMatthew G Knepley /*@C 9668e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 967ccd2543fSMatthew G Knepley 968ccd2543fSMatthew G Knepley Collective on DM 969ccd2543fSMatthew G Knepley 970ccd2543fSMatthew G Knepley Input Arguments: 971ccd2543fSMatthew G Knepley + dm - the DM 972ccd2543fSMatthew G Knepley - cell - the cell 973ccd2543fSMatthew G Knepley 974ccd2543fSMatthew G Knepley Output Arguments: 975ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 976ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 977ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 978ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 979ccd2543fSMatthew G Knepley 980ccd2543fSMatthew G Knepley Level: advanced 981ccd2543fSMatthew G Knepley 982ccd2543fSMatthew G Knepley Fortran Notes: 983ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 984ccd2543fSMatthew G Knepley include petsc.h90 in your code. 985ccd2543fSMatthew G Knepley 9868e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 987ccd2543fSMatthew G Knepley @*/ 9888e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 989ccd2543fSMatthew G Knepley { 99049dc4407SMatthew G. Knepley PetscInt depth, dim, coneSize; 991ccd2543fSMatthew G Knepley PetscErrorCode ierr; 992ccd2543fSMatthew G Knepley 993ccd2543fSMatthew G Knepley PetscFunctionBegin; 994139a35ccSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 995ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 99649dc4407SMatthew G. Knepley if (depth == 1) { 9978e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 9988e0841e0SMatthew G. Knepley } else { 9998e0841e0SMatthew G. Knepley DMLabel depth; 10008e0841e0SMatthew G. Knepley 10018e0841e0SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 10028e0841e0SMatthew G. Knepley ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 10038e0841e0SMatthew G. Knepley } 1004ccd2543fSMatthew G Knepley switch (dim) { 100517fe8556SMatthew G. Knepley case 1: 100617fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 100717fe8556SMatthew G. Knepley break; 1008ccd2543fSMatthew G Knepley case 2: 1009ccd2543fSMatthew G Knepley switch (coneSize) { 1010ccd2543fSMatthew G Knepley case 3: 1011ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1012ccd2543fSMatthew G Knepley break; 1013ccd2543fSMatthew G Knepley case 4: 1014ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1015ccd2543fSMatthew G Knepley break; 1016ccd2543fSMatthew G Knepley default: 10178e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1018ccd2543fSMatthew G Knepley } 1019ccd2543fSMatthew G Knepley break; 1020ccd2543fSMatthew G Knepley case 3: 1021ccd2543fSMatthew G Knepley switch (coneSize) { 1022ccd2543fSMatthew G Knepley case 4: 1023ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1024ccd2543fSMatthew G Knepley break; 10258e0841e0SMatthew G. Knepley case 6: /* Faces */ 10268e0841e0SMatthew G. Knepley case 8: /* Vertices */ 1027ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 1028ccd2543fSMatthew G Knepley break; 1029ccd2543fSMatthew G Knepley default: 10308e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1031ccd2543fSMatthew G Knepley } 1032ccd2543fSMatthew G Knepley break; 1033ccd2543fSMatthew G Knepley default: 1034ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1035ccd2543fSMatthew G Knepley } 10368e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 10378e0841e0SMatthew G. Knepley } 10388e0841e0SMatthew G. Knepley 10398e0841e0SMatthew G. Knepley #undef __FUNCT__ 10408e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 10418e0841e0SMatthew G. Knepley static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 10428e0841e0SMatthew G. Knepley { 10438e0841e0SMatthew G. Knepley PetscQuadrature quad; 10448e0841e0SMatthew G. Knepley PetscSection coordSection; 10458e0841e0SMatthew G. Knepley Vec coordinates; 10468e0841e0SMatthew G. Knepley PetscScalar *coords = NULL; 10478e0841e0SMatthew G. Knepley const PetscReal *quadPoints; 10488e0841e0SMatthew G. Knepley PetscReal *basisDer; 10498e0841e0SMatthew G. Knepley PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 10508e0841e0SMatthew G. Knepley PetscErrorCode ierr; 10518e0841e0SMatthew G. Knepley 10528e0841e0SMatthew G. Knepley PetscFunctionBegin; 10538e0841e0SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 10548e0841e0SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 10558e0841e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 10568e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 10578e0841e0SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 10588e0841e0SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1059954b1791SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 10608e0841e0SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 10618e0841e0SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 10628e0841e0SMatthew G. Knepley *detJ = 0.0; 10638e0841e0SMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 10648e0841e0SMatthew 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); 10658e0841e0SMatthew G. Knepley if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 10668e0841e0SMatthew G. Knepley if (J) { 10678e0841e0SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 10688e0841e0SMatthew G. Knepley PetscInt i, j, k, c, r; 10698e0841e0SMatthew G. Knepley 10708e0841e0SMatthew G. Knepley /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 10718e0841e0SMatthew G. Knepley for (k = 0; k < pdim; ++k) 10728e0841e0SMatthew G. Knepley for (j = 0; j < dim; ++j) 10738e0841e0SMatthew G. Knepley for (i = 0; i < cdim; ++i) 107471d6e60fSMatthew G. Knepley J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 10758e0841e0SMatthew G. Knepley PetscLogFlops(2.0*pdim*dim*cdim); 10768e0841e0SMatthew G. Knepley if (cdim > dim) { 10778e0841e0SMatthew G. Knepley for (c = dim; c < cdim; ++c) 10788e0841e0SMatthew G. Knepley for (r = 0; r < cdim; ++r) 10798e0841e0SMatthew G. Knepley J[r*cdim+c] = r == c ? 1.0 : 0.0; 10808e0841e0SMatthew G. Knepley } 10818e0841e0SMatthew G. Knepley switch (cdim) { 10828e0841e0SMatthew G. Knepley case 3: 10838e0841e0SMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 10848e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 108517fe8556SMatthew G. Knepley break; 108649dc4407SMatthew G. Knepley case 2: 10878e0841e0SMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 10888e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 108949dc4407SMatthew G. Knepley break; 10908e0841e0SMatthew G. Knepley case 1: 10918e0841e0SMatthew G. Knepley *detJ = J[0]; 10928e0841e0SMatthew G. Knepley if (invJ) invJ[0] = 1.0/J[0]; 109349dc4407SMatthew G. Knepley } 109449dc4407SMatthew G. Knepley } 10958e0841e0SMatthew G. Knepley } 10968e0841e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 10978e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 10988e0841e0SMatthew G. Knepley } 10998e0841e0SMatthew G. Knepley 11008e0841e0SMatthew G. Knepley #undef __FUNCT__ 11018e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 11028e0841e0SMatthew G. Knepley /*@C 11038e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 11048e0841e0SMatthew G. Knepley 11058e0841e0SMatthew G. Knepley Collective on DM 11068e0841e0SMatthew G. Knepley 11078e0841e0SMatthew G. Knepley Input Arguments: 11088e0841e0SMatthew G. Knepley + dm - the DM 11098e0841e0SMatthew G. Knepley . cell - the cell 11108e0841e0SMatthew G. Knepley - fe - the finite element containing the quadrature 11118e0841e0SMatthew G. Knepley 11128e0841e0SMatthew G. Knepley Output Arguments: 11138e0841e0SMatthew G. Knepley + v0 - the translation part of this transform 11148e0841e0SMatthew G. Knepley . J - the Jacobian of the transform from the reference element at each quadrature point 11158e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point 11168e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point 11178e0841e0SMatthew G. Knepley 11188e0841e0SMatthew G. Knepley Level: advanced 11198e0841e0SMatthew G. Knepley 11208e0841e0SMatthew G. Knepley Fortran Notes: 11218e0841e0SMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 11228e0841e0SMatthew G. Knepley include petsc.h90 in your code. 11238e0841e0SMatthew G. Knepley 11248e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 11258e0841e0SMatthew G. Knepley @*/ 11268e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 11278e0841e0SMatthew G. Knepley { 11288e0841e0SMatthew G. Knepley PetscErrorCode ierr; 11298e0841e0SMatthew G. Knepley 11308e0841e0SMatthew G. Knepley PetscFunctionBegin; 11318e0841e0SMatthew G. Knepley if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 11328e0841e0SMatthew G. Knepley else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 1133ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1134ccd2543fSMatthew G Knepley } 1135834e62ceSMatthew G. Knepley 1136834e62ceSMatthew G. Knepley #undef __FUNCT__ 1137cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 1138011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1139cc08537eSMatthew G. Knepley { 1140cc08537eSMatthew G. Knepley PetscSection coordSection; 1141cc08537eSMatthew G. Knepley Vec coordinates; 1142a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 114306e2781eSMatthew G. Knepley PetscScalar tmp[2]; 1144cc08537eSMatthew G. Knepley PetscInt coordSize; 1145cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1146cc08537eSMatthew G. Knepley 1147cc08537eSMatthew G. Knepley PetscFunctionBegin; 1148cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 114969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1150cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1151011ea5d8SMatthew G. Knepley if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 115206e2781eSMatthew G. Knepley ierr = DMPlexLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1153cc08537eSMatthew G. Knepley if (centroid) { 115406e2781eSMatthew G. Knepley centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 115506e2781eSMatthew G. Knepley centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1156cc08537eSMatthew G. Knepley } 1157cc08537eSMatthew G. Knepley if (normal) { 1158a60a936bSMatthew G. Knepley PetscReal norm; 1159a60a936bSMatthew G. Knepley 116006e2781eSMatthew G. Knepley normal[0] = -PetscRealPart(coords[1] - tmp[1]); 116106e2781eSMatthew G. Knepley normal[1] = PetscRealPart(coords[0] - tmp[0]); 1162a60a936bSMatthew G. Knepley norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1163a60a936bSMatthew G. Knepley normal[0] /= norm; 1164a60a936bSMatthew G. Knepley normal[1] /= norm; 1165cc08537eSMatthew G. Knepley } 1166cc08537eSMatthew G. Knepley if (vol) { 116706e2781eSMatthew G. Knepley *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1168cc08537eSMatthew G. Knepley } 1169cc08537eSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1170cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1171cc08537eSMatthew G. Knepley } 1172cc08537eSMatthew G. Knepley 1173cc08537eSMatthew G. Knepley #undef __FUNCT__ 1174cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 1175cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1176011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1177cc08537eSMatthew G. Knepley { 1178cc08537eSMatthew G. Knepley PetscSection coordSection; 1179cc08537eSMatthew G. Knepley Vec coordinates; 1180cc08537eSMatthew G. Knepley PetscScalar *coords = NULL; 11810a1d6728SMatthew G. Knepley PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 11820a1d6728SMatthew G. Knepley PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1183cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1184cc08537eSMatthew G. Knepley 1185cc08537eSMatthew G. Knepley PetscFunctionBegin; 1186cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 11870a1d6728SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 118869d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1189cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 11900bce18caSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1191011ea5d8SMatthew G. Knepley if (normal) { 1192011ea5d8SMatthew G. Knepley if (dim > 2) { 11931ee9d5ecSMatthew G. Knepley const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 11941ee9d5ecSMatthew G. Knepley const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 11951ee9d5ecSMatthew G. Knepley const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 11960a1d6728SMatthew G. Knepley PetscReal norm; 11970a1d6728SMatthew G. Knepley 11981ee9d5ecSMatthew G. Knepley v0[0] = PetscRealPart(coords[0]); 11991ee9d5ecSMatthew G. Knepley v0[1] = PetscRealPart(coords[1]); 12001ee9d5ecSMatthew G. Knepley v0[2] = PetscRealPart(coords[2]); 12010a1d6728SMatthew G. Knepley normal[0] = y0*z1 - z0*y1; 12020a1d6728SMatthew G. Knepley normal[1] = z0*x1 - x0*z1; 12030a1d6728SMatthew G. Knepley normal[2] = x0*y1 - y0*x1; 12048b49ba18SBarry Smith norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 12050a1d6728SMatthew G. Knepley normal[0] /= norm; 12060a1d6728SMatthew G. Knepley normal[1] /= norm; 12070a1d6728SMatthew G. Knepley normal[2] /= norm; 1208011ea5d8SMatthew G. Knepley } else { 1209011ea5d8SMatthew G. Knepley for (d = 0; d < dim; ++d) normal[d] = 0.0; 1210011ea5d8SMatthew G. Knepley } 1211011ea5d8SMatthew G. Knepley } 121299dec3a6SMatthew G. Knepley if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 12130a1d6728SMatthew G. Knepley for (p = 0; p < numCorners; ++p) { 12140a1d6728SMatthew G. Knepley /* Need to do this copy to get types right */ 12150a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12161ee9d5ecSMatthew G. Knepley ctmp[d] = PetscRealPart(coords[p*tdim+d]); 12171ee9d5ecSMatthew G. Knepley ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 12180a1d6728SMatthew G. Knepley } 12190a1d6728SMatthew G. Knepley Volume_Triangle_Origin_Internal(&vtmp, ctmp); 12200a1d6728SMatthew G. Knepley vsum += vtmp; 12210a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12220a1d6728SMatthew G. Knepley csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 12230a1d6728SMatthew G. Knepley } 12240a1d6728SMatthew G. Knepley } 12250a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 12260a1d6728SMatthew G. Knepley csum[d] /= (tdim+1)*vsum; 12270a1d6728SMatthew G. Knepley } 12280a1d6728SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1229ee6bbdb2SSatish Balay if (vol) *vol = PetscAbsReal(vsum); 12300a1d6728SMatthew G. Knepley if (centroid) { 12310a1d6728SMatthew G. Knepley if (dim > 2) { 12320a1d6728SMatthew G. Knepley for (d = 0; d < dim; ++d) { 12330a1d6728SMatthew G. Knepley centroid[d] = v0[d]; 12340a1d6728SMatthew G. Knepley for (e = 0; e < dim; ++e) { 12350a1d6728SMatthew G. Knepley centroid[d] += R[d*dim+e]*csum[e]; 12360a1d6728SMatthew G. Knepley } 12370a1d6728SMatthew G. Knepley } 12380a1d6728SMatthew G. Knepley } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 12390a1d6728SMatthew G. Knepley } 1240cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1241cc08537eSMatthew G. Knepley } 1242cc08537eSMatthew G. Knepley 1243cc08537eSMatthew G. Knepley #undef __FUNCT__ 12440ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 12450ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1246011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 12470ec8681fSMatthew G. Knepley { 12480ec8681fSMatthew G. Knepley PetscSection coordSection; 12490ec8681fSMatthew G. Knepley Vec coordinates; 12500ec8681fSMatthew G. Knepley PetscScalar *coords = NULL; 125186623015SMatthew G. Knepley PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1252a7df9edeSMatthew G. Knepley const PetscInt *faces, *facesO; 12530ec8681fSMatthew G. Knepley PetscInt numFaces, f, coordSize, numCorners, p, d; 12540ec8681fSMatthew G. Knepley PetscErrorCode ierr; 12550ec8681fSMatthew G. Knepley 12560ec8681fSMatthew G. Knepley PetscFunctionBegin; 1257f6dae198SJed Brown if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 12580ec8681fSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 125969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 12600ec8681fSMatthew G. Knepley 1261d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 12620ec8681fSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 12630ec8681fSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1264a7df9edeSMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 12650ec8681fSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 1266011ea5d8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 12670ec8681fSMatthew G. Knepley numCorners = coordSize/dim; 12680ec8681fSMatthew G. Knepley switch (numCorners) { 12690ec8681fSMatthew G. Knepley case 3: 12700ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12711ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 12721ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 12731ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 12740ec8681fSMatthew G. Knepley } 12750ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1276a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 12770ec8681fSMatthew G. Knepley vsum += vtmp; 12784f25033aSJed Brown if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 12790ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12801ee9d5ecSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 12810ec8681fSMatthew G. Knepley } 12820ec8681fSMatthew G. Knepley } 12830ec8681fSMatthew G. Knepley break; 12840ec8681fSMatthew G. Knepley case 4: 12850ec8681fSMatthew G. Knepley /* DO FOR PYRAMID */ 12860ec8681fSMatthew G. Knepley /* First tet */ 12870ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12881ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 12891ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 12901ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 12910ec8681fSMatthew G. Knepley } 12920ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1293a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 12940ec8681fSMatthew G. Knepley vsum += vtmp; 12950ec8681fSMatthew G. Knepley if (centroid) { 12960ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12970ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 12980ec8681fSMatthew G. Knepley } 12990ec8681fSMatthew G. Knepley } 13000ec8681fSMatthew G. Knepley /* Second tet */ 13010ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13021ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 13031ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 13041ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 13050ec8681fSMatthew G. Knepley } 13060ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1307a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 13080ec8681fSMatthew G. Knepley vsum += vtmp; 13090ec8681fSMatthew G. Knepley if (centroid) { 13100ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 13110ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 13120ec8681fSMatthew G. Knepley } 13130ec8681fSMatthew G. Knepley } 13140ec8681fSMatthew G. Knepley break; 13150ec8681fSMatthew G. Knepley default: 1316796f034aSJed Brown SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 13170ec8681fSMatthew G. Knepley } 13184f25033aSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 13190ec8681fSMatthew G. Knepley } 13208763be8eSMatthew G. Knepley if (vol) *vol = PetscAbsReal(vsum); 13210ec8681fSMatthew G. Knepley if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1322d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 13230ec8681fSMatthew G. Knepley PetscFunctionReturn(0); 13240ec8681fSMatthew G. Knepley } 13250ec8681fSMatthew G. Knepley 13260ec8681fSMatthew G. Knepley #undef __FUNCT__ 1327834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1328834e62ceSMatthew G. Knepley /*@C 1329834e62ceSMatthew G. Knepley DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1330834e62ceSMatthew G. Knepley 1331834e62ceSMatthew G. Knepley Collective on DM 1332834e62ceSMatthew G. Knepley 1333834e62ceSMatthew G. Knepley Input Arguments: 1334834e62ceSMatthew G. Knepley + dm - the DM 1335834e62ceSMatthew G. Knepley - cell - the cell 1336834e62ceSMatthew G. Knepley 1337834e62ceSMatthew G. Knepley Output Arguments: 1338834e62ceSMatthew G. Knepley + volume - the cell volume 1339cc08537eSMatthew G. Knepley . centroid - the cell centroid 1340cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate 1341834e62ceSMatthew G. Knepley 1342834e62ceSMatthew G. Knepley Level: advanced 1343834e62ceSMatthew G. Knepley 1344834e62ceSMatthew G. Knepley Fortran Notes: 1345834e62ceSMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1346834e62ceSMatthew G. Knepley include petsc.h90 in your code. 1347834e62ceSMatthew G. Knepley 134869d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1349834e62ceSMatthew G. Knepley @*/ 1350cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1351834e62ceSMatthew G. Knepley { 13520ec8681fSMatthew G. Knepley PetscInt depth, dim; 1353834e62ceSMatthew G. Knepley PetscErrorCode ierr; 1354834e62ceSMatthew G. Knepley 1355834e62ceSMatthew G. Knepley PetscFunctionBegin; 1356834e62ceSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1357c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1358834e62ceSMatthew G. Knepley if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1359834e62ceSMatthew G. Knepley /* We need to keep a pointer to the depth label */ 1360011ea5d8SMatthew G. Knepley ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1361834e62ceSMatthew G. Knepley /* Cone size is now the number of faces */ 1362011ea5d8SMatthew G. Knepley switch (depth) { 1363cc08537eSMatthew G. Knepley case 1: 1364011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1365cc08537eSMatthew G. Knepley break; 1366834e62ceSMatthew G. Knepley case 2: 1367011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1368834e62ceSMatthew G. Knepley break; 1369834e62ceSMatthew G. Knepley case 3: 1370011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1371834e62ceSMatthew G. Knepley break; 1372834e62ceSMatthew G. Knepley default: 1373834e62ceSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1374834e62ceSMatthew G. Knepley } 1375834e62ceSMatthew G. Knepley PetscFunctionReturn(0); 1376834e62ceSMatthew G. Knepley } 1377113c68e6SMatthew G. Knepley 1378113c68e6SMatthew G. Knepley #undef __FUNCT__ 1379c0d900a5SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFEM" 1380c0d900a5SMatthew G. Knepley /* This should also take a PetscFE argument I think */ 1381c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1382c0d900a5SMatthew G. Knepley { 1383c0d900a5SMatthew G. Knepley DM dmCell; 1384c0d900a5SMatthew G. Knepley Vec coordinates; 1385c0d900a5SMatthew G. Knepley PetscSection coordSection, sectionCell; 1386c0d900a5SMatthew G. Knepley PetscScalar *cgeom; 1387c0d900a5SMatthew G. Knepley PetscInt cStart, cEnd, cMax, c; 1388c0d900a5SMatthew G. Knepley PetscErrorCode ierr; 1389c0d900a5SMatthew G. Knepley 1390c0d900a5SMatthew G. Knepley PetscFunctionBegin; 1391c0d900a5SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1392c0d900a5SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1393c0d900a5SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1394c0d900a5SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1395c0d900a5SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1396c0d900a5SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1397c0d900a5SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1398c0d900a5SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1399c0d900a5SMatthew G. Knepley cEnd = cMax < 0 ? cEnd : cMax; 1400c0d900a5SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1401c0d900a5SMatthew G. Knepley /* TODO This needs to be multiplied by Nq for non-affine */ 14029e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1403c0d900a5SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1404c0d900a5SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1405c0d900a5SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1406c0d900a5SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1407c0d900a5SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1408c0d900a5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1409c0d900a5SMatthew G. Knepley PetscFECellGeom *cg; 1410c0d900a5SMatthew G. Knepley 1411c0d900a5SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1412c0d900a5SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1413c0d900a5SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1414c0d900a5SMatthew G. Knepley if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1415c0d900a5SMatthew G. Knepley } 1416c0d900a5SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1417c0d900a5SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1418c0d900a5SMatthew G. Knepley PetscFunctionReturn(0); 1419c0d900a5SMatthew G. Knepley } 1420c0d900a5SMatthew G. Knepley 1421c0d900a5SMatthew G. Knepley #undef __FUNCT__ 1422113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM" 1423891a9168SMatthew G. Knepley /*@ 1424891a9168SMatthew G. Knepley DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1425891a9168SMatthew G. Knepley 1426891a9168SMatthew G. Knepley Input Parameter: 1427891a9168SMatthew G. Knepley . dm - The DM 1428891a9168SMatthew G. Knepley 1429891a9168SMatthew G. Knepley Output Parameters: 1430891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data 1431891a9168SMatthew G. Knepley . facegeom - A Vec of PetscFVFaceGeom data 1432891a9168SMatthew G. Knepley 1433891a9168SMatthew G. Knepley Level: developer 1434891a9168SMatthew G. Knepley 1435891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1436891a9168SMatthew G. Knepley @*/ 1437113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1438113c68e6SMatthew G. Knepley { 1439113c68e6SMatthew G. Knepley DM dmFace, dmCell; 1440113c68e6SMatthew G. Knepley DMLabel ghostLabel; 1441113c68e6SMatthew G. Knepley PetscSection sectionFace, sectionCell; 1442113c68e6SMatthew G. Knepley PetscSection coordSection; 1443113c68e6SMatthew G. Knepley Vec coordinates; 1444113c68e6SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1445113c68e6SMatthew G. Knepley PetscReal minradius, gminradius; 1446113c68e6SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1447113c68e6SMatthew G. Knepley PetscErrorCode ierr; 1448113c68e6SMatthew G. Knepley 1449113c68e6SMatthew G. Knepley PetscFunctionBegin; 1450113c68e6SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1451113c68e6SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1452113c68e6SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1453113c68e6SMatthew G. Knepley /* Make cell centroids and volumes */ 1454113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1455113c68e6SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1456113c68e6SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1457113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1458113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1459113c68e6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1460113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 14619e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1462113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1463113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1464113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1465113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1466113c68e6SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1467113c68e6SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 1468113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1469113c68e6SMatthew G. Knepley 1470113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1471113c68e6SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1472113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1473113c68e6SMatthew G. Knepley } 1474113c68e6SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 1475113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1476113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1477113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1478113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 14799e5edeeeSMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1480113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1481113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1482113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1483113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1484113c68e6SMatthew G. Knepley ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1485113c68e6SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1486113c68e6SMatthew G. Knepley minradius = PETSC_MAX_REAL; 1487113c68e6SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 1488113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1489113c68e6SMatthew G. Knepley PetscReal area; 14909ac3fadcSMatthew G. Knepley PetscInt ghost = -1, d; 1491113c68e6SMatthew G. Knepley 14929ac3fadcSMatthew G. Knepley if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1493113c68e6SMatthew G. Knepley if (ghost >= 0) continue; 1494113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1495113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1496113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1497113c68e6SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1498113c68e6SMatthew G. Knepley { 1499113c68e6SMatthew G. Knepley PetscFVCellGeom *cL, *cR; 1500113c68e6SMatthew G. Knepley const PetscInt *cells; 1501113c68e6SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 15020453c0cdSMatthew G. Knepley PetscReal l[3], r[3], v[3]; 1503113c68e6SMatthew G. Knepley 1504113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 1505113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1506113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1507113c68e6SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 1508113c68e6SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 1509f170fd80SMatthew G. Knepley ierr = DMPlexLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 1510f170fd80SMatthew G. Knepley ierr = DMPlexLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 15110453c0cdSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1512113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1513113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1514113c68e6SMatthew G. Knepley } 1515113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1516113c68e6SMatthew 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]); 1517113c68e6SMatthew 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]); 1518113c68e6SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1519113c68e6SMatthew G. Knepley } 1520113c68e6SMatthew G. Knepley if (cells[0] < cEndInterior) { 1521113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1522113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1523113c68e6SMatthew G. Knepley } 1524113c68e6SMatthew G. Knepley if (cells[1] < cEndInterior) { 1525113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1526113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1527113c68e6SMatthew G. Knepley } 1528113c68e6SMatthew G. Knepley } 1529113c68e6SMatthew G. Knepley } 1530a9b180a6SBarry Smith ierr = MPI_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1531113c68e6SMatthew G. Knepley ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1532113c68e6SMatthew G. Knepley /* Compute centroids of ghost cells */ 1533113c68e6SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 1534113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1535113c68e6SMatthew G. Knepley const PetscInt *cone, *support; 1536113c68e6SMatthew G. Knepley PetscInt coneSize, supportSize, s; 1537113c68e6SMatthew G. Knepley 1538113c68e6SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1539113c68e6SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1540113c68e6SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1541113c68e6SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 1542113c68e6SMatthew G. Knepley if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize); 1543113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1544113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1545113c68e6SMatthew G. Knepley for (s = 0; s < 2; ++s) { 1546113c68e6SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 1547113c68e6SMatthew G. Knepley if (support[s] == c) { 1548113c68e6SMatthew G. Knepley const PetscFVCellGeom *ci; 1549113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1550113c68e6SMatthew G. Knepley PetscReal c2f[3], a; 1551113c68e6SMatthew G. Knepley 1552113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1553113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1554113c68e6SMatthew G. Knepley a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1555113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1556113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1557113c68e6SMatthew G. Knepley cg->volume = ci->volume; 1558113c68e6SMatthew G. Knepley } 1559113c68e6SMatthew G. Knepley } 1560113c68e6SMatthew G. Knepley } 1561113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1562113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1563113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1564113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1565113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1566113c68e6SMatthew G. Knepley } 1567113c68e6SMatthew G. Knepley 1568113c68e6SMatthew G. Knepley #undef __FUNCT__ 1569113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexGetMinRadius" 1570113c68e6SMatthew G. Knepley /*@C 1571113c68e6SMatthew G. Knepley DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1572113c68e6SMatthew G. Knepley 1573113c68e6SMatthew G. Knepley Not collective 1574113c68e6SMatthew G. Knepley 1575113c68e6SMatthew G. Knepley Input Argument: 1576113c68e6SMatthew G. Knepley . dm - the DM 1577113c68e6SMatthew G. Knepley 1578113c68e6SMatthew G. Knepley Output Argument: 1579113c68e6SMatthew G. Knepley . minradius - the minium cell radius 1580113c68e6SMatthew G. Knepley 1581113c68e6SMatthew G. Knepley Level: developer 1582113c68e6SMatthew G. Knepley 1583113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates() 1584113c68e6SMatthew G. Knepley @*/ 1585113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1586113c68e6SMatthew G. Knepley { 1587113c68e6SMatthew G. Knepley PetscFunctionBegin; 1588113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1589113c68e6SMatthew G. Knepley PetscValidPointer(minradius,2); 1590113c68e6SMatthew G. Knepley *minradius = ((DM_Plex*) dm->data)->minradius; 1591113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1592113c68e6SMatthew G. Knepley } 1593113c68e6SMatthew G. Knepley 1594113c68e6SMatthew G. Knepley #undef __FUNCT__ 1595113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexSetMinRadius" 1596113c68e6SMatthew G. Knepley /*@C 1597113c68e6SMatthew G. Knepley DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1598113c68e6SMatthew G. Knepley 1599113c68e6SMatthew G. Knepley Logically collective 1600113c68e6SMatthew G. Knepley 1601113c68e6SMatthew G. Knepley Input Arguments: 1602113c68e6SMatthew G. Knepley + dm - the DM 1603113c68e6SMatthew G. Knepley - minradius - the minium cell radius 1604113c68e6SMatthew G. Knepley 1605113c68e6SMatthew G. Knepley Level: developer 1606113c68e6SMatthew G. Knepley 1607113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates() 1608113c68e6SMatthew G. Knepley @*/ 1609113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1610113c68e6SMatthew G. Knepley { 1611113c68e6SMatthew G. Knepley PetscFunctionBegin; 1612113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1613113c68e6SMatthew G. Knepley ((DM_Plex*) dm->data)->minradius = minradius; 1614113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1615113c68e6SMatthew G. Knepley } 1616856ac710SMatthew G. Knepley 1617856ac710SMatthew G. Knepley #undef __FUNCT__ 1618856ac710SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction_Internal" 1619856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1620856ac710SMatthew G. Knepley { 1621856ac710SMatthew G. Knepley DMLabel ghostLabel; 1622856ac710SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 1623856ac710SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1624856ac710SMatthew G. Knepley PetscErrorCode ierr; 1625856ac710SMatthew G. Knepley 1626856ac710SMatthew G. Knepley PetscFunctionBegin; 1627856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1628856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1629856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1630856ac710SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1631856ac710SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1632856ac710SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1633856ac710SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1634856ac710SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 1635856ac710SMatthew G. Knepley const PetscInt *faces; 1636856ac710SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 1637856ac710SMatthew G. Knepley const PetscFVCellGeom *cg; 1638856ac710SMatthew G. Knepley PetscBool boundary; 1639856ac710SMatthew G. Knepley PetscInt ghost; 1640856ac710SMatthew G. Knepley 1641856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1642856ac710SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1643856ac710SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1644856ac710SMatthew 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); 1645856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1646856ac710SMatthew G. Knepley const PetscFVCellGeom *cg1; 1647856ac710SMatthew G. Knepley PetscFVFaceGeom *fg; 1648856ac710SMatthew G. Knepley const PetscInt *fcells; 1649856ac710SMatthew G. Knepley PetscInt ncell, side; 1650856ac710SMatthew G. Knepley 1651856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1652856ac710SMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1653856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1654856ac710SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1655856ac710SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1656856ac710SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 1657856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1658856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1659856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1660856ac710SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1661856ac710SMatthew G. Knepley } 1662856ac710SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1663856ac710SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1664856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1665856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1666856ac710SMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1667856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1668856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1669856ac710SMatthew G. Knepley ++usedFaces; 1670856ac710SMatthew G. Knepley } 1671856ac710SMatthew G. Knepley } 1672856ac710SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1673856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1674856ac710SMatthew G. Knepley } 1675856ac710SMatthew G. Knepley 1676856ac710SMatthew G. Knepley #undef __FUNCT__ 1677856ac710SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGradientFVM" 1678856ac710SMatthew G. Knepley /*@ 1679856ac710SMatthew G. Knepley DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 1680856ac710SMatthew G. Knepley 1681856ac710SMatthew G. Knepley Collective on DM 1682856ac710SMatthew G. Knepley 1683856ac710SMatthew G. Knepley Input Arguments: 1684856ac710SMatthew G. Knepley + dm - The DM 1685856ac710SMatthew G. Knepley . fvm - The PetscFV 1686856ac710SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 1687856ac710SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 1688856ac710SMatthew G. Knepley 1689856ac710SMatthew G. Knepley Output Parameters: 1690856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted 1691856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data 1692856ac710SMatthew G. Knepley 1693856ac710SMatthew G. Knepley Level: developer 1694856ac710SMatthew G. Knepley 1695856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 1696856ac710SMatthew G. Knepley @*/ 1697856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 1698856ac710SMatthew G. Knepley { 1699856ac710SMatthew G. Knepley DM dmFace, dmCell; 1700856ac710SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1701856ac710SMatthew G. Knepley PetscSection sectionGrad; 1702856ac710SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 1703856ac710SMatthew G. Knepley PetscErrorCode ierr; 1704856ac710SMatthew G. Knepley 1705856ac710SMatthew G. Knepley PetscFunctionBegin; 1706856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1707856ac710SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1708856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1709856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1710856ac710SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 1711856ac710SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1712856ac710SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1713856ac710SMatthew G. Knepley ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1714856ac710SMatthew G. Knepley ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1715856ac710SMatthew G. Knepley ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1716856ac710SMatthew G. Knepley ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1717856ac710SMatthew G. Knepley ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1718856ac710SMatthew G. Knepley /* Create storage for gradients */ 1719856ac710SMatthew G. Knepley ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 1720856ac710SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 1721856ac710SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 1722856ac710SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 1723856ac710SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 1724856ac710SMatthew G. Knepley ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 1725856ac710SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 1726856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1727856ac710SMatthew G. Knepley } 1728