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); 29*7a8be351SBarry Smith PetscCheck(da,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 38580bdb30SBarry Smith ierr = PetscArraycpy(mxGetPr(mat),array,n*m);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 477087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 4847c6ae99SBarry Smith { 4947c6ae99SBarry Smith PetscErrorCode ierr; 5047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith PetscFunctionBegin; 5347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5447c6ae99SBarry Smith PetscValidPointer(g,2); 5547c6ae99SBarry Smith ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr); 5647c6ae99SBarry Smith ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr); 5747c6ae99SBarry Smith ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr); 58401ddaa8SBarry Smith ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr); 5974427ab1SRichard Tran Mills if (dd->nlocal < da->bind_below) { 6074427ab1SRichard Tran Mills ierr = VecSetBindingPropagates(*g,PETSC_TRUE);CHKERRQ(ierr); 6174427ab1SRichard Tran Mills ierr = VecBindToCPU(*g,PETSC_TRUE);CHKERRQ(ierr); 6274427ab1SRichard Tran Mills } 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 6947c6ae99SBarry Smith PetscFunctionReturn(0); 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith 72939f6067SMatthew G. Knepley /*@ 73939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 74939f6067SMatthew G. Knepley 75939f6067SMatthew G. Knepley Input Parameter: 76939f6067SMatthew G. Knepley . dm - The DM object 77939f6067SMatthew G. Knepley 78939f6067SMatthew G. Knepley Output Parameters: 79939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 80939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 81939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 82939f6067SMatthew G. Knepley - numCells - The number of local cells 83939f6067SMatthew G. Knepley 84939f6067SMatthew G. Knepley Level: developer 85939f6067SMatthew G. Knepley 86939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint() 87939f6067SMatthew G. Knepley @*/ 883582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 8957459e9aSMatthew G Knepley { 9057459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 91c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 9257459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 9357459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 9457459e9aSMatthew G Knepley 9557459e9aSMatthew G Knepley PetscFunctionBegin; 96a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 973582350cSMatthew G. Knepley if (numCellsX) { 983582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 993582350cSMatthew G. Knepley *numCellsX = mx; 1003582350cSMatthew G. Knepley } 1013582350cSMatthew G. Knepley if (numCellsY) { 1028135c375SStefano Zampini PetscValidIntPointer(numCellsY,3); 1033582350cSMatthew G. Knepley *numCellsY = my; 1043582350cSMatthew G. Knepley } 1053582350cSMatthew G. Knepley if (numCellsZ) { 1068135c375SStefano Zampini PetscValidIntPointer(numCellsZ,4); 1073582350cSMatthew G. Knepley *numCellsZ = mz; 1083582350cSMatthew G. Knepley } 10957459e9aSMatthew G Knepley if (numCells) { 1103582350cSMatthew G. Knepley PetscValidIntPointer(numCells,5); 11157459e9aSMatthew G Knepley *numCells = nC; 11257459e9aSMatthew G Knepley } 11357459e9aSMatthew G Knepley PetscFunctionReturn(0); 11457459e9aSMatthew G Knepley } 11557459e9aSMatthew G Knepley 116d6401b93SMatthew G. Knepley /*@ 117d6401b93SMatthew G. Knepley DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 118d6401b93SMatthew G. Knepley 119d6401b93SMatthew G. Knepley Input Parameters: 120d6401b93SMatthew G. Knepley + dm - The DM object 121d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell 122d6401b93SMatthew G. Knepley 123d6401b93SMatthew G. Knepley Output Parameters: 124d6401b93SMatthew G. Knepley . point - The local DM point 125d6401b93SMatthew G. Knepley 126d6401b93SMatthew G. Knepley Level: developer 127d6401b93SMatthew G. Knepley 128d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells() 129d6401b93SMatthew G. Knepley @*/ 130d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) 131d6401b93SMatthew G. Knepley { 132c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 133d6401b93SMatthew G. Knepley DMDALocalInfo info; 134d6401b93SMatthew G. Knepley PetscErrorCode ierr; 135d6401b93SMatthew G. Knepley 136d6401b93SMatthew G. Knepley PetscFunctionBegin; 137a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 138d6401b93SMatthew G. Knepley PetscValidIntPointer(point,5); 139d6401b93SMatthew G. Knepley ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 1402c71b3e2SJacob Faibussowitsch if (dim > 0) PetscCheckFalse((i < info.gxs) || (i >= info.gxs+info.gxm),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %d not in [%d, %d)", i, info.gxs, info.gxs+info.gxm); 1412c71b3e2SJacob Faibussowitsch if (dim > 1) PetscCheckFalse((j < info.gys) || (j >= info.gys+info.gym),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %d not in [%d, %d)", j, info.gys, info.gys+info.gym); 1422c71b3e2SJacob Faibussowitsch if (dim > 2) PetscCheckFalse((k < info.gzs) || (k >= info.gzs+info.gzm),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %d not in [%d, %d)", k, info.gzs, info.gzs+info.gzm); 143d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0); 144d6401b93SMatthew G. Knepley PetscFunctionReturn(0); 145d6401b93SMatthew G. Knepley } 146d6401b93SMatthew G. Knepley 14757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 14857459e9aSMatthew G Knepley { 14957459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 150c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 15157459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 15257459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 15357459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 15457459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 15557459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 15657459e9aSMatthew G Knepley 15757459e9aSMatthew G Knepley PetscFunctionBegin; 15857459e9aSMatthew G Knepley if (numVerticesX) { 15957459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 16057459e9aSMatthew G Knepley *numVerticesX = nVx; 16157459e9aSMatthew G Knepley } 16257459e9aSMatthew G Knepley if (numVerticesY) { 16357459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 16457459e9aSMatthew G Knepley *numVerticesY = nVy; 16557459e9aSMatthew G Knepley } 16657459e9aSMatthew G Knepley if (numVerticesZ) { 16757459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 16857459e9aSMatthew G Knepley *numVerticesZ = nVz; 16957459e9aSMatthew G Knepley } 17057459e9aSMatthew G Knepley if (numVertices) { 17157459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 17257459e9aSMatthew G Knepley *numVertices = nV; 17357459e9aSMatthew G Knepley } 17457459e9aSMatthew G Knepley PetscFunctionReturn(0); 17557459e9aSMatthew G Knepley } 17657459e9aSMatthew G Knepley 17757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 17857459e9aSMatthew G Knepley { 17957459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 180c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 18157459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 18257459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 18357459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 18457459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 18557459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 18657459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 18757459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 18857459e9aSMatthew G Knepley 18957459e9aSMatthew G Knepley PetscFunctionBegin; 19057459e9aSMatthew G Knepley if (numXFacesX) { 19157459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 19257459e9aSMatthew G Knepley *numXFacesX = nxF; 19357459e9aSMatthew G Knepley } 19457459e9aSMatthew G Knepley if (numXFaces) { 19557459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 19657459e9aSMatthew G Knepley *numXFaces = nXF; 19757459e9aSMatthew G Knepley } 19857459e9aSMatthew G Knepley if (numYFacesY) { 19957459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 20057459e9aSMatthew G Knepley *numYFacesY = nyF; 20157459e9aSMatthew G Knepley } 20257459e9aSMatthew G Knepley if (numYFaces) { 20357459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 20457459e9aSMatthew G Knepley *numYFaces = nYF; 20557459e9aSMatthew G Knepley } 20657459e9aSMatthew G Knepley if (numZFacesZ) { 20757459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 20857459e9aSMatthew G Knepley *numZFacesZ = nzF; 20957459e9aSMatthew G Knepley } 21057459e9aSMatthew G Knepley if (numZFaces) { 21157459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 21257459e9aSMatthew G Knepley *numZFaces = nZF; 21357459e9aSMatthew G Knepley } 21457459e9aSMatthew G Knepley PetscFunctionReturn(0); 21557459e9aSMatthew G Knepley } 21657459e9aSMatthew G Knepley 21757459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 21857459e9aSMatthew G Knepley { 219c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 22057459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 22157459e9aSMatthew G Knepley PetscErrorCode ierr; 22257459e9aSMatthew G Knepley 22357459e9aSMatthew G Knepley PetscFunctionBegin; 2248865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 2258865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 2263582350cSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2270298fd71SBarry Smith ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2280298fd71SBarry Smith ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 22957459e9aSMatthew G Knepley if (height == 0) { 23057459e9aSMatthew G Knepley /* Cells */ 2318865f1eaSKarl Rupp if (pStart) *pStart = 0; 2328865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 23357459e9aSMatthew G Knepley } else if (height == 1) { 23457459e9aSMatthew G Knepley /* Faces */ 2358865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 2368865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 23757459e9aSMatthew G Knepley } else if (height == dim) { 23857459e9aSMatthew G Knepley /* Vertices */ 2398865f1eaSKarl Rupp if (pStart) *pStart = nC; 2408865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 24157459e9aSMatthew G Knepley } else if (height < 0) { 24257459e9aSMatthew G Knepley /* All points */ 2438865f1eaSKarl Rupp if (pStart) *pStart = 0; 2448865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 24598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 24657459e9aSMatthew G Knepley PetscFunctionReturn(0); 24757459e9aSMatthew G Knepley } 24857459e9aSMatthew G Knepley 2493385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) 2503385ff7eSMatthew G. Knepley { 251c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2523385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2533385ff7eSMatthew G. Knepley PetscErrorCode ierr; 2543385ff7eSMatthew G. Knepley 2553385ff7eSMatthew G. Knepley PetscFunctionBegin; 2563385ff7eSMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart,3); 2573385ff7eSMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd,4); 2583385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2593385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2603385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 2613385ff7eSMatthew G. Knepley if (depth == dim) { 2623385ff7eSMatthew G. Knepley /* Cells */ 2633385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2643385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC; 2653385ff7eSMatthew G. Knepley } else if (depth == dim-1) { 2663385ff7eSMatthew G. Knepley /* Faces */ 2673385ff7eSMatthew G. Knepley if (pStart) *pStart = nC+nV; 2683385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2693385ff7eSMatthew G. Knepley } else if (depth == 0) { 2703385ff7eSMatthew G. Knepley /* Vertices */ 2713385ff7eSMatthew G. Knepley if (pStart) *pStart = nC; 2723385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV; 2733385ff7eSMatthew G. Knepley } else if (depth < 0) { 2743385ff7eSMatthew G. Knepley /* All points */ 2753385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2763385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 27798921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth); 2783385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2793385ff7eSMatthew G. Knepley } 2803385ff7eSMatthew G. Knepley 2813385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) 2823385ff7eSMatthew G. Knepley { 283c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2843385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2853385ff7eSMatthew G. Knepley PetscErrorCode ierr; 2863385ff7eSMatthew G. Knepley 2873385ff7eSMatthew G. Knepley PetscFunctionBegin; 2883385ff7eSMatthew G. Knepley *coneSize = 0; 2893385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr); 2903385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr); 2913385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr); 2923385ff7eSMatthew G. Knepley switch (dim) { 2933385ff7eSMatthew G. Knepley case 2: 2943385ff7eSMatthew G. Knepley if (p >= 0) { 2953385ff7eSMatthew G. Knepley if (p < nC) { 2963385ff7eSMatthew G. Knepley *coneSize = 4; 2973385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 2983385ff7eSMatthew G. Knepley *coneSize = 0; 2993385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF+nZF) { 3003385ff7eSMatthew G. Knepley *coneSize = 2; 30198921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 30298921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3033385ff7eSMatthew G. Knepley break; 3043385ff7eSMatthew G. Knepley case 3: 3053385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3063385ff7eSMatthew G. Knepley } 3073385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3083385ff7eSMatthew G. Knepley } 3093385ff7eSMatthew G. Knepley 3103385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) 3113385ff7eSMatthew G. Knepley { 312c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 3133385ff7eSMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 3143385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3153385ff7eSMatthew G. Knepley 3163385ff7eSMatthew G. Knepley PetscFunctionBegin; 31769291d52SBarry Smith if (!cone) {ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr);} 3183385ff7eSMatthew G. Knepley ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr); 3193385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 3203385ff7eSMatthew G. Knepley ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 3213385ff7eSMatthew G. Knepley switch (dim) { 3223385ff7eSMatthew G. Knepley case 2: 3233385ff7eSMatthew G. Knepley if (p >= 0) { 3243385ff7eSMatthew G. Knepley if (p < nC) { 3253385ff7eSMatthew G. Knepley const PetscInt cy = p / nCx; 3263385ff7eSMatthew G. Knepley const PetscInt cx = p % nCx; 3273385ff7eSMatthew G. Knepley 3283385ff7eSMatthew G. Knepley (*cone)[0] = cy*nxF + cx + nC+nV; 3293385ff7eSMatthew G. Knepley (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF; 3303385ff7eSMatthew G. Knepley (*cone)[2] = cy*nxF + cx + nxF + nC+nV; 3313385ff7eSMatthew G. Knepley (*cone)[3] = cx*nyF + cy + nC+nV+nXF; 3323385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 3333385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 3343385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF) { 3353385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV) / nxF; 3363385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV) % nxF; 3373385ff7eSMatthew G. Knepley 3383385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3393385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + 1 + nC; 3403385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF) { 3413385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV+nXF) / nyF; 3423385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV+nXF) % nyF; 3433385ff7eSMatthew G. Knepley 3443385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3453385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + nVx + nC; 34698921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 34798921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3483385ff7eSMatthew G. Knepley break; 3493385ff7eSMatthew G. Knepley case 3: 3503385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3513385ff7eSMatthew G. Knepley } 3523385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3533385ff7eSMatthew G. Knepley } 3543385ff7eSMatthew G. Knepley 3553385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) 3563385ff7eSMatthew G. Knepley { 3573385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3583385ff7eSMatthew G. Knepley 3593385ff7eSMatthew G. Knepley PetscFunctionBegin; 36069291d52SBarry Smith ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr); 3613385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3623385ff7eSMatthew G. Knepley } 3633385ff7eSMatthew G. Knepley 3643385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) 3653385ff7eSMatthew G. Knepley { 3663385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 3673385ff7eSMatthew G. Knepley Vec coordinates; 3683385ff7eSMatthew G. Knepley PetscSection section; 3693385ff7eSMatthew G. Knepley PetscScalar *coords; 3703385ff7eSMatthew G. Knepley PetscReal h[3]; 3713385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 3723385ff7eSMatthew G. Knepley PetscErrorCode ierr; 3733385ff7eSMatthew G. Knepley 3743385ff7eSMatthew G. Knepley PetscFunctionBegin; 375a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 376ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 3772c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim > 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"The following code only works for dim <= 3"); 3783385ff7eSMatthew G. Knepley h[0] = (xu - xl)/M; 3793385ff7eSMatthew G. Knepley h[1] = (yu - yl)/N; 3803385ff7eSMatthew G. Knepley h[2] = (zu - zl)/P; 3813385ff7eSMatthew G. Knepley ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 3823385ff7eSMatthew G. Knepley ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr); 3833385ff7eSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 3843385ff7eSMatthew G. Knepley ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr); 3853385ff7eSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr); 3863385ff7eSMatthew G. Knepley ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr); 3873385ff7eSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 3883385ff7eSMatthew G. Knepley ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr); 3893385ff7eSMatthew G. Knepley } 3903385ff7eSMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 3913385ff7eSMatthew G. Knepley ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 3923385ff7eSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr); 39324640c55SToby Isaac ierr = PetscObjectSetName((PetscObject)coordinates,"coordinates");CHKERRQ(ierr); 3943385ff7eSMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 3953385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) { 396e4d69e08SBarry Smith PetscInt ind[3], d, off; 3973385ff7eSMatthew G. Knepley 398e4d69e08SBarry Smith ind[0] = 0; 399e4d69e08SBarry Smith ind[1] = 0; 400e4d69e08SBarry Smith ind[2] = k + da->zs; 4013385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) { 4023385ff7eSMatthew G. Knepley ind[1] = j + da->ys; 4033385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) { 4043385ff7eSMatthew G. Knepley const PetscInt vertex = (k*nVy + j)*nVx + i + vStart; 4053385ff7eSMatthew G. Knepley 4063385ff7eSMatthew G. Knepley ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr); 4073385ff7eSMatthew G. Knepley ind[0] = i + da->xs; 4083385ff7eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 4093385ff7eSMatthew G. Knepley coords[off+d] = h[d]*ind[d]; 4103385ff7eSMatthew G. Knepley } 4113385ff7eSMatthew G. Knepley } 4123385ff7eSMatthew G. Knepley } 4133385ff7eSMatthew G. Knepley } 4143385ff7eSMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 41546e270d4SMatthew G. Knepley ierr = DMSetCoordinateSection(dm, PETSC_DETERMINE, section);CHKERRQ(ierr); 4163385ff7eSMatthew G. Knepley ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr); 417a4b60ecfSMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 4183385ff7eSMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 4193385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 4203385ff7eSMatthew G. Knepley } 4219a800dd8SMatthew G. Knepley 42247c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 42347c6ae99SBarry Smith 42447c6ae99SBarry Smith /*@C 425aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 42647c6ae99SBarry Smith 427d8d19677SJose E. Roman Input Parameters: 42847c6ae99SBarry Smith + da - information about my local patch 42947c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 43047c6ae99SBarry Smith 43147c6ae99SBarry Smith Output Parameters: 43247c6ae99SBarry Smith . vptr - array data structured 43347c6ae99SBarry Smith 43447c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 43547c6ae99SBarry Smith to zero them. 43647c6ae99SBarry Smith 43747c6ae99SBarry Smith Level: advanced 43847c6ae99SBarry Smith 439bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 44047c6ae99SBarry Smith 44147c6ae99SBarry Smith @*/ 4427087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 44347c6ae99SBarry Smith { 44447c6ae99SBarry Smith PetscErrorCode ierr; 44547c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 44647c6ae99SBarry Smith char *iarray_start; 44747c6ae99SBarry Smith void **iptr = (void**)vptr; 44847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 44947c6ae99SBarry Smith 45047c6ae99SBarry Smith PetscFunctionBegin; 451a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 45247c6ae99SBarry Smith if (ghosted) { 453aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 45447c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 45547c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 45647c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 4570298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 4580298fd71SBarry Smith dd->startghostedin[i] = NULL; 45947c6ae99SBarry Smith 46047c6ae99SBarry Smith goto done; 46147c6ae99SBarry Smith } 46247c6ae99SBarry Smith } 46347c6ae99SBarry Smith xs = dd->Xs; 46447c6ae99SBarry Smith ys = dd->Ys; 46547c6ae99SBarry Smith zs = dd->Zs; 46647c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 46747c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 46847c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 46947c6ae99SBarry Smith } else { 470aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 47147c6ae99SBarry Smith if (dd->arrayin[i]) { 47247c6ae99SBarry Smith *iptr = dd->arrayin[i]; 47347c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 4740298fd71SBarry Smith dd->arrayin[i] = NULL; 4750298fd71SBarry Smith dd->startin[i] = NULL; 47647c6ae99SBarry Smith 47747c6ae99SBarry Smith goto done; 47847c6ae99SBarry Smith } 47947c6ae99SBarry Smith } 48047c6ae99SBarry Smith xs = dd->xs; 48147c6ae99SBarry Smith ys = dd->ys; 48247c6ae99SBarry Smith zs = dd->zs; 48347c6ae99SBarry Smith xm = dd->xe-dd->xs; 48447c6ae99SBarry Smith ym = dd->ye-dd->ys; 48547c6ae99SBarry Smith zm = dd->ze-dd->zs; 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith 488c73cfb54SMatthew G. Knepley switch (da->dim) { 48947c6ae99SBarry Smith case 1: { 49047c6ae99SBarry Smith void *ptr; 49147c6ae99SBarry Smith 492901f1932SJed Brown ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 49347c6ae99SBarry Smith 49447c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 49547c6ae99SBarry Smith *iptr = (void*)ptr; 4968865f1eaSKarl Rupp break; 4978865f1eaSKarl Rupp } 49847c6ae99SBarry Smith case 2: { 49947c6ae99SBarry Smith void **ptr; 50047c6ae99SBarry Smith 501901f1932SJed Brown ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 50247c6ae99SBarry Smith 50347c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 5048865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 50547c6ae99SBarry Smith *iptr = (void*)ptr; 5068865f1eaSKarl Rupp break; 5078865f1eaSKarl Rupp } 50847c6ae99SBarry Smith case 3: { 50947c6ae99SBarry Smith void ***ptr,**bptr; 51047c6ae99SBarry Smith 511901f1932SJed Brown ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr); 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 51447c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 5158865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 51647c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 51747c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 51847c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 51947c6ae99SBarry Smith } 52047c6ae99SBarry Smith } 52147c6ae99SBarry Smith *iptr = (void*)ptr; 5228865f1eaSKarl Rupp break; 5238865f1eaSKarl Rupp } 52447c6ae99SBarry Smith default: 52598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",da->dim); 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith 52847c6ae99SBarry Smith done: 52947c6ae99SBarry Smith /* add arrays to the checked out list */ 53047c6ae99SBarry Smith if (ghosted) { 531aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 53247c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 53347c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 53447c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 53547c6ae99SBarry Smith break; 53647c6ae99SBarry Smith } 53747c6ae99SBarry Smith } 53847c6ae99SBarry Smith } else { 539aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 54047c6ae99SBarry Smith if (!dd->arrayout[i]) { 54147c6ae99SBarry Smith dd->arrayout[i] = *iptr; 54247c6ae99SBarry Smith dd->startout[i] = iarray_start; 54347c6ae99SBarry Smith break; 54447c6ae99SBarry Smith } 54547c6ae99SBarry Smith } 54647c6ae99SBarry Smith } 54747c6ae99SBarry Smith PetscFunctionReturn(0); 54847c6ae99SBarry Smith } 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith /*@C 551aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 55247c6ae99SBarry Smith 553d8d19677SJose E. Roman Input Parameters: 55447c6ae99SBarry Smith + da - information about my local patch 55547c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 556cb004a26SBarry Smith - vptr - array data structured 55747c6ae99SBarry Smith 55847c6ae99SBarry Smith Level: advanced 55947c6ae99SBarry Smith 560bcaeba4dSBarry Smith .seealso: DMDAGetArray() 56147c6ae99SBarry Smith 56247c6ae99SBarry Smith @*/ 5637087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 56447c6ae99SBarry Smith { 56547c6ae99SBarry Smith PetscInt i; 566ea78f98cSLisandro Dalcin void **iptr = (void**)vptr,*iarray_start = NULL; 56747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 56847c6ae99SBarry Smith 56947c6ae99SBarry Smith PetscFunctionBegin; 570a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 57147c6ae99SBarry Smith if (ghosted) { 572aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 57347c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 57447c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 5750298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 5760298fd71SBarry Smith dd->startghostedout[i] = NULL; 57747c6ae99SBarry Smith break; 57847c6ae99SBarry Smith } 57947c6ae99SBarry Smith } 580aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 58147c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 58247c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 58347c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 58447c6ae99SBarry Smith break; 58547c6ae99SBarry Smith } 58647c6ae99SBarry Smith } 58747c6ae99SBarry Smith } else { 588aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 58947c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 59047c6ae99SBarry Smith iarray_start = dd->startout[i]; 5910298fd71SBarry Smith dd->arrayout[i] = NULL; 5920298fd71SBarry Smith dd->startout[i] = NULL; 59347c6ae99SBarry Smith break; 59447c6ae99SBarry Smith } 59547c6ae99SBarry Smith } 596aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 59747c6ae99SBarry Smith if (!dd->arrayin[i]) { 59847c6ae99SBarry Smith dd->arrayin[i] = *iptr; 59947c6ae99SBarry Smith dd->startin[i] = iarray_start; 60047c6ae99SBarry Smith break; 60147c6ae99SBarry Smith } 60247c6ae99SBarry Smith } 60347c6ae99SBarry Smith } 60447c6ae99SBarry Smith PetscFunctionReturn(0); 60547c6ae99SBarry Smith } 60647c6ae99SBarry Smith 607