147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 847c6ae99SBarry Smith /* 9e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions 1047c6ae99SBarry Smith */ 1147c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 12c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */ 13c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */ 1447c6ae99SBarry Smith EXTERN_C_BEGIN 1547c6ae99SBarry Smith #undef __FUNCT__ 1647c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d" 177087cfbeSBarry Smith 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); 28aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)vec)->comm,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 EXTERN_C_END 4547c6ae99SBarry Smith #endif 4647c6ae99SBarry Smith 4747c6ae99SBarry Smith 4847c6ae99SBarry Smith #undef __FUNCT__ 49564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA" 507087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec* g) 5147c6ae99SBarry Smith { 5247c6ae99SBarry Smith PetscErrorCode ierr; 5347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5447c6ae99SBarry Smith 5547c6ae99SBarry Smith PetscFunctionBegin; 5647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5747c6ae99SBarry Smith PetscValidPointer(g,2); 5811689aebSJed Brown if (da->defaultSection) { 5911689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 6011689aebSJed Brown } else { 6147c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6247c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6347c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 64401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 65c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6647c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 6747c6ae99SBarry Smith if (dd->w == 1 && dd->dim == 2) { 6847c6ae99SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 6947c6ae99SBarry Smith } 7047c6ae99SBarry Smith #endif 7111689aebSJed Brown } 7247c6ae99SBarry Smith PetscFunctionReturn(0); 7347c6ae99SBarry Smith } 7447c6ae99SBarry Smith 75a66d4d66SMatthew G Knepley #undef __FUNCT__ 7657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells" 7757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCells) 7857459e9aSMatthew G Knepley { 7957459e9aSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 8057459e9aSMatthew G Knepley const PetscInt dim = da->dim; 8157459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 8257459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 8357459e9aSMatthew G Knepley 8457459e9aSMatthew G Knepley PetscFunctionBegin; 8557459e9aSMatthew G Knepley if (numCells) { 8657459e9aSMatthew G Knepley PetscValidIntPointer(numCells,2); 8757459e9aSMatthew G Knepley *numCells = nC; 8857459e9aSMatthew G Knepley } 8957459e9aSMatthew G Knepley PetscFunctionReturn(0); 9057459e9aSMatthew G Knepley } 9157459e9aSMatthew G Knepley 9257459e9aSMatthew G Knepley #undef __FUNCT__ 9357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices" 9457459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 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 nVx = mx+1; 10057459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 10157459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 10257459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 10357459e9aSMatthew G Knepley 10457459e9aSMatthew G Knepley PetscFunctionBegin; 10557459e9aSMatthew G Knepley if (numVerticesX) { 10657459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 10757459e9aSMatthew G Knepley *numVerticesX = nVx; 10857459e9aSMatthew G Knepley } 10957459e9aSMatthew G Knepley if (numVerticesY) { 11057459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 11157459e9aSMatthew G Knepley *numVerticesY = nVy; 11257459e9aSMatthew G Knepley } 11357459e9aSMatthew G Knepley if (numVerticesZ) { 11457459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 11557459e9aSMatthew G Knepley *numVerticesZ = nVz; 11657459e9aSMatthew G Knepley } 11757459e9aSMatthew G Knepley if (numVertices) { 11857459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 11957459e9aSMatthew G Knepley *numVertices = nV; 12057459e9aSMatthew G Knepley } 12157459e9aSMatthew G Knepley PetscFunctionReturn(0); 12257459e9aSMatthew G Knepley } 12357459e9aSMatthew G Knepley 12457459e9aSMatthew G Knepley #undef __FUNCT__ 12557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces" 12657459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 12757459e9aSMatthew G Knepley { 12857459e9aSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 12957459e9aSMatthew G Knepley const PetscInt dim = da->dim; 13057459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 13157459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 13257459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 13357459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 13457459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 13557459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 13657459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 13757459e9aSMatthew G Knepley 13857459e9aSMatthew G Knepley PetscFunctionBegin; 13957459e9aSMatthew G Knepley if (numXFacesX) { 14057459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 14157459e9aSMatthew G Knepley *numXFacesX = nxF; 14257459e9aSMatthew G Knepley } 14357459e9aSMatthew G Knepley if (numXFaces) { 14457459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 14557459e9aSMatthew G Knepley *numXFaces = nXF; 14657459e9aSMatthew G Knepley } 14757459e9aSMatthew G Knepley if (numYFacesY) { 14857459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 14957459e9aSMatthew G Knepley *numYFacesY = nyF; 15057459e9aSMatthew G Knepley } 15157459e9aSMatthew G Knepley if (numYFaces) { 15257459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 15357459e9aSMatthew G Knepley *numYFaces = nYF; 15457459e9aSMatthew G Knepley } 15557459e9aSMatthew G Knepley if (numZFacesZ) { 15657459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 15757459e9aSMatthew G Knepley *numZFacesZ = nzF; 15857459e9aSMatthew G Knepley } 15957459e9aSMatthew G Knepley if (numZFaces) { 16057459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 16157459e9aSMatthew G Knepley *numZFaces = nZF; 16257459e9aSMatthew G Knepley } 16357459e9aSMatthew G Knepley PetscFunctionReturn(0); 16457459e9aSMatthew G Knepley } 16557459e9aSMatthew G Knepley 16657459e9aSMatthew G Knepley #undef __FUNCT__ 16757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum" 16857459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 16957459e9aSMatthew G Knepley { 17057459e9aSMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 17157459e9aSMatthew G Knepley const PetscInt dim = da->dim; 17257459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 17357459e9aSMatthew G Knepley PetscErrorCode ierr; 17457459e9aSMatthew G Knepley 17557459e9aSMatthew G Knepley PetscFunctionBegin; 17657459e9aSMatthew G Knepley if (pStart) {PetscValidIntPointer(pStart,3);} 17757459e9aSMatthew G Knepley if (pEnd) {PetscValidIntPointer(pEnd,4);} 17857459e9aSMatthew G Knepley ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr); 17957459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, PETSC_NULL, PETSC_NULL, PETSC_NULL, &nV);CHKERRQ(ierr); 18057459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, PETSC_NULL, &nXF, PETSC_NULL, &nYF, PETSC_NULL, &nZF);CHKERRQ(ierr); 18157459e9aSMatthew G Knepley if (height == 0) { 18257459e9aSMatthew G Knepley /* Cells */ 18357459e9aSMatthew G Knepley if (pStart) {*pStart = 0;} 18457459e9aSMatthew G Knepley if (pEnd) {*pEnd = nC;} 18557459e9aSMatthew G Knepley } else if (height == 1) { 18657459e9aSMatthew G Knepley /* Faces */ 18757459e9aSMatthew G Knepley if (pStart) {*pStart = nC+nV;} 18857459e9aSMatthew G Knepley if (pEnd) {*pEnd = nC+nV+nXF+nYF+nZF;} 18957459e9aSMatthew G Knepley } else if (height == dim) { 19057459e9aSMatthew G Knepley /* Vertices */ 19157459e9aSMatthew G Knepley if (pStart) {*pStart = nC;} 19257459e9aSMatthew G Knepley if (pEnd) {*pEnd = nC+nV;} 19357459e9aSMatthew G Knepley } else if (height < 0) { 19457459e9aSMatthew G Knepley /* All points */ 19557459e9aSMatthew G Knepley if (pStart) {*pStart = 0;} 19657459e9aSMatthew G Knepley if (pEnd) {*pEnd = nC+nV+nXF+nYF+nZF;} 197*f23aa3ddSBarry Smith } else SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 19857459e9aSMatthew G Knepley PetscFunctionReturn(0); 19957459e9aSMatthew G Knepley } 20057459e9aSMatthew G Knepley 20157459e9aSMatthew G Knepley #undef __FUNCT__ 202a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection" 203a66d4d66SMatthew G Knepley /*@C 204a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 205a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 206a66d4d66SMatthew G Knepley 207a66d4d66SMatthew G Knepley Input Parameters: 208a66d4d66SMatthew G Knepley + dm- The DMDA 209a66d4d66SMatthew G Knepley . numFields - The number of fields 210a66d4d66SMatthew G Knepley . numComp - The number of components in each field, or PETSC_NULL for 1 211a66d4d66SMatthew G Knepley . numVertexDof - The number of dofs per vertex for each field, or PETSC_NULL 212a66d4d66SMatthew G Knepley . numFaceDof - The number of dofs per face for each field and direction, or PETSC_NULL 213a66d4d66SMatthew G Knepley - numCellDof - The number of dofs per cell for each field, or PETSC_NULL 214a66d4d66SMatthew G Knepley 215a66d4d66SMatthew G Knepley Level: developer 216a66d4d66SMatthew G Knepley 217a66d4d66SMatthew G Knepley Note: 218a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 219a66d4d66SMatthew G Knepley 220a66d4d66SMatthew G Knepley - Cells: [0, nC) 221a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 22288ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 22388ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 22488ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 225a66d4d66SMatthew G Knepley 226a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 227a66d4d66SMatthew G Knepley @*/ 22880800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[]) 229a66d4d66SMatthew G Knepley { 230a66d4d66SMatthew G Knepley DM_DA *da = (DM_DA *) dm->data; 23188ed4aceSMatthew G Knepley const PetscInt dim = da->dim; 23280800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 23388ed4aceSMatthew G Knepley PetscSF sf; 234feafbc5cSMatthew G Knepley PetscMPIInt rank; 23588ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 23688ed4aceSMatthew G Knepley PetscInt *localPoints; 23788ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 238f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 23957459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 24057459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 24188ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 242a66d4d66SMatthew G Knepley PetscErrorCode ierr; 243a66d4d66SMatthew G Knepley 244a66d4d66SMatthew G Knepley PetscFunctionBegin; 245a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 246feafbc5cSMatthew G Knepley ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr); 24757459e9aSMatthew G Knepley ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr); 24857459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 24957459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 25057459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 25157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 25257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 25357459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 25457459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 25557459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 25657459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 25757459e9aSMatthew G Knepley if (zfEnd != fEnd) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 25888ed4aceSMatthew G Knepley /* Create local section */ 25980800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 260a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 261a66d4d66SMatthew G Knepley if (numVertexDof) {numVertexTotDof += numVertexDof[f];} 262a66d4d66SMatthew G Knepley if (numCellDof) {numCellTotDof += numCellDof[f];} 26388ed4aceSMatthew G Knepley if (numFaceDof) {numFaceTotDof[0] += numFaceDof[f*dim+0]; 26488ed4aceSMatthew G Knepley numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0; 2650ad7597dSKarl Rupp numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0; 2660ad7597dSKarl Rupp } 267a66d4d66SMatthew G Knepley } 26888ed4aceSMatthew G Knepley ierr = PetscSectionCreate(((PetscObject) dm)->comm, &dm->defaultSection);CHKERRQ(ierr); 269a66d4d66SMatthew G Knepley if (numFields > 1) { 27088ed4aceSMatthew G Knepley ierr = PetscSectionSetNumFields(dm->defaultSection, numFields);CHKERRQ(ierr); 271a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 27280800b1aSMatthew G Knepley const char *name; 27380800b1aSMatthew G Knepley 27480800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 27588ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldName(dm->defaultSection, f, name);CHKERRQ(ierr); 27680800b1aSMatthew G Knepley if (numComp) { 27788ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldComponents(dm->defaultSection, f, numComp[f]);CHKERRQ(ierr); 278a66d4d66SMatthew G Knepley } 279a66d4d66SMatthew G Knepley } 280a66d4d66SMatthew G Knepley } else { 281a66d4d66SMatthew G Knepley numFields = 0; 282a66d4d66SMatthew G Knepley } 28388ed4aceSMatthew G Knepley ierr = PetscSectionSetChart(dm->defaultSection, pStart, pEnd);CHKERRQ(ierr); 284a66d4d66SMatthew G Knepley if (numVertexDof) { 285a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 286a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 28788ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, v, f, numVertexDof[f]);CHKERRQ(ierr); 288a66d4d66SMatthew G Knepley } 28988ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, v, numVertexTotDof);CHKERRQ(ierr); 290a66d4d66SMatthew G Knepley } 291a66d4d66SMatthew G Knepley } 292a66d4d66SMatthew G Knepley if (numFaceDof) { 293a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 294a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 29588ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr); 296a66d4d66SMatthew G Knepley } 29788ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, xf, numFaceTotDof[0]);CHKERRQ(ierr); 298a66d4d66SMatthew G Knepley } 299a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 300a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 30188ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr); 302a66d4d66SMatthew G Knepley } 30388ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, yf, numFaceTotDof[1]);CHKERRQ(ierr); 304a66d4d66SMatthew G Knepley } 305a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 306a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 30788ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr); 308a66d4d66SMatthew G Knepley } 30988ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, zf, numFaceTotDof[2]);CHKERRQ(ierr); 310a66d4d66SMatthew G Knepley } 311a66d4d66SMatthew G Knepley } 312a66d4d66SMatthew G Knepley if (numCellDof) { 313a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 314a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 31588ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, c, f, numCellDof[f]);CHKERRQ(ierr); 316a66d4d66SMatthew G Knepley } 31788ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, c, numCellTotDof);CHKERRQ(ierr); 318a66d4d66SMatthew G Knepley } 319a66d4d66SMatthew G Knepley } 32088ed4aceSMatthew G Knepley ierr = PetscSectionSetUp(dm->defaultSection);CHKERRQ(ierr); 32188ed4aceSMatthew G Knepley /* Create mesh point SF */ 32288ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 32388ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 32488ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 32588ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3267128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 32788ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 32888ed4aceSMatthew G Knepley 3293814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 330feafbc5cSMatthew G Knepley nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */ 33188ed4aceSMatthew G Knepley if (xp && !yp && !zp) { 33288ed4aceSMatthew G Knepley nleaves += nxF; /* x faces */ 33388ed4aceSMatthew G Knepley } else if (yp && !zp && !xp) { 33488ed4aceSMatthew G Knepley nleaves += nyF; /* y faces */ 33588ed4aceSMatthew G Knepley } else if (zp && !xp && !yp) { 33688ed4aceSMatthew G Knepley nleaves += nzF; /* z faces */ 33788ed4aceSMatthew G Knepley } 33888ed4aceSMatthew G Knepley } 33988ed4aceSMatthew G Knepley } 34088ed4aceSMatthew G Knepley } 34188ed4aceSMatthew G Knepley } 34288ed4aceSMatthew G Knepley ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr); 34388ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 34488ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 34588ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3467128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 34788ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 348f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 34988ed4aceSMatthew G Knepley 3503814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 35188ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 35288ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 35388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 354f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; 355f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 356628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom back vertex */ 357f5eeb527SMatthew G Knepley 358f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 359f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 360f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 361f5eeb527SMatthew G Knepley ++nL; 36288ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 363f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; 364f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 365628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom front vertex */ 366f5eeb527SMatthew G Knepley 367f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 368f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 369f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 370f5eeb527SMatthew G Knepley ++nL; 37188ed4aceSMatthew G Knepley } else { 37288ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left bottom vertices */ 373f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 374f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; 375f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 376f5eeb527SMatthew G Knepley 377f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 378f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 379f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 380f5eeb527SMatthew G Knepley } 38188ed4aceSMatthew G Knepley } 38288ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 38388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 384f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; 385f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 386628e3b0dSSatish Balay nleavesCheck += 1; /* left top back vertex */ 387f5eeb527SMatthew G Knepley 388f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 389f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 390f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 391f5eeb527SMatthew G Knepley ++nL; 39288ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 393f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; 394f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 395628e3b0dSSatish Balay nleavesCheck += 1; /* left top front vertex */ 396f5eeb527SMatthew G Knepley 397f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 398f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 399f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 400f5eeb527SMatthew G Knepley ++nL; 40188ed4aceSMatthew G Knepley } else { 40288ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left top vertices */ 403f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 404f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; 405f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 406f5eeb527SMatthew G Knepley 407f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 408f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 409f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 410f5eeb527SMatthew G Knepley } 41188ed4aceSMatthew G Knepley } 41288ed4aceSMatthew G Knepley } else { 41388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 41488ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left back vertices */ 415f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 416f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; 417f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 418f5eeb527SMatthew G Knepley 419f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 420f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 421f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 422f5eeb527SMatthew G Knepley } 42388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 42488ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left front vertices */ 425f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 426f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; 427f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 428f5eeb527SMatthew G Knepley 429f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 430f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 431f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 432f5eeb527SMatthew G Knepley } 43388ed4aceSMatthew G Knepley } else { 43488ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* left vertices */ 435f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 436f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 437f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; 438f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 439f5eeb527SMatthew G Knepley 440f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 441f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 442f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 443f5eeb527SMatthew G Knepley } 444f5eeb527SMatthew G Knepley } 44588ed4aceSMatthew G Knepley nleavesCheck += nxF; /* left faces */ 446f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 447f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 448f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 449f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 450f5eeb527SMatthew G Knepley 451f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 452f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 453f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 454f5eeb527SMatthew G Knepley } 45588ed4aceSMatthew G Knepley } 45688ed4aceSMatthew G Knepley } 45788ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 45888ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 45988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 460f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; 461f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 462628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom back vertex */ 463f5eeb527SMatthew G Knepley 464f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 465f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 466f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 467f5eeb527SMatthew G Knepley ++nL; 46888ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 469f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; 470f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 471628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom front vertex */ 472f5eeb527SMatthew G Knepley 473f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 474f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 475f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 476f5eeb527SMatthew G Knepley ++nL; 47788ed4aceSMatthew G Knepley } else { 47888ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right bottom vertices */ 479f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 480f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; 481f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 482f5eeb527SMatthew G Knepley 483f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 484f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 485f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 486f5eeb527SMatthew G Knepley } 48788ed4aceSMatthew G Knepley } 48888ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 48988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 490f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; 491f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 492628e3b0dSSatish Balay nleavesCheck += 1; /* right top back vertex */ 493f5eeb527SMatthew G Knepley 494f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 495f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 496f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 497f5eeb527SMatthew G Knepley ++nL; 49888ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 499f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; 500f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 501628e3b0dSSatish Balay nleavesCheck += 1; /* right top front vertex */ 502f5eeb527SMatthew G Knepley 503f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 504f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 505f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 506f5eeb527SMatthew G Knepley ++nL; 50788ed4aceSMatthew G Knepley } else { 50888ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right top vertices */ 509f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 510f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; 511f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 512f5eeb527SMatthew G Knepley 513f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 514f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 515f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 516f5eeb527SMatthew G Knepley } 51788ed4aceSMatthew G Knepley } 51888ed4aceSMatthew G Knepley } else { 51988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 52088ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right back vertices */ 521f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 522f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; 523f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 524f5eeb527SMatthew G Knepley 525f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 526f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 527f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 528f5eeb527SMatthew G Knepley } 52988ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 53088ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right front vertices */ 531f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 532f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; 533f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 534f5eeb527SMatthew G Knepley 535f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 536f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 537f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 538f5eeb527SMatthew G Knepley } 53988ed4aceSMatthew G Knepley } else { 54088ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* right vertices */ 541f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 542f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 543f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; 544f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 545f5eeb527SMatthew G Knepley 546f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 547f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 548f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 549f5eeb527SMatthew G Knepley } 550f5eeb527SMatthew G Knepley } 55188ed4aceSMatthew G Knepley nleavesCheck += nxF; /* right faces */ 552f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 553f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 554f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 555f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 556f5eeb527SMatthew G Knepley 557f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 558f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 559f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 560f5eeb527SMatthew G Knepley } 56188ed4aceSMatthew G Knepley } 56288ed4aceSMatthew G Knepley } 56388ed4aceSMatthew G Knepley } else { 56488ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 56588ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 56688ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom back vertices */ 567f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 568f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; 569f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 570f5eeb527SMatthew G Knepley 571f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 572f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 573f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 574f5eeb527SMatthew G Knepley } 57588ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 57688ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom front vertices */ 577f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 578f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; 579f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 580f5eeb527SMatthew G Knepley 581f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 582f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 583f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 584f5eeb527SMatthew G Knepley } 58588ed4aceSMatthew G Knepley } else { 58688ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* bottom vertices */ 587f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 588f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 589f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; 590f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 591f5eeb527SMatthew G Knepley 592f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 593f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 594f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 595f5eeb527SMatthew G Knepley } 596f5eeb527SMatthew G Knepley } 59788ed4aceSMatthew G Knepley nleavesCheck += nyF; /* bottom faces */ 598f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 599f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 600f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 601f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 602f5eeb527SMatthew G Knepley 603f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 604f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 605f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 606f5eeb527SMatthew G Knepley } 60788ed4aceSMatthew G Knepley } 60888ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 60988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 61088ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top back vertices */ 611f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 612f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; 613f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 614f5eeb527SMatthew G Knepley 615f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 616f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 617f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 618f5eeb527SMatthew G Knepley } 61988ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 62088ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top front vertices */ 621f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 622f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; 623f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 624f5eeb527SMatthew G Knepley 625f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 626f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 627f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 628f5eeb527SMatthew G Knepley } 62988ed4aceSMatthew G Knepley } else { 63088ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* top vertices */ 631f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 632f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 633f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; 634f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 635f5eeb527SMatthew G Knepley 636f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 637f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 638f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 639f5eeb527SMatthew G Knepley } 640f5eeb527SMatthew G Knepley } 64188ed4aceSMatthew G Knepley nleavesCheck += nyF; /* top faces */ 642f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 643f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 644f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 645f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 646f5eeb527SMatthew G Knepley 647f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 648f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 649f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 650f5eeb527SMatthew G Knepley } 65188ed4aceSMatthew G Knepley } 65288ed4aceSMatthew G Knepley } else { 65388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 65488ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* back vertices */ 655f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 656f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 657f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; 658f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 659f5eeb527SMatthew G Knepley 660f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 661f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 662f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 663f5eeb527SMatthew G Knepley } 664f5eeb527SMatthew G Knepley } 66588ed4aceSMatthew G Knepley nleavesCheck += nzF; /* back faces */ 666f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 667f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 668f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 669f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 670f5eeb527SMatthew G Knepley 671f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 672f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 673f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 674f5eeb527SMatthew G Knepley } 67588ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 67688ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* front vertices */ 677f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 678f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 679f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; 680f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 681f5eeb527SMatthew G Knepley 682f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 683f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 684f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 685f5eeb527SMatthew G Knepley } 686f5eeb527SMatthew G Knepley } 68788ed4aceSMatthew G Knepley nleavesCheck += nzF; /* front faces */ 688f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 689f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 690f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 691f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 692f5eeb527SMatthew G Knepley 693f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 694f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 695f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 696f5eeb527SMatthew G Knepley } 69788ed4aceSMatthew G Knepley } else { 69888ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 69988ed4aceSMatthew G Knepley } 70088ed4aceSMatthew G Knepley } 70188ed4aceSMatthew G Knepley } 70288ed4aceSMatthew G Knepley } 70388ed4aceSMatthew G Knepley } 70488ed4aceSMatthew G Knepley } 70588ed4aceSMatthew G Knepley } 7063814d604SMatthew G Knepley /* TODO: Remove duplication in leaf determination */ 70788ed4aceSMatthew G Knepley if (nleaves != nleavesCheck) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck); 70888ed4aceSMatthew G Knepley ierr = PetscSFCreate(((PetscObject) dm)->comm, &sf);CHKERRQ(ierr); 70988ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 71088ed4aceSMatthew G Knepley /* Create global section */ 711eaf8d80aSMatthew G Knepley ierr = PetscSectionCreateGlobalSection(dm->defaultSection, sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr); 71288ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 71388ed4aceSMatthew G Knepley /* Create default SF */ 71488ed4aceSMatthew G Knepley ierr = DMCreateDefaultSF(dm, dm->defaultSection, dm->defaultGlobalSection);CHKERRQ(ierr); 715a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 716a66d4d66SMatthew G Knepley } 717a66d4d66SMatthew G Knepley 71847c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 71947c6ae99SBarry Smith 72047c6ae99SBarry Smith #undef __FUNCT__ 721aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray" 72247c6ae99SBarry Smith /*@C 723aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith Input Parameter: 72647c6ae99SBarry Smith + da - information about my local patch 72747c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith Output Parameters: 73047c6ae99SBarry Smith . vptr - array data structured 73147c6ae99SBarry Smith 73247c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 73347c6ae99SBarry Smith to zero them. 73447c6ae99SBarry Smith 73547c6ae99SBarry Smith Level: advanced 73647c6ae99SBarry Smith 737bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith @*/ 7407087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 74147c6ae99SBarry Smith { 74247c6ae99SBarry Smith PetscErrorCode ierr; 74347c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 74447c6ae99SBarry Smith char *iarray_start; 74547c6ae99SBarry Smith void **iptr = (void**)vptr; 74647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 74747c6ae99SBarry Smith 74847c6ae99SBarry Smith PetscFunctionBegin; 74947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 75047c6ae99SBarry Smith if (ghosted) { 751aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 75247c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 75347c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 75447c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 75547c6ae99SBarry Smith dd->arrayghostedin[i] = PETSC_NULL; 75647c6ae99SBarry Smith dd->startghostedin[i] = PETSC_NULL; 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith goto done; 75947c6ae99SBarry Smith } 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith xs = dd->Xs; 76247c6ae99SBarry Smith ys = dd->Ys; 76347c6ae99SBarry Smith zs = dd->Zs; 76447c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 76547c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 76647c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 76747c6ae99SBarry Smith } else { 768aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 76947c6ae99SBarry Smith if (dd->arrayin[i]) { 77047c6ae99SBarry Smith *iptr = dd->arrayin[i]; 77147c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 77247c6ae99SBarry Smith dd->arrayin[i] = PETSC_NULL; 77347c6ae99SBarry Smith dd->startin[i] = PETSC_NULL; 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith goto done; 77647c6ae99SBarry Smith } 77747c6ae99SBarry Smith } 77847c6ae99SBarry Smith xs = dd->xs; 77947c6ae99SBarry Smith ys = dd->ys; 78047c6ae99SBarry Smith zs = dd->zs; 78147c6ae99SBarry Smith xm = dd->xe-dd->xs; 78247c6ae99SBarry Smith ym = dd->ye-dd->ys; 78347c6ae99SBarry Smith zm = dd->ze-dd->zs; 78447c6ae99SBarry Smith } 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith switch (dd->dim) { 78747c6ae99SBarry Smith case 1: { 78847c6ae99SBarry Smith void *ptr; 78947c6ae99SBarry Smith 79047c6ae99SBarry Smith ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 79347c6ae99SBarry Smith *iptr = (void*)ptr; 79447c6ae99SBarry Smith break;} 79547c6ae99SBarry Smith case 2: { 79647c6ae99SBarry Smith void **ptr; 79747c6ae99SBarry Smith 79847c6ae99SBarry Smith ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 79947c6ae99SBarry Smith 80047c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 80147c6ae99SBarry Smith for (j=ys;j<ys+ym;j++) { 80247c6ae99SBarry Smith ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 80347c6ae99SBarry Smith } 80447c6ae99SBarry Smith *iptr = (void*)ptr; 80547c6ae99SBarry Smith break;} 80647c6ae99SBarry Smith case 3: { 80747c6ae99SBarry Smith void ***ptr,**bptr; 80847c6ae99SBarry Smith 80947c6ae99SBarry Smith ierr = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 81047c6ae99SBarry Smith 81147c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 81247c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 81347c6ae99SBarry Smith for (i=zs;i<zs+zm;i++) { 81447c6ae99SBarry Smith ptr[i] = bptr + ((i-zs)*ym - ys); 81547c6ae99SBarry Smith } 81647c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 81747c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 81847c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 81947c6ae99SBarry Smith } 82047c6ae99SBarry Smith } 82147c6ae99SBarry Smith 82247c6ae99SBarry Smith *iptr = (void*)ptr; 82347c6ae99SBarry Smith break;} 82447c6ae99SBarry Smith default: 82547c6ae99SBarry Smith SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith 82847c6ae99SBarry Smith done: 82947c6ae99SBarry Smith /* add arrays to the checked out list */ 83047c6ae99SBarry Smith if (ghosted) { 831aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 83247c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 83347c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr ; 83447c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 83547c6ae99SBarry Smith break; 83647c6ae99SBarry Smith } 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith } else { 839aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 84047c6ae99SBarry Smith if (!dd->arrayout[i]) { 84147c6ae99SBarry Smith dd->arrayout[i] = *iptr ; 84247c6ae99SBarry Smith dd->startout[i] = iarray_start; 84347c6ae99SBarry Smith break; 84447c6ae99SBarry Smith } 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith } 84747c6ae99SBarry Smith PetscFunctionReturn(0); 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith 85047c6ae99SBarry Smith #undef __FUNCT__ 851aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray" 85247c6ae99SBarry Smith /*@C 853aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith Input Parameter: 85647c6ae99SBarry Smith + da - information about my local patch 85747c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 85847c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith Level: advanced 86147c6ae99SBarry Smith 862bcaeba4dSBarry Smith .seealso: DMDAGetArray() 86347c6ae99SBarry Smith 86447c6ae99SBarry Smith @*/ 8657087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 86647c6ae99SBarry Smith { 86747c6ae99SBarry Smith PetscInt i; 86847c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 86947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 87047c6ae99SBarry Smith 87147c6ae99SBarry Smith PetscFunctionBegin; 87247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 87347c6ae99SBarry Smith if (ghosted) { 874aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 87547c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 87647c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 87747c6ae99SBarry Smith dd->arrayghostedout[i] = PETSC_NULL; 87847c6ae99SBarry Smith dd->startghostedout[i] = PETSC_NULL; 87947c6ae99SBarry Smith break; 88047c6ae99SBarry Smith } 88147c6ae99SBarry Smith } 882aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 88347c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 88447c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 88547c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 88647c6ae99SBarry Smith break; 88747c6ae99SBarry Smith } 88847c6ae99SBarry Smith } 88947c6ae99SBarry Smith } else { 890aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 89147c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 89247c6ae99SBarry Smith iarray_start = dd->startout[i]; 89347c6ae99SBarry Smith dd->arrayout[i] = PETSC_NULL; 89447c6ae99SBarry Smith dd->startout[i] = PETSC_NULL; 89547c6ae99SBarry Smith break; 89647c6ae99SBarry Smith } 89747c6ae99SBarry Smith } 898aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 89947c6ae99SBarry Smith if (!dd->arrayin[i]) { 90047c6ae99SBarry Smith dd->arrayin[i] = *iptr; 90147c6ae99SBarry Smith dd->startin[i] = iarray_start; 90247c6ae99SBarry Smith break; 90347c6ae99SBarry Smith } 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith } 90647c6ae99SBarry Smith PetscFunctionReturn(0); 90747c6ae99SBarry Smith } 90847c6ae99SBarry Smith 909