157459e9aSMatthew G Knepley #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 257459e9aSMatthew G Knepley 357459e9aSMatthew G Knepley #undef __FUNCT__ 457459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureArray_Private" 557459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar **array) 657459e9aSMatthew G Knepley { 757459e9aSMatthew G Knepley PetscScalar *a; 857459e9aSMatthew G Knepley PetscInt size = 0, dof, off, d, k, i; 957459e9aSMatthew G Knepley PetscErrorCode ierr; 1057459e9aSMatthew G Knepley 1157459e9aSMatthew G Knepley PetscFunctionBegin; 1257459e9aSMatthew G Knepley for(i = 0; i < nP; ++i) { 1357459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 1457459e9aSMatthew G Knepley size += dof; 1557459e9aSMatthew G Knepley } 164c06c320SMatthew G Knepley ierr = DMGetWorkArray(dm, size, &a);CHKERRQ(ierr); 1757459e9aSMatthew G Knepley for(i = 0, k = 0; i < nP; ++i) { 1857459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 1957459e9aSMatthew G Knepley ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 2057459e9aSMatthew G Knepley 2157459e9aSMatthew G Knepley for(d = 0; d < dof; ++d, ++k) { 2257459e9aSMatthew G Knepley a[k] = vArray[off+d]; 2357459e9aSMatthew G Knepley } 2457459e9aSMatthew G Knepley } 2557459e9aSMatthew G Knepley *array = a; 2657459e9aSMatthew G Knepley PetscFunctionReturn(0); 2757459e9aSMatthew G Knepley } 2857459e9aSMatthew G Knepley 2957459e9aSMatthew G Knepley #undef __FUNCT__ 3057459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private" 3157459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode) 3257459e9aSMatthew G Knepley { 3357459e9aSMatthew G Knepley PetscInt dof, off, d, k, i; 3457459e9aSMatthew G Knepley PetscErrorCode ierr; 3557459e9aSMatthew G Knepley 3657459e9aSMatthew G Knepley PetscFunctionBegin; 3757459e9aSMatthew G Knepley if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) { 3857459e9aSMatthew G Knepley for(i = 0, k = 0; i < nP; ++i) { 3957459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 4057459e9aSMatthew G Knepley ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 4157459e9aSMatthew G Knepley 4257459e9aSMatthew G Knepley for(d = 0; d < dof; ++d, ++k) { 4357459e9aSMatthew G Knepley vArray[off+d] = array[k]; 4457459e9aSMatthew G Knepley } 4557459e9aSMatthew G Knepley } 4657459e9aSMatthew G Knepley } else { 4757459e9aSMatthew G Knepley for(i = 0, k = 0; i < nP; ++i) { 4857459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 4957459e9aSMatthew G Knepley ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 5057459e9aSMatthew G Knepley 5157459e9aSMatthew G Knepley for(d = 0; d < dof; ++d, ++k) { 5257459e9aSMatthew G Knepley vArray[off+d] += array[k]; 5357459e9aSMatthew G Knepley } 5457459e9aSMatthew G Knepley } 5557459e9aSMatthew G Knepley } 5657459e9aSMatthew G Knepley PetscFunctionReturn(0); 5757459e9aSMatthew G Knepley } 5857459e9aSMatthew G Knepley 5957459e9aSMatthew G Knepley #undef __FUNCT__ 6057459e9aSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure" 6157459e9aSMatthew G Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values) 6257459e9aSMatthew G Knepley { 6357459e9aSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 6457459e9aSMatthew G Knepley PetscInt dim = da->dim; 6557459e9aSMatthew G Knepley PetscScalar *vArray; 6657459e9aSMatthew G Knepley PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 6795babc02SJed Brown PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 6857459e9aSMatthew G Knepley PetscErrorCode ierr; 6957459e9aSMatthew G Knepley 7057459e9aSMatthew G Knepley PetscFunctionBegin; 7157459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7257459e9aSMatthew G Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 7357459e9aSMatthew G Knepley PetscValidPointer(values, 5); 7457459e9aSMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 7557459e9aSMatthew G Knepley if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 7657459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 7757459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 7857459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 7957459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 8057459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 8157459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 8257459e9aSMatthew G Knepley xfStart = fStart; xfEnd = xfStart+nXF; 8357459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 8495babc02SJed Brown zfStart = yfEnd; 8557459e9aSMatthew G Knepley if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 8657459e9aSMatthew G Knepley ierr = VecGetArray(v, &vArray);CHKERRQ(ierr); 8757459e9aSMatthew G Knepley if ((p >= cStart) || (p < cEnd)) { 8857459e9aSMatthew G Knepley /* Cell */ 8957459e9aSMatthew G Knepley if (dim == 1) { 9057459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 9157459e9aSMatthew G Knepley } else if (dim == 2) { 9257459e9aSMatthew G Knepley /* 4 faces, 4 vertices 9357459e9aSMatthew G Knepley Bottom-left vertex follows same order as cells 9457459e9aSMatthew G Knepley Bottom y-face same order as cells 9557459e9aSMatthew G Knepley Left x-face follows same order as cells 9657459e9aSMatthew G Knepley We number the quad: 9757459e9aSMatthew G Knepley 9857459e9aSMatthew G Knepley 8--3--7 9957459e9aSMatthew G Knepley | | 10057459e9aSMatthew G Knepley 4 0 2 10157459e9aSMatthew G Knepley | | 10257459e9aSMatthew G Knepley 5--1--6 10357459e9aSMatthew G Knepley */ 10457459e9aSMatthew G Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 10557459e9aSMatthew G Knepley PetscInt v = cy*nVx + cx + vStart; 10657459e9aSMatthew G Knepley PetscInt xf = cy*nxF + cx + xfStart; 10757459e9aSMatthew G Knepley PetscInt yf = c + yfStart; 10857459e9aSMatthew G Knepley PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0}; 10957459e9aSMatthew G Knepley 11057459e9aSMatthew G Knepley ierr = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr); 11157459e9aSMatthew G Knepley } else { 11257459e9aSMatthew G Knepley /* 6 faces, 8 vertices 11357459e9aSMatthew G Knepley Bottom-left-back vertex follows same order as cells 11457459e9aSMatthew G Knepley Back z-face follows same order as cells 11557459e9aSMatthew G Knepley Bottom y-face follows same order as cells 11657459e9aSMatthew G Knepley Left x-face follows same order as cells 11757459e9aSMatthew G Knepley 11857459e9aSMatthew G Knepley 14-----13 11957459e9aSMatthew G Knepley /| /| 12057459e9aSMatthew G Knepley / | 2 / | 12157459e9aSMatthew G Knepley / 5| / | 12257459e9aSMatthew G Knepley 10-----9 4| 12357459e9aSMatthew G Knepley | 11-|---12 12457459e9aSMatthew G Knepley |6 / | / 12557459e9aSMatthew G Knepley | /1 3| / 12657459e9aSMatthew G Knepley |/ |/ 12757459e9aSMatthew G Knepley 7-----8 12857459e9aSMatthew G Knepley */ 12957459e9aSMatthew G Knepley PetscInt c = p - cStart; 13057459e9aSMatthew G Knepley PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart, 13157459e9aSMatthew G Knepley c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0}; 13257459e9aSMatthew G Knepley 13357459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken"); 13457459e9aSMatthew G Knepley ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr); 13557459e9aSMatthew G Knepley } 13657459e9aSMatthew G Knepley } else if ((p >= vStart) || (p < vEnd)) { 13757459e9aSMatthew G Knepley /* Vertex */ 13857459e9aSMatthew G Knepley ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr); 13957459e9aSMatthew G Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 14057459e9aSMatthew G Knepley /* X Face */ 14157459e9aSMatthew G Knepley if (dim == 1) { 14257459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 14357459e9aSMatthew G Knepley } else if (dim == 2) { 14457459e9aSMatthew G Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 14557459e9aSMatthew G Knepley PetscInt f = p - xfStart; 14657459e9aSMatthew G Knepley PetscInt points[3] = {p, f, f+nVx}; 14757459e9aSMatthew G Knepley 14857459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken"); 14957459e9aSMatthew G Knepley ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr); 15057459e9aSMatthew G Knepley } else if (dim == 3) { 15157459e9aSMatthew G Knepley /* 4 vertices */ 15257459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 15357459e9aSMatthew G Knepley } 15457459e9aSMatthew G Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 15557459e9aSMatthew G Knepley /* Y Face */ 15657459e9aSMatthew G Knepley if (dim == 1) { 15757459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 15857459e9aSMatthew G Knepley } else if (dim == 2) { 15957459e9aSMatthew G Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 16057459e9aSMatthew G Knepley PetscInt f = p - yfStart; 16157459e9aSMatthew G Knepley PetscInt points[3] = {p, f, f+1}; 16257459e9aSMatthew G Knepley 16357459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken"); 16457459e9aSMatthew G Knepley ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr); 16557459e9aSMatthew G Knepley } else if (dim == 3) { 16657459e9aSMatthew G Knepley /* 4 vertices */ 16757459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 16857459e9aSMatthew G Knepley } 16957459e9aSMatthew G Knepley } else { 17057459e9aSMatthew G Knepley /* Z Face */ 17157459e9aSMatthew G Knepley if (dim == 1) { 17257459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 17357459e9aSMatthew G Knepley } else if (dim == 2) { 17457459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D"); 17557459e9aSMatthew G Knepley } else if (dim == 3) { 17657459e9aSMatthew G Knepley /* 4 vertices */ 17757459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 17857459e9aSMatthew G Knepley } 17957459e9aSMatthew G Knepley } 18057459e9aSMatthew G Knepley ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr); 18157459e9aSMatthew G Knepley PetscFunctionReturn(0); 18257459e9aSMatthew G Knepley } 18357459e9aSMatthew G Knepley 18457459e9aSMatthew G Knepley #undef __FUNCT__ 18557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure" 18657459e9aSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode) 18757459e9aSMatthew G Knepley { 18857459e9aSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 18957459e9aSMatthew G Knepley PetscInt dim = da->dim; 19057459e9aSMatthew G Knepley PetscScalar *vArray; 19157459e9aSMatthew G Knepley PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 19295babc02SJed Brown PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 19357459e9aSMatthew G Knepley PetscErrorCode ierr; 19457459e9aSMatthew G Knepley 19557459e9aSMatthew G Knepley PetscFunctionBegin; 19657459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 19757459e9aSMatthew G Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 19857459e9aSMatthew G Knepley PetscValidPointer(values, 5); 19957459e9aSMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 20057459e9aSMatthew G Knepley if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 20157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 20257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 20357459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 20457459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 20557459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 20657459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 20757459e9aSMatthew G Knepley xfStart = fStart; xfEnd = xfStart+nXF; 20857459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 20995babc02SJed Brown zfStart = yfEnd; 21057459e9aSMatthew G Knepley if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 21157459e9aSMatthew G Knepley ierr = VecGetArray(v, &vArray);CHKERRQ(ierr); 21257459e9aSMatthew G Knepley if ((p >= cStart) || (p < cEnd)) { 21357459e9aSMatthew G Knepley /* Cell */ 21457459e9aSMatthew G Knepley if (dim == 1) { 21557459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 21657459e9aSMatthew G Knepley } else if (dim == 2) { 21757459e9aSMatthew G Knepley /* 4 faces, 4 vertices 21857459e9aSMatthew G Knepley Bottom-left vertex follows same order as cells 21957459e9aSMatthew G Knepley Bottom y-face same order as cells 22057459e9aSMatthew G Knepley Left x-face follows same order as cells 22157459e9aSMatthew G Knepley We number the quad: 22257459e9aSMatthew G Knepley 22357459e9aSMatthew G Knepley 8--3--7 22457459e9aSMatthew G Knepley | | 22557459e9aSMatthew G Knepley 4 0 2 22657459e9aSMatthew G Knepley | | 22757459e9aSMatthew G Knepley 5--1--6 22857459e9aSMatthew G Knepley */ 22957459e9aSMatthew G Knepley PetscInt c = p - cStart; 23057459e9aSMatthew G Knepley PetscInt points[9] = {p, c+yfStart, c+xfStart+1, c+yfStart+nyF, c+xfStart+0, c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0}; 23157459e9aSMatthew G Knepley 23257459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr); 23357459e9aSMatthew G Knepley } else { 23457459e9aSMatthew G Knepley /* 6 faces, 8 vertices 23557459e9aSMatthew G Knepley Bottom-left-back vertex follows same order as cells 23657459e9aSMatthew G Knepley Back z-face follows same order as cells 23757459e9aSMatthew G Knepley Bottom y-face follows same order as cells 23857459e9aSMatthew G Knepley Left x-face follows same order as cells 23957459e9aSMatthew G Knepley 24057459e9aSMatthew G Knepley 14-----13 24157459e9aSMatthew G Knepley /| /| 24257459e9aSMatthew G Knepley / | 2 / | 24357459e9aSMatthew G Knepley / 5| / | 24457459e9aSMatthew G Knepley 10-----9 4| 24557459e9aSMatthew G Knepley | 11-|---12 24657459e9aSMatthew G Knepley |6 / | / 24757459e9aSMatthew G Knepley | /1 3| / 24857459e9aSMatthew G Knepley |/ |/ 24957459e9aSMatthew G Knepley 7-----8 25057459e9aSMatthew G Knepley */ 25157459e9aSMatthew G Knepley PetscInt c = p - cStart; 25257459e9aSMatthew G Knepley PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart, 25357459e9aSMatthew G Knepley c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0}; 25457459e9aSMatthew G Knepley 25557459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr); 25657459e9aSMatthew G Knepley } 25757459e9aSMatthew G Knepley } else if ((p >= vStart) || (p < vEnd)) { 25857459e9aSMatthew G Knepley /* Vertex */ 25957459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr); 26057459e9aSMatthew G Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 26157459e9aSMatthew G Knepley /* X Face */ 26257459e9aSMatthew G Knepley if (dim == 1) { 26357459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 26457459e9aSMatthew G Knepley } else if (dim == 2) { 26557459e9aSMatthew G Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 26657459e9aSMatthew G Knepley PetscInt f = p - xfStart; 26757459e9aSMatthew G Knepley PetscInt points[3] = {p, f, f+nVx}; 26857459e9aSMatthew G Knepley 26957459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 27057459e9aSMatthew G Knepley } else if (dim == 3) { 27157459e9aSMatthew G Knepley /* 4 vertices */ 27257459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 27357459e9aSMatthew G Knepley } 27457459e9aSMatthew G Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 27557459e9aSMatthew G Knepley /* Y Face */ 27657459e9aSMatthew G Knepley if (dim == 1) { 27757459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 27857459e9aSMatthew G Knepley } else if (dim == 2) { 27957459e9aSMatthew G Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 28057459e9aSMatthew G Knepley PetscInt f = p - yfStart; 28157459e9aSMatthew G Knepley PetscInt points[3] = {p, f, f+1}; 28257459e9aSMatthew G Knepley 28357459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 28457459e9aSMatthew G Knepley } else if (dim == 3) { 28557459e9aSMatthew G Knepley /* 4 vertices */ 28657459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 28757459e9aSMatthew G Knepley } 28857459e9aSMatthew G Knepley } else { 28957459e9aSMatthew G Knepley /* Z Face */ 29057459e9aSMatthew G Knepley if (dim == 1) { 29157459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D"); 29257459e9aSMatthew G Knepley } else if (dim == 2) { 29357459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D"); 29457459e9aSMatthew G Knepley } else if (dim == 3) { 29557459e9aSMatthew G Knepley /* 4 vertices */ 29657459e9aSMatthew G Knepley SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented"); 29757459e9aSMatthew G Knepley } 29857459e9aSMatthew G Knepley } 29957459e9aSMatthew G Knepley ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr); 30057459e9aSMatthew G Knepley PetscFunctionReturn(0); 30157459e9aSMatthew G Knepley } 30257459e9aSMatthew G Knepley 303*ed4e918cSMatthew G Knepley #undef __FUNCT__ 304*ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell" 305*ed4e918cSMatthew G Knepley /*@ 306*ed4e918cSMatthew G Knepley DMDACnvertToCell - Convert (i,j,k) to local cell number 307341c9bdaSMatthew G Knepley 308*ed4e918cSMatthew G Knepley Not Collective 309*ed4e918cSMatthew G Knepley 310*ed4e918cSMatthew G Knepley Input Parameter: 311*ed4e918cSMatthew G Knepley + da - the distributed array 312*ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k) 313*ed4e918cSMatthew G Knepley 314*ed4e918cSMatthew G Knepley Output Parameter: 315*ed4e918cSMatthew G Knepley . cell - the local cell number 316*ed4e918cSMatthew G Knepley 317*ed4e918cSMatthew G Knepley Level: developer 318*ed4e918cSMatthew G Knepley 319*ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure() 320*ed4e918cSMatthew G Knepley @*/ 321*ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) 322341c9bdaSMatthew G Knepley { 323341c9bdaSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 324341c9bdaSMatthew G Knepley const PetscInt dim = da->dim; 325341c9bdaSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 326341c9bdaSMatthew G Knepley const PetscInt nC = (mx )*(dim > 1 ? (my )*(dim > 2 ? (mz ) : 1) : 1); 327*ed4e918cSMatthew G Knepley const PetscInt il = s.i - da->Xs/da->w, jl = dim > 1 ? s.j - da->Ys : 0, kl = dim > 2 ? s.k - da->Zs : 0; 328341c9bdaSMatthew G Knepley 329341c9bdaSMatthew G Knepley PetscFunctionBegin; 330*ed4e918cSMatthew G Knepley *cell = -1; 331*ed4e918cSMatthew G Knepley if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w); 332*ed4e918cSMatthew G Knepley if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye); 333*ed4e918cSMatthew G Knepley if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze); 334*ed4e918cSMatthew G Knepley *cell = (kl*my + jl)*mx + il; 335*ed4e918cSMatthew G Knepley PetscFunctionReturn(0); 336341c9bdaSMatthew G Knepley } 337341c9bdaSMatthew G Knepley 33857459e9aSMatthew G Knepley #undef __FUNCT__ 33957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D" 34015d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 34157459e9aSMatthew G Knepley { 34257459e9aSMatthew G Knepley const PetscScalar x0 = vertices[0]; 34357459e9aSMatthew G Knepley const PetscScalar y0 = vertices[1]; 34457459e9aSMatthew G Knepley const PetscScalar x1 = vertices[2]; 34557459e9aSMatthew G Knepley const PetscScalar y1 = vertices[3]; 34657459e9aSMatthew G Knepley const PetscScalar x2 = vertices[4]; 34757459e9aSMatthew G Knepley const PetscScalar y2 = vertices[5]; 34857459e9aSMatthew G Knepley const PetscScalar x3 = vertices[6]; 34957459e9aSMatthew G Knepley const PetscScalar y3 = vertices[7]; 35057459e9aSMatthew G Knepley const PetscScalar f_01 = x2 - x1 - x3 + x0; 35157459e9aSMatthew G Knepley const PetscScalar g_01 = y2 - y1 - y3 + y0; 35257459e9aSMatthew G Knepley const PetscScalar x = refPoint[0]; 35357459e9aSMatthew G Knepley const PetscScalar y = refPoint[1]; 35457459e9aSMatthew G Knepley PetscReal invDet; 35557459e9aSMatthew G Knepley PetscErrorCode ierr; 35657459e9aSMatthew G Knepley 35757459e9aSMatthew G Knepley PetscFunctionBegin; 35815d2e57cSJed Brown #if defined(PETSC_USE_DEBUG) 35915d2e57cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n", 3604c06c320SMatthew G Knepley PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr); 36115d2e57cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr); 36215d2e57cSJed Brown #endif 36315d2e57cSJed Brown J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5; 36415d2e57cSJed Brown J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5; 36557459e9aSMatthew G Knepley *detJ = J[0]*J[3] - J[1]*J[2]; 36657459e9aSMatthew G Knepley invDet = 1.0/(*detJ); 36757459e9aSMatthew G Knepley invJ[0] = invDet*J[3]; invJ[1] = -invDet*J[1]; 36857459e9aSMatthew G Knepley invJ[2] = -invDet*J[2]; invJ[3] = invDet*J[0]; 36957459e9aSMatthew G Knepley ierr = PetscLogFlops(30);CHKERRQ(ierr); 37057459e9aSMatthew G Knepley PetscFunctionReturn(0); 37157459e9aSMatthew G Knepley } 37257459e9aSMatthew G Knepley 37357459e9aSMatthew G Knepley #undef __FUNCT__ 37457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry" 37557459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[]) 37657459e9aSMatthew G Knepley { 37757459e9aSMatthew G Knepley DM cdm; 37857459e9aSMatthew G Knepley Vec coordinates; 37957459e9aSMatthew G Knepley const PetscScalar *vertices; 38057459e9aSMatthew G Knepley PetscInt dim, d, q; 38157459e9aSMatthew G Knepley PetscErrorCode ierr; 38257459e9aSMatthew G Knepley 38357459e9aSMatthew G Knepley PetscFunctionBegin; 38457459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 38557459e9aSMatthew G Knepley ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 38657459e9aSMatthew G Knepley ierr = DMDAGetCoordinates(dm, &coordinates);CHKERRQ(ierr); 38757459e9aSMatthew G Knepley ierr = DMDAGetCoordinateDA(dm, &cdm);CHKERRQ(ierr); 38857459e9aSMatthew G Knepley ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr); 38957459e9aSMatthew G Knepley for(d = 0; d < dim; ++d) { 39015d2e57cSJed Brown v0[d] = PetscRealPart(vertices[d]); 39157459e9aSMatthew G Knepley } 39257459e9aSMatthew G Knepley switch(dim) { 39357459e9aSMatthew G Knepley case 2: 39457459e9aSMatthew G Knepley for(q = 0; q < quad->numQuadPoints; ++q) { 39557459e9aSMatthew G Knepley ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr); 39657459e9aSMatthew G Knepley } 39757459e9aSMatthew G Knepley break; 39857459e9aSMatthew G Knepley default: 39957459e9aSMatthew G Knepley SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim); 40057459e9aSMatthew G Knepley } 40157459e9aSMatthew G Knepley PetscFunctionReturn(0); 40257459e9aSMatthew G Knepley } 403