PetscErrorCode DMDASubDomainIS_Private(DM dm,PetscInt n,DM *subdm,IS **iis,IS **ois) { PetscErrorCode ierr; PetscInt i; DMDALocalInfo info,subinfo; MatStencil lower,upper; PetscFunctionBegin; ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); if (iis) {ierr = PetscMalloc1(n,iis);CHKERRQ(ierr);} if (ois) {ierr = PetscMalloc1(n,ois);CHKERRQ(ierr);} for (i = 0;i < n; i++) { ierr = DMDAGetLocalInfo(subdm[i],&subinfo);CHKERRQ(ierr); if (iis) { /* create the inner IS */ lower.i = info.xs; lower.j = info.ys; lower.k = info.zs; upper.i = info.xs+info.xm; upper.j = info.ys+info.ym; upper.k = info.zs+info.zm; ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*iis)[i]);CHKERRQ(ierr); } if (ois) { /* create the outer IS */ lower.i = subinfo.xs; lower.j = subinfo.ys; lower.k = subinfo.zs; upper.i = subinfo.xs+subinfo.xm; upper.j = subinfo.ys+subinfo.ym; upper.k = subinfo.zs+subinfo.zm; ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*ois)[i]);CHKERRQ(ierr); } } PetscFunctionReturn(0); }
int main(int argc,char **argv) { PetscErrorCode ierr; DM da,*subda; PetscInt i,dim=3; PetscMPIInt size,rank; Vec v; Vec slvec,sgvec; IS *ois,*iis; VecScatter oscata; VecScatter *iscat,*oscat,*gscat; DMDALocalInfo info; ierr = PetscInitialize(&argc,&argv,(char*)0,help);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); /* Create distributed array and get vectors */ ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); if (dim == 2) { ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,&da);CHKERRQ(ierr); } else if (dim == 3) { ierr = DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,-4,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,NULL,&da);CHKERRQ(ierr); } ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); ierr = DMCreateDomainDecomposition(da,NULL,NULL,&iis,&ois,&subda);CHKERRQ(ierr); ierr = DMCreateDomainDecompositionScatters(da,1,subda,&iscat,&oscat,&gscat);CHKERRQ(ierr); { DMDALocalInfo subinfo; MatStencil lower,upper; IS patchis,subpatchis; Vec smallvec; Vec largevec; VecScatter patchscat; ierr = DMDAGetLocalInfo(subda[0],&subinfo);CHKERRQ(ierr); lower.i = info.xs; lower.j = info.ys; lower.k = info.zs; upper.i = info.xs+info.xm; upper.j = info.ys+info.ym; upper.k = info.zs+info.zm; /* test the patch IS as a thing to scatter to/from */ ierr = DMDACreatePatchIS(da,&lower,&upper,&patchis);CHKERRQ(ierr); ierr = DMGetGlobalVector(da,&largevec);CHKERRQ(ierr); ierr = VecCreate(PETSC_COMM_SELF,&smallvec);CHKERRQ(ierr); ierr = VecSetSizes(smallvec,info.dof*(upper.i - lower.i)*(upper.j - lower.j)*(upper.k - lower.k),PETSC_DECIDE);CHKERRQ(ierr); ierr = VecSetFromOptions(smallvec);CHKERRQ(ierr); ierr = VecScatterCreate(smallvec,NULL,largevec,patchis,&patchscat);CHKERRQ(ierr); ierr = FillLocalSubdomain(subda[0],smallvec);CHKERRQ(ierr); ierr = VecSet(largevec,0);CHKERRQ(ierr); ierr = VecScatterBegin(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(patchscat,smallvec,largevec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); for (i = 0; i < size; i++) { if (i == rank) { ierr = ISView(patchis,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); ierr = VecScatterView(patchscat,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); ierr = VecView(smallvec,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); } ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); } ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); ierr = VecView(largevec,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = VecDestroy(&smallvec);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(da,&largevec);CHKERRQ(ierr); ierr = ISDestroy(&patchis);CHKERRQ(ierr); ierr = VecScatterDestroy(&patchscat);CHKERRQ(ierr); } /* view the various parts */ { for (i = 0; i < size; i++) { if (i == rank) { ierr = PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i);CHKERRQ(ierr); ierr = DMView(subda[0],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); } ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); } ierr = DMGetLocalVector(subda[0],&slvec);CHKERRQ(ierr); ierr = DMGetGlobalVector(subda[0],&sgvec);CHKERRQ(ierr); ierr = DMGetGlobalVector(da,&v);CHKERRQ(ierr); /* test filling outer between the big DM and the small ones with the IS scatter*/ ierr = VecScatterCreate(v,ois[0],sgvec,NULL,&oscata);CHKERRQ(ierr); ierr = FillLocalSubdomain(subda[0],sgvec);CHKERRQ(ierr); ierr = VecScatterBegin(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterEnd(oscata,sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); /* test the local-to-local scatter */ /* fill up the local subdomain and then add them together */ ierr = FillLocalSubdomain(da,v);CHKERRQ(ierr); ierr = VecScatterBegin(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(gscat[0],v,slvec,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecView(v,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* test ghost scattering backwards */ ierr = VecSet(v,0);CHKERRQ(ierr); ierr = VecScatterBegin(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterEnd(gscat[0],slvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecView(v,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* test overlap scattering backwards */ ierr = DMLocalToGlobalBegin(subda[0],slvec,ADD_VALUES,sgvec);CHKERRQ(ierr); ierr = DMLocalToGlobalEnd(subda[0],slvec,ADD_VALUES,sgvec);CHKERRQ(ierr); ierr = VecSet(v,0);CHKERRQ(ierr); ierr = VecScatterBegin(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterEnd(oscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecView(v,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* test interior scattering backwards */ ierr = VecSet(v,0);CHKERRQ(ierr); ierr = VecScatterBegin(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterEnd(iscat[0],sgvec,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecView(v,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* test matrix allocation */ for (i = 0; i < size; i++) { if (i == rank) { Mat m; ierr = PetscPrintf(PETSC_COMM_SELF,"Processor %d: \n",i);CHKERRQ(ierr); ierr = DMSetMatType(subda[0],MATAIJ);CHKERRQ(ierr); ierr = DMCreateMatrix(subda[0],&m);CHKERRQ(ierr); ierr = MatView(m,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); ierr = MatDestroy(&m);CHKERRQ(ierr); } ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); } ierr = DMRestoreLocalVector(subda[0],&slvec);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(subda[0],&sgvec);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(da,&v);CHKERRQ(ierr); } ierr = DMDestroy(&subda[0]);CHKERRQ(ierr); ierr = ISDestroy(&ois[0]);CHKERRQ(ierr); ierr = ISDestroy(&iis[0]);CHKERRQ(ierr); ierr = VecScatterDestroy(&iscat[0]);CHKERRQ(ierr); ierr = VecScatterDestroy(&oscat[0]);CHKERRQ(ierr); ierr = VecScatterDestroy(&gscat[0]);CHKERRQ(ierr); ierr = VecScatterDestroy(&oscata);CHKERRQ(ierr); ierr = PetscFree(iscat);CHKERRQ(ierr); ierr = PetscFree(oscat);CHKERRQ(ierr); ierr = PetscFree(gscat);CHKERRQ(ierr); ierr = PetscFree(oscata);CHKERRQ(ierr); ierr = PetscFree(subda);CHKERRQ(ierr); ierr = PetscFree(ois);CHKERRQ(ierr); ierr = PetscFree(iis);CHKERRQ(ierr); ierr = DMDestroy(&da);CHKERRQ(ierr); ierr = PetscFinalize(); return 0; }
/* DMPatchZoom - Create a version of the coarse patch (identified by rank) with halo on communicator commz Collective on DM Input Parameters: + dm - the DM . rank - the rank which holds the given patch - commz - the new communicator for the patch Output Parameters: + dmz - the patch DM . sfz - the PetscSF mapping the patch+halo to the zoomed version . sfzr - the PetscSF mapping the patch to the restricted zoomed version Level: intermediate Note: All processes in commz should have the same rank (could autosplit comm) .seealso: DMPatchSolve() */ PetscErrorCode DMPatchZoom(DM dm, Vec X, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr) { DMDAStencilType st; MatStencil blower, bupper, loclower, locupper; IS is; const PetscInt *ranges, *indices; PetscInt *localPoints = NULL; PetscSFNode *remotePoints = NULL; PetscInt dim, dof; PetscInt M, N, P, rM, rN, rP, halo = 1, sxb, syb, szb, sxr, syr, szr, exr, eyr, ezr, mxb, myb, mzb, i, j, k, q; PetscMPIInt size; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRQ(ierr); /* Create patch DM */ ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0, &dof, 0,0,0,0, &st);CHKERRQ(ierr); /* Get piece for rank r, expanded by halo */ bupper.i = PetscMin(M, upper.i + halo); blower.i = PetscMax(lower.i - halo, 0); bupper.j = PetscMin(N, upper.j + halo); blower.j = PetscMax(lower.j - halo, 0); bupper.k = PetscMin(P, upper.k + halo); blower.k = PetscMax(lower.k - halo, 0); rM = bupper.i - blower.i; rN = bupper.j - blower.j; rP = bupper.k - blower.k; if (commz != MPI_COMM_NULL) { ierr = DMDACreate(commz, dmz);CHKERRQ(ierr); ierr = DMSetDimension(*dmz, dim);CHKERRQ(ierr); ierr = DMDASetSizes(*dmz, rM, rN, rP);CHKERRQ(ierr); ierr = DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr); ierr = DMDASetBoundaryType(*dmz, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);CHKERRQ(ierr); ierr = DMDASetDof(*dmz, dof);CHKERRQ(ierr); ierr = DMDASetStencilType(*dmz, st);CHKERRQ(ierr); ierr = DMDASetStencilWidth(*dmz, 0);CHKERRQ(ierr); ierr = DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL);CHKERRQ(ierr); ierr = DMSetFromOptions(*dmz);CHKERRQ(ierr); ierr = DMSetUp(*dmz);CHKERRQ(ierr); ierr = DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb);CHKERRQ(ierr); sxr = PetscMax(sxb, lower.i - blower.i); syr = PetscMax(syb, lower.j - blower.j); szr = PetscMax(szb, lower.k - blower.k); exr = PetscMin(sxb+mxb, upper.i - blower.i); eyr = PetscMin(syb+myb, upper.j - blower.j); ezr = PetscMin(szb+mzb, upper.k - blower.k); ierr = PetscMalloc2(rM*rN*rP,&localPoints,rM*rN*rP,&remotePoints);CHKERRQ(ierr); } else { sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0; } /* Create SF for restricted map */ ierr = VecGetOwnershipRanges(X,&ranges);CHKERRQ(ierr); loclower.i = blower.i + sxr; locupper.i = blower.i + exr; loclower.j = blower.j + syr; locupper.j = blower.j + eyr; loclower.k = blower.k + szr; locupper.k = blower.k + ezr; ierr = DMDACreatePatchIS(dm, &loclower, &locupper, &is);CHKERRQ(ierr); ierr = ISGetIndices(is, &indices);CHKERRQ(ierr); q = 0; for (k = szb; k < szb+mzb; ++k) { if ((k < szr) || (k >= ezr)) continue; for (j = syb; j < syb+myb; ++j) { if ((j < syr) || (j >= eyr)) continue; for (i = sxb; i < sxb+mxb; ++i) { const PetscInt lp = ((k-szb)*rN + (j-syb))*rM + i-sxb; PetscInt r; if ((i < sxr) || (i >= exr)) continue; localPoints[q] = lp; ierr = PetscFindInt(indices[q], size+1, ranges, &r);CHKERRQ(ierr); remotePoints[q].rank = r < 0 ? -(r+1) - 1 : r; remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank]; ++q; } } } ierr = ISRestoreIndices(is, &indices);CHKERRQ(ierr); ierr = ISDestroy(&is);CHKERRQ(ierr); ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) *sfzr, "Restricted Map");CHKERRQ(ierr); ierr = PetscSFSetGraph(*sfzr, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);CHKERRQ(ierr); /* Create SF for buffered map */ loclower.i = blower.i + sxb; locupper.i = blower.i + sxb+mxb; loclower.j = blower.j + syb; locupper.j = blower.j + syb+myb; loclower.k = blower.k + szb; locupper.k = blower.k + szb+mzb; ierr = DMDACreatePatchIS(dm, &loclower, &locupper, &is);CHKERRQ(ierr); ierr = ISGetIndices(is, &indices);CHKERRQ(ierr); q = 0; for (k = szb; k < szb+mzb; ++k) { for (j = syb; j < syb+myb; ++j) { for (i = sxb; i < sxb+mxb; ++i, ++q) { PetscInt r; localPoints[q] = q; ierr = PetscFindInt(indices[q], size+1, ranges, &r);CHKERRQ(ierr); remotePoints[q].rank = r < 0 ? -(r+1) - 1 : r; remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank]; } } } ierr = ISRestoreIndices(is, &indices);CHKERRQ(ierr); ierr = ISDestroy(&is);CHKERRQ(ierr); ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) *sfz, "Buffered Map");CHKERRQ(ierr); ierr = PetscSFSetGraph(*sfz, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscFree2(localPoints, remotePoints);CHKERRQ(ierr); PetscFunctionReturn(0); }
PETSC_EXTERN void PETSC_STDCALL dmdacreatepatchis_(DM da,MatStencil *lower,MatStencil *upper,IS *is, int *__ierr ){ *__ierr = DMDACreatePatchIS( (DM)PetscToPointer((da) ),lower,upper,is); }
/* Fills the local vector problem on the subdomain from the global problem. Right now this assumes one subdomain per processor. */ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) { PetscErrorCode ierr; DMDALocalInfo info,subinfo; DM subdm; MatStencil upper,lower; IS idis,isis,odis,osis,gdis; Vec svec,dvec,slvec; PetscInt xm,ym,zm,xs,ys,zs; PetscInt i; PetscFunctionBegin; /* allocate the arrays of scatters */ if (iscat) {ierr = PetscMalloc1(nsubdms,iscat);CHKERRQ(ierr);} if (oscat) {ierr = PetscMalloc1(nsubdms,oscat);CHKERRQ(ierr);} if (lscat) {ierr = PetscMalloc1(nsubdms,lscat);CHKERRQ(ierr);} ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); for (i = 0; i < nsubdms; i++) { subdm = subdms[i]; ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); ierr = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); /* create the global and subdomain index sets for the inner domain */ lower.i = xs; lower.j = ys; lower.k = zs; upper.i = xs+xm; upper.j = ys+ym; upper.k = zs+zm; ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); /* create the global and subdomain index sets for the outer subdomain */ lower.i = subinfo.xs; lower.j = subinfo.ys; lower.k = subinfo.zs; upper.i = subinfo.xs+subinfo.xm; upper.j = subinfo.ys+subinfo.ym; upper.k = subinfo.zs+subinfo.zm; ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); /* global and subdomain ISes for the local indices of the subdomain */ /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ lower.i = subinfo.gxs; lower.j = subinfo.gys; lower.k = subinfo.gzs; upper.i = subinfo.gxs+subinfo.gxm; upper.j = subinfo.gys+subinfo.gym; upper.k = subinfo.gzs+subinfo.gzm; ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); /* form the scatter */ ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[i]);CHKERRQ(ierr);} if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[i]);CHKERRQ(ierr);} if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[i]);CHKERRQ(ierr);} ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); ierr = ISDestroy(&idis);CHKERRQ(ierr); ierr = ISDestroy(&isis);CHKERRQ(ierr); ierr = ISDestroy(&odis);CHKERRQ(ierr); ierr = ISDestroy(&osis);CHKERRQ(ierr); ierr = ISDestroy(&gdis);CHKERRQ(ierr); } PetscFunctionReturn(0); }