190c77898SPeter Brune #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3e30e807fSPeter Brune #undef __FUNCT__ 4ff9846d9SPeter Brune #define __FUNCT__ "DMDACreatePatchIS" 5063654ddSPeter Brune /*@ 6063654ddSPeter Brune DMDACreatePatchIS - Creates an index set corresponding to a patch of the DA. 7ff9846d9SPeter Brune 8063654ddSPeter Brune Not Collective 9063654ddSPeter Brune 10063654ddSPeter Brune Input Parameters: 11063654ddSPeter Brune + da - the DMDA 12063654ddSPeter Brune . lower - a matstencil with i, j and k corresponding to the lower corner of the patch 13063654ddSPeter Brune - upper - a matstencil with i, j and k corresponding to the upper corner of the patch 14063654ddSPeter Brune 15063654ddSPeter Brune Output Parameters: 16063654ddSPeter Brune . is - the IS corresponding to the patch 17063654ddSPeter Brune 18063654ddSPeter Brune Level: developer 19063654ddSPeter Brune 20063654ddSPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDACreateDomainDecompositionScatters() 21063654ddSPeter Brune @*/ 22ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) 23ff9846d9SPeter Brune { 24063654ddSPeter Brune PetscInt ms=0,ns=0,ps=0; 25063654ddSPeter Brune PetscInt me=1,ne=1,pe=1; 26063654ddSPeter Brune PetscInt mr=0,nr=0,pr=0; 27ff9846d9SPeter Brune PetscInt ii,jj,kk; 28063654ddSPeter Brune PetscInt si,sj,sk; 29063654ddSPeter Brune PetscInt i,j,k,l,idx; 30ff9846d9SPeter Brune PetscInt base; 31063654ddSPeter Brune PetscInt xm=1,ym=1,zm=1; 32063654ddSPeter Brune const PetscInt *lx,*ly,*lz; 33ff9846d9SPeter Brune PetscInt ox,oy,oz; 34063654ddSPeter Brune PetscInt m,n,p,M,N,P,dof; 35063654ddSPeter Brune PetscInt nindices; 36063654ddSPeter Brune PetscInt *indices; 37063654ddSPeter Brune DM_DA *dd = (DM_DA*)da->data; 38063654ddSPeter Brune PetscErrorCode ierr; 39ff9846d9SPeter Brune 400adebc6cSBarry Smith PetscFunctionBegin; 41063654ddSPeter Brune /* need to get the sizes of the actual DM rather than the "global" space of a subdomain DM */ 42063654ddSPeter Brune M = dd->M;N = dd->N;P=dd->P; 43063654ddSPeter Brune m = dd->m;n = dd->n;p=dd->p; 44063654ddSPeter Brune dof = dd->w; 450298fd71SBarry Smith ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr); 46063654ddSPeter Brune ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 47063654ddSPeter Brune nindices = (upper->i - lower->i)*(upper->j - lower->j)*(upper->k - lower->k)*dof; 48063654ddSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*nindices,&indices);CHKERRQ(ierr); 49063654ddSPeter Brune /* start at index 0 on processor 0 */ 50063654ddSPeter Brune mr = 0; 51063654ddSPeter Brune nr = 0; 52063654ddSPeter Brune pr = 0; 53063654ddSPeter Brune ms = 0; 54063654ddSPeter Brune ns = 0; 55063654ddSPeter Brune ps = 0; 56063654ddSPeter Brune if (lx) me = lx[0]; 57063654ddSPeter Brune if (ly) ne = ly[0]; 58063654ddSPeter Brune if (lz) pe = lz[0]; 59ff9846d9SPeter Brune idx = 0; 60ff9846d9SPeter Brune for (k=lower->k-oz;k<upper->k-oz;k++) { 61ff9846d9SPeter Brune for (j=lower->j-oy;j < upper->j-oy;j++) { 62063654ddSPeter Brune for (i=lower->i-ox;i < upper->i-ox;i++) { 63063654ddSPeter Brune /* "actual" indices rather than ones outside of the domain */ 64063654ddSPeter Brune ii = i; 65063654ddSPeter Brune jj = j; 66063654ddSPeter Brune kk = k; 67063654ddSPeter Brune if (ii < 0) ii = M + ii; 68063654ddSPeter Brune if (jj < 0) jj = N + jj; 69063654ddSPeter Brune if (kk < 0) kk = P + kk; 70063654ddSPeter Brune if (ii > M-1) ii = ii - M; 71063654ddSPeter Brune if (jj > N-1) jj = jj - N; 72063654ddSPeter Brune if (kk > P-1) kk = kk - P; 73063654ddSPeter Brune /* gone out of processor range on x axis */ 74063654ddSPeter Brune while(ii > me-1 || ii < ms) { 75063654ddSPeter Brune if (mr == m-1) { 76063654ddSPeter Brune ms = 0; 77063654ddSPeter Brune me = lx[0]; 78063654ddSPeter Brune mr = 0; 79063654ddSPeter Brune } else { 80063654ddSPeter Brune mr++; 81063654ddSPeter Brune ms = me; 82063654ddSPeter Brune me += lx[mr]; 83063654ddSPeter Brune } 84063654ddSPeter Brune } 85063654ddSPeter Brune /* gone out of processor range on y axis */ 86063654ddSPeter Brune while(jj > ne-1 || jj < ns) { 87063654ddSPeter Brune if (nr == n-1) { 88063654ddSPeter Brune ns = 0; 89063654ddSPeter Brune ne = ly[0]; 90063654ddSPeter Brune nr = 0; 91063654ddSPeter Brune } else { 92063654ddSPeter Brune nr++; 93063654ddSPeter Brune ns = ne; 94063654ddSPeter Brune ne += ly[nr]; 95063654ddSPeter Brune } 96063654ddSPeter Brune } 97063654ddSPeter Brune /* gone out of processor range on z axis */ 98063654ddSPeter Brune while(kk > pe-1 || kk < ps) { 99063654ddSPeter Brune if (pr == p-1) { 100063654ddSPeter Brune ps = 0; 101063654ddSPeter Brune pe = lz[0]; 102063654ddSPeter Brune pr = 0; 103063654ddSPeter Brune } else { 104063654ddSPeter Brune pr++; 105063654ddSPeter Brune ps = pe; 106063654ddSPeter Brune pe += lz[pr]; 107063654ddSPeter Brune } 108063654ddSPeter Brune } 109063654ddSPeter Brune /* compute the vector base on owning processor */ 110063654ddSPeter Brune xm = me - ms; 111063654ddSPeter Brune ym = ne - ns; 112063654ddSPeter Brune zm = pe - ps; 113063654ddSPeter Brune base = ms*ym*zm + ns*M + ps*M*N; 114063654ddSPeter Brune /* compute the local coordinates on owning processor */ 115063654ddSPeter Brune si = ii - ms; 116063654ddSPeter Brune sj = jj - ns; 117063654ddSPeter Brune sk = kk - ps; 118063654ddSPeter Brune for (l=0;l<dof;l++) { 119063654ddSPeter Brune indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk); 120ff9846d9SPeter Brune idx++; 121ff9846d9SPeter Brune } 122ff9846d9SPeter Brune } 123ff9846d9SPeter Brune } 124063654ddSPeter Brune } 125063654ddSPeter Brune ISCreateGeneral(PETSC_COMM_SELF,idx,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr); 126ff9846d9SPeter Brune PetscFunctionReturn(0); 127ff9846d9SPeter Brune } 128ff9846d9SPeter Brune 129ff9846d9SPeter Brune #undef __FUNCT__ 130e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private" 13190c77898SPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) 1320adebc6cSBarry Smith { 13390c77898SPeter Brune DM *da; 13490c77898SPeter Brune PetscInt dim,size,i,j,k,idx; 135e30e807fSPeter Brune PetscErrorCode ierr; 136e30e807fSPeter Brune DMDALocalInfo info; 137ff9846d9SPeter Brune PetscInt xsize,ysize,zsize; 138e30e807fSPeter Brune PetscInt xo,yo,zo; 13990c77898SPeter Brune PetscInt xs,ys,zs; 14090c77898SPeter Brune PetscInt xm=1,ym=1,zm=1; 1417ddda789SPeter Brune PetscInt xol,yol,zol; 14290c77898SPeter Brune PetscInt m=1,n=1,p=1; 14390c77898SPeter Brune PetscInt M,N,P; 14490c77898SPeter Brune PetscInt pm,mtmp; 145e30e807fSPeter Brune 146e30e807fSPeter Brune PetscFunctionBegin; 147e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1487ddda789SPeter Brune ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr); 14990c77898SPeter Brune ierr = DMDAGetNumLocalSubDomains(dm,&size);CHKERRQ(ierr); 15090c77898SPeter Brune ierr = PetscMalloc(size*sizeof(DM),&da);CHKERRQ(ierr); 1517ddda789SPeter Brune 15290c77898SPeter Brune if (nlocal) *nlocal = size; 153e30e807fSPeter Brune 15490c77898SPeter Brune dim = info.dim; 155e30e807fSPeter Brune 15690c77898SPeter Brune M = info.xm; 15790c77898SPeter Brune N = info.ym; 15890c77898SPeter Brune P = info.zm; 15990c77898SPeter Brune 16090c77898SPeter Brune if (dim == 1) { 16190c77898SPeter Brune m = size; 16290c77898SPeter Brune } else if (dim == 2) { 16390c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 16490c77898SPeter Brune while (m > 0) { 16590c77898SPeter Brune n = size/m; 16690c77898SPeter Brune if (m*n*p == size) break; 16790c77898SPeter Brune m--; 16890c77898SPeter Brune } 16990c77898SPeter Brune } else if (dim == 3) { 17090c77898SPeter Brune n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); if (!n) n = 1; 17190c77898SPeter Brune while (n > 0) { 17290c77898SPeter Brune pm = size/n; 17390c77898SPeter Brune if (n*pm == size) break; 17490c77898SPeter Brune n--; 17590c77898SPeter Brune } 17690c77898SPeter Brune if (!n) n = 1; 17790c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 17890c77898SPeter Brune if (!m) m = 1; 17990c77898SPeter Brune while (m > 0) { 18090c77898SPeter Brune p = size/(m*n); 18190c77898SPeter Brune if (m*n*p == size) break; 18290c77898SPeter Brune m--; 18390c77898SPeter Brune } 18490c77898SPeter Brune if (M > P && m < p) {mtmp = m; m = p; p = mtmp;} 18590c77898SPeter Brune } 18690c77898SPeter Brune 18790c77898SPeter Brune zs = info.zs; 18890c77898SPeter Brune idx = 0; 18990c77898SPeter Brune for (k = 0; k < p; k++) { 19090c77898SPeter Brune ys = info.ys; 19190c77898SPeter Brune for (j = 0; j < n; j++) { 19290c77898SPeter Brune xs = info.xs; 19390c77898SPeter Brune for (i = 0; i < m; i++) { 19490c77898SPeter Brune if (dim == 1) { 19590c77898SPeter Brune xm = M/m + ((M % m) > i); 19690c77898SPeter Brune } else if (dim == 2) { 19790c77898SPeter Brune xm = M/m + ((M % m) > i); 19890c77898SPeter Brune ym = N/n + ((N % n) > j); 19990c77898SPeter Brune } else if (dim == 3) { 20090c77898SPeter Brune xm = M/m + ((M % m) > i); 20190c77898SPeter Brune ym = N/n + ((N % n) > j); 20290c77898SPeter Brune zm = P/p + ((P % p) > k); 20390c77898SPeter Brune } 20490c77898SPeter Brune 20590c77898SPeter Brune xsize = xm; 20690c77898SPeter Brune ysize = ym; 20790c77898SPeter Brune zsize = zm; 20890c77898SPeter Brune xo = xs; 20990c77898SPeter Brune yo = ys; 21090c77898SPeter Brune zo = zs; 21190c77898SPeter Brune 21290c77898SPeter Brune ierr = DMDACreate(PETSC_COMM_SELF,&(da[idx]));CHKERRQ(ierr); 21390c77898SPeter Brune ierr = DMSetOptionsPrefix(da[idx],"sub_");CHKERRQ(ierr); 21490c77898SPeter Brune ierr = DMDASetDim(da[idx], info.dim);CHKERRQ(ierr); 21590c77898SPeter Brune ierr = DMDASetDof(da[idx], info.dof);CHKERRQ(ierr); 21690c77898SPeter Brune 21790c77898SPeter Brune ierr = DMDASetStencilType(da[idx],info.st);CHKERRQ(ierr); 21890c77898SPeter Brune ierr = DMDASetStencilWidth(da[idx],info.sw);CHKERRQ(ierr); 21990c77898SPeter Brune 22090c77898SPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (xs != 0)) { 2217ddda789SPeter Brune xsize += xol; 2227ddda789SPeter Brune xo -= xol; 223e30e807fSPeter Brune } 22490c77898SPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (ys != 0)) { 2257ddda789SPeter Brune ysize += yol; 2267ddda789SPeter Brune yo -= yol; 227e30e807fSPeter Brune } 22890c77898SPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (zs != 0)) { 2297ddda789SPeter Brune zsize += zol; 2307ddda789SPeter Brune zo -= zol; 231e30e807fSPeter Brune } 232e30e807fSPeter Brune 23390c77898SPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (xs+xm != info.mx)) xsize += xol; 23490c77898SPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (ys+ym != info.my)) ysize += yol; 23590c77898SPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (zs+zm != info.mz)) zsize += zol; 2368865f1eaSKarl Rupp 23790c77898SPeter Brune if (info.bx != DMDA_BOUNDARY_PERIODIC) { 23890c77898SPeter Brune if (xo < 0) { 23990c77898SPeter Brune xsize += xo; 24090c77898SPeter Brune xo = 0; 24190c77898SPeter Brune } 24290c77898SPeter Brune if (xo+xsize > info.mx-1) { 24390c77898SPeter Brune xsize -= xo+xsize - info.mx; 24490c77898SPeter Brune } 24590c77898SPeter Brune } 24690c77898SPeter Brune if (info.by != DMDA_BOUNDARY_PERIODIC) { 24790c77898SPeter Brune if (yo < 0) { 24890c77898SPeter Brune ysize += yo; 24990c77898SPeter Brune yo = 0; 25090c77898SPeter Brune } 25190c77898SPeter Brune if (yo+ysize > info.my-1) { 25290c77898SPeter Brune ysize -= yo+ysize - info.my; 25390c77898SPeter Brune } 25490c77898SPeter Brune } 25590c77898SPeter Brune if (info.bz != DMDA_BOUNDARY_PERIODIC) { 25690c77898SPeter Brune if (zo < 0) { 25790c77898SPeter Brune zsize += zo; 25890c77898SPeter Brune zo = 0; 25990c77898SPeter Brune } 26090c77898SPeter Brune if (zo+zsize > info.mz-1) { 26190c77898SPeter Brune zsize -= zo+zsize - info.mz; 26290c77898SPeter Brune } 26390c77898SPeter Brune } 26490c77898SPeter Brune 26590c77898SPeter Brune ierr = DMDASetSizes(da[idx], xsize, ysize, zsize);CHKERRQ(ierr); 26690c77898SPeter Brune ierr = DMDASetNumProcs(da[idx], 1, 1, 1);CHKERRQ(ierr); 26790c77898SPeter Brune ierr = DMDASetBoundaryType(da[idx], DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 268e30e807fSPeter Brune 269e30e807fSPeter Brune /* set up as a block instead */ 27090c77898SPeter Brune ierr = DMSetUp(da[idx]);CHKERRQ(ierr); 271e30e807fSPeter Brune 27240234c92SPeter Brune /* nonoverlapping region */ 27390c77898SPeter Brune ierr = DMDASetNonOverlappingRegion(da[idx],xs,ys,zs,xm,ym,zm);CHKERRQ(ierr); 27440234c92SPeter Brune 27595c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 27690c77898SPeter Brune ierr = DMDASetOffset(da[idx],xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr); 27790c77898SPeter Brune xs += xm; 27890c77898SPeter Brune idx++; 27990c77898SPeter Brune } 28090c77898SPeter Brune ys += ym; 28190c77898SPeter Brune } 28290c77898SPeter Brune zs += zm; 28390c77898SPeter Brune } 28490c77898SPeter Brune *sdm = da; 285e30e807fSPeter Brune PetscFunctionReturn(0); 286e30e807fSPeter Brune } 287e30e807fSPeter Brune 288e30e807fSPeter Brune #undef __FUNCT__ 289e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 290e30e807fSPeter Brune /* 291e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 292e30e807fSPeter Brune 293285ae305SPeter Brune Right now this assumes one subdomain per processor. 294285ae305SPeter Brune 295e30e807fSPeter Brune */ 2960adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) 2970adebc6cSBarry Smith { 298e30e807fSPeter Brune PetscErrorCode ierr; 299285ae305SPeter Brune DMDALocalInfo info,subinfo; 300e30e807fSPeter Brune DM subdm; 301285ae305SPeter Brune MatStencil upper,lower; 302285ae305SPeter Brune IS idis,isis,odis,osis,gdis; 303285ae305SPeter Brune Vec svec,dvec,slvec; 30440234c92SPeter Brune PetscInt xm,ym,zm,xs,ys,zs; 30590c77898SPeter Brune PetscInt i; 306e30e807fSPeter Brune 307e30e807fSPeter Brune PetscFunctionBegin; 308285ae305SPeter Brune 309e30e807fSPeter Brune /* allocate the arrays of scatters */ 31090c77898SPeter Brune if (iscat) {ierr = PetscMalloc(nsubdms*sizeof(VecScatter*),iscat);CHKERRQ(ierr);} 31190c77898SPeter Brune if (oscat) {ierr = PetscMalloc(nsubdms*sizeof(VecScatter*),oscat);CHKERRQ(ierr);} 31290c77898SPeter Brune if (lscat) {ierr = PetscMalloc(nsubdms*sizeof(VecScatter*),lscat);CHKERRQ(ierr);} 313e30e807fSPeter Brune 314285ae305SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 31590c77898SPeter Brune for (i = 0; i < nsubdms; i++) { 31690c77898SPeter Brune subdm = subdms[i]; 317285ae305SPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 31890c77898SPeter Brune ierr = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 319e30e807fSPeter Brune 320285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 32140234c92SPeter Brune lower.i = xs; 32240234c92SPeter Brune lower.j = ys; 32340234c92SPeter Brune lower.k = zs; 32440234c92SPeter Brune upper.i = xs+xm; 32540234c92SPeter Brune upper.j = ys+ym; 32640234c92SPeter Brune upper.k = zs+zm; 327285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); 328285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); 329e30e807fSPeter Brune 330285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 331285ae305SPeter Brune lower.i = subinfo.xs; 332285ae305SPeter Brune lower.j = subinfo.ys; 333285ae305SPeter Brune lower.k = subinfo.zs; 334285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 335285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 336285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 337285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); 338285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); 339e30e807fSPeter Brune 340285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 341285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 342285ae305SPeter Brune lower.i = subinfo.gxs; 343285ae305SPeter Brune lower.j = subinfo.gys; 344285ae305SPeter Brune lower.k = subinfo.gzs; 345285ae305SPeter Brune upper.i = subinfo.gxs+subinfo.gxm; 346285ae305SPeter Brune upper.j = subinfo.gys+subinfo.gym; 347285ae305SPeter Brune upper.k = subinfo.gzs+subinfo.gzm; 348e30e807fSPeter Brune 349285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); 350e30e807fSPeter Brune 351e30e807fSPeter Brune /* form the scatter */ 352e30e807fSPeter Brune ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 353e30e807fSPeter Brune ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 354e30e807fSPeter Brune ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 355e30e807fSPeter Brune 35690c77898SPeter Brune if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[i]);CHKERRQ(ierr);} 35790c77898SPeter Brune if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[i]);CHKERRQ(ierr);} 35890c77898SPeter Brune if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[i]);CHKERRQ(ierr);} 359e30e807fSPeter Brune 360e30e807fSPeter Brune ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 361e30e807fSPeter Brune ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 362e30e807fSPeter Brune ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 363e30e807fSPeter Brune 364e30e807fSPeter Brune ierr = ISDestroy(&idis);CHKERRQ(ierr); 365e30e807fSPeter Brune ierr = ISDestroy(&isis);CHKERRQ(ierr); 366e30e807fSPeter Brune 367e30e807fSPeter Brune ierr = ISDestroy(&odis);CHKERRQ(ierr); 368e30e807fSPeter Brune ierr = ISDestroy(&osis);CHKERRQ(ierr); 369e30e807fSPeter Brune 370e30e807fSPeter Brune ierr = ISDestroy(&gdis);CHKERRQ(ierr); 37190c77898SPeter Brune } 372e30e807fSPeter Brune PetscFunctionReturn(0); 373e30e807fSPeter Brune } 374e30e807fSPeter Brune 375e30e807fSPeter Brune #undef __FUNCT__ 376e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private" 37790c77898SPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,PetscInt n,DM *subdm,IS **iis,IS **ois) 3780adebc6cSBarry Smith { 379e30e807fSPeter Brune PetscErrorCode ierr; 38090c77898SPeter Brune PetscInt i; 381e30e807fSPeter Brune DMDALocalInfo info,subinfo; 382285ae305SPeter Brune MatStencil lower,upper; 383e30e807fSPeter Brune 384e30e807fSPeter Brune PetscFunctionBegin; 385e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 38690c77898SPeter Brune if (iis) {ierr = PetscMalloc(n*sizeof(IS*),iis);CHKERRQ(ierr);} 38790c77898SPeter Brune if (ois) {ierr = PetscMalloc(n*sizeof(IS*),ois);CHKERRQ(ierr);} 388e30e807fSPeter Brune 38990c77898SPeter Brune for (i = 0;i < n; i++) { 39090c77898SPeter Brune ierr = DMDAGetLocalInfo(subdm[i],&subinfo);CHKERRQ(ierr); 39190c77898SPeter Brune if (iis) { 392285ae305SPeter Brune /* create the inner IS */ 393285ae305SPeter Brune lower.i = info.xs; 394285ae305SPeter Brune lower.j = info.ys; 395285ae305SPeter Brune lower.k = info.zs; 396285ae305SPeter Brune upper.i = info.xs+info.xm; 397285ae305SPeter Brune upper.j = info.ys+info.ym; 398285ae305SPeter Brune upper.k = info.zs+info.zm; 39990c77898SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*iis)[i]);CHKERRQ(ierr); 40090c77898SPeter Brune } 401e30e807fSPeter Brune 40290c77898SPeter Brune if (ois) { 403285ae305SPeter Brune /* create the outer IS */ 404285ae305SPeter Brune lower.i = subinfo.xs; 405285ae305SPeter Brune lower.j = subinfo.ys; 406285ae305SPeter Brune lower.k = subinfo.zs; 407285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 408285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 409285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 41090c77898SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*ois)[i]);CHKERRQ(ierr); 41190c77898SPeter Brune } 41290c77898SPeter Brune } 413e30e807fSPeter Brune PetscFunctionReturn(0); 414e30e807fSPeter Brune } 415e30e807fSPeter Brune 416e30e807fSPeter Brune #undef __FUNCT__ 417e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA" 4180adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) 4190adebc6cSBarry Smith { 420e30e807fSPeter Brune PetscErrorCode ierr; 42190c77898SPeter Brune DM *sdm; 42290c77898SPeter Brune PetscInt n,i; 4230a696f66SPeter Brune 4240adebc6cSBarry Smith PetscFunctionBegin; 42590c77898SPeter Brune ierr = DMDASubDomainDA_Private(dm,&n,&sdm);CHKERRQ(ierr); 42690c77898SPeter Brune if (names) { 42790c77898SPeter Brune ierr = PetscMalloc(n*sizeof(char*),names);CHKERRQ(ierr); 428*e2c616c8SPeter Brune for (i=0;i<n;i++) (*names)[i] = 0; 429e30e807fSPeter Brune } 43090c77898SPeter Brune ierr = DMDASubDomainIS_Private(dm,n,sdm,iis,ois);CHKERRQ(ierr); 43190c77898SPeter Brune if (subdm) *subdm = sdm; 432*e2c616c8SPeter Brune else { 433*e2c616c8SPeter Brune for (i=0;i<n;i++) { 434*e2c616c8SPeter Brune ierr = DMDestroy(&sdm[i]);CHKERRQ(ierr); 435*e2c616c8SPeter Brune } 436*e2c616c8SPeter Brune } 43790c77898SPeter Brune if (len) *len = n; 438e30e807fSPeter Brune PetscFunctionReturn(0); 439e30e807fSPeter Brune } 440