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__ 22428dbe442SToby Isaac #define __FUNCT__ "DMPlexComputeProjection3Dto1D_Internal" 22528dbe442SToby Isaac /* 22628dbe442SToby Isaac DMPlexComputeProjection3Dto1D_Internal - Rewrite coordinates to be the 1D projection of the 3D 22728dbe442SToby Isaac 22828dbe442SToby Isaac This uses the basis completion described by Frisvad, 22928dbe442SToby Isaac 23028dbe442SToby Isaac http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html 23128dbe442SToby Isaac DOI:10.1080/2165347X.2012.689606 23228dbe442SToby Isaac */ 23328dbe442SToby Isaac static PetscErrorCode DMPlexComputeProjection3Dto1D_Internal(PetscScalar coords[], PetscReal R[]) 23428dbe442SToby Isaac { 23528dbe442SToby Isaac PetscReal x = PetscRealPart(coords[3] - coords[0]); 23628dbe442SToby Isaac PetscReal y = PetscRealPart(coords[4] - coords[1]); 23728dbe442SToby Isaac PetscReal z = PetscRealPart(coords[5] - coords[2]); 23828dbe442SToby Isaac PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 23928dbe442SToby Isaac PetscReal rinv = 1. / r; 24028dbe442SToby Isaac PetscFunctionBegin; 24128dbe442SToby Isaac 24228dbe442SToby Isaac x *= rinv; y *= rinv; z *= rinv; 24328dbe442SToby Isaac if (x > 0.) { 24428dbe442SToby Isaac PetscReal inv1pX = 1./ (1. + x); 24528dbe442SToby Isaac 24628dbe442SToby Isaac R[0] = x; R[1] = -y; R[2] = -z; 24728dbe442SToby Isaac R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 24828dbe442SToby Isaac R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 24928dbe442SToby Isaac } 25028dbe442SToby Isaac else { 25128dbe442SToby Isaac PetscReal inv1mX = 1./ (1. - x); 25228dbe442SToby Isaac 25328dbe442SToby Isaac R[0] = x; R[1] = z; R[2] = y; 25428dbe442SToby Isaac R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 25528dbe442SToby Isaac R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 25628dbe442SToby Isaac } 25728dbe442SToby Isaac coords[0] = 0.0; 25828dbe442SToby Isaac coords[1] = r; 25928dbe442SToby Isaac PetscFunctionReturn(0); 26028dbe442SToby Isaac } 26128dbe442SToby Isaac 26228dbe442SToby Isaac #undef __FUNCT__ 263ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeProjection3Dto2D_Internal" 264ccd2543fSMatthew G Knepley /* 265ccd2543fSMatthew G Knepley DMPlexComputeProjection3Dto2D_Internal - Rewrite coordinates to be the 2D projection of the 3D 266ccd2543fSMatthew G Knepley */ 26799dec3a6SMatthew G. Knepley static PetscErrorCode DMPlexComputeProjection3Dto2D_Internal(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 268ccd2543fSMatthew G Knepley { 2691ee9d5ecSMatthew G. Knepley PetscReal x1[3], x2[3], n[3], norm; 27099dec3a6SMatthew G. Knepley PetscReal x1p[3], x2p[3], xnp[3]; 2714a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 272ccd2543fSMatthew G Knepley const PetscInt dim = 3; 27399dec3a6SMatthew G. Knepley PetscInt d, e, p; 274ccd2543fSMatthew G Knepley 275ccd2543fSMatthew G Knepley PetscFunctionBegin; 276ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 277ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 2781ee9d5ecSMatthew G. Knepley x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 2791ee9d5ecSMatthew G. Knepley x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 280ccd2543fSMatthew G Knepley } 281ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 282ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 283ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 2848b49ba18SBarry Smith norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 285ccd2543fSMatthew G Knepley n[0] /= norm; 286ccd2543fSMatthew G Knepley n[1] /= norm; 287ccd2543fSMatthew G Knepley n[2] /= norm; 288ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 289ccd2543fSMatthew G Knepley 290ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 291ccd2543fSMatthew G Knepley 292ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 293ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 294ccd2543fSMatthew G Knepley \ nx ny nz / 295ccd2543fSMatthew G Knepley 296ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 297ccd2543fSMatthew G Knepley */ 2988b49ba18SBarry Smith sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 29973868372SMatthew G. Knepley /* Check for n = z */ 30073868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 3011ee9d5ecSMatthew G. Knepley if (n[2] < 0.0) { 30299dec3a6SMatthew G. Knepley if (coordSize > 9) { 30399dec3a6SMatthew G. Knepley coords[2] = PetscRealPart(coords[3*dim+0] - coords[0*dim+0]); 30408b58242SMatthew G. Knepley coords[3] = PetscRealPart(coords[3*dim+1] - coords[0*dim+1]); 30599dec3a6SMatthew G. Knepley coords[4] = x2[0]; 30699dec3a6SMatthew G. Knepley coords[5] = x2[1]; 30799dec3a6SMatthew G. Knepley coords[6] = x1[0]; 30899dec3a6SMatthew G. Knepley coords[7] = x1[1]; 30999dec3a6SMatthew G. Knepley } else { 31073868372SMatthew G. Knepley coords[2] = x2[0]; 31173868372SMatthew G. Knepley coords[3] = x2[1]; 31273868372SMatthew G. Knepley coords[4] = x1[0]; 31373868372SMatthew G. Knepley coords[5] = x1[1]; 31499dec3a6SMatthew G. Knepley } 315b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 316b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 317b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = -1.0; 31873868372SMatthew G. Knepley } else { 31999dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 32099dec3a6SMatthew G. Knepley coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 32199dec3a6SMatthew G. Knepley coords[p*2+1] = PetscRealPart(coords[p*dim+1] - coords[0*dim+1]); 32299dec3a6SMatthew G. Knepley } 32373868372SMatthew G. Knepley coords[2] = x1[0]; 32473868372SMatthew G. Knepley coords[3] = x1[1]; 32573868372SMatthew G. Knepley coords[4] = x2[0]; 32673868372SMatthew G. Knepley coords[5] = x2[1]; 327b7ad821dSMatthew G. Knepley R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 328b7ad821dSMatthew G. Knepley R[3] = 0.0; R[4] = 1.0; R[5] = 0.0; 329b7ad821dSMatthew G. Knepley R[6] = 0.0; R[7] = 0.0; R[8] = 1.0; 33073868372SMatthew G. Knepley } 33199dec3a6SMatthew G. Knepley coords[0] = 0.0; 33299dec3a6SMatthew G. Knepley coords[1] = 0.0; 33373868372SMatthew G. Knepley PetscFunctionReturn(0); 33473868372SMatthew G. Knepley } 335da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 336ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 337ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 338ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 339ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 340ccd2543fSMatthew G Knepley x1p[d] = 0.0; 341ccd2543fSMatthew G Knepley x2p[d] = 0.0; 342ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 343ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 344ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 345ccd2543fSMatthew G Knepley } 346ccd2543fSMatthew G Knepley } 3478763be8eSMatthew G. Knepley if (PetscAbsReal(x1p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 3488763be8eSMatthew G. Knepley if (PetscAbsReal(x2p[2]) > 1.0e-9) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 349ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 35099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 35199dec3a6SMatthew G. Knepley for (d = 0; d < dim; ++d) { 35299dec3a6SMatthew G. Knepley xnp[d] = 0.0; 35399dec3a6SMatthew G. Knepley for (e = 0; e < dim; ++e) { 35499dec3a6SMatthew G. Knepley xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 35599dec3a6SMatthew G. Knepley } 35699dec3a6SMatthew G. Knepley if (d < dim-1) coords[p*2+d] = xnp[d]; 35799dec3a6SMatthew G. Knepley } 35899dec3a6SMatthew G. Knepley } 359ccd2543fSMatthew G Knepley coords[0] = 0.0; 360ccd2543fSMatthew G Knepley coords[1] = 0.0; 361ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 362ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 363ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 364ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 3657f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 3667f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 3677f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 3687f07f362SMatthew G. Knepley PetscReal tmp; 3697f07f362SMatthew G. Knepley 3707f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 3717f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 3727f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 3737f07f362SMatthew G. Knepley } 3747f07f362SMatthew G. Knepley } 375ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 376ccd2543fSMatthew G Knepley } 377ccd2543fSMatthew G Knepley 378ccd2543fSMatthew G Knepley #undef __FUNCT__ 379834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Internal" 3806322fe33SJed Brown PETSC_UNUSED 381834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 382834e62ceSMatthew G. Knepley { 383834e62ceSMatthew G. Knepley /* Signed volume is 1/2 the determinant 384834e62ceSMatthew G. Knepley 385834e62ceSMatthew G. Knepley | 1 1 1 | 386834e62ceSMatthew G. Knepley | x0 x1 x2 | 387834e62ceSMatthew G. Knepley | y0 y1 y2 | 388834e62ceSMatthew G. Knepley 389834e62ceSMatthew G. Knepley but if x0,y0 is the origin, we have 390834e62ceSMatthew G. Knepley 391834e62ceSMatthew G. Knepley | x1 x2 | 392834e62ceSMatthew G. Knepley | y1 y2 | 393834e62ceSMatthew G. Knepley */ 394834e62ceSMatthew G. Knepley const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 395834e62ceSMatthew G. Knepley const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 396834e62ceSMatthew G. Knepley PetscReal M[4], detM; 397834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; 39886623015SMatthew G. Knepley M[2] = y1; M[3] = y2; 399923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(&detM, M); 400834e62ceSMatthew G. Knepley *vol = 0.5*detM; 401834e62ceSMatthew G. Knepley PetscLogFlops(5.0); 402834e62ceSMatthew G. Knepley } 403834e62ceSMatthew G. Knepley 404834e62ceSMatthew G. Knepley #undef __FUNCT__ 405834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Triangle_Origin_Internal" 406834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 407834e62ceSMatthew G. Knepley { 408923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(vol, coords); 409834e62ceSMatthew G. Knepley *vol *= 0.5; 410834e62ceSMatthew G. Knepley } 411834e62ceSMatthew G. Knepley 412834e62ceSMatthew G. Knepley #undef __FUNCT__ 413834e62ceSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Internal" 4146322fe33SJed Brown PETSC_UNUSED 415834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 416834e62ceSMatthew G. Knepley { 417834e62ceSMatthew G. Knepley /* Signed volume is 1/6th of the determinant 418834e62ceSMatthew G. Knepley 419834e62ceSMatthew G. Knepley | 1 1 1 1 | 420834e62ceSMatthew G. Knepley | x0 x1 x2 x3 | 421834e62ceSMatthew G. Knepley | y0 y1 y2 y3 | 422834e62ceSMatthew G. Knepley | z0 z1 z2 z3 | 423834e62ceSMatthew G. Knepley 424834e62ceSMatthew G. Knepley but if x0,y0,z0 is the origin, we have 425834e62ceSMatthew G. Knepley 426834e62ceSMatthew G. Knepley | x1 x2 x3 | 427834e62ceSMatthew G. Knepley | y1 y2 y3 | 428834e62ceSMatthew G. Knepley | z1 z2 z3 | 429834e62ceSMatthew G. Knepley */ 430834e62ceSMatthew G. Knepley const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 431834e62ceSMatthew G. Knepley const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 432834e62ceSMatthew G. Knepley const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 433834e62ceSMatthew G. Knepley PetscReal M[9], detM; 434834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; M[2] = x3; 435834e62ceSMatthew G. Knepley M[3] = y1; M[4] = y2; M[5] = y3; 436834e62ceSMatthew G. Knepley M[6] = z1; M[7] = z2; M[8] = z3; 437923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(&detM, M); 438b7ad821dSMatthew G. Knepley *vol = -0.16666666666666666666666*detM; 439834e62ceSMatthew G. Knepley PetscLogFlops(10.0); 440834e62ceSMatthew G. Knepley } 441834e62ceSMatthew G. Knepley 442834e62ceSMatthew G. Knepley #undef __FUNCT__ 4430ec8681fSMatthew G. Knepley #define __FUNCT__ "Volume_Tetrahedron_Origin_Internal" 4440ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 4450ec8681fSMatthew G. Knepley { 446923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(vol, coords); 447b7ad821dSMatthew G. Knepley *vol *= -0.16666666666666666666666; 4480ec8681fSMatthew G. Knepley } 4490ec8681fSMatthew G. Knepley 4500ec8681fSMatthew G. Knepley #undef __FUNCT__ 45117fe8556SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeLineGeometry_Internal" 45217fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 45317fe8556SMatthew G. Knepley { 45417fe8556SMatthew G. Knepley PetscSection coordSection; 45517fe8556SMatthew G. Knepley Vec coordinates; 456a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 4577f07f362SMatthew G. Knepley PetscInt numCoords, d; 45817fe8556SMatthew G. Knepley PetscErrorCode ierr; 45917fe8556SMatthew G. Knepley 46017fe8556SMatthew G. Knepley PetscFunctionBegin; 46117fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 46269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 46317fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 4647f07f362SMatthew G. Knepley *detJ = 0.0; 46528dbe442SToby Isaac if (numCoords == 6) { 46628dbe442SToby Isaac const PetscInt dim = 3; 46728dbe442SToby Isaac PetscReal R[9], J0; 46828dbe442SToby Isaac 46928dbe442SToby Isaac if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 47028dbe442SToby Isaac ierr = DMPlexComputeProjection3Dto1D_Internal(coords, R);CHKERRQ(ierr); 47128dbe442SToby Isaac if (J) { 47228dbe442SToby Isaac J0 = 0.5*PetscRealPart(coords[1]); 47328dbe442SToby Isaac J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 47428dbe442SToby Isaac J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 47528dbe442SToby Isaac J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 47628dbe442SToby Isaac DMPlex_Det3D_Internal(detJ, J); 47728dbe442SToby Isaac } 47828dbe442SToby Isaac if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 47928dbe442SToby Isaac } else if (numCoords == 4) { 4807f07f362SMatthew G. Knepley const PetscInt dim = 2; 4817f07f362SMatthew G. Knepley PetscReal R[4], J0; 4827f07f362SMatthew G. Knepley 4837f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 4847f07f362SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D_Internal(coords, R);CHKERRQ(ierr); 48517fe8556SMatthew G. Knepley if (J) { 4867f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 4877f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 4887f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 489923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 49017fe8556SMatthew G. Knepley } 491923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 4927f07f362SMatthew G. Knepley } else if (numCoords == 2) { 4937f07f362SMatthew G. Knepley const PetscInt dim = 1; 4947f07f362SMatthew G. Knepley 4957f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 4967f07f362SMatthew G. Knepley if (J) { 4977f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 49817fe8556SMatthew G. Knepley *detJ = J[0]; 49917fe8556SMatthew G. Knepley PetscLogFlops(2.0); 50017fe8556SMatthew G. Knepley } 5017f07f362SMatthew G. Knepley if (invJ) {invJ[0] = 1.0/J[0]; PetscLogFlops(1.0);} 502796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 50317fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 50417fe8556SMatthew G. Knepley PetscFunctionReturn(0); 50517fe8556SMatthew G. Knepley } 50617fe8556SMatthew G. Knepley 50717fe8556SMatthew G. Knepley #undef __FUNCT__ 508ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTriangleGeometry_Internal" 509ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 510ccd2543fSMatthew G Knepley { 511ccd2543fSMatthew G Knepley PetscSection coordSection; 512ccd2543fSMatthew G Knepley Vec coordinates; 513a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 5147f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 515ccd2543fSMatthew G Knepley PetscErrorCode ierr; 516ccd2543fSMatthew G Knepley 517ccd2543fSMatthew G Knepley PetscFunctionBegin; 518ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 51969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 520ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 5217f07f362SMatthew G. Knepley *detJ = 0.0; 522ccd2543fSMatthew G Knepley if (numCoords == 9) { 5237f07f362SMatthew G. Knepley const PetscInt dim = 3; 5247f07f362SMatthew 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}; 5257f07f362SMatthew G. Knepley 5267f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 52799dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 5287f07f362SMatthew G. Knepley if (J) { 529b7ad821dSMatthew G. Knepley const PetscInt pdim = 2; 530b7ad821dSMatthew G. Knepley 531b7ad821dSMatthew G. Knepley for (d = 0; d < pdim; d++) { 532b7ad821dSMatthew G. Knepley for (f = 0; f < pdim; f++) { 533b7ad821dSMatthew G. Knepley J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 534ccd2543fSMatthew G Knepley } 5357f07f362SMatthew G. Knepley } 5367f07f362SMatthew G. Knepley PetscLogFlops(8.0); 537923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 5387f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 5397f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 5407f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 5417f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 5427f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 5437f07f362SMatthew G. Knepley } 5447f07f362SMatthew G. Knepley } 5457f07f362SMatthew G. Knepley } 5467f07f362SMatthew G. Knepley PetscLogFlops(18.0); 5477f07f362SMatthew G. Knepley } 548923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 5497f07f362SMatthew G. Knepley } else if (numCoords == 6) { 5507f07f362SMatthew G. Knepley const PetscInt dim = 2; 5517f07f362SMatthew G. Knepley 5527f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 553ccd2543fSMatthew G Knepley if (J) { 554ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 555ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 556ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 557ccd2543fSMatthew G Knepley } 558ccd2543fSMatthew G Knepley } 5597f07f362SMatthew G. Knepley PetscLogFlops(8.0); 560923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 561ccd2543fSMatthew G Knepley } 562923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 563796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 564ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 565ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 566ccd2543fSMatthew G Knepley } 567ccd2543fSMatthew G Knepley 568ccd2543fSMatthew G Knepley #undef __FUNCT__ 569ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeRectangleGeometry_Internal" 570ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 571ccd2543fSMatthew G Knepley { 572ccd2543fSMatthew G Knepley PetscSection coordSection; 573ccd2543fSMatthew G Knepley Vec coordinates; 574a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 57599dec3a6SMatthew G. Knepley PetscInt numCoords, d, f, g; 576ccd2543fSMatthew G Knepley PetscErrorCode ierr; 577ccd2543fSMatthew G Knepley 578ccd2543fSMatthew G Knepley PetscFunctionBegin; 579ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 58069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 58199dec3a6SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 5827f07f362SMatthew G. Knepley *detJ = 0.0; 58399dec3a6SMatthew G. Knepley if (numCoords == 12) { 58499dec3a6SMatthew G. Knepley const PetscInt dim = 3; 58599dec3a6SMatthew 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}; 58699dec3a6SMatthew G. Knepley 58799dec3a6SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 58899dec3a6SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D_Internal(numCoords, coords, R);CHKERRQ(ierr); 58999dec3a6SMatthew G. Knepley if (J) { 59099dec3a6SMatthew G. Knepley const PetscInt pdim = 2; 59199dec3a6SMatthew G. Knepley 59299dec3a6SMatthew G. Knepley for (d = 0; d < pdim; d++) { 59399dec3a6SMatthew G. Knepley J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 59499dec3a6SMatthew G. Knepley J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 59599dec3a6SMatthew G. Knepley } 59699dec3a6SMatthew G. Knepley PetscLogFlops(8.0); 597923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 59899dec3a6SMatthew G. Knepley for (d = 0; d < dim; d++) { 59999dec3a6SMatthew G. Knepley for (f = 0; f < dim; f++) { 60099dec3a6SMatthew G. Knepley J[d*dim+f] = 0.0; 60199dec3a6SMatthew G. Knepley for (g = 0; g < dim; g++) { 60299dec3a6SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 60399dec3a6SMatthew G. Knepley } 60499dec3a6SMatthew G. Knepley } 60599dec3a6SMatthew G. Knepley } 60699dec3a6SMatthew G. Knepley PetscLogFlops(18.0); 60799dec3a6SMatthew G. Knepley } 608923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 60997f1a218SMatthew G. Knepley } else if ((numCoords == 8) || (numCoords == 16)) { 61099dec3a6SMatthew G. Knepley const PetscInt dim = 2; 61199dec3a6SMatthew G. Knepley 6127f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 613ccd2543fSMatthew G Knepley if (J) { 614ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 61599dec3a6SMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 61699dec3a6SMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 617ccd2543fSMatthew G Knepley } 6187f07f362SMatthew G. Knepley PetscLogFlops(8.0); 619923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 620ccd2543fSMatthew G Knepley } 621923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 622796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 62399dec3a6SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 624ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 625ccd2543fSMatthew G Knepley } 626ccd2543fSMatthew G Knepley 627ccd2543fSMatthew G Knepley #undef __FUNCT__ 628ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeTetrahedronGeometry_Internal" 629ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 630ccd2543fSMatthew G Knepley { 631ccd2543fSMatthew G Knepley PetscSection coordSection; 632ccd2543fSMatthew G Knepley Vec coordinates; 633a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 634ccd2543fSMatthew G Knepley const PetscInt dim = 3; 63599dec3a6SMatthew G. Knepley PetscInt d; 636ccd2543fSMatthew G Knepley PetscErrorCode ierr; 637ccd2543fSMatthew G Knepley 638ccd2543fSMatthew G Knepley PetscFunctionBegin; 639ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 64069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 641ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 6427f07f362SMatthew G. Knepley *detJ = 0.0; 6437f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 644ccd2543fSMatthew G Knepley if (J) { 645ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 646f0df753eSMatthew G. Knepley /* I orient with outward face normals */ 647f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 648f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 649f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 650ccd2543fSMatthew G Knepley } 6517f07f362SMatthew G. Knepley PetscLogFlops(18.0); 652923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 653ccd2543fSMatthew G Knepley } 654923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 655ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 656ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 657ccd2543fSMatthew G Knepley } 658ccd2543fSMatthew G Knepley 659ccd2543fSMatthew G Knepley #undef __FUNCT__ 660ccd2543fSMatthew G Knepley #define __FUNCT__ "DMPlexComputeHexahedronGeometry_Internal" 661ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 662ccd2543fSMatthew G Knepley { 663ccd2543fSMatthew G Knepley PetscSection coordSection; 664ccd2543fSMatthew G Knepley Vec coordinates; 665a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 666ccd2543fSMatthew G Knepley const PetscInt dim = 3; 667ccd2543fSMatthew G Knepley PetscInt d; 668ccd2543fSMatthew G Knepley PetscErrorCode ierr; 669ccd2543fSMatthew G Knepley 670ccd2543fSMatthew G Knepley PetscFunctionBegin; 671ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 67269d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 673ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 6747f07f362SMatthew G. Knepley *detJ = 0.0; 6757f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 676ccd2543fSMatthew G Knepley if (J) { 677ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 678f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 679f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 680f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 681ccd2543fSMatthew G Knepley } 6827f07f362SMatthew G. Knepley PetscLogFlops(18.0); 683923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 684ccd2543fSMatthew G Knepley } 685923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 686ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 687ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 688ccd2543fSMatthew G Knepley } 689ccd2543fSMatthew G Knepley 690ccd2543fSMatthew G Knepley #undef __FUNCT__ 6918e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryAffineFEM" 692ccd2543fSMatthew G Knepley /*@C 6938e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 694ccd2543fSMatthew G Knepley 695ccd2543fSMatthew G Knepley Collective on DM 696ccd2543fSMatthew G Knepley 697ccd2543fSMatthew G Knepley Input Arguments: 698ccd2543fSMatthew G Knepley + dm - the DM 699ccd2543fSMatthew G Knepley - cell - the cell 700ccd2543fSMatthew G Knepley 701ccd2543fSMatthew G Knepley Output Arguments: 702ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 703ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 704ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 705ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 706ccd2543fSMatthew G Knepley 707ccd2543fSMatthew G Knepley Level: advanced 708ccd2543fSMatthew G Knepley 709ccd2543fSMatthew G Knepley Fortran Notes: 710ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 711ccd2543fSMatthew G Knepley include petsc.h90 in your code. 712ccd2543fSMatthew G Knepley 7138e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 714ccd2543fSMatthew G Knepley @*/ 7158e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 716ccd2543fSMatthew G Knepley { 71749dc4407SMatthew G. Knepley PetscInt depth, dim, coneSize; 718ccd2543fSMatthew G Knepley PetscErrorCode ierr; 719ccd2543fSMatthew G Knepley 720ccd2543fSMatthew G Knepley PetscFunctionBegin; 721139a35ccSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 722ccd2543fSMatthew G Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 72349dc4407SMatthew G. Knepley if (depth == 1) { 7248e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 7258e0841e0SMatthew G. Knepley } else { 7268e0841e0SMatthew G. Knepley DMLabel depth; 7278e0841e0SMatthew G. Knepley 7288e0841e0SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 7298e0841e0SMatthew G. Knepley ierr = DMLabelGetValue(depth, cell, &dim);CHKERRQ(ierr); 7308e0841e0SMatthew G. Knepley } 731ccd2543fSMatthew G Knepley switch (dim) { 73217fe8556SMatthew G. Knepley case 1: 73317fe8556SMatthew G. Knepley ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 73417fe8556SMatthew G. Knepley break; 735ccd2543fSMatthew G Knepley case 2: 736ccd2543fSMatthew G Knepley switch (coneSize) { 737ccd2543fSMatthew G Knepley case 3: 738ccd2543fSMatthew G Knepley ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 739ccd2543fSMatthew G Knepley break; 740ccd2543fSMatthew G Knepley case 4: 741ccd2543fSMatthew G Knepley ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 742ccd2543fSMatthew G Knepley break; 743ccd2543fSMatthew G Knepley default: 7448e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 745ccd2543fSMatthew G Knepley } 746ccd2543fSMatthew G Knepley break; 747ccd2543fSMatthew G Knepley case 3: 748ccd2543fSMatthew G Knepley switch (coneSize) { 749ccd2543fSMatthew G Knepley case 4: 750ccd2543fSMatthew G Knepley ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 751ccd2543fSMatthew G Knepley break; 7528e0841e0SMatthew G. Knepley case 6: /* Faces */ 7538e0841e0SMatthew G. Knepley case 8: /* Vertices */ 754ccd2543fSMatthew G Knepley ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr); 755ccd2543fSMatthew G Knepley break; 756ccd2543fSMatthew G Knepley default: 7578e0841e0SMatthew G. Knepley SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 758ccd2543fSMatthew G Knepley } 759ccd2543fSMatthew G Knepley break; 760ccd2543fSMatthew G Knepley default: 761ccd2543fSMatthew G Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 762ccd2543fSMatthew G Knepley } 7638e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 7648e0841e0SMatthew G. Knepley } 7658e0841e0SMatthew G. Knepley 7668e0841e0SMatthew G. Knepley #undef __FUNCT__ 7678e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIsoparametricGeometry_Internal" 7688e0841e0SMatthew G. Knepley static PetscErrorCode DMPlexComputeIsoparametricGeometry_Internal(DM dm, PetscFE fe, PetscInt point, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 7698e0841e0SMatthew G. Knepley { 7708e0841e0SMatthew G. Knepley PetscQuadrature quad; 7718e0841e0SMatthew G. Knepley PetscSection coordSection; 7728e0841e0SMatthew G. Knepley Vec coordinates; 7738e0841e0SMatthew G. Knepley PetscScalar *coords = NULL; 7748e0841e0SMatthew G. Knepley const PetscReal *quadPoints; 7758e0841e0SMatthew G. Knepley PetscReal *basisDer; 7768e0841e0SMatthew G. Knepley PetscInt dim, cdim, pdim, qdim, Nq, numCoords, d, q; 7778e0841e0SMatthew G. Knepley PetscErrorCode ierr; 7788e0841e0SMatthew G. Knepley 7798e0841e0SMatthew G. Knepley PetscFunctionBegin; 7808e0841e0SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 7818e0841e0SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 7828e0841e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 7838e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 7848e0841e0SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 7858e0841e0SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 786954b1791SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 7878e0841e0SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 7888e0841e0SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 7898e0841e0SMatthew G. Knepley *detJ = 0.0; 7908e0841e0SMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 7918e0841e0SMatthew 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); 7928e0841e0SMatthew G. Knepley if (v0) {for (d = 0; d < cdim; d++) v0[d] = PetscRealPart(coords[d]);} 7938e0841e0SMatthew G. Knepley if (J) { 7948e0841e0SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 7958e0841e0SMatthew G. Knepley PetscInt i, j, k, c, r; 7968e0841e0SMatthew G. Knepley 7978e0841e0SMatthew G. Knepley /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 7988e0841e0SMatthew G. Knepley for (k = 0; k < pdim; ++k) 7998e0841e0SMatthew G. Knepley for (j = 0; j < dim; ++j) 8008e0841e0SMatthew G. Knepley for (i = 0; i < cdim; ++i) 80171d6e60fSMatthew G. Knepley J[(q*cdim + i)*dim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 8028e0841e0SMatthew G. Knepley PetscLogFlops(2.0*pdim*dim*cdim); 8038e0841e0SMatthew G. Knepley if (cdim > dim) { 8048e0841e0SMatthew G. Knepley for (c = dim; c < cdim; ++c) 8058e0841e0SMatthew G. Knepley for (r = 0; r < cdim; ++r) 8068e0841e0SMatthew G. Knepley J[r*cdim+c] = r == c ? 1.0 : 0.0; 8078e0841e0SMatthew G. Knepley } 8088e0841e0SMatthew G. Knepley switch (cdim) { 8098e0841e0SMatthew G. Knepley case 3: 8108e0841e0SMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 8118e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 81217fe8556SMatthew G. Knepley break; 81349dc4407SMatthew G. Knepley case 2: 8148e0841e0SMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 8158e0841e0SMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 81649dc4407SMatthew G. Knepley break; 8178e0841e0SMatthew G. Knepley case 1: 8188e0841e0SMatthew G. Knepley *detJ = J[0]; 8198e0841e0SMatthew G. Knepley if (invJ) invJ[0] = 1.0/J[0]; 82049dc4407SMatthew G. Knepley } 82149dc4407SMatthew G. Knepley } 8228e0841e0SMatthew G. Knepley } 8238e0841e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 8248e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 8258e0841e0SMatthew G. Knepley } 8268e0841e0SMatthew G. Knepley 8278e0841e0SMatthew G. Knepley #undef __FUNCT__ 8288e0841e0SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFEM" 8298e0841e0SMatthew G. Knepley /*@C 8308e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 8318e0841e0SMatthew G. Knepley 8328e0841e0SMatthew G. Knepley Collective on DM 8338e0841e0SMatthew G. Knepley 8348e0841e0SMatthew G. Knepley Input Arguments: 8358e0841e0SMatthew G. Knepley + dm - the DM 8368e0841e0SMatthew G. Knepley . cell - the cell 8378e0841e0SMatthew G. Knepley - fe - the finite element containing the quadrature 8388e0841e0SMatthew G. Knepley 8398e0841e0SMatthew G. Knepley Output Arguments: 8408e0841e0SMatthew G. Knepley + v0 - the translation part of this transform 8418e0841e0SMatthew G. Knepley . J - the Jacobian of the transform from the reference element at each quadrature point 8428e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point 8438e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point 8448e0841e0SMatthew G. Knepley 8458e0841e0SMatthew G. Knepley Level: advanced 8468e0841e0SMatthew G. Knepley 8478e0841e0SMatthew G. Knepley Fortran Notes: 8488e0841e0SMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 8498e0841e0SMatthew G. Knepley include petsc.h90 in your code. 8508e0841e0SMatthew G. Knepley 8518e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 8528e0841e0SMatthew G. Knepley @*/ 8538e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscFE fe, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 8548e0841e0SMatthew G. Knepley { 8558e0841e0SMatthew G. Knepley PetscErrorCode ierr; 8568e0841e0SMatthew G. Knepley 8578e0841e0SMatthew G. Knepley PetscFunctionBegin; 8588e0841e0SMatthew G. Knepley if (!fe) {ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 8598e0841e0SMatthew G. Knepley else {ierr = DMPlexComputeIsoparametricGeometry_Internal(dm, fe, cell, v0, J, invJ, detJ);CHKERRQ(ierr);} 860ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 861ccd2543fSMatthew G Knepley } 862834e62ceSMatthew G. Knepley 863834e62ceSMatthew G. Knepley #undef __FUNCT__ 864cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_1D_Internal" 865011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 866cc08537eSMatthew G. Knepley { 867cc08537eSMatthew G. Knepley PetscSection coordSection; 868cc08537eSMatthew G. Knepley Vec coordinates; 869a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 870cc08537eSMatthew G. Knepley PetscInt coordSize; 871cc08537eSMatthew G. Knepley PetscErrorCode ierr; 872cc08537eSMatthew G. Knepley 873cc08537eSMatthew G. Knepley PetscFunctionBegin; 874cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 87569d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 876cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 877011ea5d8SMatthew G. Knepley if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 878cc08537eSMatthew G. Knepley if (centroid) { 8791ee9d5ecSMatthew G. Knepley centroid[0] = 0.5*PetscRealPart(coords[0] + coords[dim+0]); 8801ee9d5ecSMatthew G. Knepley centroid[1] = 0.5*PetscRealPart(coords[1] + coords[dim+1]); 881cc08537eSMatthew G. Knepley } 882cc08537eSMatthew G. Knepley if (normal) { 883a60a936bSMatthew G. Knepley PetscReal norm; 884a60a936bSMatthew G. Knepley 8850194f449SMatthew G. Knepley normal[0] = -PetscRealPart(coords[1] - coords[dim+1]); 8860194f449SMatthew G. Knepley normal[1] = PetscRealPart(coords[0] - coords[dim+0]); 887a60a936bSMatthew G. Knepley norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 888a60a936bSMatthew G. Knepley normal[0] /= norm; 889a60a936bSMatthew G. Knepley normal[1] /= norm; 890cc08537eSMatthew G. Knepley } 891cc08537eSMatthew G. Knepley if (vol) { 8928b49ba18SBarry Smith *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - coords[dim+0])) + PetscSqr(PetscRealPart(coords[1] - coords[dim+1]))); 893cc08537eSMatthew G. Knepley } 894cc08537eSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 895cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 896cc08537eSMatthew G. Knepley } 897cc08537eSMatthew G. Knepley 898cc08537eSMatthew G. Knepley #undef __FUNCT__ 899cc08537eSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_2D_Internal" 900cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 901011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 902cc08537eSMatthew G. Knepley { 903cc08537eSMatthew G. Knepley PetscSection coordSection; 904cc08537eSMatthew G. Knepley Vec coordinates; 905cc08537eSMatthew G. Knepley PetscScalar *coords = NULL; 9060a1d6728SMatthew G. Knepley PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 9070a1d6728SMatthew G. Knepley PetscInt tdim = 2, coordSize, numCorners, p, d, e; 908cc08537eSMatthew G. Knepley PetscErrorCode ierr; 909cc08537eSMatthew G. Knepley 910cc08537eSMatthew G. Knepley PetscFunctionBegin; 911cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 9120a1d6728SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 91369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 914cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 9150bce18caSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 916011ea5d8SMatthew G. Knepley if (normal) { 917011ea5d8SMatthew G. Knepley if (dim > 2) { 9181ee9d5ecSMatthew G. Knepley const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 9191ee9d5ecSMatthew G. Knepley const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 9201ee9d5ecSMatthew G. Knepley const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 9210a1d6728SMatthew G. Knepley PetscReal norm; 9220a1d6728SMatthew G. Knepley 9231ee9d5ecSMatthew G. Knepley v0[0] = PetscRealPart(coords[0]); 9241ee9d5ecSMatthew G. Knepley v0[1] = PetscRealPart(coords[1]); 9251ee9d5ecSMatthew G. Knepley v0[2] = PetscRealPart(coords[2]); 9260a1d6728SMatthew G. Knepley normal[0] = y0*z1 - z0*y1; 9270a1d6728SMatthew G. Knepley normal[1] = z0*x1 - x0*z1; 9280a1d6728SMatthew G. Knepley normal[2] = x0*y1 - y0*x1; 9298b49ba18SBarry Smith norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 9300a1d6728SMatthew G. Knepley normal[0] /= norm; 9310a1d6728SMatthew G. Knepley normal[1] /= norm; 9320a1d6728SMatthew G. Knepley normal[2] /= norm; 933011ea5d8SMatthew G. Knepley } else { 934011ea5d8SMatthew G. Knepley for (d = 0; d < dim; ++d) normal[d] = 0.0; 935011ea5d8SMatthew G. Knepley } 936011ea5d8SMatthew G. Knepley } 93799dec3a6SMatthew G. Knepley if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D_Internal(coordSize, coords, R);CHKERRQ(ierr);} 9380a1d6728SMatthew G. Knepley for (p = 0; p < numCorners; ++p) { 9390a1d6728SMatthew G. Knepley /* Need to do this copy to get types right */ 9400a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 9411ee9d5ecSMatthew G. Knepley ctmp[d] = PetscRealPart(coords[p*tdim+d]); 9421ee9d5ecSMatthew G. Knepley ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 9430a1d6728SMatthew G. Knepley } 9440a1d6728SMatthew G. Knepley Volume_Triangle_Origin_Internal(&vtmp, ctmp); 9450a1d6728SMatthew G. Knepley vsum += vtmp; 9460a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 9470a1d6728SMatthew G. Knepley csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 9480a1d6728SMatthew G. Knepley } 9490a1d6728SMatthew G. Knepley } 9500a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 9510a1d6728SMatthew G. Knepley csum[d] /= (tdim+1)*vsum; 9520a1d6728SMatthew G. Knepley } 9530a1d6728SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 954ee6bbdb2SSatish Balay if (vol) *vol = PetscAbsReal(vsum); 9550a1d6728SMatthew G. Knepley if (centroid) { 9560a1d6728SMatthew G. Knepley if (dim > 2) { 9570a1d6728SMatthew G. Knepley for (d = 0; d < dim; ++d) { 9580a1d6728SMatthew G. Knepley centroid[d] = v0[d]; 9590a1d6728SMatthew G. Knepley for (e = 0; e < dim; ++e) { 9600a1d6728SMatthew G. Knepley centroid[d] += R[d*dim+e]*csum[e]; 9610a1d6728SMatthew G. Knepley } 9620a1d6728SMatthew G. Knepley } 9630a1d6728SMatthew G. Knepley } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 9640a1d6728SMatthew G. Knepley } 965cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 966cc08537eSMatthew G. Knepley } 967cc08537eSMatthew G. Knepley 968cc08537eSMatthew G. Knepley #undef __FUNCT__ 9690ec8681fSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM_3D_Internal" 9700ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 971011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 9720ec8681fSMatthew G. Knepley { 9730ec8681fSMatthew G. Knepley PetscSection coordSection; 9740ec8681fSMatthew G. Knepley Vec coordinates; 9750ec8681fSMatthew G. Knepley PetscScalar *coords = NULL; 97686623015SMatthew G. Knepley PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 977a7df9edeSMatthew G. Knepley const PetscInt *faces, *facesO; 9780ec8681fSMatthew G. Knepley PetscInt numFaces, f, coordSize, numCorners, p, d; 9790ec8681fSMatthew G. Knepley PetscErrorCode ierr; 9800ec8681fSMatthew G. Knepley 9810ec8681fSMatthew G. Knepley PetscFunctionBegin; 9820ec8681fSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 98369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 9840ec8681fSMatthew G. Knepley 985d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 9860ec8681fSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 9870ec8681fSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 988a7df9edeSMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 9890ec8681fSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 990011ea5d8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 9910ec8681fSMatthew G. Knepley numCorners = coordSize/dim; 9920ec8681fSMatthew G. Knepley switch (numCorners) { 9930ec8681fSMatthew G. Knepley case 3: 9940ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 9951ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 9961ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 9971ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 9980ec8681fSMatthew G. Knepley } 9990ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1000a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 10010ec8681fSMatthew G. Knepley vsum += vtmp; 10024f25033aSJed Brown if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 10030ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 10041ee9d5ecSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 10050ec8681fSMatthew G. Knepley } 10060ec8681fSMatthew G. Knepley } 10070ec8681fSMatthew G. Knepley break; 10080ec8681fSMatthew G. Knepley case 4: 10090ec8681fSMatthew G. Knepley /* DO FOR PYRAMID */ 10100ec8681fSMatthew G. Knepley /* First tet */ 10110ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 10121ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 10131ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 10141ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 10150ec8681fSMatthew G. Knepley } 10160ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1017a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 10180ec8681fSMatthew G. Knepley vsum += vtmp; 10190ec8681fSMatthew G. Knepley if (centroid) { 10200ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 10210ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 10220ec8681fSMatthew G. Knepley } 10230ec8681fSMatthew G. Knepley } 10240ec8681fSMatthew G. Knepley /* Second tet */ 10250ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 10261ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 10271ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 10281ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 10290ec8681fSMatthew G. Knepley } 10300ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1031a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 10320ec8681fSMatthew G. Knepley vsum += vtmp; 10330ec8681fSMatthew G. Knepley if (centroid) { 10340ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 10350ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 10360ec8681fSMatthew G. Knepley } 10370ec8681fSMatthew G. Knepley } 10380ec8681fSMatthew G. Knepley break; 10390ec8681fSMatthew G. Knepley default: 1040796f034aSJed Brown SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 10410ec8681fSMatthew G. Knepley } 10424f25033aSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 10430ec8681fSMatthew G. Knepley } 10448763be8eSMatthew G. Knepley if (vol) *vol = PetscAbsReal(vsum); 10450ec8681fSMatthew G. Knepley if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1046d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 10470ec8681fSMatthew G. Knepley PetscFunctionReturn(0); 10480ec8681fSMatthew G. Knepley } 10490ec8681fSMatthew G. Knepley 10500ec8681fSMatthew G. Knepley #undef __FUNCT__ 1051834e62ceSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeCellGeometryFVM" 1052834e62ceSMatthew G. Knepley /*@C 1053834e62ceSMatthew G. Knepley DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1054834e62ceSMatthew G. Knepley 1055834e62ceSMatthew G. Knepley Collective on DM 1056834e62ceSMatthew G. Knepley 1057834e62ceSMatthew G. Knepley Input Arguments: 1058834e62ceSMatthew G. Knepley + dm - the DM 1059834e62ceSMatthew G. Knepley - cell - the cell 1060834e62ceSMatthew G. Knepley 1061834e62ceSMatthew G. Knepley Output Arguments: 1062834e62ceSMatthew G. Knepley + volume - the cell volume 1063cc08537eSMatthew G. Knepley . centroid - the cell centroid 1064cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate 1065834e62ceSMatthew G. Knepley 1066834e62ceSMatthew G. Knepley Level: advanced 1067834e62ceSMatthew G. Knepley 1068834e62ceSMatthew G. Knepley Fortran Notes: 1069834e62ceSMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1070834e62ceSMatthew G. Knepley include petsc.h90 in your code. 1071834e62ceSMatthew G. Knepley 107269d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1073834e62ceSMatthew G. Knepley @*/ 1074cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1075834e62ceSMatthew G. Knepley { 10760ec8681fSMatthew G. Knepley PetscInt depth, dim; 1077834e62ceSMatthew G. Knepley PetscErrorCode ierr; 1078834e62ceSMatthew G. Knepley 1079834e62ceSMatthew G. Knepley PetscFunctionBegin; 1080834e62ceSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1081c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1082834e62ceSMatthew G. Knepley if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1083834e62ceSMatthew G. Knepley /* We need to keep a pointer to the depth label */ 1084011ea5d8SMatthew G. Knepley ierr = DMPlexGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1085834e62ceSMatthew G. Knepley /* Cone size is now the number of faces */ 1086011ea5d8SMatthew G. Knepley switch (depth) { 1087cc08537eSMatthew G. Knepley case 1: 1088011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1089cc08537eSMatthew G. Knepley break; 1090834e62ceSMatthew G. Knepley case 2: 1091011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1092834e62ceSMatthew G. Knepley break; 1093834e62ceSMatthew G. Knepley case 3: 1094011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1095834e62ceSMatthew G. Knepley break; 1096834e62ceSMatthew G. Knepley default: 1097834e62ceSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1098834e62ceSMatthew G. Knepley } 1099834e62ceSMatthew G. Knepley PetscFunctionReturn(0); 1100834e62ceSMatthew G. Knepley } 1101113c68e6SMatthew G. Knepley 1102113c68e6SMatthew G. Knepley #undef __FUNCT__ 1103c0d900a5SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFEM" 1104c0d900a5SMatthew G. Knepley /* This should also take a PetscFE argument I think */ 1105c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1106c0d900a5SMatthew G. Knepley { 1107c0d900a5SMatthew G. Knepley DM dmCell; 1108c0d900a5SMatthew G. Knepley Vec coordinates; 1109c0d900a5SMatthew G. Knepley PetscSection coordSection, sectionCell; 1110c0d900a5SMatthew G. Knepley PetscScalar *cgeom; 1111c0d900a5SMatthew G. Knepley PetscInt cStart, cEnd, cMax, c; 1112c0d900a5SMatthew G. Knepley PetscErrorCode ierr; 1113c0d900a5SMatthew G. Knepley 1114c0d900a5SMatthew G. Knepley PetscFunctionBegin; 1115c0d900a5SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1116c0d900a5SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1117c0d900a5SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1118c0d900a5SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1119c0d900a5SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1120c0d900a5SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1121c0d900a5SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1122c0d900a5SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1123c0d900a5SMatthew G. Knepley cEnd = cMax < 0 ? cEnd : cMax; 1124c0d900a5SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1125c0d900a5SMatthew G. Knepley /* TODO This needs to be multiplied by Nq for non-affine */ 1126c0d900a5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(PetscFECellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 1127c0d900a5SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1128c0d900a5SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1129c0d900a5SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1130c0d900a5SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1131c0d900a5SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1132c0d900a5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1133c0d900a5SMatthew G. Knepley PetscFECellGeom *cg; 1134c0d900a5SMatthew G. Knepley 1135c0d900a5SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1136c0d900a5SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1137c0d900a5SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1138c0d900a5SMatthew G. Knepley if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1139c0d900a5SMatthew G. Knepley } 1140c0d900a5SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1141c0d900a5SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1142c0d900a5SMatthew G. Knepley PetscFunctionReturn(0); 1143c0d900a5SMatthew G. Knepley } 1144c0d900a5SMatthew G. Knepley 1145c0d900a5SMatthew G. Knepley #undef __FUNCT__ 1146113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGeometryFVM" 1147113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 1148113c68e6SMatthew G. Knepley { 1149113c68e6SMatthew G. Knepley DM dmFace, dmCell; 1150113c68e6SMatthew G. Knepley DMLabel ghostLabel; 1151113c68e6SMatthew G. Knepley PetscSection sectionFace, sectionCell; 1152113c68e6SMatthew G. Knepley PetscSection coordSection; 1153113c68e6SMatthew G. Knepley Vec coordinates; 1154113c68e6SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1155113c68e6SMatthew G. Knepley PetscReal minradius, gminradius; 1156113c68e6SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 1157113c68e6SMatthew G. Knepley PetscErrorCode ierr; 1158113c68e6SMatthew G. Knepley 1159113c68e6SMatthew G. Knepley PetscFunctionBegin; 1160113c68e6SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1161113c68e6SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1162113c68e6SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1163113c68e6SMatthew G. Knepley /* Make cell centroids and volumes */ 1164113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1165113c68e6SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1166113c68e6SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1167113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1168113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1169113c68e6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1170113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1171113c68e6SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(PetscFVCellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 1172113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1173113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1174113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1175113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1176113c68e6SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1177113c68e6SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 1178113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1179113c68e6SMatthew G. Knepley 1180113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1181113c68e6SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1182113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 1183113c68e6SMatthew G. Knepley } 1184113c68e6SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 1185113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 1186113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 1187113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1188113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 1189113c68e6SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, sizeof(PetscFVFaceGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 1190113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 1191113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 1192113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 1193113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 1194113c68e6SMatthew G. Knepley ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 1195113c68e6SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1196113c68e6SMatthew G. Knepley minradius = PETSC_MAX_REAL; 1197113c68e6SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 1198113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1199113c68e6SMatthew G. Knepley PetscReal area; 12009ac3fadcSMatthew G. Knepley PetscInt ghost = -1, d; 1201113c68e6SMatthew G. Knepley 12029ac3fadcSMatthew G. Knepley if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 1203113c68e6SMatthew G. Knepley if (ghost >= 0) continue; 1204113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 1205113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 1206113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 1207113c68e6SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 1208113c68e6SMatthew G. Knepley { 1209113c68e6SMatthew G. Knepley PetscFVCellGeom *cL, *cR; 1210113c68e6SMatthew G. Knepley const PetscInt *cells; 1211113c68e6SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 1212*0453c0cdSMatthew G. Knepley PetscReal l[3], r[3], v[3]; 1213113c68e6SMatthew G. Knepley 1214113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 1215113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 1216113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 1217113c68e6SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 1218113c68e6SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 1219*0453c0cdSMatthew G. Knepley ierr = DMPlexLocalizeCoordinate_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 1220*0453c0cdSMatthew G. Knepley ierr = DMPlexLocalizeCoordinate_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 1221*0453c0cdSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 1222113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 1223113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 1224113c68e6SMatthew G. Knepley } 1225113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 1226113c68e6SMatthew G. Knepley if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]); 1227113c68e6SMatthew G. Knepley if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]); 1228113c68e6SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 1229113c68e6SMatthew G. Knepley } 1230113c68e6SMatthew G. Knepley if (cells[0] < cEndInterior) { 1231113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 1232113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1233113c68e6SMatthew G. Knepley } 1234113c68e6SMatthew G. Knepley if (cells[1] < cEndInterior) { 1235113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 1236113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 1237113c68e6SMatthew G. Knepley } 1238113c68e6SMatthew G. Knepley } 1239113c68e6SMatthew G. Knepley } 1240113c68e6SMatthew G. Knepley ierr = MPI_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1241113c68e6SMatthew G. Knepley ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 1242113c68e6SMatthew G. Knepley /* Compute centroids of ghost cells */ 1243113c68e6SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 1244113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 1245113c68e6SMatthew G. Knepley const PetscInt *cone, *support; 1246113c68e6SMatthew G. Knepley PetscInt coneSize, supportSize, s; 1247113c68e6SMatthew G. Knepley 1248113c68e6SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 1249113c68e6SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 1250113c68e6SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 1251113c68e6SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 1252113c68e6SMatthew G. Knepley if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize); 1253113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 1254113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 1255113c68e6SMatthew G. Knepley for (s = 0; s < 2; ++s) { 1256113c68e6SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 1257113c68e6SMatthew G. Knepley if (support[s] == c) { 1258113c68e6SMatthew G. Knepley const PetscFVCellGeom *ci; 1259113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 1260113c68e6SMatthew G. Knepley PetscReal c2f[3], a; 1261113c68e6SMatthew G. Knepley 1262113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 1263113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1264113c68e6SMatthew G. Knepley a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 1265113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 1266113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 1267113c68e6SMatthew G. Knepley cg->volume = ci->volume; 1268113c68e6SMatthew G. Knepley } 1269113c68e6SMatthew G. Knepley } 1270113c68e6SMatthew G. Knepley } 1271113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 1272113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1273113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1274113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 1275113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1276113c68e6SMatthew G. Knepley } 1277113c68e6SMatthew G. Knepley 1278113c68e6SMatthew G. Knepley #undef __FUNCT__ 1279113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexGetMinRadius" 1280113c68e6SMatthew G. Knepley /*@C 1281113c68e6SMatthew G. Knepley DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 1282113c68e6SMatthew G. Knepley 1283113c68e6SMatthew G. Knepley Not collective 1284113c68e6SMatthew G. Knepley 1285113c68e6SMatthew G. Knepley Input Argument: 1286113c68e6SMatthew G. Knepley . dm - the DM 1287113c68e6SMatthew G. Knepley 1288113c68e6SMatthew G. Knepley Output Argument: 1289113c68e6SMatthew G. Knepley . minradius - the minium cell radius 1290113c68e6SMatthew G. Knepley 1291113c68e6SMatthew G. Knepley Level: developer 1292113c68e6SMatthew G. Knepley 1293113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates() 1294113c68e6SMatthew G. Knepley @*/ 1295113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 1296113c68e6SMatthew G. Knepley { 1297113c68e6SMatthew G. Knepley PetscFunctionBegin; 1298113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1299113c68e6SMatthew G. Knepley PetscValidPointer(minradius,2); 1300113c68e6SMatthew G. Knepley *minradius = ((DM_Plex*) dm->data)->minradius; 1301113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1302113c68e6SMatthew G. Knepley } 1303113c68e6SMatthew G. Knepley 1304113c68e6SMatthew G. Knepley #undef __FUNCT__ 1305113c68e6SMatthew G. Knepley #define __FUNCT__ "DMPlexSetMinRadius" 1306113c68e6SMatthew G. Knepley /*@C 1307113c68e6SMatthew G. Knepley DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 1308113c68e6SMatthew G. Knepley 1309113c68e6SMatthew G. Knepley Logically collective 1310113c68e6SMatthew G. Knepley 1311113c68e6SMatthew G. Knepley Input Arguments: 1312113c68e6SMatthew G. Knepley + dm - the DM 1313113c68e6SMatthew G. Knepley - minradius - the minium cell radius 1314113c68e6SMatthew G. Knepley 1315113c68e6SMatthew G. Knepley Level: developer 1316113c68e6SMatthew G. Knepley 1317113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates() 1318113c68e6SMatthew G. Knepley @*/ 1319113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 1320113c68e6SMatthew G. Knepley { 1321113c68e6SMatthew G. Knepley PetscFunctionBegin; 1322113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1323113c68e6SMatthew G. Knepley ((DM_Plex*) dm->data)->minradius = minradius; 1324113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 1325113c68e6SMatthew G. Knepley } 1326856ac710SMatthew G. Knepley 1327856ac710SMatthew G. Knepley #undef __FUNCT__ 1328856ac710SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction_Internal" 1329856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 1330856ac710SMatthew G. Knepley { 1331856ac710SMatthew G. Knepley DMLabel ghostLabel; 1332856ac710SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 1333856ac710SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 1334856ac710SMatthew G. Knepley PetscErrorCode ierr; 1335856ac710SMatthew G. Knepley 1336856ac710SMatthew G. Knepley PetscFunctionBegin; 1337856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1338856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1339856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1340856ac710SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 1341856ac710SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 1342856ac710SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 1343856ac710SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 1344856ac710SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 1345856ac710SMatthew G. Knepley const PetscInt *faces; 1346856ac710SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 1347856ac710SMatthew G. Knepley const PetscFVCellGeom *cg; 1348856ac710SMatthew G. Knepley PetscBool boundary; 1349856ac710SMatthew G. Knepley PetscInt ghost; 1350856ac710SMatthew G. Knepley 1351856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1352856ac710SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 1353856ac710SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 1354856ac710SMatthew G. Knepley if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 1355856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1356856ac710SMatthew G. Knepley const PetscFVCellGeom *cg1; 1357856ac710SMatthew G. Knepley PetscFVFaceGeom *fg; 1358856ac710SMatthew G. Knepley const PetscInt *fcells; 1359856ac710SMatthew G. Knepley PetscInt ncell, side; 1360856ac710SMatthew G. Knepley 1361856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1362856ac710SMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1363856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1364856ac710SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 1365856ac710SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 1366856ac710SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 1367856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 1368856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 1369856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 1370856ac710SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 1371856ac710SMatthew G. Knepley } 1372856ac710SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 1373856ac710SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 1374856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1375856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1376856ac710SMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1377856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 1378856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 1379856ac710SMatthew G. Knepley ++usedFaces; 1380856ac710SMatthew G. Knepley } 1381856ac710SMatthew G. Knepley } 1382856ac710SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 1383856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1384856ac710SMatthew G. Knepley } 1385856ac710SMatthew G. Knepley 1386856ac710SMatthew G. Knepley #undef __FUNCT__ 1387856ac710SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeGradientFVM" 1388856ac710SMatthew G. Knepley /*@ 1389856ac710SMatthew G. Knepley DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 1390856ac710SMatthew G. Knepley 1391856ac710SMatthew G. Knepley Collective on DM 1392856ac710SMatthew G. Knepley 1393856ac710SMatthew G. Knepley Input Arguments: 1394856ac710SMatthew G. Knepley + dm - The DM 1395856ac710SMatthew G. Knepley . fvm - The PetscFV 1396856ac710SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexGetFaceGeometryFVM() 1397856ac710SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexGetCellGeometryFVM() 1398856ac710SMatthew G. Knepley 1399856ac710SMatthew G. Knepley Output Parameters: 1400856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted 1401856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data 1402856ac710SMatthew G. Knepley 1403856ac710SMatthew G. Knepley Level: developer 1404856ac710SMatthew G. Knepley 1405856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 1406856ac710SMatthew G. Knepley @*/ 1407856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 1408856ac710SMatthew G. Knepley { 1409856ac710SMatthew G. Knepley DM dmFace, dmCell; 1410856ac710SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 1411856ac710SMatthew G. Knepley PetscSection sectionGrad; 1412856ac710SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 1413856ac710SMatthew G. Knepley PetscErrorCode ierr; 1414856ac710SMatthew G. Knepley 1415856ac710SMatthew G. Knepley PetscFunctionBegin; 1416856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1417856ac710SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 1418856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1419856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 1420856ac710SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 1421856ac710SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 1422856ac710SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 1423856ac710SMatthew G. Knepley ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1424856ac710SMatthew G. Knepley ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1425856ac710SMatthew G. Knepley ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 1426856ac710SMatthew G. Knepley ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 1427856ac710SMatthew G. Knepley ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 1428856ac710SMatthew G. Knepley /* Create storage for gradients */ 1429856ac710SMatthew G. Knepley ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 1430856ac710SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 1431856ac710SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 1432856ac710SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 1433856ac710SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 1434856ac710SMatthew G. Knepley ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 1435856ac710SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 1436856ac710SMatthew G. Knepley PetscFunctionReturn(0); 1437856ac710SMatthew G. Knepley } 1438