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*/ 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 #undef __FUNCT__ 1547c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d" 16*1e6b0712SBarry Smith static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 1747c6ae99SBarry Smith { 1847c6ae99SBarry Smith PetscErrorCode ierr; 1947c6ae99SBarry Smith PetscInt n,m; 2047c6ae99SBarry Smith Vec vec = (Vec)obj; 2147c6ae99SBarry Smith PetscScalar *array; 2247c6ae99SBarry Smith mxArray *mat; 239a42bb27SBarry Smith DM da; 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith PetscFunctionBegin; 26c688c046SMatthew G Knepley ierr = VecGetDM(vec, &da);CHKERRQ(ierr); 27ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 28aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 2947c6ae99SBarry Smith 3047c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 3147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3247c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3347c6ae99SBarry Smith #else 3447c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3547c6ae99SBarry Smith #endif 3647c6ae99SBarry Smith ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 3747c6ae99SBarry Smith ierr = PetscObjectName(obj);CHKERRQ(ierr); 3847c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 3947c6ae99SBarry Smith 4047c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 4147c6ae99SBarry Smith PetscFunctionReturn(0); 4247c6ae99SBarry Smith } 4347c6ae99SBarry Smith #endif 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith #undef __FUNCT__ 47564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA" 487087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 4947c6ae99SBarry Smith { 5047c6ae99SBarry Smith PetscErrorCode ierr; 5147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith PetscFunctionBegin; 5447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5547c6ae99SBarry Smith PetscValidPointer(g,2); 5611689aebSJed Brown if (da->defaultSection) { 5711689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 5811689aebSJed Brown } else { 5947c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6047c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6147c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 62401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 63c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6447c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 6547c6ae99SBarry Smith if (dd->w == 1 && dd->dim == 2) { 6600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 6747c6ae99SBarry Smith } 6847c6ae99SBarry Smith #endif 6911689aebSJed Brown } 7047c6ae99SBarry Smith PetscFunctionReturn(0); 7147c6ae99SBarry Smith } 7247c6ae99SBarry Smith 73a66d4d66SMatthew G Knepley #undef __FUNCT__ 7457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells" 7557459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCells) 7657459e9aSMatthew G Knepley { 7757459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 7857459e9aSMatthew G Knepley const PetscInt dim = da->dim; 7957459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 8057459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 8157459e9aSMatthew G Knepley 8257459e9aSMatthew G Knepley PetscFunctionBegin; 8357459e9aSMatthew G Knepley if (numCells) { 8457459e9aSMatthew G Knepley PetscValidIntPointer(numCells,2); 8557459e9aSMatthew G Knepley *numCells = nC; 8657459e9aSMatthew G Knepley } 8757459e9aSMatthew G Knepley PetscFunctionReturn(0); 8857459e9aSMatthew G Knepley } 8957459e9aSMatthew G Knepley 9057459e9aSMatthew G Knepley #undef __FUNCT__ 9157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices" 9257459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 9357459e9aSMatthew G Knepley { 9457459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 9557459e9aSMatthew G Knepley const PetscInt dim = da->dim; 9657459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 9757459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 9857459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 9957459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 10057459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 10157459e9aSMatthew G Knepley 10257459e9aSMatthew G Knepley PetscFunctionBegin; 10357459e9aSMatthew G Knepley if (numVerticesX) { 10457459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 10557459e9aSMatthew G Knepley *numVerticesX = nVx; 10657459e9aSMatthew G Knepley } 10757459e9aSMatthew G Knepley if (numVerticesY) { 10857459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 10957459e9aSMatthew G Knepley *numVerticesY = nVy; 11057459e9aSMatthew G Knepley } 11157459e9aSMatthew G Knepley if (numVerticesZ) { 11257459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 11357459e9aSMatthew G Knepley *numVerticesZ = nVz; 11457459e9aSMatthew G Knepley } 11557459e9aSMatthew G Knepley if (numVertices) { 11657459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 11757459e9aSMatthew G Knepley *numVertices = nV; 11857459e9aSMatthew G Knepley } 11957459e9aSMatthew G Knepley PetscFunctionReturn(0); 12057459e9aSMatthew G Knepley } 12157459e9aSMatthew G Knepley 12257459e9aSMatthew G Knepley #undef __FUNCT__ 12357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces" 12457459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 12557459e9aSMatthew G Knepley { 12657459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 12757459e9aSMatthew G Knepley const PetscInt dim = da->dim; 12857459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 12957459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 13057459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 13157459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 13257459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 13357459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 13457459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 13557459e9aSMatthew G Knepley 13657459e9aSMatthew G Knepley PetscFunctionBegin; 13757459e9aSMatthew G Knepley if (numXFacesX) { 13857459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 13957459e9aSMatthew G Knepley *numXFacesX = nxF; 14057459e9aSMatthew G Knepley } 14157459e9aSMatthew G Knepley if (numXFaces) { 14257459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 14357459e9aSMatthew G Knepley *numXFaces = nXF; 14457459e9aSMatthew G Knepley } 14557459e9aSMatthew G Knepley if (numYFacesY) { 14657459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 14757459e9aSMatthew G Knepley *numYFacesY = nyF; 14857459e9aSMatthew G Knepley } 14957459e9aSMatthew G Knepley if (numYFaces) { 15057459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 15157459e9aSMatthew G Knepley *numYFaces = nYF; 15257459e9aSMatthew G Knepley } 15357459e9aSMatthew G Knepley if (numZFacesZ) { 15457459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 15557459e9aSMatthew G Knepley *numZFacesZ = nzF; 15657459e9aSMatthew G Knepley } 15757459e9aSMatthew G Knepley if (numZFaces) { 15857459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 15957459e9aSMatthew G Knepley *numZFaces = nZF; 16057459e9aSMatthew G Knepley } 16157459e9aSMatthew G Knepley PetscFunctionReturn(0); 16257459e9aSMatthew G Knepley } 16357459e9aSMatthew G Knepley 16457459e9aSMatthew G Knepley #undef __FUNCT__ 16557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum" 16657459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 16757459e9aSMatthew G Knepley { 16857459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 16957459e9aSMatthew G Knepley const PetscInt dim = da->dim; 17057459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 17157459e9aSMatthew G Knepley PetscErrorCode ierr; 17257459e9aSMatthew G Knepley 17357459e9aSMatthew G Knepley PetscFunctionBegin; 1748865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 1758865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 17657459e9aSMatthew G Knepley ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr); 1770298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 1780298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 17957459e9aSMatthew G Knepley if (height == 0) { 18057459e9aSMatthew G Knepley /* Cells */ 1818865f1eaSKarl Rupp if (pStart) *pStart = 0; 1828865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 18357459e9aSMatthew G Knepley } else if (height == 1) { 18457459e9aSMatthew G Knepley /* Faces */ 1858865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 1868865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 18757459e9aSMatthew G Knepley } else if (height == dim) { 18857459e9aSMatthew G Knepley /* Vertices */ 1898865f1eaSKarl Rupp if (pStart) *pStart = nC; 1908865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 19157459e9aSMatthew G Knepley } else if (height < 0) { 19257459e9aSMatthew G Knepley /* All points */ 1938865f1eaSKarl Rupp if (pStart) *pStart = 0; 1948865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 19582f516ccSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 19657459e9aSMatthew G Knepley PetscFunctionReturn(0); 19757459e9aSMatthew G Knepley } 19857459e9aSMatthew G Knepley 19957459e9aSMatthew G Knepley #undef __FUNCT__ 200a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection" 201a66d4d66SMatthew G Knepley /*@C 202a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 203a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 204a66d4d66SMatthew G Knepley 205a66d4d66SMatthew G Knepley Input Parameters: 206a66d4d66SMatthew G Knepley + dm- The DMDA 207a66d4d66SMatthew G Knepley . numFields - The number of fields 2080298fd71SBarry Smith . numComp - The number of components in each field, or NULL for 1 2090298fd71SBarry Smith . numVertexDof - The number of dofs per vertex for each field, or NULL 2100298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL 2110298fd71SBarry Smith - numCellDof - The number of dofs per cell for each field, or NULL 212a66d4d66SMatthew G Knepley 213a66d4d66SMatthew G Knepley Level: developer 214a66d4d66SMatthew G Knepley 215a66d4d66SMatthew G Knepley Note: 216a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 217a66d4d66SMatthew G Knepley 218a66d4d66SMatthew G Knepley - Cells: [0, nC) 219a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 22088ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 22188ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 22288ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 223a66d4d66SMatthew G Knepley 224a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 225a66d4d66SMatthew G Knepley @*/ 22680800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[]) 227a66d4d66SMatthew G Knepley { 228a66d4d66SMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 22988ed4aceSMatthew G Knepley const PetscInt dim = da->dim; 23080800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 23188ed4aceSMatthew G Knepley PetscSF sf; 232feafbc5cSMatthew G Knepley PetscMPIInt rank; 23388ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 23488ed4aceSMatthew G Knepley PetscInt *localPoints; 23588ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 236f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 23757459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 23857459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 23988ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 240a66d4d66SMatthew G Knepley PetscErrorCode ierr; 241a66d4d66SMatthew G Knepley 242a66d4d66SMatthew G Knepley PetscFunctionBegin; 243a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 24482f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 24557459e9aSMatthew G Knepley ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr); 24657459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 24757459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 24857459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 24957459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 25057459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 25157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 25257459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 25357459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 25457459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 25582f516ccSBarry Smith if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 25688ed4aceSMatthew G Knepley /* Create local section */ 25780800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 258a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 2598865f1eaSKarl Rupp if (numVertexDof) numVertexTotDof += numVertexDof[f]; 2608865f1eaSKarl Rupp if (numCellDof) numCellTotDof += numCellDof[f]; 2618865f1eaSKarl Rupp if (numFaceDof) { 2628865f1eaSKarl Rupp numFaceTotDof[0] += numFaceDof[f*dim+0]; 26388ed4aceSMatthew G Knepley numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0; 2640ad7597dSKarl Rupp numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0; 2650ad7597dSKarl Rupp } 266a66d4d66SMatthew G Knepley } 26782f516ccSBarry Smith ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &dm->defaultSection);CHKERRQ(ierr); 268a66d4d66SMatthew G Knepley if (numFields > 1) { 26988ed4aceSMatthew G Knepley ierr = PetscSectionSetNumFields(dm->defaultSection, numFields);CHKERRQ(ierr); 270a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 27180800b1aSMatthew G Knepley const char *name; 27280800b1aSMatthew G Knepley 27380800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 27488ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldName(dm->defaultSection, f, name);CHKERRQ(ierr); 27580800b1aSMatthew G Knepley if (numComp) { 27688ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldComponents(dm->defaultSection, f, numComp[f]);CHKERRQ(ierr); 277a66d4d66SMatthew G Knepley } 278a66d4d66SMatthew G Knepley } 279a66d4d66SMatthew G Knepley } else { 280a66d4d66SMatthew G Knepley numFields = 0; 281a66d4d66SMatthew G Knepley } 28288ed4aceSMatthew G Knepley ierr = PetscSectionSetChart(dm->defaultSection, pStart, pEnd);CHKERRQ(ierr); 283a66d4d66SMatthew G Knepley if (numVertexDof) { 284a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 285a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 28688ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, v, f, numVertexDof[f]);CHKERRQ(ierr); 287a66d4d66SMatthew G Knepley } 28888ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, v, numVertexTotDof);CHKERRQ(ierr); 289a66d4d66SMatthew G Knepley } 290a66d4d66SMatthew G Knepley } 291a66d4d66SMatthew G Knepley if (numFaceDof) { 292a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 293a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 29488ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr); 295a66d4d66SMatthew G Knepley } 29688ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, xf, numFaceTotDof[0]);CHKERRQ(ierr); 297a66d4d66SMatthew G Knepley } 298a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 299a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 30088ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr); 301a66d4d66SMatthew G Knepley } 30288ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, yf, numFaceTotDof[1]);CHKERRQ(ierr); 303a66d4d66SMatthew G Knepley } 304a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 305a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 30688ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr); 307a66d4d66SMatthew G Knepley } 30888ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, zf, numFaceTotDof[2]);CHKERRQ(ierr); 309a66d4d66SMatthew G Knepley } 310a66d4d66SMatthew G Knepley } 311a66d4d66SMatthew G Knepley if (numCellDof) { 312a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 313a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 31488ed4aceSMatthew G Knepley ierr = PetscSectionSetFieldDof(dm->defaultSection, c, f, numCellDof[f]);CHKERRQ(ierr); 315a66d4d66SMatthew G Knepley } 31688ed4aceSMatthew G Knepley ierr = PetscSectionSetDof(dm->defaultSection, c, numCellTotDof);CHKERRQ(ierr); 317a66d4d66SMatthew G Knepley } 318a66d4d66SMatthew G Knepley } 31988ed4aceSMatthew G Knepley ierr = PetscSectionSetUp(dm->defaultSection);CHKERRQ(ierr); 32088ed4aceSMatthew G Knepley /* Create mesh point SF */ 32188ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 32288ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 32388ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 32488ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3257128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 32688ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 32788ed4aceSMatthew G Knepley 3283814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 329feafbc5cSMatthew G Knepley nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */ 33088ed4aceSMatthew G Knepley if (xp && !yp && !zp) { 33188ed4aceSMatthew G Knepley nleaves += nxF; /* x faces */ 33288ed4aceSMatthew G Knepley } else if (yp && !zp && !xp) { 33388ed4aceSMatthew G Knepley nleaves += nyF; /* y faces */ 33488ed4aceSMatthew G Knepley } else if (zp && !xp && !yp) { 33588ed4aceSMatthew G Knepley nleaves += nzF; /* z faces */ 33688ed4aceSMatthew G Knepley } 33788ed4aceSMatthew G Knepley } 33888ed4aceSMatthew G Knepley } 33988ed4aceSMatthew G Knepley } 34088ed4aceSMatthew G Knepley } 34188ed4aceSMatthew G Knepley ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr); 34288ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 34388ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 34488ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 3457128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 34688ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 347f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 34888ed4aceSMatthew G Knepley 3493814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 35088ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 35188ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 35288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 353f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; 354f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 355628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom back vertex */ 356f5eeb527SMatthew G Knepley 357f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 358f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 359f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 360f5eeb527SMatthew G Knepley ++nL; 36188ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 362f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; 363f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 364628e3b0dSSatish Balay nleavesCheck += 1; /* left bottom front vertex */ 365f5eeb527SMatthew G Knepley 366f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 367f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 368f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 369f5eeb527SMatthew G Knepley ++nL; 37088ed4aceSMatthew G Knepley } else { 37188ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left bottom vertices */ 372f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 373f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; 374f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 375f5eeb527SMatthew G Knepley 376f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 377f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 378f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 379f5eeb527SMatthew G Knepley } 38088ed4aceSMatthew G Knepley } 38188ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 38288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 383f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; 384f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 385628e3b0dSSatish Balay nleavesCheck += 1; /* left top back vertex */ 386f5eeb527SMatthew G Knepley 387f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 388f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 389f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 390f5eeb527SMatthew G Knepley ++nL; 39188ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 392f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; 393f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 394628e3b0dSSatish Balay nleavesCheck += 1; /* left top front vertex */ 395f5eeb527SMatthew G Knepley 396f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 397f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 398f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 399f5eeb527SMatthew G Knepley ++nL; 40088ed4aceSMatthew G Knepley } else { 40188ed4aceSMatthew G Knepley nleavesCheck += nVz; /* left top vertices */ 402f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 403f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; 404f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 405f5eeb527SMatthew G Knepley 406f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 407f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 408f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 409f5eeb527SMatthew G Knepley } 41088ed4aceSMatthew G Knepley } 41188ed4aceSMatthew G Knepley } else { 41288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 41388ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left back vertices */ 414f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 415f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; 416f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 417f5eeb527SMatthew G Knepley 418f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 419f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 420f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 421f5eeb527SMatthew G Knepley } 42288ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 42388ed4aceSMatthew G Knepley nleavesCheck += nVy; /* left front vertices */ 424f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 425f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; 426f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 427f5eeb527SMatthew G Knepley 428f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 429f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 430f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 431f5eeb527SMatthew G Knepley } 43288ed4aceSMatthew G Knepley } else { 43388ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* left vertices */ 434f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 435f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 436f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; 437f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 438f5eeb527SMatthew G Knepley 439f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 440f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 441f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 442f5eeb527SMatthew G Knepley } 443f5eeb527SMatthew G Knepley } 44488ed4aceSMatthew G Knepley nleavesCheck += nxF; /* left faces */ 445f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 446f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 447f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 448f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 449f5eeb527SMatthew G Knepley 450f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 451f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 452f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 453f5eeb527SMatthew G Knepley } 45488ed4aceSMatthew G Knepley } 45588ed4aceSMatthew G Knepley } 45688ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 45788ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 45888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 459f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; 460f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 461628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom back vertex */ 462f5eeb527SMatthew G Knepley 463f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 464f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 465f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 466f5eeb527SMatthew G Knepley ++nL; 46788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 468f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; 469f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 470628e3b0dSSatish Balay nleavesCheck += 1; /* right bottom front vertex */ 471f5eeb527SMatthew G Knepley 472f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 473f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 474f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 475f5eeb527SMatthew G Knepley ++nL; 47688ed4aceSMatthew G Knepley } else { 47788ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right bottom vertices */ 478f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 479f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; 480f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 481f5eeb527SMatthew G Knepley 482f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 483f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 484f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 485f5eeb527SMatthew G Knepley } 48688ed4aceSMatthew G Knepley } 48788ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 48888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 489f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; 490f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 491628e3b0dSSatish Balay nleavesCheck += 1; /* right top back vertex */ 492f5eeb527SMatthew G Knepley 493f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 494f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 495f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 496f5eeb527SMatthew G Knepley ++nL; 49788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 498f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; 499f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 500628e3b0dSSatish Balay nleavesCheck += 1; /* right top front vertex */ 501f5eeb527SMatthew G Knepley 502f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 503f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 504f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 505f5eeb527SMatthew G Knepley ++nL; 50688ed4aceSMatthew G Knepley } else { 50788ed4aceSMatthew G Knepley nleavesCheck += nVz; /* right top vertices */ 508f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv, ++nL) { 509f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; 510f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 511f5eeb527SMatthew G Knepley 512f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 513f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 514f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 515f5eeb527SMatthew G Knepley } 51688ed4aceSMatthew G Knepley } 51788ed4aceSMatthew G Knepley } else { 51888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 51988ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right back vertices */ 520f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 521f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; 522f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 523f5eeb527SMatthew G Knepley 524f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 525f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 526f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 527f5eeb527SMatthew G Knepley } 52888ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 52988ed4aceSMatthew G Knepley nleavesCheck += nVy; /* right front vertices */ 530f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 531f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; 532f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 533f5eeb527SMatthew G Knepley 534f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 535f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 536f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 537f5eeb527SMatthew G Knepley } 53888ed4aceSMatthew G Knepley } else { 53988ed4aceSMatthew G Knepley nleavesCheck += nVy*nVz; /* right vertices */ 540f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 541f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv, ++nL) { 542f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; 543f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 544f5eeb527SMatthew G Knepley 545f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 546f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 547f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 548f5eeb527SMatthew G Knepley } 549f5eeb527SMatthew G Knepley } 55088ed4aceSMatthew G Knepley nleavesCheck += nxF; /* right faces */ 551f5eeb527SMatthew G Knepley for (xf = 0; xf < nxF; ++xf, ++nL) { 552f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 553f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 554f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 555f5eeb527SMatthew G Knepley 556f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 557f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 558f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 559f5eeb527SMatthew G Knepley } 56088ed4aceSMatthew G Knepley } 56188ed4aceSMatthew G Knepley } 56288ed4aceSMatthew G Knepley } else { 56388ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 56488ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 56588ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom back vertices */ 566f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 567f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; 568f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 569f5eeb527SMatthew G Knepley 570f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 571f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 572f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 573f5eeb527SMatthew G Knepley } 57488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 57588ed4aceSMatthew G Knepley nleavesCheck += nVx; /* bottom front vertices */ 576f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 577f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; 578f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 579f5eeb527SMatthew G Knepley 580f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 581f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 582f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 583f5eeb527SMatthew G Knepley } 58488ed4aceSMatthew G Knepley } else { 58588ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* bottom vertices */ 586f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 587f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 588f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; 589f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 590f5eeb527SMatthew G Knepley 591f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 592f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 593f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 594f5eeb527SMatthew G Knepley } 595f5eeb527SMatthew G Knepley } 59688ed4aceSMatthew G Knepley nleavesCheck += nyF; /* bottom faces */ 597f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 598f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 599f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 600f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 601f5eeb527SMatthew G Knepley 602f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 603f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 604f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 605f5eeb527SMatthew G Knepley } 60688ed4aceSMatthew G Knepley } 60788ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 60888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 60988ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top back vertices */ 610f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 611f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; 612f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 613f5eeb527SMatthew G Knepley 614f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 615f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 616f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 617f5eeb527SMatthew G Knepley } 61888ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 61988ed4aceSMatthew G Knepley nleavesCheck += nVx; /* top front vertices */ 620f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 621f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; 622f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 623f5eeb527SMatthew G Knepley 624f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 625f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 626f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 627f5eeb527SMatthew G Knepley } 62888ed4aceSMatthew G Knepley } else { 62988ed4aceSMatthew G Knepley nleavesCheck += nVx*nVz; /* top vertices */ 630f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 631f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 632f5eeb527SMatthew G Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; 633f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 634f5eeb527SMatthew G Knepley 635f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 636f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 637f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 638f5eeb527SMatthew G Knepley } 639f5eeb527SMatthew G Knepley } 64088ed4aceSMatthew G Knepley nleavesCheck += nyF; /* top faces */ 641f5eeb527SMatthew G Knepley for (yf = 0; yf < nyF; ++yf, ++nL) { 642f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 643f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 644f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 645f5eeb527SMatthew G Knepley 646f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 647f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 648f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 649f5eeb527SMatthew G Knepley } 65088ed4aceSMatthew G Knepley } 65188ed4aceSMatthew G Knepley } else { 65288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 65388ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* back vertices */ 654f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 655f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 656f5eeb527SMatthew G Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; 657f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 658f5eeb527SMatthew G Knepley 659f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 660f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 661f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 662f5eeb527SMatthew G Knepley } 663f5eeb527SMatthew G Knepley } 66488ed4aceSMatthew G Knepley nleavesCheck += nzF; /* back faces */ 665f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 666f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 667f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 668f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 669f5eeb527SMatthew G Knepley 670f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 671f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 672f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 673f5eeb527SMatthew G Knepley } 67488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 67588ed4aceSMatthew G Knepley nleavesCheck += nVx*nVy; /* front vertices */ 676f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 677f5eeb527SMatthew G Knepley for (xv = 0; xv < nVx; ++xv, ++nL) { 678f5eeb527SMatthew G Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; 679f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 680f5eeb527SMatthew G Knepley 681f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 682f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 683f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 684f5eeb527SMatthew G Knepley } 685f5eeb527SMatthew G Knepley } 68688ed4aceSMatthew G Knepley nleavesCheck += nzF; /* front faces */ 687f5eeb527SMatthew G Knepley for (zf = 0; zf < nzF; ++zf, ++nL) { 688f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 689f5eeb527SMatthew G Knepley const PetscInt localFace = 0 + nC+nV; 690f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 691f5eeb527SMatthew G Knepley 692f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 693f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 694f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 695f5eeb527SMatthew G Knepley } 69688ed4aceSMatthew G Knepley } else { 69788ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 69888ed4aceSMatthew G Knepley } 69988ed4aceSMatthew G Knepley } 70088ed4aceSMatthew G Knepley } 70188ed4aceSMatthew G Knepley } 70288ed4aceSMatthew G Knepley } 70388ed4aceSMatthew G Knepley } 70488ed4aceSMatthew G Knepley } 7053814d604SMatthew G Knepley /* TODO: Remove duplication in leaf determination */ 70682f516ccSBarry 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); 70782f516ccSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr); 70888ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 70988ed4aceSMatthew G Knepley /* Create global section */ 710eaf8d80aSMatthew G Knepley ierr = PetscSectionCreateGlobalSection(dm->defaultSection, sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr); 71188ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 71288ed4aceSMatthew G Knepley /* Create default SF */ 71388ed4aceSMatthew G Knepley ierr = DMCreateDefaultSF(dm, dm->defaultSection, dm->defaultGlobalSection);CHKERRQ(ierr); 714a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 715a66d4d66SMatthew G Knepley } 716a66d4d66SMatthew G Knepley 71747c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 71847c6ae99SBarry Smith 71947c6ae99SBarry Smith #undef __FUNCT__ 720aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray" 72147c6ae99SBarry Smith /*@C 722aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 72347c6ae99SBarry Smith 72447c6ae99SBarry Smith Input Parameter: 72547c6ae99SBarry Smith + da - information about my local patch 72647c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 72747c6ae99SBarry Smith 72847c6ae99SBarry Smith Output Parameters: 72947c6ae99SBarry Smith . vptr - array data structured 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 73247c6ae99SBarry Smith to zero them. 73347c6ae99SBarry Smith 73447c6ae99SBarry Smith Level: advanced 73547c6ae99SBarry Smith 736bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 73747c6ae99SBarry Smith 73847c6ae99SBarry Smith @*/ 7397087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 74047c6ae99SBarry Smith { 74147c6ae99SBarry Smith PetscErrorCode ierr; 74247c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 74347c6ae99SBarry Smith char *iarray_start; 74447c6ae99SBarry Smith void **iptr = (void**)vptr; 74547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 74647c6ae99SBarry Smith 74747c6ae99SBarry Smith PetscFunctionBegin; 74847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 74947c6ae99SBarry Smith if (ghosted) { 750aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 75147c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 75247c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 75347c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 7540298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 7550298fd71SBarry Smith dd->startghostedin[i] = NULL; 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith goto done; 75847c6ae99SBarry Smith } 75947c6ae99SBarry Smith } 76047c6ae99SBarry Smith xs = dd->Xs; 76147c6ae99SBarry Smith ys = dd->Ys; 76247c6ae99SBarry Smith zs = dd->Zs; 76347c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 76447c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 76547c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 76647c6ae99SBarry Smith } else { 767aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 76847c6ae99SBarry Smith if (dd->arrayin[i]) { 76947c6ae99SBarry Smith *iptr = dd->arrayin[i]; 77047c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 7710298fd71SBarry Smith dd->arrayin[i] = NULL; 7720298fd71SBarry Smith dd->startin[i] = NULL; 77347c6ae99SBarry Smith 77447c6ae99SBarry Smith goto done; 77547c6ae99SBarry Smith } 77647c6ae99SBarry Smith } 77747c6ae99SBarry Smith xs = dd->xs; 77847c6ae99SBarry Smith ys = dd->ys; 77947c6ae99SBarry Smith zs = dd->zs; 78047c6ae99SBarry Smith xm = dd->xe-dd->xs; 78147c6ae99SBarry Smith ym = dd->ye-dd->ys; 78247c6ae99SBarry Smith zm = dd->ze-dd->zs; 78347c6ae99SBarry Smith } 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith switch (dd->dim) { 78647c6ae99SBarry Smith case 1: { 78747c6ae99SBarry Smith void *ptr; 78847c6ae99SBarry Smith 78947c6ae99SBarry Smith ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 79047c6ae99SBarry Smith 79147c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 79247c6ae99SBarry Smith *iptr = (void*)ptr; 7938865f1eaSKarl Rupp break; 7948865f1eaSKarl Rupp } 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*)); 8018865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 80247c6ae99SBarry Smith *iptr = (void*)ptr; 8038865f1eaSKarl Rupp break; 8048865f1eaSKarl Rupp } 80547c6ae99SBarry Smith case 3: { 80647c6ae99SBarry Smith void ***ptr,**bptr; 80747c6ae99SBarry Smith 80847c6ae99SBarry Smith ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 80947c6ae99SBarry Smith 81047c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 81147c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 8128865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 81347c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 81447c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 81547c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 81647c6ae99SBarry Smith } 81747c6ae99SBarry Smith } 81847c6ae99SBarry Smith *iptr = (void*)ptr; 8198865f1eaSKarl Rupp break; 8208865f1eaSKarl Rupp } 82147c6ae99SBarry Smith default: 822ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim); 82347c6ae99SBarry Smith } 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith done: 82647c6ae99SBarry Smith /* add arrays to the checked out list */ 82747c6ae99SBarry Smith if (ghosted) { 828aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 82947c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 83047c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 83147c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 83247c6ae99SBarry Smith break; 83347c6ae99SBarry Smith } 83447c6ae99SBarry Smith } 83547c6ae99SBarry Smith } else { 836aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 83747c6ae99SBarry Smith if (!dd->arrayout[i]) { 83847c6ae99SBarry Smith dd->arrayout[i] = *iptr; 83947c6ae99SBarry Smith dd->startout[i] = iarray_start; 84047c6ae99SBarry Smith break; 84147c6ae99SBarry Smith } 84247c6ae99SBarry Smith } 84347c6ae99SBarry Smith } 84447c6ae99SBarry Smith PetscFunctionReturn(0); 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith #undef __FUNCT__ 848aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray" 84947c6ae99SBarry Smith /*@C 850aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 85147c6ae99SBarry Smith 85247c6ae99SBarry Smith Input Parameter: 85347c6ae99SBarry Smith + da - information about my local patch 85447c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 85547c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 85647c6ae99SBarry Smith 85747c6ae99SBarry Smith Level: advanced 85847c6ae99SBarry Smith 859bcaeba4dSBarry Smith .seealso: DMDAGetArray() 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith @*/ 8627087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 86347c6ae99SBarry Smith { 86447c6ae99SBarry Smith PetscInt i; 86547c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 86647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 86747c6ae99SBarry Smith 86847c6ae99SBarry Smith PetscFunctionBegin; 86947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 87047c6ae99SBarry Smith if (ghosted) { 871aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 87247c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 87347c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 8740298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 8750298fd71SBarry Smith dd->startghostedout[i] = NULL; 87647c6ae99SBarry Smith break; 87747c6ae99SBarry Smith } 87847c6ae99SBarry Smith } 879aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 88047c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 88147c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 88247c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 88347c6ae99SBarry Smith break; 88447c6ae99SBarry Smith } 88547c6ae99SBarry Smith } 88647c6ae99SBarry Smith } else { 887aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 88847c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 88947c6ae99SBarry Smith iarray_start = dd->startout[i]; 8900298fd71SBarry Smith dd->arrayout[i] = NULL; 8910298fd71SBarry Smith dd->startout[i] = NULL; 89247c6ae99SBarry Smith break; 89347c6ae99SBarry Smith } 89447c6ae99SBarry Smith } 895aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 89647c6ae99SBarry Smith if (!dd->arrayin[i]) { 89747c6ae99SBarry Smith dd->arrayin[i] = *iptr; 89847c6ae99SBarry Smith dd->startin[i] = iarray_start; 89947c6ae99SBarry Smith break; 90047c6ae99SBarry Smith } 90147c6ae99SBarry Smith } 90247c6ae99SBarry Smith } 90347c6ae99SBarry Smith PetscFunctionReturn(0); 90447c6ae99SBarry Smith } 90547c6ae99SBarry Smith 906