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*/ 70c312b8eSJed Brown #include <petscsf.h> 847c6ae99SBarry Smith 947c6ae99SBarry Smith /* 10e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions 1147c6ae99SBarry Smith */ 1247c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 13c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */ 14c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */ 1547c6ae99SBarry Smith #undef __FUNCT__ 1647c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d" 171e6b0712SBarry Smith static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 1847c6ae99SBarry Smith { 1947c6ae99SBarry Smith PetscErrorCode ierr; 2047c6ae99SBarry Smith PetscInt n,m; 2147c6ae99SBarry Smith Vec vec = (Vec)obj; 2247c6ae99SBarry Smith PetscScalar *array; 2347c6ae99SBarry Smith mxArray *mat; 249a42bb27SBarry Smith DM da; 2547c6ae99SBarry Smith 2647c6ae99SBarry Smith PetscFunctionBegin; 27c688c046SMatthew G Knepley ierr = VecGetDM(vec, &da);CHKERRQ(ierr); 28ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 29aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 3047c6ae99SBarry Smith 3147c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 3247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3347c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3447c6ae99SBarry Smith #else 3547c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3647c6ae99SBarry Smith #endif 3747c6ae99SBarry Smith ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 3847c6ae99SBarry Smith ierr = PetscObjectName(obj);CHKERRQ(ierr); 3947c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 4047c6ae99SBarry Smith 4147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 4247c6ae99SBarry Smith PetscFunctionReturn(0); 4347c6ae99SBarry Smith } 4447c6ae99SBarry Smith #endif 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith 4747c6ae99SBarry Smith #undef __FUNCT__ 48564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA" 497087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 5047c6ae99SBarry Smith { 5147c6ae99SBarry Smith PetscErrorCode ierr; 5247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5347c6ae99SBarry Smith 5447c6ae99SBarry Smith PetscFunctionBegin; 5547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5647c6ae99SBarry Smith PetscValidPointer(g,2); 5711689aebSJed Brown if (da->defaultSection) { 5811689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 5911689aebSJed Brown } else { 6047c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6147c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6247c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 63401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 64c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 6647c6ae99SBarry Smith if (dd->w == 1 && dd->dim == 2) { 67bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith #endif 7011689aebSJed Brown } 7147c6ae99SBarry Smith PetscFunctionReturn(0); 7247c6ae99SBarry Smith } 7347c6ae99SBarry Smith 74a66d4d66SMatthew G Knepley #undef __FUNCT__ 7557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells" 76*3582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 7757459e9aSMatthew G Knepley { 7857459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 7957459e9aSMatthew G Knepley const PetscInt dim = da->dim; 8057459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 8157459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 8257459e9aSMatthew G Knepley 8357459e9aSMatthew G Knepley PetscFunctionBegin; 84*3582350cSMatthew G. Knepley if (numCellsX) { 85*3582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 86*3582350cSMatthew G. Knepley *numCellsX = mx; 87*3582350cSMatthew G. Knepley } 88*3582350cSMatthew G. Knepley if (numCellsY) { 89*3582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,3); 90*3582350cSMatthew G. Knepley *numCellsY = my; 91*3582350cSMatthew G. Knepley } 92*3582350cSMatthew G. Knepley if (numCellsZ) { 93*3582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,4); 94*3582350cSMatthew G. Knepley *numCellsZ = mz; 95*3582350cSMatthew G. Knepley } 9657459e9aSMatthew G Knepley if (numCells) { 97*3582350cSMatthew G. Knepley PetscValidIntPointer(numCells,5); 9857459e9aSMatthew G Knepley *numCells = nC; 9957459e9aSMatthew G Knepley } 10057459e9aSMatthew G Knepley PetscFunctionReturn(0); 10157459e9aSMatthew G Knepley } 10257459e9aSMatthew G Knepley 10357459e9aSMatthew G Knepley #undef __FUNCT__ 10457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices" 10557459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 10657459e9aSMatthew G Knepley { 10757459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 10857459e9aSMatthew G Knepley const PetscInt dim = da->dim; 10957459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 11057459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 11157459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 11257459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 11357459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 11457459e9aSMatthew G Knepley 11557459e9aSMatthew G Knepley PetscFunctionBegin; 11657459e9aSMatthew G Knepley if (numVerticesX) { 11757459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 11857459e9aSMatthew G Knepley *numVerticesX = nVx; 11957459e9aSMatthew G Knepley } 12057459e9aSMatthew G Knepley if (numVerticesY) { 12157459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 12257459e9aSMatthew G Knepley *numVerticesY = nVy; 12357459e9aSMatthew G Knepley } 12457459e9aSMatthew G Knepley if (numVerticesZ) { 12557459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 12657459e9aSMatthew G Knepley *numVerticesZ = nVz; 12757459e9aSMatthew G Knepley } 12857459e9aSMatthew G Knepley if (numVertices) { 12957459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 13057459e9aSMatthew G Knepley *numVertices = nV; 13157459e9aSMatthew G Knepley } 13257459e9aSMatthew G Knepley PetscFunctionReturn(0); 13357459e9aSMatthew G Knepley } 13457459e9aSMatthew G Knepley 13557459e9aSMatthew G Knepley #undef __FUNCT__ 13657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces" 13757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 13857459e9aSMatthew G Knepley { 13957459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 14057459e9aSMatthew G Knepley const PetscInt dim = da->dim; 14157459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 14257459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 14357459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 14457459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 14557459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 14657459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 14757459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 14857459e9aSMatthew G Knepley 14957459e9aSMatthew G Knepley PetscFunctionBegin; 15057459e9aSMatthew G Knepley if (numXFacesX) { 15157459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 15257459e9aSMatthew G Knepley *numXFacesX = nxF; 15357459e9aSMatthew G Knepley } 15457459e9aSMatthew G Knepley if (numXFaces) { 15557459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 15657459e9aSMatthew G Knepley *numXFaces = nXF; 15757459e9aSMatthew G Knepley } 15857459e9aSMatthew G Knepley if (numYFacesY) { 15957459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 16057459e9aSMatthew G Knepley *numYFacesY = nyF; 16157459e9aSMatthew G Knepley } 16257459e9aSMatthew G Knepley if (numYFaces) { 16357459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 16457459e9aSMatthew G Knepley *numYFaces = nYF; 16557459e9aSMatthew G Knepley } 16657459e9aSMatthew G Knepley if (numZFacesZ) { 16757459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 16857459e9aSMatthew G Knepley *numZFacesZ = nzF; 16957459e9aSMatthew G Knepley } 17057459e9aSMatthew G Knepley if (numZFaces) { 17157459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 17257459e9aSMatthew G Knepley *numZFaces = nZF; 17357459e9aSMatthew G Knepley } 17457459e9aSMatthew G Knepley PetscFunctionReturn(0); 17557459e9aSMatthew G Knepley } 17657459e9aSMatthew G Knepley 17757459e9aSMatthew G Knepley #undef __FUNCT__ 17857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum" 17957459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 18057459e9aSMatthew G Knepley { 18157459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 18257459e9aSMatthew G Knepley const PetscInt dim = da->dim; 18357459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 18457459e9aSMatthew G Knepley PetscErrorCode ierr; 18557459e9aSMatthew G Knepley 18657459e9aSMatthew G Knepley PetscFunctionBegin; 1878865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 1888865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 189*3582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 1900298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 1910298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 19257459e9aSMatthew G Knepley if (height == 0) { 19357459e9aSMatthew G Knepley /* Cells */ 1948865f1eaSKarl Rupp if (pStart) *pStart = 0; 1958865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 19657459e9aSMatthew G Knepley } else if (height == 1) { 19757459e9aSMatthew G Knepley /* Faces */ 1988865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 1998865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 20057459e9aSMatthew G Knepley } else if (height == dim) { 20157459e9aSMatthew G Knepley /* Vertices */ 2028865f1eaSKarl Rupp if (pStart) *pStart = nC; 2038865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 20457459e9aSMatthew G Knepley } else if (height < 0) { 20557459e9aSMatthew G Knepley /* All points */ 2068865f1eaSKarl Rupp if (pStart) *pStart = 0; 2078865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 20882f516ccSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 20957459e9aSMatthew G Knepley PetscFunctionReturn(0); 21057459e9aSMatthew G Knepley } 21157459e9aSMatthew G Knepley 21257459e9aSMatthew G Knepley #undef __FUNCT__ 213a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection" 214a66d4d66SMatthew G Knepley /*@C 215a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 216a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 217a66d4d66SMatthew G Knepley 218a66d4d66SMatthew G Knepley Input Parameters: 219a66d4d66SMatthew G Knepley + dm- The DMDA 220a66d4d66SMatthew G Knepley . numFields - The number of fields 2210298fd71SBarry Smith . numComp - The number of components in each field, or NULL for 1 2220298fd71SBarry Smith . numVertexDof - The number of dofs per vertex for each field, or NULL 2230298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL 2240298fd71SBarry Smith - numCellDof - The number of dofs per cell for each field, or NULL 225a66d4d66SMatthew G Knepley 226a66d4d66SMatthew G Knepley Level: developer 227a66d4d66SMatthew G Knepley 228a66d4d66SMatthew G Knepley Note: 229a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 230a66d4d66SMatthew G Knepley 231a66d4d66SMatthew G Knepley - Cells: [0, nC) 232a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 23388ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 23488ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 23588ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 236a66d4d66SMatthew G Knepley 237a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 238a66d4d66SMatthew G Knepley @*/ 23980800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[]) 240a66d4d66SMatthew G Knepley { 241a66d4d66SMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 242a4b60ecfSMatthew G. Knepley PetscSection section; 24388ed4aceSMatthew G Knepley const PetscInt dim = da->dim; 24480800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 24588ed4aceSMatthew G Knepley PetscSF sf; 246feafbc5cSMatthew G Knepley PetscMPIInt rank; 24788ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 24888ed4aceSMatthew G Knepley PetscInt *localPoints; 24988ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 250f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 25157459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 25257459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 25388ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 254a66d4d66SMatthew G Knepley PetscErrorCode ierr; 255a66d4d66SMatthew G Knepley 256a66d4d66SMatthew G Knepley PetscFunctionBegin; 257a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 258*3582350cSMatthew G. Knepley PetscValidPointer(s, 4); 25982f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 260*3582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 26157459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 26257459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 26357459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 26457459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 26557459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 26657459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 26757459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 26857459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 26957459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 27082f516ccSBarry Smith if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 27188ed4aceSMatthew G Knepley /* Create local section */ 27280800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 273a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 2748865f1eaSKarl Rupp if (numVertexDof) numVertexTotDof += numVertexDof[f]; 2758865f1eaSKarl Rupp if (numCellDof) numCellTotDof += numCellDof[f]; 2768865f1eaSKarl Rupp if (numFaceDof) { 2778865f1eaSKarl Rupp numFaceTotDof[0] += numFaceDof[f*dim+0]; 27888ed4aceSMatthew G Knepley numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0; 2790ad7597dSKarl Rupp numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0; 2800ad7597dSKarl Rupp } 281a66d4d66SMatthew G Knepley } 282a4b60ecfSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);CHKERRQ(ierr); 283a66d4d66SMatthew G Knepley if (numFields > 1) { 284a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr); 285a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 28680800b1aSMatthew G Knepley const char *name; 28780800b1aSMatthew G Knepley 28880800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 289a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr); 29080800b1aSMatthew G Knepley if (numComp) { 291a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr); 292a66d4d66SMatthew G Knepley } 293a66d4d66SMatthew G Knepley } 294a66d4d66SMatthew G Knepley } else { 295a66d4d66SMatthew G Knepley numFields = 0; 296a66d4d66SMatthew G Knepley } 297a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 298a66d4d66SMatthew G Knepley if (numVertexDof) { 299a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 300a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 301a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, v, f, numVertexDof[f]);CHKERRQ(ierr); 302a66d4d66SMatthew G Knepley } 303a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr); 304a66d4d66SMatthew G Knepley } 305a66d4d66SMatthew G Knepley } 306a66d4d66SMatthew G Knepley if (numFaceDof) { 307a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 308a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 309a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr); 310a66d4d66SMatthew G Knepley } 311a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr); 312a66d4d66SMatthew G Knepley } 313a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 314a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 315a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr); 316a66d4d66SMatthew G Knepley } 317a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr); 318a66d4d66SMatthew G Knepley } 319a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 320a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 321a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr); 322a66d4d66SMatthew G Knepley } 323a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr); 324a66d4d66SMatthew G Knepley } 325a66d4d66SMatthew G Knepley } 326a66d4d66SMatthew G Knepley if (numCellDof) { 327a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 328a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 329a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, c, f, numCellDof[f]);CHKERRQ(ierr); 330a66d4d66SMatthew G Knepley } 331a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr); 332a66d4d66SMatthew G Knepley } 333a66d4d66SMatthew G Knepley } 334a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 33588ed4aceSMatthew G Knepley /* Create mesh point SF */ 33688ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 33788ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 33888ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 33988ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3407128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 34188ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 34288ed4aceSMatthew G Knepley 3433814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 344feafbc5cSMatthew G Knepley nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */ 34588ed4aceSMatthew G Knepley if (xp && !yp && !zp) { 34688ed4aceSMatthew G Knepley nleaves += nxF; /* x faces */ 34788ed4aceSMatthew G Knepley } else if (yp && !zp && !xp) { 34888ed4aceSMatthew G Knepley nleaves += nyF; /* y faces */ 34988ed4aceSMatthew G Knepley } else if (zp && !xp && !yp) { 35088ed4aceSMatthew G Knepley nleaves += nzF; /* z faces */ 35188ed4aceSMatthew G Knepley } 35288ed4aceSMatthew G Knepley } 35388ed4aceSMatthew G Knepley } 35488ed4aceSMatthew G Knepley } 35588ed4aceSMatthew G Knepley } 356dcca6d9dSJed Brown ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr); 35788ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 35888ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 35988ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3607128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 36188ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 362f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 36388ed4aceSMatthew G Knepley 3643814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 36588ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 36688ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 36788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 368f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; 369f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 370628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom back vertex */ 371f5eeb527SMatthew G Knepley 372f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 373f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 374f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 375f5eeb527SMatthew G Knepley ++nL; 37688ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 377f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; 378f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 379628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom front vertex */ 380f5eeb527SMatthew G Knepley 381f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 382f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 383f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 384f5eeb527SMatthew G Knepley ++nL; 38588ed4aceSMatthew G Knepley } else { 38688ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left bottom vertices */ 387f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 388f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; 389f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 390f5eeb527SMatthew G Knepley 391f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 392f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 393f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 394f5eeb527SMatthew G Knepley } 39588ed4aceSMatthew G Knepley } 39688ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 39788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 398f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; 399f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 400628e3b0dSSatish Balay nleavesCheck += 1; /* left top back vertex */ 401f5eeb527SMatthew G Knepley 402f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 403f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 404f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 405f5eeb527SMatthew G Knepley ++nL; 40688ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 407f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; 408f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 409628e3b0dSSatish Balay nleavesCheck += 1; /* left top front vertex */ 410f5eeb527SMatthew G Knepley 411f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 412f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 413f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 414f5eeb527SMatthew G Knepley ++nL; 41588ed4aceSMatthew G Knepley } else { 41688ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left top vertices */ 417f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 418f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; 419f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 420f5eeb527SMatthew G Knepley 421f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 422f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 423f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 424f5eeb527SMatthew G Knepley } 42588ed4aceSMatthew G Knepley } 42688ed4aceSMatthew G Knepley } else { 42788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 42888ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left back vertices */ 429f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 430f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; 431f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 432f5eeb527SMatthew G Knepley 433f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 434f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 435f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 436f5eeb527SMatthew G Knepley } 43788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 43888ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left front vertices */ 439f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 440f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; 441f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 442f5eeb527SMatthew G Knepley 443f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 444f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 445f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 446f5eeb527SMatthew G Knepley } 44788ed4aceSMatthew G Knepley } else { 44888ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* left vertices */ 449f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 450f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 451f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; 452f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 453f5eeb527SMatthew G Knepley 454f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 455f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 456f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 457f5eeb527SMatthew G Knepley } 458f5eeb527SMatthew G Knepley } 45988ed4aceSMatthew G Knepley nleavesCheck += nxF; /* left faces */ 460f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 461f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 462f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 463f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 464f5eeb527SMatthew G Knepley 465f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 466f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 467f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 468f5eeb527SMatthew G Knepley } 46988ed4aceSMatthew G Knepley } 47088ed4aceSMatthew G Knepley } 47188ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 47288ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 47388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 474f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; 475f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 476628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom back vertex */ 477f5eeb527SMatthew G Knepley 478f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 479f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 480f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 481f5eeb527SMatthew G Knepley ++nL; 48288ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 483f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; 484f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 485628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom front vertex */ 486f5eeb527SMatthew G Knepley 487f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 488f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 489f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 490f5eeb527SMatthew G Knepley ++nL; 49188ed4aceSMatthew G Knepley } else { 49288ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right bottom vertices */ 493f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 494f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; 495f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 496f5eeb527SMatthew G Knepley 497f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 498f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 499f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 500f5eeb527SMatthew G Knepley } 50188ed4aceSMatthew G Knepley } 50288ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 50388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 504f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; 505f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 506628e3b0dSSatish Balay nleavesCheck += 1; /* right top back vertex */ 507f5eeb527SMatthew G Knepley 508f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 509f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 510f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 511f5eeb527SMatthew G Knepley ++nL; 51288ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 513f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; 514f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 515628e3b0dSSatish Balay nleavesCheck += 1; /* right top front vertex */ 516f5eeb527SMatthew G Knepley 517f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 518f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 519f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 520f5eeb527SMatthew G Knepley ++nL; 52188ed4aceSMatthew G Knepley } else { 52288ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right top vertices */ 523f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 524f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; 525f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 526f5eeb527SMatthew G Knepley 527f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 528f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 529f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 530f5eeb527SMatthew G Knepley } 53188ed4aceSMatthew G Knepley } 53288ed4aceSMatthew G Knepley } else { 53388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 53488ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right back vertices */ 535f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 536f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; 537f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 538f5eeb527SMatthew G Knepley 539f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 540f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 541f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 542f5eeb527SMatthew G Knepley } 54388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 54488ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right front vertices */ 545f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 546f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; 547f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 548f5eeb527SMatthew G Knepley 549f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 550f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 551f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 552f5eeb527SMatthew G Knepley } 55388ed4aceSMatthew G Knepley } else { 55488ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* right vertices */ 555f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 556f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 557f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; 558f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 559f5eeb527SMatthew G Knepley 560f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 561f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 562f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 563f5eeb527SMatthew G Knepley } 564f5eeb527SMatthew G Knepley } 56588ed4aceSMatthew G Knepley nleavesCheck += nxF; /* right faces */ 566f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 567f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 568f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 569f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 570f5eeb527SMatthew G Knepley 571f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 572f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 573f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 574f5eeb527SMatthew G Knepley } 57588ed4aceSMatthew G Knepley } 57688ed4aceSMatthew G Knepley } 57788ed4aceSMatthew G Knepley } else { 57888ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 57988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 58088ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom back vertices */ 581f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 582f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; 583f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 584f5eeb527SMatthew G Knepley 585f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 586f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 587f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 588f5eeb527SMatthew G Knepley } 58988ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 59088ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom front vertices */ 591f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 592f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; 593f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 594f5eeb527SMatthew G Knepley 595f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 596f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 597f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 598f5eeb527SMatthew G Knepley } 59988ed4aceSMatthew G Knepley } else { 60088ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* bottom vertices */ 601f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 602f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 603f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; 604f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 605f5eeb527SMatthew G Knepley 606f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 607f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 608f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 609f5eeb527SMatthew G Knepley } 610f5eeb527SMatthew G Knepley } 61188ed4aceSMatthew G Knepley nleavesCheck += nyF; /* bottom faces */ 612f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 613f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 614f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 615f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 616f5eeb527SMatthew G Knepley 617f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 618f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 619f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 620f5eeb527SMatthew G Knepley } 62188ed4aceSMatthew G Knepley } 62288ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 62388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 62488ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top back vertices */ 625f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 626f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; 627f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 628f5eeb527SMatthew G Knepley 629f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 630f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 631f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 632f5eeb527SMatthew G Knepley } 63388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 63488ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top front vertices */ 635f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 636f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; 637f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 638f5eeb527SMatthew G Knepley 639f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 640f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 641f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 642f5eeb527SMatthew G Knepley } 64388ed4aceSMatthew G Knepley } else { 64488ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* top vertices */ 645f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 646f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 647f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; 648f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 649f5eeb527SMatthew G Knepley 650f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 651f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 652f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 653f5eeb527SMatthew G Knepley } 654f5eeb527SMatthew G Knepley } 65588ed4aceSMatthew G Knepley nleavesCheck += nyF; /* top faces */ 656f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 657f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 658f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 659f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 660f5eeb527SMatthew G Knepley 661f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 662f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 663f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 664f5eeb527SMatthew G Knepley } 66588ed4aceSMatthew G Knepley } 66688ed4aceSMatthew G Knepley } else { 66788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 66888ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* back vertices */ 669f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 670f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 671f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; 672f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 673f5eeb527SMatthew G Knepley 674f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 675f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 676f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 677f5eeb527SMatthew G Knepley } 678f5eeb527SMatthew G Knepley } 67988ed4aceSMatthew G Knepley nleavesCheck += nzF; /* back faces */ 680f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 681f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 682f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 683f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 684f5eeb527SMatthew G Knepley 685f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 686f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 687f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 688f5eeb527SMatthew G Knepley } 68988ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 69088ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* front vertices */ 691f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 692f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 693f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; 694f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 695f5eeb527SMatthew G Knepley 696f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 697f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 698f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 699f5eeb527SMatthew G Knepley } 700f5eeb527SMatthew G Knepley } 70188ed4aceSMatthew G Knepley nleavesCheck += nzF; /* front faces */ 702f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 703f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 704f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 705f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 706f5eeb527SMatthew G Knepley 707f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 708f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 709f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 710f5eeb527SMatthew G Knepley } 71188ed4aceSMatthew G Knepley } else { 71288ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 71388ed4aceSMatthew G Knepley } 71488ed4aceSMatthew G Knepley } 71588ed4aceSMatthew G Knepley } 71688ed4aceSMatthew G Knepley } 71788ed4aceSMatthew G Knepley } 71888ed4aceSMatthew G Knepley } 71988ed4aceSMatthew G Knepley } 7203814d604SMatthew G Knepley /* TODO: Remove duplication in leaf determination */ 72182f516ccSBarry Smith if (nleaves != nleavesCheck) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck); 72282f516ccSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr); 72388ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 724a4b60ecfSMatthew G. Knepley ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 72588ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 726a4b60ecfSMatthew G. Knepley ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr); 727a4b60ecfSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 728a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 729a66d4d66SMatthew G Knepley } 730a66d4d66SMatthew G Knepley 73147c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith #undef __FUNCT__ 734aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray" 73547c6ae99SBarry Smith /*@C 736aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 73747c6ae99SBarry Smith 73847c6ae99SBarry Smith Input Parameter: 73947c6ae99SBarry Smith + da - information about my local patch 74047c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 74147c6ae99SBarry Smith 74247c6ae99SBarry Smith Output Parameters: 74347c6ae99SBarry Smith . vptr - array data structured 74447c6ae99SBarry Smith 74547c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 74647c6ae99SBarry Smith to zero them. 74747c6ae99SBarry Smith 74847c6ae99SBarry Smith Level: advanced 74947c6ae99SBarry Smith 750bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 75147c6ae99SBarry Smith 75247c6ae99SBarry Smith @*/ 7537087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 75447c6ae99SBarry Smith { 75547c6ae99SBarry Smith PetscErrorCode ierr; 75647c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 75747c6ae99SBarry Smith char *iarray_start; 75847c6ae99SBarry Smith void **iptr = (void**)vptr; 75947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith PetscFunctionBegin; 76247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 76347c6ae99SBarry Smith if (ghosted) { 764aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 76547c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 76647c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 76747c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 7680298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 7690298fd71SBarry Smith dd->startghostedin[i] = NULL; 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith goto done; 77247c6ae99SBarry Smith } 77347c6ae99SBarry Smith } 77447c6ae99SBarry Smith xs = dd->Xs; 77547c6ae99SBarry Smith ys = dd->Ys; 77647c6ae99SBarry Smith zs = dd->Zs; 77747c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 77847c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 77947c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 78047c6ae99SBarry Smith } else { 781aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 78247c6ae99SBarry Smith if (dd->arrayin[i]) { 78347c6ae99SBarry Smith *iptr = dd->arrayin[i]; 78447c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 7850298fd71SBarry Smith dd->arrayin[i] = NULL; 7860298fd71SBarry Smith dd->startin[i] = NULL; 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith goto done; 78947c6ae99SBarry Smith } 79047c6ae99SBarry Smith } 79147c6ae99SBarry Smith xs = dd->xs; 79247c6ae99SBarry Smith ys = dd->ys; 79347c6ae99SBarry Smith zs = dd->zs; 79447c6ae99SBarry Smith xm = dd->xe-dd->xs; 79547c6ae99SBarry Smith ym = dd->ye-dd->ys; 79647c6ae99SBarry Smith zm = dd->ze-dd->zs; 79747c6ae99SBarry Smith } 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith switch (dd->dim) { 80047c6ae99SBarry Smith case 1: { 80147c6ae99SBarry Smith void *ptr; 80247c6ae99SBarry Smith 803901f1932SJed Brown ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 80447c6ae99SBarry Smith 80547c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 80647c6ae99SBarry Smith *iptr = (void*)ptr; 8078865f1eaSKarl Rupp break; 8088865f1eaSKarl Rupp } 80947c6ae99SBarry Smith case 2: { 81047c6ae99SBarry Smith void **ptr; 81147c6ae99SBarry Smith 812901f1932SJed Brown ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 81347c6ae99SBarry Smith 81447c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 8158865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 81647c6ae99SBarry Smith *iptr = (void*)ptr; 8178865f1eaSKarl Rupp break; 8188865f1eaSKarl Rupp } 81947c6ae99SBarry Smith case 3: { 82047c6ae99SBarry Smith void ***ptr,**bptr; 82147c6ae99SBarry Smith 822901f1932SJed Brown ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 82547c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 8268865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 82747c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 82847c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 82947c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith } 83247c6ae99SBarry Smith *iptr = (void*)ptr; 8338865f1eaSKarl Rupp break; 8348865f1eaSKarl Rupp } 83547c6ae99SBarry Smith default: 836ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith 83947c6ae99SBarry Smith done: 84047c6ae99SBarry Smith /* add arrays to the checked out list */ 84147c6ae99SBarry Smith if (ghosted) { 842aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 84347c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 84447c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 84547c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 84647c6ae99SBarry Smith break; 84747c6ae99SBarry Smith } 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith } else { 850aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 85147c6ae99SBarry Smith if (!dd->arrayout[i]) { 85247c6ae99SBarry Smith dd->arrayout[i] = *iptr; 85347c6ae99SBarry Smith dd->startout[i] = iarray_start; 85447c6ae99SBarry Smith break; 85547c6ae99SBarry Smith } 85647c6ae99SBarry Smith } 85747c6ae99SBarry Smith } 85847c6ae99SBarry Smith PetscFunctionReturn(0); 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith #undef __FUNCT__ 862aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray" 86347c6ae99SBarry Smith /*@C 864aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 86547c6ae99SBarry Smith 86647c6ae99SBarry Smith Input Parameter: 86747c6ae99SBarry Smith + da - information about my local patch 86847c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 86947c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 87047c6ae99SBarry Smith 87147c6ae99SBarry Smith Level: advanced 87247c6ae99SBarry Smith 873bcaeba4dSBarry Smith .seealso: DMDAGetArray() 87447c6ae99SBarry Smith 87547c6ae99SBarry Smith @*/ 8767087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 87747c6ae99SBarry Smith { 87847c6ae99SBarry Smith PetscInt i; 87947c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 88047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 88147c6ae99SBarry Smith 88247c6ae99SBarry Smith PetscFunctionBegin; 88347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 88447c6ae99SBarry Smith if (ghosted) { 885aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 88647c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 88747c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 8880298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 8890298fd71SBarry Smith dd->startghostedout[i] = NULL; 89047c6ae99SBarry Smith break; 89147c6ae99SBarry Smith } 89247c6ae99SBarry Smith } 893aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 89447c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 89547c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 89647c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 89747c6ae99SBarry Smith break; 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith } 90047c6ae99SBarry Smith } else { 901aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 90247c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 90347c6ae99SBarry Smith iarray_start = dd->startout[i]; 9040298fd71SBarry Smith dd->arrayout[i] = NULL; 9050298fd71SBarry Smith dd->startout[i] = NULL; 90647c6ae99SBarry Smith break; 90747c6ae99SBarry Smith } 90847c6ae99SBarry Smith } 909aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 91047c6ae99SBarry Smith if (!dd->arrayin[i]) { 91147c6ae99SBarry Smith dd->arrayin[i] = *iptr; 91247c6ae99SBarry Smith dd->startin[i] = iarray_start; 91347c6ae99SBarry Smith break; 91447c6ae99SBarry Smith } 91547c6ae99SBarry Smith } 91647c6ae99SBarry Smith } 91747c6ae99SBarry Smith PetscFunctionReturn(0); 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith 920