147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7b2fff234SMatthew G. Knepley #include <petscbt.h> 80c312b8eSJed Brown #include <petscsf.h> 937b16e26SMatthew G. Knepley #include <petscfe.h> 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith /* 12e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions 1347c6ae99SBarry Smith */ 1447c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 15c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */ 16c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */ 1747c6ae99SBarry Smith #undef __FUNCT__ 1847c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d" 191e6b0712SBarry Smith static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 2047c6ae99SBarry Smith { 2147c6ae99SBarry Smith PetscErrorCode ierr; 2247c6ae99SBarry Smith PetscInt n,m; 2347c6ae99SBarry Smith Vec vec = (Vec)obj; 2447c6ae99SBarry Smith PetscScalar *array; 2547c6ae99SBarry Smith mxArray *mat; 269a42bb27SBarry Smith DM da; 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith PetscFunctionBegin; 29c688c046SMatthew G Knepley ierr = VecGetDM(vec, &da);CHKERRQ(ierr); 30ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 31aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 3247c6ae99SBarry Smith 3347c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 3447c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3547c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3647c6ae99SBarry Smith #else 3747c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3847c6ae99SBarry Smith #endif 3947c6ae99SBarry Smith ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 4047c6ae99SBarry Smith ierr = PetscObjectName(obj);CHKERRQ(ierr); 4147c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 4247c6ae99SBarry Smith 4347c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 4447c6ae99SBarry Smith PetscFunctionReturn(0); 4547c6ae99SBarry Smith } 4647c6ae99SBarry Smith #endif 4747c6ae99SBarry Smith 4847c6ae99SBarry Smith 4947c6ae99SBarry Smith #undef __FUNCT__ 50564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA" 517087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 5247c6ae99SBarry Smith { 5347c6ae99SBarry Smith PetscErrorCode ierr; 5447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5547c6ae99SBarry Smith 5647c6ae99SBarry Smith PetscFunctionBegin; 5747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5847c6ae99SBarry Smith PetscValidPointer(g,2); 5911689aebSJed Brown if (da->defaultSection) { 6011689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 6111689aebSJed Brown } else { 6247c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6347c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6447c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 65401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 66c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6747c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 6847c6ae99SBarry Smith if (dd->w == 1 && dd->dim == 2) { 69bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith #endif 7211689aebSJed Brown } 7347c6ae99SBarry Smith PetscFunctionReturn(0); 7447c6ae99SBarry Smith } 7547c6ae99SBarry Smith 76a66d4d66SMatthew G Knepley #undef __FUNCT__ 7757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells" 78*939f6067SMatthew G. Knepley /*@ 79*939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 80*939f6067SMatthew G. Knepley 81*939f6067SMatthew G. Knepley Input Parameter: 82*939f6067SMatthew G. Knepley . dm - The DM object 83*939f6067SMatthew G. Knepley 84*939f6067SMatthew G. Knepley Output Parameters: 85*939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 86*939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 87*939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 88*939f6067SMatthew G. Knepley - numCells - The number of local cells 89*939f6067SMatthew G. Knepley 90*939f6067SMatthew G. Knepley Level: developer 91*939f6067SMatthew G. Knepley 92*939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint() 93*939f6067SMatthew G. Knepley @*/ 94e42e3c58SMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 9557459e9aSMatthew G Knepley { 9657459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 9757459e9aSMatthew G Knepley const PetscInt dim = da->dim; 9857459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 9957459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 10057459e9aSMatthew G Knepley 10157459e9aSMatthew G Knepley PetscFunctionBegin; 102e42e3c58SMatthew G. Knepley if (numCellsX) { 103e42e3c58SMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 104e42e3c58SMatthew G. Knepley *numCellsX = mx; 105e42e3c58SMatthew G. Knepley } 106e42e3c58SMatthew G. Knepley if (numCellsY) { 107e42e3c58SMatthew G. Knepley PetscValidIntPointer(numCellsX,3); 108e42e3c58SMatthew G. Knepley *numCellsY = my; 109e42e3c58SMatthew G. Knepley } 110e42e3c58SMatthew G. Knepley if (numCellsZ) { 111e42e3c58SMatthew G. Knepley PetscValidIntPointer(numCellsX,4); 112e42e3c58SMatthew G. Knepley *numCellsZ = mz; 113e42e3c58SMatthew G. Knepley } 11457459e9aSMatthew G Knepley if (numCells) { 115e42e3c58SMatthew G. Knepley PetscValidIntPointer(numCells,5); 11657459e9aSMatthew G Knepley *numCells = nC; 11757459e9aSMatthew G Knepley } 11857459e9aSMatthew G Knepley PetscFunctionReturn(0); 11957459e9aSMatthew G Knepley } 12057459e9aSMatthew G Knepley 12157459e9aSMatthew G Knepley #undef __FUNCT__ 12257459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices" 12357459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 12457459e9aSMatthew G Knepley { 12557459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 12657459e9aSMatthew G Knepley const PetscInt dim = da->dim; 12757459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 12857459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 12957459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 13057459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 13157459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 13257459e9aSMatthew G Knepley 13357459e9aSMatthew G Knepley PetscFunctionBegin; 13457459e9aSMatthew G Knepley if (numVerticesX) { 13557459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 13657459e9aSMatthew G Knepley *numVerticesX = nVx; 13757459e9aSMatthew G Knepley } 13857459e9aSMatthew G Knepley if (numVerticesY) { 13957459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 14057459e9aSMatthew G Knepley *numVerticesY = nVy; 14157459e9aSMatthew G Knepley } 14257459e9aSMatthew G Knepley if (numVerticesZ) { 14357459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 14457459e9aSMatthew G Knepley *numVerticesZ = nVz; 14557459e9aSMatthew G Knepley } 14657459e9aSMatthew G Knepley if (numVertices) { 14757459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 14857459e9aSMatthew G Knepley *numVertices = nV; 14957459e9aSMatthew G Knepley } 15057459e9aSMatthew G Knepley PetscFunctionReturn(0); 15157459e9aSMatthew G Knepley } 15257459e9aSMatthew G Knepley 15357459e9aSMatthew G Knepley #undef __FUNCT__ 15457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces" 15557459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 15657459e9aSMatthew G Knepley { 15757459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 15857459e9aSMatthew G Knepley const PetscInt dim = da->dim; 15957459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 16057459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 16157459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 16257459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 16357459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 16457459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 16557459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 16657459e9aSMatthew G Knepley 16757459e9aSMatthew G Knepley PetscFunctionBegin; 16857459e9aSMatthew G Knepley if (numXFacesX) { 16957459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 17057459e9aSMatthew G Knepley *numXFacesX = nxF; 17157459e9aSMatthew G Knepley } 17257459e9aSMatthew G Knepley if (numXFaces) { 17357459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 17457459e9aSMatthew G Knepley *numXFaces = nXF; 17557459e9aSMatthew G Knepley } 17657459e9aSMatthew G Knepley if (numYFacesY) { 17757459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 17857459e9aSMatthew G Knepley *numYFacesY = nyF; 17957459e9aSMatthew G Knepley } 18057459e9aSMatthew G Knepley if (numYFaces) { 18157459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 18257459e9aSMatthew G Knepley *numYFaces = nYF; 18357459e9aSMatthew G Knepley } 18457459e9aSMatthew G Knepley if (numZFacesZ) { 18557459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 18657459e9aSMatthew G Knepley *numZFacesZ = nzF; 18757459e9aSMatthew G Knepley } 18857459e9aSMatthew G Knepley if (numZFaces) { 18957459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 19057459e9aSMatthew G Knepley *numZFaces = nZF; 19157459e9aSMatthew G Knepley } 19257459e9aSMatthew G Knepley PetscFunctionReturn(0); 19357459e9aSMatthew G Knepley } 19457459e9aSMatthew G Knepley 19557459e9aSMatthew G Knepley #undef __FUNCT__ 19657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum" 19757459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 19857459e9aSMatthew G Knepley { 19957459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 20057459e9aSMatthew G Knepley const PetscInt dim = da->dim; 20157459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 20257459e9aSMatthew G Knepley PetscErrorCode ierr; 20357459e9aSMatthew G Knepley 20457459e9aSMatthew G Knepley PetscFunctionBegin; 2058865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 2068865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 207e42e3c58SMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2080298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2090298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 21057459e9aSMatthew G Knepley if (height == 0) { 21157459e9aSMatthew G Knepley /* Cells */ 2128865f1eaSKarl Rupp if (pStart) *pStart = 0; 2138865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 21457459e9aSMatthew G Knepley } else if (height == 1) { 21557459e9aSMatthew G Knepley /* Faces */ 2168865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 2178865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 21857459e9aSMatthew G Knepley } else if (height == dim) { 21957459e9aSMatthew G Knepley /* Vertices */ 2208865f1eaSKarl Rupp if (pStart) *pStart = nC; 2218865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 22257459e9aSMatthew G Knepley } else if (height < 0) { 22357459e9aSMatthew G Knepley /* All points */ 2248865f1eaSKarl Rupp if (pStart) *pStart = 0; 2258865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 22682f516ccSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 22757459e9aSMatthew G Knepley PetscFunctionReturn(0); 22857459e9aSMatthew G Knepley } 22957459e9aSMatthew G Knepley 23057459e9aSMatthew G Knepley #undef __FUNCT__ 231d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum" 232d0892670SMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) 233d0892670SMatthew G. Knepley { 234d0892670SMatthew G. Knepley DM_DA *da = (DM_DA*) dm->data; 235d0892670SMatthew G. Knepley const PetscInt dim = da->dim; 236d0892670SMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 237d0892670SMatthew G. Knepley PetscErrorCode ierr; 238d0892670SMatthew G. Knepley 239d0892670SMatthew G. Knepley PetscFunctionBegin; 240d0892670SMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart,3); 241d0892670SMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd,4); 242d0892670SMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 243d0892670SMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 244d0892670SMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 245d0892670SMatthew G. Knepley if (depth == dim) { 246d0892670SMatthew G. Knepley /* Cells */ 247d0892670SMatthew G. Knepley if (pStart) *pStart = 0; 248d0892670SMatthew G. Knepley if (pEnd) *pEnd = nC; 249d0892670SMatthew G. Knepley } else if (depth == dim-1) { 250d0892670SMatthew G. Knepley /* Faces */ 251d0892670SMatthew G. Knepley if (pStart) *pStart = nC+nV; 252d0892670SMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 253d0892670SMatthew G. Knepley } else if (depth == 0) { 254d0892670SMatthew G. Knepley /* Vertices */ 255d0892670SMatthew G. Knepley if (pStart) *pStart = nC; 256d0892670SMatthew G. Knepley if (pEnd) *pEnd = nC+nV; 257d0892670SMatthew G. Knepley } else if (depth < 0) { 258d0892670SMatthew G. Knepley /* All points */ 259d0892670SMatthew G. Knepley if (pStart) *pStart = 0; 260d0892670SMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 261d0892670SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth); 262d0892670SMatthew G. Knepley PetscFunctionReturn(0); 263d0892670SMatthew G. Knepley } 264d0892670SMatthew G. Knepley 265d0892670SMatthew G. Knepley #undef __FUNCT__ 266d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize" 267d0892670SMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) 268d0892670SMatthew G. Knepley { 269d0892670SMatthew G. Knepley DM_DA *da = (DM_DA*) dm->data; 270d0892670SMatthew G. Knepley const PetscInt dim = da->dim; 271d0892670SMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 272d0892670SMatthew G. Knepley PetscErrorCode ierr; 273d0892670SMatthew G. Knepley 274d0892670SMatthew G. Knepley PetscFunctionBegin; 275d0892670SMatthew G. Knepley *coneSize = 0; 276d0892670SMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 277d0892670SMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 278d0892670SMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 279d0892670SMatthew G. Knepley switch (dim) { 280d0892670SMatthew G. Knepley case 2: 281d0892670SMatthew G. Knepley if (p >= 0) { 282d0892670SMatthew G. Knepley if (p < nC) { 283d0892670SMatthew G. Knepley *coneSize = 4; 284d0892670SMatthew G. Knepley } else if (p < nC+nV) { 285d0892670SMatthew G. Knepley *coneSize = 0; 286d0892670SMatthew G. Knepley } else if (p < nC+nV+nXF+nYF+nZF) { 287d0892670SMatthew G. Knepley *coneSize = 2; 288d0892670SMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 289d0892670SMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 290d0892670SMatthew G. Knepley break; 291d0892670SMatthew G. Knepley case 3: 292d0892670SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 293d0892670SMatthew G. Knepley break; 294d0892670SMatthew G. Knepley } 295d0892670SMatthew G. Knepley PetscFunctionReturn(0); 296d0892670SMatthew G. Knepley } 297d0892670SMatthew G. Knepley 298d0892670SMatthew G. Knepley #undef __FUNCT__ 299d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetCone" 300d0892670SMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) 301d0892670SMatthew G. Knepley { 302d0892670SMatthew G. Knepley DM_DA *da = (DM_DA*) dm->data; 303d0892670SMatthew G. Knepley const PetscInt dim = da->dim; 304d0892670SMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 305d0892670SMatthew G. Knepley PetscErrorCode ierr; 306d0892670SMatthew G. Knepley 307d0892670SMatthew G. Knepley PetscFunctionBegin; 308d0892670SMatthew G. Knepley if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);} 309d0892670SMatthew G. Knepley ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr); 310d0892670SMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 311d0892670SMatthew G. Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 312d0892670SMatthew G. Knepley switch (dim) { 313d0892670SMatthew G. Knepley case 2: 314d0892670SMatthew G. Knepley if (p >= 0) { 315d0892670SMatthew G. Knepley if (p < nC) { 316d0892670SMatthew G. Knepley const PetscInt cy = p / nCx; 317d0892670SMatthew G. Knepley const PetscInt cx = p % nCx; 318d0892670SMatthew G. Knepley 319d0892670SMatthew G. Knepley (*cone)[0] = cy*nxF + cx + nC+nV; 320d0892670SMatthew G. Knepley (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF; 321d0892670SMatthew G. Knepley (*cone)[2] = cy*nxF + cx + nxF + nC+nV; 322d0892670SMatthew G. Knepley (*cone)[3] = cx*nyF + cy + nC+nV+nXF; 323d0892670SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 324d0892670SMatthew G. Knepley } else if (p < nC+nV) { 325d0892670SMatthew G. Knepley } else if (p < nC+nV+nXF) { 326d0892670SMatthew G. Knepley const PetscInt fy = (p - nC+nV) / nxF; 327d0892670SMatthew G. Knepley const PetscInt fx = (p - nC+nV) % nxF; 328d0892670SMatthew G. Knepley 329d0892670SMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 330d0892670SMatthew G. Knepley (*cone)[1] = fy*nVx + fx + 1 + nC; 331d0892670SMatthew G. Knepley } else if (p < nC+nV+nXF+nYF) { 332d0892670SMatthew G. Knepley const PetscInt fx = (p - nC+nV+nXF) / nyF; 333d0892670SMatthew G. Knepley const PetscInt fy = (p - nC+nV+nXF) % nyF; 334d0892670SMatthew G. Knepley 335d0892670SMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 336d0892670SMatthew G. Knepley (*cone)[1] = fy*nVx + fx + nVx + nC; 337d0892670SMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 338d0892670SMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 339d0892670SMatthew G. Knepley break; 340d0892670SMatthew G. Knepley case 3: 341d0892670SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 342d0892670SMatthew G. Knepley break; 343d0892670SMatthew G. Knepley } 344d0892670SMatthew G. Knepley PetscFunctionReturn(0); 345d0892670SMatthew G. Knepley } 346d0892670SMatthew G. Knepley 347d0892670SMatthew G. Knepley #undef __FUNCT__ 348d0892670SMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone" 349d0892670SMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) 350d0892670SMatthew G. Knepley { 351d0892670SMatthew G. Knepley PetscErrorCode ierr; 352d0892670SMatthew G. Knepley 353d0892670SMatthew G. Knepley PetscFunctionBegin; 354d0892670SMatthew G. Knepley ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr); 355d0892670SMatthew G. Knepley PetscFunctionReturn(0); 356d0892670SMatthew G. Knepley } 357d0892670SMatthew G. Knepley 358d0892670SMatthew G. Knepley #undef __FUNCT__ 359a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection" 360a66d4d66SMatthew G Knepley /*@C 361a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 362a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 363a66d4d66SMatthew G Knepley 364a66d4d66SMatthew G Knepley Input Parameters: 365a66d4d66SMatthew G Knepley + dm- The DMDA 366a66d4d66SMatthew G Knepley . numFields - The number of fields 367a4294c51SMatthew G. Knepley . numComp - The number of components in each field 368a4294c51SMatthew G. Knepley . numDof - The number of dofs per dimension for each field 3690298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL 370a66d4d66SMatthew G Knepley 371a66d4d66SMatthew G Knepley Level: developer 372a66d4d66SMatthew G Knepley 373a66d4d66SMatthew G Knepley Note: 374a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 375a66d4d66SMatthew G Knepley 376a66d4d66SMatthew G Knepley - Cells: [0, nC) 377a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 37888ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 37988ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 38088ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 381a66d4d66SMatthew G Knepley 382a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 383a66d4d66SMatthew G Knepley @*/ 384a4294c51SMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numDof[], PetscInt numFaceDof[], PetscSection *s) 385a66d4d66SMatthew G Knepley { 386a66d4d66SMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 387a4b60ecfSMatthew G. Knepley PetscSection section; 38888ed4aceSMatthew G Knepley const PetscInt dim = da->dim; 38980800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 390b2fff234SMatthew G. Knepley PetscBT isLeaf; 39188ed4aceSMatthew G Knepley PetscSF sf; 392feafbc5cSMatthew G Knepley PetscMPIInt rank; 39388ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 39488ed4aceSMatthew G Knepley PetscInt *localPoints; 39588ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 396f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 39757459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 39857459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 39988ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 400a66d4d66SMatthew G Knepley PetscErrorCode ierr; 401a66d4d66SMatthew G Knepley 402a66d4d66SMatthew G Knepley PetscFunctionBegin; 403a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 404e42e3c58SMatthew G. Knepley PetscValidPointer(s, 4); 40582f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 406e42e3c58SMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 40757459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 40857459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 40957459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 41057459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 41157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 41257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 41357459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 41457459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 41557459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 41682f516ccSBarry Smith if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 41788ed4aceSMatthew G Knepley /* Create local section */ 41880800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 419a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 420a4294c51SMatthew G. Knepley numVertexTotDof += numDof[f*(dim+1)+0]; 421a4294c51SMatthew G. Knepley numCellTotDof += numDof[f*(dim+1)+dim]; 422a4294c51SMatthew G. Knepley numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0; 423a4294c51SMatthew G. Knepley numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0; 424a4294c51SMatthew G. Knepley numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0; 425a66d4d66SMatthew G Knepley } 426a4b60ecfSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);CHKERRQ(ierr); 427a4294c51SMatthew G. Knepley if (numFields > 0) { 428a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr); 429a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 43080800b1aSMatthew G Knepley const char *name; 43180800b1aSMatthew G Knepley 43280800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 433a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr); 43480800b1aSMatthew G Knepley if (numComp) { 435a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr); 436a66d4d66SMatthew G Knepley } 437a66d4d66SMatthew G Knepley } 438a66d4d66SMatthew G Knepley } 439a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 440a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 441a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 442a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr); 443a66d4d66SMatthew G Knepley } 444a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr); 445a66d4d66SMatthew G Knepley } 446a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 447a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 448a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 449a66d4d66SMatthew G Knepley } 450a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr); 451a66d4d66SMatthew G Knepley } 452a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 453a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 454a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 455a66d4d66SMatthew G Knepley } 456a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr); 457a66d4d66SMatthew G Knepley } 458a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 459a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 460a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 461a66d4d66SMatthew G Knepley } 462a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr); 463a66d4d66SMatthew G Knepley } 464a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 465a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 466a4294c51SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr); 467a66d4d66SMatthew G Knepley } 468a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr); 469a66d4d66SMatthew G Knepley } 470a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 47188ed4aceSMatthew G Knepley /* Create mesh point SF */ 472b2fff234SMatthew G. Knepley ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr); 47388ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 47488ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 47588ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 47688ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 4777128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 47888ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 479b2fff234SMatthew G. Knepley PetscInt xv, yv, zv; 48088ed4aceSMatthew G Knepley 4813814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 482b2fff234SMatthew G. Knepley if (xp < 0) { /* left */ 483b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 484b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 485b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 486b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 487b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 488b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 489b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 490b2fff234SMatthew G. Knepley } else { 491b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 492b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 493b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 494b2fff234SMatthew G. Knepley } 495b2fff234SMatthew G. Knepley } 496b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 497b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 498b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 499b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 500b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 501b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 502b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 503b2fff234SMatthew G. Knepley } else { 504b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 505b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 506b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 507b2fff234SMatthew G. Knepley } 508b2fff234SMatthew G. Knepley } 509b2fff234SMatthew G. Knepley } else { 510b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 511b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 512b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 513b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 514b2fff234SMatthew G. Knepley } 515b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 516b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 517b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 518b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 519b2fff234SMatthew G. Knepley } 520b2fff234SMatthew G. Knepley } else { 521b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 522b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 523b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 524b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 525b2fff234SMatthew G. Knepley } 526b2fff234SMatthew G. Knepley } 527b2fff234SMatthew G. Knepley #if 0 528b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 529b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 530b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 531b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 532b2fff234SMatthew G. Knepley } 533b2fff234SMatthew G. Knepley #endif 534b2fff234SMatthew G. Knepley } 535b2fff234SMatthew G. Knepley } 536b2fff234SMatthew G. Knepley } else if (xp > 0) { /* right */ 537b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 538b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 539b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 540b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 541b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 542b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 543b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 544b2fff234SMatthew G. Knepley } else { 545b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 546b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 547b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 548b2fff234SMatthew G. Knepley } 549b2fff234SMatthew G. Knepley } 550b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 551b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 552b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 553b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 554b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 555b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 556b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 557b2fff234SMatthew G. Knepley } else { 558b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 559b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 560b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 561b2fff234SMatthew G. Knepley } 562b2fff234SMatthew G. Knepley } 563b2fff234SMatthew G. Knepley } else { 564b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 565b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 566b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 567b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 568b2fff234SMatthew G. Knepley } 569b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 570b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 571b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 572b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 573b2fff234SMatthew G. Knepley } 574b2fff234SMatthew G. Knepley } else { 575b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 576b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 577b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 578b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 579b2fff234SMatthew G. Knepley } 580b2fff234SMatthew G. Knepley } 581b2fff234SMatthew G. Knepley #if 0 582b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 583b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 584b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 585b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 586b2fff234SMatthew G. Knepley } 587b2fff234SMatthew G. Knepley #endif 588b2fff234SMatthew G. Knepley } 589b2fff234SMatthew G. Knepley } 590b2fff234SMatthew G. Knepley } else { 591b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 592b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 593b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 594b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 595b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 596b2fff234SMatthew G. Knepley } 597b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 598b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 599b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 600b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 601b2fff234SMatthew G. Knepley } 602b2fff234SMatthew G. Knepley } else { 603b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 604b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 605b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 606b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 607b2fff234SMatthew G. Knepley } 608b2fff234SMatthew G. Knepley } 609b2fff234SMatthew G. Knepley #if 0 610b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 611b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 612b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 613b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 614b2fff234SMatthew G. Knepley } 615b2fff234SMatthew G. Knepley #endif 616b2fff234SMatthew G. Knepley } 617b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 618b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 619b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 620b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 621b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 622b2fff234SMatthew G. Knepley } 623b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 624b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 625b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 626b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 627b2fff234SMatthew G. Knepley } 628b2fff234SMatthew G. Knepley } else { 629b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 630b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 631b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 632b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 633b2fff234SMatthew G. Knepley } 634b2fff234SMatthew G. Knepley } 635b2fff234SMatthew G. Knepley #if 0 636b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 637b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 638b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 639b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 640b2fff234SMatthew G. Knepley } 641b2fff234SMatthew G. Knepley #endif 642b2fff234SMatthew G. Knepley } 643b2fff234SMatthew G. Knepley } else { 644b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 645b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 646b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 647b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 648b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 649b2fff234SMatthew G. Knepley } 650b2fff234SMatthew G. Knepley } 651b2fff234SMatthew G. Knepley #if 0 652b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 653b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 654b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 655b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 656b2fff234SMatthew G. Knepley } 657b2fff234SMatthew G. Knepley #endif 658b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 659b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 660b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 661b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 662b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 663b2fff234SMatthew G. Knepley } 664b2fff234SMatthew G. Knepley } 665b2fff234SMatthew G. Knepley #if 0 666b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 667b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 668b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 669b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 670b2fff234SMatthew G. Knepley } 671b2fff234SMatthew G. Knepley #endif 672b2fff234SMatthew G. Knepley } else { 673b2fff234SMatthew G. Knepley /* Nothing is shared from the interior */ 67488ed4aceSMatthew G Knepley } 67588ed4aceSMatthew G Knepley } 67688ed4aceSMatthew G Knepley } 67788ed4aceSMatthew G Knepley } 67888ed4aceSMatthew G Knepley } 679b2fff234SMatthew G. Knepley } 680b2fff234SMatthew G. Knepley } 681b2fff234SMatthew G. Knepley ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr); 68288ed4aceSMatthew G Knepley ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr); 68388ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 68488ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 68588ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 6867128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 68788ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 688f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 68988ed4aceSMatthew G Knepley 6903814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 69188ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 69288ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 69388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 694b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 695f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 696f5eeb527SMatthew G Knepley 697b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 698f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 699f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 700f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 701f5eeb527SMatthew G Knepley ++nL; 702b2fff234SMatthew G. Knepley } 70388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 704b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 705f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 706f5eeb527SMatthew G Knepley 707b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 708f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 709f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 710f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 711f5eeb527SMatthew G Knepley ++nL; 712b2fff234SMatthew G. Knepley } 71388ed4aceSMatthew G Knepley } else { 714b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 715b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 716f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 717f5eeb527SMatthew G Knepley 718b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 719f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 720f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 721f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 722b2fff234SMatthew G. Knepley ++nL; 723b2fff234SMatthew G. Knepley } 724f5eeb527SMatthew G Knepley } 72588ed4aceSMatthew G Knepley } 72688ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 72788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 728b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 729f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 730f5eeb527SMatthew G Knepley 731b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 732f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 733f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 734f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 735f5eeb527SMatthew G Knepley ++nL; 736b2fff234SMatthew G. Knepley } 73788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 738b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 739f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 740f5eeb527SMatthew G Knepley 741b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 742f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 743f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 744f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 745f5eeb527SMatthew G Knepley ++nL; 746b2fff234SMatthew G. Knepley } 74788ed4aceSMatthew G Knepley } else { 748b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 749b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 750f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 751f5eeb527SMatthew G Knepley 752b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 753f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 754f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 755f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 756b2fff234SMatthew G. Knepley ++nL; 757b2fff234SMatthew G. Knepley } 758f5eeb527SMatthew G Knepley } 75988ed4aceSMatthew G Knepley } 76088ed4aceSMatthew G Knepley } else { 76188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 762b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 763b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 764f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 765f5eeb527SMatthew G Knepley 766b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 767f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 768f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 769f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 770b2fff234SMatthew G. Knepley ++nL; 771b2fff234SMatthew G. Knepley } 772f5eeb527SMatthew G Knepley } 77388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 774b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 775b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 776f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 777f5eeb527SMatthew G Knepley 778b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 779f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 780f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 781f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 782b2fff234SMatthew G. Knepley ++nL; 783b2fff234SMatthew G. Knepley } 784f5eeb527SMatthew G Knepley } 78588ed4aceSMatthew G Knepley } else { 786f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 787b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 788b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 789f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 790f5eeb527SMatthew G Knepley 791b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 792f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 793f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 794f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 795b2fff234SMatthew G. Knepley ++nL; 796f5eeb527SMatthew G Knepley } 797f5eeb527SMatthew G Knepley } 798b2fff234SMatthew G. Knepley } 799b2fff234SMatthew G. Knepley #if 0 800b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 801f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 802b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 803f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 804f5eeb527SMatthew G Knepley 805b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 806f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 807f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 808f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 809f5eeb527SMatthew G Knepley } 81088ed4aceSMatthew G Knepley } 811b2fff234SMatthew G. Knepley #endif 812b2fff234SMatthew G. Knepley } 81388ed4aceSMatthew G Knepley } 81488ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 81588ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 81688ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 817b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 818f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 819f5eeb527SMatthew G Knepley 820b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 821f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 822f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 823f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 824f5eeb527SMatthew G Knepley ++nL; 825b2fff234SMatthew G. Knepley } 82688ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 827b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 828f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 829f5eeb527SMatthew G Knepley 830b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 831f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 832f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 833f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 834f5eeb527SMatthew G Knepley ++nL; 835b2fff234SMatthew G. Knepley } 83688ed4aceSMatthew G Knepley } else { 837b2fff234SMatthew G. Knepley nleavesCheck += nVz; 838b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 839b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 840f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 841f5eeb527SMatthew G Knepley 842b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 843f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 844f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 845f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 846b2fff234SMatthew G. Knepley ++nL; 847b2fff234SMatthew G. Knepley } 848f5eeb527SMatthew G Knepley } 84988ed4aceSMatthew G Knepley } 85088ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 85188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 852b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 853f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 854f5eeb527SMatthew G Knepley 855b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 856f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 857f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 858f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 859f5eeb527SMatthew G Knepley ++nL; 860b2fff234SMatthew G. Knepley } 86188ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 862b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 863f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 864f5eeb527SMatthew G Knepley 865b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 866f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 867f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 868f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 869f5eeb527SMatthew G Knepley ++nL; 870b2fff234SMatthew G. Knepley } 87188ed4aceSMatthew G Knepley } else { 872b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 873b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 874f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 875f5eeb527SMatthew G Knepley 876b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 877f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 878f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 879f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 880b2fff234SMatthew G. Knepley ++nL; 881b2fff234SMatthew G. Knepley } 882f5eeb527SMatthew G Knepley } 88388ed4aceSMatthew G Knepley } 88488ed4aceSMatthew G Knepley } else { 88588ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 886b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 887b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 888f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 889f5eeb527SMatthew G Knepley 890b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 891f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 892f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 893f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 894b2fff234SMatthew G. Knepley ++nL; 895b2fff234SMatthew G. Knepley } 896f5eeb527SMatthew G Knepley } 89788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 898b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 899b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 900f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 901f5eeb527SMatthew G Knepley 902b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 903f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 904f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 905f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 906b2fff234SMatthew G. Knepley ++nL; 907b2fff234SMatthew G. Knepley } 908f5eeb527SMatthew G Knepley } 90988ed4aceSMatthew G Knepley } else { 910f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 911b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 912b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 913f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 914f5eeb527SMatthew G Knepley 915b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 916f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 917f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 918f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 919b2fff234SMatthew G. Knepley ++nL; 920f5eeb527SMatthew G Knepley } 921f5eeb527SMatthew G Knepley } 922b2fff234SMatthew G. Knepley } 923b2fff234SMatthew G. Knepley #if 0 924b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 925f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 926b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 927f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 928f5eeb527SMatthew G Knepley 929b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 930f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 931f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 932f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 933b2fff234SMatthew G. Knepley ++nL; 934f5eeb527SMatthew G Knepley } 93588ed4aceSMatthew G Knepley } 936b2fff234SMatthew G. Knepley #endif 937b2fff234SMatthew G. Knepley } 93888ed4aceSMatthew G Knepley } 93988ed4aceSMatthew G Knepley } else { 94088ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 94188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 942b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 943b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 944f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 945f5eeb527SMatthew G Knepley 946b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 947f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 948f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 949f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 950b2fff234SMatthew G. Knepley ++nL; 951b2fff234SMatthew G. Knepley } 952f5eeb527SMatthew G Knepley } 95388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 954b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 955b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 956f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 957f5eeb527SMatthew G Knepley 958b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 959f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 960f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 961f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 962b2fff234SMatthew G. Knepley ++nL; 963b2fff234SMatthew G. Knepley } 964f5eeb527SMatthew G Knepley } 96588ed4aceSMatthew G Knepley } else { 966f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 967b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 968b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 969f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 970f5eeb527SMatthew G Knepley 971b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 972f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 973f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 974f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 975b2fff234SMatthew G. Knepley ++nL; 976f5eeb527SMatthew G Knepley } 977f5eeb527SMatthew G Knepley } 978b2fff234SMatthew G. Knepley } 979b2fff234SMatthew G. Knepley #if 0 980b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 981f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 982b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 983f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 984f5eeb527SMatthew G Knepley 985b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 986f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 987f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 988f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 989b2fff234SMatthew G. Knepley ++nL; 990f5eeb527SMatthew G Knepley } 99188ed4aceSMatthew G Knepley } 992b2fff234SMatthew G. Knepley #endif 993b2fff234SMatthew G. Knepley } 99488ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 99588ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 996b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 997b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 998f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 999f5eeb527SMatthew G Knepley 1000b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1001f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1002f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1003f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1004b2fff234SMatthew G. Knepley ++nL; 1005b2fff234SMatthew G. Knepley } 1006f5eeb527SMatthew G Knepley } 100788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1008b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1009b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 1010f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1011f5eeb527SMatthew G Knepley 1012b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1013f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1014f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1015f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1016b2fff234SMatthew G. Knepley ++nL; 1017b2fff234SMatthew G. Knepley } 1018f5eeb527SMatthew G Knepley } 101988ed4aceSMatthew G Knepley } else { 1020f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 1021b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1022b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 1023f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1024f5eeb527SMatthew G Knepley 1025b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1026f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1027f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1028f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1029b2fff234SMatthew G. Knepley ++nL; 1030f5eeb527SMatthew G Knepley } 1031f5eeb527SMatthew G Knepley } 1032b2fff234SMatthew G. Knepley } 1033b2fff234SMatthew G. Knepley #if 0 1034b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 1035f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1036b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 1037f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1038f5eeb527SMatthew G Knepley 1039b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1040f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1041f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1042f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1043b2fff234SMatthew G. Knepley ++nL; 1044f5eeb527SMatthew G Knepley } 104588ed4aceSMatthew G Knepley } 1046b2fff234SMatthew G. Knepley #endif 1047b2fff234SMatthew G. Knepley } 104888ed4aceSMatthew G Knepley } else { 104988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 1050f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1051b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1052b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 1053f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1054f5eeb527SMatthew G Knepley 1055b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1056f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1057f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1058f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1059b2fff234SMatthew G. Knepley ++nL; 1060f5eeb527SMatthew G Knepley } 1061f5eeb527SMatthew G Knepley } 1062b2fff234SMatthew G. Knepley } 1063b2fff234SMatthew G. Knepley #if 0 1064b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1065f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1066b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 1067f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1068f5eeb527SMatthew G Knepley 1069b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1070f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1071f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1072f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1073b2fff234SMatthew G. Knepley ++nL; 1074f5eeb527SMatthew G Knepley } 1075b2fff234SMatthew G. Knepley } 1076b2fff234SMatthew G. Knepley #endif 107788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1078f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1079b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1080b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 1081f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1082f5eeb527SMatthew G Knepley 1083b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1084f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1085f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1086f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1087b2fff234SMatthew G. Knepley ++nL; 1088f5eeb527SMatthew G Knepley } 1089f5eeb527SMatthew G Knepley } 1090b2fff234SMatthew G. Knepley } 1091b2fff234SMatthew G. Knepley #if 0 1092b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1093f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1094b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 1095f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1096f5eeb527SMatthew G Knepley 1097b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1098f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1099f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1100f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1101b2fff234SMatthew G. Knepley ++nL; 1102f5eeb527SMatthew G Knepley } 1103b2fff234SMatthew G. Knepley } 1104b2fff234SMatthew G. Knepley #endif 110588ed4aceSMatthew G Knepley } else { 110688ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 110788ed4aceSMatthew G Knepley } 110888ed4aceSMatthew G Knepley } 110988ed4aceSMatthew G Knepley } 111088ed4aceSMatthew G Knepley } 111188ed4aceSMatthew G Knepley } 111288ed4aceSMatthew G Knepley } 111388ed4aceSMatthew G Knepley } 1114b2fff234SMatthew G. Knepley ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr); 1115b2fff234SMatthew G. Knepley /* Remove duplication in leaf determination */ 1116b2fff234SMatthew G. Knepley if (nleaves != nL) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck); 111782f516ccSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr); 111888ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1119a4b60ecfSMatthew G. Knepley ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 112088ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1121a4294c51SMatthew G. Knepley *s = section; 1122a4294c51SMatthew G. Knepley PetscFunctionReturn(0); 1123a4294c51SMatthew G. Knepley } 1124a4294c51SMatthew G. Knepley 1125d0892670SMatthew G. Knepley #undef __FUNCT__ 1126d0892670SMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates" 1127d0892670SMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) 1128d0892670SMatthew G. Knepley { 1129d0892670SMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 1130d0892670SMatthew G. Knepley Vec coordinates; 1131d0892670SMatthew G. Knepley PetscSection section; 1132d0892670SMatthew G. Knepley PetscScalar *coords; 1133d0892670SMatthew G. Knepley PetscReal h[3]; 1134d0892670SMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 1135d0892670SMatthew G. Knepley PetscErrorCode ierr; 1136d0892670SMatthew G. Knepley 1137d0892670SMatthew G. Knepley PetscFunctionBegin; 1138d0892670SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1139d0892670SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 1140d0892670SMatthew G. Knepley h[0] = (xu - xl)/M; 1141d0892670SMatthew G. Knepley h[1] = (yu - yl)/N; 1142d0892670SMatthew G. Knepley h[2] = (zu - zl)/P; 1143d0892670SMatthew G. Knepley ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1144d0892670SMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 1145d0892670SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 1146d0892670SMatthew G. Knepley ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr); 1147d0892670SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr); 1148d0892670SMatthew G. Knepley ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr); 1149d0892670SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 1150d0892670SMatthew G. Knepley ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr); 1151d0892670SMatthew G. Knepley } 1152d0892670SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 1153d0892670SMatthew G. Knepley ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 1154d0892670SMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr); 1155d0892670SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1156d0892670SMatthew G. Knepley for (k = 0; k < nVz; ++k) { 1157d0892670SMatthew G. Knepley PetscInt ind[3] = {0, 0, k + da->zs}, d, off; 1158d0892670SMatthew G. Knepley 1159d0892670SMatthew G. Knepley for (j = 0; j < nVy; ++j) { 1160d0892670SMatthew G. Knepley ind[1] = j + da->ys; 1161d0892670SMatthew G. Knepley for (i = 0; i < nVx; ++i) { 1162d0892670SMatthew G. Knepley const PetscInt vertex = (k*nVy + j)*nVx + i + vStart; 1163d0892670SMatthew G. Knepley 1164d0892670SMatthew G. Knepley ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr); 1165d0892670SMatthew G. Knepley ind[0] = i + da->xs; 1166d0892670SMatthew G. Knepley for (d = 0; d < dim; ++d) { 1167d0892670SMatthew G. Knepley coords[off+d] = h[d]*ind[d]; 1168d0892670SMatthew G. Knepley } 1169d0892670SMatthew G. Knepley } 1170d0892670SMatthew G. Knepley } 1171d0892670SMatthew G. Knepley } 1172d0892670SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1173d0892670SMatthew G. Knepley ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr); 1174d0892670SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1175a4b60ecfSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 1176d0892670SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1177d0892670SMatthew G. Knepley PetscFunctionReturn(0); 1178d0892670SMatthew G. Knepley } 117937b16e26SMatthew G. Knepley 118037b16e26SMatthew G. Knepley #undef __FUNCT__ 118137b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunctionLocal" 118237b16e26SMatthew G. Knepley PetscErrorCode DMDAProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX) 118337b16e26SMatthew G. Knepley { 118437b16e26SMatthew G. Knepley PetscDualSpace *sp; 118537b16e26SMatthew G. Knepley PetscQuadrature q; 118637b16e26SMatthew G. Knepley PetscSection section; 118737b16e26SMatthew G. Knepley PetscScalar *values; 118837b16e26SMatthew G. Knepley PetscReal *v0, *J, *detJ; 118937b16e26SMatthew G. Knepley PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, c, v, d; 119037b16e26SMatthew G. Knepley PetscErrorCode ierr; 119137b16e26SMatthew G. Knepley 119237b16e26SMatthew G. Knepley PetscFunctionBegin; 119337b16e26SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 119437b16e26SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr); 119537b16e26SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 119637b16e26SMatthew G. Knepley ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &sp);CHKERRQ(ierr); 119737b16e26SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 119837b16e26SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 119937b16e26SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 120037b16e26SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 120137b16e26SMatthew G. Knepley totDim += spDim*numComp; 120237b16e26SMatthew G. Knepley } 120337b16e26SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 120437b16e26SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 120537b16e26SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 120637b16e26SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 120737b16e26SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 120837b16e26SMatthew G. Knepley ierr = PetscMalloc3(dim*q.numPoints,PetscReal,&v0,dim*dim*q.numPoints,PetscReal,&J,q.numPoints,PetscReal,&detJ);CHKERRQ(ierr); 120937b16e26SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 121037b16e26SMatthew G. Knepley PetscCellGeometry geom; 121137b16e26SMatthew G. Knepley 121237b16e26SMatthew G. Knepley ierr = DMDAComputeCellGeometry(dm, c, &q, v0, J, NULL, detJ);CHKERRQ(ierr); 121337b16e26SMatthew G. Knepley geom.v0 = v0; 121437b16e26SMatthew G. Knepley geom.J = J; 121537b16e26SMatthew G. Knepley geom.detJ = detJ; 121637b16e26SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 121737b16e26SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 121837b16e26SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 121937b16e26SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 122037b16e26SMatthew G. Knepley ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr); 122137b16e26SMatthew G. Knepley v += numComp; 122237b16e26SMatthew G. Knepley } 122337b16e26SMatthew G. Knepley } 122437b16e26SMatthew G. Knepley ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 122537b16e26SMatthew G. Knepley } 122637b16e26SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 122737b16e26SMatthew G. Knepley ierr = PetscFree3(v0,J,detJ);CHKERRQ(ierr); 122837b16e26SMatthew G. Knepley ierr = PetscFree(sp);CHKERRQ(ierr); 122937b16e26SMatthew G. Knepley PetscFunctionReturn(0); 123037b16e26SMatthew G. Knepley } 123137b16e26SMatthew G. Knepley 123237b16e26SMatthew G. Knepley #undef __FUNCT__ 123337b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunction" 123437b16e26SMatthew G. Knepley /*@C 123537b16e26SMatthew G. Knepley DMDAProjectFunction - This projects the given function into the function space provided. 123637b16e26SMatthew G. Knepley 123737b16e26SMatthew G. Knepley Input Parameters: 123837b16e26SMatthew G. Knepley + dm - The DM 123937b16e26SMatthew G. Knepley . fe - The PetscFE associated with the field 124037b16e26SMatthew G. Knepley . funcs - The coordinate functions to evaluate 124137b16e26SMatthew G. Knepley - mode - The insertion mode for values 124237b16e26SMatthew G. Knepley 124337b16e26SMatthew G. Knepley Output Parameter: 124437b16e26SMatthew G. Knepley . X - vector 124537b16e26SMatthew G. Knepley 124637b16e26SMatthew G. Knepley Level: developer 124737b16e26SMatthew G. Knepley 124837b16e26SMatthew G. Knepley .seealso: DMDAComputeL2Diff() 124937b16e26SMatthew G. Knepley @*/ 125037b16e26SMatthew G. Knepley PetscErrorCode DMDAProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X) 125137b16e26SMatthew G. Knepley { 125237b16e26SMatthew G. Knepley Vec localX; 125337b16e26SMatthew G. Knepley PetscErrorCode ierr; 125437b16e26SMatthew G. Knepley 125537b16e26SMatthew G. Knepley PetscFunctionBegin; 125637b16e26SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 125737b16e26SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 125837b16e26SMatthew G. Knepley ierr = DMDAProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr); 125937b16e26SMatthew G. Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 126037b16e26SMatthew G. Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 126137b16e26SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 126237b16e26SMatthew G. Knepley PetscFunctionReturn(0); 126337b16e26SMatthew G. Knepley } 126437b16e26SMatthew G. Knepley 126537b16e26SMatthew G. Knepley #undef __FUNCT__ 126637b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2Diff" 126737b16e26SMatthew G. Knepley /*@C 126837b16e26SMatthew G. Knepley DMDAComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 126937b16e26SMatthew G. Knepley 127037b16e26SMatthew G. Knepley Input Parameters: 127137b16e26SMatthew G. Knepley + dm - The DM 127237b16e26SMatthew G. Knepley . fe - The PetscFE object for each field 127337b16e26SMatthew G. Knepley . funcs - The functions to evaluate for each field component 127437b16e26SMatthew G. Knepley - X - The coefficient vector u_h 127537b16e26SMatthew G. Knepley 127637b16e26SMatthew G. Knepley Output Parameter: 127737b16e26SMatthew G. Knepley . diff - The diff ||u - u_h||_2 127837b16e26SMatthew G. Knepley 127937b16e26SMatthew G. Knepley Level: developer 128037b16e26SMatthew G. Knepley 128137b16e26SMatthew G. Knepley .seealso: DMDAProjectFunction() 128237b16e26SMatthew G. Knepley @*/ 128337b16e26SMatthew G. Knepley PetscErrorCode DMDAComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff) 128437b16e26SMatthew G. Knepley { 128537b16e26SMatthew G. Knepley const PetscInt debug = 0; 128637b16e26SMatthew G. Knepley PetscSection section; 128737b16e26SMatthew G. Knepley PetscQuadrature quad; 128837b16e26SMatthew G. Knepley Vec localX; 128937b16e26SMatthew G. Knepley PetscScalar *funcVal; 129037b16e26SMatthew G. Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 129137b16e26SMatthew G. Knepley PetscReal localDiff = 0.0; 129237b16e26SMatthew G. Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 129337b16e26SMatthew G. Knepley PetscErrorCode ierr; 129437b16e26SMatthew G. Knepley 129537b16e26SMatthew G. Knepley PetscFunctionBegin; 129637b16e26SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 129737b16e26SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 129837b16e26SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 129937b16e26SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 130037b16e26SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 130137b16e26SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 130237b16e26SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 130337b16e26SMatthew G. Knepley PetscInt Nc; 130437b16e26SMatthew G. Knepley 130537b16e26SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 130637b16e26SMatthew G. Knepley numComponents += Nc; 130737b16e26SMatthew G. Knepley } 130837b16e26SMatthew G. Knepley /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 130937b16e26SMatthew G. Knepley ierr = PetscMalloc5(numComponents,PetscScalar,&funcVal,dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr); 131037b16e26SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 131137b16e26SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 131237b16e26SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1313df66330eSJed Brown const PetscScalar *x = NULL; 131437b16e26SMatthew G. Knepley PetscReal elemDiff = 0.0; 131537b16e26SMatthew G. Knepley 131637b16e26SMatthew G. Knepley ierr = DMDAComputeCellGeometry(dm, c, &quad, v0, J, invJ, &detJ);CHKERRQ(ierr); 131737b16e26SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 131837b16e26SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 131937b16e26SMatthew G. Knepley 132037b16e26SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 132137b16e26SMatthew G. Knepley const PetscInt numQuadPoints = quad.numPoints; 132237b16e26SMatthew G. Knepley const PetscReal *quadPoints = quad.points; 132337b16e26SMatthew G. Knepley const PetscReal *quadWeights = quad.weights; 132437b16e26SMatthew G. Knepley PetscReal *basis; 132537b16e26SMatthew G. Knepley PetscInt numBasisFuncs, numBasisComps, q, d, e, fc, f; 132637b16e26SMatthew G. Knepley 132737b16e26SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 132837b16e26SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 132937b16e26SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 133037b16e26SMatthew G. Knepley if (debug) { 133137b16e26SMatthew G. Knepley char title[1024]; 133237b16e26SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 133337b16e26SMatthew G. Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 133437b16e26SMatthew G. Knepley } 133537b16e26SMatthew G. Knepley for (q = 0; q < numQuadPoints; ++q) { 133637b16e26SMatthew G. Knepley for (d = 0; d < dim; d++) { 133737b16e26SMatthew G. Knepley coords[d] = v0[d]; 133837b16e26SMatthew G. Knepley for (e = 0; e < dim; e++) { 133937b16e26SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 134037b16e26SMatthew G. Knepley } 134137b16e26SMatthew G. Knepley } 134237b16e26SMatthew G. Knepley (*funcs[field])(coords, funcVal); 134337b16e26SMatthew G. Knepley for (fc = 0; fc < numBasisComps; ++fc) { 134437b16e26SMatthew G. Knepley PetscScalar interpolant = 0.0; 134537b16e26SMatthew G. Knepley 134637b16e26SMatthew G. Knepley for (f = 0; f < numBasisFuncs; ++f) { 134737b16e26SMatthew G. Knepley const PetscInt fidx = f*numBasisComps+fc; 134837b16e26SMatthew G. Knepley interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 134937b16e26SMatthew G. Knepley } 135037b16e26SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 135137b16e26SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 135237b16e26SMatthew G. Knepley } 135337b16e26SMatthew G. Knepley } 135437b16e26SMatthew G. Knepley comp += numBasisComps; 135537b16e26SMatthew G. Knepley fieldOffset += numBasisFuncs*numBasisComps; 135637b16e26SMatthew G. Knepley } 135737b16e26SMatthew G. Knepley ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 135837b16e26SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 135937b16e26SMatthew G. Knepley localDiff += elemDiff; 136037b16e26SMatthew G. Knepley } 136137b16e26SMatthew G. Knepley ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 136237b16e26SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 136337b16e26SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 136437b16e26SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 1365a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 1366a66d4d66SMatthew G Knepley } 1367a66d4d66SMatthew G Knepley 136847c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 136947c6ae99SBarry Smith 137047c6ae99SBarry Smith #undef __FUNCT__ 1371aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray" 137247c6ae99SBarry Smith /*@C 1373aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 137447c6ae99SBarry Smith 137547c6ae99SBarry Smith Input Parameter: 137647c6ae99SBarry Smith + da - information about my local patch 137747c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 137847c6ae99SBarry Smith 137947c6ae99SBarry Smith Output Parameters: 138047c6ae99SBarry Smith . vptr - array data structured 138147c6ae99SBarry Smith 138247c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 138347c6ae99SBarry Smith to zero them. 138447c6ae99SBarry Smith 138547c6ae99SBarry Smith Level: advanced 138647c6ae99SBarry Smith 1387bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 138847c6ae99SBarry Smith 138947c6ae99SBarry Smith @*/ 13907087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 139147c6ae99SBarry Smith { 139247c6ae99SBarry Smith PetscErrorCode ierr; 139347c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 139447c6ae99SBarry Smith char *iarray_start; 139547c6ae99SBarry Smith void **iptr = (void**)vptr; 139647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 139747c6ae99SBarry Smith 139847c6ae99SBarry Smith PetscFunctionBegin; 139947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 140047c6ae99SBarry Smith if (ghosted) { 1401aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 140247c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 140347c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 140447c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 14050298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 14060298fd71SBarry Smith dd->startghostedin[i] = NULL; 140747c6ae99SBarry Smith 140847c6ae99SBarry Smith goto done; 140947c6ae99SBarry Smith } 141047c6ae99SBarry Smith } 141147c6ae99SBarry Smith xs = dd->Xs; 141247c6ae99SBarry Smith ys = dd->Ys; 141347c6ae99SBarry Smith zs = dd->Zs; 141447c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 141547c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 141647c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 141747c6ae99SBarry Smith } else { 1418aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 141947c6ae99SBarry Smith if (dd->arrayin[i]) { 142047c6ae99SBarry Smith *iptr = dd->arrayin[i]; 142147c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 14220298fd71SBarry Smith dd->arrayin[i] = NULL; 14230298fd71SBarry Smith dd->startin[i] = NULL; 142447c6ae99SBarry Smith 142547c6ae99SBarry Smith goto done; 142647c6ae99SBarry Smith } 142747c6ae99SBarry Smith } 142847c6ae99SBarry Smith xs = dd->xs; 142947c6ae99SBarry Smith ys = dd->ys; 143047c6ae99SBarry Smith zs = dd->zs; 143147c6ae99SBarry Smith xm = dd->xe-dd->xs; 143247c6ae99SBarry Smith ym = dd->ye-dd->ys; 143347c6ae99SBarry Smith zm = dd->ze-dd->zs; 143447c6ae99SBarry Smith } 143547c6ae99SBarry Smith 143647c6ae99SBarry Smith switch (dd->dim) { 143747c6ae99SBarry Smith case 1: { 143847c6ae99SBarry Smith void *ptr; 143947c6ae99SBarry Smith 144047c6ae99SBarry Smith ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 144147c6ae99SBarry Smith 144247c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 144347c6ae99SBarry Smith *iptr = (void*)ptr; 14448865f1eaSKarl Rupp break; 14458865f1eaSKarl Rupp } 144647c6ae99SBarry Smith case 2: { 144747c6ae99SBarry Smith void **ptr; 144847c6ae99SBarry Smith 144947c6ae99SBarry Smith ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 145047c6ae99SBarry Smith 145147c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 14528865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 145347c6ae99SBarry Smith *iptr = (void*)ptr; 14548865f1eaSKarl Rupp break; 14558865f1eaSKarl Rupp } 145647c6ae99SBarry Smith case 3: { 145747c6ae99SBarry Smith void ***ptr,**bptr; 145847c6ae99SBarry Smith 145947c6ae99SBarry Smith ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 146047c6ae99SBarry Smith 146147c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 146247c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 14638865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 146447c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 146547c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 146647c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 146747c6ae99SBarry Smith } 146847c6ae99SBarry Smith } 146947c6ae99SBarry Smith *iptr = (void*)ptr; 14708865f1eaSKarl Rupp break; 14718865f1eaSKarl Rupp } 147247c6ae99SBarry Smith default: 1473ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 147447c6ae99SBarry Smith } 147547c6ae99SBarry Smith 147647c6ae99SBarry Smith done: 147747c6ae99SBarry Smith /* add arrays to the checked out list */ 147847c6ae99SBarry Smith if (ghosted) { 1479aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 148047c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 148147c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 148247c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 148347c6ae99SBarry Smith break; 148447c6ae99SBarry Smith } 148547c6ae99SBarry Smith } 148647c6ae99SBarry Smith } else { 1487aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 148847c6ae99SBarry Smith if (!dd->arrayout[i]) { 148947c6ae99SBarry Smith dd->arrayout[i] = *iptr; 149047c6ae99SBarry Smith dd->startout[i] = iarray_start; 149147c6ae99SBarry Smith break; 149247c6ae99SBarry Smith } 149347c6ae99SBarry Smith } 149447c6ae99SBarry Smith } 149547c6ae99SBarry Smith PetscFunctionReturn(0); 149647c6ae99SBarry Smith } 149747c6ae99SBarry Smith 149847c6ae99SBarry Smith #undef __FUNCT__ 1499aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray" 150047c6ae99SBarry Smith /*@C 1501aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 150247c6ae99SBarry Smith 150347c6ae99SBarry Smith Input Parameter: 150447c6ae99SBarry Smith + da - information about my local patch 150547c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 150647c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 150747c6ae99SBarry Smith 150847c6ae99SBarry Smith Level: advanced 150947c6ae99SBarry Smith 1510bcaeba4dSBarry Smith .seealso: DMDAGetArray() 151147c6ae99SBarry Smith 151247c6ae99SBarry Smith @*/ 15137087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 151447c6ae99SBarry Smith { 151547c6ae99SBarry Smith PetscInt i; 151647c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 151747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 151847c6ae99SBarry Smith 151947c6ae99SBarry Smith PetscFunctionBegin; 152047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 152147c6ae99SBarry Smith if (ghosted) { 1522aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 152347c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 152447c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 15250298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 15260298fd71SBarry Smith dd->startghostedout[i] = NULL; 152747c6ae99SBarry Smith break; 152847c6ae99SBarry Smith } 152947c6ae99SBarry Smith } 1530aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 153147c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 153247c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 153347c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 153447c6ae99SBarry Smith break; 153547c6ae99SBarry Smith } 153647c6ae99SBarry Smith } 153747c6ae99SBarry Smith } else { 1538aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 153947c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 154047c6ae99SBarry Smith iarray_start = dd->startout[i]; 15410298fd71SBarry Smith dd->arrayout[i] = NULL; 15420298fd71SBarry Smith dd->startout[i] = NULL; 154347c6ae99SBarry Smith break; 154447c6ae99SBarry Smith } 154547c6ae99SBarry Smith } 1546aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 154747c6ae99SBarry Smith if (!dd->arrayin[i]) { 154847c6ae99SBarry Smith dd->arrayin[i] = *iptr; 154947c6ae99SBarry Smith dd->startin[i] = iarray_start; 155047c6ae99SBarry Smith break; 155147c6ae99SBarry Smith } 155247c6ae99SBarry Smith } 155347c6ae99SBarry Smith } 155447c6ae99SBarry Smith PetscFunctionReturn(0); 155547c6ae99SBarry Smith } 155647c6ae99SBarry Smith 1557