1ccd2543fSMatthew G Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2ccd2543fSMatthew G Knepley 3ccd2543fSMatthew G Knepley #undef __FUNCT__ 4ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_2D_Internal" 5ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 6ccd2543fSMatthew G Knepley { 7ccd2543fSMatthew G Knepley const PetscInt embedDim = 2; 8ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 9ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 10ccd2543fSMatthew G Knepley PetscReal v0[2], J[4], invJ[4], detJ; 11ccd2543fSMatthew G Knepley PetscReal xi, eta; 12ccd2543fSMatthew G Knepley PetscErrorCode ierr; 13ccd2543fSMatthew G Knepley 14ccd2543fSMatthew G Knepley PetscFunctionBegin; 158e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 16ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 17ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 18ccd2543fSMatthew G Knepley 19ccd2543fSMatthew G Knepley if ((xi >= 0.0) && (eta >= 0.0) && (xi + eta <= 2.0)) *cell = c; 20ccd2543fSMatthew G Knepley else *cell = -1; 21ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 22ccd2543fSMatthew G Knepley } 23ccd2543fSMatthew G Knepley 24ccd2543fSMatthew G Knepley #undef __FUNCT__ 25ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_2D_Internal" 26ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 27ccd2543fSMatthew G Knepley { 28ccd2543fSMatthew G Knepley PetscSection coordSection; 29ccd2543fSMatthew G Knepley Vec coordsLocal; 30a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 31ccd2543fSMatthew G Knepley const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 32ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 33ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 34ccd2543fSMatthew G Knepley PetscInt crossings = 0, f; 35ccd2543fSMatthew G Knepley PetscErrorCode ierr; 36ccd2543fSMatthew G Knepley 37ccd2543fSMatthew G Knepley PetscFunctionBegin; 38ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 3969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 40ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 41ccd2543fSMatthew G Knepley for (f = 0; f < 4; ++f) { 42ccd2543fSMatthew G Knepley PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 43ccd2543fSMatthew G Knepley PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 44ccd2543fSMatthew G Knepley PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 45ccd2543fSMatthew G Knepley PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 46ccd2543fSMatthew G Knepley PetscReal slope = (y_j - y_i) / (x_j - x_i); 47ccd2543fSMatthew G Knepley PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 48ccd2543fSMatthew G Knepley PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 49ccd2543fSMatthew G Knepley PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 50ccd2543fSMatthew G Knepley if ((cond1 || cond2) && above) ++crossings; 51ccd2543fSMatthew G Knepley } 52ccd2543fSMatthew G Knepley if (crossings % 2) *cell = c; 53ccd2543fSMatthew G Knepley else *cell = -1; 54ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 55ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 56ccd2543fSMatthew G Knepley } 57ccd2543fSMatthew G Knepley 58ccd2543fSMatthew G Knepley #undef __FUNCT__ 59ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_Simplex_3D_Internal" 60ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 61ccd2543fSMatthew G Knepley { 62ccd2543fSMatthew G Knepley const PetscInt embedDim = 3; 63ccd2543fSMatthew G Knepley PetscReal v0[3], J[9], invJ[9], detJ; 64ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 65ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 66ccd2543fSMatthew G Knepley PetscReal z = PetscRealPart(point[2]); 67ccd2543fSMatthew G Knepley PetscReal xi, eta, zeta; 68ccd2543fSMatthew G Knepley PetscErrorCode ierr; 69ccd2543fSMatthew G Knepley 70ccd2543fSMatthew G Knepley PetscFunctionBegin; 718e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 72ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 73ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 74ccd2543fSMatthew G Knepley zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 75ccd2543fSMatthew G Knepley 76ccd2543fSMatthew G Knepley if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 77ccd2543fSMatthew G Knepley else *cell = -1; 78ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 79ccd2543fSMatthew G Knepley } 80ccd2543fSMatthew G Knepley 81ccd2543fSMatthew G Knepley #undef __FUNCT__ 82ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexLocatePoint_General_3D_Internal" 83ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 84ccd2543fSMatthew G Knepley { 85ccd2543fSMatthew G Knepley PetscSection coordSection; 86ccd2543fSMatthew G Knepley Vec coordsLocal; 877c1f9639SMatthew G Knepley PetscScalar *coords; 88fb150da6SMatthew G. Knepley const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 89fb150da6SMatthew G. Knepley 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 90ccd2543fSMatthew G Knepley PetscBool found = PETSC_TRUE; 91ccd2543fSMatthew G Knepley PetscInt f; 92ccd2543fSMatthew G Knepley PetscErrorCode ierr; 93ccd2543fSMatthew G Knepley 94ccd2543fSMatthew G Knepley PetscFunctionBegin; 95ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 9669d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 97ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 98ccd2543fSMatthew G Knepley for (f = 0; f < 6; ++f) { 99ccd2543fSMatthew G Knepley /* Check the point is under plane */ 100ccd2543fSMatthew G Knepley /* Get face normal */ 101ccd2543fSMatthew G Knepley PetscReal v_i[3]; 102ccd2543fSMatthew G Knepley PetscReal v_j[3]; 103ccd2543fSMatthew G Knepley PetscReal normal[3]; 104ccd2543fSMatthew G Knepley PetscReal pp[3]; 105ccd2543fSMatthew G Knepley PetscReal dot; 106ccd2543fSMatthew G Knepley 107ccd2543fSMatthew G Knepley v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 108ccd2543fSMatthew G Knepley v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 109ccd2543fSMatthew G Knepley v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 110ccd2543fSMatthew G Knepley v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 111ccd2543fSMatthew G Knepley v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 112ccd2543fSMatthew G Knepley v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 113ccd2543fSMatthew G Knepley normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 114ccd2543fSMatthew G Knepley normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 115ccd2543fSMatthew G Knepley normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 116ccd2543fSMatthew G Knepley pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 117ccd2543fSMatthew G Knepley pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 118ccd2543fSMatthew G Knepley pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 119ccd2543fSMatthew G Knepley dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 120ccd2543fSMatthew G Knepley 121ccd2543fSMatthew G Knepley /* Check that projected point is in face (2D location problem) */ 122ccd2543fSMatthew G Knepley if (dot < 0.0) { 123ccd2543fSMatthew G Knepley found = PETSC_FALSE; 124ccd2543fSMatthew G Knepley break; 125ccd2543fSMatthew G Knepley } 126ccd2543fSMatthew G Knepley } 127ccd2543fSMatthew G Knepley if (found) *cell = c; 128ccd2543fSMatthew G Knepley else *cell = -1; 129ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 130ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 131ccd2543fSMatthew G Knepley } 132ccd2543fSMatthew G Knepley 133ccd2543fSMatthew G Knepley #undef __FUNCT__ 134ccd2543fSMatthew G Knepley #define __FUNCT__ "DMLocatePoints_Plex" 135ccd2543fSMatthew G Knepley /* 136ccd2543fSMatthew G Knepley Need to implement using the guess 137ccd2543fSMatthew G Knepley */ 138ccd2543fSMatthew G Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, IS *cellIS) 139ccd2543fSMatthew G Knepley { 140ccd2543fSMatthew G Knepley PetscInt cell = -1 /*, guess = -1*/; 141ccd2543fSMatthew G Knepley PetscInt bs, numPoints, p; 142ccd2543fSMatthew G Knepley PetscInt dim, cStart, cEnd, cMax, c, coneSize; 143ccd2543fSMatthew G Knepley PetscInt *cells; 144ccd2543fSMatthew G Knepley PetscScalar *a; 145ccd2543fSMatthew G Knepley PetscErrorCode ierr; 146ccd2543fSMatthew G Knepley 147ccd2543fSMatthew G Knepley PetscFunctionBegin; 148c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 149ccd2543fSMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 150ccd2543fSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 151ccd2543fSMatthew G Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 152ccd2543fSMatthew G Knepley ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 153ccd2543fSMatthew G Knepley ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 154ccd2543fSMatthew G Knepley ierr = VecGetArray(v, &a);CHKERRQ(ierr); 155796f034aSJed Brown 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); 156ccd2543fSMatthew G Knepley numPoints /= bs; 157785e854fSJed Brown ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 158ccd2543fSMatthew G Knepley for (p = 0; p < numPoints; ++p) { 159ccd2543fSMatthew G Knepley const PetscScalar *point = &a[p*bs]; 160ccd2543fSMatthew G Knepley 161ccd2543fSMatthew G Knepley switch (dim) { 162ccd2543fSMatthew G Knepley case 2: 163ccd2543fSMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 164ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 165ccd2543fSMatthew G Knepley switch (coneSize) { 166ccd2543fSMatthew G Knepley case 3: 167ccd2543fSMatthew G Knepley ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 168ccd2543fSMatthew G Knepley break; 169ccd2543fSMatthew G Knepley case 4: 170ccd2543fSMatthew G Knepley ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 171ccd2543fSMatthew G Knepley break; 172ccd2543fSMatthew G Knepley default: 173796f034aSJed Brown SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 174ccd2543fSMatthew G Knepley } 175ccd2543fSMatthew G Knepley if (cell >= 0) break; 176ccd2543fSMatthew G Knepley } 177ccd2543fSMatthew G Knepley break; 178ccd2543fSMatthew G Knepley case 3: 179ccd2543fSMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 180ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr); 181ccd2543fSMatthew G Knepley switch (coneSize) { 182ccd2543fSMatthew G Knepley case 4: 183ccd2543fSMatthew G Knepley ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 184ccd2543fSMatthew G Knepley break; 18535953a02SMatthew G. Knepley case 6: 186ccd2543fSMatthew G Knepley ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 187ccd2543fSMatthew G Knepley break; 188ccd2543fSMatthew G Knepley default: 189796f034aSJed Brown SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 190ccd2543fSMatthew G Knepley } 191ccd2543fSMatthew G Knepley if (cell >= 0) break; 192ccd2543fSMatthew G Knepley } 193ccd2543fSMatthew G Knepley break; 194ccd2543fSMatthew G Knepley default: 195796f034aSJed Brown SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 196ccd2543fSMatthew G Knepley } 197ccd2543fSMatthew G Knepley cells[p] = cell; 198ccd2543fSMatthew G Knepley } 199ccd2543fSMatthew G Knepley ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 200ccd2543fSMatthew G Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, cells, PETSC_OWN_POINTER, cellIS);CHKERRQ(ierr); 201ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 202ccd2543fSMatthew G Knepley } 203ccd2543fSMatthew G Knepley 204ccd2543fSMatthew G Knepley #undef __FUNCT__ 20517fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeProjection2Dto1D_Internal" 20617fe8556SMatthew G. Knepley /* 20717fe8556SMatthew G. Knepley DMPlexComputeProjection2Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 2D 20817fe8556SMatthew G. Knepley */ 2097f07f362SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection2Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 21017fe8556SMatthew G. Knepley { 21117fe8556SMatthew G. Knepley const PetscReal x = PetscRealPart(coords[2] - coords[0]); 21217fe8556SMatthew G. Knepley const PetscReal y = PetscRealPart(coords[3] - coords[1]); 2138b49ba18SBarry Smith const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 21417fe8556SMatthew G. Knepley 21517fe8556SMatthew G. Knepley PetscFunctionBegin; 2161c99cf0cSGeoffrey Irving R[0] = c; R[1] = -s; 2171c99cf0cSGeoffrey Irving R[2] = s; R[3] = c; 21817fe8556SMatthew G. Knepley coords[0] = 0.0; 2197f07f362SMatthew G. Knepley coords[1] = r; 22017fe8556SMatthew G. Knepley PetscFunctionReturn(0); 22117fe8556SMatthew G. Knepley } 22217fe8556SMatthew G. Knepley 22317fe8556SMatthew G. Knepley #undef __FUNCT__ 224ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 225ccd2543fSMatthew G Knepley /* 226ccd2543fSMatthew G Knepley DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 227ccd2543fSMatthew G Knepley */ 22899dec3a6SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 229ccd2543fSMatthew G Knepley { 2301ee9d5ecSMatthew G. Knepley PetscReal x1[3], x2[3], n[3], norm; 23199dec3a6SMatthew G. Knepley PetscReal x1p[3], x2p[3], xnp[3]; 2324a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 233ccd2543fSMatthew G Knepley const PetscInt dim = 3; 23499dec3a6SMatthew G. Knepley PetscInt d, e, p; 235ccd2543fSMatthew G Knepley 236ccd2543fSMatthew G Knepley PetscFunctionBegin; 237ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 238ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 2391ee9d5ecSMatthew G. Knepley x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 2401ee9d5ecSMatthew G. Knepley x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 241ccd2543fSMatthew G Knepley } 242ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 243ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 244ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 2458b49ba18SBarry Smith norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 246ccd2543fSMatthew G Knepley n[0] /= norm; 247ccd2543fSMatthew G Knepley n[1] /= norm; 248ccd2543fSMatthew G Knepley n[2] /= norm; 249ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 250ccd2543fSMatthew G Knepley 251ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 252ccd2543fSMatthew G Knepley 253ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 254ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 255ccd2543fSMatthew G Knepley \ nx ny nz / 256ccd2543fSMatthew G Knepley 257ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 258ccd2543fSMatthew G Knepley */ 2598b49ba18SBarry Smith sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 26073868372SMatthew G. Knepley /* Check for n = z */ 26173868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 2621ee9d5ecSMatthew G. Knepley if (n[2] < 0.0) { 26399dec3a6SMatthew G. Knepley if (coordSize > 9) { 26499dec3a6SMatthew G. Knepley coords[2] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]); 26508b58242SMatthew G. Knepley coords[3] = PetscRealPart(coords[3*dim+1] - coords[0*dim+1]); 26699dec3a6SMatthew G. Knepley coords[4] = x2[0]; 26799dec3a6SMatthew G. Knepley coords[5] = x2[1]; 26899dec3a6SMatthew G. Knepley coords[6] = x1[0]; 26999dec3a6SMatthew G. Knepley coords[7] = x1[1]; 27099dec3a6SMatthew G. Knepley } else { 27173868372SMatthew G. Knepley coords[2] = x2[0]; 27273868372SMatthew G. Knepley coords[3] = x2[1]; 27373868372SMatthew G. Knepley coords[4] = x1[0]; 27473868372SMatthew G. Knepley coords[5] = x1[1]; 27599dec3a6SMatthew G. Knepley } 276b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 277b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 278b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = -1.0; 27973868372SMatthew G. Knepley } else { 28099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 28199dec3a6SMatthew G. Knepley coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 28299dec3a6SMatthew G. Knepley coords[p*2+1] = PetscRealPart(coords[p*dim+1] - coords[0*dim+1]); 28399dec3a6SMatthew G. Knepley } 28473868372SMatthew G. Knepley coords[2] = x1[0]; 28573868372SMatthew G. Knepley coords[3] = x1[1]; 28673868372SMatthew G. Knepley coords[4] = x2[0]; 28773868372SMatthew G. Knepley coords[5] = x2[1]; 288b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 289b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 290b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = 1.0; 29173868372SMatthew G. Knepley } 29299dec3a6SMatthew G. Knepley coords[0] = 0.0; 29399dec3a6SMatthew G. Knepley coords[1] = 0.0; 29473868372SMatthew G. Knepley PetscFunctionReturn(0); 29573868372SMatthew G. Knepley } 296da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 297ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 298ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 299ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 300ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 301ccd2543fSMatthew G Knepley x1p[d] = 0.0; 302ccd2543fSMatthew G Knepley x2p[d] = 0.0; 303ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 304ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 305ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 306ccd2543fSMatthew G Knepley } 307ccd2543fSMatthew G Knepley } 3088763be8eSMatthew G. Knepley if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 3098763be8eSMatthew G. Knepley if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 310ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 31199dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 31299dec3a6SMatthew G. Knepley for (d = 0; d < dim; ++d) { 31399dec3a6SMatthew G. Knepley xnp[d] = 0.0; 31499dec3a6SMatthew G. Knepley for (e = 0; e < dim; ++e) { 31599dec3a6SMatthew G. Knepley xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 31699dec3a6SMatthew G. Knepley } 31799dec3a6SMatthew G. Knepley if (d < dim-1) coords[p*2+d] = xnp[d]; 31899dec3a6SMatthew G. Knepley } 31999dec3a6SMatthew G. Knepley } 320ccd2543fSMatthew G Knepley coords[0] = 0.0; 321ccd2543fSMatthew G Knepley coords[1] = 0.0; 322ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 323ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 324ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 325ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 3267f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 3277f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 3287f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 3297f07f362SMatthew G. Knepley PetscReal tmp; 3307f07f362SMatthew G. Knepley 3317f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 3327f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 3337f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 3347f07f362SMatthew G. Knepley } 3357f07f362SMatthew G. Knepley } 336ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 337ccd2543fSMatthew G Knepley } 338ccd2543fSMatthew G Knepley 339ccd2543fSMatthew G Knepley #undef __FUNCT__ 340834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal" 3416322fe33SJed Brown PETSC_UNUSED 342834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 343834e62ceSMatthew G. Knepley { 344834e62ceSMatthew G. Knepley /* Signed volume is 1/2 the determinant 345834e62ceSMatthew G. Knepley 346834e62ceSMatthew G. Knepley | 1 1 1 | 347834e62ceSMatthew G. Knepley | x0 x1 x2 | 348834e62ceSMatthew G. Knepley | y0 y1 y2 | 349834e62ceSMatthew G. Knepley 350834e62ceSMatthew G. Knepley but if x0,y0 is the origin, we have 351834e62ceSMatthew G. Knepley 352834e62ceSMatthew G. Knepley | x1 x2 | 353834e62ceSMatthew G. Knepley | y1 y2 | 354834e62ceSMatthew G. Knepley */ 355834e62ceSMatthew G. Knepley const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 356834e62ceSMatthew G. Knepley const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 357834e62ceSMatthew G. Knepley PetscReal M[4], detM; 358834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; 35986623015SMatthew G. Knepley M[2] = y1; M[3] = y2; 360923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(&detM, M); 361834e62ceSMatthew G. Knepley *vol = 0.5*detM; 362834e62ceSMatthew G. Knepley PetscLogFlops(5.0); 363834e62ceSMatthew G. Knepley } 364834e62ceSMatthew G. Knepley 365834e62ceSMatthew G. Knepley #undef __FUNCT__ 366834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal" 367834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 368834e62ceSMatthew G. Knepley { 369923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(vol, coords); 370834e62ceSMatthew G. Knepley *vol *= 0.5; 371834e62ceSMatthew G. Knepley } 372834e62ceSMatthew G. Knepley 373834e62ceSMatthew G. Knepley #undef __FUNCT__ 374834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal" 3756322fe33SJed Brown PETSC_UNUSED 376834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 377834e62ceSMatthew G. Knepley { 378834e62ceSMatthew G. Knepley /* Signed volume is 1/6th of the determinant 379834e62ceSMatthew G. Knepley 380834e62ceSMatthew G. Knepley | 1 1 1 1 | 381834e62ceSMatthew G. Knepley | x0 x1 x2 x3 | 382834e62ceSMatthew G. Knepley | y0 y1 y2 y3 | 383834e62ceSMatthew G. Knepley | z0 z1 z2 z3 | 384834e62ceSMatthew G. Knepley 385834e62ceSMatthew G. Knepley but if x0,y0,z0 is the origin, we have 386834e62ceSMatthew G. Knepley 387834e62ceSMatthew G. Knepley | x1 x2 x3 | 388834e62ceSMatthew G. Knepley | y1 y2 y3 | 389834e62ceSMatthew G. Knepley | z1 z2 z3 | 390834e62ceSMatthew G. Knepley */ 391834e62ceSMatthew G. Knepley const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 392834e62ceSMatthew G. Knepley const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 393834e62ceSMatthew G. Knepley const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 394834e62ceSMatthew G. Knepley PetscReal M[9], detM; 395834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; M[2] = x3; 396834e62ceSMatthew G. Knepley M[3] = y1; M[4] = y2; M[5] = y3; 397834e62ceSMatthew G. Knepley M[6] = z1; M[7] = z2; M[8] = z3; 398923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(&detM, M); 399b7ad821dSMatthew G. Knepley *vol = -0.16666666666666666666666*detM; 400834e62ceSMatthew G. Knepley PetscLogFlops(10.0); 401834e62ceSMatthew G. Knepley } 402834e62ceSMatthew G. Knepley 403834e62ceSMatthew G. Knepley #undef __FUNCT__ 4040ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 4050ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 4060ec8681fSMatthew G. Knepley { 407923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(vol, coords); 408b7ad821dSMatthew G. Knepley *vol *= -0.16666666666666666666666; 4090ec8681fSMatthew G. Knepley } 4100ec8681fSMatthew G. Knepley 4110ec8681fSMatthew G. Knepley #undef __FUNCT__ 41217fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 41317fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 41417fe8556SMatthew G. Knepley { 41517fe8556SMatthew G. Knepley PetscSection coordSection; 41617fe8556SMatthew G. Knepley Vec coordinates; 417a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 4187f07f362SMatthew G. Knepley PetscInt numCoords, d; 41917fe8556SMatthew G. Knepley PetscErrorCode ierr; 42017fe8556SMatthew G. Knepley 42117fe8556SMatthew G. Knepley PetscFunctionBegin; 42217fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 42369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 42417fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 4257f07f362SMatthew G. Knepley *detJ = 0.0; 42617fe8556SMatthew G. Knepley if (numCoords == 4) { 4277f07f362SMatthew G. Knepley const PetscInt dim = 2; 4287f07f362SMatthew G. Knepley PetscReal R[4], J0; 4297f07f362SMatthew G. Knepley 4307f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 4317f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 43217fe8556SMatthew G. Knepley if (J) { 4337f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 4347f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 4357f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 436923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 43717fe8556SMatthew G. Knepley } 438923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 4397f07f362SMatthew G. Knepley } else if (numCoords == 2) { 4407f07f362SMatthew G. Knepley const PetscInt dim = 1; 4417f07f362SMatthew G. Knepley 4427f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 4437f07f362SMatthew G. Knepley if (J) { 4447f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 44517fe8556SMatthew G. Knepley *detJ = J[0]; 44617fe8556SMatthew G. Knepley PetscLogFlops(2.0); 44717fe8556SMatthew G. Knepley } 4487f07f362SMatthew G. Knepley if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 449796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 45017fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 45117fe8556SMatthew G. Knepley PetscFunctionReturn(0); 45217fe8556SMatthew G. Knepley } 45317fe8556SMatthew G. Knepley 45417fe8556SMatthew G. Knepley #undef __FUNCT__ 455ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 456ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 457ccd2543fSMatthew G Knepley { 458ccd2543fSMatthew G Knepley PetscSection coordSection; 459ccd2543fSMatthew G Knepley Vec coordinates; 460a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 4617f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 462ccd2543fSMatthew G Knepley PetscErrorCode ierr; 463ccd2543fSMatthew G Knepley 464ccd2543fSMatthew G Knepley PetscFunctionBegin; 465ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 46669d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 467ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 4687f07f362SMatthew G. Knepley *detJ = 0.0; 469ccd2543fSMatthew G Knepley if (numCoords == 9) { 4707f07f362SMatthew G. Knepley const PetscInt dim = 3; 4717f07f362SMatthew 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}; 4727f07f362SMatthew G. Knepley 4737f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 47499dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 4757f07f362SMatthew G. Knepley if (J) { 476b7ad821dSMatthew G. Knepley const PetscInt pdim = 2; 477b7ad821dSMatthew G. Knepley 478b7ad821dSMatthew G. Knepley for (d = 0; d < pdim; d++) { 479b7ad821dSMatthew G. Knepley for (f = 0; f < pdim; f++) { 480b7ad821dSMatthew G. Knepley J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 481ccd2543fSMatthew G Knepley } 4827f07f362SMatthew G. Knepley } 4837f07f362SMatthew G. Knepley PetscLogFlops(8.0); 484923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 4857f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 4867f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 4877f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 4887f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 4897f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 4907f07f362SMatthew G. Knepley } 4917f07f362SMatthew G. Knepley } 4927f07f362SMatthew G. Knepley } 4937f07f362SMatthew G. Knepley PetscLogFlops(18.0); 4947f07f362SMatthew G. Knepley } 495923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 4967f07f362SMatthew G. Knepley } else if (numCoords == 6) { 4977f07f362SMatthew G. Knepley const PetscInt dim = 2; 4987f07f362SMatthew G. Knepley 4997f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 500ccd2543fSMatthew G Knepley if (J) { 501ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 502ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 503ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 504ccd2543fSMatthew G Knepley } 505ccd2543fSMatthew G Knepley } 5067f07f362SMatthew G. Knepley PetscLogFlops(8.0); 507923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 508ccd2543fSMatthew G Knepley } 509923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 510796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 511ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 512ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 513ccd2543fSMatthew G Knepley } 514ccd2543fSMatthew G Knepley 515ccd2543fSMatthew G Knepley #undef __FUNCT__ 516ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 517ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 518ccd2543fSMatthew G Knepley { 519ccd2543fSMatthew G Knepley PetscSection coordSection; 520ccd2543fSMatthew G Knepley Vec coordinates; 521a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 52299dec3a6SMatthew G. Knepley PetscInt numCoords, d, f, g; 523ccd2543fSMatthew G Knepley PetscErrorCode ierr; 524ccd2543fSMatthew G Knepley 525ccd2543fSMatthew G Knepley PetscFunctionBegin; 526ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 52769d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 52899dec3a6SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 5297f07f362SMatthew G. Knepley *detJ = 0.0; 53099dec3a6SMatthew G. Knepley if (numCoords == 12) { 53199dec3a6SMatthew G. Knepley const PetscInt dim = 3; 53299dec3a6SMatthew 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}; 53399dec3a6SMatthew G. Knepley 53499dec3a6SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 53599dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 53699dec3a6SMatthew G. Knepley if (J) { 53799dec3a6SMatthew G. Knepley const PetscInt pdim = 2; 53899dec3a6SMatthew G. Knepley 53999dec3a6SMatthew G. Knepley for (d = 0; d < pdim; d++) { 54099dec3a6SMatthew G. Knepley J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 54199dec3a6SMatthew G. Knepley J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 54299dec3a6SMatthew G. Knepley } 54399dec3a6SMatthew G. Knepley PetscLogFlops(8.0); 544923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 54599dec3a6SMatthew G. Knepley for (d = 0; d < dim; d++) { 54699dec3a6SMatthew G. Knepley for (f = 0; f < dim; f++) { 54799dec3a6SMatthew G. Knepley J[d*dim+f] = 0.0; 54899dec3a6SMatthew G. Knepley for (g = 0; g < dim; g++) { 54999dec3a6SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 55099dec3a6SMatthew G. Knepley } 55199dec3a6SMatthew G. Knepley } 55299dec3a6SMatthew G. Knepley } 55399dec3a6SMatthew G. Knepley PetscLogFlops(18.0); 55499dec3a6SMatthew G. Knepley } 555923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 55697f1a218SMatthew G. Knepley } else if ((numCoords == 8) || (numCoords == 16)) { 55799dec3a6SMatthew G. Knepley const PetscInt dim = 2; 55899dec3a6SMatthew G. Knepley 5597f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 560ccd2543fSMatthew G Knepley if (J) { 561ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 56299dec3a6SMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 56399dec3a6SMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 564ccd2543fSMatthew G Knepley } 5657f07f362SMatthew G. Knepley PetscLogFlops(8.0); 566923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 567ccd2543fSMatthew G Knepley } 568923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 569796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 57099dec3a6SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 571ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 572ccd2543fSMatthew G Knepley } 573ccd2543fSMatthew G Knepley 574ccd2543fSMatthew G Knepley #undef __FUNCT__ 575ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 576ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 577ccd2543fSMatthew G Knepley { 578ccd2543fSMatthew G Knepley PetscSection coordSection; 579ccd2543fSMatthew G Knepley Vec coordinates; 580a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 581ccd2543fSMatthew G Knepley const PetscInt dim = 3; 58299dec3a6SMatthew G. Knepley PetscInt d; 583ccd2543fSMatthew G Knepley PetscErrorCode ierr; 584ccd2543fSMatthew G Knepley 585ccd2543fSMatthew G Knepley PetscFunctionBegin; 586ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 58769d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 588ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 5897f07f362SMatthew G. Knepley *detJ = 0.0; 5907f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 591ccd2543fSMatthew G Knepley if (J) { 592ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 593f0df753eSMatthew G. Knepley /* I orient with outward face normals */ 594f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 595f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 596f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 597ccd2543fSMatthew G Knepley } 5987f07f362SMatthew G. Knepley PetscLogFlops(18.0); 599923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 600ccd2543fSMatthew G Knepley } 601923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 602ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 603ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 604ccd2543fSMatthew G Knepley } 605ccd2543fSMatthew G Knepley 606ccd2543fSMatthew G Knepley #undef __FUNCT__ 607ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 608ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 609ccd2543fSMatthew G Knepley { 610ccd2543fSMatthew G Knepley PetscSection coordSection; 611ccd2543fSMatthew G Knepley Vec coordinates; 612a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 613ccd2543fSMatthew G Knepley const PetscInt dim = 3; 614ccd2543fSMatthew G Knepley PetscInt d; 615ccd2543fSMatthew G Knepley PetscErrorCode ierr; 616ccd2543fSMatthew G Knepley 617ccd2543fSMatthew G Knepley PetscFunctionBegin; 618ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 61969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 620ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 6217f07f362SMatthew G. Knepley *detJ = 0.0; 6227f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 623ccd2543fSMatthew G Knepley if (J) { 624ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 625f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 626f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 627f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 628ccd2543fSMatthew G Knepley } 6297f07f362SMatthew G. Knepley PetscLogFlops(18.0); 630923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 631ccd2543fSMatthew G Knepley } 632923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 633ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 634ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 635ccd2543fSMatthew G Knepley } 636ccd2543fSMatthew G Knepley 637ccd2543fSMatthew G Knepley #undef __FUNCT__ 6388e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 639ccd2543fSMatthew G Knepley /*@C 6408e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 641ccd2543fSMatthew G Knepley 642ccd2543fSMatthew G Knepley Collective on DM 643ccd2543fSMatthew G Knepley 644ccd2543fSMatthew G Knepley Input Arguments: 645ccd2543fSMatthew G Knepley + dm - the DM 646ccd2543fSMatthew G Knepley - cell - the cell 647ccd2543fSMatthew G Knepley 648ccd2543fSMatthew G Knepley Output Arguments: 649ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 650ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 651ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 652ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 653ccd2543fSMatthew G Knepley 654ccd2543fSMatthew G Knepley Level: advanced 655ccd2543fSMatthew G Knepley 656ccd2543fSMatthew G Knepley Fortran Notes: 657ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 658ccd2543fSMatthew G Knepley include petsc.h90 in your code. 659ccd2543fSMatthew G Knepley 6608e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 661ccd2543fSMatthew G Knepley @*/ 6628e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 663ccd2543fSMatthew G Knepley { 66449dc4407SMatthew G. Knepley PetscInt depth, dim, coneSize; 665ccd2543fSMatthew G Knepley PetscErrorCode ierr; 666ccd2543fSMatthew G Knepley 667ccd2543fSMatthew G Knepley PetscFunctionBegin; 668139a35ccSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 669ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 67049dc4407SMatthew G. Knepley if (depth == 1) { 6718e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 6728e0841e0SMatthew G. Knepley } else { 6738e0841e0SMatthew G. Knepley DMLabel depth; 6748e0841e0SMatthew G. Knepley 6758e0841e0SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 6768e0841e0SMatthew G. Knepley ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 6778e0841e0SMatthew G. Knepley } 678ccd2543fSMatthew G Knepley switch (dim) { 67917fe8556SMatthew G. Knepley case 1: 68017fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 68117fe8556SMatthew G. Knepley break; 682ccd2543fSMatthew G Knepley case 2: 683ccd2543fSMatthew G Knepley switch (coneSize) { 684ccd2543fSMatthew G Knepley case 3: 685ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 686ccd2543fSMatthew G Knepley break; 687ccd2543fSMatthew G Knepley case 4: 688ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 689ccd2543fSMatthew G Knepley break; 690ccd2543fSMatthew G Knepley default: 6918e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 692ccd2543fSMatthew G Knepley } 693ccd2543fSMatthew G Knepley break; 694ccd2543fSMatthew G Knepley case 3: 695ccd2543fSMatthew G Knepley switch (coneSize) { 696ccd2543fSMatthew G Knepley case 4: 697ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 698ccd2543fSMatthew G Knepley break; 6998e0841e0SMatthew G. Knepley case 6: /* Faces */ 7008e0841e0SMatthew G. Knepley case 8: /* Vertices */ 701ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 702ccd2543fSMatthew G Knepley break; 703ccd2543fSMatthew G Knepley default: 7048e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 705ccd2543fSMatthew G Knepley } 706ccd2543fSMatthew G Knepley break; 707ccd2543fSMatthew G Knepley default: 708ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 709ccd2543fSMatthew G Knepley } 7108e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 7118e0841e0SMatthew G. Knepley } 7128e0841e0SMatthew G. Knepley 7138e0841e0SMatthew G. Knepley #undef __FUNCT__ 7148e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 7158e0841e0SMatthew G. Knepley static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 7168e0841e0SMatthew G. Knepley { 7178e0841e0SMatthew G. Knepley PetscQuadrature quad; 7188e0841e0SMatthew G. Knepley PetscSection coordSection; 7198e0841e0SMatthew G. Knepley Vec coordinates; 7208e0841e0SMatthew G. Knepley PetscScalar *coords = NULL; 7218e0841e0SMatthew G. Knepley const PetscReal *quadPoints; 7228e0841e0SMatthew G. Knepley PetscReal *basisDer; 7238e0841e0SMatthew G. Knepley PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 7248e0841e0SMatthew G. Knepley PetscErrorCode ierr; 7258e0841e0SMatthew G. Knepley 7268e0841e0SMatthew G. Knepley PetscFunctionBegin; 7278e0841e0SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 7288e0841e0SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 7298e0841e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 7308e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 7318e0841e0SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 7328e0841e0SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 733*954b1791SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 7348e0841e0SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 7358e0841e0SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 7368e0841e0SMatthew G. Knepley *detJ = 0.0; 7378e0841e0SMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 7388e0841e0SMatthew 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); 7398e0841e0SMatthew G. Knepley if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 7408e0841e0SMatthew G. Knepley if (J) { 7418e0841e0SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 7428e0841e0SMatthew G. Knepley PetscInt i, j, k, c, r; 7438e0841e0SMatthew G. Knepley 7448e0841e0SMatthew G. Knepley /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 7458e0841e0SMatthew G. Knepley for (k = 0; k < pdim; ++k) 7468e0841e0SMatthew G. Knepley for (j = 0; j < dim; ++j) 7478e0841e0SMatthew G. Knepley for (i = 0; i < cdim; ++i) 74871d6e60fSMatthew G. Knepley J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 7498e0841e0SMatthew G. Knepley PetscLogFlops(2.0*pdim*dim*cdim); 7508e0841e0SMatthew G. Knepley if (cdim > dim) { 7518e0841e0SMatthew G. Knepley for (c = dim; c < cdim; ++c) 7528e0841e0SMatthew G. Knepley for (r = 0; r < cdim; ++r) 7538e0841e0SMatthew G. Knepley J[r*cdim+c] = r == c ? 1.0 : 0.0; 7548e0841e0SMatthew G. Knepley } 7558e0841e0SMatthew G. Knepley switch (cdim) { 7568e0841e0SMatthew G. Knepley case 3: 7578e0841e0SMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 7588e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 75917fe8556SMatthew G. Knepley break; 76049dc4407SMatthew G. Knepley case 2: 7618e0841e0SMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 7628e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 76349dc4407SMatthew G. Knepley break; 7648e0841e0SMatthew G. Knepley case 1: 7658e0841e0SMatthew G. Knepley *detJ = J[0]; 7668e0841e0SMatthew G. Knepley if (invJ) invJ[0] = 1.0/J[0]; 76749dc4407SMatthew G. Knepley } 76849dc4407SMatthew G. Knepley } 7698e0841e0SMatthew G. Knepley } 7708e0841e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 7718e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 7728e0841e0SMatthew G. Knepley } 7738e0841e0SMatthew G. Knepley 7748e0841e0SMatthew G. Knepley #undef __FUNCT__ 7758e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 7768e0841e0SMatthew G. Knepley /*@C 7778e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 7788e0841e0SMatthew G. Knepley 7798e0841e0SMatthew G. Knepley Collective on DM 7808e0841e0SMatthew G. Knepley 7818e0841e0SMatthew G. Knepley Input Arguments: 7828e0841e0SMatthew G. Knepley + dm - the DM 7838e0841e0SMatthew G. Knepley . cell - the cell 7848e0841e0SMatthew G. Knepley - fe - the finite element containing the quadrature 7858e0841e0SMatthew G. Knepley 7868e0841e0SMatthew G. Knepley Output Arguments: 7878e0841e0SMatthew G. Knepley + v0 - the translation part of this transform 7888e0841e0SMatthew G. Knepley . J - the Jacobian of the transform from the reference element at each quadrature point 7898e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point 7908e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point 7918e0841e0SMatthew G. Knepley 7928e0841e0SMatthew G. Knepley Level: advanced 7938e0841e0SMatthew G. Knepley 7948e0841e0SMatthew G. Knepley Fortran Notes: 7958e0841e0SMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 7968e0841e0SMatthew G. Knepley include petsc.h90 in your code. 7978e0841e0SMatthew G. Knepley 7988e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 7998e0841e0SMatthew G. Knepley @*/ 8008e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 8018e0841e0SMatthew G. Knepley { 8028e0841e0SMatthew G. Knepley PetscErrorCode ierr; 8038e0841e0SMatthew G. Knepley 8048e0841e0SMatthew G. Knepley PetscFunctionBegin; 8058e0841e0SMatthew G. Knepley if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 8068e0841e0SMatthew G. Knepley else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 807ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 808ccd2543fSMatthew G Knepley } 809834e62ceSMatthew G. Knepley 810834e62ceSMatthew G. Knepley #undef __FUNCT__ 811cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 812011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 813cc08537eSMatthew G. Knepley { 814cc08537eSMatthew G. Knepley PetscSection coordSection; 815cc08537eSMatthew G. Knepley Vec coordinates; 816a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 817cc08537eSMatthew G. Knepley PetscInt coordSize; 818cc08537eSMatthew G. Knepley PetscErrorCode ierr; 819cc08537eSMatthew G. Knepley 820cc08537eSMatthew G. Knepley PetscFunctionBegin; 821cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 82269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 823cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 824011ea5d8SMatthew G. Knepley if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 825cc08537eSMatthew G. Knepley if (centroid) { 8261ee9d5ecSMatthew G. Knepley centroid[0] = 0.5*PetscRealPart(coords[0] + coords[dim+0]); 8271ee9d5ecSMatthew G. Knepley centroid[1] = 0.5*PetscRealPart(coords[1] + coords[dim+1]); 828cc08537eSMatthew G. Knepley } 829cc08537eSMatthew G. Knepley if (normal) { 830a60a936bSMatthew G. Knepley PetscReal norm; 831a60a936bSMatthew G. Knepley 8320194f449SMatthew G. Knepley normal[0] = -PetscRealPart(coords[1] - coords[dim+1]); 8330194f449SMatthew G. Knepley normal[1] = PetscRealPart(coords[0] - coords[dim+0]); 834a60a936bSMatthew G. Knepley norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 835a60a936bSMatthew G. Knepley normal[0] /= norm; 836a60a936bSMatthew G. Knepley normal[1] /= norm; 837cc08537eSMatthew G. Knepley } 838cc08537eSMatthew G. Knepley if (vol) { 8398b49ba18SBarry Smith *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - coords[dim+0])) + PetscSqr(PetscRealPart(coords[1] - coords[dim+1]))); 840cc08537eSMatthew G. Knepley } 841cc08537eSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 842cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 843cc08537eSMatthew G. Knepley } 844cc08537eSMatthew G. Knepley 845cc08537eSMatthew G. Knepley #undef __FUNCT__ 846cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 847cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 848011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 849cc08537eSMatthew G. Knepley { 850cc08537eSMatthew G. Knepley PetscSection coordSection; 851cc08537eSMatthew G. Knepley Vec coordinates; 852cc08537eSMatthew G. Knepley PetscScalar *coords = NULL; 8530a1d6728SMatthew G. Knepley PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 8540a1d6728SMatthew G. Knepley PetscInt tdim = 2, coordSize, numCorners, p, d, e; 855cc08537eSMatthew G. Knepley PetscErrorCode ierr; 856cc08537eSMatthew G. Knepley 857cc08537eSMatthew G. Knepley PetscFunctionBegin; 858cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 8590a1d6728SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 86069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 861cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 8620a1d6728SMatthew G. Knepley dim = coordSize/numCorners; 863011ea5d8SMatthew G. Knepley if (normal) { 864011ea5d8SMatthew G. Knepley if (dim > 2) { 8651ee9d5ecSMatthew G. Knepley const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 8661ee9d5ecSMatthew G. Knepley const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 8671ee9d5ecSMatthew G. Knepley const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 8680a1d6728SMatthew G. Knepley PetscReal norm; 8690a1d6728SMatthew G. Knepley 8701ee9d5ecSMatthew G. Knepley v0[0] = PetscRealPart(coords[0]); 8711ee9d5ecSMatthew G. Knepley v0[1] = PetscRealPart(coords[1]); 8721ee9d5ecSMatthew G. Knepley v0[2] = PetscRealPart(coords[2]); 8730a1d6728SMatthew G. Knepley normal[0] = y0*z1 - z0*y1; 8740a1d6728SMatthew G. Knepley normal[1] = z0*x1 - x0*z1; 8750a1d6728SMatthew G. Knepley normal[2] = x0*y1 - y0*x1; 8768b49ba18SBarry Smith norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 8770a1d6728SMatthew G. Knepley normal[0] /= norm; 8780a1d6728SMatthew G. Knepley normal[1] /= norm; 8790a1d6728SMatthew G. Knepley normal[2] /= norm; 880011ea5d8SMatthew G. Knepley } else { 881011ea5d8SMatthew G. Knepley for (d = 0; d < dim; ++d) normal[d] = 0.0; 882011ea5d8SMatthew G. Knepley } 883011ea5d8SMatthew G. Knepley } 88499dec3a6SMatthew G. Knepley if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 8850a1d6728SMatthew G. Knepley for (p = 0; p < numCorners; ++p) { 8860a1d6728SMatthew G. Knepley /* Need to do this copy to get types right */ 8870a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 8881ee9d5ecSMatthew G. Knepley ctmp[d] = PetscRealPart(coords[p*tdim+d]); 8891ee9d5ecSMatthew G. Knepley ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 8900a1d6728SMatthew G. Knepley } 8910a1d6728SMatthew G. Knepley Volume_Triangle_Origin_Internal(&vtmp, ctmp); 8920a1d6728SMatthew G. Knepley vsum += vtmp; 8930a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 8940a1d6728SMatthew G. Knepley csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 8950a1d6728SMatthew G. Knepley } 8960a1d6728SMatthew G. Knepley } 8970a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 8980a1d6728SMatthew G. Knepley csum[d] /= (tdim+1)*vsum; 8990a1d6728SMatthew G. Knepley } 9000a1d6728SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 901ee6bbdb2SSatish Balay if (vol) *vol = PetscAbsReal(vsum); 9020a1d6728SMatthew G. Knepley if (centroid) { 9030a1d6728SMatthew G. Knepley if (dim > 2) { 9040a1d6728SMatthew G. Knepley for (d = 0; d < dim; ++d) { 9050a1d6728SMatthew G. Knepley centroid[d] = v0[d]; 9060a1d6728SMatthew G. Knepley for (e = 0; e < dim; ++e) { 9070a1d6728SMatthew G. Knepley centroid[d] += R[d*dim+e]*csum[e]; 9080a1d6728SMatthew G. Knepley } 9090a1d6728SMatthew G. Knepley } 9100a1d6728SMatthew G. Knepley } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 9110a1d6728SMatthew G. Knepley } 912cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 913cc08537eSMatthew G. Knepley } 914cc08537eSMatthew G. Knepley 915cc08537eSMatthew G. Knepley #undef __FUNCT__ 9160ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 9170ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 918011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 9190ec8681fSMatthew G. Knepley { 9200ec8681fSMatthew G. Knepley PetscSection coordSection; 9210ec8681fSMatthew G. Knepley Vec coordinates; 9220ec8681fSMatthew G. Knepley PetscScalar *coords = NULL; 92386623015SMatthew G. Knepley PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 924a7df9edeSMatthew G. Knepley const PetscInt *faces, *facesO; 9250ec8681fSMatthew G. Knepley PetscInt numFaces, f, coordSize, numCorners, p, d; 9260ec8681fSMatthew G. Knepley PetscErrorCode ierr; 9270ec8681fSMatthew G. Knepley 9280ec8681fSMatthew G. Knepley PetscFunctionBegin; 9290ec8681fSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 93069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 9310ec8681fSMatthew G. Knepley 932d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 9330ec8681fSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 9340ec8681fSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 935a7df9edeSMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 9360ec8681fSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 937011ea5d8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 9380ec8681fSMatthew G. Knepley numCorners = coordSize/dim; 9390ec8681fSMatthew G. Knepley switch (numCorners) { 9400ec8681fSMatthew G. Knepley case 3: 9410ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9421ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 9431ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 9441ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 9450ec8681fSMatthew G. Knepley } 9460ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 947a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 9480ec8681fSMatthew G. Knepley vsum += vtmp; 9494f25033aSJed Brown if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 9500ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9511ee9d5ecSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 9520ec8681fSMatthew G. Knepley } 9530ec8681fSMatthew G. Knepley } 9540ec8681fSMatthew G. Knepley break; 9550ec8681fSMatthew G. Knepley case 4: 9560ec8681fSMatthew G. Knepley /* DO FOR PYRAMID */ 9570ec8681fSMatthew G. Knepley /* First tet */ 9580ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9591ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 9601ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 9611ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 9620ec8681fSMatthew G. Knepley } 9630ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 964a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 9650ec8681fSMatthew G. Knepley vsum += vtmp; 9660ec8681fSMatthew G. Knepley if (centroid) { 9670ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9680ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 9690ec8681fSMatthew G. Knepley } 9700ec8681fSMatthew G. Knepley } 9710ec8681fSMatthew G. Knepley /* Second tet */ 9720ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9731ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 9741ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 9751ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 9760ec8681fSMatthew G. Knepley } 9770ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 978a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 9790ec8681fSMatthew G. Knepley vsum += vtmp; 9800ec8681fSMatthew G. Knepley if (centroid) { 9810ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9820ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 9830ec8681fSMatthew G. Knepley } 9840ec8681fSMatthew G. Knepley } 9850ec8681fSMatthew G. Knepley break; 9860ec8681fSMatthew G. Knepley default: 987796f034aSJed Brown SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 9880ec8681fSMatthew G. Knepley } 9894f25033aSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 9900ec8681fSMatthew G. Knepley } 9918763be8eSMatthew G. Knepley if (vol) *vol = PetscAbsReal(vsum); 9920ec8681fSMatthew G. Knepley if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 993d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 9940ec8681fSMatthew G. Knepley PetscFunctionReturn(0); 9950ec8681fSMatthew G. Knepley } 9960ec8681fSMatthew G. Knepley 9970ec8681fSMatthew G. Knepley #undef __FUNCT__ 998834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 999834e62ceSMatthew G. Knepley /*@C 1000834e62ceSMatthew G. Knepley DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1001834e62ceSMatthew G. Knepley 1002834e62ceSMatthew G. Knepley Collective on DM 1003834e62ceSMatthew G. Knepley 1004834e62ceSMatthew G. Knepley Input Arguments: 1005834e62ceSMatthew G. Knepley + dm - the DM 1006834e62ceSMatthew G. Knepley - cell - the cell 1007834e62ceSMatthew G. Knepley 1008834e62ceSMatthew G. Knepley Output Arguments: 1009834e62ceSMatthew G. Knepley + volume - the cell volume 1010cc08537eSMatthew G. Knepley . centroid - the cell centroid 1011cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate 1012834e62ceSMatthew G. Knepley 1013834e62ceSMatthew G. Knepley Level: advanced 1014834e62ceSMatthew G. Knepley 1015834e62ceSMatthew G. Knepley Fortran Notes: 1016834e62ceSMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1017834e62ceSMatthew G. Knepley include petsc.h90 in your code. 1018834e62ceSMatthew G. Knepley 101969d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1020834e62ceSMatthew G. Knepley @*/ 1021cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1022834e62ceSMatthew G. Knepley { 10230ec8681fSMatthew G. Knepley PetscInt depth, dim; 1024834e62ceSMatthew G. Knepley PetscErrorCode ierr; 1025834e62ceSMatthew G. Knepley 1026834e62ceSMatthew G. Knepley PetscFunctionBegin; 1027834e62ceSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1028c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1029834e62ceSMatthew G. Knepley if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1030834e62ceSMatthew G. Knepley /* We need to keep a pointer to the depth label */ 1031011ea5d8SMatthew G. Knepley ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1032834e62ceSMatthew G. Knepley /* Cone size is now the number of faces */ 1033011ea5d8SMatthew G. Knepley switch (depth) { 1034cc08537eSMatthew G. Knepley case 1: 1035011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1036cc08537eSMatthew G. Knepley break; 1037834e62ceSMatthew G. Knepley case 2: 1038011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1039834e62ceSMatthew G. Knepley break; 1040834e62ceSMatthew G. Knepley case 3: 1041011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1042834e62ceSMatthew G. Knepley break; 1043834e62ceSMatthew G. Knepley default: 1044834e62ceSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1045834e62ceSMatthew G. Knepley } 1046834e62ceSMatthew G. Knepley PetscFunctionReturn(0); 1047834e62ceSMatthew G. Knepley } 1048