147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7b2fff234SMatthew G. Knepley #include <petscbt.h> 80c312b8eSJed Brown #include <petscsf.h> 92764a2aaSMatthew G. Knepley #include <petscds.h> 109a800dd8SMatthew G. Knepley #include <petscfe.h> 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith /* 13e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions 1447c6ae99SBarry Smith */ 1547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 16c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */ 17c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */ 1847c6ae99SBarry Smith #undef __FUNCT__ 1947c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d" 201e6b0712SBarry Smith static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 2147c6ae99SBarry Smith { 2247c6ae99SBarry Smith PetscErrorCode ierr; 2347c6ae99SBarry Smith PetscInt n,m; 2447c6ae99SBarry Smith Vec vec = (Vec)obj; 2547c6ae99SBarry Smith PetscScalar *array; 2647c6ae99SBarry Smith mxArray *mat; 279a42bb27SBarry Smith DM da; 2847c6ae99SBarry Smith 2947c6ae99SBarry Smith PetscFunctionBegin; 30c688c046SMatthew G Knepley ierr = VecGetDM(vec, &da);CHKERRQ(ierr); 31ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 32aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 3347c6ae99SBarry Smith 3447c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 3547c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3647c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3747c6ae99SBarry Smith #else 3847c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3947c6ae99SBarry Smith #endif 4047c6ae99SBarry Smith ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 4147c6ae99SBarry Smith ierr = PetscObjectName(obj);CHKERRQ(ierr); 4247c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 4347c6ae99SBarry Smith 4447c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 4547c6ae99SBarry Smith PetscFunctionReturn(0); 4647c6ae99SBarry Smith } 4747c6ae99SBarry Smith #endif 4847c6ae99SBarry Smith 4947c6ae99SBarry Smith 5047c6ae99SBarry Smith #undef __FUNCT__ 51564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA" 527087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 5347c6ae99SBarry Smith { 5447c6ae99SBarry Smith PetscErrorCode ierr; 5547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5647c6ae99SBarry Smith 5747c6ae99SBarry Smith PetscFunctionBegin; 5847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5947c6ae99SBarry Smith PetscValidPointer(g,2); 6011689aebSJed Brown if (da->defaultSection) { 6111689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 6211689aebSJed Brown } else { 6347c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6447c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6547c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 66401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 67c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6847c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 6947c6ae99SBarry Smith if (dd->w == 1 && dd->dim == 2) { 70bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 7147c6ae99SBarry Smith } 7247c6ae99SBarry Smith #endif 7311689aebSJed Brown } 7447c6ae99SBarry Smith PetscFunctionReturn(0); 7547c6ae99SBarry Smith } 7647c6ae99SBarry Smith 77a66d4d66SMatthew G Knepley #undef __FUNCT__ 7857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells" 79939f6067SMatthew G. Knepley /*@ 80939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 81939f6067SMatthew G. Knepley 82939f6067SMatthew G. Knepley Input Parameter: 83939f6067SMatthew G. Knepley . dm - The DM object 84939f6067SMatthew G. Knepley 85939f6067SMatthew G. Knepley Output Parameters: 86939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 87939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 88939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 89939f6067SMatthew G. Knepley - numCells - The number of local cells 90939f6067SMatthew G. Knepley 91939f6067SMatthew G. Knepley Level: developer 92939f6067SMatthew G. Knepley 93939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint() 94939f6067SMatthew G. Knepley @*/ 953582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 9657459e9aSMatthew G Knepley { 9757459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 98*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 9957459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 10057459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 10157459e9aSMatthew G Knepley 10257459e9aSMatthew G Knepley PetscFunctionBegin; 103d6401b93SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1043582350cSMatthew G. Knepley if (numCellsX) { 1053582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 1063582350cSMatthew G. Knepley *numCellsX = mx; 1073582350cSMatthew G. Knepley } 1083582350cSMatthew G. Knepley if (numCellsY) { 1093582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,3); 1103582350cSMatthew G. Knepley *numCellsY = my; 1113582350cSMatthew G. Knepley } 1123582350cSMatthew G. Knepley if (numCellsZ) { 1133582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,4); 1143582350cSMatthew G. Knepley *numCellsZ = mz; 1153582350cSMatthew G. Knepley } 11657459e9aSMatthew G Knepley if (numCells) { 1173582350cSMatthew G. Knepley PetscValidIntPointer(numCells,5); 11857459e9aSMatthew G Knepley *numCells = nC; 11957459e9aSMatthew G Knepley } 12057459e9aSMatthew G Knepley PetscFunctionReturn(0); 12157459e9aSMatthew G Knepley } 12257459e9aSMatthew G Knepley 12357459e9aSMatthew G Knepley #undef __FUNCT__ 124d6401b93SMatthew G. Knepley #define __FUNCT__ "DMDAGetCellPoint" 125d6401b93SMatthew G. Knepley /*@ 126d6401b93SMatthew G. Knepley DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 127d6401b93SMatthew G. Knepley 128d6401b93SMatthew G. Knepley Input Parameters: 129d6401b93SMatthew G. Knepley + dm - The DM object 130d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell 131d6401b93SMatthew G. Knepley 132d6401b93SMatthew G. Knepley Output Parameters: 133d6401b93SMatthew G. Knepley . point - The local DM point 134d6401b93SMatthew G. Knepley 135d6401b93SMatthew G. Knepley Level: developer 136d6401b93SMatthew G. Knepley 137d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells() 138d6401b93SMatthew G. Knepley @*/ 139d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) 140d6401b93SMatthew G. Knepley { 141*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 142d6401b93SMatthew G. Knepley DMDALocalInfo info; 143d6401b93SMatthew G. Knepley PetscErrorCode ierr; 144d6401b93SMatthew G. Knepley 145d6401b93SMatthew G. Knepley PetscFunctionBegin; 146d6401b93SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 147d6401b93SMatthew G. Knepley PetscValidIntPointer(point,5); 148d6401b93SMatthew G. Knepley ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 149d6401b93SMatthew G. Knepley if (dim > 0) {if ((i < info.gxs) || (i >= info.gxs+info.gxm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %d not in [%d, %d)", i, info.gxs, info.gxs+info.gxm);} 150d6401b93SMatthew G. Knepley if (dim > 1) {if ((i < info.gys) || (i >= info.gys+info.gym)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %d not in [%d, %d)", i, info.gys, info.gys+info.gym);} 151d6401b93SMatthew G. Knepley if (dim > 2) {if ((i < info.gzs) || (i >= info.gzs+info.gzm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %d not in [%d, %d)", i, info.gzs, info.gzs+info.gzm);} 152d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0); 153d6401b93SMatthew G. Knepley PetscFunctionReturn(0); 154d6401b93SMatthew G. Knepley } 155d6401b93SMatthew G. Knepley 156d6401b93SMatthew G. Knepley #undef __FUNCT__ 15757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices" 15857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 15957459e9aSMatthew G Knepley { 16057459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 161*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 16257459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 16357459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 16457459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 16557459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 16657459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 16757459e9aSMatthew G Knepley 16857459e9aSMatthew G Knepley PetscFunctionBegin; 16957459e9aSMatthew G Knepley if (numVerticesX) { 17057459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 17157459e9aSMatthew G Knepley *numVerticesX = nVx; 17257459e9aSMatthew G Knepley } 17357459e9aSMatthew G Knepley if (numVerticesY) { 17457459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 17557459e9aSMatthew G Knepley *numVerticesY = nVy; 17657459e9aSMatthew G Knepley } 17757459e9aSMatthew G Knepley if (numVerticesZ) { 17857459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 17957459e9aSMatthew G Knepley *numVerticesZ = nVz; 18057459e9aSMatthew G Knepley } 18157459e9aSMatthew G Knepley if (numVertices) { 18257459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 18357459e9aSMatthew G Knepley *numVertices = nV; 18457459e9aSMatthew G Knepley } 18557459e9aSMatthew G Knepley PetscFunctionReturn(0); 18657459e9aSMatthew G Knepley } 18757459e9aSMatthew G Knepley 18857459e9aSMatthew G Knepley #undef __FUNCT__ 18957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces" 19057459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 19157459e9aSMatthew G Knepley { 19257459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 193*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 19457459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 19557459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 19657459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 19757459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 19857459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 19957459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 20057459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 20157459e9aSMatthew G Knepley 20257459e9aSMatthew G Knepley PetscFunctionBegin; 20357459e9aSMatthew G Knepley if (numXFacesX) { 20457459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 20557459e9aSMatthew G Knepley *numXFacesX = nxF; 20657459e9aSMatthew G Knepley } 20757459e9aSMatthew G Knepley if (numXFaces) { 20857459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 20957459e9aSMatthew G Knepley *numXFaces = nXF; 21057459e9aSMatthew G Knepley } 21157459e9aSMatthew G Knepley if (numYFacesY) { 21257459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 21357459e9aSMatthew G Knepley *numYFacesY = nyF; 21457459e9aSMatthew G Knepley } 21557459e9aSMatthew G Knepley if (numYFaces) { 21657459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 21757459e9aSMatthew G Knepley *numYFaces = nYF; 21857459e9aSMatthew G Knepley } 21957459e9aSMatthew G Knepley if (numZFacesZ) { 22057459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 22157459e9aSMatthew G Knepley *numZFacesZ = nzF; 22257459e9aSMatthew G Knepley } 22357459e9aSMatthew G Knepley if (numZFaces) { 22457459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 22557459e9aSMatthew G Knepley *numZFaces = nZF; 22657459e9aSMatthew G Knepley } 22757459e9aSMatthew G Knepley PetscFunctionReturn(0); 22857459e9aSMatthew G Knepley } 22957459e9aSMatthew G Knepley 23057459e9aSMatthew G Knepley #undef __FUNCT__ 23157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum" 23257459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 23357459e9aSMatthew G Knepley { 234*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 23557459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 23657459e9aSMatthew G Knepley PetscErrorCode ierr; 23757459e9aSMatthew G Knepley 23857459e9aSMatthew G Knepley PetscFunctionBegin; 2398865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 2408865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 2413582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2420298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2430298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 24457459e9aSMatthew G Knepley if (height == 0) { 24557459e9aSMatthew G Knepley /* Cells */ 2468865f1eaSKarl Rupp if (pStart) *pStart = 0; 2478865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 24857459e9aSMatthew G Knepley } else if (height == 1) { 24957459e9aSMatthew G Knepley /* Faces */ 2508865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 2518865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 25257459e9aSMatthew G Knepley } else if (height == dim) { 25357459e9aSMatthew G Knepley /* Vertices */ 2548865f1eaSKarl Rupp if (pStart) *pStart = nC; 2558865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 25657459e9aSMatthew G Knepley } else if (height < 0) { 25757459e9aSMatthew G Knepley /* All points */ 2588865f1eaSKarl Rupp if (pStart) *pStart = 0; 2598865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 26082f516ccSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 26157459e9aSMatthew G Knepley PetscFunctionReturn(0); 26257459e9aSMatthew G Knepley } 26357459e9aSMatthew G Knepley 26457459e9aSMatthew G Knepley #undef __FUNCT__ 2653385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum" 2663385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) 2673385ff7eSMatthew G. Knepley { 268*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2693385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2703385ff7eSMatthew G. Knepley PetscErrorCode ierr; 2713385ff7eSMatthew G. Knepley 2723385ff7eSMatthew G. Knepley PetscFunctionBegin; 2733385ff7eSMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart,3); 2743385ff7eSMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd,4); 2753385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2763385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2773385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 2783385ff7eSMatthew G. Knepley if (depth == dim) { 2793385ff7eSMatthew G. Knepley /* Cells */ 2803385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2813385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC; 2823385ff7eSMatthew G. Knepley } else if (depth == dim-1) { 2833385ff7eSMatthew G. Knepley /* Faces */ 2843385ff7eSMatthew G. Knepley if (pStart) *pStart = nC+nV; 2853385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2863385ff7eSMatthew G. Knepley } else if (depth == 0) { 2873385ff7eSMatthew G. Knepley /* Vertices */ 2883385ff7eSMatthew G. Knepley if (pStart) *pStart = nC; 2893385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV; 2903385ff7eSMatthew G. Knepley } else if (depth < 0) { 2913385ff7eSMatthew G. Knepley /* All points */ 2923385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2933385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2943385ff7eSMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth); 2953385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2963385ff7eSMatthew G. Knepley } 2973385ff7eSMatthew G. Knepley 2983385ff7eSMatthew G. Knepley #undef __FUNCT__ 2993385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize" 3003385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) 3013385ff7eSMatthew G. Knepley { 302*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 3033385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 3043385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3053385ff7eSMatthew G. Knepley 3063385ff7eSMatthew G. Knepley PetscFunctionBegin; 3073385ff7eSMatthew G. Knepley *coneSize = 0; 3083385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 3093385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 3103385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 3113385ff7eSMatthew G. Knepley switch (dim) { 3123385ff7eSMatthew G. Knepley case 2: 3133385ff7eSMatthew G. Knepley if (p >= 0) { 3143385ff7eSMatthew G. Knepley if (p < nC) { 3153385ff7eSMatthew G. Knepley *coneSize = 4; 3163385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 3173385ff7eSMatthew G. Knepley *coneSize = 0; 3183385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF+nZF) { 3193385ff7eSMatthew G. Knepley *coneSize = 2; 3203385ff7eSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 3213385ff7eSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3223385ff7eSMatthew G. Knepley break; 3233385ff7eSMatthew G. Knepley case 3: 3243385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3253385ff7eSMatthew G. Knepley break; 3263385ff7eSMatthew G. Knepley } 3273385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3283385ff7eSMatthew G. Knepley } 3293385ff7eSMatthew G. Knepley 3303385ff7eSMatthew G. Knepley #undef __FUNCT__ 3313385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetCone" 3323385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) 3333385ff7eSMatthew G. Knepley { 334*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 3353385ff7eSMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 3363385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3373385ff7eSMatthew G. Knepley 3383385ff7eSMatthew G. Knepley PetscFunctionBegin; 3393385ff7eSMatthew G. Knepley if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);} 3403385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr); 3413385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 3423385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 3433385ff7eSMatthew G. Knepley switch (dim) { 3443385ff7eSMatthew G. Knepley case 2: 3453385ff7eSMatthew G. Knepley if (p >= 0) { 3463385ff7eSMatthew G. Knepley if (p < nC) { 3473385ff7eSMatthew G. Knepley const PetscInt cy = p / nCx; 3483385ff7eSMatthew G. Knepley const PetscInt cx = p % nCx; 3493385ff7eSMatthew G. Knepley 3503385ff7eSMatthew G. Knepley (*cone)[0] = cy*nxF + cx + nC+nV; 3513385ff7eSMatthew G. Knepley (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF; 3523385ff7eSMatthew G. Knepley (*cone)[2] = cy*nxF + cx + nxF + nC+nV; 3533385ff7eSMatthew G. Knepley (*cone)[3] = cx*nyF + cy + nC+nV+nXF; 3543385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 3553385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 3563385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF) { 3573385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV) / nxF; 3583385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV) % nxF; 3593385ff7eSMatthew G. Knepley 3603385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3613385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + 1 + nC; 3623385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF) { 3633385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV+nXF) / nyF; 3643385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV+nXF) % nyF; 3653385ff7eSMatthew G. Knepley 3663385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3673385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + nVx + nC; 3683385ff7eSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 3693385ff7eSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3703385ff7eSMatthew G. Knepley break; 3713385ff7eSMatthew G. Knepley case 3: 3723385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3733385ff7eSMatthew G. Knepley break; 3743385ff7eSMatthew G. Knepley } 3753385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3763385ff7eSMatthew G. Knepley } 3773385ff7eSMatthew G. Knepley 3783385ff7eSMatthew G. Knepley #undef __FUNCT__ 3793385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone" 3803385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) 3813385ff7eSMatthew G. Knepley { 3823385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3833385ff7eSMatthew G. Knepley 3843385ff7eSMatthew G. Knepley PetscFunctionBegin; 3853385ff7eSMatthew G. Knepley ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr); 3863385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3873385ff7eSMatthew G. Knepley } 3883385ff7eSMatthew G. Knepley 3893385ff7eSMatthew G. Knepley #undef __FUNCT__ 390a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection" 391a66d4d66SMatthew G Knepley /*@C 392a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 393a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 394a66d4d66SMatthew G Knepley 395a66d4d66SMatthew G Knepley Input Parameters: 396a66d4d66SMatthew G Knepley + dm- The DMDA 397a66d4d66SMatthew G Knepley . numFields - The number of fields 398affa55c7SMatthew G. Knepley . numComp - The number of components in each field 399affa55c7SMatthew G. Knepley . numDof - The number of dofs per dimension for each field 4000298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL 401a66d4d66SMatthew G Knepley 402a66d4d66SMatthew G Knepley Level: developer 403a66d4d66SMatthew G Knepley 404a66d4d66SMatthew G Knepley Note: 405a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 406a66d4d66SMatthew G Knepley 407a66d4d66SMatthew G Knepley - Cells: [0, nC) 408a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 40988ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 41088ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 41188ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 412a66d4d66SMatthew G Knepley 413a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 414a66d4d66SMatthew G Knepley @*/ 41544e1f9abSMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, const PetscInt numComp[], const PetscInt numDof[], const PetscInt numFaceDof[], PetscSection *s) 416a66d4d66SMatthew G Knepley { 417a4b60ecfSMatthew G. Knepley PetscSection section; 418*c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 41980800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 420b2fff234SMatthew G. Knepley PetscBT isLeaf; 42188ed4aceSMatthew G Knepley PetscSF sf; 422feafbc5cSMatthew G Knepley PetscMPIInt rank; 42388ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 42488ed4aceSMatthew G Knepley PetscInt *localPoints; 42588ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 426f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 42757459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 42857459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 42988ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 430a66d4d66SMatthew G Knepley PetscErrorCode ierr; 431a66d4d66SMatthew G Knepley 432a66d4d66SMatthew G Knepley PetscFunctionBegin; 433a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4343582350cSMatthew G. Knepley PetscValidPointer(s, 4); 43582f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 4363582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 43757459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 43857459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 43957459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 44057459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 44157459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 44257459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 44357459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 44457459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 44557459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 44682f516ccSBarry Smith if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 44788ed4aceSMatthew G Knepley /* Create local section */ 44880800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 449a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 450affa55c7SMatthew G. Knepley numVertexTotDof += numDof[f*(dim+1)+0]; 451affa55c7SMatthew G. Knepley numCellTotDof += numDof[f*(dim+1)+dim]; 452affa55c7SMatthew G. Knepley numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0; 453affa55c7SMatthew G. Knepley numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0; 454affa55c7SMatthew G. Knepley numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0; 455a66d4d66SMatthew G Knepley } 456a4b60ecfSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);CHKERRQ(ierr); 457affa55c7SMatthew G. Knepley if (numFields > 0) { 458a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr); 459a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 46080800b1aSMatthew G Knepley const char *name; 46180800b1aSMatthew G Knepley 46280800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 463affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr); 46480800b1aSMatthew G Knepley if (numComp) { 465a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr); 466a66d4d66SMatthew G Knepley } 467a66d4d66SMatthew G Knepley } 468a66d4d66SMatthew G Knepley } 469a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 470a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 471a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 472affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr); 473a66d4d66SMatthew G Knepley } 474a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr); 475a66d4d66SMatthew G Knepley } 476a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 477a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 478affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 479a66d4d66SMatthew G Knepley } 480a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr); 481a66d4d66SMatthew G Knepley } 482a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 483a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 484affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 485a66d4d66SMatthew G Knepley } 486a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr); 487a66d4d66SMatthew G Knepley } 488a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 489a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 490affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 491a66d4d66SMatthew G Knepley } 492a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr); 493a66d4d66SMatthew G Knepley } 494a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 495a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 496affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr); 497a66d4d66SMatthew G Knepley } 498a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr); 499a66d4d66SMatthew G Knepley } 500a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 50188ed4aceSMatthew G Knepley /* Create mesh point SF */ 502b2fff234SMatthew G. Knepley ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr); 50388ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 50488ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 50588ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 50688ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 5077128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 50888ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 509b2fff234SMatthew G. Knepley PetscInt xv, yv, zv; 51088ed4aceSMatthew G Knepley 5113814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 512b2fff234SMatthew G. Knepley if (xp < 0) { /* left */ 513b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 514b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 515b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 516b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 517b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 518b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 519b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 520b2fff234SMatthew G. Knepley } else { 521b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 522b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 523b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 524b2fff234SMatthew G. Knepley } 525b2fff234SMatthew G. Knepley } 526b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 527b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 528b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 529b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 530b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 531b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 532b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 533b2fff234SMatthew G. Knepley } else { 534b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 535b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 536b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 537b2fff234SMatthew G. Knepley } 538b2fff234SMatthew G. Knepley } 539b2fff234SMatthew G. Knepley } else { 540b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 541b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 542b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 543b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 544b2fff234SMatthew G. Knepley } 545b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 546b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 547b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 548b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 549b2fff234SMatthew G. Knepley } 550b2fff234SMatthew G. Knepley } else { 551b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 552b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 553b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 554b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 555b2fff234SMatthew G. Knepley } 556b2fff234SMatthew G. Knepley } 557b2fff234SMatthew G. Knepley #if 0 558b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 559b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 560b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 561b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 562b2fff234SMatthew G. Knepley } 563b2fff234SMatthew G. Knepley #endif 564b2fff234SMatthew G. Knepley } 565b2fff234SMatthew G. Knepley } 566b2fff234SMatthew G. Knepley } else if (xp > 0) { /* right */ 567b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 568b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 569b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 570b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 571b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 572b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 573b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 574b2fff234SMatthew G. Knepley } else { 575b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 576b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 577b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 578b2fff234SMatthew G. Knepley } 579b2fff234SMatthew G. Knepley } 580b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 581b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 582b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 583b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 584b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 585b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 586b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 587b2fff234SMatthew G. Knepley } else { 588b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 589b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 590b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 591b2fff234SMatthew G. Knepley } 592b2fff234SMatthew G. Knepley } 593b2fff234SMatthew G. Knepley } else { 594b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 595b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 596b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 597b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 598b2fff234SMatthew G. Knepley } 599b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 600b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 601b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 602b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 603b2fff234SMatthew G. Knepley } 604b2fff234SMatthew G. Knepley } else { 605b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 606b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 607b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 608b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 609b2fff234SMatthew G. Knepley } 610b2fff234SMatthew G. Knepley } 611b2fff234SMatthew G. Knepley #if 0 612b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 613b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 614b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 615b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 616b2fff234SMatthew G. Knepley } 617b2fff234SMatthew G. Knepley #endif 618b2fff234SMatthew G. Knepley } 619b2fff234SMatthew G. Knepley } 620b2fff234SMatthew G. Knepley } else { 621b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 622b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 623b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 624b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 625b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 626b2fff234SMatthew G. Knepley } 627b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 628b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 629b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 630b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 631b2fff234SMatthew G. Knepley } 632b2fff234SMatthew G. Knepley } else { 633b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 634b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 635b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 636b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 637b2fff234SMatthew G. Knepley } 638b2fff234SMatthew G. Knepley } 639b2fff234SMatthew G. Knepley #if 0 640b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 641b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 642b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 643b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 644b2fff234SMatthew G. Knepley } 645b2fff234SMatthew G. Knepley #endif 646b2fff234SMatthew G. Knepley } 647b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 648b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 649b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 650b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 651b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 652b2fff234SMatthew G. Knepley } 653b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 654b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 655b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 656b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 657b2fff234SMatthew G. Knepley } 658b2fff234SMatthew G. Knepley } else { 659b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 660b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 661b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 662b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 663b2fff234SMatthew G. Knepley } 664b2fff234SMatthew G. Knepley } 665b2fff234SMatthew G. Knepley #if 0 666b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 667b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 668b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 669b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 670b2fff234SMatthew G. Knepley } 671b2fff234SMatthew G. Knepley #endif 672b2fff234SMatthew G. Knepley } 673b2fff234SMatthew G. Knepley } else { 674b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 675b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 676b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 677b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 678b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 679b2fff234SMatthew G. Knepley } 680b2fff234SMatthew G. Knepley } 681b2fff234SMatthew G. Knepley #if 0 682b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 683b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 684b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 685b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 686b2fff234SMatthew G. Knepley } 687b2fff234SMatthew G. Knepley #endif 688b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 689b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 690b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 691b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 692b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 693b2fff234SMatthew G. Knepley } 694b2fff234SMatthew G. Knepley } 695b2fff234SMatthew G. Knepley #if 0 696b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 697b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 698b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 699b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 700b2fff234SMatthew G. Knepley } 701b2fff234SMatthew G. Knepley #endif 702b2fff234SMatthew G. Knepley } else { 703b2fff234SMatthew G. Knepley /* Nothing is shared from the interior */ 70488ed4aceSMatthew G Knepley } 70588ed4aceSMatthew G Knepley } 70688ed4aceSMatthew G Knepley } 70788ed4aceSMatthew G Knepley } 70888ed4aceSMatthew G Knepley } 709b2fff234SMatthew G. Knepley } 710b2fff234SMatthew G. Knepley } 711b2fff234SMatthew G. Knepley ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr); 712dcca6d9dSJed Brown ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr); 71388ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 71488ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 71588ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 7167128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 71788ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 718f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 71988ed4aceSMatthew G Knepley 7203814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 72188ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 72288ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 72388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 724b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 725f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 726f5eeb527SMatthew G Knepley 727b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 728f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 729f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 730f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 731f5eeb527SMatthew G Knepley ++nL; 732b2fff234SMatthew G. Knepley } 73388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 734b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 735f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 736f5eeb527SMatthew G Knepley 737b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 738f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 739f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 740f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 741f5eeb527SMatthew G Knepley ++nL; 742b2fff234SMatthew G. Knepley } 74388ed4aceSMatthew G Knepley } else { 744b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 745b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 746f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 747f5eeb527SMatthew G Knepley 748b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 749f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 750f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 751f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 752b2fff234SMatthew G. Knepley ++nL; 753b2fff234SMatthew G. Knepley } 754f5eeb527SMatthew G Knepley } 75588ed4aceSMatthew G Knepley } 75688ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 75788ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 758b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 759f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 760f5eeb527SMatthew G Knepley 761b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 762f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 763f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 764f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 765f5eeb527SMatthew G Knepley ++nL; 766b2fff234SMatthew G. Knepley } 76788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 768b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 769f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 770f5eeb527SMatthew G Knepley 771b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 772f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 773f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 774f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 775f5eeb527SMatthew G Knepley ++nL; 776b2fff234SMatthew G. Knepley } 77788ed4aceSMatthew G Knepley } else { 778b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 779b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 780f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 781f5eeb527SMatthew G Knepley 782b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 783f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 784f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 785f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 786b2fff234SMatthew G. Knepley ++nL; 787b2fff234SMatthew G. Knepley } 788f5eeb527SMatthew G Knepley } 78988ed4aceSMatthew G Knepley } 79088ed4aceSMatthew G Knepley } else { 79188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 792b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 793b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 794f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 795f5eeb527SMatthew G Knepley 796b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 797f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 798f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 799f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 800b2fff234SMatthew G. Knepley ++nL; 801b2fff234SMatthew G. Knepley } 802f5eeb527SMatthew G Knepley } 80388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 804b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 805b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 806f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 807f5eeb527SMatthew G Knepley 808b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 809f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 810f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 811f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 812b2fff234SMatthew G. Knepley ++nL; 813b2fff234SMatthew G. Knepley } 814f5eeb527SMatthew G Knepley } 81588ed4aceSMatthew G Knepley } else { 816f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 817b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 818b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 819f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 820f5eeb527SMatthew G Knepley 821b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 822f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 823f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 824f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 825b2fff234SMatthew G. Knepley ++nL; 826f5eeb527SMatthew G Knepley } 827f5eeb527SMatthew G Knepley } 828b2fff234SMatthew G. Knepley } 829b2fff234SMatthew G. Knepley #if 0 830b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 831f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 832b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 833f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 834f5eeb527SMatthew G Knepley 835b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 836f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 837f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 838f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 839f5eeb527SMatthew G Knepley } 84088ed4aceSMatthew G Knepley } 841b2fff234SMatthew G. Knepley #endif 842b2fff234SMatthew G. Knepley } 84388ed4aceSMatthew G Knepley } 84488ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 84588ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 84688ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 847b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 848f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 849f5eeb527SMatthew G Knepley 850b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 851f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 852f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 853f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 854f5eeb527SMatthew G Knepley ++nL; 855b2fff234SMatthew G. Knepley } 85688ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 857b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 858f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 859f5eeb527SMatthew G Knepley 860b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 861f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 862f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 863f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 864f5eeb527SMatthew G Knepley ++nL; 865b2fff234SMatthew G. Knepley } 86688ed4aceSMatthew G Knepley } else { 867b2fff234SMatthew G. Knepley nleavesCheck += nVz; 868b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 869b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 870f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 871f5eeb527SMatthew G Knepley 872b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 873f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 874f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 875f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 876b2fff234SMatthew G. Knepley ++nL; 877b2fff234SMatthew G. Knepley } 878f5eeb527SMatthew G Knepley } 87988ed4aceSMatthew G Knepley } 88088ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 88188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 882b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 883f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 884f5eeb527SMatthew G Knepley 885b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 886f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 887f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 888f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 889f5eeb527SMatthew G Knepley ++nL; 890b2fff234SMatthew G. Knepley } 89188ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 892b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 893f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 894f5eeb527SMatthew G Knepley 895b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 896f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 897f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 898f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 899f5eeb527SMatthew G Knepley ++nL; 900b2fff234SMatthew G. Knepley } 90188ed4aceSMatthew G Knepley } else { 902b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 903b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 904f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 905f5eeb527SMatthew G Knepley 906b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 907f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 908f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 909f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 910b2fff234SMatthew G. Knepley ++nL; 911b2fff234SMatthew G. Knepley } 912f5eeb527SMatthew G Knepley } 91388ed4aceSMatthew G Knepley } 91488ed4aceSMatthew G Knepley } else { 91588ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 916b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 917b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 918f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 919f5eeb527SMatthew G Knepley 920b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 921f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 922f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 923f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 924b2fff234SMatthew G. Knepley ++nL; 925b2fff234SMatthew G. Knepley } 926f5eeb527SMatthew G Knepley } 92788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 928b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 929b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 930f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 931f5eeb527SMatthew G Knepley 932b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 933f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 934f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 935f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 936b2fff234SMatthew G. Knepley ++nL; 937b2fff234SMatthew G. Knepley } 938f5eeb527SMatthew G Knepley } 93988ed4aceSMatthew G Knepley } else { 940f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 941b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 942b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 943f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 944f5eeb527SMatthew G Knepley 945b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 946f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 947f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 948f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 949b2fff234SMatthew G. Knepley ++nL; 950f5eeb527SMatthew G Knepley } 951f5eeb527SMatthew G Knepley } 952b2fff234SMatthew G. Knepley } 953b2fff234SMatthew G. Knepley #if 0 954b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 955f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 956b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 957f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 958f5eeb527SMatthew G Knepley 959b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 960f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 961f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 962f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 963b2fff234SMatthew G. Knepley ++nL; 964f5eeb527SMatthew G Knepley } 96588ed4aceSMatthew G Knepley } 966b2fff234SMatthew G. Knepley #endif 967b2fff234SMatthew G. Knepley } 96888ed4aceSMatthew G Knepley } 96988ed4aceSMatthew G Knepley } else { 97088ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 97188ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 972b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 973b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 974f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 975f5eeb527SMatthew G Knepley 976b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 977f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 978f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 979f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 980b2fff234SMatthew G. Knepley ++nL; 981b2fff234SMatthew G. Knepley } 982f5eeb527SMatthew G Knepley } 98388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 984b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 985b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 986f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 987f5eeb527SMatthew G Knepley 988b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 989f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 990f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 991f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 992b2fff234SMatthew G. Knepley ++nL; 993b2fff234SMatthew G. Knepley } 994f5eeb527SMatthew G Knepley } 99588ed4aceSMatthew G Knepley } else { 996f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 997b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 998b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 999f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1000f5eeb527SMatthew G Knepley 1001b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1002f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1003f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1004f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1005b2fff234SMatthew G. Knepley ++nL; 1006f5eeb527SMatthew G Knepley } 1007f5eeb527SMatthew G Knepley } 1008b2fff234SMatthew G. Knepley } 1009b2fff234SMatthew G. Knepley #if 0 1010b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 1011f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1012b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 1013f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1014f5eeb527SMatthew G Knepley 1015b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1016f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1017f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1018f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1019b2fff234SMatthew G. Knepley ++nL; 1020f5eeb527SMatthew G Knepley } 102188ed4aceSMatthew G Knepley } 1022b2fff234SMatthew G. Knepley #endif 1023b2fff234SMatthew G. Knepley } 102488ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 102588ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 1026b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1027b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 1028f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1029f5eeb527SMatthew G Knepley 1030b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1031f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1032f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1033f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1034b2fff234SMatthew G. Knepley ++nL; 1035b2fff234SMatthew G. Knepley } 1036f5eeb527SMatthew G Knepley } 103788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1038b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1039b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 1040f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1041f5eeb527SMatthew G Knepley 1042b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1043f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1044f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1045f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1046b2fff234SMatthew G. Knepley ++nL; 1047b2fff234SMatthew G. Knepley } 1048f5eeb527SMatthew G Knepley } 104988ed4aceSMatthew G Knepley } else { 1050f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 1051b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1052b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 1053f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1054f5eeb527SMatthew G Knepley 1055b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1056f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1057f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1058f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1059b2fff234SMatthew G. Knepley ++nL; 1060f5eeb527SMatthew G Knepley } 1061f5eeb527SMatthew G Knepley } 1062b2fff234SMatthew G. Knepley } 1063b2fff234SMatthew G. Knepley #if 0 1064b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 1065f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1066b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 1067f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1068f5eeb527SMatthew G Knepley 1069b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1070f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1071f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1072f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1073b2fff234SMatthew G. Knepley ++nL; 1074f5eeb527SMatthew G Knepley } 107588ed4aceSMatthew G Knepley } 1076b2fff234SMatthew G. Knepley #endif 1077b2fff234SMatthew G. Knepley } 107888ed4aceSMatthew G Knepley } else { 107988ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 1080f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1081b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1082b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 1083f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1084f5eeb527SMatthew G Knepley 1085b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1086f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1087f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1088f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1089b2fff234SMatthew G. Knepley ++nL; 1090f5eeb527SMatthew G Knepley } 1091f5eeb527SMatthew G Knepley } 1092b2fff234SMatthew G. Knepley } 1093b2fff234SMatthew G. Knepley #if 0 1094b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1095f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1096b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 1097f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1098f5eeb527SMatthew G Knepley 1099b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1100f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1101f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1102f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1103b2fff234SMatthew G. Knepley ++nL; 1104f5eeb527SMatthew G Knepley } 1105b2fff234SMatthew G. Knepley } 1106b2fff234SMatthew G. Knepley #endif 110788ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1108f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1109b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1110b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 1111f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1112f5eeb527SMatthew G Knepley 1113b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1114f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1115f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1116f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1117b2fff234SMatthew G. Knepley ++nL; 1118f5eeb527SMatthew G Knepley } 1119f5eeb527SMatthew G Knepley } 1120b2fff234SMatthew G. Knepley } 1121b2fff234SMatthew G. Knepley #if 0 1122b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1123f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1124b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 1125f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1126f5eeb527SMatthew G Knepley 1127b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1128f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1129f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1130f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1131b2fff234SMatthew G. Knepley ++nL; 1132f5eeb527SMatthew G Knepley } 1133b2fff234SMatthew G. Knepley } 1134b2fff234SMatthew G. Knepley #endif 113588ed4aceSMatthew G Knepley } else { 113688ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 113788ed4aceSMatthew G Knepley } 113888ed4aceSMatthew G Knepley } 113988ed4aceSMatthew G Knepley } 114088ed4aceSMatthew G Knepley } 114188ed4aceSMatthew G Knepley } 114288ed4aceSMatthew G Knepley } 114388ed4aceSMatthew G Knepley } 1144b2fff234SMatthew G. Knepley ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr); 1145b2fff234SMatthew G. Knepley /* Remove duplication in leaf determination */ 1146b2fff234SMatthew G. Knepley if (nleaves != nL) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck); 114782f516ccSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr); 114888ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1149a4b60ecfSMatthew G. Knepley ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 115088ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1151affa55c7SMatthew G. Knepley *s = section; 1152affa55c7SMatthew G. Knepley PetscFunctionReturn(0); 1153affa55c7SMatthew G. Knepley } 1154affa55c7SMatthew G. Knepley 11553385ff7eSMatthew G. Knepley #undef __FUNCT__ 11563385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates" 11573385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) 11583385ff7eSMatthew G. Knepley { 11593385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 11603385ff7eSMatthew G. Knepley Vec coordinates; 11613385ff7eSMatthew G. Knepley PetscSection section; 11623385ff7eSMatthew G. Knepley PetscScalar *coords; 11633385ff7eSMatthew G. Knepley PetscReal h[3]; 11643385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 11653385ff7eSMatthew G. Knepley PetscErrorCode ierr; 11663385ff7eSMatthew G. Knepley 11673385ff7eSMatthew G. Knepley PetscFunctionBegin; 11683385ff7eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11693385ff7eSMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 11703385ff7eSMatthew G. Knepley h[0] = (xu - xl)/M; 11713385ff7eSMatthew G. Knepley h[1] = (yu - yl)/N; 11723385ff7eSMatthew G. Knepley h[2] = (zu - zl)/P; 11733385ff7eSMatthew G. Knepley ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 11743385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 11753385ff7eSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 11763385ff7eSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr); 11773385ff7eSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr); 11783385ff7eSMatthew G. Knepley ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr); 11793385ff7eSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 11803385ff7eSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr); 11813385ff7eSMatthew G. Knepley } 11823385ff7eSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 11833385ff7eSMatthew G. Knepley ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 11843385ff7eSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr); 11853385ff7eSMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 11863385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) { 1187e4d69e08SBarry Smith PetscInt ind[3], d, off; 11883385ff7eSMatthew G. Knepley 1189e4d69e08SBarry Smith ind[0] = 0; 1190e4d69e08SBarry Smith ind[1] = 0; 1191e4d69e08SBarry Smith ind[2] = k + da->zs; 11923385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) { 11933385ff7eSMatthew G. Knepley ind[1] = j + da->ys; 11943385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) { 11953385ff7eSMatthew G. Knepley const PetscInt vertex = (k*nVy + j)*nVx + i + vStart; 11963385ff7eSMatthew G. Knepley 11973385ff7eSMatthew G. Knepley ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr); 11983385ff7eSMatthew G. Knepley ind[0] = i + da->xs; 11993385ff7eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 12003385ff7eSMatthew G. Knepley coords[off+d] = h[d]*ind[d]; 12013385ff7eSMatthew G. Knepley } 12023385ff7eSMatthew G. Knepley } 12033385ff7eSMatthew G. Knepley } 12043385ff7eSMatthew G. Knepley } 12053385ff7eSMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 12063385ff7eSMatthew G. Knepley ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr); 12073385ff7eSMatthew G. Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1208a4b60ecfSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 12093385ff7eSMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 12103385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 12113385ff7eSMatthew G. Knepley } 12129a800dd8SMatthew G. Knepley 12139a800dd8SMatthew G. Knepley #undef __FUNCT__ 12149a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunctionLocal" 1215b7a4d31fSMatthew G. Knepley PetscErrorCode DMDAProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 12169a800dd8SMatthew G. Knepley { 12172764a2aaSMatthew G. Knepley PetscDS prob; 1218b7a4d31fSMatthew G. Knepley PetscFE fe; 1219b7a4d31fSMatthew G. Knepley PetscDualSpace sp; 12209a800dd8SMatthew G. Knepley PetscQuadrature q; 12219a800dd8SMatthew G. Knepley PetscSection section; 12229a800dd8SMatthew G. Knepley PetscScalar *values; 12239a800dd8SMatthew G. Knepley PetscReal *v0, *J, *detJ; 1224b7a4d31fSMatthew G. Knepley PetscInt numFields, numComp, numPoints, dim, spDim, totDim, numValues, cStart, cEnd, f, c, v, d; 12259a800dd8SMatthew G. Knepley PetscErrorCode ierr; 12269a800dd8SMatthew G. Knepley 12279a800dd8SMatthew G. Knepley PetscFunctionBegin; 12282764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 12292764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 12302764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1231b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 12329a800dd8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 12339a800dd8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 12349a800dd8SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 12359a800dd8SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 12369a800dd8SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 12379a800dd8SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 12389a800dd8SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 123921454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, &numPoints, NULL, NULL);CHKERRQ(ierr); 124021454ff5SMatthew G. Knepley ierr = PetscMalloc3(dim*numPoints,&v0,dim*dim*numPoints,&J,numPoints,&detJ);CHKERRQ(ierr); 12419a800dd8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 12429a800dd8SMatthew G. Knepley PetscCellGeometry geom; 12439a800dd8SMatthew G. Knepley 124421454ff5SMatthew G. Knepley ierr = DMDAComputeCellGeometry(dm, c, q, v0, J, NULL, detJ);CHKERRQ(ierr); 12459a800dd8SMatthew G. Knepley geom.v0 = v0; 12469a800dd8SMatthew G. Knepley geom.J = J; 12479a800dd8SMatthew G. Knepley geom.detJ = detJ; 12489a800dd8SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 1249c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[f] : NULL; 1250b7a4d31fSMatthew G. Knepley 12512764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1252b7a4d31fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 1253b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 1254b7a4d31fSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 12559a800dd8SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 1256b7a4d31fSMatthew G. Knepley ierr = PetscDualSpaceApply(sp, d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 12579a800dd8SMatthew G. Knepley v += numComp; 12589a800dd8SMatthew G. Knepley } 12599a800dd8SMatthew G. Knepley } 12609a800dd8SMatthew G. Knepley ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 12619a800dd8SMatthew G. Knepley } 12629a800dd8SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 12639a800dd8SMatthew G. Knepley ierr = PetscFree3(v0,J,detJ);CHKERRQ(ierr); 12649a800dd8SMatthew G. Knepley PetscFunctionReturn(0); 12659a800dd8SMatthew G. Knepley } 12669a800dd8SMatthew G. Knepley 12679a800dd8SMatthew G. Knepley #undef __FUNCT__ 12689a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunction" 12699a800dd8SMatthew G. Knepley /*@C 12709a800dd8SMatthew G. Knepley DMDAProjectFunction - This projects the given function into the function space provided. 12719a800dd8SMatthew G. Knepley 12729a800dd8SMatthew G. Knepley Input Parameters: 12739a800dd8SMatthew G. Knepley + dm - The DM 12749a800dd8SMatthew G. Knepley . funcs - The coordinate functions to evaluate 1275c110b1eeSGeoffrey Irving . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 12769a800dd8SMatthew G. Knepley - mode - The insertion mode for values 12779a800dd8SMatthew G. Knepley 12789a800dd8SMatthew G. Knepley Output Parameter: 12799a800dd8SMatthew G. Knepley . X - vector 12809a800dd8SMatthew G. Knepley 12819a800dd8SMatthew G. Knepley Level: developer 12829a800dd8SMatthew G. Knepley 12839a800dd8SMatthew G. Knepley .seealso: DMDAComputeL2Diff() 12849a800dd8SMatthew G. Knepley @*/ 1285b7a4d31fSMatthew G. Knepley PetscErrorCode DMDAProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 12869a800dd8SMatthew G. Knepley { 12879a800dd8SMatthew G. Knepley Vec localX; 12889a800dd8SMatthew G. Knepley PetscErrorCode ierr; 12899a800dd8SMatthew G. Knepley 12909a800dd8SMatthew G. Knepley PetscFunctionBegin; 12919a800dd8SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12929a800dd8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 1293b7a4d31fSMatthew G. Knepley ierr = DMDAProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr); 12949a800dd8SMatthew G. Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 12959a800dd8SMatthew G. Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 12969a800dd8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 12979a800dd8SMatthew G. Knepley PetscFunctionReturn(0); 12989a800dd8SMatthew G. Knepley } 12999a800dd8SMatthew G. Knepley 13009a800dd8SMatthew G. Knepley #undef __FUNCT__ 13019a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2Diff" 13029a800dd8SMatthew G. Knepley /*@C 13039a800dd8SMatthew G. Knepley DMDAComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 13049a800dd8SMatthew G. Knepley 13059a800dd8SMatthew G. Knepley Input Parameters: 13069a800dd8SMatthew G. Knepley + dm - The DM 13079a800dd8SMatthew G. Knepley . funcs - The functions to evaluate for each field component 1308c110b1eeSGeoffrey Irving . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 13099a800dd8SMatthew G. Knepley - X - The coefficient vector u_h 13109a800dd8SMatthew G. Knepley 13119a800dd8SMatthew G. Knepley Output Parameter: 13129a800dd8SMatthew G. Knepley . diff - The diff ||u - u_h||_2 13139a800dd8SMatthew G. Knepley 13149a800dd8SMatthew G. Knepley Level: developer 13159a800dd8SMatthew G. Knepley 131623d86601SMatthew G. Knepley .seealso: DMDAProjectFunction(), DMDAComputeL2GradientDiff() 13179a800dd8SMatthew G. Knepley @*/ 1318b7a4d31fSMatthew G. Knepley PetscErrorCode DMDAComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 13199a800dd8SMatthew G. Knepley { 13209a800dd8SMatthew G. Knepley const PetscInt debug = 0; 13212764a2aaSMatthew G. Knepley PetscDS prob; 1322b7a4d31fSMatthew G. Knepley PetscFE fe; 13239a800dd8SMatthew G. Knepley PetscQuadrature quad; 1324b7a4d31fSMatthew G. Knepley PetscSection section; 13259a800dd8SMatthew G. Knepley Vec localX; 13269a800dd8SMatthew G. Knepley PetscScalar *funcVal; 13279a800dd8SMatthew G. Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 13289a800dd8SMatthew G. Knepley PetscReal localDiff = 0.0; 1329b7a4d31fSMatthew G. Knepley PetscInt dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp; 13309a800dd8SMatthew G. Knepley PetscErrorCode ierr; 13319a800dd8SMatthew G. Knepley 13329a800dd8SMatthew G. Knepley PetscFunctionBegin; 13339a800dd8SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 13342764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 13352764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 13362764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1337b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 13389a800dd8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 13399a800dd8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 13409a800dd8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 13419a800dd8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 13429a800dd8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 13439a800dd8SMatthew G. Knepley /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 1344b7a4d31fSMatthew G. Knepley ierr = PetscMalloc5(Nc,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 13459a800dd8SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 13469a800dd8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 13479a800dd8SMatthew G. Knepley PetscScalar *x = NULL; 13489a800dd8SMatthew G. Knepley PetscReal elemDiff = 0.0; 13499a800dd8SMatthew G. Knepley 135021454ff5SMatthew G. Knepley ierr = DMDAComputeCellGeometry(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr); 13519a800dd8SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 13529a800dd8SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 13539a800dd8SMatthew G. Knepley 13549a800dd8SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 1355c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[field] : NULL; 135621454ff5SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 13579a800dd8SMatthew G. Knepley PetscReal *basis; 135821454ff5SMatthew G. Knepley PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 13599a800dd8SMatthew G. Knepley 13602764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 136121454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 1362b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr); 1363b7a4d31fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr); 1364b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr); 13659a800dd8SMatthew G. Knepley if (debug) { 13669a800dd8SMatthew G. Knepley char title[1024]; 13679a800dd8SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 13689a800dd8SMatthew G. Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 13699a800dd8SMatthew G. Knepley } 13709a800dd8SMatthew G. Knepley for (q = 0; q < numQuadPoints; ++q) { 13719a800dd8SMatthew G. Knepley for (d = 0; d < dim; d++) { 13729a800dd8SMatthew G. Knepley coords[d] = v0[d]; 13739a800dd8SMatthew G. Knepley for (e = 0; e < dim; e++) { 13749a800dd8SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 13759a800dd8SMatthew G. Knepley } 13769a800dd8SMatthew G. Knepley } 1377c110b1eeSGeoffrey Irving (*funcs[field])(coords, funcVal, ctx); 13789a800dd8SMatthew G. Knepley for (fc = 0; fc < numBasisComps; ++fc) { 13799a800dd8SMatthew G. Knepley PetscScalar interpolant = 0.0; 13809a800dd8SMatthew G. Knepley 13819a800dd8SMatthew G. Knepley for (f = 0; f < numBasisFuncs; ++f) { 13829a800dd8SMatthew G. Knepley const PetscInt fidx = f*numBasisComps+fc; 13839a800dd8SMatthew G. Knepley interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 13849a800dd8SMatthew G. Knepley } 13859a800dd8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 13869a800dd8SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 13879a800dd8SMatthew G. Knepley } 13889a800dd8SMatthew G. Knepley } 13899a800dd8SMatthew G. Knepley comp += numBasisComps; 13909a800dd8SMatthew G. Knepley fieldOffset += numBasisFuncs*numBasisComps; 13919a800dd8SMatthew G. Knepley } 13929a800dd8SMatthew G. Knepley ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 13939a800dd8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 13949a800dd8SMatthew G. Knepley localDiff += elemDiff; 13959a800dd8SMatthew G. Knepley } 13969a800dd8SMatthew G. Knepley ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 13979a800dd8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 13989a800dd8SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 13999a800dd8SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 1400a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 1401a66d4d66SMatthew G Knepley } 1402a66d4d66SMatthew G Knepley 140323d86601SMatthew G. Knepley #undef __FUNCT__ 140423d86601SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2GradientDiff" 140523d86601SMatthew G. Knepley /*@C 140623d86601SMatthew G. Knepley DMDAComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 140723d86601SMatthew G. Knepley 140823d86601SMatthew G. Knepley Input Parameters: 140923d86601SMatthew G. Knepley + dm - The DM 141023d86601SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component 141123d86601SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 141223d86601SMatthew G. Knepley . X - The coefficient vector u_h 141323d86601SMatthew G. Knepley - n - The vector to project along 141423d86601SMatthew G. Knepley 141523d86601SMatthew G. Knepley Output Parameter: 141623d86601SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2 141723d86601SMatthew G. Knepley 141823d86601SMatthew G. Knepley Level: developer 141923d86601SMatthew G. Knepley 142023d86601SMatthew G. Knepley .seealso: DMDAProjectFunction(), DMPlexComputeL2Diff() 142123d86601SMatthew G. Knepley @*/ 1422b7a4d31fSMatthew G. Knepley PetscErrorCode DMDAComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 142323d86601SMatthew G. Knepley { 142423d86601SMatthew G. Knepley const PetscInt debug = 0; 14252764a2aaSMatthew G. Knepley PetscDS prob; 1426b7a4d31fSMatthew G. Knepley PetscFE fe; 142723d86601SMatthew G. Knepley PetscQuadrature quad; 1428b7a4d31fSMatthew G. Knepley PetscSection section; 142923d86601SMatthew G. Knepley Vec localX; 143023d86601SMatthew G. Knepley PetscScalar *funcVal, *interpolantVec; 143123d86601SMatthew G. Knepley PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 143223d86601SMatthew G. Knepley PetscReal localDiff = 0.0; 1433b7a4d31fSMatthew G. Knepley PetscInt dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp; 143423d86601SMatthew G. Knepley PetscErrorCode ierr; 143523d86601SMatthew G. Knepley 143623d86601SMatthew G. Knepley PetscFunctionBegin; 143723d86601SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 14382764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 14392764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 14402764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1441b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 144223d86601SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 144323d86601SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 144423d86601SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 144523d86601SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 144623d86601SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 144723d86601SMatthew G. Knepley /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 1448b7a4d31fSMatthew G. Knepley ierr = PetscMalloc7(Nc,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 144923d86601SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 145023d86601SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 145123d86601SMatthew G. Knepley PetscScalar *x = NULL; 145223d86601SMatthew G. Knepley PetscReal elemDiff = 0.0; 145323d86601SMatthew G. Knepley 145423d86601SMatthew G. Knepley ierr = DMDAComputeCellGeometry(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr); 145523d86601SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 145623d86601SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 145723d86601SMatthew G. Knepley 145823d86601SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 145923d86601SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 146023d86601SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 146123d86601SMatthew G. Knepley PetscReal *basisDer; 146223d86601SMatthew G. Knepley PetscInt Nq, Nb, Nc, q, d, e, fc, f, g; 146323d86601SMatthew G. Knepley 14642764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 146523d86601SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1466b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1467b7a4d31fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1468b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 146923d86601SMatthew G. Knepley if (debug) { 147023d86601SMatthew G. Knepley char title[1024]; 147123d86601SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 147223d86601SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 147323d86601SMatthew G. Knepley } 147423d86601SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 147523d86601SMatthew G. Knepley for (d = 0; d < dim; d++) { 147623d86601SMatthew G. Knepley coords[d] = v0[d]; 147723d86601SMatthew G. Knepley for (e = 0; e < dim; e++) { 147823d86601SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 147923d86601SMatthew G. Knepley } 148023d86601SMatthew G. Knepley } 148123d86601SMatthew G. Knepley (*funcs[field])(coords, n, funcVal, ctx); 148223d86601SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 148323d86601SMatthew G. Knepley PetscScalar interpolant = 0.0; 148423d86601SMatthew G. Knepley 148523d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 148623d86601SMatthew G. Knepley for (f = 0; f < Nb; ++f) { 148723d86601SMatthew G. Knepley const PetscInt fidx = f*Nc+fc; 148823d86601SMatthew G. Knepley 148923d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) { 149023d86601SMatthew G. Knepley realSpaceDer[d] = 0.0; 149123d86601SMatthew G. Knepley for (g = 0; g < dim; ++g) { 149223d86601SMatthew G. Knepley realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Nc+fidx)*dim+g]; 149323d86601SMatthew G. Knepley } 149423d86601SMatthew G. Knepley interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 149523d86601SMatthew G. Knepley } 149623d86601SMatthew G. Knepley } 149723d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 149823d86601SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 149923d86601SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 150023d86601SMatthew G. Knepley } 150123d86601SMatthew G. Knepley } 150223d86601SMatthew G. Knepley comp += Nc; 150323d86601SMatthew G. Knepley fieldOffset += Nb*Nc; 150423d86601SMatthew G. Knepley } 150523d86601SMatthew G. Knepley ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 150623d86601SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 150723d86601SMatthew G. Knepley localDiff += elemDiff; 150823d86601SMatthew G. Knepley } 150923d86601SMatthew G. Knepley ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 151023d86601SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 151123d86601SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 151223d86601SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 151323d86601SMatthew G. Knepley PetscFunctionReturn(0); 151423d86601SMatthew G. Knepley } 151523d86601SMatthew G. Knepley 151647c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 151747c6ae99SBarry Smith 151847c6ae99SBarry Smith #undef __FUNCT__ 1519aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray" 152047c6ae99SBarry Smith /*@C 1521aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 152247c6ae99SBarry Smith 152347c6ae99SBarry Smith Input Parameter: 152447c6ae99SBarry Smith + da - information about my local patch 152547c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 152647c6ae99SBarry Smith 152747c6ae99SBarry Smith Output Parameters: 152847c6ae99SBarry Smith . vptr - array data structured 152947c6ae99SBarry Smith 153047c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 153147c6ae99SBarry Smith to zero them. 153247c6ae99SBarry Smith 153347c6ae99SBarry Smith Level: advanced 153447c6ae99SBarry Smith 1535bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 153647c6ae99SBarry Smith 153747c6ae99SBarry Smith @*/ 15387087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 153947c6ae99SBarry Smith { 154047c6ae99SBarry Smith PetscErrorCode ierr; 154147c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 154247c6ae99SBarry Smith char *iarray_start; 154347c6ae99SBarry Smith void **iptr = (void**)vptr; 154447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 154547c6ae99SBarry Smith 154647c6ae99SBarry Smith PetscFunctionBegin; 154747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 154847c6ae99SBarry Smith if (ghosted) { 1549aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 155047c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 155147c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 155247c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 15530298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 15540298fd71SBarry Smith dd->startghostedin[i] = NULL; 155547c6ae99SBarry Smith 155647c6ae99SBarry Smith goto done; 155747c6ae99SBarry Smith } 155847c6ae99SBarry Smith } 155947c6ae99SBarry Smith xs = dd->Xs; 156047c6ae99SBarry Smith ys = dd->Ys; 156147c6ae99SBarry Smith zs = dd->Zs; 156247c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 156347c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 156447c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 156547c6ae99SBarry Smith } else { 1566aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 156747c6ae99SBarry Smith if (dd->arrayin[i]) { 156847c6ae99SBarry Smith *iptr = dd->arrayin[i]; 156947c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 15700298fd71SBarry Smith dd->arrayin[i] = NULL; 15710298fd71SBarry Smith dd->startin[i] = NULL; 157247c6ae99SBarry Smith 157347c6ae99SBarry Smith goto done; 157447c6ae99SBarry Smith } 157547c6ae99SBarry Smith } 157647c6ae99SBarry Smith xs = dd->xs; 157747c6ae99SBarry Smith ys = dd->ys; 157847c6ae99SBarry Smith zs = dd->zs; 157947c6ae99SBarry Smith xm = dd->xe-dd->xs; 158047c6ae99SBarry Smith ym = dd->ye-dd->ys; 158147c6ae99SBarry Smith zm = dd->ze-dd->zs; 158247c6ae99SBarry Smith } 158347c6ae99SBarry Smith 1584*c73cfb54SMatthew G. Knepley switch (da->dim) { 158547c6ae99SBarry Smith case 1: { 158647c6ae99SBarry Smith void *ptr; 158747c6ae99SBarry Smith 1588901f1932SJed Brown ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 158947c6ae99SBarry Smith 159047c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 159147c6ae99SBarry Smith *iptr = (void*)ptr; 15928865f1eaSKarl Rupp break; 15938865f1eaSKarl Rupp } 159447c6ae99SBarry Smith case 2: { 159547c6ae99SBarry Smith void **ptr; 159647c6ae99SBarry Smith 1597901f1932SJed Brown ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 159847c6ae99SBarry Smith 159947c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 16008865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 160147c6ae99SBarry Smith *iptr = (void*)ptr; 16028865f1eaSKarl Rupp break; 16038865f1eaSKarl Rupp } 160447c6ae99SBarry Smith case 3: { 160547c6ae99SBarry Smith void ***ptr,**bptr; 160647c6ae99SBarry Smith 1607901f1932SJed Brown ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 160847c6ae99SBarry Smith 160947c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 161047c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 16118865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 161247c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 161347c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 161447c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 161547c6ae99SBarry Smith } 161647c6ae99SBarry Smith } 161747c6ae99SBarry Smith *iptr = (void*)ptr; 16188865f1eaSKarl Rupp break; 16198865f1eaSKarl Rupp } 162047c6ae99SBarry Smith default: 1621*c73cfb54SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",da->dim); 162247c6ae99SBarry Smith } 162347c6ae99SBarry Smith 162447c6ae99SBarry Smith done: 162547c6ae99SBarry Smith /* add arrays to the checked out list */ 162647c6ae99SBarry Smith if (ghosted) { 1627aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 162847c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 162947c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 163047c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 163147c6ae99SBarry Smith break; 163247c6ae99SBarry Smith } 163347c6ae99SBarry Smith } 163447c6ae99SBarry Smith } else { 1635aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 163647c6ae99SBarry Smith if (!dd->arrayout[i]) { 163747c6ae99SBarry Smith dd->arrayout[i] = *iptr; 163847c6ae99SBarry Smith dd->startout[i] = iarray_start; 163947c6ae99SBarry Smith break; 164047c6ae99SBarry Smith } 164147c6ae99SBarry Smith } 164247c6ae99SBarry Smith } 164347c6ae99SBarry Smith PetscFunctionReturn(0); 164447c6ae99SBarry Smith } 164547c6ae99SBarry Smith 164647c6ae99SBarry Smith #undef __FUNCT__ 1647aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray" 164847c6ae99SBarry Smith /*@C 1649aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 165047c6ae99SBarry Smith 165147c6ae99SBarry Smith Input Parameter: 165247c6ae99SBarry Smith + da - information about my local patch 165347c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 165447c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 165547c6ae99SBarry Smith 165647c6ae99SBarry Smith Level: advanced 165747c6ae99SBarry Smith 1658bcaeba4dSBarry Smith .seealso: DMDAGetArray() 165947c6ae99SBarry Smith 166047c6ae99SBarry Smith @*/ 16617087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 166247c6ae99SBarry Smith { 166347c6ae99SBarry Smith PetscInt i; 166447c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 166547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 166647c6ae99SBarry Smith 166747c6ae99SBarry Smith PetscFunctionBegin; 166847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 166947c6ae99SBarry Smith if (ghosted) { 1670aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 167147c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 167247c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 16730298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 16740298fd71SBarry Smith dd->startghostedout[i] = NULL; 167547c6ae99SBarry Smith break; 167647c6ae99SBarry Smith } 167747c6ae99SBarry Smith } 1678aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 167947c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 168047c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 168147c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 168247c6ae99SBarry Smith break; 168347c6ae99SBarry Smith } 168447c6ae99SBarry Smith } 168547c6ae99SBarry Smith } else { 1686aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 168747c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 168847c6ae99SBarry Smith iarray_start = dd->startout[i]; 16890298fd71SBarry Smith dd->arrayout[i] = NULL; 16900298fd71SBarry Smith dd->startout[i] = NULL; 169147c6ae99SBarry Smith break; 169247c6ae99SBarry Smith } 169347c6ae99SBarry Smith } 1694aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 169547c6ae99SBarry Smith if (!dd->arrayin[i]) { 169647c6ae99SBarry Smith dd->arrayin[i] = *iptr; 169747c6ae99SBarry Smith dd->startin[i] = iarray_start; 169847c6ae99SBarry Smith break; 169947c6ae99SBarry Smith } 170047c6ae99SBarry Smith } 170147c6ae99SBarry Smith } 170247c6ae99SBarry Smith PetscFunctionReturn(0); 170347c6ae99SBarry Smith } 170447c6ae99SBarry Smith 1705