147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry 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 */ 181e6b0712SBarry Smith static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine) 1947c6ae99SBarry Smith { 2047c6ae99SBarry Smith PetscErrorCode ierr; 2147c6ae99SBarry Smith PetscInt n,m; 2247c6ae99SBarry Smith Vec vec = (Vec)obj; 2347c6ae99SBarry Smith PetscScalar *array; 2447c6ae99SBarry Smith mxArray *mat; 259a42bb27SBarry Smith DM da; 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith PetscFunctionBegin; 28c688c046SMatthew G Knepley ierr = VecGetDM(vec, &da);CHKERRQ(ierr); 29ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 30aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr); 3147c6ae99SBarry Smith 3247c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 3347c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3447c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3547c6ae99SBarry Smith #else 3647c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3747c6ae99SBarry Smith #endif 3847c6ae99SBarry Smith ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr); 3947c6ae99SBarry Smith ierr = PetscObjectName(obj);CHKERRQ(ierr); 4047c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 4147c6ae99SBarry Smith 4247c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 4347c6ae99SBarry Smith PetscFunctionReturn(0); 4447c6ae99SBarry Smith } 4547c6ae99SBarry Smith #endif 4647c6ae99SBarry Smith 4747c6ae99SBarry Smith 487087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 4947c6ae99SBarry Smith { 5047c6ae99SBarry Smith PetscErrorCode ierr; 5147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith PetscFunctionBegin; 5447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5547c6ae99SBarry Smith PetscValidPointer(g,2); 5611689aebSJed Brown if (da->defaultSection) { 5711689aebSJed Brown ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr); 5811689aebSJed Brown } else { 5947c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 6047c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 6147c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 62401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 63c688c046SMatthew G Knepley ierr = VecSetDM(*g, da);CHKERRQ(ierr); 6447c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 657bc9226cSBarry Smith if (dd->w == 1 && da->dim == 2) { 66bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr); 6747c6ae99SBarry Smith } 6847c6ae99SBarry Smith #endif 6911689aebSJed Brown } 7047c6ae99SBarry Smith PetscFunctionReturn(0); 7147c6ae99SBarry Smith } 7247c6ae99SBarry Smith 73939f6067SMatthew G. Knepley /*@ 74939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 75939f6067SMatthew G. Knepley 76939f6067SMatthew G. Knepley Input Parameter: 77939f6067SMatthew G. Knepley . dm - The DM object 78939f6067SMatthew G. Knepley 79939f6067SMatthew G. Knepley Output Parameters: 80939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 81939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 82939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 83939f6067SMatthew G. Knepley - numCells - The number of local cells 84939f6067SMatthew G. Knepley 85939f6067SMatthew G. Knepley Level: developer 86939f6067SMatthew G. Knepley 87939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint() 88939f6067SMatthew G. Knepley @*/ 893582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 9057459e9aSMatthew G Knepley { 9157459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 92c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 9357459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 9457459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 9557459e9aSMatthew G Knepley 9657459e9aSMatthew G Knepley PetscFunctionBegin; 97d6401b93SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 983582350cSMatthew G. Knepley if (numCellsX) { 993582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 1003582350cSMatthew G. Knepley *numCellsX = mx; 1013582350cSMatthew G. Knepley } 1023582350cSMatthew G. Knepley if (numCellsY) { 1038135c375SStefano Zampini PetscValidIntPointer(numCellsY,3); 1043582350cSMatthew G. Knepley *numCellsY = my; 1053582350cSMatthew G. Knepley } 1063582350cSMatthew G. Knepley if (numCellsZ) { 1078135c375SStefano Zampini PetscValidIntPointer(numCellsZ,4); 1083582350cSMatthew G. Knepley *numCellsZ = mz; 1093582350cSMatthew G. Knepley } 11057459e9aSMatthew G Knepley if (numCells) { 1113582350cSMatthew G. Knepley PetscValidIntPointer(numCells,5); 11257459e9aSMatthew G Knepley *numCells = nC; 11357459e9aSMatthew G Knepley } 11457459e9aSMatthew G Knepley PetscFunctionReturn(0); 11557459e9aSMatthew G Knepley } 11657459e9aSMatthew G Knepley 117d6401b93SMatthew G. Knepley /*@ 118d6401b93SMatthew G. Knepley DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 119d6401b93SMatthew G. Knepley 120d6401b93SMatthew G. Knepley Input Parameters: 121d6401b93SMatthew G. Knepley + dm - The DM object 122d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell 123d6401b93SMatthew G. Knepley 124d6401b93SMatthew G. Knepley Output Parameters: 125d6401b93SMatthew G. Knepley . point - The local DM point 126d6401b93SMatthew G. Knepley 127d6401b93SMatthew G. Knepley Level: developer 128d6401b93SMatthew G. Knepley 129d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells() 130d6401b93SMatthew G. Knepley @*/ 131d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) 132d6401b93SMatthew G. Knepley { 133c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 134d6401b93SMatthew G. Knepley DMDALocalInfo info; 135d6401b93SMatthew G. Knepley PetscErrorCode ierr; 136d6401b93SMatthew G. Knepley 137d6401b93SMatthew G. Knepley PetscFunctionBegin; 138d6401b93SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 139d6401b93SMatthew G. Knepley PetscValidIntPointer(point,5); 140d6401b93SMatthew G. Knepley ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 141d6401b93SMatthew 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);} 14263b6bffdSJed Brown if (dim > 1) {if ((j < info.gys) || (j >= info.gys+info.gym)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %d not in [%d, %d)", j, info.gys, info.gys+info.gym);} 14363b6bffdSJed Brown if (dim > 2) {if ((k < info.gzs) || (k >= info.gzs+info.gzm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %d not in [%d, %d)", k, info.gzs, info.gzs+info.gzm);} 144d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0); 145d6401b93SMatthew G. Knepley PetscFunctionReturn(0); 146d6401b93SMatthew G. Knepley } 147d6401b93SMatthew G. Knepley 14857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 14957459e9aSMatthew G Knepley { 15057459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 151c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 15257459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 15357459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 15457459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 15557459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 15657459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 15757459e9aSMatthew G Knepley 15857459e9aSMatthew G Knepley PetscFunctionBegin; 15957459e9aSMatthew G Knepley if (numVerticesX) { 16057459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 16157459e9aSMatthew G Knepley *numVerticesX = nVx; 16257459e9aSMatthew G Knepley } 16357459e9aSMatthew G Knepley if (numVerticesY) { 16457459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 16557459e9aSMatthew G Knepley *numVerticesY = nVy; 16657459e9aSMatthew G Knepley } 16757459e9aSMatthew G Knepley if (numVerticesZ) { 16857459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 16957459e9aSMatthew G Knepley *numVerticesZ = nVz; 17057459e9aSMatthew G Knepley } 17157459e9aSMatthew G Knepley if (numVertices) { 17257459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 17357459e9aSMatthew G Knepley *numVertices = nV; 17457459e9aSMatthew G Knepley } 17557459e9aSMatthew G Knepley PetscFunctionReturn(0); 17657459e9aSMatthew G Knepley } 17757459e9aSMatthew G Knepley 17857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 17957459e9aSMatthew G Knepley { 18057459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 181c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 18257459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 18357459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 18457459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 18557459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 18657459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 18757459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 18857459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 18957459e9aSMatthew G Knepley 19057459e9aSMatthew G Knepley PetscFunctionBegin; 19157459e9aSMatthew G Knepley if (numXFacesX) { 19257459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 19357459e9aSMatthew G Knepley *numXFacesX = nxF; 19457459e9aSMatthew G Knepley } 19557459e9aSMatthew G Knepley if (numXFaces) { 19657459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 19757459e9aSMatthew G Knepley *numXFaces = nXF; 19857459e9aSMatthew G Knepley } 19957459e9aSMatthew G Knepley if (numYFacesY) { 20057459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 20157459e9aSMatthew G Knepley *numYFacesY = nyF; 20257459e9aSMatthew G Knepley } 20357459e9aSMatthew G Knepley if (numYFaces) { 20457459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 20557459e9aSMatthew G Knepley *numYFaces = nYF; 20657459e9aSMatthew G Knepley } 20757459e9aSMatthew G Knepley if (numZFacesZ) { 20857459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 20957459e9aSMatthew G Knepley *numZFacesZ = nzF; 21057459e9aSMatthew G Knepley } 21157459e9aSMatthew G Knepley if (numZFaces) { 21257459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 21357459e9aSMatthew G Knepley *numZFaces = nZF; 21457459e9aSMatthew G Knepley } 21557459e9aSMatthew G Knepley PetscFunctionReturn(0); 21657459e9aSMatthew G Knepley } 21757459e9aSMatthew G Knepley 21857459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 21957459e9aSMatthew G Knepley { 220c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 22157459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 22257459e9aSMatthew G Knepley PetscErrorCode ierr; 22357459e9aSMatthew G Knepley 22457459e9aSMatthew G Knepley PetscFunctionBegin; 2258865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 2268865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 2273582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2280298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2290298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 23057459e9aSMatthew G Knepley if (height == 0) { 23157459e9aSMatthew G Knepley /* Cells */ 2328865f1eaSKarl Rupp if (pStart) *pStart = 0; 2338865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 23457459e9aSMatthew G Knepley } else if (height == 1) { 23557459e9aSMatthew G Knepley /* Faces */ 2368865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 2378865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 23857459e9aSMatthew G Knepley } else if (height == dim) { 23957459e9aSMatthew G Knepley /* Vertices */ 2408865f1eaSKarl Rupp if (pStart) *pStart = nC; 2418865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 24257459e9aSMatthew G Knepley } else if (height < 0) { 24357459e9aSMatthew G Knepley /* All points */ 2448865f1eaSKarl Rupp if (pStart) *pStart = 0; 2458865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 24682f516ccSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 24757459e9aSMatthew G Knepley PetscFunctionReturn(0); 24857459e9aSMatthew G Knepley } 24957459e9aSMatthew G Knepley 2503385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) 2513385ff7eSMatthew G. Knepley { 252c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2533385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2543385ff7eSMatthew G. Knepley PetscErrorCode ierr; 2553385ff7eSMatthew G. Knepley 2563385ff7eSMatthew G. Knepley PetscFunctionBegin; 2573385ff7eSMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart,3); 2583385ff7eSMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd,4); 2593385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2603385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2613385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 2623385ff7eSMatthew G. Knepley if (depth == dim) { 2633385ff7eSMatthew G. Knepley /* Cells */ 2643385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2653385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC; 2663385ff7eSMatthew G. Knepley } else if (depth == dim-1) { 2673385ff7eSMatthew G. Knepley /* Faces */ 2683385ff7eSMatthew G. Knepley if (pStart) *pStart = nC+nV; 2693385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2703385ff7eSMatthew G. Knepley } else if (depth == 0) { 2713385ff7eSMatthew G. Knepley /* Vertices */ 2723385ff7eSMatthew G. Knepley if (pStart) *pStart = nC; 2733385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV; 2743385ff7eSMatthew G. Knepley } else if (depth < 0) { 2753385ff7eSMatthew G. Knepley /* All points */ 2763385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2773385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2783385ff7eSMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth); 2793385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2803385ff7eSMatthew G. Knepley } 2813385ff7eSMatthew G. Knepley 2823385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) 2833385ff7eSMatthew G. Knepley { 284c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2853385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2863385ff7eSMatthew G. Knepley PetscErrorCode ierr; 2873385ff7eSMatthew G. Knepley 2883385ff7eSMatthew G. Knepley PetscFunctionBegin; 2893385ff7eSMatthew G. Knepley *coneSize = 0; 2903385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2913385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2923385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 2933385ff7eSMatthew G. Knepley switch (dim) { 2943385ff7eSMatthew G. Knepley case 2: 2953385ff7eSMatthew G. Knepley if (p >= 0) { 2963385ff7eSMatthew G. Knepley if (p < nC) { 2973385ff7eSMatthew G. Knepley *coneSize = 4; 2983385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 2993385ff7eSMatthew G. Knepley *coneSize = 0; 3003385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF+nZF) { 3013385ff7eSMatthew G. Knepley *coneSize = 2; 3023385ff7eSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 3033385ff7eSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3043385ff7eSMatthew G. Knepley break; 3053385ff7eSMatthew G. Knepley case 3: 3063385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3073385ff7eSMatthew G. Knepley break; 3083385ff7eSMatthew G. Knepley } 3093385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3103385ff7eSMatthew G. Knepley } 3113385ff7eSMatthew G. Knepley 3123385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) 3133385ff7eSMatthew G. Knepley { 314c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 3153385ff7eSMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 3163385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3173385ff7eSMatthew G. Knepley 3183385ff7eSMatthew G. Knepley PetscFunctionBegin; 319*69291d52SBarry Smith if (!cone) {ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr);} 3203385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr); 3213385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 3223385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 3233385ff7eSMatthew G. Knepley switch (dim) { 3243385ff7eSMatthew G. Knepley case 2: 3253385ff7eSMatthew G. Knepley if (p >= 0) { 3263385ff7eSMatthew G. Knepley if (p < nC) { 3273385ff7eSMatthew G. Knepley const PetscInt cy = p / nCx; 3283385ff7eSMatthew G. Knepley const PetscInt cx = p % nCx; 3293385ff7eSMatthew G. Knepley 3303385ff7eSMatthew G. Knepley (*cone)[0] = cy*nxF + cx + nC+nV; 3313385ff7eSMatthew G. Knepley (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF; 3323385ff7eSMatthew G. Knepley (*cone)[2] = cy*nxF + cx + nxF + nC+nV; 3333385ff7eSMatthew G. Knepley (*cone)[3] = cx*nyF + cy + nC+nV+nXF; 3343385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 3353385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 3363385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF) { 3373385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV) / nxF; 3383385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV) % nxF; 3393385ff7eSMatthew G. Knepley 3403385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3413385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + 1 + nC; 3423385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF) { 3433385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV+nXF) / nyF; 3443385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV+nXF) % nyF; 3453385ff7eSMatthew G. Knepley 3463385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3473385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + nVx + nC; 3483385ff7eSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 3493385ff7eSMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3503385ff7eSMatthew G. Knepley break; 3513385ff7eSMatthew G. Knepley case 3: 3523385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3533385ff7eSMatthew G. Knepley break; 3543385ff7eSMatthew G. Knepley } 3553385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3563385ff7eSMatthew G. Knepley } 3573385ff7eSMatthew G. Knepley 3583385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) 3593385ff7eSMatthew G. Knepley { 3603385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3613385ff7eSMatthew G. Knepley 3623385ff7eSMatthew G. Knepley PetscFunctionBegin; 363*69291d52SBarry Smith ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr); 3643385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3653385ff7eSMatthew G. Knepley } 3663385ff7eSMatthew G. Knepley 367a66d4d66SMatthew G Knepley /*@C 368a66d4d66SMatthew G Knepley DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with 369a66d4d66SMatthew G Knepley different numbers of dofs on vertices, cells, and faces in each direction. 370a66d4d66SMatthew G Knepley 371a66d4d66SMatthew G Knepley Input Parameters: 372a66d4d66SMatthew G Knepley + dm- The DMDA 373a66d4d66SMatthew G Knepley . numFields - The number of fields 374affa55c7SMatthew G. Knepley . numComp - The number of components in each field 375affa55c7SMatthew G. Knepley . numDof - The number of dofs per dimension for each field 3760298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL 377a66d4d66SMatthew G Knepley 378a66d4d66SMatthew G Knepley Level: developer 379a66d4d66SMatthew G Knepley 380a66d4d66SMatthew G Knepley Note: 381a66d4d66SMatthew G Knepley The default DMDA numbering is as follows: 382a66d4d66SMatthew G Knepley 383a66d4d66SMatthew G Knepley - Cells: [0, nC) 384a66d4d66SMatthew G Knepley - Vertices: [nC, nC+nV) 38588ed4aceSMatthew G Knepley - X-Faces: [nC+nV, nC+nV+nXF) normal is +- x-dir 38688ed4aceSMatthew G Knepley - Y-Faces: [nC+nV+nXF, nC+nV+nXF+nYF) normal is +- y-dir 38788ed4aceSMatthew G Knepley - Z-Faces: [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir 388a66d4d66SMatthew G Knepley 389a66d4d66SMatthew G Knepley We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment. 390a66d4d66SMatthew G Knepley @*/ 39144e1f9abSMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, const PetscInt numComp[], const PetscInt numDof[], const PetscInt numFaceDof[], PetscSection *s) 392a66d4d66SMatthew G Knepley { 393a4b60ecfSMatthew G. Knepley PetscSection section; 394c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 39580800b1aSMatthew G Knepley PetscInt numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0}; 396b2fff234SMatthew G. Knepley PetscBT isLeaf; 39788ed4aceSMatthew G Knepley PetscSF sf; 398feafbc5cSMatthew G Knepley PetscMPIInt rank; 39988ed4aceSMatthew G Knepley const PetscMPIInt *neighbors; 40088ed4aceSMatthew G Knepley PetscInt *localPoints; 40188ed4aceSMatthew G Knepley PetscSFNode *remotePoints; 402f5eeb527SMatthew G Knepley PetscInt nleaves = 0, nleavesCheck = 0, nL = 0; 40357459e9aSMatthew G Knepley PetscInt nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF; 40457459e9aSMatthew G Knepley PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd; 40588ed4aceSMatthew G Knepley PetscInt f, v, c, xf, yf, zf, xn, yn, zn; 406a66d4d66SMatthew G Knepley PetscErrorCode ierr; 407a66d4d66SMatthew G Knepley 408a66d4d66SMatthew G Knepley PetscFunctionBegin; 409a66d4d66SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4103582350cSMatthew G. Knepley PetscValidPointer(s, 4); 41182f516ccSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 4123582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 41357459e9aSMatthew G Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 41457459e9aSMatthew G Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 41557459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 41657459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 41757459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 41857459e9aSMatthew G Knepley ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 41957459e9aSMatthew G Knepley xfStart = vEnd; xfEnd = xfStart+nXF; 42057459e9aSMatthew G Knepley yfStart = xfEnd; yfEnd = yfStart+nYF; 42157459e9aSMatthew G Knepley zfStart = yfEnd; zfEnd = zfStart+nZF; 42282f516ccSBarry Smith if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd); 42388ed4aceSMatthew G Knepley /* Create local section */ 42480800b1aSMatthew G Knepley ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr); 425a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 426affa55c7SMatthew G. Knepley numVertexTotDof += numDof[f*(dim+1)+0]; 427affa55c7SMatthew G. Knepley numCellTotDof += numDof[f*(dim+1)+dim]; 428affa55c7SMatthew G. Knepley numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0; 429affa55c7SMatthew G. Knepley numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0; 430affa55c7SMatthew G. Knepley numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0; 431a66d4d66SMatthew G Knepley } 432a4b60ecfSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);CHKERRQ(ierr); 433affa55c7SMatthew G. Knepley if (numFields > 0) { 434a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr); 435a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 43680800b1aSMatthew G Knepley const char *name; 43780800b1aSMatthew G Knepley 43880800b1aSMatthew G Knepley ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr); 439affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr); 44080800b1aSMatthew G Knepley if (numComp) { 441a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr); 442a66d4d66SMatthew G Knepley } 443a66d4d66SMatthew G Knepley } 444a66d4d66SMatthew G Knepley } 445a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 446a66d4d66SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 447a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 448affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr); 449a66d4d66SMatthew G Knepley } 450a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr); 451a66d4d66SMatthew G Knepley } 452a66d4d66SMatthew G Knepley for (xf = xfStart; xf < xfEnd; ++xf) { 453a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 454affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 455a66d4d66SMatthew G Knepley } 456a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr); 457a66d4d66SMatthew G Knepley } 458a66d4d66SMatthew G Knepley for (yf = yfStart; yf < yfEnd; ++yf) { 459a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 460affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 461a66d4d66SMatthew G Knepley } 462a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr); 463a66d4d66SMatthew G Knepley } 464a66d4d66SMatthew G Knepley for (zf = zfStart; zf < zfEnd; ++zf) { 465a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 466affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr); 467a66d4d66SMatthew G Knepley } 468a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr); 469a66d4d66SMatthew G Knepley } 470a66d4d66SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 471a66d4d66SMatthew G Knepley for (f = 0; f < numFields; ++f) { 472affa55c7SMatthew G. Knepley ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr); 473a66d4d66SMatthew G Knepley } 474a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr); 475a66d4d66SMatthew G Knepley } 476a4b60ecfSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 47788ed4aceSMatthew G Knepley /* Create mesh point SF */ 478b2fff234SMatthew G. Knepley ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr); 47988ed4aceSMatthew G Knepley ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr); 48088ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 48188ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 48288ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 4837128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 48488ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 485b2fff234SMatthew G. Knepley PetscInt xv, yv, zv; 48688ed4aceSMatthew G Knepley 4873814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 488b2fff234SMatthew G. Knepley if (xp < 0) { /* left */ 489b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 490b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 491b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 492b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 493b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 494b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 495b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 496b2fff234SMatthew G. Knepley } else { 497b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 498b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 499b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 500b2fff234SMatthew G. Knepley } 501b2fff234SMatthew G. Knepley } 502b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 503b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 504b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 505b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 506b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 507b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 508b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 509b2fff234SMatthew G. Knepley } else { 510b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 511b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 512b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 513b2fff234SMatthew G. Knepley } 514b2fff234SMatthew G. Knepley } 515b2fff234SMatthew G. Knepley } else { 516b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 517b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 518b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 519b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 520b2fff234SMatthew G. Knepley } 521b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 522b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 523b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 524b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 525b2fff234SMatthew G. Knepley } 526b2fff234SMatthew G. Knepley } else { 527b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 528b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 529b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 530b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 531b2fff234SMatthew G. Knepley } 532b2fff234SMatthew G. Knepley } 533b2fff234SMatthew G. Knepley #if 0 534b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 535b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 536b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 537b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 538b2fff234SMatthew G. Knepley } 539b2fff234SMatthew G. Knepley #endif 540b2fff234SMatthew G. Knepley } 541b2fff234SMatthew G. Knepley } 542b2fff234SMatthew G. Knepley } else if (xp > 0) { /* right */ 543b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 544b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 545b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 546b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 547b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 548b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 549b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 550b2fff234SMatthew G. Knepley } else { 551b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 552b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 553b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 554b2fff234SMatthew G. Knepley } 555b2fff234SMatthew G. Knepley } 556b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 557b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 558b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 559b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 560b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 561b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 562b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 563b2fff234SMatthew G. Knepley } else { 564b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 565b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 566b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 567b2fff234SMatthew G. Knepley } 568b2fff234SMatthew G. Knepley } 569b2fff234SMatthew G. Knepley } else { 570b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 571b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 572b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 573b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 574b2fff234SMatthew G. Knepley } 575b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 576b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 577b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 578b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 579b2fff234SMatthew G. Knepley } 580b2fff234SMatthew G. Knepley } else { 581b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 582b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 583b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 584b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 585b2fff234SMatthew G. Knepley } 586b2fff234SMatthew G. Knepley } 587b2fff234SMatthew G. Knepley #if 0 588b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 589b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 590b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 591b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 592b2fff234SMatthew G. Knepley } 593b2fff234SMatthew G. Knepley #endif 594b2fff234SMatthew G. Knepley } 595b2fff234SMatthew G. Knepley } 596b2fff234SMatthew G. Knepley } else { 597b2fff234SMatthew G. Knepley if (yp < 0) { /* bottom */ 598b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 599b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 600b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 601b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 602b2fff234SMatthew G. Knepley } 603b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 604b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 605b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 606b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 607b2fff234SMatthew G. Knepley } 608b2fff234SMatthew G. Knepley } else { 609b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 610b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 611b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 612b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 613b2fff234SMatthew G. Knepley } 614b2fff234SMatthew G. Knepley } 615b2fff234SMatthew G. Knepley #if 0 616b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 617b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 618b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 619b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 620b2fff234SMatthew G. Knepley } 621b2fff234SMatthew G. Knepley #endif 622b2fff234SMatthew G. Knepley } 623b2fff234SMatthew G. Knepley } else if (yp > 0) { /* top */ 624b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 625b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 626b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 627b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 628b2fff234SMatthew G. Knepley } 629b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 630b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 631b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 632b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 633b2fff234SMatthew G. Knepley } 634b2fff234SMatthew G. Knepley } else { 635b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 636b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 637b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 638b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 639b2fff234SMatthew G. Knepley } 640b2fff234SMatthew G. Knepley } 641b2fff234SMatthew G. Knepley #if 0 642b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 643b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 644b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 645b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves; 646b2fff234SMatthew G. Knepley } 647b2fff234SMatthew G. Knepley #endif 648b2fff234SMatthew G. Knepley } 649b2fff234SMatthew G. Knepley } else { 650b2fff234SMatthew G. Knepley if (zp < 0) { /* back */ 651b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 652b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 653b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 654b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 655b2fff234SMatthew G. Knepley } 656b2fff234SMatthew G. Knepley } 657b2fff234SMatthew G. Knepley #if 0 658b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 659b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 660b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 661b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 662b2fff234SMatthew G. Knepley } 663b2fff234SMatthew G. Knepley #endif 664b2fff234SMatthew G. Knepley } else if (zp > 0) { /* front */ 665b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 666b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 667b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 668b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves; 669b2fff234SMatthew G. Knepley } 670b2fff234SMatthew G. Knepley } 671b2fff234SMatthew G. Knepley #if 0 672b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 673b2fff234SMatthew G. Knepley /* THIS IS WRONG */ 674b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 675b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves; 676b2fff234SMatthew G. Knepley } 677b2fff234SMatthew G. Knepley #endif 678b2fff234SMatthew G. Knepley } else { 679b2fff234SMatthew G. Knepley /* Nothing is shared from the interior */ 68088ed4aceSMatthew G Knepley } 68188ed4aceSMatthew G Knepley } 68288ed4aceSMatthew G Knepley } 68388ed4aceSMatthew G Knepley } 68488ed4aceSMatthew G Knepley } 685b2fff234SMatthew G. Knepley } 686b2fff234SMatthew G. Knepley } 687b2fff234SMatthew G. Knepley ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr); 6880d6c776aSBarry Smith ierr = PetscMalloc1(nleaves,&localPoints);CHKERRQ(ierr); 6890d6c776aSBarry Smith ierr = PetscMalloc1(nleaves,&remotePoints);CHKERRQ(ierr); 69088ed4aceSMatthew G Knepley for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) { 69188ed4aceSMatthew G Knepley for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) { 69288ed4aceSMatthew G Knepley for (xn = 0; xn < 3; ++xn) { 6937128ae9fSMatthew G Knepley const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0; 69488ed4aceSMatthew G Knepley const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn]; 695f5eeb527SMatthew G Knepley PetscInt xv, yv, zv; 69688ed4aceSMatthew G Knepley 6973814d604SMatthew G Knepley if (neighbor >= 0 && neighbor < rank) { 69888ed4aceSMatthew G Knepley if (xp < 0) { /* left */ 69988ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 70088ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 701b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + 0 + nC; /* left bottom back vertex */ 702f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 703f5eeb527SMatthew G Knepley 704b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 705f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 706f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 707f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 708f5eeb527SMatthew G Knepley ++nL; 709b2fff234SMatthew G. Knepley } 71088ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 711b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* left bottom front vertex */ 712f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 713f5eeb527SMatthew G Knepley 714b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 715f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 716f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 717f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 718f5eeb527SMatthew G Knepley ++nL; 719b2fff234SMatthew G. Knepley } 72088ed4aceSMatthew G Knepley } else { 721b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 722b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + 0 + nC; /* left bottom vertices */ 723f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 724f5eeb527SMatthew G Knepley 725b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 726f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 727f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 728f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 729b2fff234SMatthew G. Knepley ++nL; 730b2fff234SMatthew G. Knepley } 731f5eeb527SMatthew G Knepley } 73288ed4aceSMatthew G Knepley } 73388ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 73488ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 735b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* left top back vertex */ 736f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 737f5eeb527SMatthew G Knepley 738b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 739f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 740f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 741f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 742f5eeb527SMatthew G Knepley ++nL; 743b2fff234SMatthew G. Knepley } 74488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 745b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* left top front vertex */ 746f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*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; 752f5eeb527SMatthew G Knepley ++nL; 753b2fff234SMatthew G. Knepley } 75488ed4aceSMatthew G Knepley } else { 755b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 756b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* left top vertices */ 757f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 758f5eeb527SMatthew G Knepley 759b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 760f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 761f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 762f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 763b2fff234SMatthew G. Knepley ++nL; 764b2fff234SMatthew G. Knepley } 765f5eeb527SMatthew G Knepley } 76688ed4aceSMatthew G Knepley } 76788ed4aceSMatthew G Knepley } else { 76888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 769b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 770b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* left back vertices */ 771f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 772f5eeb527SMatthew G Knepley 773b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 774f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 775f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 776f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 777b2fff234SMatthew G. Knepley ++nL; 778b2fff234SMatthew G. Knepley } 779f5eeb527SMatthew G Knepley } 78088ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 781b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 782b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* left front vertices */ 783f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 784f5eeb527SMatthew G Knepley 785b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 786f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 787f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 788f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 789b2fff234SMatthew G. Knepley ++nL; 790b2fff234SMatthew G. Knepley } 791f5eeb527SMatthew G Knepley } 79288ed4aceSMatthew G Knepley } else { 793f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 794b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 795b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + 0 + nC; /* left vertices */ 796f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */ 797f5eeb527SMatthew G Knepley 798b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 799f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 800f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 801f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 802b2fff234SMatthew G. Knepley ++nL; 803f5eeb527SMatthew G Knepley } 804f5eeb527SMatthew G Knepley } 805b2fff234SMatthew G. Knepley } 806b2fff234SMatthew G. Knepley #if 0 807b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 808f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 809b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* left faces */ 810f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 811f5eeb527SMatthew G Knepley 812b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 813f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 814f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 815f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 816f5eeb527SMatthew G Knepley } 81788ed4aceSMatthew G Knepley } 818b2fff234SMatthew G. Knepley #endif 819b2fff234SMatthew G. Knepley } 82088ed4aceSMatthew G Knepley } 82188ed4aceSMatthew G Knepley } else if (xp > 0) { /* right */ 82288ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 82388ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 824b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + nVx-1 + nC; /* right bottom back vertex */ 825f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 826f5eeb527SMatthew G Knepley 827b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 828f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 829f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 830f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 831f5eeb527SMatthew G Knepley ++nL; 832b2fff234SMatthew G. Knepley } 83388ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 834b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + nVx-1 + nC; /* right bottom front vertex */ 835f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 836f5eeb527SMatthew G Knepley 837b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 838f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 839f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 840f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 841f5eeb527SMatthew G Knepley ++nL; 842b2fff234SMatthew G. Knepley } 84388ed4aceSMatthew G Knepley } else { 844b2fff234SMatthew G. Knepley nleavesCheck += nVz; 845b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 846b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + nVx-1 + nC; /* right bottom vertices */ 847f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 848f5eeb527SMatthew G Knepley 849b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 850f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 851f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 852f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 853b2fff234SMatthew G. Knepley ++nL; 854b2fff234SMatthew G. Knepley } 855f5eeb527SMatthew G Knepley } 85688ed4aceSMatthew G Knepley } 85788ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 85888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 859b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */ 860f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 861f5eeb527SMatthew G Knepley 862b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 863f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 864f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 865f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 866f5eeb527SMatthew G Knepley ++nL; 867b2fff234SMatthew G. Knepley } 86888ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 869b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */ 870f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*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; 876f5eeb527SMatthew G Knepley ++nL; 877b2fff234SMatthew G. Knepley } 87888ed4aceSMatthew G Knepley } else { 879b2fff234SMatthew G. Knepley for (zv = 0; zv < nVz; ++zv) { 880b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */ 881f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 882f5eeb527SMatthew G Knepley 883b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 884f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 885f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 886f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 887b2fff234SMatthew G. Knepley ++nL; 888b2fff234SMatthew G. Knepley } 889f5eeb527SMatthew G Knepley } 89088ed4aceSMatthew G Knepley } 89188ed4aceSMatthew G Knepley } else { 89288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 893b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 894b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */ 895f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 896f5eeb527SMatthew G Knepley 897b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 898f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 899f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 900f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 901b2fff234SMatthew G. Knepley ++nL; 902b2fff234SMatthew G. Knepley } 903f5eeb527SMatthew G Knepley } 90488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 905b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 906b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */ 907f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 908f5eeb527SMatthew G Knepley 909b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 910f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 911f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 912f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 913b2fff234SMatthew G. Knepley ++nL; 914b2fff234SMatthew G. Knepley } 915f5eeb527SMatthew G Knepley } 91688ed4aceSMatthew G Knepley } else { 917f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 918b2fff234SMatthew G. Knepley for (yv = 0; yv < nVy; ++yv) { 919b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */ 920f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0 + nC; /* TODO: Correct this for neighbor sizes */ 921f5eeb527SMatthew G Knepley 922b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 923f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 924f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 925f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 926b2fff234SMatthew G. Knepley ++nL; 927f5eeb527SMatthew G Knepley } 928f5eeb527SMatthew G Knepley } 929b2fff234SMatthew G. Knepley } 930b2fff234SMatthew G. Knepley #if 0 931b2fff234SMatthew G. Knepley for (xf = 0; xf < nxF; ++xf) { 932f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 933b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* right faces */ 934f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 935f5eeb527SMatthew G Knepley 936b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 937f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 938f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 939f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 940b2fff234SMatthew G. Knepley ++nL; 941f5eeb527SMatthew G Knepley } 94288ed4aceSMatthew G Knepley } 943b2fff234SMatthew G. Knepley #endif 944b2fff234SMatthew G. Knepley } 94588ed4aceSMatthew G Knepley } 94688ed4aceSMatthew G Knepley } else { 94788ed4aceSMatthew G Knepley if (yp < 0) { /* bottom */ 94888ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 949b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 950b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + 0)*nVx + xv + nC; /* bottom back vertices */ 951f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 952f5eeb527SMatthew G Knepley 953b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 954f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 955f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 956f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 957b2fff234SMatthew G. Knepley ++nL; 958b2fff234SMatthew G. Knepley } 959f5eeb527SMatthew G Knepley } 96088ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 961b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 962b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* bottom front vertices */ 963f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 964f5eeb527SMatthew G Knepley 965b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 966f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 967f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 968f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 969b2fff234SMatthew G. Knepley ++nL; 970b2fff234SMatthew G. Knepley } 971f5eeb527SMatthew G Knepley } 97288ed4aceSMatthew G Knepley } else { 973f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 974b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 975b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + 0)*nVx + xv + nC; /* bottom vertices */ 976f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 977f5eeb527SMatthew G Knepley 978b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 979f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 980f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 981f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 982b2fff234SMatthew G. Knepley ++nL; 983f5eeb527SMatthew G Knepley } 984f5eeb527SMatthew G Knepley } 985b2fff234SMatthew G. Knepley } 986b2fff234SMatthew G. Knepley #if 0 987b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 988f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 989b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* bottom faces */ 990f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 991f5eeb527SMatthew G Knepley 992b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 993f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 994f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 995f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 996b2fff234SMatthew G. Knepley ++nL; 997f5eeb527SMatthew G Knepley } 99888ed4aceSMatthew G Knepley } 999b2fff234SMatthew G. Knepley #endif 1000b2fff234SMatthew G. Knepley } 100188ed4aceSMatthew G Knepley } else if (yp > 0) { /* top */ 100288ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 1003b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1004b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */ 1005f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1006f5eeb527SMatthew G Knepley 1007b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1008f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1009f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1010f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1011b2fff234SMatthew G. Knepley ++nL; 1012b2fff234SMatthew G. Knepley } 1013f5eeb527SMatthew G Knepley } 101488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1015b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1016b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */ 1017f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1018f5eeb527SMatthew G Knepley 1019b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1020f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1021f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1022f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1023b2fff234SMatthew G. Knepley ++nL; 1024b2fff234SMatthew G. Knepley } 1025f5eeb527SMatthew G Knepley } 102688ed4aceSMatthew G Knepley } else { 1027f5eeb527SMatthew G Knepley for (zv = 0; zv < nVz; ++zv) { 1028b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1029b2fff234SMatthew G. Knepley const PetscInt localVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */ 1030f5eeb527SMatthew G Knepley const PetscInt remoteVertex = (zv*nVy + 0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1031f5eeb527SMatthew G Knepley 1032b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1033f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1034f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1035f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1036b2fff234SMatthew G. Knepley ++nL; 1037f5eeb527SMatthew G Knepley } 1038f5eeb527SMatthew G Knepley } 1039b2fff234SMatthew G. Knepley } 1040b2fff234SMatthew G. Knepley #if 0 1041b2fff234SMatthew G. Knepley for (yf = 0; yf < nyF; ++yf) { 1042f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1043b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* top faces */ 1044f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1045f5eeb527SMatthew G Knepley 1046b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1047f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1048f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1049f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1050b2fff234SMatthew G. Knepley ++nL; 1051f5eeb527SMatthew G Knepley } 105288ed4aceSMatthew G Knepley } 1053b2fff234SMatthew G. Knepley #endif 1054b2fff234SMatthew G. Knepley } 105588ed4aceSMatthew G Knepley } else { 105688ed4aceSMatthew G Knepley if (zp < 0) { /* back */ 1057f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1058b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1059b2fff234SMatthew G. Knepley const PetscInt localVertex = ( 0*nVy + yv)*nVx + xv + nC; /* back vertices */ 1060f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1061f5eeb527SMatthew G Knepley 1062b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1063f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1064f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1065f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1066b2fff234SMatthew G. Knepley ++nL; 1067f5eeb527SMatthew G Knepley } 1068f5eeb527SMatthew G Knepley } 1069b2fff234SMatthew G. Knepley } 1070b2fff234SMatthew G. Knepley #if 0 1071b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1072f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1073b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* back faces */ 1074f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1075f5eeb527SMatthew G Knepley 1076b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1077f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1078f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1079f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1080b2fff234SMatthew G. Knepley ++nL; 1081f5eeb527SMatthew G Knepley } 1082b2fff234SMatthew G. Knepley } 1083b2fff234SMatthew G. Knepley #endif 108488ed4aceSMatthew G Knepley } else if (zp > 0) { /* front */ 1085f5eeb527SMatthew G Knepley for (yv = 0; yv < nVy; ++yv) { 1086b2fff234SMatthew G. Knepley for (xv = 0; xv < nVx; ++xv) { 1087b2fff234SMatthew G. Knepley const PetscInt localVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */ 1088f5eeb527SMatthew G Knepley const PetscInt remoteVertex = ( 0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */ 1089f5eeb527SMatthew G Knepley 1090b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localVertex)) { 1091f5eeb527SMatthew G Knepley localPoints[nL] = localVertex; 1092f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1093f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteVertex; 1094b2fff234SMatthew G. Knepley ++nL; 1095f5eeb527SMatthew G Knepley } 1096f5eeb527SMatthew G Knepley } 1097b2fff234SMatthew G. Knepley } 1098b2fff234SMatthew G. Knepley #if 0 1099b2fff234SMatthew G. Knepley for (zf = 0; zf < nzF; ++zf) { 1100f5eeb527SMatthew G Knepley /* THIS IS WRONG */ 1101b2fff234SMatthew G. Knepley const PetscInt localFace = 0 + nC+nV; /* front faces */ 1102f5eeb527SMatthew G Knepley const PetscInt remoteFace = 0 + nC+nV; 1103f5eeb527SMatthew G Knepley 1104b2fff234SMatthew G. Knepley if (!PetscBTLookupSet(isLeaf, localFace)) { 1105f5eeb527SMatthew G Knepley localPoints[nL] = localFace; 1106f5eeb527SMatthew G Knepley remotePoints[nL].rank = neighbor; 1107f5eeb527SMatthew G Knepley remotePoints[nL].index = remoteFace; 1108b2fff234SMatthew G. Knepley ++nL; 1109f5eeb527SMatthew G Knepley } 1110b2fff234SMatthew G. Knepley } 1111b2fff234SMatthew G. Knepley #endif 111288ed4aceSMatthew G Knepley } else { 111388ed4aceSMatthew G Knepley /* Nothing is shared from the interior */ 111488ed4aceSMatthew G Knepley } 111588ed4aceSMatthew G Knepley } 111688ed4aceSMatthew G Knepley } 111788ed4aceSMatthew G Knepley } 111888ed4aceSMatthew G Knepley } 111988ed4aceSMatthew G Knepley } 112088ed4aceSMatthew G Knepley } 1121b2fff234SMatthew G. Knepley ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr); 1122b2fff234SMatthew G. Knepley /* Remove duplication in leaf determination */ 1123b2fff234SMatthew 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); 112482f516ccSBarry Smith ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr); 112588ed4aceSMatthew G Knepley ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1126a4b60ecfSMatthew G. Knepley ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr); 112788ed4aceSMatthew G Knepley ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1128affa55c7SMatthew G. Knepley *s = section; 1129affa55c7SMatthew G. Knepley PetscFunctionReturn(0); 1130affa55c7SMatthew G. Knepley } 1131affa55c7SMatthew G. Knepley 11323385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) 11333385ff7eSMatthew G. Knepley { 11343385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 11353385ff7eSMatthew G. Knepley Vec coordinates; 11363385ff7eSMatthew G. Knepley PetscSection section; 11373385ff7eSMatthew G. Knepley PetscScalar *coords; 11383385ff7eSMatthew G. Knepley PetscReal h[3]; 11393385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 11403385ff7eSMatthew G. Knepley PetscErrorCode ierr; 11413385ff7eSMatthew G. Knepley 11423385ff7eSMatthew G. Knepley PetscFunctionBegin; 11433385ff7eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11443385ff7eSMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 11454a2f8832SBarry Smith if (dim > 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"The following code only works for dim <= 3"); 11463385ff7eSMatthew G. Knepley h[0] = (xu - xl)/M; 11473385ff7eSMatthew G. Knepley h[1] = (yu - yl)/N; 11483385ff7eSMatthew G. Knepley h[2] = (zu - zl)/P; 11493385ff7eSMatthew G. Knepley ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 11503385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 11513385ff7eSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 11523385ff7eSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr); 11533385ff7eSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr); 11543385ff7eSMatthew G. Knepley ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr); 11553385ff7eSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 11563385ff7eSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr); 11573385ff7eSMatthew G. Knepley } 11583385ff7eSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 11593385ff7eSMatthew G. Knepley ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 11603385ff7eSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr); 116124640c55SToby Isaac ierr = PetscObjectSetName((PetscObject)coordinates,"coordinates");CHKERRQ(ierr); 11623385ff7eSMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 11633385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) { 1164e4d69e08SBarry Smith PetscInt ind[3], d, off; 11653385ff7eSMatthew G. Knepley 1166e4d69e08SBarry Smith ind[0] = 0; 1167e4d69e08SBarry Smith ind[1] = 0; 1168e4d69e08SBarry Smith ind[2] = k + da->zs; 11693385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) { 11703385ff7eSMatthew G. Knepley ind[1] = j + da->ys; 11713385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) { 11723385ff7eSMatthew G. Knepley const PetscInt vertex = (k*nVy + j)*nVx + i + vStart; 11733385ff7eSMatthew G. Knepley 11743385ff7eSMatthew G. Knepley ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr); 11753385ff7eSMatthew G. Knepley ind[0] = i + da->xs; 11763385ff7eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 11773385ff7eSMatthew G. Knepley coords[off+d] = h[d]*ind[d]; 11783385ff7eSMatthew G. Knepley } 11793385ff7eSMatthew G. Knepley } 11803385ff7eSMatthew G. Knepley } 11813385ff7eSMatthew G. Knepley } 11823385ff7eSMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 118346e270d4SMatthew G. Knepley ierr = DMSetCoordinateSection(dm, PETSC_DETERMINE, section);CHKERRQ(ierr); 11843385ff7eSMatthew G. Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 1185a4b60ecfSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 11863385ff7eSMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 11873385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 11883385ff7eSMatthew G. Knepley } 11899a800dd8SMatthew G. Knepley 1190c970b36aSToby Isaac PetscErrorCode DMProjectFunctionLocal_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 11919a800dd8SMatthew G. Knepley { 11922764a2aaSMatthew G. Knepley PetscDS prob; 1193b7a4d31fSMatthew G. Knepley PetscFE fe; 1194b7a4d31fSMatthew G. Knepley PetscDualSpace sp; 11959a800dd8SMatthew G. Knepley PetscQuadrature q; 11969a800dd8SMatthew G. Knepley PetscSection section; 11979a800dd8SMatthew G. Knepley PetscScalar *values; 11989c3cf19fSMatthew G. Knepley PetscInt numFields, Nc, numPoints, dim, dimEmbed, spDim, totDim, numValues, cStart, cEnd, f, c, v, d; 11999a800dd8SMatthew G. Knepley PetscErrorCode ierr; 12009a800dd8SMatthew G. Knepley 12019a800dd8SMatthew G. Knepley PetscFunctionBegin; 12022764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 12032764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 12042764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1205b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 12069a800dd8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 12079a800dd8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 12089a800dd8SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 12097a1a1bd4SToby Isaac ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 12109a800dd8SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 12119a800dd8SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 12129a800dd8SMatthew 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); 1213*69291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 12149c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, NULL, &numPoints, NULL, NULL);CHKERRQ(ierr); 12159a800dd8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1216de2092bfSMatthew G. Knepley PetscFECellGeom geom; 12179a800dd8SMatthew G. Knepley 1218de2092bfSMatthew G. Knepley ierr = DMDAComputeCellGeometryFEM(dm, c, q, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 12197a1a1bd4SToby Isaac geom.dim = dim; 12207a1a1bd4SToby Isaac geom.dimEmbed = dimEmbed; 12219a800dd8SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 1222c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[f] : NULL; 1223b7a4d31fSMatthew G. Knepley 12242764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 12259c3cf19fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1226b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 1227b7a4d31fSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 12289a800dd8SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 12299c3cf19fSMatthew G. Knepley ierr = PetscDualSpaceApply(sp, d, time, &geom, Nc, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 12309a800dd8SMatthew G. Knepley } 12319a800dd8SMatthew G. Knepley } 12329a800dd8SMatthew G. Knepley ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 12339a800dd8SMatthew G. Knepley } 1234*69291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 12359a800dd8SMatthew G. Knepley PetscFunctionReturn(0); 12369a800dd8SMatthew G. Knepley } 12379a800dd8SMatthew G. Knepley 1238c970b36aSToby Isaac PetscErrorCode DMComputeL2Diff_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 12399a800dd8SMatthew G. Knepley { 12409a800dd8SMatthew G. Knepley const PetscInt debug = 0; 12412764a2aaSMatthew G. Knepley PetscDS prob; 1242b7a4d31fSMatthew G. Knepley PetscFE fe; 12439a800dd8SMatthew G. Knepley PetscQuadrature quad; 1244b7a4d31fSMatthew G. Knepley PetscSection section; 12459a800dd8SMatthew G. Knepley Vec localX; 12469a800dd8SMatthew G. Knepley PetscScalar *funcVal; 12479a800dd8SMatthew G. Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 12489a800dd8SMatthew G. Knepley PetscReal localDiff = 0.0; 1249b7a4d31fSMatthew G. Knepley PetscInt dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp; 12509a800dd8SMatthew G. Knepley PetscErrorCode ierr; 12519a800dd8SMatthew G. Knepley 12529a800dd8SMatthew G. Knepley PetscFunctionBegin; 12539a800dd8SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 12542764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 12552764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 12562764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1257b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 12589a800dd8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 12599a800dd8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 12609a800dd8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 12619a800dd8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 12629a800dd8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 12639a800dd8SMatthew G. Knepley /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 1264b7a4d31fSMatthew G. Knepley ierr = PetscMalloc5(Nc,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 12659a800dd8SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 12669a800dd8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 12679a800dd8SMatthew G. Knepley PetscScalar *x = NULL; 12689a800dd8SMatthew G. Knepley PetscReal elemDiff = 0.0; 12699a800dd8SMatthew G. Knepley 12708e0841e0SMatthew G. Knepley ierr = DMDAComputeCellGeometryFEM(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr); 12719a800dd8SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 12729a800dd8SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 12739a800dd8SMatthew G. Knepley 12749a800dd8SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 1275c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[field] : NULL; 127621454ff5SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 12779a800dd8SMatthew G. Knepley PetscReal *basis; 12789c3cf19fSMatthew G. Knepley PetscInt qNc, Nq, Nb, Nc, q, d, e, fc, f; 12799a800dd8SMatthew G. Knepley 12802764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 12819c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 12829c3cf19fSMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 12839c3cf19fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 12849c3cf19fSMatthew G. Knepley if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Quadrature components %D != %D field components\n", qNc, Nc); 1285b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr); 12869a800dd8SMatthew G. Knepley if (debug) { 12879a800dd8SMatthew G. Knepley char title[1024]; 12889a800dd8SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 12899c3cf19fSMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr); 12909a800dd8SMatthew G. Knepley } 12919c3cf19fSMatthew G. Knepley for (q = 0; q < Nq; ++q) { 12929a800dd8SMatthew G. Knepley for (d = 0; d < dim; d++) { 12939a800dd8SMatthew G. Knepley coords[d] = v0[d]; 12949a800dd8SMatthew G. Knepley for (e = 0; e < dim; e++) { 12959a800dd8SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 12969a800dd8SMatthew G. Knepley } 12979a800dd8SMatthew G. Knepley } 1298a98e6df5SMatthew G. Knepley ierr = (*funcs[field])(dim, time, coords, numFields, funcVal, ctx);CHKERRQ(ierr); 12999c3cf19fSMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 13009a800dd8SMatthew G. Knepley PetscScalar interpolant = 0.0; 13019a800dd8SMatthew G. Knepley 13029c3cf19fSMatthew G. Knepley for (f = 0; f < Nb; ++f) { 13039c3cf19fSMatthew G. Knepley interpolant += x[fieldOffset+f]*basis[(q*Nb+f)*Nc+fc]; 13049a800dd8SMatthew G. Knepley } 13059c3cf19fSMatthew 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*Nc+fc]*detJ);CHKERRQ(ierr);} 13069c3cf19fSMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q*Nc+c]*detJ; 13079a800dd8SMatthew G. Knepley } 13089a800dd8SMatthew G. Knepley } 13099c3cf19fSMatthew G. Knepley comp += Nc; 13109c3cf19fSMatthew G. Knepley fieldOffset += Nb; 13119a800dd8SMatthew G. Knepley } 13129a800dd8SMatthew G. Knepley ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 13139a800dd8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 13149a800dd8SMatthew G. Knepley localDiff += elemDiff; 13159a800dd8SMatthew G. Knepley } 13169a800dd8SMatthew G. Knepley ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 13179a800dd8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1318b2566f29SBarry Smith ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 13199a800dd8SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 1320a66d4d66SMatthew G Knepley PetscFunctionReturn(0); 1321a66d4d66SMatthew G Knepley } 1322a66d4d66SMatthew G Knepley 1323b698f381SToby Isaac PetscErrorCode DMComputeL2GradientDiff_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 132423d86601SMatthew G. Knepley { 132523d86601SMatthew G. Knepley const PetscInt debug = 0; 13262764a2aaSMatthew G. Knepley PetscDS prob; 1327b7a4d31fSMatthew G. Knepley PetscFE fe; 132823d86601SMatthew G. Knepley PetscQuadrature quad; 1329b7a4d31fSMatthew G. Knepley PetscSection section; 133023d86601SMatthew G. Knepley Vec localX; 133123d86601SMatthew G. Knepley PetscScalar *funcVal, *interpolantVec; 133223d86601SMatthew G. Knepley PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 133323d86601SMatthew G. Knepley PetscReal localDiff = 0.0; 1334b7a4d31fSMatthew G. Knepley PetscInt dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp; 133523d86601SMatthew G. Knepley PetscErrorCode ierr; 133623d86601SMatthew G. Knepley 133723d86601SMatthew G. Knepley PetscFunctionBegin; 133823d86601SMatthew G. Knepley ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 13392764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 13402764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr); 13412764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr); 1342b7a4d31fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 134323d86601SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 134423d86601SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 134523d86601SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 134623d86601SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 134723d86601SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 134823d86601SMatthew G. Knepley /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 1349b7a4d31fSMatthew G. Knepley ierr = PetscMalloc7(Nc,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 135023d86601SMatthew G. Knepley ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 135123d86601SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 135223d86601SMatthew G. Knepley PetscScalar *x = NULL; 135323d86601SMatthew G. Knepley PetscReal elemDiff = 0.0; 135423d86601SMatthew G. Knepley 13558e0841e0SMatthew G. Knepley ierr = DMDAComputeCellGeometryFEM(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr); 135623d86601SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 135723d86601SMatthew G. Knepley ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 135823d86601SMatthew G. Knepley 135923d86601SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 136023d86601SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 136123d86601SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 136223d86601SMatthew G. Knepley PetscReal *basisDer; 13639c3cf19fSMatthew G. Knepley PetscInt qNc, Nq, Nb, Nc, q, d, e, fc, f, g; 136423d86601SMatthew G. Knepley 13652764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr); 13669c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1367b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 1368b7a4d31fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 13699c3cf19fSMatthew G. Knepley if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Quadrature components %D != %D field components\n", qNc, Nc); 1370b7a4d31fSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 137123d86601SMatthew G. Knepley if (debug) { 137223d86601SMatthew G. Knepley char title[1024]; 137323d86601SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 137423d86601SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 137523d86601SMatthew G. Knepley } 137623d86601SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 137723d86601SMatthew G. Knepley for (d = 0; d < dim; d++) { 137823d86601SMatthew G. Knepley coords[d] = v0[d]; 137923d86601SMatthew G. Knepley for (e = 0; e < dim; e++) { 138023d86601SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 138123d86601SMatthew G. Knepley } 138223d86601SMatthew G. Knepley } 1383a98e6df5SMatthew G. Knepley ierr = (*funcs[field])(dim, time, coords, n, numFields, funcVal, ctx);CHKERRQ(ierr); 138423d86601SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 138523d86601SMatthew G. Knepley PetscScalar interpolant = 0.0; 138623d86601SMatthew G. Knepley 138723d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 138823d86601SMatthew G. Knepley for (f = 0; f < Nb; ++f) { 138923d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) { 139023d86601SMatthew G. Knepley realSpaceDer[d] = 0.0; 139123d86601SMatthew G. Knepley for (g = 0; g < dim; ++g) { 13929c3cf19fSMatthew G. Knepley realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb+f)*dim+g]; 139323d86601SMatthew G. Knepley } 13949c3cf19fSMatthew G. Knepley interpolantVec[d] += x[fieldOffset+f]*realSpaceDer[d]; 139523d86601SMatthew G. Knepley } 139623d86601SMatthew G. Knepley } 139723d86601SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 13989c3cf19fSMatthew 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*Nc+c]*detJ);CHKERRQ(ierr);} 13999c3cf19fSMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q*Nc+c]*detJ; 140023d86601SMatthew G. Knepley } 140123d86601SMatthew G. Knepley } 140223d86601SMatthew G. Knepley comp += Nc; 140323d86601SMatthew G. Knepley fieldOffset += Nb*Nc; 140423d86601SMatthew G. Knepley } 140523d86601SMatthew G. Knepley ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 140623d86601SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 140723d86601SMatthew G. Knepley localDiff += elemDiff; 140823d86601SMatthew G. Knepley } 140923d86601SMatthew G. Knepley ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 141023d86601SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1411b2566f29SBarry Smith ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 141223d86601SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 141323d86601SMatthew G. Knepley PetscFunctionReturn(0); 141423d86601SMatthew G. Knepley } 141523d86601SMatthew G. Knepley 141647c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 141747c6ae99SBarry Smith 141847c6ae99SBarry Smith /*@C 1419aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 142047c6ae99SBarry Smith 142147c6ae99SBarry Smith Input Parameter: 142247c6ae99SBarry Smith + da - information about my local patch 142347c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 142447c6ae99SBarry Smith 142547c6ae99SBarry Smith Output Parameters: 142647c6ae99SBarry Smith . vptr - array data structured 142747c6ae99SBarry Smith 142847c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 142947c6ae99SBarry Smith to zero them. 143047c6ae99SBarry Smith 143147c6ae99SBarry Smith Level: advanced 143247c6ae99SBarry Smith 1433bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 143447c6ae99SBarry Smith 143547c6ae99SBarry Smith @*/ 14367087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 143747c6ae99SBarry Smith { 143847c6ae99SBarry Smith PetscErrorCode ierr; 143947c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 144047c6ae99SBarry Smith char *iarray_start; 144147c6ae99SBarry Smith void **iptr = (void**)vptr; 144247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 144347c6ae99SBarry Smith 144447c6ae99SBarry Smith PetscFunctionBegin; 144547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 144647c6ae99SBarry Smith if (ghosted) { 1447aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 144847c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 144947c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 145047c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 14510298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 14520298fd71SBarry Smith dd->startghostedin[i] = NULL; 145347c6ae99SBarry Smith 145447c6ae99SBarry Smith goto done; 145547c6ae99SBarry Smith } 145647c6ae99SBarry Smith } 145747c6ae99SBarry Smith xs = dd->Xs; 145847c6ae99SBarry Smith ys = dd->Ys; 145947c6ae99SBarry Smith zs = dd->Zs; 146047c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 146147c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 146247c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 146347c6ae99SBarry Smith } else { 1464aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 146547c6ae99SBarry Smith if (dd->arrayin[i]) { 146647c6ae99SBarry Smith *iptr = dd->arrayin[i]; 146747c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 14680298fd71SBarry Smith dd->arrayin[i] = NULL; 14690298fd71SBarry Smith dd->startin[i] = NULL; 147047c6ae99SBarry Smith 147147c6ae99SBarry Smith goto done; 147247c6ae99SBarry Smith } 147347c6ae99SBarry Smith } 147447c6ae99SBarry Smith xs = dd->xs; 147547c6ae99SBarry Smith ys = dd->ys; 147647c6ae99SBarry Smith zs = dd->zs; 147747c6ae99SBarry Smith xm = dd->xe-dd->xs; 147847c6ae99SBarry Smith ym = dd->ye-dd->ys; 147947c6ae99SBarry Smith zm = dd->ze-dd->zs; 148047c6ae99SBarry Smith } 148147c6ae99SBarry Smith 1482c73cfb54SMatthew G. Knepley switch (da->dim) { 148347c6ae99SBarry Smith case 1: { 148447c6ae99SBarry Smith void *ptr; 148547c6ae99SBarry Smith 1486901f1932SJed Brown ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 148747c6ae99SBarry Smith 148847c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 148947c6ae99SBarry Smith *iptr = (void*)ptr; 14908865f1eaSKarl Rupp break; 14918865f1eaSKarl Rupp } 149247c6ae99SBarry Smith case 2: { 149347c6ae99SBarry Smith void **ptr; 149447c6ae99SBarry Smith 1495901f1932SJed Brown ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 149647c6ae99SBarry Smith 149747c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 14988865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 149947c6ae99SBarry Smith *iptr = (void*)ptr; 15008865f1eaSKarl Rupp break; 15018865f1eaSKarl Rupp } 150247c6ae99SBarry Smith case 3: { 150347c6ae99SBarry Smith void ***ptr,**bptr; 150447c6ae99SBarry Smith 1505901f1932SJed Brown ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 150647c6ae99SBarry Smith 150747c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 150847c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 15098865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 151047c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 151147c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 151247c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 151347c6ae99SBarry Smith } 151447c6ae99SBarry Smith } 151547c6ae99SBarry Smith *iptr = (void*)ptr; 15168865f1eaSKarl Rupp break; 15178865f1eaSKarl Rupp } 151847c6ae99SBarry Smith default: 1519c73cfb54SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",da->dim); 152047c6ae99SBarry Smith } 152147c6ae99SBarry Smith 152247c6ae99SBarry Smith done: 152347c6ae99SBarry Smith /* add arrays to the checked out list */ 152447c6ae99SBarry Smith if (ghosted) { 1525aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 152647c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 152747c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 152847c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 152947c6ae99SBarry Smith break; 153047c6ae99SBarry Smith } 153147c6ae99SBarry Smith } 153247c6ae99SBarry Smith } else { 1533aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 153447c6ae99SBarry Smith if (!dd->arrayout[i]) { 153547c6ae99SBarry Smith dd->arrayout[i] = *iptr; 153647c6ae99SBarry Smith dd->startout[i] = iarray_start; 153747c6ae99SBarry Smith break; 153847c6ae99SBarry Smith } 153947c6ae99SBarry Smith } 154047c6ae99SBarry Smith } 154147c6ae99SBarry Smith PetscFunctionReturn(0); 154247c6ae99SBarry Smith } 154347c6ae99SBarry Smith 154447c6ae99SBarry Smith /*@C 1545aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 154647c6ae99SBarry Smith 154747c6ae99SBarry Smith Input Parameter: 154847c6ae99SBarry Smith + da - information about my local patch 154947c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 155047c6ae99SBarry Smith - vptr - array data structured to be passed to ad_FormFunctionLocal() 155147c6ae99SBarry Smith 155247c6ae99SBarry Smith Level: advanced 155347c6ae99SBarry Smith 1554bcaeba4dSBarry Smith .seealso: DMDAGetArray() 155547c6ae99SBarry Smith 155647c6ae99SBarry Smith @*/ 15577087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 155847c6ae99SBarry Smith { 155947c6ae99SBarry Smith PetscInt i; 156047c6ae99SBarry Smith void **iptr = (void**)vptr,*iarray_start = 0; 156147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 156247c6ae99SBarry Smith 156347c6ae99SBarry Smith PetscFunctionBegin; 156447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 156547c6ae99SBarry Smith if (ghosted) { 1566aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 156747c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 156847c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 15690298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 15700298fd71SBarry Smith dd->startghostedout[i] = NULL; 157147c6ae99SBarry Smith break; 157247c6ae99SBarry Smith } 157347c6ae99SBarry Smith } 1574aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 157547c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 157647c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 157747c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 157847c6ae99SBarry Smith break; 157947c6ae99SBarry Smith } 158047c6ae99SBarry Smith } 158147c6ae99SBarry Smith } else { 1582aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 158347c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 158447c6ae99SBarry Smith iarray_start = dd->startout[i]; 15850298fd71SBarry Smith dd->arrayout[i] = NULL; 15860298fd71SBarry Smith dd->startout[i] = NULL; 158747c6ae99SBarry Smith break; 158847c6ae99SBarry Smith } 158947c6ae99SBarry Smith } 1590aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 159147c6ae99SBarry Smith if (!dd->arrayin[i]) { 159247c6ae99SBarry Smith dd->arrayin[i] = *iptr; 159347c6ae99SBarry Smith dd->startin[i] = iarray_start; 159447c6ae99SBarry Smith break; 159547c6ae99SBarry Smith } 159647c6ae99SBarry Smith } 159747c6ae99SBarry Smith } 159847c6ae99SBarry Smith PetscFunctionReturn(0); 159947c6ae99SBarry Smith } 160047c6ae99SBarry Smith 1601