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; 15ccd2543fSMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, 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; 307c1f9639SMatthew G Knepley PetscScalar *coords; 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); 39ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(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; 71ccd2543fSMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, 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; 88ccd2543fSMatthew G Knepley const PetscInt faces[24] = {0, 1, 2, 3, 5, 4, 7, 6, 1, 0, 4, 5, 89ccd2543fSMatthew G Knepley 3, 2, 6, 7, 1, 5, 6, 2, 0, 3, 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); 96ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(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; 148ccd2543fSMatthew G Knepley ierr = DMPlexGetDimension(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); 155ccd2543fSMatthew 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); 156ccd2543fSMatthew G Knepley numPoints /= bs; 157ccd2543fSMatthew G Knepley ierr = PetscMalloc(numPoints * sizeof(PetscInt), &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: 173ccd2543fSMatthew G Knepley 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; 185ccd2543fSMatthew G Knepley case 8: 186ccd2543fSMatthew G Knepley ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, c, &cell);CHKERRQ(ierr); 187ccd2543fSMatthew G Knepley break; 188ccd2543fSMatthew G Knepley default: 189ccd2543fSMatthew G Knepley 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: 195ccd2543fSMatthew G Knepley 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]); 2137f07f362SMatthew G. Knepley const PetscReal r = sqrt(x*x + y*y), c = x/r, s = y/r; 21417fe8556SMatthew G. Knepley 21517fe8556SMatthew G. Knepley PetscFunctionBegin; 2167f07f362SMatthew G. Knepley R[0] = c; R[1] = s; 2177f07f362SMatthew G. Knepley 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 */ 2287f07f362SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[], PetscReal R[]) 229ccd2543fSMatthew G Knepley { 230ccd2543fSMatthew G Knepley PetscScalar x1[3], x2[3], n[3], norm; 2317f07f362SMatthew G. Knepley PetscScalar x1p[3], x2p[3]; 2324a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 233ccd2543fSMatthew G Knepley const PetscInt dim = 3; 234ccd2543fSMatthew G Knepley PetscInt d, e; 235ccd2543fSMatthew G Knepley 236ccd2543fSMatthew G Knepley PetscFunctionBegin; 237ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 238ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 239ccd2543fSMatthew G Knepley x1[d] = coords[1*dim+d] - coords[0*dim+d]; 240ccd2543fSMatthew G Knepley x2[d] = 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]; 245ccd2543fSMatthew G Knepley norm = sqrt(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 */ 2594a217a95SMatthew G. Knepley sqrtz = sqrt(1.0 - PetscAbsScalar(n[2]*n[2])); 26073868372SMatthew G. Knepley /* Check for n = z */ 26173868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 26273868372SMatthew G. Knepley coords[0] = 0.0; 26373868372SMatthew G. Knepley coords[1] = 0.0; 2644a217a95SMatthew G. Knepley if (PetscRealPart(n[2]) < 0.0) { 26573868372SMatthew G. Knepley coords[2] = x2[0]; 26673868372SMatthew G. Knepley coords[3] = x2[1]; 26773868372SMatthew G. Knepley coords[4] = x1[0]; 26873868372SMatthew G. Knepley coords[5] = x1[1]; 26973868372SMatthew G. Knepley } else { 27073868372SMatthew G. Knepley coords[2] = x1[0]; 27173868372SMatthew G. Knepley coords[3] = x1[1]; 27273868372SMatthew G. Knepley coords[4] = x2[0]; 27373868372SMatthew G. Knepley coords[5] = x2[1]; 27473868372SMatthew G. Knepley } 27573868372SMatthew G. Knepley PetscFunctionReturn(0); 27673868372SMatthew G. Knepley } 277da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 278ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 279ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 280ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 281ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 282ccd2543fSMatthew G Knepley x1p[d] = 0.0; 283ccd2543fSMatthew G Knepley x2p[d] = 0.0; 284ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 285ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 286ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 287ccd2543fSMatthew G Knepley } 288ccd2543fSMatthew G Knepley } 28988790e04SMatthew G Knepley if (PetscAbsScalar(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 29088790e04SMatthew G Knepley if (PetscAbsScalar(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 291ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 292ccd2543fSMatthew G Knepley coords[0] = 0.0; 293ccd2543fSMatthew G Knepley coords[1] = 0.0; 294ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 295ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 296ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 297ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 2987f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 2997f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 3007f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 3017f07f362SMatthew G. Knepley PetscReal tmp; 3027f07f362SMatthew G. Knepley 3037f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 3047f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 3057f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 3067f07f362SMatthew G. Knepley } 3077f07f362SMatthew G. Knepley } 308ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 309ccd2543fSMatthew G Knepley } 310ccd2543fSMatthew G Knepley 311ccd2543fSMatthew G Knepley #undef __FUNCT__ 3127f07f362SMatthew G. Knepley #define __FUNCT__ "Invert2D_Internal" 3137f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 3147f07f362SMatthew G. Knepley { 3157f07f362SMatthew G. Knepley const PetscReal invDet = 1.0/detJ; 3167f07f362SMatthew G. Knepley 3177f07f362SMatthew G. Knepley invJ[0] = invDet*J[3]; 3187f07f362SMatthew G. Knepley invJ[1] = -invDet*J[1]; 3197f07f362SMatthew G. Knepley invJ[2] = -invDet*J[2]; 3207f07f362SMatthew G. Knepley invJ[3] = invDet*J[0]; 3217f07f362SMatthew G. Knepley PetscLogFlops(5.0); 3227f07f362SMatthew G. Knepley } 3237f07f362SMatthew G. Knepley 3247f07f362SMatthew G. Knepley #undef __FUNCT__ 3257f07f362SMatthew G. Knepley #define __FUNCT__ "Invert3D_Internal" 3267f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ) 3277f07f362SMatthew G. Knepley { 3287f07f362SMatthew G. Knepley const PetscReal invDet = 1.0/detJ; 3297f07f362SMatthew G. Knepley 3307f07f362SMatthew G. Knepley invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 3317f07f362SMatthew G. Knepley invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 3327f07f362SMatthew G. Knepley invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 3337f07f362SMatthew G. Knepley invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 3347f07f362SMatthew G. Knepley invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 3357f07f362SMatthew G. Knepley invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 3367f07f362SMatthew G. Knepley invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 3377f07f362SMatthew G. Knepley invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 3387f07f362SMatthew G. Knepley invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 3397f07f362SMatthew G. Knepley PetscLogFlops(37.0); 3407f07f362SMatthew G. Knepley } 3417f07f362SMatthew G. Knepley 3427f07f362SMatthew G. Knepley #undef __FUNCT__ 3437f07f362SMatthew G. Knepley #define __FUNCT__ "Det2D_Internal" 3447f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det2D_Internal(PetscReal *detJ, PetscReal J[]) 3457f07f362SMatthew G. Knepley { 3467f07f362SMatthew G. Knepley *detJ = J[0]*J[3] - J[1]*J[2]; 3477f07f362SMatthew G. Knepley PetscLogFlops(3.0); 3487f07f362SMatthew G. Knepley } 3497f07f362SMatthew G. Knepley 3507f07f362SMatthew G. Knepley #undef __FUNCT__ 3517f07f362SMatthew G. Knepley #define __FUNCT__ "Det3D_Internal" 3527f07f362SMatthew G. Knepley PETSC_STATIC_INLINE void Det3D_Internal(PetscReal *detJ, PetscReal J[]) 3537f07f362SMatthew G. Knepley { 3547f07f362SMatthew G. Knepley *detJ = (J[0*3+0]*(J[1*3+2]*J[2*3+1] - J[1*3+1]*J[2*3+2]) + 3557f07f362SMatthew G. Knepley J[0*3+1]*(J[1*3+0]*J[2*3+2] - J[1*3+2]*J[2*3+0]) + 3567f07f362SMatthew G. Knepley J[0*3+2]*(J[1*3+1]*J[2*3+0] - J[1*3+0]*J[2*3+1])); 3577f07f362SMatthew G. Knepley PetscLogFlops(12.0); 3587f07f362SMatthew G. Knepley } 3597f07f362SMatthew G. Knepley 3607f07f362SMatthew G. Knepley #undef __FUNCT__ 36117fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 36217fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 36317fe8556SMatthew G. Knepley { 36417fe8556SMatthew G. Knepley PetscSection coordSection; 36517fe8556SMatthew G. Knepley Vec coordinates; 36617fe8556SMatthew G. Knepley PetscScalar *coords; 3677f07f362SMatthew G. Knepley PetscInt numCoords, d; 36817fe8556SMatthew G. Knepley PetscErrorCode ierr; 36917fe8556SMatthew G. Knepley 37017fe8556SMatthew G. Knepley PetscFunctionBegin; 37117fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 37217fe8556SMatthew G. Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 37317fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 3747f07f362SMatthew G. Knepley *detJ = 0.0; 37517fe8556SMatthew G. Knepley if (numCoords == 4) { 3767f07f362SMatthew G. Knepley const PetscInt dim = 2; 3777f07f362SMatthew G. Knepley PetscReal R[4], J0; 3787f07f362SMatthew G. Knepley 3797f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 3807f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 38117fe8556SMatthew G. Knepley if (J) { 3827f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 3837f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 3847f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 3857f07f362SMatthew G. Knepley Det2D_Internal(detJ, J); 38617fe8556SMatthew G. Knepley } 3877f07f362SMatthew G. Knepley if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 3887f07f362SMatthew G. Knepley } else if (numCoords == 2) { 3897f07f362SMatthew G. Knepley const PetscInt dim = 1; 3907f07f362SMatthew G. Knepley 3917f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 3927f07f362SMatthew G. Knepley if (J) { 3937f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 39417fe8556SMatthew G. Knepley *detJ = J[0]; 39517fe8556SMatthew G. Knepley PetscLogFlops(2.0); 39617fe8556SMatthew G. Knepley } 3977f07f362SMatthew G. Knepley if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 3987f07f362SMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %d != 2", numCoords); 39917fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 40017fe8556SMatthew G. Knepley PetscFunctionReturn(0); 40117fe8556SMatthew G. Knepley } 40217fe8556SMatthew G. Knepley 40317fe8556SMatthew G. Knepley #undef __FUNCT__ 404ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 405ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 406ccd2543fSMatthew G Knepley { 407ccd2543fSMatthew G Knepley PetscSection coordSection; 408ccd2543fSMatthew G Knepley Vec coordinates; 4097c1f9639SMatthew G Knepley PetscScalar *coords; 4107f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 411ccd2543fSMatthew G Knepley PetscErrorCode ierr; 412ccd2543fSMatthew G Knepley 413ccd2543fSMatthew G Knepley PetscFunctionBegin; 414ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 415ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 416ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 4177f07f362SMatthew G. Knepley *detJ = 0.0; 418ccd2543fSMatthew G Knepley if (numCoords == 9) { 4197f07f362SMatthew G. Knepley const PetscInt dim = 3; 4207f07f362SMatthew 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}; 4217f07f362SMatthew G. Knepley 4227f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 4237f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(coords, R);CHKERRQ(ierr); 4247f07f362SMatthew G. Knepley if (J) { 4257f07f362SMatthew G. Knepley for (d = 0; d < dim-1; d++) { 4267f07f362SMatthew G. Knepley for (f = 0; f < dim-1; f++) { 4277f07f362SMatthew G. Knepley J0[d*(dim+1)+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 428ccd2543fSMatthew G Knepley } 4297f07f362SMatthew G. Knepley } 4307f07f362SMatthew G. Knepley PetscLogFlops(8.0); 4317f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 4327f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 4337f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 4347f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 4357f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 4367f07f362SMatthew G. Knepley } 4377f07f362SMatthew G. Knepley } 4387f07f362SMatthew G. Knepley } 4397f07f362SMatthew G. Knepley PetscLogFlops(18.0); 4407f07f362SMatthew G. Knepley Det3D_Internal(detJ, J); 4417f07f362SMatthew G. Knepley } 4427f07f362SMatthew G. Knepley if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 4437f07f362SMatthew G. Knepley } else if (numCoords == 6) { 4447f07f362SMatthew G. Knepley const PetscInt dim = 2; 4457f07f362SMatthew G. Knepley 4467f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 447ccd2543fSMatthew G Knepley if (J) { 448ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 449ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 450ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 451ccd2543fSMatthew G Knepley } 452ccd2543fSMatthew G Knepley } 4537f07f362SMatthew G. Knepley PetscLogFlops(8.0); 4547f07f362SMatthew G. Knepley Det2D_Internal(detJ, J); 455ccd2543fSMatthew G Knepley } 4567f07f362SMatthew G. Knepley if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 4577f07f362SMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords); 458ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 459ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 460ccd2543fSMatthew G Knepley } 461ccd2543fSMatthew G Knepley 462ccd2543fSMatthew G Knepley #undef __FUNCT__ 463ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 464ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 465ccd2543fSMatthew G Knepley { 466ccd2543fSMatthew G Knepley PetscSection coordSection; 467ccd2543fSMatthew G Knepley Vec coordinates; 4687c1f9639SMatthew G Knepley PetscScalar *coords; 469ccd2543fSMatthew G Knepley const PetscInt dim = 2; 470ccd2543fSMatthew G Knepley PetscInt d, f; 471ccd2543fSMatthew G Knepley PetscErrorCode ierr; 472ccd2543fSMatthew G Knepley 473ccd2543fSMatthew G Knepley PetscFunctionBegin; 474ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 475ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 476ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 4777f07f362SMatthew G. Knepley *detJ = 0.0; 4787f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 479ccd2543fSMatthew G Knepley if (J) { 480ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 481ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 482ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 483ccd2543fSMatthew G Knepley } 484ccd2543fSMatthew G Knepley } 4857f07f362SMatthew G. Knepley PetscLogFlops(8.0); 4867f07f362SMatthew G. Knepley Det2D_Internal(detJ, J); 487ccd2543fSMatthew G Knepley } 4887f07f362SMatthew G. Knepley if (invJ) {Invert2D_Internal(invJ, J, *detJ);} 489ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 490ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 491ccd2543fSMatthew G Knepley } 492ccd2543fSMatthew G Knepley 493ccd2543fSMatthew G Knepley #undef __FUNCT__ 494ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 495ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 496ccd2543fSMatthew G Knepley { 497ccd2543fSMatthew G Knepley PetscSection coordSection; 498ccd2543fSMatthew G Knepley Vec coordinates; 4997c1f9639SMatthew G Knepley PetscScalar *coords; 500ccd2543fSMatthew G Knepley const PetscInt dim = 3; 501ccd2543fSMatthew G Knepley PetscInt d, f; 502ccd2543fSMatthew G Knepley PetscErrorCode ierr; 503ccd2543fSMatthew G Knepley 504ccd2543fSMatthew G Knepley PetscFunctionBegin; 505ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 506ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 507ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 5087f07f362SMatthew G. Knepley *detJ = 0.0; 5097f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 510ccd2543fSMatthew G Knepley if (J) { 511ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 512ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 513ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 514ccd2543fSMatthew G Knepley } 515ccd2543fSMatthew G Knepley } 5167f07f362SMatthew G. Knepley PetscLogFlops(18.0); 5177f07f362SMatthew G. Knepley Det3D_Internal(detJ, J); 518ccd2543fSMatthew G Knepley } 5197f07f362SMatthew G. Knepley if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 520ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 521ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 522ccd2543fSMatthew G Knepley } 523ccd2543fSMatthew G Knepley 524ccd2543fSMatthew G Knepley #undef __FUNCT__ 525ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 526ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 527ccd2543fSMatthew G Knepley { 528ccd2543fSMatthew G Knepley PetscSection coordSection; 529ccd2543fSMatthew G Knepley Vec coordinates; 5307c1f9639SMatthew G Knepley PetscScalar *coords; 531ccd2543fSMatthew G Knepley const PetscInt dim = 3; 532ccd2543fSMatthew G Knepley PetscInt d; 533ccd2543fSMatthew G Knepley PetscErrorCode ierr; 534ccd2543fSMatthew G Knepley 535ccd2543fSMatthew G Knepley PetscFunctionBegin; 536ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 537ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 538ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 5397f07f362SMatthew G. Knepley *detJ = 0.0; 5407f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 541ccd2543fSMatthew G Knepley if (J) { 542ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 543*de36ddddSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[(2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 544ccd2543fSMatthew G Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 545ccd2543fSMatthew G Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 546ccd2543fSMatthew G Knepley } 5477f07f362SMatthew G. Knepley PetscLogFlops(18.0); 5487f07f362SMatthew G. Knepley Det3D_Internal(detJ, J); 549ccd2543fSMatthew G Knepley } 5507f07f362SMatthew G. Knepley if (invJ) {Invert3D_Internal(invJ, J, *detJ);} 5517f07f362SMatthew G. Knepley *detJ *= -8.0; 552ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 553ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 554ccd2543fSMatthew G Knepley } 555ccd2543fSMatthew G Knepley 556ccd2543fSMatthew G Knepley #undef __FUNCT__ 557ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry" 558ccd2543fSMatthew G Knepley /*@C 559ccd2543fSMatthew G Knepley DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 560ccd2543fSMatthew G Knepley 561ccd2543fSMatthew G Knepley Collective on DM 562ccd2543fSMatthew G Knepley 563ccd2543fSMatthew G Knepley Input Arguments: 564ccd2543fSMatthew G Knepley + dm - the DM 565ccd2543fSMatthew G Knepley - cell - the cell 566ccd2543fSMatthew G Knepley 567ccd2543fSMatthew G Knepley Output Arguments: 568ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 569ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 570ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 571ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 572ccd2543fSMatthew G Knepley 573ccd2543fSMatthew G Knepley Level: advanced 574ccd2543fSMatthew G Knepley 575ccd2543fSMatthew G Knepley Fortran Notes: 576ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 577ccd2543fSMatthew G Knepley include petsc.h90 in your code. 578ccd2543fSMatthew G Knepley 579ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 580ccd2543fSMatthew G Knepley @*/ 581ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 582ccd2543fSMatthew G Knepley { 58349dc4407SMatthew G. Knepley PetscInt depth, dim, coneSize; 584ccd2543fSMatthew G Knepley PetscErrorCode ierr; 585ccd2543fSMatthew G Knepley 586ccd2543fSMatthew G Knepley PetscFunctionBegin; 587139a35ccSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 588ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 58949dc4407SMatthew G. Knepley if (depth == 1) { 5905743f1d7SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 591ccd2543fSMatthew G Knepley switch (dim) { 59217fe8556SMatthew G. Knepley case 1: 59317fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 59417fe8556SMatthew G. Knepley break; 595ccd2543fSMatthew G Knepley case 2: 596ccd2543fSMatthew G Knepley switch (coneSize) { 597ccd2543fSMatthew G Knepley case 3: 598ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 599ccd2543fSMatthew G Knepley break; 600ccd2543fSMatthew G Knepley case 4: 601ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 602ccd2543fSMatthew G Knepley break; 603ccd2543fSMatthew G Knepley default: 604ccd2543fSMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 605ccd2543fSMatthew G Knepley } 606ccd2543fSMatthew G Knepley break; 607ccd2543fSMatthew G Knepley case 3: 608ccd2543fSMatthew G Knepley switch (coneSize) { 609ccd2543fSMatthew G Knepley case 4: 610ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 611ccd2543fSMatthew G Knepley break; 612ccd2543fSMatthew G Knepley case 8: 613ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 614ccd2543fSMatthew G Knepley break; 615ccd2543fSMatthew G Knepley default: 616ccd2543fSMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 617ccd2543fSMatthew G Knepley } 618ccd2543fSMatthew G Knepley break; 619ccd2543fSMatthew G Knepley default: 620ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 621ccd2543fSMatthew G Knepley } 62249dc4407SMatthew G. Knepley } else { 6235743f1d7SMatthew G. Knepley /* We need to keep a pointer to the depth label */ 6245743f1d7SMatthew G. Knepley ierr = DMPlexGetLabelValue(dm, "depth", cell, &dim);CHKERRQ(ierr); 62549dc4407SMatthew G. Knepley /* Cone size is now the number of faces */ 62649dc4407SMatthew G. Knepley switch (dim) { 62717fe8556SMatthew G. Knepley case 1: 62817fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 62917fe8556SMatthew G. Knepley break; 63049dc4407SMatthew G. Knepley case 2: 63149dc4407SMatthew G. Knepley switch (coneSize) { 63249dc4407SMatthew G. Knepley case 3: 63349dc4407SMatthew G. Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 63449dc4407SMatthew G. Knepley break; 63549dc4407SMatthew G. Knepley case 4: 63649dc4407SMatthew G. Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 63749dc4407SMatthew G. Knepley break; 63849dc4407SMatthew G. Knepley default: 63949dc4407SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 64049dc4407SMatthew G. Knepley } 64149dc4407SMatthew G. Knepley break; 64249dc4407SMatthew G. Knepley case 3: 64349dc4407SMatthew G. Knepley switch (coneSize) { 64449dc4407SMatthew G. Knepley case 4: 64549dc4407SMatthew G. Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 64649dc4407SMatthew G. Knepley break; 64749dc4407SMatthew G. Knepley case 6: 64849dc4407SMatthew G. Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 64949dc4407SMatthew G. Knepley break; 65049dc4407SMatthew G. Knepley default: 65149dc4407SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 65249dc4407SMatthew G. Knepley } 65349dc4407SMatthew G. Knepley break; 65449dc4407SMatthew G. Knepley default: 65549dc4407SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 65649dc4407SMatthew G. Knepley } 65749dc4407SMatthew G. Knepley } 658ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 659ccd2543fSMatthew G Knepley } 660