1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 29d150b73SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 39d150b73SToby Isaac #include <petscblaslapack.h> 4*af74b616SDave May #include <petsctime.h> 5ccd2543fSMatthew G Knepley 6fea14342SMatthew G. Knepley static PetscErrorCode DMPlexGetLineIntersection_2D_Internal(const PetscReal segmentA[], const PetscReal segmentB[], PetscReal intersection[], PetscBool *hasIntersection) 7fea14342SMatthew G. Knepley { 8fea14342SMatthew G. Knepley const PetscReal p0_x = segmentA[0*2+0]; 9fea14342SMatthew G. Knepley const PetscReal p0_y = segmentA[0*2+1]; 10fea14342SMatthew G. Knepley const PetscReal p1_x = segmentA[1*2+0]; 11fea14342SMatthew G. Knepley const PetscReal p1_y = segmentA[1*2+1]; 12fea14342SMatthew G. Knepley const PetscReal p2_x = segmentB[0*2+0]; 13fea14342SMatthew G. Knepley const PetscReal p2_y = segmentB[0*2+1]; 14fea14342SMatthew G. Knepley const PetscReal p3_x = segmentB[1*2+0]; 15fea14342SMatthew G. Knepley const PetscReal p3_y = segmentB[1*2+1]; 16fea14342SMatthew G. Knepley const PetscReal s1_x = p1_x - p0_x; 17fea14342SMatthew G. Knepley const PetscReal s1_y = p1_y - p0_y; 18fea14342SMatthew G. Knepley const PetscReal s2_x = p3_x - p2_x; 19fea14342SMatthew G. Knepley const PetscReal s2_y = p3_y - p2_y; 20fea14342SMatthew G. Knepley const PetscReal denom = (-s2_x * s1_y + s1_x * s2_y); 21fea14342SMatthew G. Knepley 22fea14342SMatthew G. Knepley PetscFunctionBegin; 23fea14342SMatthew G. Knepley *hasIntersection = PETSC_FALSE; 24fea14342SMatthew G. Knepley /* Non-parallel lines */ 25fea14342SMatthew G. Knepley if (denom != 0.0) { 26fea14342SMatthew G. Knepley const PetscReal s = (-s1_y * (p0_x - p2_x) + s1_x * (p0_y - p2_y)) / denom; 27fea14342SMatthew G. Knepley const PetscReal t = ( s2_x * (p0_y - p2_y) - s2_y * (p0_x - p2_x)) / denom; 28fea14342SMatthew G. Knepley 29fea14342SMatthew G. Knepley if (s >= 0 && s <= 1 && t >= 0 && t <= 1) { 30fea14342SMatthew G. Knepley *hasIntersection = PETSC_TRUE; 31fea14342SMatthew G. Knepley if (intersection) { 32fea14342SMatthew G. Knepley intersection[0] = p0_x + (t * s1_x); 33fea14342SMatthew G. Knepley intersection[1] = p0_y + (t * s1_y); 34fea14342SMatthew G. Knepley } 35fea14342SMatthew G. Knepley } 36fea14342SMatthew G. Knepley } 37fea14342SMatthew G. Knepley PetscFunctionReturn(0); 38fea14342SMatthew G. Knepley } 39fea14342SMatthew G. Knepley 40ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 41ccd2543fSMatthew G Knepley { 42ccd2543fSMatthew G Knepley const PetscInt embedDim = 2; 43f5ebc837SMatthew G. Knepley const PetscReal eps = PETSC_SQRT_MACHINE_EPSILON; 44ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 45ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 46ccd2543fSMatthew G Knepley PetscReal v0[2], J[4], invJ[4], detJ; 47ccd2543fSMatthew G Knepley PetscReal xi, eta; 48ccd2543fSMatthew G Knepley PetscErrorCode ierr; 49ccd2543fSMatthew G Knepley 50ccd2543fSMatthew G Knepley PetscFunctionBegin; 518e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 52ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 53ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 54ccd2543fSMatthew G Knepley 55f5ebc837SMatthew G. Knepley if ((xi >= -eps) && (eta >= -eps) && (xi + eta <= 2.0+eps)) *cell = c; 56c1496c66SMatthew G. Knepley else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 57ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 58ccd2543fSMatthew G Knepley } 59ccd2543fSMatthew G Knepley 6062a38674SMatthew G. Knepley static PetscErrorCode DMPlexClosestPoint_Simplex_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscReal cpoint[]) 6162a38674SMatthew G. Knepley { 6262a38674SMatthew G. Knepley const PetscInt embedDim = 2; 6362a38674SMatthew G. Knepley PetscReal x = PetscRealPart(point[0]); 6462a38674SMatthew G. Knepley PetscReal y = PetscRealPart(point[1]); 6562a38674SMatthew G. Knepley PetscReal v0[2], J[4], invJ[4], detJ; 6662a38674SMatthew G. Knepley PetscReal xi, eta, r; 6762a38674SMatthew G. Knepley PetscErrorCode ierr; 6862a38674SMatthew G. Knepley 6962a38674SMatthew G. Knepley PetscFunctionBegin; 7062a38674SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 7162a38674SMatthew G. Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]); 7262a38674SMatthew G. Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]); 7362a38674SMatthew G. Knepley 7462a38674SMatthew G. Knepley xi = PetscMax(xi, 0.0); 7562a38674SMatthew G. Knepley eta = PetscMax(eta, 0.0); 7662a38674SMatthew G. Knepley r = (xi + eta)/2.0; 7762a38674SMatthew G. Knepley if (xi + eta > 2.0) { 7862a38674SMatthew G. Knepley r = (xi + eta)/2.0; 7962a38674SMatthew G. Knepley xi /= r; 8062a38674SMatthew G. Knepley eta /= r; 8162a38674SMatthew G. Knepley } 8262a38674SMatthew G. Knepley cpoint[0] = J[0*embedDim+0]*xi + J[0*embedDim+1]*eta + v0[0]; 8362a38674SMatthew G. Knepley cpoint[1] = J[1*embedDim+0]*xi + J[1*embedDim+1]*eta + v0[1]; 8462a38674SMatthew G. Knepley PetscFunctionReturn(0); 8562a38674SMatthew G. Knepley } 8662a38674SMatthew G. Knepley 87ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_2D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 88ccd2543fSMatthew G Knepley { 89ccd2543fSMatthew G Knepley PetscSection coordSection; 90ccd2543fSMatthew G Knepley Vec coordsLocal; 91a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 92ccd2543fSMatthew G Knepley const PetscInt faces[8] = {0, 1, 1, 2, 2, 3, 3, 0}; 93ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 94ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 95ccd2543fSMatthew G Knepley PetscInt crossings = 0, f; 96ccd2543fSMatthew G Knepley PetscErrorCode ierr; 97ccd2543fSMatthew G Knepley 98ccd2543fSMatthew G Knepley PetscFunctionBegin; 99ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 10069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 101ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 102ccd2543fSMatthew G Knepley for (f = 0; f < 4; ++f) { 103ccd2543fSMatthew G Knepley PetscReal x_i = PetscRealPart(coords[faces[2*f+0]*2+0]); 104ccd2543fSMatthew G Knepley PetscReal y_i = PetscRealPart(coords[faces[2*f+0]*2+1]); 105ccd2543fSMatthew G Knepley PetscReal x_j = PetscRealPart(coords[faces[2*f+1]*2+0]); 106ccd2543fSMatthew G Knepley PetscReal y_j = PetscRealPart(coords[faces[2*f+1]*2+1]); 107ccd2543fSMatthew G Knepley PetscReal slope = (y_j - y_i) / (x_j - x_i); 108ccd2543fSMatthew G Knepley PetscBool cond1 = (x_i <= x) && (x < x_j) ? PETSC_TRUE : PETSC_FALSE; 109ccd2543fSMatthew G Knepley PetscBool cond2 = (x_j <= x) && (x < x_i) ? PETSC_TRUE : PETSC_FALSE; 110ccd2543fSMatthew G Knepley PetscBool above = (y < slope * (x - x_i) + y_i) ? PETSC_TRUE : PETSC_FALSE; 111ccd2543fSMatthew G Knepley if ((cond1 || cond2) && above) ++crossings; 112ccd2543fSMatthew G Knepley } 113ccd2543fSMatthew G Knepley if (crossings % 2) *cell = c; 114c1496c66SMatthew G. Knepley else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 115ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 116ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 117ccd2543fSMatthew G Knepley } 118ccd2543fSMatthew G Knepley 119ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_Simplex_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 120ccd2543fSMatthew G Knepley { 121ccd2543fSMatthew G Knepley const PetscInt embedDim = 3; 122ccd2543fSMatthew G Knepley PetscReal v0[3], J[9], invJ[9], detJ; 123ccd2543fSMatthew G Knepley PetscReal x = PetscRealPart(point[0]); 124ccd2543fSMatthew G Knepley PetscReal y = PetscRealPart(point[1]); 125ccd2543fSMatthew G Knepley PetscReal z = PetscRealPart(point[2]); 126ccd2543fSMatthew G Knepley PetscReal xi, eta, zeta; 127ccd2543fSMatthew G Knepley PetscErrorCode ierr; 128ccd2543fSMatthew G Knepley 129ccd2543fSMatthew G Knepley PetscFunctionBegin; 1308e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 131ccd2543fSMatthew G Knepley xi = invJ[0*embedDim+0]*(x - v0[0]) + invJ[0*embedDim+1]*(y - v0[1]) + invJ[0*embedDim+2]*(z - v0[2]); 132ccd2543fSMatthew G Knepley eta = invJ[1*embedDim+0]*(x - v0[0]) + invJ[1*embedDim+1]*(y - v0[1]) + invJ[1*embedDim+2]*(z - v0[2]); 133ccd2543fSMatthew G Knepley zeta = invJ[2*embedDim+0]*(x - v0[0]) + invJ[2*embedDim+1]*(y - v0[1]) + invJ[2*embedDim+2]*(z - v0[2]); 134ccd2543fSMatthew G Knepley 135ccd2543fSMatthew G Knepley if ((xi >= 0.0) && (eta >= 0.0) && (zeta >= 0.0) && (xi + eta + zeta <= 2.0)) *cell = c; 136c1496c66SMatthew G. Knepley else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 137ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 138ccd2543fSMatthew G Knepley } 139ccd2543fSMatthew G Knepley 140ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexLocatePoint_General_3D_Internal(DM dm, const PetscScalar point[], PetscInt c, PetscInt *cell) 141ccd2543fSMatthew G Knepley { 142ccd2543fSMatthew G Knepley PetscSection coordSection; 143ccd2543fSMatthew G Knepley Vec coordsLocal; 1447c1f9639SMatthew G Knepley PetscScalar *coords; 145fb150da6SMatthew G. Knepley const PetscInt faces[24] = {0, 3, 2, 1, 5, 4, 7, 6, 3, 0, 4, 5, 146fb150da6SMatthew G. Knepley 1, 2, 6, 7, 3, 5, 6, 2, 0, 1, 7, 4}; 147ccd2543fSMatthew G Knepley PetscBool found = PETSC_TRUE; 148ccd2543fSMatthew G Knepley PetscInt f; 149ccd2543fSMatthew G Knepley PetscErrorCode ierr; 150ccd2543fSMatthew G Knepley 151ccd2543fSMatthew G Knepley PetscFunctionBegin; 152ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 15369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 154ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 155ccd2543fSMatthew G Knepley for (f = 0; f < 6; ++f) { 156ccd2543fSMatthew G Knepley /* Check the point is under plane */ 157ccd2543fSMatthew G Knepley /* Get face normal */ 158ccd2543fSMatthew G Knepley PetscReal v_i[3]; 159ccd2543fSMatthew G Knepley PetscReal v_j[3]; 160ccd2543fSMatthew G Knepley PetscReal normal[3]; 161ccd2543fSMatthew G Knepley PetscReal pp[3]; 162ccd2543fSMatthew G Knepley PetscReal dot; 163ccd2543fSMatthew G Knepley 164ccd2543fSMatthew G Knepley v_i[0] = PetscRealPart(coords[faces[f*4+3]*3+0]-coords[faces[f*4+0]*3+0]); 165ccd2543fSMatthew G Knepley v_i[1] = PetscRealPart(coords[faces[f*4+3]*3+1]-coords[faces[f*4+0]*3+1]); 166ccd2543fSMatthew G Knepley v_i[2] = PetscRealPart(coords[faces[f*4+3]*3+2]-coords[faces[f*4+0]*3+2]); 167ccd2543fSMatthew G Knepley v_j[0] = PetscRealPart(coords[faces[f*4+1]*3+0]-coords[faces[f*4+0]*3+0]); 168ccd2543fSMatthew G Knepley v_j[1] = PetscRealPart(coords[faces[f*4+1]*3+1]-coords[faces[f*4+0]*3+1]); 169ccd2543fSMatthew G Knepley v_j[2] = PetscRealPart(coords[faces[f*4+1]*3+2]-coords[faces[f*4+0]*3+2]); 170ccd2543fSMatthew G Knepley normal[0] = v_i[1]*v_j[2] - v_i[2]*v_j[1]; 171ccd2543fSMatthew G Knepley normal[1] = v_i[2]*v_j[0] - v_i[0]*v_j[2]; 172ccd2543fSMatthew G Knepley normal[2] = v_i[0]*v_j[1] - v_i[1]*v_j[0]; 173ccd2543fSMatthew G Knepley pp[0] = PetscRealPart(coords[faces[f*4+0]*3+0] - point[0]); 174ccd2543fSMatthew G Knepley pp[1] = PetscRealPart(coords[faces[f*4+0]*3+1] - point[1]); 175ccd2543fSMatthew G Knepley pp[2] = PetscRealPart(coords[faces[f*4+0]*3+2] - point[2]); 176ccd2543fSMatthew G Knepley dot = normal[0]*pp[0] + normal[1]*pp[1] + normal[2]*pp[2]; 177ccd2543fSMatthew G Knepley 178ccd2543fSMatthew G Knepley /* Check that projected point is in face (2D location problem) */ 179ccd2543fSMatthew G Knepley if (dot < 0.0) { 180ccd2543fSMatthew G Knepley found = PETSC_FALSE; 181ccd2543fSMatthew G Knepley break; 182ccd2543fSMatthew G Knepley } 183ccd2543fSMatthew G Knepley } 184ccd2543fSMatthew G Knepley if (found) *cell = c; 185c1496c66SMatthew G. Knepley else *cell = DMLOCATEPOINT_POINT_NOT_FOUND; 186ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &coords);CHKERRQ(ierr); 187ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 188ccd2543fSMatthew G Knepley } 189ccd2543fSMatthew G Knepley 190c4eade1cSMatthew G. Knepley static PetscErrorCode PetscGridHashInitialize_Internal(PetscGridHash box, PetscInt dim, const PetscScalar point[]) 191c4eade1cSMatthew G. Knepley { 192c4eade1cSMatthew G. Knepley PetscInt d; 193c4eade1cSMatthew G. Knepley 194c4eade1cSMatthew G. Knepley PetscFunctionBegin; 195c4eade1cSMatthew G. Knepley box->dim = dim; 196c4eade1cSMatthew G. Knepley for (d = 0; d < dim; ++d) box->lower[d] = box->upper[d] = PetscRealPart(point[d]); 197c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 198c4eade1cSMatthew G. Knepley } 199c4eade1cSMatthew G. Knepley 200c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashCreate(MPI_Comm comm, PetscInt dim, const PetscScalar point[], PetscGridHash *box) 201c4eade1cSMatthew G. Knepley { 202c4eade1cSMatthew G. Knepley PetscErrorCode ierr; 203c4eade1cSMatthew G. Knepley 204c4eade1cSMatthew G. Knepley PetscFunctionBegin; 205c4eade1cSMatthew G. Knepley ierr = PetscMalloc1(1, box);CHKERRQ(ierr); 206c4eade1cSMatthew G. Knepley ierr = PetscGridHashInitialize_Internal(*box, dim, point);CHKERRQ(ierr); 207c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 208c4eade1cSMatthew G. Knepley } 209c4eade1cSMatthew G. Knepley 210c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashEnlarge(PetscGridHash box, const PetscScalar point[]) 211c4eade1cSMatthew G. Knepley { 212c4eade1cSMatthew G. Knepley PetscInt d; 213c4eade1cSMatthew G. Knepley 214c4eade1cSMatthew G. Knepley PetscFunctionBegin; 215c4eade1cSMatthew G. Knepley for (d = 0; d < box->dim; ++d) { 216c4eade1cSMatthew G. Knepley box->lower[d] = PetscMin(box->lower[d], PetscRealPart(point[d])); 217c4eade1cSMatthew G. Knepley box->upper[d] = PetscMax(box->upper[d], PetscRealPart(point[d])); 218c4eade1cSMatthew G. Knepley } 219c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 220c4eade1cSMatthew G. Knepley } 221c4eade1cSMatthew G. Knepley 22262a38674SMatthew G. Knepley /* 22362a38674SMatthew G. Knepley PetscGridHashSetGrid - Divide the grid into boxes 22462a38674SMatthew G. Knepley 22562a38674SMatthew G. Knepley Not collective 22662a38674SMatthew G. Knepley 22762a38674SMatthew G. Knepley Input Parameters: 22862a38674SMatthew G. Knepley + box - The grid hash object 22962a38674SMatthew G. Knepley . n - The number of boxes in each dimension, or PETSC_DETERMINE 23062a38674SMatthew G. Knepley - h - The box size in each dimension, only used if n[d] == PETSC_DETERMINE 23162a38674SMatthew G. Knepley 23262a38674SMatthew G. Knepley Level: developer 23362a38674SMatthew G. Knepley 23462a38674SMatthew G. Knepley .seealso: PetscGridHashCreate() 23562a38674SMatthew G. Knepley */ 236c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashSetGrid(PetscGridHash box, const PetscInt n[], const PetscReal h[]) 237c4eade1cSMatthew G. Knepley { 238c4eade1cSMatthew G. Knepley PetscInt d; 239c4eade1cSMatthew G. Knepley 240c4eade1cSMatthew G. Knepley PetscFunctionBegin; 241c4eade1cSMatthew G. Knepley for (d = 0; d < box->dim; ++d) { 242c4eade1cSMatthew G. Knepley box->extent[d] = box->upper[d] - box->lower[d]; 243c4eade1cSMatthew G. Knepley if (n[d] == PETSC_DETERMINE) { 244c4eade1cSMatthew G. Knepley box->h[d] = h[d]; 245c4eade1cSMatthew G. Knepley box->n[d] = PetscCeilReal(box->extent[d]/h[d]); 246c4eade1cSMatthew G. Knepley } else { 247c4eade1cSMatthew G. Knepley box->n[d] = n[d]; 248c4eade1cSMatthew G. Knepley box->h[d] = box->extent[d]/n[d]; 249c4eade1cSMatthew G. Knepley } 250c4eade1cSMatthew G. Knepley } 251c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 252c4eade1cSMatthew G. Knepley } 253c4eade1cSMatthew G. Knepley 25462a38674SMatthew G. Knepley /* 25562a38674SMatthew G. Knepley PetscGridHashGetEnclosingBox - Find the grid boxes containing each input point 25662a38674SMatthew G. Knepley 25762a38674SMatthew G. Knepley Not collective 25862a38674SMatthew G. Knepley 25962a38674SMatthew G. Knepley Input Parameters: 26062a38674SMatthew G. Knepley + box - The grid hash object 26162a38674SMatthew G. Knepley . numPoints - The number of input points 26262a38674SMatthew G. Knepley - points - The input point coordinates 26362a38674SMatthew G. Knepley 26462a38674SMatthew G. Knepley Output Parameters: 26562a38674SMatthew G. Knepley + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 26662a38674SMatthew G. Knepley - boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL 26762a38674SMatthew G. Knepley 26862a38674SMatthew G. Knepley Level: developer 26962a38674SMatthew G. Knepley 27062a38674SMatthew G. Knepley .seealso: PetscGridHashCreate() 27162a38674SMatthew G. Knepley */ 2721c6dfc3eSMatthew G. Knepley PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[]) 273c4eade1cSMatthew G. Knepley { 274c4eade1cSMatthew G. Knepley const PetscReal *lower = box->lower; 275c4eade1cSMatthew G. Knepley const PetscReal *upper = box->upper; 276c4eade1cSMatthew G. Knepley const PetscReal *h = box->h; 277c4eade1cSMatthew G. Knepley const PetscInt *n = box->n; 278c4eade1cSMatthew G. Knepley const PetscInt dim = box->dim; 279c4eade1cSMatthew G. Knepley PetscInt d, p; 280c4eade1cSMatthew G. Knepley 281c4eade1cSMatthew G. Knepley PetscFunctionBegin; 282c4eade1cSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 283c4eade1cSMatthew G. Knepley for (d = 0; d < dim; ++d) { 2841c6dfc3eSMatthew G. Knepley PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 285c4eade1cSMatthew G. Knepley 2861c6dfc3eSMatthew G. Knepley if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 287c4eade1cSMatthew G. Knepley if (dbox < 0 || dbox >= n[d]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input point %d (%g, %g, %g) is outside of our bounding box", 2881c6dfc3eSMatthew G. Knepley p, PetscRealPart(points[p*dim+0]), dim > 1 ? PetscRealPart(points[p*dim+1]) : 0.0, dim > 2 ? PetscRealPart(points[p*dim+2]) : 0.0); 289c4eade1cSMatthew G. Knepley dboxes[p*dim+d] = dbox; 290c4eade1cSMatthew G. Knepley } 291c4eade1cSMatthew G. Knepley if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 292c4eade1cSMatthew G. Knepley } 293c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 294c4eade1cSMatthew G. Knepley } 295c4eade1cSMatthew G. Knepley 296*af74b616SDave May /* 297*af74b616SDave May PetscGridHashGetEnclosingBoxQuery - Find the grid boxes containing each input point 298*af74b616SDave May 299*af74b616SDave May Not collective 300*af74b616SDave May 301*af74b616SDave May Input Parameters: 302*af74b616SDave May + box - The grid hash object 303*af74b616SDave May . numPoints - The number of input points 304*af74b616SDave May - points - The input point coordinates 305*af74b616SDave May 306*af74b616SDave May Output Parameters: 307*af74b616SDave May + dboxes - An array of numPoints*dim integers expressing the enclosing box as (i_0, i_1, ..., i_dim) 308*af74b616SDave May . boxes - An array of numPoints integers expressing the enclosing box as single number, or NULL 309*af74b616SDave May - found - Flag indicating if point was located within a box 310*af74b616SDave May 311*af74b616SDave May Level: developer 312*af74b616SDave May 313*af74b616SDave May .seealso: PetscGridHashGetEnclosingBox() 314*af74b616SDave May */ 315*af74b616SDave May PetscErrorCode PetscGridHashGetEnclosingBoxQuery(PetscGridHash box, PetscInt numPoints, const PetscScalar points[], PetscInt dboxes[], PetscInt boxes[],PetscBool *found) 316*af74b616SDave May { 317*af74b616SDave May const PetscReal *lower = box->lower; 318*af74b616SDave May const PetscReal *upper = box->upper; 319*af74b616SDave May const PetscReal *h = box->h; 320*af74b616SDave May const PetscInt *n = box->n; 321*af74b616SDave May const PetscInt dim = box->dim; 322*af74b616SDave May PetscInt d, p; 323*af74b616SDave May 324*af74b616SDave May PetscFunctionBegin; 325*af74b616SDave May *found = PETSC_FALSE; 326*af74b616SDave May for (p = 0; p < numPoints; ++p) { 327*af74b616SDave May for (d = 0; d < dim; ++d) { 328*af74b616SDave May PetscInt dbox = PetscFloorReal((PetscRealPart(points[p*dim+d]) - lower[d])/h[d]); 329*af74b616SDave May 330*af74b616SDave May if (dbox == n[d] && PetscAbsReal(PetscRealPart(points[p*dim+d]) - upper[d]) < 1.0e-9) dbox = n[d]-1; 331*af74b616SDave May if (dbox < 0 || dbox >= n[d]) { 332*af74b616SDave May PetscFunctionReturn(0); 333*af74b616SDave May } 334*af74b616SDave May dboxes[p*dim+d] = dbox; 335*af74b616SDave May } 336*af74b616SDave May if (boxes) for (d = 1, boxes[p] = dboxes[p*dim]; d < dim; ++d) boxes[p] += dboxes[p*dim+d]*n[d-1]; 337*af74b616SDave May } 338*af74b616SDave May *found = PETSC_TRUE; 339*af74b616SDave May PetscFunctionReturn(0); 340*af74b616SDave May } 341*af74b616SDave May 342c4eade1cSMatthew G. Knepley PetscErrorCode PetscGridHashDestroy(PetscGridHash *box) 343c4eade1cSMatthew G. Knepley { 344c4eade1cSMatthew G. Knepley PetscErrorCode ierr; 345c4eade1cSMatthew G. Knepley 346c4eade1cSMatthew G. Knepley PetscFunctionBegin; 347c4eade1cSMatthew G. Knepley if (*box) { 348c4eade1cSMatthew G. Knepley ierr = PetscSectionDestroy(&(*box)->cellSection);CHKERRQ(ierr); 349c4eade1cSMatthew G. Knepley ierr = ISDestroy(&(*box)->cells);CHKERRQ(ierr); 350c4eade1cSMatthew G. Knepley ierr = DMLabelDestroy(&(*box)->cellsSparse);CHKERRQ(ierr); 351c4eade1cSMatthew G. Knepley } 352c4eade1cSMatthew G. Knepley ierr = PetscFree(*box);CHKERRQ(ierr); 353c4eade1cSMatthew G. Knepley PetscFunctionReturn(0); 354c4eade1cSMatthew G. Knepley } 355c4eade1cSMatthew G. Knepley 356cafe43deSMatthew G. Knepley PetscErrorCode DMPlexLocatePoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cellStart, PetscInt *cell) 357cafe43deSMatthew G. Knepley { 358cafe43deSMatthew G. Knepley PetscInt coneSize; 359cafe43deSMatthew G. Knepley PetscErrorCode ierr; 360cafe43deSMatthew G. Knepley 361cafe43deSMatthew G. Knepley PetscFunctionBegin; 362cafe43deSMatthew G. Knepley switch (dim) { 363cafe43deSMatthew G. Knepley case 2: 364cafe43deSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 365cafe43deSMatthew G. Knepley switch (coneSize) { 366cafe43deSMatthew G. Knepley case 3: 367cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Simplex_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 368cafe43deSMatthew G. Knepley break; 369cafe43deSMatthew G. Knepley case 4: 370cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_General_2D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 371cafe43deSMatthew G. Knepley break; 372cafe43deSMatthew G. Knepley default: 373cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 374cafe43deSMatthew G. Knepley } 375cafe43deSMatthew G. Knepley break; 376cafe43deSMatthew G. Knepley case 3: 377cafe43deSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cellStart, &coneSize);CHKERRQ(ierr); 378cafe43deSMatthew G. Knepley switch (coneSize) { 379cafe43deSMatthew G. Knepley case 4: 380cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Simplex_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 381cafe43deSMatthew G. Knepley break; 382cafe43deSMatthew G. Knepley case 6: 383cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_General_3D_Internal(dm, point, cellStart, cell);CHKERRQ(ierr); 384cafe43deSMatthew G. Knepley break; 385cafe43deSMatthew G. Knepley default: 386cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for cell with cone size %D", coneSize); 387cafe43deSMatthew G. Knepley } 388cafe43deSMatthew G. Knepley break; 389cafe43deSMatthew G. Knepley default: 390cafe43deSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No point location for mesh dimension %D", dim); 391cafe43deSMatthew G. Knepley } 392cafe43deSMatthew G. Knepley PetscFunctionReturn(0); 393cafe43deSMatthew G. Knepley } 394cafe43deSMatthew G. Knepley 39562a38674SMatthew G. Knepley /* 39662a38674SMatthew G. Knepley DMPlexClosestPoint_Internal - Returns the closest point in the cell to the given point 39762a38674SMatthew G. Knepley */ 39862a38674SMatthew G. Knepley PetscErrorCode DMPlexClosestPoint_Internal(DM dm, PetscInt dim, const PetscScalar point[], PetscInt cell, PetscReal cpoint[]) 39962a38674SMatthew G. Knepley { 40062a38674SMatthew G. Knepley PetscInt coneSize; 40162a38674SMatthew G. Knepley PetscErrorCode ierr; 40262a38674SMatthew G. Knepley 40362a38674SMatthew G. Knepley PetscFunctionBegin; 40462a38674SMatthew G. Knepley switch (dim) { 40562a38674SMatthew G. Knepley case 2: 40662a38674SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 40762a38674SMatthew G. Knepley switch (coneSize) { 40862a38674SMatthew G. Knepley case 3: 40962a38674SMatthew G. Knepley ierr = DMPlexClosestPoint_Simplex_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 41062a38674SMatthew G. Knepley break; 41162a38674SMatthew G. Knepley #if 0 41262a38674SMatthew G. Knepley case 4: 41362a38674SMatthew G. Knepley ierr = DMPlexClosestPoint_General_2D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 41462a38674SMatthew G. Knepley break; 41562a38674SMatthew G. Knepley #endif 41662a38674SMatthew G. Knepley default: 41762a38674SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 41862a38674SMatthew G. Knepley } 41962a38674SMatthew G. Knepley break; 42062a38674SMatthew G. Knepley #if 0 42162a38674SMatthew G. Knepley case 3: 42262a38674SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 42362a38674SMatthew G. Knepley switch (coneSize) { 42462a38674SMatthew G. Knepley case 4: 42562a38674SMatthew G. Knepley ierr = DMPlexClosestPoint_Simplex_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 42662a38674SMatthew G. Knepley break; 42762a38674SMatthew G. Knepley case 6: 42862a38674SMatthew G. Knepley ierr = DMPlexClosestPoint_General_3D_Internal(dm, point, cell, cpoint);CHKERRQ(ierr); 42962a38674SMatthew G. Knepley break; 43062a38674SMatthew G. Knepley default: 43162a38674SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for cell with cone size %D", coneSize); 43262a38674SMatthew G. Knepley } 43362a38674SMatthew G. Knepley break; 43462a38674SMatthew G. Knepley #endif 43562a38674SMatthew G. Knepley default: 43662a38674SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No closest point location for mesh dimension %D", dim); 43762a38674SMatthew G. Knepley } 43862a38674SMatthew G. Knepley PetscFunctionReturn(0); 43962a38674SMatthew G. Knepley } 44062a38674SMatthew G. Knepley 44162a38674SMatthew G. Knepley /* 44262a38674SMatthew G. Knepley DMPlexComputeGridHash_Internal - Create a grid hash structure covering the Plex 44362a38674SMatthew G. Knepley 44462a38674SMatthew G. Knepley Collective on DM 44562a38674SMatthew G. Knepley 44662a38674SMatthew G. Knepley Input Parameter: 44762a38674SMatthew G. Knepley . dm - The Plex 44862a38674SMatthew G. Knepley 44962a38674SMatthew G. Knepley Output Parameter: 45062a38674SMatthew G. Knepley . localBox - The grid hash object 45162a38674SMatthew G. Knepley 45262a38674SMatthew G. Knepley Level: developer 45362a38674SMatthew G. Knepley 45462a38674SMatthew G. Knepley .seealso: PetscGridHashCreate(), PetscGridHashGetEnclosingBox() 45562a38674SMatthew G. Knepley */ 456cafe43deSMatthew G. Knepley PetscErrorCode DMPlexComputeGridHash_Internal(DM dm, PetscGridHash *localBox) 457cafe43deSMatthew G. Knepley { 458cafe43deSMatthew G. Knepley MPI_Comm comm; 459cafe43deSMatthew G. Knepley PetscGridHash lbox; 460cafe43deSMatthew G. Knepley Vec coordinates; 461cafe43deSMatthew G. Knepley PetscSection coordSection; 462cafe43deSMatthew G. Knepley Vec coordsLocal; 463cafe43deSMatthew G. Knepley const PetscScalar *coords; 464722d0f5cSMatthew G. Knepley PetscInt *dboxes, *boxes; 465cafe43deSMatthew G. Knepley PetscInt n[3] = {10, 10, 10}; 4661d0c6c94SMatthew G. Knepley PetscInt dim, N, cStart, cEnd, cMax, c, i; 467cafe43deSMatthew G. Knepley PetscErrorCode ierr; 468cafe43deSMatthew G. Knepley 469cafe43deSMatthew G. Knepley PetscFunctionBegin; 470cafe43deSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 471cafe43deSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 472cafe43deSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 4735b3353d8SMatthew G. Knepley if (dim != 2) SETERRQ(comm, PETSC_ERR_SUP, "I have only coded this for 2D"); 474cafe43deSMatthew G. Knepley ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 475cafe43deSMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 476cafe43deSMatthew G. Knepley ierr = PetscGridHashCreate(comm, dim, coords, &lbox);CHKERRQ(ierr); 477cafe43deSMatthew G. Knepley for (i = 0; i < N; i += dim) {ierr = PetscGridHashEnlarge(lbox, &coords[i]);CHKERRQ(ierr);} 478cafe43deSMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 479cafe43deSMatthew G. Knepley ierr = PetscGridHashSetGrid(lbox, n, NULL);CHKERRQ(ierr); 480cafe43deSMatthew G. Knepley #if 0 481cafe43deSMatthew G. Knepley /* Could define a custom reduction to merge these */ 482b2566f29SBarry Smith ierr = MPIU_Allreduce(lbox->lower, gbox->lower, 3, MPIU_REAL, MPI_MIN, comm);CHKERRQ(ierr); 483b2566f29SBarry Smith ierr = MPIU_Allreduce(lbox->upper, gbox->upper, 3, MPIU_REAL, MPI_MAX, comm);CHKERRQ(ierr); 484cafe43deSMatthew G. Knepley #endif 485cafe43deSMatthew G. Knepley /* Is there a reason to snap the local bounding box to a division of the global box? */ 486cafe43deSMatthew G. Knepley /* Should we compute all overlaps of local boxes? We could do this with a rendevouz scheme partitioning the global box */ 487cafe43deSMatthew G. Knepley /* Create label */ 488cafe43deSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 4891d0c6c94SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 4901d0c6c94SMatthew G. Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 491cafe43deSMatthew G. Knepley ierr = DMLabelCreate("cells", &lbox->cellsSparse);CHKERRQ(ierr); 492cafe43deSMatthew G. Knepley ierr = DMLabelCreateIndex(lbox->cellsSparse, cStart, cEnd);CHKERRQ(ierr); 493722d0f5cSMatthew G. Knepley /* Compute boxes which overlap each cell: http://stackoverflow.com/questions/13790208/triangle-square-intersection-test-in-2d */ 494cafe43deSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 495cafe43deSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 49638353de4SMatthew G. Knepley ierr = PetscCalloc2(16 * dim, &dboxes, 16, &boxes);CHKERRQ(ierr); 497cafe43deSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 498cafe43deSMatthew G. Knepley const PetscReal *h = lbox->h; 499cafe43deSMatthew G. Knepley PetscScalar *ccoords = NULL; 50038353de4SMatthew G. Knepley PetscInt csize = 0; 501cafe43deSMatthew G. Knepley PetscScalar point[3]; 502cafe43deSMatthew G. Knepley PetscInt dlim[6], d, e, i, j, k; 503cafe43deSMatthew G. Knepley 504cafe43deSMatthew G. Knepley /* Find boxes enclosing each vertex */ 50538353de4SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &csize, &ccoords);CHKERRQ(ierr); 50638353de4SMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(lbox, csize/dim, ccoords, dboxes, boxes);CHKERRQ(ierr); 507722d0f5cSMatthew G. Knepley /* Mark cells containing the vertices */ 50838353de4SMatthew G. Knepley for (e = 0; e < csize/dim; ++e) {ierr = DMLabelSetValue(lbox->cellsSparse, c, boxes[e]);CHKERRQ(ierr);} 509cafe43deSMatthew G. Knepley /* Get grid of boxes containing these */ 510cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) {dlim[d*2+0] = dlim[d*2+1] = dboxes[d];} 5112291669eSMatthew G. Knepley for (d = dim; d < 3; ++d) {dlim[d*2+0] = dlim[d*2+1] = 0;} 512cafe43deSMatthew G. Knepley for (e = 1; e < dim+1; ++e) { 513cafe43deSMatthew G. Knepley for (d = 0; d < dim; ++d) { 514cafe43deSMatthew G. Knepley dlim[d*2+0] = PetscMin(dlim[d*2+0], dboxes[e*dim+d]); 515cafe43deSMatthew G. Knepley dlim[d*2+1] = PetscMax(dlim[d*2+1], dboxes[e*dim+d]); 516cafe43deSMatthew G. Knepley } 517cafe43deSMatthew G. Knepley } 518fea14342SMatthew G. Knepley /* Check for intersection of box with cell */ 519cafe43deSMatthew G. Knepley for (k = dlim[2*2+0], point[2] = lbox->lower[2] + k*h[2]; k <= dlim[2*2+1]; ++k, point[2] += h[2]) { 520cafe43deSMatthew G. Knepley for (j = dlim[1*2+0], point[1] = lbox->lower[1] + j*h[1]; j <= dlim[1*2+1]; ++j, point[1] += h[1]) { 521cafe43deSMatthew G. Knepley for (i = dlim[0*2+0], point[0] = lbox->lower[0] + i*h[0]; i <= dlim[0*2+1]; ++i, point[0] += h[0]) { 522cafe43deSMatthew G. Knepley const PetscInt box = (k*lbox->n[1] + j)*lbox->n[0] + i; 523cafe43deSMatthew G. Knepley PetscScalar cpoint[3]; 524fea14342SMatthew G. Knepley PetscInt cell, edge, ii, jj, kk; 525cafe43deSMatthew G. Knepley 526fea14342SMatthew G. Knepley /* Check whether cell contains any vertex of these subboxes TODO vectorize this */ 527cafe43deSMatthew G. Knepley for (kk = 0, cpoint[2] = point[2]; kk < (dim > 2 ? 2 : 1); ++kk, cpoint[2] += h[2]) { 528cafe43deSMatthew G. Knepley for (jj = 0, cpoint[1] = point[1]; jj < (dim > 1 ? 2 : 1); ++jj, cpoint[1] += h[1]) { 529cafe43deSMatthew G. Knepley for (ii = 0, cpoint[0] = point[0]; ii < 2; ++ii, cpoint[0] += h[0]) { 530cafe43deSMatthew G. Knepley 531cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, cpoint, c, &cell);CHKERRQ(ierr); 532cafe43deSMatthew G. Knepley if (cell >= 0) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); ii = jj = kk = 2;} 533cafe43deSMatthew G. Knepley } 534cafe43deSMatthew G. Knepley } 535cafe43deSMatthew G. Knepley } 536fea14342SMatthew G. Knepley /* Check whether cell edge intersects any edge of these subboxes TODO vectorize this */ 537fea14342SMatthew G. Knepley for (edge = 0; edge < dim+1; ++edge) { 538fea14342SMatthew G. Knepley PetscReal segA[6], segB[6]; 539fea14342SMatthew G. Knepley 540fea14342SMatthew G. Knepley for (d = 0; d < dim; ++d) {segA[d] = PetscRealPart(ccoords[edge*dim+d]); segA[dim+d] = PetscRealPart(ccoords[((edge+1)%(dim+1))*dim+d]);} 541fea14342SMatthew G. Knepley for (kk = 0; kk < (dim > 2 ? 2 : 1); ++kk) { 5429a128ed2SMatthew G. Knepley if (dim > 2) {segB[2] = PetscRealPart(point[2]); 5439a128ed2SMatthew G. Knepley segB[dim+2] = PetscRealPart(point[2]) + kk*h[2];} 544fea14342SMatthew G. Knepley for (jj = 0; jj < (dim > 1 ? 2 : 1); ++jj) { 5459a128ed2SMatthew G. Knepley if (dim > 1) {segB[1] = PetscRealPart(point[1]); 5469a128ed2SMatthew G. Knepley segB[dim+1] = PetscRealPart(point[1]) + jj*h[1];} 547fea14342SMatthew G. Knepley for (ii = 0; ii < 2; ++ii) { 548fea14342SMatthew G. Knepley PetscBool intersects; 549fea14342SMatthew G. Knepley 5509a128ed2SMatthew G. Knepley segB[0] = PetscRealPart(point[0]); 5519a128ed2SMatthew G. Knepley segB[dim+0] = PetscRealPart(point[0]) + ii*h[0]; 552fea14342SMatthew G. Knepley ierr = DMPlexGetLineIntersection_2D_Internal(segA, segB, NULL, &intersects);CHKERRQ(ierr); 553fea14342SMatthew G. Knepley if (intersects) {DMLabelSetValue(lbox->cellsSparse, c, box);CHKERRQ(ierr); edge = ii = jj = kk = dim+1;} 554cafe43deSMatthew G. Knepley } 555cafe43deSMatthew G. Knepley } 556cafe43deSMatthew G. Knepley } 557cafe43deSMatthew G. Knepley } 558fea14342SMatthew G. Knepley } 559fea14342SMatthew G. Knepley } 560fea14342SMatthew G. Knepley } 561fea14342SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &ccoords);CHKERRQ(ierr); 562fea14342SMatthew G. Knepley } 563722d0f5cSMatthew G. Knepley ierr = PetscFree2(dboxes, boxes);CHKERRQ(ierr); 564cafe43deSMatthew G. Knepley ierr = DMLabelConvertToSection(lbox->cellsSparse, &lbox->cellSection, &lbox->cells);CHKERRQ(ierr); 565cafe43deSMatthew G. Knepley ierr = DMLabelDestroy(&lbox->cellsSparse);CHKERRQ(ierr); 566cafe43deSMatthew G. Knepley *localBox = lbox; 567cafe43deSMatthew G. Knepley PetscFunctionReturn(0); 568cafe43deSMatthew G. Knepley } 569cafe43deSMatthew G. Knepley 57062a38674SMatthew G. Knepley PetscErrorCode DMLocatePoints_Plex(DM dm, Vec v, DMPointLocationType ltype, PetscSF cellSF) 571ccd2543fSMatthew G Knepley { 572cafe43deSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 573*af74b616SDave May PetscBool hash = mesh->useHashLocation, reuse = PETSC_FALSE; 5743a93e3b7SToby Isaac PetscInt bs, numPoints, p, numFound, *found = NULL; 5751318edbeSMatthew G. Knepley PetscInt dim, cStart, cEnd, cMax, numCells, c; 576cafe43deSMatthew G. Knepley const PetscInt *boxCells; 5773a93e3b7SToby Isaac PetscSFNode *cells; 578ccd2543fSMatthew G Knepley PetscScalar *a; 5793a93e3b7SToby Isaac PetscMPIInt result; 580*af74b616SDave May PetscLogDouble t0,t1; 581ccd2543fSMatthew G Knepley PetscErrorCode ierr; 582ccd2543fSMatthew G Knepley 583ccd2543fSMatthew G Knepley PetscFunctionBegin; 584*af74b616SDave May ierr = PetscTime(&t0);CHKERRQ(ierr); 585080342d1SMatthew G. Knepley if (ltype == DM_POINTLOCATION_NEAREST && !hash) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Nearest point location only supported with grid hashing. Use -dm_plex_hash_location to enable it."); 586cafe43deSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 587cafe43deSMatthew G. Knepley ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 5883a93e3b7SToby Isaac ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)cellSF),PETSC_COMM_SELF,&result);CHKERRQ(ierr); 5893a93e3b7SToby Isaac if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)cellSF),PETSC_ERR_SUP, "Trying parallel point location: only local point location supported"); 590cafe43deSMatthew G. Knepley if (bs != dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Block size for point vector %D must be the mesh coordinate dimension %D", bs, dim); 591ccd2543fSMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 592ccd2543fSMatthew G Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 593ccd2543fSMatthew G Knepley if (cMax >= 0) cEnd = PetscMin(cEnd, cMax); 594ccd2543fSMatthew G Knepley ierr = VecGetLocalSize(v, &numPoints);CHKERRQ(ierr); 595ccd2543fSMatthew G Knepley ierr = VecGetArray(v, &a);CHKERRQ(ierr); 596ccd2543fSMatthew G Knepley numPoints /= bs; 597*af74b616SDave May { 598*af74b616SDave May const PetscSFNode *sf_cells; 599*af74b616SDave May 600*af74b616SDave May ierr = PetscSFGetGraph(cellSF,NULL,NULL,NULL,&sf_cells);CHKERRQ(ierr); 601*af74b616SDave May if (sf_cells) { 602*af74b616SDave May ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Re-using existing StarForest node list\n");CHKERRQ(ierr); 603*af74b616SDave May cells = (PetscSFNode*)sf_cells; 604*af74b616SDave May reuse = PETSC_TRUE; 605*af74b616SDave May } else { 606*af74b616SDave May ierr = PetscInfo(dm,"[DMLocatePoints_Plex] Creating and initializing new StarForest node list\n");CHKERRQ(ierr); 607785e854fSJed Brown ierr = PetscMalloc1(numPoints, &cells);CHKERRQ(ierr); 608*af74b616SDave May /* initialize cells if created */ 609*af74b616SDave May for (p=0; p<numPoints; p++) { 610*af74b616SDave May cells[p].rank = 0; 611*af74b616SDave May cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 612*af74b616SDave May } 613*af74b616SDave May } 614*af74b616SDave May } 615953fc75cSMatthew G. Knepley if (hash) { 616ac6ec2abSMatthew G. Knepley if (!mesh->lbox) {ierr = PetscInfo(dm, "Initializing grid hashing");CHKERRQ(ierr);ierr = DMPlexComputeGridHash_Internal(dm, &mesh->lbox);CHKERRQ(ierr);} 617cafe43deSMatthew G. Knepley /* Designate the local box for each point */ 618cafe43deSMatthew G. Knepley /* Send points to correct process */ 619cafe43deSMatthew G. Knepley /* Search cells that lie in each subbox */ 620cafe43deSMatthew G. Knepley /* Should we bin points before doing search? */ 621cafe43deSMatthew G. Knepley ierr = ISGetIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr); 622953fc75cSMatthew G. Knepley } 6233a93e3b7SToby Isaac for (p = 0, numFound = 0; p < numPoints; ++p) { 624ccd2543fSMatthew G Knepley const PetscScalar *point = &a[p*bs]; 625953fc75cSMatthew G. Knepley PetscInt dbin[3], bin, cell = -1, cellOffset; 626ccd2543fSMatthew G Knepley 627*af74b616SDave May /* check initial values in cells[].index - abort early if found */ 628*af74b616SDave May if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 629*af74b616SDave May c = cells[p].index; 630e9b685f5SMatthew G. Knepley cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 631*af74b616SDave May ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 632*af74b616SDave May if (cell >= 0) { 633*af74b616SDave May cells[p].rank = 0; 634*af74b616SDave May cells[p].index = cell; 635*af74b616SDave May numFound++; 636*af74b616SDave May } 637*af74b616SDave May } 638*af74b616SDave May if (cells[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { continue; } 639*af74b616SDave May 640953fc75cSMatthew G. Knepley if (hash) { 641*af74b616SDave May PetscBool found_box; 642*af74b616SDave May 643*af74b616SDave May //ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 644*af74b616SDave May /* allow for case that point is outside box - abort early */ 645*af74b616SDave May ierr = PetscGridHashGetEnclosingBoxQuery(mesh->lbox, 1, point, dbin, &bin,&found_box);CHKERRQ(ierr); 646*af74b616SDave May if (found_box) { 647cafe43deSMatthew G. Knepley /* TODO Lay an interface over this so we can switch between Section (dense) and Label (sparse) */ 648cafe43deSMatthew G. Knepley ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 649cafe43deSMatthew G. Knepley ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 650cafe43deSMatthew G. Knepley for (c = cellOffset; c < cellOffset + numCells; ++c) { 651cafe43deSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, point, boxCells[c], &cell);CHKERRQ(ierr); 6523a93e3b7SToby Isaac if (cell >= 0) { 6533a93e3b7SToby Isaac cells[p].rank = 0; 6543a93e3b7SToby Isaac cells[p].index = cell; 6553a93e3b7SToby Isaac numFound++; 6563a93e3b7SToby Isaac break; 657ccd2543fSMatthew G Knepley } 6583a93e3b7SToby Isaac } 659*af74b616SDave May } 660953fc75cSMatthew G. Knepley } else { 661953fc75cSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 662953fc75cSMatthew G. Knepley ierr = DMPlexLocatePoint_Internal(dm, dim, point, c, &cell);CHKERRQ(ierr); 6633a93e3b7SToby Isaac if (cell >= 0) { 6643a93e3b7SToby Isaac cells[p].rank = 0; 6653a93e3b7SToby Isaac cells[p].index = cell; 6663a93e3b7SToby Isaac numFound++; 6673a93e3b7SToby Isaac break; 668953fc75cSMatthew G. Knepley } 669953fc75cSMatthew G. Knepley } 6703a93e3b7SToby Isaac } 671ccd2543fSMatthew G Knepley } 672953fc75cSMatthew G. Knepley if (hash) {ierr = ISRestoreIndices(mesh->lbox->cells, &boxCells);CHKERRQ(ierr);} 67362a38674SMatthew G. Knepley if (ltype == DM_POINTLOCATION_NEAREST && hash && numFound < numPoints) { 67462a38674SMatthew G. Knepley for (p = 0; p < numPoints; p++) { 67562a38674SMatthew G. Knepley const PetscScalar *point = &a[p*bs]; 67662a38674SMatthew G. Knepley PetscReal cpoint[3], diff[3], dist, distMax = PETSC_MAX_REAL; 677b716b415SMatthew G. Knepley PetscInt dbin[3], bin, cellOffset, d; 67862a38674SMatthew G. Knepley 679e9b685f5SMatthew G. Knepley if (cells[p].index < 0) { 68062a38674SMatthew G. Knepley ++numFound; 68162a38674SMatthew G. Knepley ierr = PetscGridHashGetEnclosingBox(mesh->lbox, 1, point, dbin, &bin);CHKERRQ(ierr); 68262a38674SMatthew G. Knepley ierr = PetscSectionGetDof(mesh->lbox->cellSection, bin, &numCells);CHKERRQ(ierr); 68362a38674SMatthew G. Knepley ierr = PetscSectionGetOffset(mesh->lbox->cellSection, bin, &cellOffset);CHKERRQ(ierr); 68462a38674SMatthew G. Knepley for (c = cellOffset; c < cellOffset + numCells; ++c) { 68562a38674SMatthew G. Knepley ierr = DMPlexClosestPoint_Internal(dm, dim, point, boxCells[c], cpoint);CHKERRQ(ierr); 686b716b415SMatthew G. Knepley for (d = 0; d < dim; ++d) diff[d] = cpoint[d] - PetscRealPart(point[d]); 68762a38674SMatthew G. Knepley dist = DMPlex_NormD_Internal(dim, diff); 68862a38674SMatthew G. Knepley if (dist < distMax) { 68962a38674SMatthew G. Knepley for (d = 0; d < dim; ++d) a[p*bs+d] = cpoint[d]; 69062a38674SMatthew G. Knepley cells[p].rank = 0; 69162a38674SMatthew G. Knepley cells[p].index = boxCells[c]; 69262a38674SMatthew G. Knepley distMax = dist; 69362a38674SMatthew G. Knepley break; 69462a38674SMatthew G. Knepley } 69562a38674SMatthew G. Knepley } 69662a38674SMatthew G. Knepley } 69762a38674SMatthew G. Knepley } 69862a38674SMatthew G. Knepley } 69962a38674SMatthew G. Knepley /* This code is only be relevant when interfaced to parallel point location */ 700cafe43deSMatthew G. Knepley /* Check for highest numbered proc that claims a point (do we care?) */ 7012d1fa6caSMatthew G. Knepley if (ltype == DM_POINTLOCATION_REMOVE && numFound < numPoints) { 7023a93e3b7SToby Isaac ierr = PetscMalloc1(numFound,&found);CHKERRQ(ierr); 7033a93e3b7SToby Isaac for (p = 0, numFound = 0; p < numPoints; p++) { 7043a93e3b7SToby Isaac if (cells[p].rank >= 0 && cells[p].index >= 0) { 7053a93e3b7SToby Isaac if (numFound < p) { 7063a93e3b7SToby Isaac cells[numFound] = cells[p]; 7073a93e3b7SToby Isaac } 7083a93e3b7SToby Isaac found[numFound++] = p; 7093a93e3b7SToby Isaac } 7103a93e3b7SToby Isaac } 7113a93e3b7SToby Isaac } 71262a38674SMatthew G. Knepley ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 713*af74b616SDave May if (!reuse) { 7143a93e3b7SToby Isaac ierr = PetscSFSetGraph(cellSF, cEnd - cStart, numFound, found, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 715*af74b616SDave May } 716*af74b616SDave May ierr = PetscTime(&t1);CHKERRQ(ierr); 717*af74b616SDave May ierr = PetscInfo3(dm,"[DMLocatePoints_Plex] npoints %D : time(rank0) %1.2e (sec): points/sec %1.4e\n",numPoints,t1-t0,(double)((double)numPoints/(t1-t0)));CHKERRQ(ierr); 718ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 719ccd2543fSMatthew G Knepley } 720ccd2543fSMatthew G Knepley 721741bfc07SMatthew G. Knepley /*@C 722741bfc07SMatthew G. Knepley DMPlexComputeProjection2Dto1D - Rewrite coordinates to be the 1D projection of the 2D coordinates 723741bfc07SMatthew G. Knepley 724741bfc07SMatthew G. Knepley Not collective 725741bfc07SMatthew G. Knepley 726741bfc07SMatthew G. Knepley Input Parameter: 727741bfc07SMatthew G. Knepley . coords - The coordinates of a segment 728741bfc07SMatthew G. Knepley 729741bfc07SMatthew G. Knepley Output Parameters: 730741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x 731741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection 732741bfc07SMatthew G. Knepley 733741bfc07SMatthew G. Knepley Level: developer 734741bfc07SMatthew G. Knepley 735741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection3Dto1D(), DMPlexComputeProjection3Dto2D() 736741bfc07SMatthew G. Knepley @*/ 737741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar coords[], PetscReal R[]) 73817fe8556SMatthew G. Knepley { 73917fe8556SMatthew G. Knepley const PetscReal x = PetscRealPart(coords[2] - coords[0]); 74017fe8556SMatthew G. Knepley const PetscReal y = PetscRealPart(coords[3] - coords[1]); 7418b49ba18SBarry Smith const PetscReal r = PetscSqrtReal(x*x + y*y), c = x/r, s = y/r; 74217fe8556SMatthew G. Knepley 74317fe8556SMatthew G. Knepley PetscFunctionBegin; 7441c99cf0cSGeoffrey Irving R[0] = c; R[1] = -s; 7451c99cf0cSGeoffrey Irving R[2] = s; R[3] = c; 74617fe8556SMatthew G. Knepley coords[0] = 0.0; 7477f07f362SMatthew G. Knepley coords[1] = r; 74817fe8556SMatthew G. Knepley PetscFunctionReturn(0); 74917fe8556SMatthew G. Knepley } 75017fe8556SMatthew G. Knepley 751741bfc07SMatthew G. Knepley /*@C 752741bfc07SMatthew G. Knepley DMPlexComputeProjection3Dto1D - Rewrite coordinates to be the 1D projection of the 3D coordinates 75328dbe442SToby Isaac 754741bfc07SMatthew G. Knepley Not collective 75528dbe442SToby Isaac 756741bfc07SMatthew G. Knepley Input Parameter: 757741bfc07SMatthew G. Knepley . coords - The coordinates of a segment 758741bfc07SMatthew G. Knepley 759741bfc07SMatthew G. Knepley Output Parameters: 760741bfc07SMatthew G. Knepley + coords - The new y-coordinate, and 0 for x and z 761741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection 762741bfc07SMatthew G. Knepley 763741bfc07SMatthew G. Knepley Note: This uses the basis completion described by Frisvad in http://www.imm.dtu.dk/~jerf/papers/abstracts/onb.html, DOI:10.1080/2165347X.2012.689606 764741bfc07SMatthew G. Knepley 765741bfc07SMatthew G. Knepley Level: developer 766741bfc07SMatthew G. Knepley 767741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto2D() 768741bfc07SMatthew G. Knepley @*/ 769741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar coords[], PetscReal R[]) 77028dbe442SToby Isaac { 77128dbe442SToby Isaac PetscReal x = PetscRealPart(coords[3] - coords[0]); 77228dbe442SToby Isaac PetscReal y = PetscRealPart(coords[4] - coords[1]); 77328dbe442SToby Isaac PetscReal z = PetscRealPart(coords[5] - coords[2]); 77428dbe442SToby Isaac PetscReal r = PetscSqrtReal(x*x + y*y + z*z); 77528dbe442SToby Isaac PetscReal rinv = 1. / r; 77628dbe442SToby Isaac PetscFunctionBegin; 77728dbe442SToby Isaac 77828dbe442SToby Isaac x *= rinv; y *= rinv; z *= rinv; 77928dbe442SToby Isaac if (x > 0.) { 78028dbe442SToby Isaac PetscReal inv1pX = 1./ (1. + x); 78128dbe442SToby Isaac 78228dbe442SToby Isaac R[0] = x; R[1] = -y; R[2] = -z; 78328dbe442SToby Isaac R[3] = y; R[4] = 1. - y*y*inv1pX; R[5] = -y*z*inv1pX; 78428dbe442SToby Isaac R[6] = z; R[7] = -y*z*inv1pX; R[8] = 1. - z*z*inv1pX; 78528dbe442SToby Isaac } 78628dbe442SToby Isaac else { 78728dbe442SToby Isaac PetscReal inv1mX = 1./ (1. - x); 78828dbe442SToby Isaac 78928dbe442SToby Isaac R[0] = x; R[1] = z; R[2] = y; 79028dbe442SToby Isaac R[3] = y; R[4] = -y*z*inv1mX; R[5] = 1. - y*y*inv1mX; 79128dbe442SToby Isaac R[6] = z; R[7] = 1. - z*z*inv1mX; R[8] = -y*z*inv1mX; 79228dbe442SToby Isaac } 79328dbe442SToby Isaac coords[0] = 0.0; 79428dbe442SToby Isaac coords[1] = r; 79528dbe442SToby Isaac PetscFunctionReturn(0); 79628dbe442SToby Isaac } 79728dbe442SToby Isaac 798741bfc07SMatthew G. Knepley /*@ 799741bfc07SMatthew G. Knepley DMPlexComputeProjection3Dto2D - Rewrite coordinates to be the 2D projection of the 3D coordinates 800741bfc07SMatthew G. Knepley 801741bfc07SMatthew G. Knepley Not collective 802741bfc07SMatthew G. Knepley 803741bfc07SMatthew G. Knepley Input Parameter: 804741bfc07SMatthew G. Knepley . coords - The coordinates of a segment 805741bfc07SMatthew G. Knepley 806741bfc07SMatthew G. Knepley Output Parameters: 807741bfc07SMatthew G. Knepley + coords - The new y- and z-coordinates, and 0 for x 808741bfc07SMatthew G. Knepley - R - The rotation which accomplishes the projection 809741bfc07SMatthew G. Knepley 810741bfc07SMatthew G. Knepley Level: developer 811741bfc07SMatthew G. Knepley 812741bfc07SMatthew G. Knepley .seealso: DMPlexComputeProjection2Dto1D(), DMPlexComputeProjection3Dto1D() 813741bfc07SMatthew G. Knepley @*/ 814741bfc07SMatthew G. Knepley PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt coordSize, PetscScalar coords[], PetscReal R[]) 815ccd2543fSMatthew G Knepley { 8161ee9d5ecSMatthew G. Knepley PetscReal x1[3], x2[3], n[3], norm; 81799dec3a6SMatthew G. Knepley PetscReal x1p[3], x2p[3], xnp[3]; 8184a217a95SMatthew G. Knepley PetscReal sqrtz, alpha; 819ccd2543fSMatthew G Knepley const PetscInt dim = 3; 82099dec3a6SMatthew G. Knepley PetscInt d, e, p; 821ccd2543fSMatthew G Knepley 822ccd2543fSMatthew G Knepley PetscFunctionBegin; 823ccd2543fSMatthew G Knepley /* 0) Calculate normal vector */ 824ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 8251ee9d5ecSMatthew G. Knepley x1[d] = PetscRealPart(coords[1*dim+d] - coords[0*dim+d]); 8261ee9d5ecSMatthew G. Knepley x2[d] = PetscRealPart(coords[2*dim+d] - coords[0*dim+d]); 827ccd2543fSMatthew G Knepley } 828ccd2543fSMatthew G Knepley n[0] = x1[1]*x2[2] - x1[2]*x2[1]; 829ccd2543fSMatthew G Knepley n[1] = x1[2]*x2[0] - x1[0]*x2[2]; 830ccd2543fSMatthew G Knepley n[2] = x1[0]*x2[1] - x1[1]*x2[0]; 8318b49ba18SBarry Smith norm = PetscSqrtReal(n[0]*n[0] + n[1]*n[1] + n[2]*n[2]); 832ccd2543fSMatthew G Knepley n[0] /= norm; 833ccd2543fSMatthew G Knepley n[1] /= norm; 834ccd2543fSMatthew G Knepley n[2] /= norm; 835ccd2543fSMatthew G Knepley /* 1) Take the normal vector and rotate until it is \hat z 836ccd2543fSMatthew G Knepley 837ccd2543fSMatthew G Knepley Let the normal vector be <nx, ny, nz> and alpha = 1/sqrt(1 - nz^2), then 838ccd2543fSMatthew G Knepley 839ccd2543fSMatthew G Knepley R = / alpha nx nz alpha ny nz -1/alpha \ 840ccd2543fSMatthew G Knepley | -alpha ny alpha nx 0 | 841ccd2543fSMatthew G Knepley \ nx ny nz / 842ccd2543fSMatthew G Knepley 843ccd2543fSMatthew G Knepley will rotate the normal vector to \hat z 844ccd2543fSMatthew G Knepley */ 8458b49ba18SBarry Smith sqrtz = PetscSqrtReal(1.0 - n[2]*n[2]); 84673868372SMatthew G. Knepley /* Check for n = z */ 84773868372SMatthew G. Knepley if (sqrtz < 1.0e-10) { 8487df32b8bSSanderA const PetscInt s = PetscSign(n[2]); 8497df32b8bSSanderA /* If nz < 0, rotate 180 degrees around x-axis */ 85099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 85199dec3a6SMatthew G. Knepley coords[p*2+0] = PetscRealPart(coords[p*dim+0] - coords[0*dim+0]); 8527df32b8bSSanderA coords[p*2+1] = (PetscRealPart(coords[p*dim+1] - coords[0*dim+1])) * s; 85373868372SMatthew G. Knepley } 85499dec3a6SMatthew G. Knepley coords[0] = 0.0; 85599dec3a6SMatthew G. Knepley coords[1] = 0.0; 8567df32b8bSSanderA coords[2] = x1[0]; 8577df32b8bSSanderA coords[3] = x1[1] * s; 8587df32b8bSSanderA coords[4] = x2[0]; 8597df32b8bSSanderA coords[5] = x2[1] * s; 8607df32b8bSSanderA R[0] = 1.0; R[1] = 0.0; R[2] = 0.0; 8617df32b8bSSanderA R[3] = 0.0; R[4] = 1.0 * s; R[5] = 0.0; 8627df32b8bSSanderA R[6] = 0.0; R[7] = 0.0; R[8] = 1.0 * s; 86373868372SMatthew G. Knepley PetscFunctionReturn(0); 86473868372SMatthew G. Knepley } 865da18b5e6SMatthew G Knepley alpha = 1.0/sqrtz; 866ccd2543fSMatthew G Knepley R[0] = alpha*n[0]*n[2]; R[1] = alpha*n[1]*n[2]; R[2] = -sqrtz; 867ccd2543fSMatthew G Knepley R[3] = -alpha*n[1]; R[4] = alpha*n[0]; R[5] = 0.0; 868ccd2543fSMatthew G Knepley R[6] = n[0]; R[7] = n[1]; R[8] = n[2]; 869ccd2543fSMatthew G Knepley for (d = 0; d < dim; ++d) { 870ccd2543fSMatthew G Knepley x1p[d] = 0.0; 871ccd2543fSMatthew G Knepley x2p[d] = 0.0; 872ccd2543fSMatthew G Knepley for (e = 0; e < dim; ++e) { 873ccd2543fSMatthew G Knepley x1p[d] += R[d*dim+e]*x1[e]; 874ccd2543fSMatthew G Knepley x2p[d] += R[d*dim+e]*x2[e]; 875ccd2543fSMatthew G Knepley } 876ccd2543fSMatthew G Knepley } 87748120919SToby Isaac if (PetscAbsReal(x1p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 87848120919SToby Isaac if (PetscAbsReal(x2p[2]) > 10. * PETSC_SMALL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid rotation calculated"); 879ccd2543fSMatthew G Knepley /* 2) Project to (x, y) */ 88099dec3a6SMatthew G. Knepley for (p = 3; p < coordSize/3; ++p) { 88199dec3a6SMatthew G. Knepley for (d = 0; d < dim; ++d) { 88299dec3a6SMatthew G. Knepley xnp[d] = 0.0; 88399dec3a6SMatthew G. Knepley for (e = 0; e < dim; ++e) { 88499dec3a6SMatthew G. Knepley xnp[d] += R[d*dim+e]*PetscRealPart(coords[p*dim+e] - coords[0*dim+e]); 88599dec3a6SMatthew G. Knepley } 88699dec3a6SMatthew G. Knepley if (d < dim-1) coords[p*2+d] = xnp[d]; 88799dec3a6SMatthew G. Knepley } 88899dec3a6SMatthew G. Knepley } 889ccd2543fSMatthew G Knepley coords[0] = 0.0; 890ccd2543fSMatthew G Knepley coords[1] = 0.0; 891ccd2543fSMatthew G Knepley coords[2] = x1p[0]; 892ccd2543fSMatthew G Knepley coords[3] = x1p[1]; 893ccd2543fSMatthew G Knepley coords[4] = x2p[0]; 894ccd2543fSMatthew G Knepley coords[5] = x2p[1]; 8957f07f362SMatthew G. Knepley /* Output R^T which rotates \hat z to the input normal */ 8967f07f362SMatthew G. Knepley for (d = 0; d < dim; ++d) { 8977f07f362SMatthew G. Knepley for (e = d+1; e < dim; ++e) { 8987f07f362SMatthew G. Knepley PetscReal tmp; 8997f07f362SMatthew G. Knepley 9007f07f362SMatthew G. Knepley tmp = R[d*dim+e]; 9017f07f362SMatthew G. Knepley R[d*dim+e] = R[e*dim+d]; 9027f07f362SMatthew G. Knepley R[e*dim+d] = tmp; 9037f07f362SMatthew G. Knepley } 9047f07f362SMatthew G. Knepley } 905ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 906ccd2543fSMatthew G Knepley } 907ccd2543fSMatthew G Knepley 9086322fe33SJed Brown PETSC_UNUSED 909834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Internal(PetscReal *vol, PetscReal coords[]) 910834e62ceSMatthew G. Knepley { 911834e62ceSMatthew G. Knepley /* Signed volume is 1/2 the determinant 912834e62ceSMatthew G. Knepley 913834e62ceSMatthew G. Knepley | 1 1 1 | 914834e62ceSMatthew G. Knepley | x0 x1 x2 | 915834e62ceSMatthew G. Knepley | y0 y1 y2 | 916834e62ceSMatthew G. Knepley 917834e62ceSMatthew G. Knepley but if x0,y0 is the origin, we have 918834e62ceSMatthew G. Knepley 919834e62ceSMatthew G. Knepley | x1 x2 | 920834e62ceSMatthew G. Knepley | y1 y2 | 921834e62ceSMatthew G. Knepley */ 922834e62ceSMatthew G. Knepley const PetscReal x1 = coords[2] - coords[0], y1 = coords[3] - coords[1]; 923834e62ceSMatthew G. Knepley const PetscReal x2 = coords[4] - coords[0], y2 = coords[5] - coords[1]; 924834e62ceSMatthew G. Knepley PetscReal M[4], detM; 925834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; 92686623015SMatthew G. Knepley M[2] = y1; M[3] = y2; 927923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(&detM, M); 928834e62ceSMatthew G. Knepley *vol = 0.5*detM; 9293bc0b13bSBarry Smith (void)PetscLogFlops(5.0); 930834e62ceSMatthew G. Knepley } 931834e62ceSMatthew G. Knepley 932834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Triangle_Origin_Internal(PetscReal *vol, PetscReal coords[]) 933834e62ceSMatthew G. Knepley { 934923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(vol, coords); 935834e62ceSMatthew G. Knepley *vol *= 0.5; 936834e62ceSMatthew G. Knepley } 937834e62ceSMatthew G. Knepley 9386322fe33SJed Brown PETSC_UNUSED 939834e62ceSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Internal(PetscReal *vol, PetscReal coords[]) 940834e62ceSMatthew G. Knepley { 941834e62ceSMatthew G. Knepley /* Signed volume is 1/6th of the determinant 942834e62ceSMatthew G. Knepley 943834e62ceSMatthew G. Knepley | 1 1 1 1 | 944834e62ceSMatthew G. Knepley | x0 x1 x2 x3 | 945834e62ceSMatthew G. Knepley | y0 y1 y2 y3 | 946834e62ceSMatthew G. Knepley | z0 z1 z2 z3 | 947834e62ceSMatthew G. Knepley 948834e62ceSMatthew G. Knepley but if x0,y0,z0 is the origin, we have 949834e62ceSMatthew G. Knepley 950834e62ceSMatthew G. Knepley | x1 x2 x3 | 951834e62ceSMatthew G. Knepley | y1 y2 y3 | 952834e62ceSMatthew G. Knepley | z1 z2 z3 | 953834e62ceSMatthew G. Knepley */ 954834e62ceSMatthew G. Knepley const PetscReal x1 = coords[3] - coords[0], y1 = coords[4] - coords[1], z1 = coords[5] - coords[2]; 955834e62ceSMatthew G. Knepley const PetscReal x2 = coords[6] - coords[0], y2 = coords[7] - coords[1], z2 = coords[8] - coords[2]; 956834e62ceSMatthew G. Knepley const PetscReal x3 = coords[9] - coords[0], y3 = coords[10] - coords[1], z3 = coords[11] - coords[2]; 957834e62ceSMatthew G. Knepley PetscReal M[9], detM; 958834e62ceSMatthew G. Knepley M[0] = x1; M[1] = x2; M[2] = x3; 959834e62ceSMatthew G. Knepley M[3] = y1; M[4] = y2; M[5] = y3; 960834e62ceSMatthew G. Knepley M[6] = z1; M[7] = z2; M[8] = z3; 961923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(&detM, M); 962b7ad821dSMatthew G. Knepley *vol = -0.16666666666666666666666*detM; 9633bc0b13bSBarry Smith (void)PetscLogFlops(10.0); 964834e62ceSMatthew G. Knepley } 965834e62ceSMatthew G. Knepley 9660ec8681fSMatthew G. Knepley PETSC_STATIC_INLINE void Volume_Tetrahedron_Origin_Internal(PetscReal *vol, PetscReal coords[]) 9670ec8681fSMatthew G. Knepley { 968923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(vol, coords); 969b7ad821dSMatthew G. Knepley *vol *= -0.16666666666666666666666; 9700ec8681fSMatthew G. Knepley } 9710ec8681fSMatthew G. Knepley 972cb92db44SToby Isaac static PetscErrorCode DMPlexComputePointGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 973cb92db44SToby Isaac { 974cb92db44SToby Isaac PetscSection coordSection; 975cb92db44SToby Isaac Vec coordinates; 976cb92db44SToby Isaac const PetscScalar *coords; 977cb92db44SToby Isaac PetscInt dim, d, off; 978cb92db44SToby Isaac PetscErrorCode ierr; 979cb92db44SToby Isaac 980cb92db44SToby Isaac PetscFunctionBegin; 981cb92db44SToby Isaac ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 982cb92db44SToby Isaac ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 983cb92db44SToby Isaac ierr = PetscSectionGetDof(coordSection,e,&dim);CHKERRQ(ierr); 984cb92db44SToby Isaac if (!dim) PetscFunctionReturn(0); 985cb92db44SToby Isaac ierr = PetscSectionGetOffset(coordSection,e,&off);CHKERRQ(ierr); 986cb92db44SToby Isaac ierr = VecGetArrayRead(coordinates,&coords);CHKERRQ(ierr); 987cb92db44SToby Isaac if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[off + d]);} 988cb92db44SToby Isaac ierr = VecRestoreArrayRead(coordinates,&coords);CHKERRQ(ierr); 989cb92db44SToby Isaac *detJ = 1.; 990cb92db44SToby Isaac if (J) { 991cb92db44SToby Isaac for (d = 0; d < dim * dim; d++) J[d] = 0.; 992cb92db44SToby Isaac for (d = 0; d < dim; d++) J[d * dim + d] = 1.; 993cb92db44SToby Isaac if (invJ) { 994cb92db44SToby Isaac for (d = 0; d < dim * dim; d++) invJ[d] = 0.; 995cb92db44SToby Isaac for (d = 0; d < dim; d++) invJ[d * dim + d] = 1.; 996cb92db44SToby Isaac } 997cb92db44SToby Isaac } 998cb92db44SToby Isaac PetscFunctionReturn(0); 999cb92db44SToby Isaac } 1000cb92db44SToby Isaac 100117fe8556SMatthew G. Knepley static PetscErrorCode DMPlexComputeLineGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 100217fe8556SMatthew G. Knepley { 100317fe8556SMatthew G. Knepley PetscSection coordSection; 100417fe8556SMatthew G. Knepley Vec coordinates; 1005a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 10068bf5c034SToby Isaac PetscInt numCoords, d, pStart, pEnd, numSelfCoords = 0; 100717fe8556SMatthew G. Knepley PetscErrorCode ierr; 100817fe8556SMatthew G. Knepley 100917fe8556SMatthew G. Knepley PetscFunctionBegin; 101017fe8556SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 101169d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 10128bf5c034SToby Isaac ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 10138bf5c034SToby Isaac if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 101417fe8556SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 10158bf5c034SToby Isaac numCoords = numSelfCoords ? numSelfCoords : numCoords; 1016adac9986SMatthew G. Knepley if (invJ && !J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In order to compute invJ, J must not be NULL"); 10177f07f362SMatthew G. Knepley *detJ = 0.0; 101828dbe442SToby Isaac if (numCoords == 6) { 101928dbe442SToby Isaac const PetscInt dim = 3; 102028dbe442SToby Isaac PetscReal R[9], J0; 102128dbe442SToby Isaac 102228dbe442SToby Isaac if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1023741bfc07SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto1D(coords, R);CHKERRQ(ierr); 102428dbe442SToby Isaac if (J) { 102528dbe442SToby Isaac J0 = 0.5*PetscRealPart(coords[1]); 102628dbe442SToby Isaac J[0] = R[0]*J0; J[1] = R[1]; J[2] = R[2]; 102728dbe442SToby Isaac J[3] = R[3]*J0; J[4] = R[4]; J[5] = R[5]; 102828dbe442SToby Isaac J[6] = R[6]*J0; J[7] = R[7]; J[8] = R[8]; 102928dbe442SToby Isaac DMPlex_Det3D_Internal(detJ, J); 103028dbe442SToby Isaac if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1031adac9986SMatthew G. Knepley } 103228dbe442SToby Isaac } else if (numCoords == 4) { 10337f07f362SMatthew G. Knepley const PetscInt dim = 2; 10347f07f362SMatthew G. Knepley PetscReal R[4], J0; 10357f07f362SMatthew G. Knepley 10367f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1037741bfc07SMatthew G. Knepley ierr = DMPlexComputeProjection2Dto1D(coords, R);CHKERRQ(ierr); 103817fe8556SMatthew G. Knepley if (J) { 10397f07f362SMatthew G. Knepley J0 = 0.5*PetscRealPart(coords[1]); 10407f07f362SMatthew G. Knepley J[0] = R[0]*J0; J[1] = R[1]; 10417f07f362SMatthew G. Knepley J[2] = R[2]*J0; J[3] = R[3]; 1042923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 1043923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1044adac9986SMatthew G. Knepley } 10457f07f362SMatthew G. Knepley } else if (numCoords == 2) { 10467f07f362SMatthew G. Knepley const PetscInt dim = 1; 10477f07f362SMatthew G. Knepley 10487f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 10497f07f362SMatthew G. Knepley if (J) { 10507f07f362SMatthew G. Knepley J[0] = 0.5*(PetscRealPart(coords[1]) - PetscRealPart(coords[0])); 105117fe8556SMatthew G. Knepley *detJ = J[0]; 10523bc0b13bSBarry Smith ierr = PetscLogFlops(2.0);CHKERRQ(ierr); 10533bc0b13bSBarry Smith if (invJ) {invJ[0] = 1.0/J[0]; ierr = PetscLogFlops(1.0);CHKERRQ(ierr);} 1054adac9986SMatthew G. Knepley } 1055796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this segment is %D != 2", numCoords); 105617fe8556SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 105717fe8556SMatthew G. Knepley PetscFunctionReturn(0); 105817fe8556SMatthew G. Knepley } 105917fe8556SMatthew G. Knepley 1060ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTriangleGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1061ccd2543fSMatthew G Knepley { 1062ccd2543fSMatthew G Knepley PetscSection coordSection; 1063ccd2543fSMatthew G Knepley Vec coordinates; 1064a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 10657f07f362SMatthew G. Knepley PetscInt numCoords, d, f, g; 1066ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1067ccd2543fSMatthew G Knepley 1068ccd2543fSMatthew G Knepley PetscFunctionBegin; 1069ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 107069d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1071ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 10727f07f362SMatthew G. Knepley *detJ = 0.0; 1073ccd2543fSMatthew G Knepley if (numCoords == 9) { 10747f07f362SMatthew G. Knepley const PetscInt dim = 3; 10757f07f362SMatthew 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}; 10767f07f362SMatthew G. Knepley 10777f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1078741bfc07SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 10797f07f362SMatthew G. Knepley if (J) { 1080b7ad821dSMatthew G. Knepley const PetscInt pdim = 2; 1081b7ad821dSMatthew G. Knepley 1082b7ad821dSMatthew G. Knepley for (d = 0; d < pdim; d++) { 1083b7ad821dSMatthew G. Knepley for (f = 0; f < pdim; f++) { 1084b7ad821dSMatthew G. Knepley J0[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 1085ccd2543fSMatthew G Knepley } 10867f07f362SMatthew G. Knepley } 10873bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1088923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 10897f07f362SMatthew G. Knepley for (d = 0; d < dim; d++) { 10907f07f362SMatthew G. Knepley for (f = 0; f < dim; f++) { 10917f07f362SMatthew G. Knepley J[d*dim+f] = 0.0; 10927f07f362SMatthew G. Knepley for (g = 0; g < dim; g++) { 10937f07f362SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 10947f07f362SMatthew G. Knepley } 10957f07f362SMatthew G. Knepley } 10967f07f362SMatthew G. Knepley } 10973bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 10987f07f362SMatthew G. Knepley } 1099923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 11007f07f362SMatthew G. Knepley } else if (numCoords == 6) { 11017f07f362SMatthew G. Knepley const PetscInt dim = 2; 11027f07f362SMatthew G. Knepley 11037f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1104ccd2543fSMatthew G Knepley if (J) { 1105ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 1106ccd2543fSMatthew G Knepley for (f = 0; f < dim; f++) { 1107ccd2543fSMatthew G Knepley J[d*dim+f] = 0.5*(PetscRealPart(coords[(f+1)*dim+d]) - PetscRealPart(coords[0*dim+d])); 1108ccd2543fSMatthew G Knepley } 1109ccd2543fSMatthew G Knepley } 11103bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1111923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 1112ccd2543fSMatthew G Knepley } 1113923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1114796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this triangle is %D != 6 or 9", numCoords); 1115ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1116ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1117ccd2543fSMatthew G Knepley } 1118ccd2543fSMatthew G Knepley 1119dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeRectangleGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1120ccd2543fSMatthew G Knepley { 1121ccd2543fSMatthew G Knepley PetscSection coordSection; 1122ccd2543fSMatthew G Knepley Vec coordinates; 1123a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 11240d29256aSToby Isaac PetscInt numCoords, numSelfCoords = 0, d, f, g, pStart, pEnd; 1125ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1126ccd2543fSMatthew G Knepley 1127ccd2543fSMatthew G Knepley PetscFunctionBegin; 1128ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 112969d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 11300d29256aSToby Isaac ierr = PetscSectionGetChart(coordSection,&pStart,&pEnd);CHKERRQ(ierr); 11310d29256aSToby Isaac if (e >= pStart && e < pEnd) {ierr = PetscSectionGetDof(coordSection,e,&numSelfCoords);CHKERRQ(ierr);} 113299dec3a6SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 113371f58de1SToby Isaac numCoords = numSelfCoords ? numSelfCoords : numCoords; 1134dfccc68fSToby Isaac if (!Nq) { 11357f07f362SMatthew G. Knepley *detJ = 0.0; 113699dec3a6SMatthew G. Knepley if (numCoords == 12) { 113799dec3a6SMatthew G. Knepley const PetscInt dim = 3; 113899dec3a6SMatthew 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}; 113999dec3a6SMatthew G. Knepley 1140dfccc68fSToby Isaac if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1141741bfc07SMatthew G. Knepley ierr = DMPlexComputeProjection3Dto2D(numCoords, coords, R);CHKERRQ(ierr); 114299dec3a6SMatthew G. Knepley if (J) { 114399dec3a6SMatthew G. Knepley const PetscInt pdim = 2; 114499dec3a6SMatthew G. Knepley 114599dec3a6SMatthew G. Knepley for (d = 0; d < pdim; d++) { 114699dec3a6SMatthew G. Knepley J0[d*dim+0] = 0.5*(PetscRealPart(coords[1*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 114799dec3a6SMatthew G. Knepley J0[d*dim+1] = 0.5*(PetscRealPart(coords[3*pdim+d]) - PetscRealPart(coords[0*pdim+d])); 114899dec3a6SMatthew G. Knepley } 11493bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1150923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J0); 115199dec3a6SMatthew G. Knepley for (d = 0; d < dim; d++) { 115299dec3a6SMatthew G. Knepley for (f = 0; f < dim; f++) { 115399dec3a6SMatthew G. Knepley J[d*dim+f] = 0.0; 115499dec3a6SMatthew G. Knepley for (g = 0; g < dim; g++) { 115599dec3a6SMatthew G. Knepley J[d*dim+f] += R[d*dim+g]*J0[g*dim+f]; 115699dec3a6SMatthew G. Knepley } 115799dec3a6SMatthew G. Knepley } 115899dec3a6SMatthew G. Knepley } 11593bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 116099dec3a6SMatthew G. Knepley } 1161923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 116271f58de1SToby Isaac } else if (numCoords == 8) { 116399dec3a6SMatthew G. Knepley const PetscInt dim = 2; 116499dec3a6SMatthew G. Knepley 1165dfccc68fSToby Isaac if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1166ccd2543fSMatthew G Knepley if (J) { 1167ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 116899dec3a6SMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 116999dec3a6SMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1170ccd2543fSMatthew G Knepley } 11713bc0b13bSBarry Smith ierr = PetscLogFlops(8.0);CHKERRQ(ierr); 1172923591dfSMatthew G. Knepley DMPlex_Det2D_Internal(detJ, J); 1173ccd2543fSMatthew G Knepley } 1174923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert2D_Internal(invJ, J, *detJ);} 1175796f034aSJed Brown } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 1176dfccc68fSToby Isaac } else { 1177dfccc68fSToby Isaac const PetscInt Nv = 4; 1178dfccc68fSToby Isaac const PetscInt dimR = 2; 1179dfccc68fSToby Isaac const PetscInt zToPlex[4] = {0, 1, 3, 2}; 1180dfccc68fSToby Isaac PetscReal zOrder[12]; 1181dfccc68fSToby Isaac PetscReal zCoeff[12]; 1182dfccc68fSToby Isaac PetscInt i, j, k, l, dim; 1183dfccc68fSToby Isaac 1184dfccc68fSToby Isaac if (numCoords == 12) { 1185dfccc68fSToby Isaac dim = 3; 1186dfccc68fSToby Isaac } else if (numCoords == 8) { 1187dfccc68fSToby Isaac dim = 2; 1188dfccc68fSToby Isaac } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of coordinates for this quadrilateral is %D != 8 or 12", numCoords); 1189dfccc68fSToby Isaac for (i = 0; i < Nv; i++) { 1190dfccc68fSToby Isaac PetscInt zi = zToPlex[i]; 1191dfccc68fSToby Isaac 1192dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1193dfccc68fSToby Isaac zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 1194dfccc68fSToby Isaac } 1195dfccc68fSToby Isaac } 1196dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1197dfccc68fSToby Isaac zCoeff[dim * 0 + j] = 0.25 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1198dfccc68fSToby Isaac zCoeff[dim * 1 + j] = 0.25 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1199dfccc68fSToby Isaac zCoeff[dim * 2 + j] = 0.25 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1200dfccc68fSToby Isaac zCoeff[dim * 3 + j] = 0.25 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j]); 1201dfccc68fSToby Isaac } 1202dfccc68fSToby Isaac for (i = 0; i < Nq; i++) { 1203dfccc68fSToby Isaac PetscReal xi = points[dimR * i], eta = points[dimR * i + 1]; 1204dfccc68fSToby Isaac 1205dfccc68fSToby Isaac if (v) { 1206dfccc68fSToby Isaac PetscReal extPoint[4]; 1207dfccc68fSToby Isaac 1208dfccc68fSToby Isaac extPoint[0] = 1.; 1209dfccc68fSToby Isaac extPoint[1] = xi; 1210dfccc68fSToby Isaac extPoint[2] = eta; 1211dfccc68fSToby Isaac extPoint[3] = xi * eta; 1212dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1213dfccc68fSToby Isaac PetscReal val = 0.; 1214dfccc68fSToby Isaac 1215dfccc68fSToby Isaac for (k = 0; k < Nv; k++) { 1216dfccc68fSToby Isaac val += extPoint[k] * zCoeff[dim * k + j]; 1217dfccc68fSToby Isaac } 1218dfccc68fSToby Isaac v[i * dim + j] = val; 1219dfccc68fSToby Isaac } 1220dfccc68fSToby Isaac } 1221dfccc68fSToby Isaac if (J) { 1222dfccc68fSToby Isaac PetscReal extJ[8]; 1223dfccc68fSToby Isaac 1224dfccc68fSToby Isaac extJ[0] = 0.; 1225dfccc68fSToby Isaac extJ[1] = 0.; 1226dfccc68fSToby Isaac extJ[2] = 1.; 1227dfccc68fSToby Isaac extJ[3] = 0.; 1228dfccc68fSToby Isaac extJ[4] = 0.; 1229dfccc68fSToby Isaac extJ[5] = 1.; 1230dfccc68fSToby Isaac extJ[6] = eta; 1231dfccc68fSToby Isaac extJ[7] = xi; 1232dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1233dfccc68fSToby Isaac for (k = 0; k < dimR; k++) { 1234dfccc68fSToby Isaac PetscReal val = 0.; 1235dfccc68fSToby Isaac 1236dfccc68fSToby Isaac for (l = 0; l < Nv; l++) { 1237dfccc68fSToby Isaac val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 1238dfccc68fSToby Isaac } 1239dfccc68fSToby Isaac J[i * dim * dim + dim * j + k] = val; 1240dfccc68fSToby Isaac } 1241dfccc68fSToby Isaac } 1242dfccc68fSToby Isaac if (dim == 3) { /* put the cross product in the third component of the Jacobian */ 1243dfccc68fSToby Isaac PetscReal x, y, z; 1244dfccc68fSToby Isaac PetscReal *iJ = &J[i * dim * dim]; 1245dfccc68fSToby Isaac PetscReal norm; 1246dfccc68fSToby Isaac 1247dfccc68fSToby Isaac x = iJ[1 * dim + 0] * iJ[2 * dim + 1] - iJ[1 * dim + 1] * iJ[2 * dim + 0]; 1248dfccc68fSToby Isaac y = iJ[0 * dim + 1] * iJ[2 * dim + 0] - iJ[0 * dim + 0] * iJ[2 * dim + 1]; 1249dfccc68fSToby Isaac z = iJ[0 * dim + 0] * iJ[1 * dim + 1] - iJ[0 * dim + 1] * iJ[1 * dim + 0]; 1250dfccc68fSToby Isaac norm = PetscSqrtReal(x * x + y * y + z * z); 1251dfccc68fSToby Isaac iJ[2] = x / norm; 1252dfccc68fSToby Isaac iJ[5] = y / norm; 1253dfccc68fSToby Isaac iJ[8] = z / norm; 1254dfccc68fSToby Isaac DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 1255dfccc68fSToby Isaac if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1256dfccc68fSToby Isaac } else { 1257dfccc68fSToby Isaac DMPlex_Det2D_Internal(&detJ[i], &J[i * dim * dim]); 1258dfccc68fSToby Isaac if (invJ) {DMPlex_Invert2D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1259dfccc68fSToby Isaac } 1260dfccc68fSToby Isaac } 1261dfccc68fSToby Isaac } 1262dfccc68fSToby Isaac } 126399dec3a6SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, &numCoords, &coords);CHKERRQ(ierr); 1264ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1265ccd2543fSMatthew G Knepley } 1266ccd2543fSMatthew G Knepley 1267ccd2543fSMatthew G Knepley static PetscErrorCode DMPlexComputeTetrahedronGeometry_Internal(DM dm, PetscInt e, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1268ccd2543fSMatthew G Knepley { 1269ccd2543fSMatthew G Knepley PetscSection coordSection; 1270ccd2543fSMatthew G Knepley Vec coordinates; 1271a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 1272ccd2543fSMatthew G Knepley const PetscInt dim = 3; 127399dec3a6SMatthew G. Knepley PetscInt d; 1274ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1275ccd2543fSMatthew G Knepley 1276ccd2543fSMatthew G Knepley PetscFunctionBegin; 1277ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 127869d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1279ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 12807f07f362SMatthew G. Knepley *detJ = 0.0; 12817f07f362SMatthew G. Knepley if (v0) {for (d = 0; d < dim; d++) v0[d] = PetscRealPart(coords[d]);} 1282ccd2543fSMatthew G Knepley if (J) { 1283ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 1284f0df753eSMatthew G. Knepley /* I orient with outward face normals */ 1285f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[2*dim+d]) - PetscRealPart(coords[0*dim+d])); 1286f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1287f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1288ccd2543fSMatthew G Knepley } 12893bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1290923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 1291ccd2543fSMatthew G Knepley } 1292923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1293ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1294ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1295ccd2543fSMatthew G Knepley } 1296ccd2543fSMatthew G Knepley 1297dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeHexahedronGeometry_Internal(DM dm, PetscInt e, PetscInt Nq, const PetscReal points[], PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 1298ccd2543fSMatthew G Knepley { 1299ccd2543fSMatthew G Knepley PetscSection coordSection; 1300ccd2543fSMatthew G Knepley Vec coordinates; 1301a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 1302ccd2543fSMatthew G Knepley const PetscInt dim = 3; 1303ccd2543fSMatthew G Knepley PetscInt d; 1304ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1305ccd2543fSMatthew G Knepley 1306ccd2543fSMatthew G Knepley PetscFunctionBegin; 1307ccd2543fSMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 130869d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1309ccd2543fSMatthew G Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1310dfccc68fSToby Isaac if (!Nq) { 13117f07f362SMatthew G. Knepley *detJ = 0.0; 1312dfccc68fSToby Isaac if (v) {for (d = 0; d < dim; d++) v[d] = PetscRealPart(coords[d]);} 1313ccd2543fSMatthew G Knepley if (J) { 1314ccd2543fSMatthew G Knepley for (d = 0; d < dim; d++) { 1315f0df753eSMatthew G. Knepley J[d*dim+0] = 0.5*(PetscRealPart(coords[3*dim+d]) - PetscRealPart(coords[0*dim+d])); 1316f0df753eSMatthew G. Knepley J[d*dim+1] = 0.5*(PetscRealPart(coords[1*dim+d]) - PetscRealPart(coords[0*dim+d])); 1317f0df753eSMatthew G. Knepley J[d*dim+2] = 0.5*(PetscRealPart(coords[4*dim+d]) - PetscRealPart(coords[0*dim+d])); 1318ccd2543fSMatthew G Knepley } 13193bc0b13bSBarry Smith ierr = PetscLogFlops(18.0);CHKERRQ(ierr); 1320923591dfSMatthew G. Knepley DMPlex_Det3D_Internal(detJ, J); 1321ccd2543fSMatthew G Knepley } 1322923591dfSMatthew G. Knepley if (invJ) {DMPlex_Invert3D_Internal(invJ, J, *detJ);} 1323dfccc68fSToby Isaac } else { 1324dfccc68fSToby Isaac const PetscInt Nv = 8; 1325dfccc68fSToby Isaac const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 1326dfccc68fSToby Isaac const PetscInt dim = 3; 1327dfccc68fSToby Isaac const PetscInt dimR = 3; 1328dfccc68fSToby Isaac PetscReal zOrder[24]; 1329dfccc68fSToby Isaac PetscReal zCoeff[24]; 1330dfccc68fSToby Isaac PetscInt i, j, k, l; 1331dfccc68fSToby Isaac 1332dfccc68fSToby Isaac for (i = 0; i < Nv; i++) { 1333dfccc68fSToby Isaac PetscInt zi = zToPlex[i]; 1334dfccc68fSToby Isaac 1335dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1336dfccc68fSToby Isaac zOrder[dim * i + j] = PetscRealPart(coords[dim * zi + j]); 1337dfccc68fSToby Isaac } 1338dfccc68fSToby Isaac } 1339dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1340dfccc68fSToby Isaac zCoeff[dim * 0 + j] = 0.125 * ( zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1341dfccc68fSToby Isaac zCoeff[dim * 1 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1342dfccc68fSToby Isaac zCoeff[dim * 2 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] + zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1343dfccc68fSToby Isaac zCoeff[dim * 3 + j] = 0.125 * ( zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] + zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1344dfccc68fSToby Isaac zCoeff[dim * 4 + j] = 0.125 * (- zOrder[dim * 0 + j] - zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] + zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1345dfccc68fSToby Isaac zCoeff[dim * 5 + j] = 0.125 * (+ zOrder[dim * 0 + j] - zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] + zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1346dfccc68fSToby Isaac zCoeff[dim * 6 + j] = 0.125 * (+ zOrder[dim * 0 + j] + zOrder[dim * 1 + j] - zOrder[dim * 2 + j] - zOrder[dim * 3 + j] - zOrder[dim * 4 + j] - zOrder[dim * 5 + j] + zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1347dfccc68fSToby Isaac zCoeff[dim * 7 + j] = 0.125 * (- zOrder[dim * 0 + j] + zOrder[dim * 1 + j] + zOrder[dim * 2 + j] - zOrder[dim * 3 + j] + zOrder[dim * 4 + j] - zOrder[dim * 5 + j] - zOrder[dim * 6 + j] + zOrder[dim * 7 + j]); 1348dfccc68fSToby Isaac } 1349dfccc68fSToby Isaac for (i = 0; i < Nq; i++) { 1350dfccc68fSToby Isaac PetscReal xi = points[dimR * i], eta = points[dimR * i + 1], theta = points[dimR * i + 2]; 1351dfccc68fSToby Isaac 1352dfccc68fSToby Isaac if (v) { 135391d2b7ceSToby Isaac PetscReal extPoint[8]; 1354dfccc68fSToby Isaac 1355dfccc68fSToby Isaac extPoint[0] = 1.; 1356dfccc68fSToby Isaac extPoint[1] = xi; 1357dfccc68fSToby Isaac extPoint[2] = eta; 1358dfccc68fSToby Isaac extPoint[3] = xi * eta; 1359dfccc68fSToby Isaac extPoint[4] = theta; 1360dfccc68fSToby Isaac extPoint[5] = theta * xi; 1361dfccc68fSToby Isaac extPoint[6] = theta * eta; 1362dfccc68fSToby Isaac extPoint[7] = theta * eta * xi; 1363dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1364dfccc68fSToby Isaac PetscReal val = 0.; 1365dfccc68fSToby Isaac 1366dfccc68fSToby Isaac for (k = 0; k < Nv; k++) { 1367dfccc68fSToby Isaac val += extPoint[k] * zCoeff[dim * k + j]; 1368dfccc68fSToby Isaac } 1369dfccc68fSToby Isaac v[i * dim + j] = val; 1370dfccc68fSToby Isaac } 1371dfccc68fSToby Isaac } 1372dfccc68fSToby Isaac if (J) { 1373dfccc68fSToby Isaac PetscReal extJ[24]; 1374dfccc68fSToby Isaac 1375dfccc68fSToby Isaac extJ[0] = 0. ; extJ[1] = 0. ; extJ[2] = 0. ; 1376dfccc68fSToby Isaac extJ[3] = 1. ; extJ[4] = 0. ; extJ[5] = 0. ; 1377dfccc68fSToby Isaac extJ[6] = 0. ; extJ[7] = 1. ; extJ[8] = 0. ; 1378dfccc68fSToby Isaac extJ[9] = eta ; extJ[10] = xi ; extJ[11] = 0. ; 1379dfccc68fSToby Isaac extJ[12] = 0. ; extJ[13] = 0. ; extJ[14] = 1. ; 1380dfccc68fSToby Isaac extJ[15] = theta ; extJ[16] = 0. ; extJ[17] = xi ; 1381dfccc68fSToby Isaac extJ[18] = 0. ; extJ[19] = theta ; extJ[20] = eta ; 1382dfccc68fSToby Isaac extJ[21] = theta * eta; extJ[22] = theta * xi; extJ[23] = eta * xi; 1383dfccc68fSToby Isaac 1384dfccc68fSToby Isaac for (j = 0; j < dim; j++) { 1385dfccc68fSToby Isaac for (k = 0; k < dimR; k++) { 1386dfccc68fSToby Isaac PetscReal val = 0.; 1387dfccc68fSToby Isaac 1388dfccc68fSToby Isaac for (l = 0; l < Nv; l++) { 1389dfccc68fSToby Isaac val += zCoeff[dim * l + j] * extJ[dimR * l + k]; 1390dfccc68fSToby Isaac } 1391dfccc68fSToby Isaac J[i * dim * dim + dim * j + k] = val; 1392dfccc68fSToby Isaac } 1393dfccc68fSToby Isaac } 1394dfccc68fSToby Isaac DMPlex_Det3D_Internal(&detJ[i], &J[i * dim * dim]); 1395dfccc68fSToby Isaac if (invJ) {DMPlex_Invert3D_Internal(&invJ[i * dim * dim], &J[i * dim * dim], detJ[i]);} 1396dfccc68fSToby Isaac } 1397dfccc68fSToby Isaac } 1398dfccc68fSToby Isaac } 1399ccd2543fSMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, e, NULL, &coords);CHKERRQ(ierr); 1400ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1401ccd2543fSMatthew G Knepley } 1402ccd2543fSMatthew G Knepley 1403dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_Implicit(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1404dfccc68fSToby Isaac { 1405dfccc68fSToby Isaac PetscInt depth, dim, coordDim, coneSize, i; 1406dfccc68fSToby Isaac PetscInt Nq = 0; 1407dfccc68fSToby Isaac const PetscReal *points = NULL; 1408dfccc68fSToby Isaac DMLabel depthLabel; 14097318780aSToby Isaac PetscReal v0[3] = {-1.}, J0[9] = {-1.}, detJ0 = -1.; 1410dfccc68fSToby Isaac PetscBool isAffine = PETSC_TRUE; 1411dfccc68fSToby Isaac PetscErrorCode ierr; 1412dfccc68fSToby Isaac 1413dfccc68fSToby Isaac PetscFunctionBegin; 1414dfccc68fSToby Isaac ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1415dfccc68fSToby Isaac ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 1416dfccc68fSToby Isaac ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 1417dfccc68fSToby Isaac ierr = DMLabelGetValue(depthLabel, cell, &dim);CHKERRQ(ierr); 1418dfccc68fSToby Isaac if (depth == 1 && dim == 1) { 1419dfccc68fSToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1420dfccc68fSToby Isaac } 1421dfccc68fSToby Isaac ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 1422dfccc68fSToby Isaac if (coordDim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coordinate dimension %D > 3", coordDim); 1423dfccc68fSToby Isaac if (quad) {ierr = PetscQuadratureGetData(quad, NULL, &Nq, &points, NULL);CHKERRQ(ierr);} 1424dfccc68fSToby Isaac switch (dim) { 1425dfccc68fSToby Isaac case 0: 1426dfccc68fSToby Isaac ierr = DMPlexComputePointGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 1427dfccc68fSToby Isaac isAffine = PETSC_FALSE; 1428dfccc68fSToby Isaac break; 1429dfccc68fSToby Isaac case 1: 14307318780aSToby Isaac if (Nq) { 14317318780aSToby Isaac ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 14327318780aSToby Isaac } else { 14337318780aSToby Isaac ierr = DMPlexComputeLineGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 14347318780aSToby Isaac } 1435dfccc68fSToby Isaac break; 1436dfccc68fSToby Isaac case 2: 1437dfccc68fSToby Isaac switch (coneSize) { 1438dfccc68fSToby Isaac case 3: 14397318780aSToby Isaac if (Nq) { 14407318780aSToby Isaac ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 14417318780aSToby Isaac } else { 14427318780aSToby Isaac ierr = DMPlexComputeTriangleGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 14437318780aSToby Isaac } 1444dfccc68fSToby Isaac break; 1445dfccc68fSToby Isaac case 4: 1446dfccc68fSToby Isaac ierr = DMPlexComputeRectangleGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr); 1447dfccc68fSToby Isaac isAffine = PETSC_FALSE; 1448dfccc68fSToby Isaac break; 1449dfccc68fSToby Isaac default: 1450dfccc68fSToby Isaac SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1451dfccc68fSToby Isaac } 1452dfccc68fSToby Isaac break; 1453dfccc68fSToby Isaac case 3: 1454dfccc68fSToby Isaac switch (coneSize) { 1455dfccc68fSToby Isaac case 4: 14567318780aSToby Isaac if (Nq) { 14577318780aSToby Isaac ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v0, J0, NULL, &detJ0);CHKERRQ(ierr); 14587318780aSToby Isaac } else { 14597318780aSToby Isaac ierr = DMPlexComputeTetrahedronGeometry_Internal(dm, cell, v, J, invJ, detJ);CHKERRQ(ierr); 14607318780aSToby Isaac } 1461dfccc68fSToby Isaac break; 1462dfccc68fSToby Isaac case 6: /* Faces */ 1463dfccc68fSToby Isaac case 8: /* Vertices */ 1464dfccc68fSToby Isaac ierr = DMPlexComputeHexahedronGeometry_Internal(dm, cell, Nq, points, v, J, invJ, detJ);CHKERRQ(ierr); 1465dfccc68fSToby Isaac isAffine = PETSC_FALSE; 1466dfccc68fSToby Isaac break; 1467dfccc68fSToby Isaac default: 1468dfccc68fSToby Isaac SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported number of faces %D in cell %D for element geometry computation", coneSize, cell); 1469dfccc68fSToby Isaac } 1470dfccc68fSToby Isaac break; 1471dfccc68fSToby Isaac default: 1472dfccc68fSToby Isaac SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1473dfccc68fSToby Isaac } 14747318780aSToby Isaac if (isAffine && Nq) { 1475dfccc68fSToby Isaac if (v) { 1476dfccc68fSToby Isaac for (i = 0; i < Nq; i++) { 14777318780aSToby Isaac CoordinatesRefToReal(coordDim,dim,v0, J0, &points[dim * i], &v[coordDim * i]); 1478dfccc68fSToby Isaac } 1479dfccc68fSToby Isaac } 14807318780aSToby Isaac if (detJ) { 14817318780aSToby Isaac for (i = 0; i < Nq; i++) { 14827318780aSToby Isaac detJ[i] = detJ0; 1483dfccc68fSToby Isaac } 14847318780aSToby Isaac } 14857318780aSToby Isaac if (J) { 14867318780aSToby Isaac PetscInt k; 14877318780aSToby Isaac 14887318780aSToby Isaac for (i = 0, k = 0; i < Nq; i++) { 1489dfccc68fSToby Isaac PetscInt j; 1490dfccc68fSToby Isaac 14917318780aSToby Isaac for (j = 0; j < coordDim * coordDim; j++, k++) { 14927318780aSToby Isaac J[k] = J0[j]; 14937318780aSToby Isaac } 14947318780aSToby Isaac } 14957318780aSToby Isaac } 14967318780aSToby Isaac if (invJ) { 14977318780aSToby Isaac PetscInt k; 14987318780aSToby Isaac switch (coordDim) { 14997318780aSToby Isaac case 0: 15007318780aSToby Isaac break; 15017318780aSToby Isaac case 1: 15027318780aSToby Isaac invJ[0] = 1./J0[0]; 15037318780aSToby Isaac break; 15047318780aSToby Isaac case 2: 15057318780aSToby Isaac DMPlex_Invert2D_Internal(invJ, J0, detJ0); 15067318780aSToby Isaac break; 15077318780aSToby Isaac case 3: 15087318780aSToby Isaac DMPlex_Invert3D_Internal(invJ, J0, detJ0); 15097318780aSToby Isaac break; 15107318780aSToby Isaac } 15117318780aSToby Isaac for (i = 1, k = coordDim * coordDim; i < Nq; i++) { 15127318780aSToby Isaac PetscInt j; 15137318780aSToby Isaac 15147318780aSToby Isaac for (j = 0; j < coordDim * coordDim; j++, k++) { 15157318780aSToby Isaac invJ[k] = invJ[j]; 15167318780aSToby Isaac } 1517dfccc68fSToby Isaac } 1518dfccc68fSToby Isaac } 1519dfccc68fSToby Isaac } 1520dfccc68fSToby Isaac PetscFunctionReturn(0); 1521dfccc68fSToby Isaac } 1522dfccc68fSToby Isaac 1523ccd2543fSMatthew G Knepley /*@C 15248e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryAffineFEM - Assuming an affine map, compute the Jacobian, inverse Jacobian, and Jacobian determinant for a given cell 1525ccd2543fSMatthew G Knepley 1526ccd2543fSMatthew G Knepley Collective on DM 1527ccd2543fSMatthew G Knepley 1528ccd2543fSMatthew G Knepley Input Arguments: 1529ccd2543fSMatthew G Knepley + dm - the DM 1530ccd2543fSMatthew G Knepley - cell - the cell 1531ccd2543fSMatthew G Knepley 1532ccd2543fSMatthew G Knepley Output Arguments: 1533ccd2543fSMatthew G Knepley + v0 - the translation part of this affine transform 1534ccd2543fSMatthew G Knepley . J - the Jacobian of the transform from the reference element 1535ccd2543fSMatthew G Knepley . invJ - the inverse of the Jacobian 1536ccd2543fSMatthew G Knepley - detJ - the Jacobian determinant 1537ccd2543fSMatthew G Knepley 1538ccd2543fSMatthew G Knepley Level: advanced 1539ccd2543fSMatthew G Knepley 1540ccd2543fSMatthew G Knepley Fortran Notes: 1541ccd2543fSMatthew G Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1542ccd2543fSMatthew G Knepley include petsc.h90 in your code. 1543ccd2543fSMatthew G Knepley 15448e0841e0SMatthew G. Knepley .seealso: DMPlexComputeCellGeometryFEM(), DMGetCoordinateSection(), DMGetCoordinateVec() 1545ccd2543fSMatthew G Knepley @*/ 15468e0841e0SMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM dm, PetscInt cell, PetscReal *v0, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 1547ccd2543fSMatthew G Knepley { 1548ccd2543fSMatthew G Knepley PetscErrorCode ierr; 1549ccd2543fSMatthew G Knepley 1550ccd2543fSMatthew G Knepley PetscFunctionBegin; 1551dfccc68fSToby Isaac ierr = DMPlexComputeCellGeometryFEM_Implicit(dm,cell,NULL,v0,J,invJ,detJ);CHKERRQ(ierr); 15528e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 15538e0841e0SMatthew G. Knepley } 15548e0841e0SMatthew G. Knepley 1555dfccc68fSToby Isaac static PetscErrorCode DMPlexComputeCellGeometryFEM_FE(DM dm, PetscFE fe, PetscInt point, PetscQuadrature quad, PetscReal v[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 15568e0841e0SMatthew G. Knepley { 1557dfccc68fSToby Isaac PetscQuadrature feQuad; 15588e0841e0SMatthew G. Knepley PetscSection coordSection; 15598e0841e0SMatthew G. Knepley Vec coordinates; 15608e0841e0SMatthew G. Knepley PetscScalar *coords = NULL; 15618e0841e0SMatthew G. Knepley const PetscReal *quadPoints; 1562f960e424SToby Isaac PetscReal *basisDer, *basis, detJt; 1563f960e424SToby Isaac PetscInt dim, cdim, pdim, qdim, Nq, numCoords, q; 15648e0841e0SMatthew G. Knepley PetscErrorCode ierr; 15658e0841e0SMatthew G. Knepley 15668e0841e0SMatthew G. Knepley PetscFunctionBegin; 15678e0841e0SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 15688e0841e0SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 15698e0841e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 15708e0841e0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 15718e0841e0SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 1572dfccc68fSToby Isaac if (!quad) { /* use the first point of the first functional of the dual space */ 1573dfccc68fSToby Isaac PetscDualSpace dsp; 1574dfccc68fSToby Isaac 1575dfccc68fSToby Isaac ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr); 1576dfccc68fSToby Isaac ierr = PetscDualSpaceGetFunctional(dsp, 0, &quad);CHKERRQ(ierr); 15778e0841e0SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1578dfccc68fSToby Isaac Nq = 1; 1579dfccc68fSToby Isaac } else { 1580dfccc68fSToby Isaac ierr = PetscQuadratureGetData(quad, &qdim, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 1581dfccc68fSToby Isaac } 158291d2b7ceSToby Isaac ierr = PetscFEGetDimension(fe, &pdim);CHKERRQ(ierr); 1583dfccc68fSToby Isaac ierr = PetscFEGetQuadrature(fe, &feQuad);CHKERRQ(ierr); 1584dfccc68fSToby Isaac if (feQuad == quad) { 1585f960e424SToby Isaac ierr = PetscFEGetDefaultTabulation(fe, &basis, &basisDer, NULL);CHKERRQ(ierr); 15868e0841e0SMatthew 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); 1587dfccc68fSToby Isaac } else { 1588dfccc68fSToby Isaac ierr = PetscFEGetTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr); 1589dfccc68fSToby Isaac } 1590dfccc68fSToby Isaac if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Point dimension %d != quadrature dimension %d", dim, qdim); 1591dfccc68fSToby Isaac if (v) { 1592dfccc68fSToby Isaac ierr = PetscMemzero(v, Nq*cdim*sizeof(PetscReal));CHKERRQ(ierr); 1593f960e424SToby Isaac for (q = 0; q < Nq; ++q) { 1594f960e424SToby Isaac PetscInt i, k; 1595f960e424SToby Isaac 1596f960e424SToby Isaac for (k = 0; k < pdim; ++k) 1597f960e424SToby Isaac for (i = 0; i < cdim; ++i) 1598dfccc68fSToby Isaac v[q*cdim + i] += basis[q*pdim + k] * PetscRealPart(coords[k*cdim + i]); 1599f960e424SToby Isaac ierr = PetscLogFlops(2.0*pdim*cdim);CHKERRQ(ierr); 1600f960e424SToby Isaac } 1601f960e424SToby Isaac } 16028e0841e0SMatthew G. Knepley if (J) { 1603dfccc68fSToby Isaac ierr = PetscMemzero(J, Nq*cdim*cdim*sizeof(PetscReal));CHKERRQ(ierr); 16048e0841e0SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 16058e0841e0SMatthew G. Knepley PetscInt i, j, k, c, r; 16068e0841e0SMatthew G. Knepley 16078e0841e0SMatthew G. Knepley /* J = dx_i/d\xi_j = sum[k=0,n-1] dN_k/d\xi_j * x_i(k) */ 16088e0841e0SMatthew G. Knepley for (k = 0; k < pdim; ++k) 16098e0841e0SMatthew G. Knepley for (j = 0; j < dim; ++j) 16108e0841e0SMatthew G. Knepley for (i = 0; i < cdim; ++i) 1611dfccc68fSToby Isaac J[(q*cdim + i)*cdim + j] += basisDer[(q*pdim + k)*dim + j] * PetscRealPart(coords[k*cdim + i]); 16123bc0b13bSBarry Smith ierr = PetscLogFlops(2.0*pdim*dim*cdim);CHKERRQ(ierr); 16138e0841e0SMatthew G. Knepley if (cdim > dim) { 16148e0841e0SMatthew G. Knepley for (c = dim; c < cdim; ++c) 16158e0841e0SMatthew G. Knepley for (r = 0; r < cdim; ++r) 16168e0841e0SMatthew G. Knepley J[r*cdim+c] = r == c ? 1.0 : 0.0; 16178e0841e0SMatthew G. Knepley } 1618f960e424SToby Isaac if (!detJ && !invJ) continue; 1619a63b72c6SToby Isaac detJt = 0.; 16208e0841e0SMatthew G. Knepley switch (cdim) { 16218e0841e0SMatthew G. Knepley case 3: 1622037dc194SToby Isaac DMPlex_Det3D_Internal(&detJt, &J[q*cdim*dim]); 1623037dc194SToby Isaac if (invJ) {DMPlex_Invert3D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);} 162417fe8556SMatthew G. Knepley break; 162549dc4407SMatthew G. Knepley case 2: 16269f328543SToby Isaac DMPlex_Det2D_Internal(&detJt, &J[q*cdim*dim]); 1627037dc194SToby Isaac if (invJ) {DMPlex_Invert2D_Internal(&invJ[q*cdim*dim], &J[q*cdim*dim], detJt);} 162849dc4407SMatthew G. Knepley break; 16298e0841e0SMatthew G. Knepley case 1: 1630037dc194SToby Isaac detJt = J[q*cdim*dim]; 1631037dc194SToby Isaac if (invJ) invJ[q*cdim*dim] = 1.0/detJt; 163249dc4407SMatthew G. Knepley } 1633f960e424SToby Isaac if (detJ) detJ[q] = detJt; 163449dc4407SMatthew G. Knepley } 163549dc4407SMatthew G. Knepley } 1636037dc194SToby Isaac else if (detJ || invJ) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Need J to compute invJ or detJ"); 1637dfccc68fSToby Isaac if (feQuad != quad) { 1638dfccc68fSToby Isaac ierr = PetscFERestoreTabulation(fe, Nq, quadPoints, &basis, &basisDer, NULL);CHKERRQ(ierr); 1639dfccc68fSToby Isaac } 16408e0841e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, point, &numCoords, &coords);CHKERRQ(ierr); 16418e0841e0SMatthew G. Knepley PetscFunctionReturn(0); 16428e0841e0SMatthew G. Knepley } 16438e0841e0SMatthew G. Knepley 16448e0841e0SMatthew G. Knepley /*@C 16458e0841e0SMatthew G. Knepley DMPlexComputeCellGeometryFEM - Compute the Jacobian, inverse Jacobian, and Jacobian determinant at each quadrature point in the given cell 16468e0841e0SMatthew G. Knepley 16478e0841e0SMatthew G. Knepley Collective on DM 16488e0841e0SMatthew G. Knepley 16498e0841e0SMatthew G. Knepley Input Arguments: 16508e0841e0SMatthew G. Knepley + dm - the DM 16518e0841e0SMatthew G. Knepley . cell - the cell 1652dfccc68fSToby Isaac - quad - the quadrature containing the points in the reference element where the geometry will be evaluated. If quad == NULL, geometry will be 1653dfccc68fSToby Isaac evaluated at the first vertex of the reference element 16548e0841e0SMatthew G. Knepley 16558e0841e0SMatthew G. Knepley Output Arguments: 1656dfccc68fSToby Isaac + v - the image of the transformed quadrature points, otherwise the image of the first vertex in the closure of the reference element 16578e0841e0SMatthew G. Knepley . J - the Jacobian of the transform from the reference element at each quadrature point 16588e0841e0SMatthew G. Knepley . invJ - the inverse of the Jacobian at each quadrature point 16598e0841e0SMatthew G. Knepley - detJ - the Jacobian determinant at each quadrature point 16608e0841e0SMatthew G. Knepley 16618e0841e0SMatthew G. Knepley Level: advanced 16628e0841e0SMatthew G. Knepley 16638e0841e0SMatthew G. Knepley Fortran Notes: 16648e0841e0SMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 16658e0841e0SMatthew G. Knepley include petsc.h90 in your code. 16668e0841e0SMatthew G. Knepley 16678e0841e0SMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 16688e0841e0SMatthew G. Knepley @*/ 1669dfccc68fSToby Isaac PetscErrorCode DMPlexComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal *v, PetscReal *J, PetscReal *invJ, PetscReal *detJ) 16708e0841e0SMatthew G. Knepley { 1671dfccc68fSToby Isaac PetscFE fe = NULL; 16728e0841e0SMatthew G. Knepley PetscErrorCode ierr; 16738e0841e0SMatthew G. Knepley 16748e0841e0SMatthew G. Knepley PetscFunctionBegin; 1675dfccc68fSToby Isaac if (dm->coordinateDM) { 1676dfccc68fSToby Isaac PetscClassId id; 1677dfccc68fSToby Isaac PetscInt numFields; 1678dfccc68fSToby Isaac PetscDS prob = dm->coordinateDM->prob; 1679dfccc68fSToby Isaac PetscObject disc; 1680dfccc68fSToby Isaac 1681dfccc68fSToby Isaac ierr = PetscDSGetNumFields(prob, &numFields);CHKERRQ(ierr); 1682dfccc68fSToby Isaac if (numFields) { 1683dfccc68fSToby Isaac ierr = PetscDSGetDiscretization(prob,0,&disc);CHKERRQ(ierr); 1684dfccc68fSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 1685dfccc68fSToby Isaac if (id == PETSCFE_CLASSID) { 1686dfccc68fSToby Isaac fe = (PetscFE) disc; 1687dfccc68fSToby Isaac } 1688dfccc68fSToby Isaac } 1689dfccc68fSToby Isaac } 1690dfccc68fSToby Isaac if (!fe) {ierr = DMPlexComputeCellGeometryFEM_Implicit(dm, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);} 1691dfccc68fSToby Isaac else {ierr = DMPlexComputeCellGeometryFEM_FE(dm, fe, cell, quad, v, J, invJ, detJ);CHKERRQ(ierr);} 1692ccd2543fSMatthew G Knepley PetscFunctionReturn(0); 1693ccd2543fSMatthew G Knepley } 1694834e62ceSMatthew G. Knepley 1695011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_1D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1696cc08537eSMatthew G. Knepley { 1697cc08537eSMatthew G. Knepley PetscSection coordSection; 1698cc08537eSMatthew G. Knepley Vec coordinates; 1699a1e44745SMatthew G. Knepley PetscScalar *coords = NULL; 170006e2781eSMatthew G. Knepley PetscScalar tmp[2]; 1701cc08537eSMatthew G. Knepley PetscInt coordSize; 1702cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1703cc08537eSMatthew G. Knepley 1704cc08537eSMatthew G. Knepley PetscFunctionBegin; 1705cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 170669d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1707cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1708011ea5d8SMatthew G. Knepley if (dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support 2D edges right now"); 17092e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinate_Internal(dm, dim, coords, &coords[dim], tmp);CHKERRQ(ierr); 1710cc08537eSMatthew G. Knepley if (centroid) { 171106e2781eSMatthew G. Knepley centroid[0] = 0.5*PetscRealPart(coords[0] + tmp[0]); 171206e2781eSMatthew G. Knepley centroid[1] = 0.5*PetscRealPart(coords[1] + tmp[1]); 1713cc08537eSMatthew G. Knepley } 1714cc08537eSMatthew G. Knepley if (normal) { 1715a60a936bSMatthew G. Knepley PetscReal norm; 1716a60a936bSMatthew G. Knepley 171706e2781eSMatthew G. Knepley normal[0] = -PetscRealPart(coords[1] - tmp[1]); 171806e2781eSMatthew G. Knepley normal[1] = PetscRealPart(coords[0] - tmp[0]); 1719a60a936bSMatthew G. Knepley norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1]); 1720a60a936bSMatthew G. Knepley normal[0] /= norm; 1721a60a936bSMatthew G. Knepley normal[1] /= norm; 1722cc08537eSMatthew G. Knepley } 1723cc08537eSMatthew G. Knepley if (vol) { 172406e2781eSMatthew G. Knepley *vol = PetscSqrtReal(PetscSqr(PetscRealPart(coords[0] - tmp[0])) + PetscSqr(PetscRealPart(coords[1] - tmp[1]))); 1725cc08537eSMatthew G. Knepley } 1726cc08537eSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1727cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1728cc08537eSMatthew G. Knepley } 1729cc08537eSMatthew G. Knepley 1730cc08537eSMatthew G. Knepley /* Centroid_i = (\sum_n A_n Cn_i ) / A */ 1731011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_2D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1732cc08537eSMatthew G. Knepley { 1733cc08537eSMatthew G. Knepley PetscSection coordSection; 1734cc08537eSMatthew G. Knepley Vec coordinates; 1735cc08537eSMatthew G. Knepley PetscScalar *coords = NULL; 17360a1d6728SMatthew G. Knepley PetscReal vsum = 0.0, csum[3] = {0.0, 0.0, 0.0}, vtmp, ctmp[4], v0[3], R[9]; 17370a1d6728SMatthew G. Knepley PetscInt tdim = 2, coordSize, numCorners, p, d, e; 1738cc08537eSMatthew G. Knepley PetscErrorCode ierr; 1739cc08537eSMatthew G. Knepley 1740cc08537eSMatthew G. Knepley PetscFunctionBegin; 1741cc08537eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 17420a1d6728SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numCorners);CHKERRQ(ierr); 174369d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1744cc08537eSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 17450bce18caSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 1746ceee4971SMatthew G. Knepley if (dim > 2 && centroid) { 1747ceee4971SMatthew G. Knepley v0[0] = PetscRealPart(coords[0]); 1748ceee4971SMatthew G. Knepley v0[1] = PetscRealPart(coords[1]); 1749ceee4971SMatthew G. Knepley v0[2] = PetscRealPart(coords[2]); 1750ceee4971SMatthew G. Knepley } 1751011ea5d8SMatthew G. Knepley if (normal) { 1752011ea5d8SMatthew G. Knepley if (dim > 2) { 17531ee9d5ecSMatthew G. Knepley const PetscReal x0 = PetscRealPart(coords[dim+0] - coords[0]), x1 = PetscRealPart(coords[dim*2+0] - coords[0]); 17541ee9d5ecSMatthew G. Knepley const PetscReal y0 = PetscRealPart(coords[dim+1] - coords[1]), y1 = PetscRealPart(coords[dim*2+1] - coords[1]); 17551ee9d5ecSMatthew G. Knepley const PetscReal z0 = PetscRealPart(coords[dim+2] - coords[2]), z1 = PetscRealPart(coords[dim*2+2] - coords[2]); 17560a1d6728SMatthew G. Knepley PetscReal norm; 17570a1d6728SMatthew G. Knepley 17580a1d6728SMatthew G. Knepley normal[0] = y0*z1 - z0*y1; 17590a1d6728SMatthew G. Knepley normal[1] = z0*x1 - x0*z1; 17600a1d6728SMatthew G. Knepley normal[2] = x0*y1 - y0*x1; 17618b49ba18SBarry Smith norm = PetscSqrtReal(normal[0]*normal[0] + normal[1]*normal[1] + normal[2]*normal[2]); 17620a1d6728SMatthew G. Knepley normal[0] /= norm; 17630a1d6728SMatthew G. Knepley normal[1] /= norm; 17640a1d6728SMatthew G. Knepley normal[2] /= norm; 1765011ea5d8SMatthew G. Knepley } else { 1766011ea5d8SMatthew G. Knepley for (d = 0; d < dim; ++d) normal[d] = 0.0; 1767011ea5d8SMatthew G. Knepley } 1768011ea5d8SMatthew G. Knepley } 1769741bfc07SMatthew G. Knepley if (dim == 3) {ierr = DMPlexComputeProjection3Dto2D(coordSize, coords, R);CHKERRQ(ierr);} 17700a1d6728SMatthew G. Knepley for (p = 0; p < numCorners; ++p) { 17710a1d6728SMatthew G. Knepley /* Need to do this copy to get types right */ 17720a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 17731ee9d5ecSMatthew G. Knepley ctmp[d] = PetscRealPart(coords[p*tdim+d]); 17741ee9d5ecSMatthew G. Knepley ctmp[tdim+d] = PetscRealPart(coords[((p+1)%numCorners)*tdim+d]); 17750a1d6728SMatthew G. Knepley } 17760a1d6728SMatthew G. Knepley Volume_Triangle_Origin_Internal(&vtmp, ctmp); 17770a1d6728SMatthew G. Knepley vsum += vtmp; 17780a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 17790a1d6728SMatthew G. Knepley csum[d] += (ctmp[d] + ctmp[tdim+d])*vtmp; 17800a1d6728SMatthew G. Knepley } 17810a1d6728SMatthew G. Knepley } 17820a1d6728SMatthew G. Knepley for (d = 0; d < tdim; ++d) { 17830a1d6728SMatthew G. Knepley csum[d] /= (tdim+1)*vsum; 17840a1d6728SMatthew G. Knepley } 17850a1d6728SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords);CHKERRQ(ierr); 1786ee6bbdb2SSatish Balay if (vol) *vol = PetscAbsReal(vsum); 17870a1d6728SMatthew G. Knepley if (centroid) { 17880a1d6728SMatthew G. Knepley if (dim > 2) { 17890a1d6728SMatthew G. Knepley for (d = 0; d < dim; ++d) { 17900a1d6728SMatthew G. Knepley centroid[d] = v0[d]; 17910a1d6728SMatthew G. Knepley for (e = 0; e < dim; ++e) { 17920a1d6728SMatthew G. Knepley centroid[d] += R[d*dim+e]*csum[e]; 17930a1d6728SMatthew G. Knepley } 17940a1d6728SMatthew G. Knepley } 17950a1d6728SMatthew G. Knepley } else for (d = 0; d < dim; ++d) centroid[d] = csum[d]; 17960a1d6728SMatthew G. Knepley } 1797cc08537eSMatthew G. Knepley PetscFunctionReturn(0); 1798cc08537eSMatthew G. Knepley } 1799cc08537eSMatthew G. Knepley 18000ec8681fSMatthew G. Knepley /* Centroid_i = (\sum_n V_n Cn_i ) / V */ 1801011ea5d8SMatthew G. Knepley static PetscErrorCode DMPlexComputeGeometryFVM_3D_Internal(DM dm, PetscInt dim, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 18020ec8681fSMatthew G. Knepley { 18030ec8681fSMatthew G. Knepley PetscSection coordSection; 18040ec8681fSMatthew G. Knepley Vec coordinates; 18050ec8681fSMatthew G. Knepley PetscScalar *coords = NULL; 180686623015SMatthew G. Knepley PetscReal vsum = 0.0, vtmp, coordsTmp[3*3]; 1807a7df9edeSMatthew G. Knepley const PetscInt *faces, *facesO; 18080ec8681fSMatthew G. Knepley PetscInt numFaces, f, coordSize, numCorners, p, d; 18090ec8681fSMatthew G. Knepley PetscErrorCode ierr; 18100ec8681fSMatthew G. Knepley 18110ec8681fSMatthew G. Knepley PetscFunctionBegin; 1812f6dae198SJed Brown if (PetscUnlikely(dim > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No support for dim %D > 3",dim); 18130ec8681fSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 181469d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 18150ec8681fSMatthew G. Knepley 1816d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] = 0.0; 18170ec8681fSMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr); 18180ec8681fSMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 1819a7df9edeSMatthew G. Knepley ierr = DMPlexGetConeOrientation(dm, cell, &facesO);CHKERRQ(ierr); 18200ec8681fSMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 1821011ea5d8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 18220ec8681fSMatthew G. Knepley numCorners = coordSize/dim; 18230ec8681fSMatthew G. Knepley switch (numCorners) { 18240ec8681fSMatthew G. Knepley case 3: 18250ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18261ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 18271ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 18281ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[2*dim+d]); 18290ec8681fSMatthew G. Knepley } 18300ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1831a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 18320ec8681fSMatthew G. Knepley vsum += vtmp; 18334f25033aSJed Brown if (centroid) { /* Centroid of OABC = (a+b+c)/4 */ 18340ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18351ee9d5ecSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 18360ec8681fSMatthew G. Knepley } 18370ec8681fSMatthew G. Knepley } 18380ec8681fSMatthew G. Knepley break; 18390ec8681fSMatthew G. Knepley case 4: 18400ec8681fSMatthew G. Knepley /* DO FOR PYRAMID */ 18410ec8681fSMatthew G. Knepley /* First tet */ 18420ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18431ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[0*dim+d]); 18441ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[1*dim+d]); 18451ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 18460ec8681fSMatthew G. Knepley } 18470ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1848a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 18490ec8681fSMatthew G. Knepley vsum += vtmp; 18500ec8681fSMatthew G. Knepley if (centroid) { 18510ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18520ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 18530ec8681fSMatthew G. Knepley } 18540ec8681fSMatthew G. Knepley } 18550ec8681fSMatthew G. Knepley /* Second tet */ 18560ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18571ee9d5ecSMatthew G. Knepley coordsTmp[0*dim+d] = PetscRealPart(coords[1*dim+d]); 18581ee9d5ecSMatthew G. Knepley coordsTmp[1*dim+d] = PetscRealPart(coords[2*dim+d]); 18591ee9d5ecSMatthew G. Knepley coordsTmp[2*dim+d] = PetscRealPart(coords[3*dim+d]); 18600ec8681fSMatthew G. Knepley } 18610ec8681fSMatthew G. Knepley Volume_Tetrahedron_Origin_Internal(&vtmp, coordsTmp); 1862a7df9edeSMatthew G. Knepley if (facesO[f] < 0) vtmp = -vtmp; 18630ec8681fSMatthew G. Knepley vsum += vtmp; 18640ec8681fSMatthew G. Knepley if (centroid) { 18650ec8681fSMatthew G. Knepley for (d = 0; d < dim; ++d) { 18660ec8681fSMatthew G. Knepley for (p = 0; p < 3; ++p) centroid[d] += coordsTmp[p*dim+d]*vtmp; 18670ec8681fSMatthew G. Knepley } 18680ec8681fSMatthew G. Knepley } 18690ec8681fSMatthew G. Knepley break; 18700ec8681fSMatthew G. Knepley default: 1871796f034aSJed Brown SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle faces with %D vertices", numCorners); 18720ec8681fSMatthew G. Knepley } 18734f25033aSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, faces[f], &coordSize, &coords);CHKERRQ(ierr); 18740ec8681fSMatthew G. Knepley } 18758763be8eSMatthew G. Knepley if (vol) *vol = PetscAbsReal(vsum); 18760ec8681fSMatthew G. Knepley if (normal) for (d = 0; d < dim; ++d) normal[d] = 0.0; 1877d9a81ebdSMatthew G. Knepley if (centroid) for (d = 0; d < dim; ++d) centroid[d] /= (vsum*4); 18780ec8681fSMatthew G. Knepley PetscFunctionReturn(0); 18790ec8681fSMatthew G. Knepley } 18800ec8681fSMatthew G. Knepley 1881834e62ceSMatthew G. Knepley /*@C 1882834e62ceSMatthew G. Knepley DMPlexComputeCellGeometryFVM - Compute the volume for a given cell 1883834e62ceSMatthew G. Knepley 1884834e62ceSMatthew G. Knepley Collective on DM 1885834e62ceSMatthew G. Knepley 1886834e62ceSMatthew G. Knepley Input Arguments: 1887834e62ceSMatthew G. Knepley + dm - the DM 1888834e62ceSMatthew G. Knepley - cell - the cell 1889834e62ceSMatthew G. Knepley 1890834e62ceSMatthew G. Knepley Output Arguments: 1891834e62ceSMatthew G. Knepley + volume - the cell volume 1892cc08537eSMatthew G. Knepley . centroid - the cell centroid 1893cc08537eSMatthew G. Knepley - normal - the cell normal, if appropriate 1894834e62ceSMatthew G. Knepley 1895834e62ceSMatthew G. Knepley Level: advanced 1896834e62ceSMatthew G. Knepley 1897834e62ceSMatthew G. Knepley Fortran Notes: 1898834e62ceSMatthew G. Knepley Since it returns arrays, this routine is only available in Fortran 90, and you must 1899834e62ceSMatthew G. Knepley include petsc.h90 in your code. 1900834e62ceSMatthew G. Knepley 190169d8a9ceSMatthew G. Knepley .seealso: DMGetCoordinateSection(), DMGetCoordinateVec() 1902834e62ceSMatthew G. Knepley @*/ 1903cc08537eSMatthew G. Knepley PetscErrorCode DMPlexComputeCellGeometryFVM(DM dm, PetscInt cell, PetscReal *vol, PetscReal centroid[], PetscReal normal[]) 1904834e62ceSMatthew G. Knepley { 19050ec8681fSMatthew G. Knepley PetscInt depth, dim; 1906834e62ceSMatthew G. Knepley PetscErrorCode ierr; 1907834e62ceSMatthew G. Knepley 1908834e62ceSMatthew G. Knepley PetscFunctionBegin; 1909834e62ceSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1910c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1911834e62ceSMatthew G. Knepley if (depth != dim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 1912834e62ceSMatthew G. Knepley /* We need to keep a pointer to the depth label */ 1913c58f1c22SToby Isaac ierr = DMGetLabelValue(dm, "depth", cell, &depth);CHKERRQ(ierr); 1914834e62ceSMatthew G. Knepley /* Cone size is now the number of faces */ 1915011ea5d8SMatthew G. Knepley switch (depth) { 1916cc08537eSMatthew G. Knepley case 1: 1917011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_1D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1918cc08537eSMatthew G. Knepley break; 1919834e62ceSMatthew G. Knepley case 2: 1920011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_2D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1921834e62ceSMatthew G. Knepley break; 1922834e62ceSMatthew G. Knepley case 3: 1923011ea5d8SMatthew G. Knepley ierr = DMPlexComputeGeometryFVM_3D_Internal(dm, dim, cell, vol, centroid, normal);CHKERRQ(ierr); 1924834e62ceSMatthew G. Knepley break; 1925834e62ceSMatthew G. Knepley default: 1926834e62ceSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension %D for element geometry computation", dim); 1927834e62ceSMatthew G. Knepley } 1928834e62ceSMatthew G. Knepley PetscFunctionReturn(0); 1929834e62ceSMatthew G. Knepley } 1930113c68e6SMatthew G. Knepley 1931c501906fSMatthew G. Knepley /*@ 1932c501906fSMatthew G. Knepley DMPlexComputeGeometryFEM - Precompute cell geometry for the entire mesh 1933c501906fSMatthew G. Knepley 1934c501906fSMatthew G. Knepley Collective on dm 1935c501906fSMatthew G. Knepley 1936c501906fSMatthew G. Knepley Input Parameter: 1937c501906fSMatthew G. Knepley . dm - The DMPlex 1938c501906fSMatthew G. Knepley 1939c501906fSMatthew G. Knepley Output Parameter: 1940c501906fSMatthew G. Knepley . cellgeom - A vector with the cell geometry data for each cell 1941c501906fSMatthew G. Knepley 1942c501906fSMatthew G. Knepley Level: beginner 1943c501906fSMatthew G. Knepley 1944c501906fSMatthew G. Knepley .keywords: DMPlexComputeCellGeometryFEM() 1945c501906fSMatthew G. Knepley @*/ 1946c0d900a5SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFEM(DM dm, Vec *cellgeom) 1947c0d900a5SMatthew G. Knepley { 1948c0d900a5SMatthew G. Knepley DM dmCell; 1949c0d900a5SMatthew G. Knepley Vec coordinates; 1950c0d900a5SMatthew G. Knepley PetscSection coordSection, sectionCell; 1951c0d900a5SMatthew G. Knepley PetscScalar *cgeom; 1952c0d900a5SMatthew G. Knepley PetscInt cStart, cEnd, cMax, c; 1953c0d900a5SMatthew G. Knepley PetscErrorCode ierr; 1954c0d900a5SMatthew G. Knepley 1955c0d900a5SMatthew G. Knepley PetscFunctionBegin; 1956c0d900a5SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 1957c0d900a5SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 1958c0d900a5SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1959c0d900a5SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 1960c0d900a5SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 1961c0d900a5SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 1962c0d900a5SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1963c0d900a5SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);CHKERRQ(ierr); 1964c0d900a5SMatthew G. Knepley cEnd = cMax < 0 ? cEnd : cMax; 1965c0d900a5SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 1966c0d900a5SMatthew G. Knepley /* TODO This needs to be multiplied by Nq for non-affine */ 19679e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFECellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 1968c0d900a5SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 1969c0d900a5SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 1970c0d900a5SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 1971c0d900a5SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 1972c0d900a5SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1973c0d900a5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1974c0d900a5SMatthew G. Knepley PetscFECellGeom *cg; 1975c0d900a5SMatthew G. Knepley 1976c0d900a5SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 1977c0d900a5SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 1978c0d900a5SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmCell, c, NULL, cg->v0, cg->J, cg->invJ, &cg->detJ);CHKERRQ(ierr); 1979c0d900a5SMatthew G. Knepley if (cg->detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cg->detJ, c); 1980c0d900a5SMatthew G. Knepley } 1981c0d900a5SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 1982c0d900a5SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 1983c0d900a5SMatthew G. Knepley PetscFunctionReturn(0); 1984c0d900a5SMatthew G. Knepley } 1985c0d900a5SMatthew G. Knepley 1986891a9168SMatthew G. Knepley /*@ 1987891a9168SMatthew G. Knepley DMPlexComputeGeometryFVM - Computes the cell and face geometry for a finite volume method 1988891a9168SMatthew G. Knepley 1989891a9168SMatthew G. Knepley Input Parameter: 1990891a9168SMatthew G. Knepley . dm - The DM 1991891a9168SMatthew G. Knepley 1992891a9168SMatthew G. Knepley Output Parameters: 1993891a9168SMatthew G. Knepley + cellgeom - A Vec of PetscFVCellGeom data 1994891a9168SMatthew G. Knepley . facegeom - A Vec of PetscFVFaceGeom data 1995891a9168SMatthew G. Knepley 1996891a9168SMatthew G. Knepley Level: developer 1997891a9168SMatthew G. Knepley 1998891a9168SMatthew G. Knepley .seealso: PetscFVFaceGeom, PetscFVCellGeom, DMPlexComputeGeometryFEM() 1999891a9168SMatthew G. Knepley @*/ 2000113c68e6SMatthew G. Knepley PetscErrorCode DMPlexComputeGeometryFVM(DM dm, Vec *cellgeom, Vec *facegeom) 2001113c68e6SMatthew G. Knepley { 2002113c68e6SMatthew G. Knepley DM dmFace, dmCell; 2003113c68e6SMatthew G. Knepley DMLabel ghostLabel; 2004113c68e6SMatthew G. Knepley PetscSection sectionFace, sectionCell; 2005113c68e6SMatthew G. Knepley PetscSection coordSection; 2006113c68e6SMatthew G. Knepley Vec coordinates; 2007113c68e6SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 2008113c68e6SMatthew G. Knepley PetscReal minradius, gminradius; 2009113c68e6SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 2010113c68e6SMatthew G. Knepley PetscErrorCode ierr; 2011113c68e6SMatthew G. Knepley 2012113c68e6SMatthew G. Knepley PetscFunctionBegin; 2013113c68e6SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2014113c68e6SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 2015113c68e6SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 2016113c68e6SMatthew G. Knepley /* Make cell centroids and volumes */ 2017113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 2018113c68e6SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 2019113c68e6SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 2020113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 2021113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2022113c68e6SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2023113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 20249e5edeeeSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVCellGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 2025113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 2026113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 2027113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 2028113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, cellgeom);CHKERRQ(ierr); 202906348e87SToby Isaac if (cEndInterior < 0) { 203006348e87SToby Isaac cEndInterior = cEnd; 203106348e87SToby Isaac } 2032113c68e6SMatthew G. Knepley ierr = VecGetArray(*cellgeom, &cgeom);CHKERRQ(ierr); 2033113c68e6SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 2034113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 2035113c68e6SMatthew G. Knepley 2036113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2037113c68e6SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 2038113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 2039113c68e6SMatthew G. Knepley } 2040113c68e6SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 2041113c68e6SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 2042113c68e6SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 2043113c68e6SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2044113c68e6SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 20459e5edeeeSMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, (PetscInt) PetscCeilReal(((PetscReal) sizeof(PetscFVFaceGeom))/sizeof(PetscScalar)));CHKERRQ(ierr);} 2046113c68e6SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 2047113c68e6SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 2048113c68e6SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 2049113c68e6SMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, facegeom);CHKERRQ(ierr); 2050113c68e6SMatthew G. Knepley ierr = VecGetArray(*facegeom, &fgeom);CHKERRQ(ierr); 2051c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2052113c68e6SMatthew G. Knepley minradius = PETSC_MAX_REAL; 2053113c68e6SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 2054113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 2055113c68e6SMatthew G. Knepley PetscReal area; 205650d63984SToby Isaac PetscInt ghost = -1, d, numChildren; 2057113c68e6SMatthew G. Knepley 20589ac3fadcSMatthew G. Knepley if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 205950d63984SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 206050d63984SToby Isaac if (ghost >= 0 || numChildren) continue; 2061113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 2062113c68e6SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 2063113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 2064113c68e6SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 2065113c68e6SMatthew G. Knepley { 2066113c68e6SMatthew G. Knepley PetscFVCellGeom *cL, *cR; 206706348e87SToby Isaac PetscInt ncells; 2068113c68e6SMatthew G. Knepley const PetscInt *cells; 2069113c68e6SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 20700453c0cdSMatthew G. Knepley PetscReal l[3], r[3], v[3]; 2071113c68e6SMatthew G. Knepley 2072113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 207306348e87SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &ncells);CHKERRQ(ierr); 2074113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 2075113c68e6SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 207606348e87SToby Isaac if (ncells > 1) { 207706348e87SToby Isaac ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 2078113c68e6SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 207906348e87SToby Isaac } 208006348e87SToby Isaac else { 208106348e87SToby Isaac rcentroid = fg->centroid; 208206348e87SToby Isaac } 20832e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, lcentroid, l);CHKERRQ(ierr); 20842e17dfb7SMatthew G. Knepley ierr = DMLocalizeCoordinateReal_Internal(dm, dim, fg->centroid, rcentroid, r);CHKERRQ(ierr); 20850453c0cdSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, l, r, v); 2086113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) < 0) { 2087113c68e6SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 2088113c68e6SMatthew G. Knepley } 2089113c68e6SMatthew G. Knepley if (DMPlex_DotRealD_Internal(dim, fg->normal, v) <= 0) { 2090113c68e6SMatthew 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]); 2091113c68e6SMatthew 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]); 2092113c68e6SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 2093113c68e6SMatthew G. Knepley } 2094113c68e6SMatthew G. Knepley if (cells[0] < cEndInterior) { 2095113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cL->centroid, v); 2096113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 2097113c68e6SMatthew G. Knepley } 209806348e87SToby Isaac if (ncells > 1 && cells[1] < cEndInterior) { 2099113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, fg->centroid, cR->centroid, v); 2100113c68e6SMatthew G. Knepley minradius = PetscMin(minradius, DMPlex_NormD_Internal(dim, v)); 2101113c68e6SMatthew G. Knepley } 2102113c68e6SMatthew G. Knepley } 2103113c68e6SMatthew G. Knepley } 2104b2566f29SBarry Smith ierr = MPIU_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 2105113c68e6SMatthew G. Knepley ierr = DMPlexSetMinRadius(dm, gminradius);CHKERRQ(ierr); 2106113c68e6SMatthew G. Knepley /* Compute centroids of ghost cells */ 2107113c68e6SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 2108113c68e6SMatthew G. Knepley PetscFVFaceGeom *fg; 2109113c68e6SMatthew G. Knepley const PetscInt *cone, *support; 2110113c68e6SMatthew G. Knepley PetscInt coneSize, supportSize, s; 2111113c68e6SMatthew G. Knepley 2112113c68e6SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 2113113c68e6SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 2114113c68e6SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 2115113c68e6SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 211650d63984SToby Isaac if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 2", cone[0], supportSize); 2117113c68e6SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 2118113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 2119113c68e6SMatthew G. Knepley for (s = 0; s < 2; ++s) { 2120113c68e6SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 2121113c68e6SMatthew G. Knepley if (support[s] == c) { 2122640bce14SSatish Balay PetscFVCellGeom *ci; 2123113c68e6SMatthew G. Knepley PetscFVCellGeom *cg; 2124113c68e6SMatthew G. Knepley PetscReal c2f[3], a; 2125113c68e6SMatthew G. Knepley 2126113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 2127113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 2128113c68e6SMatthew G. Knepley a = DMPlex_DotRealD_Internal(dim, c2f, fg->normal)/DMPlex_DotRealD_Internal(dim, fg->normal, fg->normal); 2129113c68e6SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 2130113c68e6SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 2131113c68e6SMatthew G. Knepley cg->volume = ci->volume; 2132113c68e6SMatthew G. Knepley } 2133113c68e6SMatthew G. Knepley } 2134113c68e6SMatthew G. Knepley } 2135113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*facegeom, &fgeom);CHKERRQ(ierr); 2136113c68e6SMatthew G. Knepley ierr = VecRestoreArray(*cellgeom, &cgeom);CHKERRQ(ierr); 2137113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 2138113c68e6SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 2139113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 2140113c68e6SMatthew G. Knepley } 2141113c68e6SMatthew G. Knepley 2142113c68e6SMatthew G. Knepley /*@C 2143113c68e6SMatthew G. Knepley DMPlexGetMinRadius - Returns the minimum distance from any cell centroid to a face 2144113c68e6SMatthew G. Knepley 2145113c68e6SMatthew G. Knepley Not collective 2146113c68e6SMatthew G. Knepley 2147113c68e6SMatthew G. Knepley Input Argument: 2148113c68e6SMatthew G. Knepley . dm - the DM 2149113c68e6SMatthew G. Knepley 2150113c68e6SMatthew G. Knepley Output Argument: 2151113c68e6SMatthew G. Knepley . minradius - the minium cell radius 2152113c68e6SMatthew G. Knepley 2153113c68e6SMatthew G. Knepley Level: developer 2154113c68e6SMatthew G. Knepley 2155113c68e6SMatthew G. Knepley .seealso: DMGetCoordinates() 2156113c68e6SMatthew G. Knepley @*/ 2157113c68e6SMatthew G. Knepley PetscErrorCode DMPlexGetMinRadius(DM dm, PetscReal *minradius) 2158113c68e6SMatthew G. Knepley { 2159113c68e6SMatthew G. Knepley PetscFunctionBegin; 2160113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2161113c68e6SMatthew G. Knepley PetscValidPointer(minradius,2); 2162113c68e6SMatthew G. Knepley *minradius = ((DM_Plex*) dm->data)->minradius; 2163113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 2164113c68e6SMatthew G. Knepley } 2165113c68e6SMatthew G. Knepley 2166113c68e6SMatthew G. Knepley /*@C 2167113c68e6SMatthew G. Knepley DMPlexSetMinRadius - Sets the minimum distance from the cell centroid to a face 2168113c68e6SMatthew G. Knepley 2169113c68e6SMatthew G. Knepley Logically collective 2170113c68e6SMatthew G. Knepley 2171113c68e6SMatthew G. Knepley Input Arguments: 2172113c68e6SMatthew G. Knepley + dm - the DM 2173113c68e6SMatthew G. Knepley - minradius - the minium cell radius 2174113c68e6SMatthew G. Knepley 2175113c68e6SMatthew G. Knepley Level: developer 2176113c68e6SMatthew G. Knepley 2177113c68e6SMatthew G. Knepley .seealso: DMSetCoordinates() 2178113c68e6SMatthew G. Knepley @*/ 2179113c68e6SMatthew G. Knepley PetscErrorCode DMPlexSetMinRadius(DM dm, PetscReal minradius) 2180113c68e6SMatthew G. Knepley { 2181113c68e6SMatthew G. Knepley PetscFunctionBegin; 2182113c68e6SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2183113c68e6SMatthew G. Knepley ((DM_Plex*) dm->data)->minradius = minradius; 2184113c68e6SMatthew G. Knepley PetscFunctionReturn(0); 2185113c68e6SMatthew G. Knepley } 2186856ac710SMatthew G. Knepley 2187856ac710SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction_Internal(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2188856ac710SMatthew G. Knepley { 2189856ac710SMatthew G. Knepley DMLabel ghostLabel; 2190856ac710SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 2191856ac710SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 2192856ac710SMatthew G. Knepley PetscErrorCode ierr; 2193856ac710SMatthew G. Knepley 2194856ac710SMatthew G. Knepley PetscFunctionBegin; 2195856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2196856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2197856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2198856ac710SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 2199856ac710SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 2200c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2201856ac710SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 2202856ac710SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 2203856ac710SMatthew G. Knepley const PetscInt *faces; 2204856ac710SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 2205640bce14SSatish Balay PetscFVCellGeom *cg; 2206856ac710SMatthew G. Knepley PetscBool boundary; 2207856ac710SMatthew G. Knepley PetscInt ghost; 2208856ac710SMatthew G. Knepley 2209856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2210856ac710SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 2211856ac710SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 2212856ac710SMatthew 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); 2213856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 2214640bce14SSatish Balay PetscFVCellGeom *cg1; 2215856ac710SMatthew G. Knepley PetscFVFaceGeom *fg; 2216856ac710SMatthew G. Knepley const PetscInt *fcells; 2217856ac710SMatthew G. Knepley PetscInt ncell, side; 2218856ac710SMatthew G. Knepley 2219856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2220a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2221856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 2222856ac710SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 2223856ac710SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 2224856ac710SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 2225856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 2226856ac710SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 2227856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 2228856ac710SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 2229856ac710SMatthew G. Knepley } 2230856ac710SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 2231856ac710SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 2232856ac710SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 2233856ac710SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2234a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2235856ac710SMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 2236856ac710SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 2237856ac710SMatthew G. Knepley ++usedFaces; 2238856ac710SMatthew G. Knepley } 2239856ac710SMatthew G. Knepley } 2240856ac710SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 2241856ac710SMatthew G. Knepley PetscFunctionReturn(0); 2242856ac710SMatthew G. Knepley } 2243856ac710SMatthew G. Knepley 2244b81db932SToby Isaac static PetscErrorCode BuildGradientReconstruction_Internal_Tree(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2245b81db932SToby Isaac { 2246b81db932SToby Isaac DMLabel ghostLabel; 2247b81db932SToby Isaac PetscScalar *dx, *grad, **gref; 2248b81db932SToby Isaac PetscInt dim, cStart, cEnd, c, cEndInterior, fStart, fEnd, f, nStart, nEnd, maxNumFaces = 0; 2249b81db932SToby Isaac PetscSection neighSec; 2250b81db932SToby Isaac PetscInt (*neighbors)[2]; 2251b81db932SToby Isaac PetscInt *counter; 2252b81db932SToby Isaac PetscErrorCode ierr; 2253b81db932SToby Isaac 2254b81db932SToby Isaac PetscFunctionBegin; 2255b81db932SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2256b81db932SToby Isaac ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2257b81db932SToby Isaac ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 22585bc680faSToby Isaac if (cEndInterior < 0) { 22595bc680faSToby Isaac cEndInterior = cEnd; 22605bc680faSToby Isaac } 2261b81db932SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&neighSec);CHKERRQ(ierr); 2262b81db932SToby Isaac ierr = PetscSectionSetChart(neighSec,cStart,cEndInterior);CHKERRQ(ierr); 2263b81db932SToby Isaac ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2264c58f1c22SToby Isaac ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 2265b81db932SToby Isaac for (f = fStart; f < fEnd; f++) { 2266b81db932SToby Isaac const PetscInt *fcells; 2267b81db932SToby Isaac PetscBool boundary; 22685bc680faSToby Isaac PetscInt ghost = -1; 2269b81db932SToby Isaac PetscInt numChildren, numCells, c; 2270b81db932SToby Isaac 227106348e87SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 2272a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 2273b81db932SToby Isaac ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 2274b81db932SToby Isaac if ((ghost >= 0) || boundary || numChildren) continue; 2275b81db932SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 227606348e87SToby Isaac if (numCells == 2) { 2277b81db932SToby Isaac ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 2278b81db932SToby Isaac for (c = 0; c < 2; c++) { 2279b81db932SToby Isaac PetscInt cell = fcells[c]; 2280b81db932SToby Isaac 2281e6885bbbSToby Isaac if (cell >= cStart && cell < cEndInterior) { 2282b81db932SToby Isaac ierr = PetscSectionAddDof(neighSec,cell,1);CHKERRQ(ierr); 2283b81db932SToby Isaac } 2284b81db932SToby Isaac } 2285b81db932SToby Isaac } 228606348e87SToby Isaac } 2287b81db932SToby Isaac ierr = PetscSectionSetUp(neighSec);CHKERRQ(ierr); 2288b81db932SToby Isaac ierr = PetscSectionGetMaxDof(neighSec,&maxNumFaces);CHKERRQ(ierr); 2289b81db932SToby Isaac ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 2290b81db932SToby Isaac nStart = 0; 2291b81db932SToby Isaac ierr = PetscSectionGetStorageSize(neighSec,&nEnd);CHKERRQ(ierr); 2292b81db932SToby Isaac ierr = PetscMalloc1((nEnd-nStart),&neighbors);CHKERRQ(ierr); 2293b81db932SToby Isaac ierr = PetscCalloc1((cEndInterior-cStart),&counter);CHKERRQ(ierr); 2294b81db932SToby Isaac for (f = fStart; f < fEnd; f++) { 2295b81db932SToby Isaac const PetscInt *fcells; 2296b81db932SToby Isaac PetscBool boundary; 22975bc680faSToby Isaac PetscInt ghost = -1; 2298b81db932SToby Isaac PetscInt numChildren, numCells, c; 2299b81db932SToby Isaac 230006348e87SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);} 2301a6ba4734SToby Isaac ierr = DMIsBoundaryPoint(dm, f, &boundary);CHKERRQ(ierr); 2302b81db932SToby Isaac ierr = DMPlexGetTreeChildren(dm, f, &numChildren, NULL);CHKERRQ(ierr); 2303b81db932SToby Isaac if ((ghost >= 0) || boundary || numChildren) continue; 2304b81db932SToby Isaac ierr = DMPlexGetSupportSize(dm, f, &numCells);CHKERRQ(ierr); 230506348e87SToby Isaac if (numCells == 2) { 2306b81db932SToby Isaac ierr = DMPlexGetSupport(dm, f, &fcells);CHKERRQ(ierr); 2307b81db932SToby Isaac for (c = 0; c < 2; c++) { 2308b81db932SToby Isaac PetscInt cell = fcells[c], off; 2309b81db932SToby Isaac 2310e6885bbbSToby Isaac if (cell >= cStart && cell < cEndInterior) { 2311b81db932SToby Isaac ierr = PetscSectionGetOffset(neighSec,cell,&off);CHKERRQ(ierr); 2312b81db932SToby Isaac off += counter[cell - cStart]++; 2313b81db932SToby Isaac neighbors[off][0] = f; 2314b81db932SToby Isaac neighbors[off][1] = fcells[1 - c]; 2315b81db932SToby Isaac } 2316b81db932SToby Isaac } 2317b81db932SToby Isaac } 231806348e87SToby Isaac } 2319b81db932SToby Isaac ierr = PetscFree(counter);CHKERRQ(ierr); 2320b81db932SToby Isaac ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 2321b81db932SToby Isaac for (c = cStart; c < cEndInterior; c++) { 2322317218b9SToby Isaac PetscInt numFaces, f, d, off, ghost = -1; 2323640bce14SSatish Balay PetscFVCellGeom *cg; 2324b81db932SToby Isaac 2325b81db932SToby Isaac ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 2326b81db932SToby Isaac ierr = PetscSectionGetDof(neighSec, c, &numFaces);CHKERRQ(ierr); 2327b81db932SToby Isaac ierr = PetscSectionGetOffset(neighSec, c, &off);CHKERRQ(ierr); 2328317218b9SToby Isaac if (ghostLabel) {ierr = DMLabelGetValue(ghostLabel, c, &ghost);CHKERRQ(ierr);} 2329317218b9SToby Isaac if (ghost < 0 && numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 2330b81db932SToby Isaac for (f = 0; f < numFaces; ++f) { 2331640bce14SSatish Balay PetscFVCellGeom *cg1; 2332b81db932SToby Isaac PetscFVFaceGeom *fg; 2333b81db932SToby Isaac const PetscInt *fcells; 2334b81db932SToby Isaac PetscInt ncell, side, nface; 2335b81db932SToby Isaac 2336b81db932SToby Isaac nface = neighbors[off + f][0]; 2337b81db932SToby Isaac ncell = neighbors[off + f][1]; 2338b81db932SToby Isaac ierr = DMPlexGetSupport(dm,nface,&fcells);CHKERRQ(ierr); 2339b81db932SToby Isaac side = (c != fcells[0]); 2340b81db932SToby Isaac ierr = DMPlexPointLocalRef(dmFace, nface, fgeom, &fg);CHKERRQ(ierr); 2341b81db932SToby Isaac ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 2342b81db932SToby Isaac for (d = 0; d < dim; ++d) dx[f*dim+d] = cg1->centroid[d] - cg->centroid[d]; 2343b81db932SToby Isaac gref[f] = fg->grad[side]; /* Gradient reconstruction term will go here */ 2344b81db932SToby Isaac } 2345b81db932SToby Isaac ierr = PetscFVComputeGradient(fvm, numFaces, dx, grad);CHKERRQ(ierr); 2346b81db932SToby Isaac for (f = 0; f < numFaces; ++f) { 2347b81db932SToby Isaac for (d = 0; d < dim; ++d) gref[f][d] = grad[f*dim+d]; 2348b81db932SToby Isaac } 2349b81db932SToby Isaac } 2350b81db932SToby Isaac ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 23515fe94518SToby Isaac ierr = PetscSectionDestroy(&neighSec);CHKERRQ(ierr); 2352b81db932SToby Isaac ierr = PetscFree(neighbors);CHKERRQ(ierr); 2353b81db932SToby Isaac PetscFunctionReturn(0); 2354b81db932SToby Isaac } 2355b81db932SToby Isaac 2356856ac710SMatthew G. Knepley /*@ 2357856ac710SMatthew G. Knepley DMPlexComputeGradientFVM - Compute geometric factors for gradient reconstruction, which are stored in the geometry data, and compute layout for gradient data 2358856ac710SMatthew G. Knepley 2359856ac710SMatthew G. Knepley Collective on DM 2360856ac710SMatthew G. Knepley 2361856ac710SMatthew G. Knepley Input Arguments: 2362856ac710SMatthew G. Knepley + dm - The DM 2363856ac710SMatthew G. Knepley . fvm - The PetscFV 23648f9f38e3SMatthew G. Knepley . faceGeometry - The face geometry from DMPlexComputeFaceGeometryFVM() 23658f9f38e3SMatthew G. Knepley - cellGeometry - The face geometry from DMPlexComputeCellGeometryFVM() 2366856ac710SMatthew G. Knepley 2367856ac710SMatthew G. Knepley Output Parameters: 2368856ac710SMatthew G. Knepley + faceGeometry - The geometric factors for gradient calculation are inserted 2369856ac710SMatthew G. Knepley - dmGrad - The DM describing the layout of gradient data 2370856ac710SMatthew G. Knepley 2371856ac710SMatthew G. Knepley Level: developer 2372856ac710SMatthew G. Knepley 2373856ac710SMatthew G. Knepley .seealso: DMPlexGetFaceGeometryFVM(), DMPlexGetCellGeometryFVM() 2374856ac710SMatthew G. Knepley @*/ 2375856ac710SMatthew G. Knepley PetscErrorCode DMPlexComputeGradientFVM(DM dm, PetscFV fvm, Vec faceGeometry, Vec cellGeometry, DM *dmGrad) 2376856ac710SMatthew G. Knepley { 2377856ac710SMatthew G. Knepley DM dmFace, dmCell; 2378856ac710SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 2379b81db932SToby Isaac PetscSection sectionGrad, parentSection; 2380856ac710SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 2381856ac710SMatthew G. Knepley PetscErrorCode ierr; 2382856ac710SMatthew G. Knepley 2383856ac710SMatthew G. Knepley PetscFunctionBegin; 2384856ac710SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2385856ac710SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 2386856ac710SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2387856ac710SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 2388856ac710SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */ 2389856ac710SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 2390856ac710SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 2391856ac710SMatthew G. Knepley ierr = VecGetArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2392856ac710SMatthew G. Knepley ierr = VecGetArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2393b81db932SToby Isaac ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2394b81db932SToby Isaac if (!parentSection) { 2395856ac710SMatthew G. Knepley ierr = BuildGradientReconstruction_Internal(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2396b5a3613cSMatthew G. Knepley } else { 2397b81db932SToby Isaac ierr = BuildGradientReconstruction_Internal_Tree(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 2398b81db932SToby Isaac } 2399856ac710SMatthew G. Knepley ierr = VecRestoreArray(faceGeometry, &fgeom);CHKERRQ(ierr); 2400856ac710SMatthew G. Knepley ierr = VecRestoreArray(cellGeometry, &cgeom);CHKERRQ(ierr); 2401856ac710SMatthew G. Knepley /* Create storage for gradients */ 2402856ac710SMatthew G. Knepley ierr = DMClone(dm, dmGrad);CHKERRQ(ierr); 2403856ac710SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 2404856ac710SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 2405856ac710SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 2406856ac710SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 2407856ac710SMatthew G. Knepley ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr); 2408856ac710SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 2409856ac710SMatthew G. Knepley PetscFunctionReturn(0); 2410856ac710SMatthew G. Knepley } 2411b27d5b9eSToby Isaac 2412c501906fSMatthew G. Knepley /*@ 2413c501906fSMatthew G. Knepley DMPlexGetDataFVM - Retrieve precomputed cell geometry 2414c501906fSMatthew G. Knepley 2415c501906fSMatthew G. Knepley Collective on DM 2416c501906fSMatthew G. Knepley 2417c501906fSMatthew G. Knepley Input Arguments: 2418c501906fSMatthew G. Knepley + dm - The DM 2419c501906fSMatthew G. Knepley - fvm - The PetscFV 2420c501906fSMatthew G. Knepley 2421c501906fSMatthew G. Knepley Output Parameters: 2422c501906fSMatthew G. Knepley + cellGeometry - The cell geometry 2423c501906fSMatthew G. Knepley . faceGeometry - The face geometry 2424c501906fSMatthew G. Knepley - dmGrad - The gradient matrices 2425c501906fSMatthew G. Knepley 2426c501906fSMatthew G. Knepley Level: developer 2427c501906fSMatthew G. Knepley 2428c501906fSMatthew G. Knepley .seealso: DMPlexComputeGeometryFVM() 2429c501906fSMatthew G. Knepley @*/ 2430b27d5b9eSToby Isaac PetscErrorCode DMPlexGetDataFVM(DM dm, PetscFV fv, Vec *cellgeom, Vec *facegeom, DM *gradDM) 2431b27d5b9eSToby Isaac { 2432b27d5b9eSToby Isaac PetscObject cellgeomobj, facegeomobj; 2433b27d5b9eSToby Isaac PetscErrorCode ierr; 2434b27d5b9eSToby Isaac 2435b27d5b9eSToby Isaac PetscFunctionBegin; 2436b27d5b9eSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2437b27d5b9eSToby Isaac if (!cellgeomobj) { 2438b27d5b9eSToby Isaac Vec cellgeomInt, facegeomInt; 2439b27d5b9eSToby Isaac 2440b27d5b9eSToby Isaac ierr = DMPlexComputeGeometryFVM(dm, &cellgeomInt, &facegeomInt);CHKERRQ(ierr); 2441b27d5b9eSToby Isaac ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_cellgeom_fvm",(PetscObject)cellgeomInt);CHKERRQ(ierr); 2442b27d5b9eSToby Isaac ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_facegeom_fvm",(PetscObject)facegeomInt);CHKERRQ(ierr); 2443b27d5b9eSToby Isaac ierr = VecDestroy(&cellgeomInt);CHKERRQ(ierr); 2444b27d5b9eSToby Isaac ierr = VecDestroy(&facegeomInt);CHKERRQ(ierr); 2445b27d5b9eSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_cellgeom_fvm", &cellgeomobj);CHKERRQ(ierr); 2446b27d5b9eSToby Isaac } 2447b27d5b9eSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_facegeom_fvm", &facegeomobj);CHKERRQ(ierr); 2448b27d5b9eSToby Isaac if (cellgeom) *cellgeom = (Vec) cellgeomobj; 2449b27d5b9eSToby Isaac if (facegeom) *facegeom = (Vec) facegeomobj; 2450b27d5b9eSToby Isaac if (gradDM) { 2451b27d5b9eSToby Isaac PetscObject gradobj; 2452b27d5b9eSToby Isaac PetscBool computeGradients; 2453b27d5b9eSToby Isaac 2454b27d5b9eSToby Isaac ierr = PetscFVGetComputeGradients(fv,&computeGradients);CHKERRQ(ierr); 2455b27d5b9eSToby Isaac if (!computeGradients) { 2456b27d5b9eSToby Isaac *gradDM = NULL; 2457b27d5b9eSToby Isaac PetscFunctionReturn(0); 2458b27d5b9eSToby Isaac } 2459b27d5b9eSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2460b27d5b9eSToby Isaac if (!gradobj) { 2461b27d5b9eSToby Isaac DM dmGradInt; 2462b27d5b9eSToby Isaac 2463b27d5b9eSToby Isaac ierr = DMPlexComputeGradientFVM(dm,fv,(Vec) facegeomobj,(Vec) cellgeomobj,&dmGradInt);CHKERRQ(ierr); 2464b27d5b9eSToby Isaac ierr = PetscObjectCompose((PetscObject) dm, "DMPlex_dmgrad_fvm", (PetscObject)dmGradInt);CHKERRQ(ierr); 2465b27d5b9eSToby Isaac ierr = DMDestroy(&dmGradInt);CHKERRQ(ierr); 2466b27d5b9eSToby Isaac ierr = PetscObjectQuery((PetscObject) dm, "DMPlex_dmgrad_fvm", &gradobj);CHKERRQ(ierr); 2467b27d5b9eSToby Isaac } 2468b27d5b9eSToby Isaac *gradDM = (DM) gradobj; 2469b27d5b9eSToby Isaac } 2470b27d5b9eSToby Isaac PetscFunctionReturn(0); 2471b27d5b9eSToby Isaac } 2472d6143a4eSToby Isaac 24739d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_NewtonUpdate(PetscInt dimC, PetscInt dimR, PetscScalar *J, PetscScalar *invJ, PetscScalar *work, PetscReal *resNeg, PetscReal *guess) 24749d150b73SToby Isaac { 24759d150b73SToby Isaac PetscInt l, m; 24769d150b73SToby Isaac 2477cd345991SToby Isaac PetscFunctionBeginHot; 24789d150b73SToby Isaac if (dimC == dimR && dimR <= 3) { 24799d150b73SToby Isaac /* invert Jacobian, multiply */ 24809d150b73SToby Isaac PetscScalar det, idet; 24819d150b73SToby Isaac 24829d150b73SToby Isaac switch (dimR) { 24839d150b73SToby Isaac case 1: 24849d150b73SToby Isaac invJ[0] = 1./ J[0]; 24859d150b73SToby Isaac break; 24869d150b73SToby Isaac case 2: 24879d150b73SToby Isaac det = J[0] * J[3] - J[1] * J[2]; 24889d150b73SToby Isaac idet = 1./det; 24899d150b73SToby Isaac invJ[0] = J[3] * idet; 24909d150b73SToby Isaac invJ[1] = -J[1] * idet; 24919d150b73SToby Isaac invJ[2] = -J[2] * idet; 24929d150b73SToby Isaac invJ[3] = J[0] * idet; 24939d150b73SToby Isaac break; 24949d150b73SToby Isaac case 3: 24959d150b73SToby Isaac { 24969d150b73SToby Isaac invJ[0] = J[4] * J[8] - J[5] * J[7]; 24979d150b73SToby Isaac invJ[1] = J[2] * J[7] - J[1] * J[8]; 24989d150b73SToby Isaac invJ[2] = J[1] * J[5] - J[2] * J[4]; 24999d150b73SToby Isaac det = invJ[0] * J[0] + invJ[1] * J[3] + invJ[2] * J[6]; 25009d150b73SToby Isaac idet = 1./det; 25019d150b73SToby Isaac invJ[0] *= idet; 25029d150b73SToby Isaac invJ[1] *= idet; 25039d150b73SToby Isaac invJ[2] *= idet; 25049d150b73SToby Isaac invJ[3] = idet * (J[5] * J[6] - J[3] * J[8]); 25059d150b73SToby Isaac invJ[4] = idet * (J[0] * J[8] - J[2] * J[6]); 25069d150b73SToby Isaac invJ[5] = idet * (J[2] * J[3] - J[0] * J[5]); 25079d150b73SToby Isaac invJ[6] = idet * (J[3] * J[7] - J[4] * J[6]); 25089d150b73SToby Isaac invJ[7] = idet * (J[1] * J[6] - J[0] * J[7]); 25099d150b73SToby Isaac invJ[8] = idet * (J[0] * J[4] - J[1] * J[3]); 25109d150b73SToby Isaac } 25119d150b73SToby Isaac break; 25129d150b73SToby Isaac } 25139d150b73SToby Isaac for (l = 0; l < dimR; l++) { 25149d150b73SToby Isaac for (m = 0; m < dimC; m++) { 2515c6e120d1SToby Isaac guess[l] += PetscRealPart(invJ[l * dimC + m]) * resNeg[m]; 25169d150b73SToby Isaac } 25179d150b73SToby Isaac } 25189d150b73SToby Isaac } else { 25199d150b73SToby Isaac #if defined(PETSC_USE_COMPLEX) 25209d150b73SToby Isaac char transpose = 'C'; 25219d150b73SToby Isaac #else 25229d150b73SToby Isaac char transpose = 'T'; 25239d150b73SToby Isaac #endif 25249d150b73SToby Isaac PetscBLASInt m = dimR; 25259d150b73SToby Isaac PetscBLASInt n = dimC; 25269d150b73SToby Isaac PetscBLASInt one = 1; 25279d150b73SToby Isaac PetscBLASInt worksize = dimR * dimC, info; 25289d150b73SToby Isaac 25299d150b73SToby Isaac for (l = 0; l < dimC; l++) {invJ[l] = resNeg[l];} 25309d150b73SToby Isaac 25319d150b73SToby Isaac PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&m,&n,&one,J,&m,invJ,&n,work,&worksize, &info)); 25329d150b73SToby Isaac if (info != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS"); 25339d150b73SToby Isaac 2534c6e120d1SToby Isaac for (l = 0; l < dimR; l++) {guess[l] += PetscRealPart(invJ[l]);} 25359d150b73SToby Isaac } 25369d150b73SToby Isaac PetscFunctionReturn(0); 25379d150b73SToby Isaac } 25389d150b73SToby Isaac 25399d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 25409d150b73SToby Isaac { 2541c0cbe899SToby Isaac PetscInt coordSize, i, j, k, l, m, maxIts = 7, numV = (1 << dimR); 25429d150b73SToby Isaac PetscScalar *coordsScalar = NULL; 25439d150b73SToby Isaac PetscReal *cellData, *cellCoords, *cellCoeffs, *extJ, *resNeg; 25449d150b73SToby Isaac PetscScalar *J, *invJ, *work; 25459d150b73SToby Isaac PetscErrorCode ierr; 25469d150b73SToby Isaac 25479d150b73SToby Isaac PetscFunctionBegin; 25489d150b73SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 25499d150b73SToby Isaac ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 25509d150b73SToby Isaac if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr); 25519d150b73SToby Isaac ierr = DMGetWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr); 25529d150b73SToby Isaac ierr = DMGetWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr); 25539d150b73SToby Isaac cellCoords = &cellData[0]; 25549d150b73SToby Isaac cellCoeffs = &cellData[coordSize]; 25559d150b73SToby Isaac extJ = &cellData[2 * coordSize]; 25569d150b73SToby Isaac resNeg = &cellData[2 * coordSize + dimR]; 25579d150b73SToby Isaac invJ = &J[dimR * dimC]; 25589d150b73SToby Isaac work = &J[2 * dimR * dimC]; 25599d150b73SToby Isaac if (dimR == 2) { 25609d150b73SToby Isaac const PetscInt zToPlex[4] = {0, 1, 3, 2}; 25619d150b73SToby Isaac 25629d150b73SToby Isaac for (i = 0; i < 4; i++) { 25639d150b73SToby Isaac PetscInt plexI = zToPlex[i]; 25649d150b73SToby Isaac 25659d150b73SToby Isaac for (j = 0; j < dimC; j++) { 25669d150b73SToby Isaac cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 25679d150b73SToby Isaac } 25689d150b73SToby Isaac } 25699d150b73SToby Isaac } else if (dimR == 3) { 25709d150b73SToby Isaac const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 25719d150b73SToby Isaac 25729d150b73SToby Isaac for (i = 0; i < 8; i++) { 25739d150b73SToby Isaac PetscInt plexI = zToPlex[i]; 25749d150b73SToby Isaac 25759d150b73SToby Isaac for (j = 0; j < dimC; j++) { 25769d150b73SToby Isaac cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 25779d150b73SToby Isaac } 25789d150b73SToby Isaac } 25799d150b73SToby Isaac } else { 25809d150b73SToby Isaac for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);} 25819d150b73SToby Isaac } 25829d150b73SToby Isaac /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 25839d150b73SToby Isaac for (i = 0; i < dimR; i++) { 25849d150b73SToby Isaac PetscReal *swap; 25859d150b73SToby Isaac 25869d150b73SToby Isaac for (j = 0; j < (numV / 2); j++) { 25879d150b73SToby Isaac for (k = 0; k < dimC; k++) { 25889d150b73SToby Isaac cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 25899d150b73SToby Isaac cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 25909d150b73SToby Isaac } 25919d150b73SToby Isaac } 25929d150b73SToby Isaac 25939d150b73SToby Isaac if (i < dimR - 1) { 25949d150b73SToby Isaac swap = cellCoeffs; 25959d150b73SToby Isaac cellCoeffs = cellCoords; 25969d150b73SToby Isaac cellCoords = swap; 25979d150b73SToby Isaac } 25989d150b73SToby Isaac } 25999d150b73SToby Isaac ierr = PetscMemzero(refCoords,numPoints * dimR * sizeof (PetscReal));CHKERRQ(ierr); 26009d150b73SToby Isaac for (j = 0; j < numPoints; j++) { 26019d150b73SToby Isaac for (i = 0; i < maxIts; i++) { 26029d150b73SToby Isaac PetscReal *guess = &refCoords[dimR * j]; 26039d150b73SToby Isaac 26049d150b73SToby Isaac /* compute -residual and Jacobian */ 26059d150b73SToby Isaac for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[dimC * j + k];} 26069d150b73SToby Isaac for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;} 26079d150b73SToby Isaac for (k = 0; k < numV; k++) { 26089d150b73SToby Isaac PetscReal extCoord = 1.; 26099d150b73SToby Isaac for (l = 0; l < dimR; l++) { 26109d150b73SToby Isaac PetscReal coord = guess[l]; 26119d150b73SToby Isaac PetscInt dep = (k & (1 << l)) >> l; 26129d150b73SToby Isaac 26139d150b73SToby Isaac extCoord *= dep * coord + !dep; 26149d150b73SToby Isaac extJ[l] = dep; 26159d150b73SToby Isaac 26169d150b73SToby Isaac for (m = 0; m < dimR; m++) { 26179d150b73SToby Isaac PetscReal coord = guess[m]; 26189d150b73SToby Isaac PetscInt dep = ((k & (1 << m)) >> m) && (m != l); 26199d150b73SToby Isaac PetscReal mult = dep * coord + !dep; 26209d150b73SToby Isaac 26219d150b73SToby Isaac extJ[l] *= mult; 26229d150b73SToby Isaac } 26239d150b73SToby Isaac } 26249d150b73SToby Isaac for (l = 0; l < dimC; l++) { 26259d150b73SToby Isaac PetscReal coeff = cellCoeffs[dimC * k + l]; 26269d150b73SToby Isaac 26279d150b73SToby Isaac resNeg[l] -= coeff * extCoord; 26289d150b73SToby Isaac for (m = 0; m < dimR; m++) { 26299d150b73SToby Isaac J[dimR * l + m] += coeff * extJ[m]; 26309d150b73SToby Isaac } 26319d150b73SToby Isaac } 26329d150b73SToby Isaac } 26330611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG) 26340611203eSToby Isaac { 26350611203eSToby Isaac PetscReal maxAbs = 0.; 26360611203eSToby Isaac 26370611203eSToby Isaac for (l = 0; l < dimC; l++) { 26380611203eSToby Isaac maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l])); 26390611203eSToby Isaac } 26400611203eSToby Isaac ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr); 26410611203eSToby Isaac } 26420611203eSToby Isaac #endif 26439d150b73SToby Isaac 26449d150b73SToby Isaac ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr); 26459d150b73SToby Isaac } 26469d150b73SToby Isaac } 26479d150b73SToby Isaac ierr = DMRestoreWorkArray(dm, 3 * dimR * dimC, PETSC_SCALAR, &J);CHKERRQ(ierr); 26489d150b73SToby Isaac ierr = DMRestoreWorkArray(dm, 2 * coordSize + dimR + dimC, PETSC_REAL, &cellData);CHKERRQ(ierr); 26499d150b73SToby Isaac ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 26509d150b73SToby Isaac PetscFunctionReturn(0); 26519d150b73SToby Isaac } 26529d150b73SToby Isaac 26539d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_Tensor(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 26549d150b73SToby Isaac { 26559d150b73SToby Isaac PetscInt coordSize, i, j, k, l, numV = (1 << dimR); 26569d150b73SToby Isaac PetscScalar *coordsScalar = NULL; 26579d150b73SToby Isaac PetscReal *cellData, *cellCoords, *cellCoeffs; 26589d150b73SToby Isaac PetscErrorCode ierr; 26599d150b73SToby Isaac 26609d150b73SToby Isaac PetscFunctionBegin; 26619d150b73SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 26629d150b73SToby Isaac ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 26639d150b73SToby Isaac if (coordSize < dimC * numV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expecting at least %D coordinates, got %D",dimC * (1 << dimR), coordSize);CHKERRQ(ierr); 26649d150b73SToby Isaac ierr = DMGetWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr); 26659d150b73SToby Isaac cellCoords = &cellData[0]; 26669d150b73SToby Isaac cellCoeffs = &cellData[coordSize]; 26679d150b73SToby Isaac if (dimR == 2) { 26689d150b73SToby Isaac const PetscInt zToPlex[4] = {0, 1, 3, 2}; 26699d150b73SToby Isaac 26709d150b73SToby Isaac for (i = 0; i < 4; i++) { 26719d150b73SToby Isaac PetscInt plexI = zToPlex[i]; 26729d150b73SToby Isaac 26739d150b73SToby Isaac for (j = 0; j < dimC; j++) { 26749d150b73SToby Isaac cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 26759d150b73SToby Isaac } 26769d150b73SToby Isaac } 26779d150b73SToby Isaac } else if (dimR == 3) { 26789d150b73SToby Isaac const PetscInt zToPlex[8] = {0, 3, 1, 2, 4, 5, 7, 6}; 26799d150b73SToby Isaac 26809d150b73SToby Isaac for (i = 0; i < 8; i++) { 26819d150b73SToby Isaac PetscInt plexI = zToPlex[i]; 26829d150b73SToby Isaac 26839d150b73SToby Isaac for (j = 0; j < dimC; j++) { 26849d150b73SToby Isaac cellCoords[dimC * i + j] = PetscRealPart(coordsScalar[dimC * plexI + j]); 26859d150b73SToby Isaac } 26869d150b73SToby Isaac } 26879d150b73SToby Isaac } else { 26889d150b73SToby Isaac for (i = 0; i < coordSize; i++) {cellCoords[i] = PetscRealPart(coordsScalar[i]);} 26899d150b73SToby Isaac } 26909d150b73SToby Isaac /* Perform the shuffling transform that converts values at the corners of [-1,1]^d to coefficients */ 26919d150b73SToby Isaac for (i = 0; i < dimR; i++) { 26929d150b73SToby Isaac PetscReal *swap; 26939d150b73SToby Isaac 26949d150b73SToby Isaac for (j = 0; j < (numV / 2); j++) { 26959d150b73SToby Isaac for (k = 0; k < dimC; k++) { 26969d150b73SToby Isaac cellCoeffs[dimC * j + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] + cellCoords[dimC * 2 * j + k]); 26979d150b73SToby Isaac cellCoeffs[dimC * (j + (numV / 2)) + k] = 0.5 * (cellCoords[dimC * (2 * j + 1) + k] - cellCoords[dimC * 2 * j + k]); 26989d150b73SToby Isaac } 26999d150b73SToby Isaac } 27009d150b73SToby Isaac 27019d150b73SToby Isaac if (i < dimR - 1) { 27029d150b73SToby Isaac swap = cellCoeffs; 27039d150b73SToby Isaac cellCoeffs = cellCoords; 27049d150b73SToby Isaac cellCoords = swap; 27059d150b73SToby Isaac } 27069d150b73SToby Isaac } 27079d150b73SToby Isaac ierr = PetscMemzero(realCoords,numPoints * dimC * sizeof (PetscReal));CHKERRQ(ierr); 27089d150b73SToby Isaac for (j = 0; j < numPoints; j++) { 27099d150b73SToby Isaac const PetscReal *guess = &refCoords[dimR * j]; 27109d150b73SToby Isaac PetscReal *mapped = &realCoords[dimC * j]; 27119d150b73SToby Isaac 27129d150b73SToby Isaac for (k = 0; k < numV; k++) { 27139d150b73SToby Isaac PetscReal extCoord = 1.; 27149d150b73SToby Isaac for (l = 0; l < dimR; l++) { 27159d150b73SToby Isaac PetscReal coord = guess[l]; 27169d150b73SToby Isaac PetscInt dep = (k & (1 << l)) >> l; 27179d150b73SToby Isaac 27189d150b73SToby Isaac extCoord *= dep * coord + !dep; 27199d150b73SToby Isaac } 27209d150b73SToby Isaac for (l = 0; l < dimC; l++) { 27219d150b73SToby Isaac PetscReal coeff = cellCoeffs[dimC * k + l]; 27229d150b73SToby Isaac 27239d150b73SToby Isaac mapped[l] += coeff * extCoord; 27249d150b73SToby Isaac } 27259d150b73SToby Isaac } 27269d150b73SToby Isaac } 27279d150b73SToby Isaac ierr = DMRestoreWorkArray(dm, 2 * coordSize, PETSC_REAL, &cellData);CHKERRQ(ierr); 27289d150b73SToby Isaac ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &coordsScalar);CHKERRQ(ierr); 27299d150b73SToby Isaac PetscFunctionReturn(0); 27309d150b73SToby Isaac } 27319d150b73SToby Isaac 27329d150b73SToby Isaac static PetscErrorCode DMPlexCoordinatesToReference_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 27339d150b73SToby Isaac { 2734c0cbe899SToby Isaac PetscInt numComp, numDof, i, j, k, l, m, maxIter = 7, coordSize; 2735c6e120d1SToby Isaac PetscScalar *nodes = NULL; 2736c6e120d1SToby Isaac PetscReal *invV, *modes; 2737c6e120d1SToby Isaac PetscReal *B, *D, *resNeg; 2738c6e120d1SToby Isaac PetscScalar *J, *invJ, *work; 27399d150b73SToby Isaac PetscErrorCode ierr; 27409d150b73SToby Isaac 27419d150b73SToby Isaac PetscFunctionBegin; 27429d150b73SToby Isaac ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr); 27439d150b73SToby Isaac ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 27449d150b73SToby Isaac if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC); 27459d150b73SToby Isaac ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 27469d150b73SToby Isaac /* convert nodes to values in the stable evaluation basis */ 2747c6e120d1SToby Isaac ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 27489d150b73SToby Isaac invV = fe->invV; 27499d150b73SToby Isaac for (i = 0; i < numDof; i++) { 27509d150b73SToby Isaac for (j = 0; j < dimC; j++) { 27519d150b73SToby Isaac modes[i * dimC + j] = 0.; 27529d150b73SToby Isaac for (k = 0; k < numDof; k++) { 2753c6e120d1SToby Isaac modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]); 27549d150b73SToby Isaac } 27559d150b73SToby Isaac } 27569d150b73SToby Isaac } 2757c6e120d1SToby Isaac ierr = DMGetWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr); 2758c6e120d1SToby Isaac D = &B[numDof]; 2759c6e120d1SToby Isaac resNeg = &D[numDof * dimR]; 2760c6e120d1SToby Isaac ierr = DMGetWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr); 27619d150b73SToby Isaac invJ = &J[dimC * dimR]; 27629d150b73SToby Isaac work = &invJ[dimC * dimR]; 27639d150b73SToby Isaac for (i = 0; i < numPoints * dimR; i++) {refCoords[i] = 0.;} 27649d150b73SToby Isaac for (j = 0; j < numPoints; j++) { 27659b1f03cbSToby Isaac for (i = 0; i < maxIter; i++) { /* we could batch this so that we're not making big B and D arrays all the time */ 27669d150b73SToby Isaac PetscReal *guess = &refCoords[j * dimR]; 27679d150b73SToby Isaac ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, D, NULL);CHKERRQ(ierr); 27689d150b73SToby Isaac for (k = 0; k < dimC; k++) {resNeg[k] = realCoords[j * dimC + k];} 27699d150b73SToby Isaac for (k = 0; k < dimC * dimR; k++) {J[k] = 0.;} 27709d150b73SToby Isaac for (k = 0; k < numDof; k++) { 27719d150b73SToby Isaac for (l = 0; l < dimC; l++) { 27729d150b73SToby Isaac resNeg[l] -= modes[k * dimC + l] * B[k]; 27739d150b73SToby Isaac for (m = 0; m < dimR; m++) { 27749d150b73SToby Isaac J[l * dimR + m] += modes[k * dimC + l] * D[k * dimR + m]; 27759d150b73SToby Isaac } 27769d150b73SToby Isaac } 27779d150b73SToby Isaac } 27780611203eSToby Isaac #if 0 && defined(PETSC_USE_DEBUG) 27790611203eSToby Isaac { 27800611203eSToby Isaac PetscReal maxAbs = 0.; 27810611203eSToby Isaac 27820611203eSToby Isaac for (l = 0; l < dimC; l++) { 27830611203eSToby Isaac maxAbs = PetscMax(maxAbs,PetscAbsReal(resNeg[l])); 27840611203eSToby Isaac } 27850611203eSToby Isaac ierr = PetscInfo4(dm,"cell %D, point %D, iter %D: res %g\n",cell,j,i,maxAbs);CHKERRQ(ierr); 27860611203eSToby Isaac } 27870611203eSToby Isaac #endif 27889d150b73SToby Isaac ierr = DMPlexCoordinatesToReference_NewtonUpdate(dimC,dimR,J,invJ,work,resNeg,guess);CHKERRQ(ierr); 27899d150b73SToby Isaac } 27909d150b73SToby Isaac } 2791c6e120d1SToby Isaac ierr = DMRestoreWorkArray(dm,3 * dimC * dimR,PETSC_SCALAR,&J);CHKERRQ(ierr); 2792c6e120d1SToby Isaac ierr = DMRestoreWorkArray(dm,numDof + numDof * dimR + dimC,PETSC_REAL,&B);CHKERRQ(ierr); 2793c6e120d1SToby Isaac ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 27949d150b73SToby Isaac ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 27959d150b73SToby Isaac PetscFunctionReturn(0); 27969d150b73SToby Isaac } 27979d150b73SToby Isaac 27989d150b73SToby Isaac static PetscErrorCode DMPlexReferenceToCoordinates_FE(DM dm, PetscFE fe, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[], Vec coords, PetscInt dimC, PetscInt dimR) 27999d150b73SToby Isaac { 28009d150b73SToby Isaac PetscInt numComp, numDof, i, j, k, l, coordSize; 2801c6e120d1SToby Isaac PetscScalar *nodes = NULL; 2802c6e120d1SToby Isaac PetscReal *invV, *modes; 28039d150b73SToby Isaac PetscReal *B; 28049d150b73SToby Isaac PetscErrorCode ierr; 28059d150b73SToby Isaac 28069d150b73SToby Isaac PetscFunctionBegin; 28079d150b73SToby Isaac ierr = PetscFEGetDimension(fe, &numDof);CHKERRQ(ierr); 28089d150b73SToby Isaac ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 28099d150b73SToby Isaac if (numComp != dimC) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"coordinate discretization must have as many components (%D) as embedding dimension (!= %D)",numComp,dimC); 28109d150b73SToby Isaac ierr = DMPlexVecGetClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 28119d150b73SToby Isaac /* convert nodes to values in the stable evaluation basis */ 2812c6e120d1SToby Isaac ierr = DMGetWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 28139d150b73SToby Isaac invV = fe->invV; 28149d150b73SToby Isaac for (i = 0; i < numDof; i++) { 28159d150b73SToby Isaac for (j = 0; j < dimC; j++) { 28169d150b73SToby Isaac modes[i * dimC + j] = 0.; 28179d150b73SToby Isaac for (k = 0; k < numDof; k++) { 2818c6e120d1SToby Isaac modes[i * dimC + j] += invV[i * numDof + k] * PetscRealPart(nodes[k * dimC + j]); 28199d150b73SToby Isaac } 28209d150b73SToby Isaac } 28219d150b73SToby Isaac } 28229b1f03cbSToby Isaac ierr = DMGetWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr); 28239d150b73SToby Isaac for (i = 0; i < numPoints * dimC; i++) {realCoords[i] = 0.;} 28249d150b73SToby Isaac for (j = 0; j < numPoints; j++) { 28259d150b73SToby Isaac const PetscReal *guess = &refCoords[j * dimR]; 28269d150b73SToby Isaac PetscReal *mapped = &realCoords[j * dimC]; 28279d150b73SToby Isaac 28289d150b73SToby Isaac ierr = PetscSpaceEvaluate(fe->basisSpace, 1, guess, B, NULL, NULL);CHKERRQ(ierr); 28299d150b73SToby Isaac for (k = 0; k < numDof; k++) { 28309d150b73SToby Isaac for (l = 0; l < dimC; l++) { 28319d150b73SToby Isaac mapped[l] += modes[k * dimC + l] * B[k]; 28329d150b73SToby Isaac } 28339d150b73SToby Isaac } 28349d150b73SToby Isaac } 28359b1f03cbSToby Isaac ierr = DMRestoreWorkArray(dm,numDof,PETSC_REAL,&B);CHKERRQ(ierr); 2836c6e120d1SToby Isaac ierr = DMRestoreWorkArray(dm,dimC * numDof,PETSC_REAL,&modes);CHKERRQ(ierr); 28379d150b73SToby Isaac ierr = DMPlexVecRestoreClosure(dm, NULL, coords, cell, &coordSize, &nodes);CHKERRQ(ierr); 28389d150b73SToby Isaac PetscFunctionReturn(0); 28399d150b73SToby Isaac } 28409d150b73SToby Isaac 2841d6143a4eSToby Isaac /*@ 2842d6143a4eSToby Isaac DMPlexCoordinatesToReference - Pull coordinates back from the mesh to the reference element using a single element 2843d6143a4eSToby Isaac map. This inversion will be accurate inside the reference element, but may be inaccurate for mappings that do not 2844d6143a4eSToby Isaac extend uniquely outside the reference cell (e.g, most non-affine maps) 2845d6143a4eSToby Isaac 2846d6143a4eSToby Isaac Not collective 2847d6143a4eSToby Isaac 2848d6143a4eSToby Isaac Input Parameters: 2849d6143a4eSToby Isaac + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or 2850d6143a4eSToby Isaac implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 2851d6143a4eSToby Isaac as a multilinear map for tensor-product elements 2852d6143a4eSToby Isaac . cell - the cell whose map is used. 2853d6143a4eSToby Isaac . numPoints - the number of points to locate 28541b266c99SBarry Smith - realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim()) 2855d6143a4eSToby Isaac 2856d6143a4eSToby Isaac Output Parameters: 2857d6143a4eSToby Isaac . refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension()) 28581b266c99SBarry Smith 28591b266c99SBarry Smith Level: intermediate 2860d6143a4eSToby Isaac @*/ 2861d6143a4eSToby Isaac PetscErrorCode DMPlexCoordinatesToReference(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal realCoords[], PetscReal refCoords[]) 2862d6143a4eSToby Isaac { 28639d150b73SToby Isaac PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i; 28649d150b73SToby Isaac DM coordDM = NULL; 28659d150b73SToby Isaac Vec coords; 28669d150b73SToby Isaac PetscFE fe = NULL; 28679d150b73SToby Isaac PetscErrorCode ierr; 28689d150b73SToby Isaac 2869d6143a4eSToby Isaac PetscFunctionBegin; 28709d150b73SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 28719d150b73SToby Isaac ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr); 28729d150b73SToby Isaac ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr); 28739d150b73SToby Isaac if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0); 28749d150b73SToby Isaac ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 28759d150b73SToby Isaac ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 28769d150b73SToby Isaac ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 28779d150b73SToby Isaac if (coordDM) { 28789d150b73SToby Isaac PetscInt coordFields; 28799d150b73SToby Isaac 28809d150b73SToby Isaac ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr); 28819d150b73SToby Isaac if (coordFields) { 28829d150b73SToby Isaac PetscClassId id; 28839d150b73SToby Isaac PetscObject disc; 28849d150b73SToby Isaac 28859d150b73SToby Isaac ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr); 28869d150b73SToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 28879d150b73SToby Isaac if (id == PETSCFE_CLASSID) { 28889d150b73SToby Isaac fe = (PetscFE) disc; 28899d150b73SToby Isaac } 28909d150b73SToby Isaac } 28919d150b73SToby Isaac } 28929d150b73SToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 28939d150b73SToby Isaac ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr); 28949d150b73SToby Isaac cEnd = cEndInterior > 0 ? cEndInterior : cEnd; 28959d150b73SToby Isaac if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr); 28969d150b73SToby Isaac if (!fe) { /* implicit discretization: affine or multilinear */ 28979d150b73SToby Isaac PetscInt coneSize; 28989d150b73SToby Isaac PetscBool isSimplex, isTensor; 28999d150b73SToby Isaac 29009d150b73SToby Isaac ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 29019d150b73SToby Isaac isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 29029d150b73SToby Isaac isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 29039d150b73SToby Isaac if (isSimplex) { 29049d150b73SToby Isaac PetscReal detJ, *v0, *J, *invJ; 29059d150b73SToby Isaac 29069d150b73SToby Isaac ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 29079d150b73SToby Isaac J = &v0[dimC]; 29089d150b73SToby Isaac invJ = &J[dimC * dimC]; 29099d150b73SToby Isaac ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, invJ, &detJ);CHKERRQ(ierr); 29109d150b73SToby Isaac for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */ 29119d150b73SToby Isaac CoordinatesRealToRef(dimC, dimR, v0, invJ, &realCoords[dimC * i], &refCoords[dimR * i]); 29129d150b73SToby Isaac } 29139d150b73SToby Isaac ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 29149d150b73SToby Isaac } else if (isTensor) { 29159d150b73SToby Isaac ierr = DMPlexCoordinatesToReference_Tensor(coordDM, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr); 29169d150b73SToby Isaac } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize); 29179d150b73SToby Isaac } else { 29189d150b73SToby Isaac ierr = DMPlexCoordinatesToReference_FE(coordDM, fe, cell, numPoints, realCoords, refCoords, coords, dimC, dimR);CHKERRQ(ierr); 29199d150b73SToby Isaac } 29209d150b73SToby Isaac PetscFunctionReturn(0); 29219d150b73SToby Isaac } 29229d150b73SToby Isaac 29239d150b73SToby Isaac /*@ 29249d150b73SToby Isaac DMPlexReferenceToCoordinates - Map references coordinates to coordinates in the the mesh for a single element map. 29259d150b73SToby Isaac 29269d150b73SToby Isaac Not collective 29279d150b73SToby Isaac 29289d150b73SToby Isaac Input Parameters: 29299d150b73SToby Isaac + dm - The mesh, with coordinate maps defined either by a PetscDS for the coordinate DM (see DMGetCoordinateDM()) or 29309d150b73SToby Isaac implicitly by the coordinates of the corner vertices of the cell: as an affine map for simplicial elements, or 29319d150b73SToby Isaac as a multilinear map for tensor-product elements 29329d150b73SToby Isaac . cell - the cell whose map is used. 29339d150b73SToby Isaac . numPoints - the number of points to locate 29349d150b73SToby Isaac + refCoords - (numPoints x dimension) array of reference coordinates (see DMGetDimension()) 29359d150b73SToby Isaac 29369d150b73SToby Isaac Output Parameters: 29379d150b73SToby Isaac . realCoords - (numPoints x coordinate dimension) array of coordinates (see DMGetCoordinateDim()) 29381b266c99SBarry Smith 29391b266c99SBarry Smith Level: intermediate 29409d150b73SToby Isaac @*/ 29419d150b73SToby Isaac PetscErrorCode DMPlexReferenceToCoordinates(DM dm, PetscInt cell, PetscInt numPoints, const PetscReal refCoords[], PetscReal realCoords[]) 29429d150b73SToby Isaac { 29439d150b73SToby Isaac PetscInt dimC, dimR, depth, cStart, cEnd, cEndInterior, i; 29449d150b73SToby Isaac DM coordDM = NULL; 29459d150b73SToby Isaac Vec coords; 29469d150b73SToby Isaac PetscFE fe = NULL; 29479d150b73SToby Isaac PetscErrorCode ierr; 29489d150b73SToby Isaac 29499d150b73SToby Isaac PetscFunctionBegin; 29509d150b73SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 29519d150b73SToby Isaac ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr); 29529d150b73SToby Isaac ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr); 29539d150b73SToby Isaac if (dimR <= 0 || dimC <= 0 || numPoints <= 0) PetscFunctionReturn(0); 29549d150b73SToby Isaac ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 29559d150b73SToby Isaac ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 29569d150b73SToby Isaac ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 29579d150b73SToby Isaac if (coordDM) { 29589d150b73SToby Isaac PetscInt coordFields; 29599d150b73SToby Isaac 29609d150b73SToby Isaac ierr = DMGetNumFields(coordDM,&coordFields);CHKERRQ(ierr); 29619d150b73SToby Isaac if (coordFields) { 29629d150b73SToby Isaac PetscClassId id; 29639d150b73SToby Isaac PetscObject disc; 29649d150b73SToby Isaac 29659d150b73SToby Isaac ierr = DMGetField(coordDM,0,&disc);CHKERRQ(ierr); 29669d150b73SToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 29679d150b73SToby Isaac if (id == PETSCFE_CLASSID) { 29689d150b73SToby Isaac fe = (PetscFE) disc; 29699d150b73SToby Isaac } 29709d150b73SToby Isaac } 29719d150b73SToby Isaac } 29729d150b73SToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 29739d150b73SToby Isaac ierr = DMPlexGetHybridBounds(dm,&cEndInterior,NULL,NULL,NULL);CHKERRQ(ierr); 29749d150b73SToby Isaac cEnd = cEndInterior > 0 ? cEndInterior : cEnd; 29759d150b73SToby Isaac if (cell < cStart || cell >= cEnd) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"point %D not in cell range [%D,%D)",cell,cStart,cEnd);CHKERRQ(ierr); 29769d150b73SToby Isaac if (!fe) { /* implicit discretization: affine or multilinear */ 29779d150b73SToby Isaac PetscInt coneSize; 29789d150b73SToby Isaac PetscBool isSimplex, isTensor; 29799d150b73SToby Isaac 29809d150b73SToby Isaac ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 29819d150b73SToby Isaac isSimplex = (coneSize == (dimR + 1)) ? PETSC_TRUE : PETSC_FALSE; 29829d150b73SToby Isaac isTensor = (coneSize == ((depth == 1) ? (1 << dimR) : (2 * dimR))) ? PETSC_TRUE : PETSC_FALSE; 29839d150b73SToby Isaac if (isSimplex) { 29849d150b73SToby Isaac PetscReal detJ, *v0, *J; 29859d150b73SToby Isaac 29869d150b73SToby Isaac ierr = DMGetWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 29879d150b73SToby Isaac J = &v0[dimC]; 29889d150b73SToby Isaac ierr = DMPlexComputeCellGeometryAffineFEM(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr); 29899d150b73SToby Isaac for (i = 0; i < numPoints; i++) { /* Apply the inverse affine transformation for each point */ 29909d150b73SToby Isaac CoordinatesRefToReal(dimC, dimR, v0, J, &refCoords[dimR * i], &realCoords[dimC * i]); 29919d150b73SToby Isaac } 29929d150b73SToby Isaac ierr = DMRestoreWorkArray(dm,dimC + 2 * dimC * dimC, PETSC_REAL, &v0);CHKERRQ(ierr); 29939d150b73SToby Isaac } else if (isTensor) { 29949d150b73SToby Isaac ierr = DMPlexReferenceToCoordinates_Tensor(coordDM, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr); 29959d150b73SToby Isaac } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unrecognized cone size %D",coneSize); 29969d150b73SToby Isaac } else { 29979d150b73SToby Isaac ierr = DMPlexReferenceToCoordinates_FE(coordDM, fe, cell, numPoints, refCoords, realCoords, coords, dimC, dimR);CHKERRQ(ierr); 29989d150b73SToby Isaac } 2999d6143a4eSToby Isaac PetscFunctionReturn(0); 3000d6143a4eSToby Isaac } 3001