105264a50SDave May #include <petscsf.h> 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 357459e9aSMatthew G Knepley 4d7a12f93SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Static(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, PetscInt *csize, PetscScalar **array) 557459e9aSMatthew G Knepley { 657459e9aSMatthew G Knepley PetscScalar *a; 70a3ada39SMatthew G. Knepley PetscInt pStart, pEnd, size = 0, dof, off, d, k, i; 857459e9aSMatthew G Knepley PetscErrorCode ierr; 957459e9aSMatthew G Knepley 1057459e9aSMatthew G Knepley PetscFunctionBegin; 110a3ada39SMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 1257459e9aSMatthew G Knepley for (i = 0; i < nP; ++i) { 130a3ada39SMatthew G. Knepley const PetscInt p = points[i]; 140a3ada39SMatthew G. Knepley 150a3ada39SMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 160a3ada39SMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 1757459e9aSMatthew G Knepley size += dof; 1857459e9aSMatthew G Knepley } 190a3ada39SMatthew G. Knepley if (csize) *csize = size; 200a3ada39SMatthew G. Knepley if (array) { 21aa1993deSMatthew G Knepley ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr); 2257459e9aSMatthew G Knepley for (i = 0, k = 0; i < nP; ++i) { 230a3ada39SMatthew G. Knepley const PetscInt p = points[i]; 240a3ada39SMatthew G. Knepley 250a3ada39SMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 260a3ada39SMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 270a3ada39SMatthew G. Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 2857459e9aSMatthew G Knepley 298865f1eaSKarl Rupp for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d]; 3057459e9aSMatthew G Knepley } 3157459e9aSMatthew G Knepley *array = a; 320a3ada39SMatthew G. Knepley } 3357459e9aSMatthew G Knepley PetscFunctionReturn(0); 3457459e9aSMatthew G Knepley } 3557459e9aSMatthew G Knepley 3657459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode) 3757459e9aSMatthew G Knepley { 3857459e9aSMatthew G Knepley PetscInt dof, off, d, k, i; 3957459e9aSMatthew G Knepley PetscErrorCode ierr; 4057459e9aSMatthew G Knepley 4157459e9aSMatthew G Knepley PetscFunctionBegin; 4257459e9aSMatthew G Knepley if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) { 4357459e9aSMatthew G Knepley for (i = 0, k = 0; i < nP; ++i) { 4457459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 4557459e9aSMatthew G Knepley ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 4657459e9aSMatthew G Knepley 478865f1eaSKarl Rupp for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k]; 4857459e9aSMatthew G Knepley } 4957459e9aSMatthew G Knepley } else { 5057459e9aSMatthew G Knepley for (i = 0, k = 0; i < nP; ++i) { 5157459e9aSMatthew G Knepley ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 5257459e9aSMatthew G Knepley ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 5357459e9aSMatthew G Knepley 548865f1eaSKarl Rupp for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k]; 5557459e9aSMatthew G Knepley } 5657459e9aSMatthew G Knepley } 5757459e9aSMatthew G Knepley PetscFunctionReturn(0); 5857459e9aSMatthew G Knepley } 5957459e9aSMatthew G Knepley 60aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints) 61aa1993deSMatthew G Knepley { 62aa1993deSMatthew G Knepley PetscErrorCode ierr; 63aa1993deSMatthew G Knepley PetscInt *work; 64aa1993deSMatthew G Knepley 65aa1993deSMatthew G Knepley PetscFunctionBegin; 66aa1993deSMatthew G Knepley if (rn) *rn = n; 67aa1993deSMatthew G Knepley if (rpoints) { 68aa1993deSMatthew G Knepley ierr = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr); 69aa1993deSMatthew G Knepley ierr = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr); 70aa1993deSMatthew G Knepley *rpoints = work; 71aa1993deSMatthew G Knepley } 72aa1993deSMatthew G Knepley PetscFunctionReturn(0); 73aa1993deSMatthew G Knepley } 74aa1993deSMatthew G Knepley 75aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints) 76aa1993deSMatthew G Knepley { 77aa1993deSMatthew G Knepley PetscErrorCode ierr; 78aa1993deSMatthew G Knepley 79aa1993deSMatthew G Knepley PetscFunctionBegin; 80aa1993deSMatthew G Knepley if (rn) *rn = 0; 818865f1eaSKarl Rupp if (rpoints) { 822a808120SBarry Smith /* note the second argument to DMRestoreWorkArray() is always ignored */ 832a808120SBarry Smith ierr = DMRestoreWorkArray(dm,0, PETSC_INT, (void*) rpoints);CHKERRQ(ierr); 848865f1eaSKarl Rupp } 85aa1993deSMatthew G Knepley PetscFunctionReturn(0); 86aa1993deSMatthew G Knepley } 87aa1993deSMatthew G Knepley 886693a731SMatthew G. Knepley PetscErrorCode DMDAGetTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure) 896693a731SMatthew G. Knepley { 90c73cfb54SMatthew G. Knepley PetscInt dim = dm->dim; 916693a731SMatthew G. Knepley PetscInt nVx, nVy, nVz, nxF, nXF, nyF, nYF, nzF, nZF; 926693a731SMatthew G. Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart; 936693a731SMatthew G. Knepley PetscErrorCode ierr; 946693a731SMatthew G. Knepley 956693a731SMatthew G. Knepley PetscFunctionBegin; 966693a731SMatthew G. Knepley if (!useClosure) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Star operation is not yet supported"); 976693a731SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 982a808120SBarry Smith PetscValidIntPointer(closureSize,4); 996693a731SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 1006693a731SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1016693a731SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1026693a731SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 1036693a731SMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, NULL);CHKERRQ(ierr); 1046693a731SMatthew G. Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 1056693a731SMatthew G. Knepley xfStart = fStart; xfEnd = xfStart+nXF; 1066693a731SMatthew G. Knepley yfStart = xfEnd; 1076693a731SMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 1086693a731SMatthew G. Knepley if ((p >= cStart) || (p < cEnd)) { 1096693a731SMatthew G. Knepley /* Cell */ 1106693a731SMatthew G. Knepley if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 1116693a731SMatthew G. Knepley else if (dim == 2) { 1126693a731SMatthew G. Knepley /* 4 faces, 4 vertices 1136693a731SMatthew G. Knepley Bottom-left vertex follows same order as cells 1146693a731SMatthew G. Knepley Bottom y-face same order as cells 1156693a731SMatthew G. Knepley Left x-face follows same order as cells 1166693a731SMatthew G. Knepley We number the quad: 1176693a731SMatthew G. Knepley 1186693a731SMatthew G. Knepley 8--3--7 1196693a731SMatthew G. Knepley | | 1206693a731SMatthew G. Knepley 4 0 2 1216693a731SMatthew G. Knepley | | 1226693a731SMatthew G. Knepley 5--1--6 1236693a731SMatthew G. Knepley */ 1246693a731SMatthew G. Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 1256693a731SMatthew G. Knepley PetscInt v = cy*nVx + cx + vStart; 1266693a731SMatthew G. Knepley PetscInt xf = cy*nxF + cx + xfStart; 1276693a731SMatthew G. Knepley PetscInt yf = c + yfStart; 1286693a731SMatthew G. Knepley 1292a808120SBarry Smith *closureSize = 9; 1306693a731SMatthew G. Knepley if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 1316693a731SMatthew G. Knepley (*closure)[0] = p; (*closure)[1] = yf; (*closure)[2] = xf+1; (*closure)[3] = yf+nyF; (*closure)[4] = xf+0; (*closure)[5] = v+0; (*closure)[6]= v+1; (*closure)[7] = v+nVx+1; (*closure)[8] = v+nVx+0; 1326693a731SMatthew G. Knepley } else { 1336693a731SMatthew G. Knepley /* 6 faces, 12 edges, 8 vertices 1346693a731SMatthew G. Knepley Bottom-left vertex follows same order as cells 1356693a731SMatthew G. Knepley Bottom y-face same order as cells 1366693a731SMatthew G. Knepley Left x-face follows same order as cells 1376693a731SMatthew G. Knepley We number the quad: 1386693a731SMatthew G. Knepley 1396693a731SMatthew G. Knepley 8--3--7 1406693a731SMatthew G. Knepley | | 1416693a731SMatthew G. Knepley 4 0 2 1426693a731SMatthew G. Knepley | | 1436693a731SMatthew G. Knepley 5--1--6 1446693a731SMatthew G. Knepley */ 145be8e0784SMatthew G. Knepley #if 0 1466693a731SMatthew G. Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1), cz = c / ((nVx-1)*(nVy-1)); 1476693a731SMatthew G. Knepley PetscInt v = cy*nVx + cx + vStart; 1486693a731SMatthew G. Knepley PetscInt xf = cy*nxF + cx + xfStart; 1496693a731SMatthew G. Knepley PetscInt yf = c + yfStart; 1506693a731SMatthew G. Knepley 1512a808120SBarry Smith *closureSize = 26; 1526693a731SMatthew G. Knepley if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 153be8e0784SMatthew G. Knepley #endif 154be8e0784SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 1556693a731SMatthew G. Knepley } 1566693a731SMatthew G. Knepley } else if ((p >= vStart) || (p < vEnd)) { 1576693a731SMatthew G. Knepley /* Vertex */ 1582a808120SBarry Smith *closureSize = 1; 1596693a731SMatthew G. Knepley if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 1606693a731SMatthew G. Knepley (*closure)[0] = p; 1616693a731SMatthew G. Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 1626693a731SMatthew G. Knepley /* X Face */ 1636693a731SMatthew G. Knepley if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 1646693a731SMatthew G. Knepley else if (dim == 2) { 1656693a731SMatthew G. Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 1666693a731SMatthew G. Knepley PetscInt f = p - xfStart; 1676693a731SMatthew G. Knepley 1682a808120SBarry Smith *closureSize = 3; 1696693a731SMatthew G. Knepley if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 1706693a731SMatthew G. Knepley (*closure)[0] = p; (*closure)[1] = f; (*closure)[2] = f+nVx; 1716693a731SMatthew G. Knepley } else if (dim == 3) { 1726693a731SMatthew G. Knepley /* 4 vertices */ 1736693a731SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 1746693a731SMatthew G. Knepley } 1756693a731SMatthew G. Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 1766693a731SMatthew G. Knepley /* Y Face */ 1776693a731SMatthew G. Knepley if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 1786693a731SMatthew G. Knepley else if (dim == 2) { 1796693a731SMatthew G. Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 1806693a731SMatthew G. Knepley PetscInt f = p - yfStart; 1816693a731SMatthew G. Knepley 1822a808120SBarry Smith *closureSize = 3; 1836693a731SMatthew G. Knepley if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 1846693a731SMatthew G. Knepley (*closure)[0] = p; (*closure)[1] = f; (*closure)[2]= f+1; 1856693a731SMatthew G. Knepley } else if (dim == 3) { 1866693a731SMatthew G. Knepley /* 4 vertices */ 1876693a731SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 1886693a731SMatthew G. Knepley } 1896693a731SMatthew G. Knepley } else { 1906693a731SMatthew G. Knepley /* Z Face */ 1916693a731SMatthew G. Knepley if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 1926693a731SMatthew G. Knepley else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 1936693a731SMatthew G. Knepley else if (dim == 3) { 1946693a731SMatthew G. Knepley /* 4 vertices */ 1956693a731SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 1966693a731SMatthew G. Knepley } 1976693a731SMatthew G. Knepley } 1986693a731SMatthew G. Knepley PetscFunctionReturn(0); 1996693a731SMatthew G. Knepley } 2006693a731SMatthew G. Knepley 2016693a731SMatthew G. Knepley PetscErrorCode DMDARestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure) 2026693a731SMatthew G. Knepley { 2036693a731SMatthew G. Knepley PetscErrorCode ierr; 2046693a731SMatthew G. Knepley 2056693a731SMatthew G. Knepley PetscFunctionBegin; 2066693a731SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, closure);CHKERRQ(ierr); 2076693a731SMatthew G. Knepley PetscFunctionReturn(0); 2086693a731SMatthew G. Knepley } 2096693a731SMatthew G. Knepley 210aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure) 21157459e9aSMatthew G Knepley { 212c73cfb54SMatthew G. Knepley PetscInt dim = dm->dim; 213aa1993deSMatthew G Knepley PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 214aa1993deSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart; 215aa1993deSMatthew G Knepley PetscErrorCode ierr; 216aa1993deSMatthew G Knepley 217aa1993deSMatthew G Knepley PetscFunctionBegin; 218aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 219aa1993deSMatthew G Knepley if (n) PetscValidIntPointer(n,4); 220aa1993deSMatthew G Knepley if (closure) PetscValidPointer(closure, 5); 221aa1993deSMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 22282f516ccSBarry Smith if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 223aa1993deSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 224aa1993deSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 225aa1993deSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 226aa1993deSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 2270298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 228aa1993deSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 229aa1993deSMatthew G Knepley xfStart = fStart; xfEnd = xfStart+nXF; 230aa1993deSMatthew G Knepley yfStart = xfEnd; 23182f516ccSBarry Smith if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 232aa1993deSMatthew G Knepley if ((p >= cStart) || (p < cEnd)) { 233aa1993deSMatthew G Knepley /* Cell */ 23482f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 235201c01f8SBarry Smith else if (dim == 2) { 236aa1993deSMatthew G Knepley /* 4 faces, 4 vertices 237aa1993deSMatthew G Knepley Bottom-left vertex follows same order as cells 238aa1993deSMatthew G Knepley Bottom y-face same order as cells 239aa1993deSMatthew G Knepley Left x-face follows same order as cells 240aa1993deSMatthew G Knepley We number the quad: 241aa1993deSMatthew G Knepley 242aa1993deSMatthew G Knepley 8--3--7 243aa1993deSMatthew G Knepley | | 244aa1993deSMatthew G Knepley 4 0 2 245aa1993deSMatthew G Knepley | | 246aa1993deSMatthew G Knepley 5--1--6 247aa1993deSMatthew G Knepley */ 248aa1993deSMatthew G Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 249aa1993deSMatthew G Knepley PetscInt v = cy*nVx + cx + vStart; 250aa1993deSMatthew G Knepley PetscInt xf = cy*nxF + cx + xfStart; 251aa1993deSMatthew G Knepley PetscInt yf = c + yfStart; 252da80777bSKarl Rupp PetscInt points[9]; 253da80777bSKarl Rupp 254da80777bSKarl Rupp /* Note: initializer list is not computable at compile time, hence we fill the array manually */ 255da80777bSKarl Rupp points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6]= v+1; points[7] = v+nVx+1; points[8] = v+nVx+0; 256aa1993deSMatthew G Knepley 257aa1993deSMatthew G Knepley ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr); 25882f516ccSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 259aa1993deSMatthew G Knepley } else if ((p >= vStart) || (p < vEnd)) { 260aa1993deSMatthew G Knepley /* Vertex */ 261aa1993deSMatthew G Knepley ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr); 262aa1993deSMatthew G Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 263aa1993deSMatthew G Knepley /* X Face */ 26482f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 265201c01f8SBarry Smith else if (dim == 2) { 266aa1993deSMatthew G Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 267aa1993deSMatthew G Knepley PetscInt f = p - xfStart; 268da80777bSKarl Rupp PetscInt points[3]; 269aa1993deSMatthew G Knepley 270da80777bSKarl Rupp points[0] = p; points[1] = f; points[2] = f+nVx; 27182f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 272aa1993deSMatthew G Knepley ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr); 273aa1993deSMatthew G Knepley } else if (dim == 3) { 274aa1993deSMatthew G Knepley /* 4 vertices */ 27582f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 276aa1993deSMatthew G Knepley } 277aa1993deSMatthew G Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 278aa1993deSMatthew G Knepley /* Y Face */ 27982f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 280201c01f8SBarry Smith else if (dim == 2) { 281aa1993deSMatthew G Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 282aa1993deSMatthew G Knepley PetscInt f = p - yfStart; 283da80777bSKarl Rupp PetscInt points[3]; 284aa1993deSMatthew G Knepley 285da80777bSKarl Rupp points[0] = p; points[1] = f; points[2]= f+1; 28682f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 287aa1993deSMatthew G Knepley ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr); 288aa1993deSMatthew G Knepley } else if (dim == 3) { 289aa1993deSMatthew G Knepley /* 4 vertices */ 29082f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 291aa1993deSMatthew G Knepley } 292aa1993deSMatthew G Knepley } else { 293aa1993deSMatthew G Knepley /* Z Face */ 29482f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 29582f516ccSBarry Smith else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 296201c01f8SBarry Smith else if (dim == 3) { 297aa1993deSMatthew G Knepley /* 4 vertices */ 29882f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 299aa1993deSMatthew G Knepley } 300aa1993deSMatthew G Knepley } 301aa1993deSMatthew G Knepley PetscFunctionReturn(0); 302aa1993deSMatthew G Knepley } 303aa1993deSMatthew G Knepley 304e5c84f05SJed Brown /* Zeros n and closure. */ 305aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure) 306aa1993deSMatthew G Knepley { 307aa1993deSMatthew G Knepley PetscErrorCode ierr; 308aa1993deSMatthew G Knepley 309aa1993deSMatthew G Knepley PetscFunctionBegin; 310aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 311e5c84f05SJed Brown if (n) PetscValidIntPointer(n,4); 312e5c84f05SJed Brown if (closure) PetscValidPointer(closure, 5); 313aa1993deSMatthew G Knepley ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr); 314aa1993deSMatthew G Knepley PetscFunctionReturn(0); 315aa1993deSMatthew G Knepley } 316aa1993deSMatthew G Knepley 317d7a12f93SMatthew G. Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, PetscScalar **values) 318aa1993deSMatthew G Knepley { 319c73cfb54SMatthew G. Knepley PetscInt dim = dm->dim; 32057459e9aSMatthew G Knepley PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 32195babc02SJed Brown PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 32257459e9aSMatthew G Knepley PetscErrorCode ierr; 32357459e9aSMatthew G Knepley 32457459e9aSMatthew G Knepley PetscFunctionBegin; 32557459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 326aa1993deSMatthew G Knepley PetscValidScalarPointer(vArray, 4); 3270a3ada39SMatthew G. Knepley if (values) PetscValidPointer(values, 6); 32857459e9aSMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 32982f516ccSBarry Smith if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 33057459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 33157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 33257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 33357459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 3340298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 33557459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 33657459e9aSMatthew G Knepley xfStart = fStart; xfEnd = xfStart+nXF; 33757459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 33895babc02SJed Brown zfStart = yfEnd; 33982f516ccSBarry Smith if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 34057459e9aSMatthew G Knepley if ((p >= cStart) || (p < cEnd)) { 34157459e9aSMatthew G Knepley /* Cell */ 34282f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 343201c01f8SBarry Smith else if (dim == 2) { 34457459e9aSMatthew G Knepley /* 4 faces, 4 vertices 34557459e9aSMatthew G Knepley Bottom-left vertex follows same order as cells 34657459e9aSMatthew G Knepley Bottom y-face same order as cells 34757459e9aSMatthew G Knepley Left x-face follows same order as cells 34857459e9aSMatthew G Knepley We number the quad: 34957459e9aSMatthew G Knepley 35057459e9aSMatthew G Knepley 8--3--7 35157459e9aSMatthew G Knepley | | 35257459e9aSMatthew G Knepley 4 0 2 35357459e9aSMatthew G Knepley | | 35457459e9aSMatthew G Knepley 5--1--6 35557459e9aSMatthew G Knepley */ 35657459e9aSMatthew G Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 35757459e9aSMatthew G Knepley PetscInt v = cy*nVx + cx + vStart; 358312b75c6SMatthew G. Knepley PetscInt xf = cx*nxF + cy + xfStart; 35957459e9aSMatthew G Knepley PetscInt yf = c + yfStart; 360da80777bSKarl Rupp PetscInt points[9]; 36157459e9aSMatthew G Knepley 362312b75c6SMatthew G. Knepley points[0] = p; 363312b75c6SMatthew G. Knepley points[1] = yf; points[2] = xf+nxF; points[3] = yf+nyF; points[4] = xf; 364312b75c6SMatthew G. Knepley points[5] = v+0; points[6] = v+1; points[7] = v+nVx+1; points[8] = v+nVx+0; 3650a3ada39SMatthew G. Knepley ierr = FillClosureArray_Static(dm, section, 9, points, vArray, csize, values);CHKERRQ(ierr); 36657459e9aSMatthew G Knepley } else { 36757459e9aSMatthew G Knepley /* 6 faces, 8 vertices 36857459e9aSMatthew G Knepley Bottom-left-back vertex follows same order as cells 36957459e9aSMatthew G Knepley Back z-face follows same order as cells 37057459e9aSMatthew G Knepley Bottom y-face follows same order as cells 37157459e9aSMatthew G Knepley Left x-face follows same order as cells 37257459e9aSMatthew G Knepley 37357459e9aSMatthew G Knepley 14-----13 37457459e9aSMatthew G Knepley /| /| 37557459e9aSMatthew G Knepley / | 2 / | 37657459e9aSMatthew G Knepley / 5| / | 37757459e9aSMatthew G Knepley 10-----9 4| 37857459e9aSMatthew G Knepley | 11-|---12 37957459e9aSMatthew G Knepley |6 / | / 38057459e9aSMatthew G Knepley | /1 3| / 38157459e9aSMatthew G Knepley |/ |/ 38257459e9aSMatthew G Knepley 7-----8 38357459e9aSMatthew G Knepley */ 38457459e9aSMatthew G Knepley PetscInt c = p - cStart; 385da80777bSKarl Rupp PetscInt points[15]; 386da80777bSKarl Rupp 387da80777bSKarl Rupp points[0] = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart; 388da80777bSKarl Rupp points[7] = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; 389da80777bSKarl Rupp points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0; 39057459e9aSMatthew G Knepley 39182f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 3920a3ada39SMatthew G. Knepley ierr = FillClosureArray_Static(dm, section, 15, points, vArray, csize, values);CHKERRQ(ierr); 39357459e9aSMatthew G Knepley } 39457459e9aSMatthew G Knepley } else if ((p >= vStart) || (p < vEnd)) { 39557459e9aSMatthew G Knepley /* Vertex */ 3960a3ada39SMatthew G. Knepley ierr = FillClosureArray_Static(dm, section, 1, &p, vArray, csize, values);CHKERRQ(ierr); 39757459e9aSMatthew G Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 39857459e9aSMatthew G Knepley /* X Face */ 39982f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 400201c01f8SBarry Smith else if (dim == 2) { 40157459e9aSMatthew G Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 40257459e9aSMatthew G Knepley PetscInt f = p - xfStart; 403da80777bSKarl Rupp PetscInt points[3]; 40457459e9aSMatthew G Knepley 405da80777bSKarl Rupp points[0] = p; points[1] = f; points[2] = f+nVx; 40682f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 4070a3ada39SMatthew G. Knepley ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr); 40882f516ccSBarry Smith } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 40957459e9aSMatthew G Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 41057459e9aSMatthew G Knepley /* Y Face */ 41182f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 412201c01f8SBarry Smith else if (dim == 2) { 41357459e9aSMatthew G Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 41457459e9aSMatthew G Knepley PetscInt f = p - yfStart; 415da80777bSKarl Rupp PetscInt points[3]; 41657459e9aSMatthew G Knepley 417da80777bSKarl Rupp points[0] = p; points[1] = f; points[2] = f+1; 41882f516ccSBarry Smith SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 4190a3ada39SMatthew G. Knepley ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr); 42082f516ccSBarry Smith } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 42157459e9aSMatthew G Knepley } else { 42257459e9aSMatthew G Knepley /* Z Face */ 42382f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 42482f516ccSBarry Smith else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 42582f516ccSBarry Smith else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 42657459e9aSMatthew G Knepley } 427aa1993deSMatthew G Knepley PetscFunctionReturn(0); 428aa1993deSMatthew G Knepley } 429aa1993deSMatthew G Knepley 430d7a12f93SMatthew G. Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, PetscScalar **values) 431aa1993deSMatthew G Knepley { 432aa1993deSMatthew G Knepley PetscScalar *vArray; 433aa1993deSMatthew G Knepley PetscErrorCode ierr; 434aa1993deSMatthew G Knepley 435aa1993deSMatthew G Knepley PetscFunctionBegin; 436aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 437aa1993deSMatthew G Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 4380a3ada39SMatthew G. Knepley if (values) PetscValidPointer(values, 6); 439aa1993deSMatthew G Knepley ierr = VecGetArray(v, &vArray);CHKERRQ(ierr); 4400a3ada39SMatthew G. Knepley ierr = DMDAGetClosureScalar(dm, section, p, vArray, csize, values);CHKERRQ(ierr); 44157459e9aSMatthew G Knepley ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr); 44257459e9aSMatthew G Knepley PetscFunctionReturn(0); 44357459e9aSMatthew G Knepley } 44457459e9aSMatthew G Knepley 445d7a12f93SMatthew G. Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, PetscScalar **values) 446aa1993deSMatthew G Knepley { 447aa1993deSMatthew G Knepley PetscErrorCode ierr; 448aa1993deSMatthew G Knepley 449aa1993deSMatthew G Knepley PetscFunctionBegin; 450aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4510a3ada39SMatthew G. Knepley PetscValidPointer(values, 6); 4520a3ada39SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, csize ? *csize : 0, PETSC_SCALAR, (void*) values);CHKERRQ(ierr); 453aa1993deSMatthew G Knepley PetscFunctionReturn(0); 454aa1993deSMatthew G Knepley } 455aa1993deSMatthew G Knepley 456d7a12f93SMatthew G. Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, PetscScalar **values) 457aa1993deSMatthew G Knepley { 458aa1993deSMatthew G Knepley PetscErrorCode ierr; 459aa1993deSMatthew G Knepley 460aa1993deSMatthew G Knepley PetscFunctionBegin; 461aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 462aa1993deSMatthew G Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 463aa1993deSMatthew G Knepley PetscValidPointer(values, 5); 4640a3ada39SMatthew G. Knepley ierr = DMDARestoreClosureScalar(dm, section, p, NULL, csize, values);CHKERRQ(ierr); 465aa1993deSMatthew G Knepley PetscFunctionReturn(0); 466aa1993deSMatthew G Knepley } 467aa1993deSMatthew G Knepley 468aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode) 46957459e9aSMatthew G Knepley { 470c73cfb54SMatthew G. Knepley PetscInt dim = dm->dim; 4716a7fb054SMatthew G. Knepley PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF, nCx, nCy; 47295babc02SJed Brown PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 47357459e9aSMatthew G Knepley PetscErrorCode ierr; 47457459e9aSMatthew G Knepley 47557459e9aSMatthew G Knepley PetscFunctionBegin; 47657459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 477aa1993deSMatthew G Knepley PetscValidScalarPointer(values, 4); 47857459e9aSMatthew G Knepley PetscValidPointer(values, 5); 47957459e9aSMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 48082f516ccSBarry Smith if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 48157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 48257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 48357459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 48457459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 4856a7fb054SMatthew G. Knepley ierr = DMDAGetNumCells(dm, &nCx, &nCy, NULL, NULL);CHKERRQ(ierr); 4860298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 48757459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 48857459e9aSMatthew G Knepley xfStart = fStart; xfEnd = xfStart+nXF; 48957459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 49095babc02SJed Brown zfStart = yfEnd; 49182f516ccSBarry Smith if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 49257459e9aSMatthew G Knepley if ((p >= cStart) || (p < cEnd)) { 49357459e9aSMatthew G Knepley /* Cell */ 49482f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 495201c01f8SBarry Smith else if (dim == 2) { 49657459e9aSMatthew G Knepley /* 4 faces, 4 vertices 49757459e9aSMatthew G Knepley Bottom-left vertex follows same order as cells 49857459e9aSMatthew G Knepley Bottom y-face same order as cells 49957459e9aSMatthew G Knepley Left x-face follows same order as cells 50057459e9aSMatthew G Knepley We number the quad: 50157459e9aSMatthew G Knepley 50257459e9aSMatthew G Knepley 8--3--7 50357459e9aSMatthew G Knepley | | 50457459e9aSMatthew G Knepley 4 0 2 50557459e9aSMatthew G Knepley | | 50657459e9aSMatthew G Knepley 5--1--6 50757459e9aSMatthew G Knepley */ 508312b75c6SMatthew G. Knepley PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 509312b75c6SMatthew G. Knepley PetscInt v = cy*nVx + cx + vStart; 510312b75c6SMatthew G. Knepley PetscInt xf = cx*nxF + cy + xfStart; 511312b75c6SMatthew G. Knepley PetscInt yf = c + yfStart; 512da80777bSKarl Rupp PetscInt points[9]; 51357459e9aSMatthew G Knepley 5146a7fb054SMatthew G. Knepley points[0] = p; 515312b75c6SMatthew G. Knepley points[1] = yf; points[2] = xf+nxF; points[3] = yf+nyF; points[4] = xf; 516312b75c6SMatthew G. Knepley points[5] = v+0; points[6] = v+1; points[7] = v+nVx+1; points[8] = v+nVx+0; 51757459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr); 51857459e9aSMatthew G Knepley } else { 51957459e9aSMatthew G Knepley /* 6 faces, 8 vertices 52057459e9aSMatthew G Knepley Bottom-left-back vertex follows same order as cells 52157459e9aSMatthew G Knepley Back z-face follows same order as cells 52257459e9aSMatthew G Knepley Bottom y-face follows same order as cells 52357459e9aSMatthew G Knepley Left x-face follows same order as cells 52457459e9aSMatthew G Knepley 52557459e9aSMatthew G Knepley 14-----13 52657459e9aSMatthew G Knepley /| /| 52757459e9aSMatthew G Knepley / | 2 / | 52857459e9aSMatthew G Knepley / 5| / | 52957459e9aSMatthew G Knepley 10-----9 4| 53057459e9aSMatthew G Knepley | 11-|---12 53157459e9aSMatthew G Knepley |6 / | / 53257459e9aSMatthew G Knepley | /1 3| / 53357459e9aSMatthew G Knepley |/ |/ 53457459e9aSMatthew G Knepley 7-----8 53557459e9aSMatthew G Knepley */ 53657459e9aSMatthew G Knepley PetscInt c = p - cStart; 537da80777bSKarl Rupp PetscInt points[15]; 53857459e9aSMatthew G Knepley 539da80777bSKarl Rupp points[0] = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart; 540da80777bSKarl Rupp points[7] = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; points[12] = c+vStart+nVx*nVy+1; 541da80777bSKarl Rupp points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0; 54257459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr); 54357459e9aSMatthew G Knepley } 54457459e9aSMatthew G Knepley } else if ((p >= vStart) || (p < vEnd)) { 54557459e9aSMatthew G Knepley /* Vertex */ 54657459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr); 54757459e9aSMatthew G Knepley } else if ((p >= fStart) || (p < fStart + nXF)) { 54857459e9aSMatthew G Knepley /* X Face */ 54982f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 550201c01f8SBarry Smith else if (dim == 2) { 55157459e9aSMatthew G Knepley /* 2 vertices: The bottom vertex has the same numbering as the face */ 55257459e9aSMatthew G Knepley PetscInt f = p - xfStart; 553da80777bSKarl Rupp PetscInt points[3]; 55457459e9aSMatthew G Knepley 555da80777bSKarl Rupp points[0] = p; points[1] = f; points[2] = f+nVx; 55657459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 55782f516ccSBarry Smith } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 55857459e9aSMatthew G Knepley } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 55957459e9aSMatthew G Knepley /* Y Face */ 56082f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 561201c01f8SBarry Smith else if (dim == 2) { 56257459e9aSMatthew G Knepley /* 2 vertices: The left vertex has the same numbering as the face */ 56357459e9aSMatthew G Knepley PetscInt f = p - yfStart; 564da80777bSKarl Rupp PetscInt points[3]; 56557459e9aSMatthew G Knepley 566da80777bSKarl Rupp points[0] = p; points[1] = f; points[2] = f+1; 56757459e9aSMatthew G Knepley ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 56882f516ccSBarry Smith } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 56957459e9aSMatthew G Knepley } else { 57057459e9aSMatthew G Knepley /* Z Face */ 57182f516ccSBarry Smith if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 57282f516ccSBarry Smith else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 57382f516ccSBarry Smith else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 57457459e9aSMatthew G Knepley } 575aa1993deSMatthew G Knepley PetscFunctionReturn(0); 576aa1993deSMatthew G Knepley } 577aa1993deSMatthew G Knepley 578aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode) 579aa1993deSMatthew G Knepley { 580aa1993deSMatthew G Knepley PetscScalar *vArray; 581aa1993deSMatthew G Knepley PetscErrorCode ierr; 582aa1993deSMatthew G Knepley 583aa1993deSMatthew G Knepley PetscFunctionBegin; 584aa1993deSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 585aa1993deSMatthew G Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 586aa1993deSMatthew G Knepley PetscValidPointer(values, 5); 587aa1993deSMatthew G Knepley ierr = VecGetArray(v,&vArray);CHKERRQ(ierr); 588aa1993deSMatthew G Knepley ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr); 58957459e9aSMatthew G Knepley ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr); 59057459e9aSMatthew G Knepley PetscFunctionReturn(0); 59157459e9aSMatthew G Knepley } 59257459e9aSMatthew G Knepley 593ed4e918cSMatthew G Knepley /*@ 594f0e3914cSSatish Balay DMDAConvertToCell - Convert (i,j,k) to local cell number 595341c9bdaSMatthew G Knepley 596ed4e918cSMatthew G Knepley Not Collective 597ed4e918cSMatthew G Knepley 598ed4e918cSMatthew G Knepley Input Parameter: 599ed4e918cSMatthew G Knepley + da - the distributed array 600907376e6SBarry Smith - s - A MatStencil giving (i,j,k) 601ed4e918cSMatthew G Knepley 602ed4e918cSMatthew G Knepley Output Parameter: 603ed4e918cSMatthew G Knepley . cell - the local cell number 604ed4e918cSMatthew G Knepley 605ed4e918cSMatthew G Knepley Level: developer 606ed4e918cSMatthew G Knepley 607ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure() 608ed4e918cSMatthew G Knepley @*/ 609ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) 610341c9bdaSMatthew G Knepley { 611341c9bdaSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 612c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 6134e846693SMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/; 614ed4e918cSMatthew 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; 615341c9bdaSMatthew G Knepley 616341c9bdaSMatthew G Knepley PetscFunctionBegin; 617ed4e918cSMatthew G Knepley *cell = -1; 61882f516ccSBarry Smith if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w); 61982f516ccSBarry Smith if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye); 62082f516ccSBarry Smith if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze); 621ed4e918cSMatthew G Knepley *cell = (kl*my + jl)*mx + il; 622ed4e918cSMatthew G Knepley PetscFunctionReturn(0); 623341c9bdaSMatthew G Knepley } 624341c9bdaSMatthew G Knepley 62515d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 62657459e9aSMatthew G Knepley { 62757459e9aSMatthew G Knepley const PetscScalar x0 = vertices[0]; 62857459e9aSMatthew G Knepley const PetscScalar y0 = vertices[1]; 62957459e9aSMatthew G Knepley const PetscScalar x1 = vertices[2]; 63057459e9aSMatthew G Knepley const PetscScalar y1 = vertices[3]; 63157459e9aSMatthew G Knepley const PetscScalar x2 = vertices[4]; 63257459e9aSMatthew G Knepley const PetscScalar y2 = vertices[5]; 63357459e9aSMatthew G Knepley const PetscScalar x3 = vertices[6]; 63457459e9aSMatthew G Knepley const PetscScalar y3 = vertices[7]; 63557459e9aSMatthew G Knepley const PetscScalar f_01 = x2 - x1 - x3 + x0; 63657459e9aSMatthew G Knepley const PetscScalar g_01 = y2 - y1 - y3 + y0; 63757459e9aSMatthew G Knepley const PetscScalar x = refPoint[0]; 63857459e9aSMatthew G Knepley const PetscScalar y = refPoint[1]; 63957459e9aSMatthew G Knepley PetscReal invDet; 64057459e9aSMatthew G Knepley PetscErrorCode ierr; 64157459e9aSMatthew G Knepley 64257459e9aSMatthew G Knepley PetscFunctionBegin; 643e477d84eSMatthew G. Knepley #if 0 64415d2e57cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n", 6454c06c320SMatthew G Knepley PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr); 64615d2e57cSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr); 64715d2e57cSJed Brown #endif 64815d2e57cSJed Brown J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5; 64915d2e57cSJed Brown J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5; 65057459e9aSMatthew G Knepley *detJ = J[0]*J[3] - J[1]*J[2]; 65157459e9aSMatthew G Knepley invDet = 1.0/(*detJ); 652e477d84eSMatthew G. Knepley if (invJ) { 65357459e9aSMatthew G Knepley invJ[0] = invDet*J[3]; invJ[1] = -invDet*J[1]; 65457459e9aSMatthew G Knepley invJ[2] = -invDet*J[2]; invJ[3] = invDet*J[0]; 655e477d84eSMatthew G. Knepley } 65657459e9aSMatthew G Knepley ierr = PetscLogFlops(30);CHKERRQ(ierr); 65757459e9aSMatthew G Knepley PetscFunctionReturn(0); 65857459e9aSMatthew G Knepley } 65957459e9aSMatthew G Knepley 6608e0841e0SMatthew G. Knepley PetscErrorCode DMDAComputeCellGeometryFEM(DM dm, PetscInt cell, PetscQuadrature quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[]) 66157459e9aSMatthew G Knepley { 66257459e9aSMatthew G Knepley DM cdm; 66357459e9aSMatthew G Knepley Vec coordinates; 66421454ff5SMatthew G. Knepley const PetscReal *quadPoints; 665d7a12f93SMatthew G. Knepley PetscScalar *vertices = NULL; 6669c3cf19fSMatthew G. Knepley PetscInt Nq, csize, dim, d, q; 66757459e9aSMatthew G Knepley PetscErrorCode ierr; 66857459e9aSMatthew G Knepley 66957459e9aSMatthew G Knepley PetscFunctionBegin; 67057459e9aSMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 67157459e9aSMatthew G Knepley ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 672e477d84eSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 6736636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 6740a3ada39SMatthew G. Knepley ierr = DMDAVecGetClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr); 6758865f1eaSKarl Rupp for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]); 67657459e9aSMatthew G Knepley switch (dim) { 67757459e9aSMatthew G Knepley case 2: 6789c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &quadPoints, NULL);CHKERRQ(ierr); 6799c3cf19fSMatthew G. Knepley for (q = 0; q < Nq; ++q) { 68021454ff5SMatthew G. Knepley ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr); 68157459e9aSMatthew G Knepley } 68257459e9aSMatthew G Knepley break; 68357459e9aSMatthew G Knepley default: 68482f516ccSBarry Smith SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension %d not supported", dim); 68557459e9aSMatthew G Knepley } 6860a3ada39SMatthew G. Knepley ierr = DMDAVecRestoreClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr); 68757459e9aSMatthew G Knepley PetscFunctionReturn(0); 68857459e9aSMatthew G Knepley } 68905264a50SDave May 69005264a50SDave May PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular,Vec pos,IS *iscell) 69105264a50SDave May { 69205264a50SDave May PetscInt n,bs,p,npoints; 69305264a50SDave May PetscInt xs,xe,Xs,Xe,mxlocal; 69405264a50SDave May PetscInt ys,ye,Ys,Ye,mylocal; 695*aab5bcd8SJed Brown PetscInt d,c0,c1; 69605264a50SDave May PetscReal gmin_l[2],gmax_l[2],dx[2]; 69705264a50SDave May PetscReal gmin[2],gmax[2]; 69805264a50SDave May PetscInt *cellidx; 69905264a50SDave May Vec coor; 70005264a50SDave May const PetscScalar *_coor; 70105264a50SDave May PetscErrorCode ierr; 70205264a50SDave May 703*aab5bcd8SJed Brown PetscFunctionBegin; 70405264a50SDave May ierr = DMDAGetCorners(dmregular,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr); 70505264a50SDave May ierr = DMDAGetGhostCorners(dmregular,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr); 70605264a50SDave May xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 70705264a50SDave May ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 70805264a50SDave May 70905264a50SDave May mxlocal = xe - xs - 1; 71005264a50SDave May mylocal = ye - ys - 1; 71105264a50SDave May 71205264a50SDave May ierr = DMGetCoordinatesLocal(dmregular,&coor);CHKERRQ(ierr); 71305264a50SDave May ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr); 71405264a50SDave May c0 = (xs-Xs) + (ys-Ys)*(Xe-Xs); 71505264a50SDave May c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs); 71605264a50SDave May 717a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[2*c0+0]); 718a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[2*c0+1]); 71905264a50SDave May 720a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[2*c1+0]); 721a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[2*c1+1]); 72205264a50SDave May 72305264a50SDave May dx[0] = (gmax_l[0]-gmin_l[0])/((PetscReal)mxlocal); 72405264a50SDave May dx[1] = (gmax_l[1]-gmin_l[1])/((PetscReal)mylocal); 72505264a50SDave May 72605264a50SDave May ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr); 72705264a50SDave May 72805264a50SDave May ierr = DMDAGetBoundingBox(dmregular,gmin,gmax);CHKERRQ(ierr); 72905264a50SDave May 73005264a50SDave May ierr = VecGetLocalSize(pos,&n);CHKERRQ(ierr); 73105264a50SDave May ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr); 73205264a50SDave May npoints = n/bs; 73305264a50SDave May 73405264a50SDave May ierr = PetscMalloc1(npoints,&cellidx);CHKERRQ(ierr); 73505264a50SDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 73605264a50SDave May for (p=0; p<npoints; p++) { 7370ca99263SDave May PetscReal coor_p[2]; 73805264a50SDave May PetscInt mi[2]; 73905264a50SDave May 7400ca99263SDave May coor_p[0] = PetscRealPart(_coor[2*p]); 7410ca99263SDave May coor_p[1] = PetscRealPart(_coor[2*p+1]); 74205264a50SDave May 74305264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 74405264a50SDave May 7450ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 7460ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 7470ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 7480ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 74905264a50SDave May 750*aab5bcd8SJed Brown for (d=0; d<2; d++) { 75105264a50SDave May mi[d] = (PetscInt)( (coor_p[d] - gmin[d])/dx[d] ); 75205264a50SDave May } 75305264a50SDave May 75405264a50SDave May if (mi[0] < xs) { continue; } 75505264a50SDave May if (mi[0] > (xe-1)) { continue; } 75605264a50SDave May if (mi[1] < ys) { continue; } 75705264a50SDave May if (mi[1] > (ye-1)) { continue; } 75805264a50SDave May 75905264a50SDave May if (mi[0] == (xe-1)) { mi[0]--; } 76005264a50SDave May if (mi[1] == (ye-1)) { mi[1]--; } 76105264a50SDave May 76205264a50SDave May cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal; 76305264a50SDave May } 76405264a50SDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 76505264a50SDave May ierr = ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell);CHKERRQ(ierr); 76605264a50SDave May PetscFunctionReturn(0); 76705264a50SDave May } 76805264a50SDave May 76905264a50SDave May PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular,Vec pos,IS *iscell) 77005264a50SDave May { 77105264a50SDave May PetscInt n,bs,p,npoints; 77205264a50SDave May PetscInt xs,xe,Xs,Xe,mxlocal; 77305264a50SDave May PetscInt ys,ye,Ys,Ye,mylocal; 77405264a50SDave May PetscInt zs,ze,Zs,Ze,mzlocal; 775*aab5bcd8SJed Brown PetscInt d,c0,c1; 77605264a50SDave May PetscReal gmin_l[3],gmax_l[3],dx[3]; 77705264a50SDave May PetscReal gmin[3],gmax[3]; 77805264a50SDave May PetscInt *cellidx; 77905264a50SDave May Vec coor; 78005264a50SDave May const PetscScalar *_coor; 78105264a50SDave May PetscErrorCode ierr; 78205264a50SDave May 783*aab5bcd8SJed Brown PetscFunctionBegin; 78405264a50SDave May ierr = DMDAGetCorners(dmregular,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr); 78505264a50SDave May ierr = DMDAGetGhostCorners(dmregular,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr); 78605264a50SDave May xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 78705264a50SDave May ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 78805264a50SDave May ze += zs; Ze += Zs; if (zs != Zs) zs -= 1; 78905264a50SDave May 79005264a50SDave May mxlocal = xe - xs - 1; 79105264a50SDave May mylocal = ye - ys - 1; 79205264a50SDave May mzlocal = ze - zs - 1; 79305264a50SDave May 79405264a50SDave May ierr = DMGetCoordinatesLocal(dmregular,&coor);CHKERRQ(ierr); 79505264a50SDave May ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr); 79605264a50SDave May c0 = (xs-Xs) + (ys-Ys) *(Xe-Xs) + (zs-Zs) *(Xe-Xs)*(Ye-Ys); 79705264a50SDave May c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs) + (ze-2-Zs+1)*(Xe-Xs)*(Ye-Ys); 79805264a50SDave May 799a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[3*c0+0]); 800a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[3*c0+1]); 801a5f152d1SDave May gmin_l[2] = PetscRealPart(_coor[3*c0+2]); 80205264a50SDave May 803a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[3*c1+0]); 804a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[3*c1+1]); 805a5f152d1SDave May gmax_l[2] = PetscRealPart(_coor[3*c1+2]); 80605264a50SDave May 80705264a50SDave May dx[0] = (gmax_l[0]-gmin_l[0])/((PetscReal)mxlocal); 80805264a50SDave May dx[1] = (gmax_l[1]-gmin_l[1])/((PetscReal)mylocal); 80905264a50SDave May dx[2] = (gmax_l[2]-gmin_l[2])/((PetscReal)mzlocal); 81005264a50SDave May 81105264a50SDave May ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr); 81205264a50SDave May 81305264a50SDave May ierr = DMDAGetBoundingBox(dmregular,gmin,gmax);CHKERRQ(ierr); 81405264a50SDave May 81505264a50SDave May ierr = VecGetLocalSize(pos,&n);CHKERRQ(ierr); 81605264a50SDave May ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr); 81705264a50SDave May npoints = n/bs; 81805264a50SDave May 81905264a50SDave May ierr = PetscMalloc1(npoints,&cellidx);CHKERRQ(ierr); 82005264a50SDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 82105264a50SDave May for (p=0; p<npoints; p++) { 8220ca99263SDave May PetscReal coor_p[3]; 82305264a50SDave May PetscInt mi[3]; 82405264a50SDave May 8250ca99263SDave May coor_p[0] = PetscRealPart(_coor[3*p]); 8260ca99263SDave May coor_p[1] = PetscRealPart(_coor[3*p+1]); 8270ca99263SDave May coor_p[2] = PetscRealPart(_coor[3*p+2]); 82805264a50SDave May 82905264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 83005264a50SDave May 8310ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 8320ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 8330ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 8340ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 8350ca99263SDave May if (coor_p[2] < gmin_l[2]) { continue; } 8360ca99263SDave May if (coor_p[2] > gmax_l[2]) { continue; } 83705264a50SDave May 838*aab5bcd8SJed Brown for (d=0; d<3; d++) { 83905264a50SDave May mi[d] = (PetscInt)( (coor_p[d] - gmin[d])/dx[d] ); 84005264a50SDave May } 84105264a50SDave May 84205264a50SDave May if (mi[0] < xs) { continue; } 84305264a50SDave May if (mi[0] > (xe-1)) { continue; } 84405264a50SDave May if (mi[1] < ys) { continue; } 84505264a50SDave May if (mi[1] > (ye-1)) { continue; } 84605264a50SDave May if (mi[2] < zs) { continue; } 84705264a50SDave May if (mi[2] > (ze-1)) { continue; } 84805264a50SDave May 84905264a50SDave May if (mi[0] == (xe-1)) { mi[0]--; } 85005264a50SDave May if (mi[1] == (ye-1)) { mi[1]--; } 85105264a50SDave May if (mi[2] == (ze-1)) { mi[2]--; } 85205264a50SDave May 85305264a50SDave May cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal + (mi[2]-zs) * mxlocal * mylocal; 85405264a50SDave May } 85505264a50SDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 85605264a50SDave May ierr = ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell);CHKERRQ(ierr); 85705264a50SDave May PetscFunctionReturn(0); 85805264a50SDave May } 85905264a50SDave May 86005264a50SDave May PetscErrorCode DMLocatePoints_DA_Regular(DM dm,Vec pos,DMPointLocationType ltype,PetscSF cellSF) 86105264a50SDave May { 86205264a50SDave May IS iscell; 86305264a50SDave May PetscSFNode *cells; 86405264a50SDave May PetscInt p,bs,dim,npoints,nfound; 86505264a50SDave May const PetscInt *boxCells; 86605264a50SDave May PetscErrorCode ierr; 86705264a50SDave May 86805264a50SDave May PetscFunctionBegin; 86905264a50SDave May ierr = VecGetBlockSize(pos,&dim);CHKERRQ(ierr); 87005264a50SDave May switch (dim) { 87105264a50SDave May case 1: 87205264a50SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support not provided for 1D"); 87305264a50SDave May break; 87405264a50SDave May case 2: 87505264a50SDave May ierr = private_DMDALocatePointsIS_2D_Regular(dm,pos,&iscell);CHKERRQ(ierr); 87605264a50SDave May break; 87705264a50SDave May case 3: 87805264a50SDave May ierr = private_DMDALocatePointsIS_3D_Regular(dm,pos,&iscell);CHKERRQ(ierr); 87905264a50SDave May break; 88005264a50SDave May default: 88105264a50SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupport spatial dimension"); 88205264a50SDave May break; 88305264a50SDave May } 88405264a50SDave May 88505264a50SDave May ierr = VecGetLocalSize(pos,&npoints);CHKERRQ(ierr); 88605264a50SDave May ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr); 88705264a50SDave May npoints = npoints / bs; 88805264a50SDave May 88905264a50SDave May ierr = PetscMalloc1(npoints, &cells);CHKERRQ(ierr); 89005264a50SDave May ierr = ISGetIndices(iscell, &boxCells);CHKERRQ(ierr); 89105264a50SDave May 89205264a50SDave May for (p=0; p<npoints; p++) { 89305264a50SDave May cells[p].rank = 0; 89405264a50SDave May cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND; 89505264a50SDave May 89605264a50SDave May cells[p].index = boxCells[p]; 89705264a50SDave May } 89805264a50SDave May ierr = ISRestoreIndices(iscell, &boxCells);CHKERRQ(ierr); 89905264a50SDave May 90005264a50SDave May nfound = npoints; 90105264a50SDave May ierr = PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr); 90205264a50SDave May ierr = ISDestroy(&iscell);CHKERRQ(ierr); 90305264a50SDave May 90405264a50SDave May PetscFunctionReturn(0); 90505264a50SDave May } 906