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; 30*7c1f9639SMatthew 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; 87*7c1f9639SMatthew 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__ 205ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 206ccd2543fSMatthew G Knepley /* 207ccd2543fSMatthew G Knepley DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 208ccd2543fSMatthew G Knepley */ 209ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscScalar coords[]) 210ccd2543fSMatthew G Knepley { 211ccd2543fSMatthew G Knepley PetscScalar x1[3], x2[3], n[3], norm; 212ccd2543fSMatthew G Knepley const PetscInt dim = 3; 213ccd2543fSMatthew G Knepley PetscInt d, e; 214ccd2543fSMatthew G Knepley 215ccd2543fSMatthew G Knepley PetscFunctionBegin; 216ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 217ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 218ccd2543fSMatthew G Knepley x1[d] = coords[1*dim+d] - coords[0*dim+d]; 219ccd2543fSMatthew G Knepley x2[d] = coords[2*dim+d] - coords[0*dim+d]; 220ccd2543fSMatthew G Knepley } 221ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 222ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 223ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 224ccd2543fSMatthew G Knepley norm = sqrt(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 225ccd2543fSMatthew G Knepley n[0] /= norm; 226ccd2543fSMatthew G Knepley n[1] /= norm; 227ccd2543fSMatthew G Knepley n[2] /= norm; 228ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 229ccd2543fSMatthew G Knepley 230ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 231ccd2543fSMatthew G Knepley 232ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 233ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 234ccd2543fSMatthew G Knepley \ nx ny nz / 235ccd2543fSMatthew G Knepley 236ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 237ccd2543fSMatthew G Knepley */ 238ccd2543fSMatthew G Knepley PetscScalar R[9], x1p[3], x2p[3]; 239ccd2543fSMatthew G Knepley PetscScalar sqrtz = sqrt(1.0 - n[2]*n[2]), alpha = 1.0/sqrtz; 240ccd2543fSMatthew G Knepley 241ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 242ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 243ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 244ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 245ccd2543fSMatthew G Knepley x1p[d] = 0.0; 246ccd2543fSMatthew G Knepley x2p[d] = 0.0; 247ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 248ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 249ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 250ccd2543fSMatthew G Knepley } 251ccd2543fSMatthew G Knepley } 252ccd2543fSMatthew G Knepley if (fabs(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 253ccd2543fSMatthew G Knepley if (fabs(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 254ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 255ccd2543fSMatthew G Knepley coords[0] = 0.0; 256ccd2543fSMatthew G Knepley coords[1] = 0.0; 257ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 258ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 259ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 260ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 261ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 262ccd2543fSMatthew G Knepley } 263ccd2543fSMatthew G Knepley 264ccd2543fSMatthew G Knepley #undef __FUNCT__ 265ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 266ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 267ccd2543fSMatthew G Knepley { 268ccd2543fSMatthew G Knepley PetscSection coordSection; 269ccd2543fSMatthew G Knepley Vec coordinates; 270*7c1f9639SMatthew G Knepley PetscScalar *coords; 271ccd2543fSMatthew G Knepley const PetscInt dim = 2; 272ccd2543fSMatthew G Knepley PetscInt numCoords, d, f; 273ccd2543fSMatthew G Knepley PetscErrorCode ierr; 274ccd2543fSMatthew G Knepley 275ccd2543fSMatthew G Knepley PetscFunctionBegin; 276ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 277ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 278ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 279ccd2543fSMatthew G Knepley if (numCoords == 9) { 280ccd2543fSMatthew G Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(coords);CHKERRQ(ierr); 281ccd2543fSMatthew G Knepley } else if (numCoords != 6) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %d != 6", numCoords); 282ccd2543fSMatthew G Knepley if (v0) { 283ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 284ccd2543fSMatthew G Knepley } 285ccd2543fSMatthew G Knepley if (J) { 286ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 287ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 288ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 289ccd2543fSMatthew G Knepley } 290ccd2543fSMatthew G Knepley } 291ccd2543fSMatthew G Knepley *detJ = J[0]*J[3] - J[1]*J[2]; 292ccd2543fSMatthew G Knepley #if 0 293ccd2543fSMatthew G Knepley if (detJ < 0.0) { 294ccd2543fSMatthew G Knepley const PetscReal xLength = mesh->periodicity[0]; 295ccd2543fSMatthew G Knepley 296ccd2543fSMatthew G Knepley if (xLength != 0.0) { 297ccd2543fSMatthew G Knepley PetscReal v0x = coords[0*dim+0]; 298ccd2543fSMatthew G Knepley 299ccd2543fSMatthew G Knepley if (v0x == 0.0) v0x = v0[0] = xLength; 300ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 301ccd2543fSMatthew G Knepley const PetscReal px = coords[(f+1)*dim+0] == 0.0 ? xLength : coords[(f+1)*dim+0]; 302ccd2543fSMatthew G Knepley 303ccd2543fSMatthew G Knepley J[0*dim+f] = 0.5*(px - v0x); 304ccd2543fSMatthew G Knepley } 305ccd2543fSMatthew G Knepley } 306ccd2543fSMatthew G Knepley detJ = J[0]*J[3] - J[1]*J[2]; 307ccd2543fSMatthew G Knepley } 308ccd2543fSMatthew G Knepley #endif 309ccd2543fSMatthew G Knepley PetscLogFlops(8.0 + 3.0); 310ccd2543fSMatthew G Knepley } 311ccd2543fSMatthew G Knepley if (invJ) { 312ccd2543fSMatthew G Knepley const PetscReal invDet = 1.0/(*detJ); 313ccd2543fSMatthew G Knepley 314ccd2543fSMatthew G Knepley invJ[0] = invDet*J[3]; 315ccd2543fSMatthew G Knepley invJ[1] = -invDet*J[1]; 316ccd2543fSMatthew G Knepley invJ[2] = -invDet*J[2]; 317ccd2543fSMatthew G Knepley invJ[3] = invDet*J[0]; 318ccd2543fSMatthew G Knepley PetscLogFlops(5.0); 319ccd2543fSMatthew G Knepley } 320ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 321ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 322ccd2543fSMatthew G Knepley } 323ccd2543fSMatthew G Knepley 324ccd2543fSMatthew G Knepley #undef __FUNCT__ 325ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 326ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 327ccd2543fSMatthew G Knepley { 328ccd2543fSMatthew G Knepley PetscSection coordSection; 329ccd2543fSMatthew G Knepley Vec coordinates; 330*7c1f9639SMatthew G Knepley PetscScalar *coords; 331ccd2543fSMatthew G Knepley const PetscInt dim = 2; 332ccd2543fSMatthew G Knepley PetscInt d, f; 333ccd2543fSMatthew G Knepley PetscErrorCode ierr; 334ccd2543fSMatthew G Knepley 335ccd2543fSMatthew G Knepley PetscFunctionBegin; 336ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 337ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 338ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 339ccd2543fSMatthew G Knepley if (v0) { 340ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 341ccd2543fSMatthew G Knepley } 342ccd2543fSMatthew G Knepley if (J) { 343ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 344ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 345ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f*2+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 346ccd2543fSMatthew G Knepley } 347ccd2543fSMatthew G Knepley } 348ccd2543fSMatthew G Knepley *detJ = J[0]*J[3] - J[1]*J[2]; 349ccd2543fSMatthew G Knepley PetscLogFlops(8.0 + 3.0); 350ccd2543fSMatthew G Knepley } 351ccd2543fSMatthew G Knepley if (invJ) { 352ccd2543fSMatthew G Knepley const PetscReal invDet = 1.0/(*detJ); 353ccd2543fSMatthew G Knepley 354ccd2543fSMatthew G Knepley invJ[0] = invDet*J[3]; 355ccd2543fSMatthew G Knepley invJ[1] = -invDet*J[1]; 356ccd2543fSMatthew G Knepley invJ[2] = -invDet*J[2]; 357ccd2543fSMatthew G Knepley invJ[3] = invDet*J[0]; 358ccd2543fSMatthew G Knepley PetscLogFlops(5.0); 359ccd2543fSMatthew G Knepley } 360ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 361ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 362ccd2543fSMatthew G Knepley } 363ccd2543fSMatthew G Knepley 364ccd2543fSMatthew G Knepley #undef __FUNCT__ 365ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 366ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 367ccd2543fSMatthew G Knepley { 368ccd2543fSMatthew G Knepley PetscSection coordSection; 369ccd2543fSMatthew G Knepley Vec coordinates; 370*7c1f9639SMatthew G Knepley PetscScalar *coords; 371ccd2543fSMatthew G Knepley const PetscInt dim = 3; 372ccd2543fSMatthew G Knepley PetscInt d, f; 373ccd2543fSMatthew G Knepley PetscErrorCode ierr; 374ccd2543fSMatthew G Knepley 375ccd2543fSMatthew G Knepley PetscFunctionBegin; 376ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 377ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 378ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 379ccd2543fSMatthew G Knepley if (v0) { 380ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 381ccd2543fSMatthew G Knepley } 382ccd2543fSMatthew G Knepley if (J) { 383ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 384ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 385ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 386ccd2543fSMatthew G Knepley } 387ccd2543fSMatthew G Knepley } 388ccd2543fSMatthew G Knepley /* ??? This does not work with CTetGen: The minus sign is here since I orient the first face to get the outward normal */ 389ccd2543fSMatthew G Knepley *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 390ccd2543fSMatthew G Knepley J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 391ccd2543fSMatthew G Knepley J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 392ccd2543fSMatthew G Knepley PetscLogFlops(18.0 + 12.0); 393ccd2543fSMatthew G Knepley } 394ccd2543fSMatthew G Knepley if (invJ) { 395ccd2543fSMatthew G Knepley const PetscReal invDet = 1.0/(*detJ); 396ccd2543fSMatthew G Knepley 397ccd2543fSMatthew G Knepley invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 398ccd2543fSMatthew G Knepley invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 399ccd2543fSMatthew G Knepley invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 400ccd2543fSMatthew G Knepley invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 401ccd2543fSMatthew G Knepley invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 402ccd2543fSMatthew G Knepley invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 403ccd2543fSMatthew G Knepley invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 404ccd2543fSMatthew G Knepley invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 405ccd2543fSMatthew G Knepley invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 406ccd2543fSMatthew G Knepley PetscLogFlops(37.0); 407ccd2543fSMatthew G Knepley } 408ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 409ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 410ccd2543fSMatthew G Knepley } 411ccd2543fSMatthew G Knepley 412ccd2543fSMatthew G Knepley #undef __FUNCT__ 413ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 414ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 415ccd2543fSMatthew G Knepley { 416ccd2543fSMatthew G Knepley PetscSection coordSection; 417ccd2543fSMatthew G Knepley Vec coordinates; 418*7c1f9639SMatthew G Knepley PetscScalar *coords; 419ccd2543fSMatthew G Knepley const PetscInt dim = 3; 420ccd2543fSMatthew G Knepley PetscInt d; 421ccd2543fSMatthew G Knepley PetscErrorCode ierr; 422ccd2543fSMatthew G Knepley 423ccd2543fSMatthew G Knepley PetscFunctionBegin; 424ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 425ccd2543fSMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 426ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 427ccd2543fSMatthew G Knepley if (v0) { 428ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]); 429ccd2543fSMatthew G Knepley } 430ccd2543fSMatthew G Knepley if (J) { 431ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 432ccd2543fSMatthew G Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[(0+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 433ccd2543fSMatthew G Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[(1+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 434ccd2543fSMatthew G Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[(3+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 435ccd2543fSMatthew G Knepley } 436ccd2543fSMatthew G Knepley *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 437ccd2543fSMatthew G Knepley J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 438ccd2543fSMatthew G Knepley J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 439ccd2543fSMatthew G Knepley PetscLogFlops(18.0 + 12.0); 440ccd2543fSMatthew G Knepley } 441ccd2543fSMatthew G Knepley if (invJ) { 442ccd2543fSMatthew G Knepley const PetscReal invDet = -1.0/(*detJ); 443ccd2543fSMatthew G Knepley 444ccd2543fSMatthew G Knepley invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]); 445ccd2543fSMatthew G Knepley invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]); 446ccd2543fSMatthew G Knepley invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]); 447ccd2543fSMatthew G Knepley invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]); 448ccd2543fSMatthew G Knepley invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]); 449ccd2543fSMatthew G Knepley invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]); 450ccd2543fSMatthew G Knepley invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]); 451ccd2543fSMatthew G Knepley invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]); 452ccd2543fSMatthew G Knepley invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]); 453ccd2543fSMatthew G Knepley PetscLogFlops(37.0); 454ccd2543fSMatthew G Knepley } 455ccd2543fSMatthew G Knepley *detJ *= 8.0; 456ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 457ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 458ccd2543fSMatthew G Knepley } 459ccd2543fSMatthew G Knepley 460ccd2543fSMatthew G Knepley #undef __FUNCT__ 461ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeCellGeometry" 462ccd2543fSMatthew G Knepley /*@C 463ccd2543fSMatthew G Knepley DMPlexComputeCellGeometry - Compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 464ccd2543fSMatthew G Knepley 465ccd2543fSMatthew G Knepley Collective on DM 466ccd2543fSMatthew G Knepley 467ccd2543fSMatthew G Knepley Input Arguments: 468ccd2543fSMatthew G Knepley + dm - the DM 469ccd2543fSMatthew G Knepley - cell - the cell 470ccd2543fSMatthew G Knepley 471ccd2543fSMatthew G Knepley Output Arguments: 472ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 473ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 474ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 475ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 476ccd2543fSMatthew G Knepley 477ccd2543fSMatthew G Knepley Level: advanced 478ccd2543fSMatthew G Knepley 479ccd2543fSMatthew G Knepley Fortran Notes: 480ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 481ccd2543fSMatthew G Knepley include petsc.h90 in your code. 482ccd2543fSMatthew G Knepley 483ccd2543fSMatthew G Knepley .seealso: DMPlexGetCoordinateSection(), DMPlexGetCoordinateVec() 484ccd2543fSMatthew G Knepley @*/ 485ccd2543fSMatthew G Knepley PetscErrorCode DMPlexComputeCellGeometry(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 486ccd2543fSMatthew G Knepley { 487ccd2543fSMatthew G Knepley PetscInt dim, coneSize; 488ccd2543fSMatthew G Knepley PetscErrorCode ierr; 489ccd2543fSMatthew G Knepley 490ccd2543fSMatthew G Knepley PetscFunctionBegin; 491ccd2543fSMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 492ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 493ccd2543fSMatthew G Knepley switch (dim) { 494ccd2543fSMatthew G Knepley case 2: 495ccd2543fSMatthew G Knepley switch (coneSize) { 496ccd2543fSMatthew G Knepley case 3: 497ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 498ccd2543fSMatthew G Knepley break; 499ccd2543fSMatthew G Knepley case 4: 500ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 501ccd2543fSMatthew G Knepley break; 502ccd2543fSMatthew G Knepley default: 503ccd2543fSMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 504ccd2543fSMatthew G Knepley } 505ccd2543fSMatthew G Knepley break; 506ccd2543fSMatthew G Knepley case 3: 507ccd2543fSMatthew G Knepley switch (coneSize) { 508ccd2543fSMatthew G Knepley case 4: 509ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 510ccd2543fSMatthew G Knepley break; 511ccd2543fSMatthew G Knepley case 8: 512ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 513ccd2543fSMatthew G Knepley break; 514ccd2543fSMatthew G Knepley default: 515ccd2543fSMatthew G Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of vertices %D in cell %D for element geometry computation", coneSize, cell); 516ccd2543fSMatthew G Knepley } 517ccd2543fSMatthew G Knepley break; 518ccd2543fSMatthew G Knepley default: 519ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 520ccd2543fSMatthew G Knepley } 521ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 522ccd2543fSMatthew G Knepley } 523