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 PetscInt n,m; 2147c6ae99SBarry Smith Vec vec = (Vec)obj; 2247c6ae99SBarry Smith PetscScalar *array; 2347c6ae99SBarry Smith mxArray *mat; 249a42bb27SBarry Smith DM da; 2547c6ae99SBarry Smith 2647c6ae99SBarry Smith PetscFunctionBegin; 279566063dSJacob Faibussowitsch PetscCall(VecGetDM(vec, &da)); 287a8be351SBarry Smith PetscCheck(da,PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA"); 299566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,0,0,0,&m,&n,0)); 3047c6ae99SBarry Smith 319566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec,&array)); 3247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3347c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxREAL); 3447c6ae99SBarry Smith #else 3547c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX); 3647c6ae99SBarry Smith #endif 379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mxGetPr(mat),array,n*m)); 389566063dSJacob Faibussowitsch PetscCall(PetscObjectName(obj)); 3947c6ae99SBarry Smith engPutVariable((Engine*)mengine,obj->name,mat); 4047c6ae99SBarry Smith 419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec,&array)); 4247c6ae99SBarry Smith PetscFunctionReturn(0); 4347c6ae99SBarry Smith } 4447c6ae99SBarry Smith #endif 4547c6ae99SBarry Smith 467087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_DA(DM da,Vec *g) 4747c6ae99SBarry Smith { 4847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 4947c6ae99SBarry Smith 5047c6ae99SBarry Smith PetscFunctionBegin; 5147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 5247c6ae99SBarry Smith PetscValidPointer(g,2); 539566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,g)); 549566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE)); 559566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*g,dd->w)); 569566063dSJacob Faibussowitsch PetscCall(VecSetType(*g,da->vectype)); 5774427ab1SRichard Tran Mills if (dd->nlocal < da->bind_below) { 589566063dSJacob Faibussowitsch PetscCall(VecSetBindingPropagates(*g,PETSC_TRUE)); 599566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(*g,PETSC_TRUE)); 6074427ab1SRichard Tran Mills } 619566063dSJacob Faibussowitsch PetscCall(VecSetDM(*g, da)); 6247c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 637bc9226cSBarry Smith if (dd->w == 1 && da->dim == 2) { 649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d)); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith #endif 6747c6ae99SBarry Smith PetscFunctionReturn(0); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith 70939f6067SMatthew G. Knepley /*@ 71939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 72939f6067SMatthew G. Knepley 73939f6067SMatthew G. Knepley Input Parameter: 74939f6067SMatthew G. Knepley . dm - The DM object 75939f6067SMatthew G. Knepley 76939f6067SMatthew G. Knepley Output Parameters: 77939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 78939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 79939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 80939f6067SMatthew G. Knepley - numCells - The number of local cells 81939f6067SMatthew G. Knepley 82939f6067SMatthew G. Knepley Level: developer 83939f6067SMatthew G. Knepley 84939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint() 85939f6067SMatthew G. Knepley @*/ 863582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) 8757459e9aSMatthew G Knepley { 8857459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 89c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 9057459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 9157459e9aSMatthew G Knepley const PetscInt nC = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 9257459e9aSMatthew G Knepley 9357459e9aSMatthew G Knepley PetscFunctionBegin; 94a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 953582350cSMatthew G. Knepley if (numCellsX) { 963582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX,2); 973582350cSMatthew G. Knepley *numCellsX = mx; 983582350cSMatthew G. Knepley } 993582350cSMatthew G. Knepley if (numCellsY) { 1008135c375SStefano Zampini PetscValidIntPointer(numCellsY,3); 1013582350cSMatthew G. Knepley *numCellsY = my; 1023582350cSMatthew G. Knepley } 1033582350cSMatthew G. Knepley if (numCellsZ) { 1048135c375SStefano Zampini PetscValidIntPointer(numCellsZ,4); 1053582350cSMatthew G. Knepley *numCellsZ = mz; 1063582350cSMatthew G. Knepley } 10757459e9aSMatthew G Knepley if (numCells) { 1083582350cSMatthew G. Knepley PetscValidIntPointer(numCells,5); 10957459e9aSMatthew G Knepley *numCells = nC; 11057459e9aSMatthew G Knepley } 11157459e9aSMatthew G Knepley PetscFunctionReturn(0); 11257459e9aSMatthew G Knepley } 11357459e9aSMatthew G Knepley 114d6401b93SMatthew G. Knepley /*@ 115d6401b93SMatthew G. Knepley DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 116d6401b93SMatthew G. Knepley 117d6401b93SMatthew G. Knepley Input Parameters: 118d6401b93SMatthew G. Knepley + dm - The DM object 119d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell 120d6401b93SMatthew G. Knepley 121d6401b93SMatthew G. Knepley Output Parameters: 122d6401b93SMatthew G. Knepley . point - The local DM point 123d6401b93SMatthew G. Knepley 124d6401b93SMatthew G. Knepley Level: developer 125d6401b93SMatthew G. Knepley 126d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells() 127d6401b93SMatthew G. Knepley @*/ 128d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) 129d6401b93SMatthew G. Knepley { 130c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 131d6401b93SMatthew G. Knepley DMDALocalInfo info; 132d6401b93SMatthew G. Knepley 133d6401b93SMatthew G. Knepley PetscFunctionBegin; 134a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 135d6401b93SMatthew G. Knepley PetscValidIntPointer(point,5); 1369566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 137*08401ef6SPierre Jolivet if (dim > 0) PetscCheck(!(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); 138*08401ef6SPierre Jolivet if (dim > 1) PetscCheck(!(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); 139*08401ef6SPierre Jolivet if (dim > 2) PetscCheck(!(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); 140d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0); 141d6401b93SMatthew G. Knepley PetscFunctionReturn(0); 142d6401b93SMatthew G. Knepley } 143d6401b93SMatthew G. Knepley 14457459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) 14557459e9aSMatthew G Knepley { 14657459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 147c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 14857459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 14957459e9aSMatthew G Knepley const PetscInt nVx = mx+1; 15057459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my+1) : 1; 15157459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz+1) : 1; 15257459e9aSMatthew G Knepley const PetscInt nV = nVx*nVy*nVz; 15357459e9aSMatthew G Knepley 15457459e9aSMatthew G Knepley PetscFunctionBegin; 15557459e9aSMatthew G Knepley if (numVerticesX) { 15657459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX,2); 15757459e9aSMatthew G Knepley *numVerticesX = nVx; 15857459e9aSMatthew G Knepley } 15957459e9aSMatthew G Knepley if (numVerticesY) { 16057459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY,3); 16157459e9aSMatthew G Knepley *numVerticesY = nVy; 16257459e9aSMatthew G Knepley } 16357459e9aSMatthew G Knepley if (numVerticesZ) { 16457459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ,4); 16557459e9aSMatthew G Knepley *numVerticesZ = nVz; 16657459e9aSMatthew G Knepley } 16757459e9aSMatthew G Knepley if (numVertices) { 16857459e9aSMatthew G Knepley PetscValidIntPointer(numVertices,5); 16957459e9aSMatthew G Knepley *numVertices = nV; 17057459e9aSMatthew G Knepley } 17157459e9aSMatthew G Knepley PetscFunctionReturn(0); 17257459e9aSMatthew G Knepley } 17357459e9aSMatthew G Knepley 17457459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) 17557459e9aSMatthew G Knepley { 17657459e9aSMatthew G Knepley DM_DA *da = (DM_DA*) dm->data; 177c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 17857459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 17957459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1); 18057459e9aSMatthew G Knepley const PetscInt nXF = (mx+1)*nxF; 18157459e9aSMatthew G Knepley const PetscInt nyF = mx*(dim > 2 ? mz : 1); 18257459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0; 18357459e9aSMatthew G Knepley const PetscInt nzF = mx*(dim > 1 ? my : 0); 18457459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0; 18557459e9aSMatthew G Knepley 18657459e9aSMatthew G Knepley PetscFunctionBegin; 18757459e9aSMatthew G Knepley if (numXFacesX) { 18857459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX,2); 18957459e9aSMatthew G Knepley *numXFacesX = nxF; 19057459e9aSMatthew G Knepley } 19157459e9aSMatthew G Knepley if (numXFaces) { 19257459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces,3); 19357459e9aSMatthew G Knepley *numXFaces = nXF; 19457459e9aSMatthew G Knepley } 19557459e9aSMatthew G Knepley if (numYFacesY) { 19657459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY,4); 19757459e9aSMatthew G Knepley *numYFacesY = nyF; 19857459e9aSMatthew G Knepley } 19957459e9aSMatthew G Knepley if (numYFaces) { 20057459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces,5); 20157459e9aSMatthew G Knepley *numYFaces = nYF; 20257459e9aSMatthew G Knepley } 20357459e9aSMatthew G Knepley if (numZFacesZ) { 20457459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ,6); 20557459e9aSMatthew G Knepley *numZFacesZ = nzF; 20657459e9aSMatthew G Knepley } 20757459e9aSMatthew G Knepley if (numZFaces) { 20857459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces,7); 20957459e9aSMatthew G Knepley *numZFaces = nZF; 21057459e9aSMatthew G Knepley } 21157459e9aSMatthew G Knepley PetscFunctionReturn(0); 21257459e9aSMatthew G Knepley } 21357459e9aSMatthew G Knepley 21457459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) 21557459e9aSMatthew G Knepley { 216c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 21757459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 21857459e9aSMatthew G Knepley 21957459e9aSMatthew G Knepley PetscFunctionBegin; 2208865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart,3); 2218865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd,4); 2229566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2239566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2249566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 22557459e9aSMatthew G Knepley if (height == 0) { 22657459e9aSMatthew G Knepley /* Cells */ 2278865f1eaSKarl Rupp if (pStart) *pStart = 0; 2288865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 22957459e9aSMatthew G Knepley } else if (height == 1) { 23057459e9aSMatthew G Knepley /* Faces */ 2318865f1eaSKarl Rupp if (pStart) *pStart = nC+nV; 2328865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 23357459e9aSMatthew G Knepley } else if (height == dim) { 23457459e9aSMatthew G Knepley /* Vertices */ 2358865f1eaSKarl Rupp if (pStart) *pStart = nC; 2368865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV; 23757459e9aSMatthew G Knepley } else if (height < 0) { 23857459e9aSMatthew G Knepley /* All points */ 2398865f1eaSKarl Rupp if (pStart) *pStart = 0; 2408865f1eaSKarl Rupp if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 24198921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height); 24257459e9aSMatthew G Knepley PetscFunctionReturn(0); 24357459e9aSMatthew G Knepley } 24457459e9aSMatthew G Knepley 2453385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) 2463385ff7eSMatthew G. Knepley { 247c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2483385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2493385ff7eSMatthew G. Knepley 2503385ff7eSMatthew G. Knepley PetscFunctionBegin; 2513385ff7eSMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart,3); 2523385ff7eSMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd,4); 2539566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2549566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2559566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 2563385ff7eSMatthew G. Knepley if (depth == dim) { 2573385ff7eSMatthew G. Knepley /* Cells */ 2583385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2593385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC; 2603385ff7eSMatthew G. Knepley } else if (depth == dim-1) { 2613385ff7eSMatthew G. Knepley /* Faces */ 2623385ff7eSMatthew G. Knepley if (pStart) *pStart = nC+nV; 2633385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 2643385ff7eSMatthew G. Knepley } else if (depth == 0) { 2653385ff7eSMatthew G. Knepley /* Vertices */ 2663385ff7eSMatthew G. Knepley if (pStart) *pStart = nC; 2673385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV; 2683385ff7eSMatthew G. Knepley } else if (depth < 0) { 2693385ff7eSMatthew G. Knepley /* All points */ 2703385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2713385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC+nV+nXF+nYF+nZF; 27298921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth); 2733385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2743385ff7eSMatthew G. Knepley } 2753385ff7eSMatthew G. Knepley 2763385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) 2773385ff7eSMatthew G. Knepley { 278c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2793385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2803385ff7eSMatthew G. Knepley 2813385ff7eSMatthew G. Knepley PetscFunctionBegin; 2823385ff7eSMatthew G. Knepley *coneSize = 0; 2839566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2849566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2859566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 2863385ff7eSMatthew G. Knepley switch (dim) { 2873385ff7eSMatthew G. Knepley case 2: 2883385ff7eSMatthew G. Knepley if (p >= 0) { 2893385ff7eSMatthew G. Knepley if (p < nC) { 2903385ff7eSMatthew G. Knepley *coneSize = 4; 2913385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 2923385ff7eSMatthew G. Knepley *coneSize = 0; 2933385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF+nZF) { 2943385ff7eSMatthew G. Knepley *coneSize = 2; 29598921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 29698921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 2973385ff7eSMatthew G. Knepley break; 2983385ff7eSMatthew G. Knepley case 3: 2993385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3003385ff7eSMatthew G. Knepley } 3013385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3023385ff7eSMatthew G. Knepley } 3033385ff7eSMatthew G. Knepley 3043385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) 3053385ff7eSMatthew G. Knepley { 306c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 3073385ff7eSMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 3083385ff7eSMatthew G. Knepley 3093385ff7eSMatthew G. Knepley PetscFunctionBegin; 3109566063dSJacob Faibussowitsch if (!cone) PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 3119566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC)); 3129566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 3139566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF)); 3143385ff7eSMatthew G. Knepley switch (dim) { 3153385ff7eSMatthew G. Knepley case 2: 3163385ff7eSMatthew G. Knepley if (p >= 0) { 3173385ff7eSMatthew G. Knepley if (p < nC) { 3183385ff7eSMatthew G. Knepley const PetscInt cy = p / nCx; 3193385ff7eSMatthew G. Knepley const PetscInt cx = p % nCx; 3203385ff7eSMatthew G. Knepley 3213385ff7eSMatthew G. Knepley (*cone)[0] = cy*nxF + cx + nC+nV; 3223385ff7eSMatthew G. Knepley (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF; 3233385ff7eSMatthew G. Knepley (*cone)[2] = cy*nxF + cx + nxF + nC+nV; 3243385ff7eSMatthew G. Knepley (*cone)[3] = cx*nyF + cy + nC+nV+nXF; 3253385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 3263385ff7eSMatthew G. Knepley } else if (p < nC+nV) { 3273385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF) { 3283385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV) / nxF; 3293385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV) % nxF; 3303385ff7eSMatthew G. Knepley 3313385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3323385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + 1 + nC; 3333385ff7eSMatthew G. Knepley } else if (p < nC+nV+nXF+nYF) { 3343385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC+nV+nXF) / nyF; 3353385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC+nV+nXF) % nyF; 3363385ff7eSMatthew G. Knepley 3373385ff7eSMatthew G. Knepley (*cone)[0] = fy*nVx + fx + nC; 3383385ff7eSMatthew G. Knepley (*cone)[1] = fy*nVx + fx + nVx + nC; 33998921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF); 34098921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p); 3413385ff7eSMatthew G. Knepley break; 3423385ff7eSMatthew G. Knepley case 3: 3433385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3443385ff7eSMatthew G. Knepley } 3453385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3463385ff7eSMatthew G. Knepley } 3473385ff7eSMatthew G. Knepley 3483385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) 3493385ff7eSMatthew G. Knepley { 3503385ff7eSMatthew G. Knepley PetscFunctionBegin; 3519566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 3523385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3533385ff7eSMatthew G. Knepley } 3543385ff7eSMatthew G. Knepley 3553385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) 3563385ff7eSMatthew G. Knepley { 3573385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 3583385ff7eSMatthew G. Knepley Vec coordinates; 3593385ff7eSMatthew G. Knepley PetscSection section; 3603385ff7eSMatthew G. Knepley PetscScalar *coords; 3613385ff7eSMatthew G. Knepley PetscReal h[3]; 3623385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 3633385ff7eSMatthew G. Knepley 3643385ff7eSMatthew G. Knepley PetscFunctionBegin; 365a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA); 3669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 367*08401ef6SPierre Jolivet PetscCheck(dim <= 3,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"The following code only works for dim <= 3"); 3683385ff7eSMatthew G. Knepley h[0] = (xu - xl)/M; 3693385ff7eSMatthew G. Knepley h[1] = (yu - yl)/N; 3703385ff7eSMatthew G. Knepley h[2] = (zu - zl)/P; 3719566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd)); 3729566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 3739566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion)); 3749566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 1)); 3759566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, 0, dim)); 3769566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, vStart, vEnd)); 3773385ff7eSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 3789566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, v, dim)); 3793385ff7eSMatthew G. Knepley } 3809566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(section)); 3819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &size)); 3829566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates)); 3839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates,"coordinates")); 3849566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 3853385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) { 386e4d69e08SBarry Smith PetscInt ind[3], d, off; 3873385ff7eSMatthew G. Knepley 388e4d69e08SBarry Smith ind[0] = 0; 389e4d69e08SBarry Smith ind[1] = 0; 390e4d69e08SBarry Smith ind[2] = k + da->zs; 3913385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) { 3923385ff7eSMatthew G. Knepley ind[1] = j + da->ys; 3933385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) { 3943385ff7eSMatthew G. Knepley const PetscInt vertex = (k*nVy + j)*nVx + i + vStart; 3953385ff7eSMatthew G. Knepley 3969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, vertex, &off)); 3973385ff7eSMatthew G. Knepley ind[0] = i + da->xs; 3983385ff7eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 3993385ff7eSMatthew G. Knepley coords[off+d] = h[d]*ind[d]; 4003385ff7eSMatthew G. Knepley } 4013385ff7eSMatthew G. Knepley } 4023385ff7eSMatthew G. Knepley } 4033385ff7eSMatthew G. Knepley } 4049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 4059566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section)); 4069566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 4079566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 4089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 4093385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 4103385ff7eSMatthew G. Knepley } 4119a800dd8SMatthew G. Knepley 41247c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 41347c6ae99SBarry Smith 41447c6ae99SBarry Smith /*@C 415aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 41647c6ae99SBarry Smith 417d8d19677SJose E. Roman Input Parameters: 41847c6ae99SBarry Smith + da - information about my local patch 41947c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith Output Parameters: 42247c6ae99SBarry Smith . vptr - array data structured 42347c6ae99SBarry Smith 42447c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 42547c6ae99SBarry Smith to zero them. 42647c6ae99SBarry Smith 42747c6ae99SBarry Smith Level: advanced 42847c6ae99SBarry Smith 429bcaeba4dSBarry Smith .seealso: DMDARestoreArray() 43047c6ae99SBarry Smith 43147c6ae99SBarry Smith @*/ 4327087cfbeSBarry Smith PetscErrorCode DMDAGetArray(DM da,PetscBool ghosted,void *vptr) 43347c6ae99SBarry Smith { 43447c6ae99SBarry Smith PetscInt j,i,xs,ys,xm,ym,zs,zm; 43547c6ae99SBarry Smith char *iarray_start; 43647c6ae99SBarry Smith void **iptr = (void**)vptr; 43747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 43847c6ae99SBarry Smith 43947c6ae99SBarry Smith PetscFunctionBegin; 440a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 44147c6ae99SBarry Smith if (ghosted) { 442aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 44347c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 44447c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 44547c6ae99SBarry Smith iarray_start = (char*)dd->startghostedin[i]; 4460298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 4470298fd71SBarry Smith dd->startghostedin[i] = NULL; 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith goto done; 45047c6ae99SBarry Smith } 45147c6ae99SBarry Smith } 45247c6ae99SBarry Smith xs = dd->Xs; 45347c6ae99SBarry Smith ys = dd->Ys; 45447c6ae99SBarry Smith zs = dd->Zs; 45547c6ae99SBarry Smith xm = dd->Xe-dd->Xs; 45647c6ae99SBarry Smith ym = dd->Ye-dd->Ys; 45747c6ae99SBarry Smith zm = dd->Ze-dd->Zs; 45847c6ae99SBarry Smith } else { 459aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 46047c6ae99SBarry Smith if (dd->arrayin[i]) { 46147c6ae99SBarry Smith *iptr = dd->arrayin[i]; 46247c6ae99SBarry Smith iarray_start = (char*)dd->startin[i]; 4630298fd71SBarry Smith dd->arrayin[i] = NULL; 4640298fd71SBarry Smith dd->startin[i] = NULL; 46547c6ae99SBarry Smith 46647c6ae99SBarry Smith goto done; 46747c6ae99SBarry Smith } 46847c6ae99SBarry Smith } 46947c6ae99SBarry Smith xs = dd->xs; 47047c6ae99SBarry Smith ys = dd->ys; 47147c6ae99SBarry Smith zs = dd->zs; 47247c6ae99SBarry Smith xm = dd->xe-dd->xs; 47347c6ae99SBarry Smith ym = dd->ye-dd->ys; 47447c6ae99SBarry Smith zm = dd->ze-dd->zs; 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith 477c73cfb54SMatthew G. Knepley switch (da->dim) { 47847c6ae99SBarry Smith case 1: { 47947c6ae99SBarry Smith void *ptr; 48047c6ae99SBarry Smith 4819566063dSJacob Faibussowitsch PetscCall(PetscMalloc(xm*sizeof(PetscScalar),&iarray_start)); 48247c6ae99SBarry Smith 48347c6ae99SBarry Smith ptr = (void*)(iarray_start - xs*sizeof(PetscScalar)); 48447c6ae99SBarry Smith *iptr = (void*)ptr; 4858865f1eaSKarl Rupp break; 4868865f1eaSKarl Rupp } 48747c6ae99SBarry Smith case 2: { 48847c6ae99SBarry Smith void **ptr; 48947c6ae99SBarry Smith 4909566063dSJacob Faibussowitsch PetscCall(PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start)); 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*)); 4938865f1eaSKarl Rupp for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs); 49447c6ae99SBarry Smith *iptr = (void*)ptr; 4958865f1eaSKarl Rupp break; 4968865f1eaSKarl Rupp } 49747c6ae99SBarry Smith case 3: { 49847c6ae99SBarry Smith void ***ptr,**bptr; 49947c6ae99SBarry Smith 5009566063dSJacob Faibussowitsch PetscCall(PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start)); 50147c6ae99SBarry Smith 50247c6ae99SBarry Smith ptr = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*)); 50347c6ae99SBarry Smith bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**)); 5048865f1eaSKarl Rupp for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys); 50547c6ae99SBarry Smith for (i=zs; i<zs+zm; i++) { 50647c6ae99SBarry Smith for (j=ys; j<ys+ym; j++) { 50747c6ae99SBarry Smith ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs); 50847c6ae99SBarry Smith } 50947c6ae99SBarry Smith } 51047c6ae99SBarry Smith *iptr = (void*)ptr; 5118865f1eaSKarl Rupp break; 5128865f1eaSKarl Rupp } 51347c6ae99SBarry Smith default: 51498921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",da->dim); 51547c6ae99SBarry Smith } 51647c6ae99SBarry Smith 51747c6ae99SBarry Smith done: 51847c6ae99SBarry Smith /* add arrays to the checked out list */ 51947c6ae99SBarry Smith if (ghosted) { 520aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 52147c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 52247c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 52347c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 52447c6ae99SBarry Smith break; 52547c6ae99SBarry Smith } 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith } else { 528aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 52947c6ae99SBarry Smith if (!dd->arrayout[i]) { 53047c6ae99SBarry Smith dd->arrayout[i] = *iptr; 53147c6ae99SBarry Smith dd->startout[i] = iarray_start; 53247c6ae99SBarry Smith break; 53347c6ae99SBarry Smith } 53447c6ae99SBarry Smith } 53547c6ae99SBarry Smith } 53647c6ae99SBarry Smith PetscFunctionReturn(0); 53747c6ae99SBarry Smith } 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith /*@C 540aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 54147c6ae99SBarry Smith 542d8d19677SJose E. Roman Input Parameters: 54347c6ae99SBarry Smith + da - information about my local patch 54447c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 545cb004a26SBarry Smith - vptr - array data structured 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith Level: advanced 54847c6ae99SBarry Smith 549bcaeba4dSBarry Smith .seealso: DMDAGetArray() 55047c6ae99SBarry Smith 55147c6ae99SBarry Smith @*/ 5527087cfbeSBarry Smith PetscErrorCode DMDARestoreArray(DM da,PetscBool ghosted,void *vptr) 55347c6ae99SBarry Smith { 55447c6ae99SBarry Smith PetscInt i; 555ea78f98cSLisandro Dalcin void **iptr = (void**)vptr,*iarray_start = NULL; 55647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 55747c6ae99SBarry Smith 55847c6ae99SBarry Smith PetscFunctionBegin; 559a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 56047c6ae99SBarry Smith if (ghosted) { 561aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 56247c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 56347c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 5640298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 5650298fd71SBarry Smith dd->startghostedout[i] = NULL; 56647c6ae99SBarry Smith break; 56747c6ae99SBarry Smith } 56847c6ae99SBarry Smith } 569aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 57047c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 57147c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 57247c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 57347c6ae99SBarry Smith break; 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith } 57647c6ae99SBarry Smith } else { 577aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 57847c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 57947c6ae99SBarry Smith iarray_start = dd->startout[i]; 5800298fd71SBarry Smith dd->arrayout[i] = NULL; 5810298fd71SBarry Smith dd->startout[i] = NULL; 58247c6ae99SBarry Smith break; 58347c6ae99SBarry Smith } 58447c6ae99SBarry Smith } 585aa219208SBarry Smith for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) { 58647c6ae99SBarry Smith if (!dd->arrayin[i]) { 58747c6ae99SBarry Smith dd->arrayin[i] = *iptr; 58847c6ae99SBarry Smith dd->startin[i] = iarray_start; 58947c6ae99SBarry Smith break; 59047c6ae99SBarry Smith } 59147c6ae99SBarry Smith } 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith PetscFunctionReturn(0); 59447c6ae99SBarry Smith } 595