1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3063654ddSPeter Brune /*@ 4063654ddSPeter Brune DMDACreatePatchIS - Creates an index set corresponding to a patch of the DA. 5ff9846d9SPeter Brune 6063654ddSPeter Brune Not Collective 7063654ddSPeter Brune 8063654ddSPeter Brune Input Parameters: 9063654ddSPeter Brune + da - the DMDA 10063654ddSPeter Brune . lower - a matstencil with i, j and k corresponding to the lower corner of the patch 11063654ddSPeter Brune - upper - a matstencil with i, j and k corresponding to the upper corner of the patch 12063654ddSPeter Brune 13063654ddSPeter Brune Output Parameters: 14063654ddSPeter Brune . is - the IS corresponding to the patch 15063654ddSPeter Brune 16063654ddSPeter Brune Level: developer 17063654ddSPeter Brune 18063654ddSPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDACreateDomainDecompositionScatters() 19063654ddSPeter Brune @*/ 20ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) 21ff9846d9SPeter Brune { 22063654ddSPeter Brune PetscInt ms=0,ns=0,ps=0; 23063654ddSPeter Brune PetscInt me=1,ne=1,pe=1; 24063654ddSPeter Brune PetscInt mr=0,nr=0,pr=0; 25ff9846d9SPeter Brune PetscInt ii,jj,kk; 26063654ddSPeter Brune PetscInt si,sj,sk; 27063654ddSPeter Brune PetscInt i,j,k,l,idx; 28ff9846d9SPeter Brune PetscInt base; 29063654ddSPeter Brune PetscInt xm=1,ym=1,zm=1; 30063654ddSPeter Brune const PetscInt *lx,*ly,*lz; 31ff9846d9SPeter Brune PetscInt ox,oy,oz; 32063654ddSPeter Brune PetscInt m,n,p,M,N,P,dof; 33063654ddSPeter Brune PetscInt nindices; 34063654ddSPeter Brune PetscInt *indices; 35063654ddSPeter Brune DM_DA *dd = (DM_DA*)da->data; 36063654ddSPeter Brune PetscErrorCode ierr; 37ff9846d9SPeter Brune 380adebc6cSBarry Smith PetscFunctionBegin; 39063654ddSPeter Brune /* need to get the sizes of the actual DM rather than the "global" space of a subdomain DM */ 40063654ddSPeter Brune M = dd->M;N = dd->N;P=dd->P; 41063654ddSPeter Brune m = dd->m;n = dd->n;p=dd->p; 42063654ddSPeter Brune dof = dd->w; 430298fd71SBarry Smith ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr); 44063654ddSPeter Brune ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 45063654ddSPeter Brune nindices = (upper->i - lower->i)*(upper->j - lower->j)*(upper->k - lower->k)*dof; 46854ce69bSBarry Smith ierr = PetscMalloc1(nindices,&indices);CHKERRQ(ierr); 47063654ddSPeter Brune /* start at index 0 on processor 0 */ 48063654ddSPeter Brune mr = 0; 49063654ddSPeter Brune nr = 0; 50063654ddSPeter Brune pr = 0; 51063654ddSPeter Brune ms = 0; 52063654ddSPeter Brune ns = 0; 53063654ddSPeter Brune ps = 0; 54063654ddSPeter Brune if (lx) me = lx[0]; 55063654ddSPeter Brune if (ly) ne = ly[0]; 56063654ddSPeter Brune if (lz) pe = lz[0]; 57ff9846d9SPeter Brune idx = 0; 58ff9846d9SPeter Brune for (k=lower->k-oz;k<upper->k-oz;k++) { 59ff9846d9SPeter Brune for (j=lower->j-oy;j < upper->j-oy;j++) { 60063654ddSPeter Brune for (i=lower->i-ox;i < upper->i-ox;i++) { 61063654ddSPeter Brune /* "actual" indices rather than ones outside of the domain */ 62063654ddSPeter Brune ii = i; 63063654ddSPeter Brune jj = j; 64063654ddSPeter Brune kk = k; 65063654ddSPeter Brune if (ii < 0) ii = M + ii; 66063654ddSPeter Brune if (jj < 0) jj = N + jj; 67063654ddSPeter Brune if (kk < 0) kk = P + kk; 68063654ddSPeter Brune if (ii > M-1) ii = ii - M; 69063654ddSPeter Brune if (jj > N-1) jj = jj - N; 70063654ddSPeter Brune if (kk > P-1) kk = kk - P; 71063654ddSPeter Brune /* gone out of processor range on x axis */ 72063654ddSPeter Brune while(ii > me-1 || ii < ms) { 73063654ddSPeter Brune if (mr == m-1) { 74063654ddSPeter Brune ms = 0; 75063654ddSPeter Brune me = lx[0]; 76063654ddSPeter Brune mr = 0; 77063654ddSPeter Brune } else { 78063654ddSPeter Brune mr++; 79063654ddSPeter Brune ms = me; 80063654ddSPeter Brune me += lx[mr]; 81063654ddSPeter Brune } 82063654ddSPeter Brune } 83063654ddSPeter Brune /* gone out of processor range on y axis */ 84063654ddSPeter Brune while(jj > ne-1 || jj < ns) { 85063654ddSPeter Brune if (nr == n-1) { 86063654ddSPeter Brune ns = 0; 87063654ddSPeter Brune ne = ly[0]; 88063654ddSPeter Brune nr = 0; 89063654ddSPeter Brune } else { 90063654ddSPeter Brune nr++; 91063654ddSPeter Brune ns = ne; 92063654ddSPeter Brune ne += ly[nr]; 93063654ddSPeter Brune } 94063654ddSPeter Brune } 95063654ddSPeter Brune /* gone out of processor range on z axis */ 96063654ddSPeter Brune while(kk > pe-1 || kk < ps) { 97063654ddSPeter Brune if (pr == p-1) { 98063654ddSPeter Brune ps = 0; 99063654ddSPeter Brune pe = lz[0]; 100063654ddSPeter Brune pr = 0; 101063654ddSPeter Brune } else { 102063654ddSPeter Brune pr++; 103063654ddSPeter Brune ps = pe; 104063654ddSPeter Brune pe += lz[pr]; 105063654ddSPeter Brune } 106063654ddSPeter Brune } 107063654ddSPeter Brune /* compute the vector base on owning processor */ 108063654ddSPeter Brune xm = me - ms; 109063654ddSPeter Brune ym = ne - ns; 110063654ddSPeter Brune zm = pe - ps; 111b0bff951SPeter Brune base = ms*ym*zm + ns*M*zm + ps*M*N; 112063654ddSPeter Brune /* compute the local coordinates on owning processor */ 113063654ddSPeter Brune si = ii - ms; 114063654ddSPeter Brune sj = jj - ns; 115063654ddSPeter Brune sk = kk - ps; 116063654ddSPeter Brune for (l=0;l<dof;l++) { 117063654ddSPeter Brune indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk); 118ff9846d9SPeter Brune idx++; 119ff9846d9SPeter Brune } 120ff9846d9SPeter Brune } 121ff9846d9SPeter Brune } 122063654ddSPeter Brune } 123063654ddSPeter Brune ISCreateGeneral(PETSC_COMM_SELF,idx,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr); 124ff9846d9SPeter Brune PetscFunctionReturn(0); 125ff9846d9SPeter Brune } 126ff9846d9SPeter Brune 12790c77898SPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) 1280adebc6cSBarry Smith { 12990c77898SPeter Brune DM *da; 13090c77898SPeter Brune PetscInt dim,size,i,j,k,idx; 131e30e807fSPeter Brune PetscErrorCode ierr; 132e30e807fSPeter Brune DMDALocalInfo info; 133ff9846d9SPeter Brune PetscInt xsize,ysize,zsize; 134e30e807fSPeter Brune PetscInt xo,yo,zo; 13590c77898SPeter Brune PetscInt xs,ys,zs; 13690c77898SPeter Brune PetscInt xm=1,ym=1,zm=1; 1377ddda789SPeter Brune PetscInt xol,yol,zol; 13890c77898SPeter Brune PetscInt m=1,n=1,p=1; 13990c77898SPeter Brune PetscInt M,N,P; 14090c77898SPeter Brune PetscInt pm,mtmp; 141e30e807fSPeter Brune 142e30e807fSPeter Brune PetscFunctionBegin; 143e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1447ddda789SPeter Brune ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr); 14590c77898SPeter Brune ierr = DMDAGetNumLocalSubDomains(dm,&size);CHKERRQ(ierr); 146785e854fSJed Brown ierr = PetscMalloc1(size,&da);CHKERRQ(ierr); 1477ddda789SPeter Brune 14890c77898SPeter Brune if (nlocal) *nlocal = size; 149e30e807fSPeter Brune 15090c77898SPeter Brune dim = info.dim; 151e30e807fSPeter Brune 15290c77898SPeter Brune M = info.xm; 15390c77898SPeter Brune N = info.ym; 15490c77898SPeter Brune P = info.zm; 15590c77898SPeter Brune 15690c77898SPeter Brune if (dim == 1) { 15790c77898SPeter Brune m = size; 15890c77898SPeter Brune } else if (dim == 2) { 15990c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); 16090c77898SPeter Brune while (m > 0) { 16190c77898SPeter Brune n = size/m; 16290c77898SPeter Brune if (m*n*p == size) break; 16390c77898SPeter Brune m--; 16490c77898SPeter Brune } 16590c77898SPeter Brune } else if (dim == 3) { 16690c77898SPeter Brune n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); if (!n) n = 1; 16790c77898SPeter Brune while (n > 0) { 16890c77898SPeter Brune pm = size/n; 16990c77898SPeter Brune if (n*pm == size) break; 17090c77898SPeter Brune n--; 17190c77898SPeter Brune } 17290c77898SPeter Brune if (!n) n = 1; 17390c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); 17490c77898SPeter Brune if (!m) m = 1; 17590c77898SPeter Brune while (m > 0) { 17690c77898SPeter Brune p = size/(m*n); 17790c77898SPeter Brune if (m*n*p == size) break; 17890c77898SPeter Brune m--; 17990c77898SPeter Brune } 18090c77898SPeter Brune if (M > P && m < p) {mtmp = m; m = p; p = mtmp;} 18190c77898SPeter Brune } 18290c77898SPeter Brune 18390c77898SPeter Brune zs = info.zs; 18490c77898SPeter Brune idx = 0; 18590c77898SPeter Brune for (k = 0; k < p; k++) { 18690c77898SPeter Brune ys = info.ys; 18790c77898SPeter Brune for (j = 0; j < n; j++) { 18890c77898SPeter Brune xs = info.xs; 18990c77898SPeter Brune for (i = 0; i < m; i++) { 19090c77898SPeter Brune if (dim == 1) { 19190c77898SPeter Brune xm = M/m + ((M % m) > i); 19290c77898SPeter Brune } else if (dim == 2) { 19390c77898SPeter Brune xm = M/m + ((M % m) > i); 19490c77898SPeter Brune ym = N/n + ((N % n) > j); 19590c77898SPeter Brune } else if (dim == 3) { 19690c77898SPeter Brune xm = M/m + ((M % m) > i); 19790c77898SPeter Brune ym = N/n + ((N % n) > j); 19890c77898SPeter Brune zm = P/p + ((P % p) > k); 19990c77898SPeter Brune } 20090c77898SPeter Brune 20190c77898SPeter Brune xsize = xm; 20290c77898SPeter Brune ysize = ym; 20390c77898SPeter Brune zsize = zm; 20490c77898SPeter Brune xo = xs; 20590c77898SPeter Brune yo = ys; 20690c77898SPeter Brune zo = zs; 20790c77898SPeter Brune 20890c77898SPeter Brune ierr = DMDACreate(PETSC_COMM_SELF,&(da[idx]));CHKERRQ(ierr); 20990c77898SPeter Brune ierr = DMSetOptionsPrefix(da[idx],"sub_");CHKERRQ(ierr); 210c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da[idx], info.dim);CHKERRQ(ierr); 21190c77898SPeter Brune ierr = DMDASetDof(da[idx], info.dof);CHKERRQ(ierr); 21290c77898SPeter Brune 21390c77898SPeter Brune ierr = DMDASetStencilType(da[idx],info.st);CHKERRQ(ierr); 21490c77898SPeter Brune ierr = DMDASetStencilWidth(da[idx],info.sw);CHKERRQ(ierr); 21590c77898SPeter Brune 216bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) { 2177ddda789SPeter Brune xsize += xol; 2187ddda789SPeter Brune xo -= xol; 219e30e807fSPeter Brune } 220bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) { 2217ddda789SPeter Brune ysize += yol; 2227ddda789SPeter Brune yo -= yol; 223e30e807fSPeter Brune } 224bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) { 2257ddda789SPeter Brune zsize += zol; 2267ddda789SPeter Brune zo -= zol; 227e30e807fSPeter Brune } 228e30e807fSPeter Brune 229bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs+xm != info.mx)) xsize += xol; 230bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys+ym != info.my)) ysize += yol; 231bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs+zm != info.mz)) zsize += zol; 2328865f1eaSKarl Rupp 233bff4a2f0SMatthew G. Knepley if (info.bx != DM_BOUNDARY_PERIODIC) { 23490c77898SPeter Brune if (xo < 0) { 23590c77898SPeter Brune xsize += xo; 23690c77898SPeter Brune xo = 0; 23790c77898SPeter Brune } 23890c77898SPeter Brune if (xo+xsize > info.mx-1) { 23990c77898SPeter Brune xsize -= xo+xsize - info.mx; 24090c77898SPeter Brune } 24190c77898SPeter Brune } 242bff4a2f0SMatthew G. Knepley if (info.by != DM_BOUNDARY_PERIODIC) { 24390c77898SPeter Brune if (yo < 0) { 24490c77898SPeter Brune ysize += yo; 24590c77898SPeter Brune yo = 0; 24690c77898SPeter Brune } 24790c77898SPeter Brune if (yo+ysize > info.my-1) { 24890c77898SPeter Brune ysize -= yo+ysize - info.my; 24990c77898SPeter Brune } 25090c77898SPeter Brune } 251bff4a2f0SMatthew G. Knepley if (info.bz != DM_BOUNDARY_PERIODIC) { 25290c77898SPeter Brune if (zo < 0) { 25390c77898SPeter Brune zsize += zo; 25490c77898SPeter Brune zo = 0; 25590c77898SPeter Brune } 25690c77898SPeter Brune if (zo+zsize > info.mz-1) { 25790c77898SPeter Brune zsize -= zo+zsize - info.mz; 25890c77898SPeter Brune } 25990c77898SPeter Brune } 26090c77898SPeter Brune 26190c77898SPeter Brune ierr = DMDASetSizes(da[idx], xsize, ysize, zsize);CHKERRQ(ierr); 26290c77898SPeter Brune ierr = DMDASetNumProcs(da[idx], 1, 1, 1);CHKERRQ(ierr); 263bff4a2f0SMatthew G. Knepley ierr = DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED);CHKERRQ(ierr); 264e30e807fSPeter Brune 265e30e807fSPeter Brune /* set up as a block instead */ 26690c77898SPeter Brune ierr = DMSetUp(da[idx]);CHKERRQ(ierr); 267e30e807fSPeter Brune 26840234c92SPeter Brune /* nonoverlapping region */ 26990c77898SPeter Brune ierr = DMDASetNonOverlappingRegion(da[idx],xs,ys,zs,xm,ym,zm);CHKERRQ(ierr); 27040234c92SPeter Brune 27195c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 27290c77898SPeter Brune ierr = DMDASetOffset(da[idx],xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr); 27390c77898SPeter Brune xs += xm; 27490c77898SPeter Brune idx++; 27590c77898SPeter Brune } 27690c77898SPeter Brune ys += ym; 27790c77898SPeter Brune } 27890c77898SPeter Brune zs += zm; 27990c77898SPeter Brune } 28090c77898SPeter Brune *sdm = da; 281e30e807fSPeter Brune PetscFunctionReturn(0); 282e30e807fSPeter Brune } 283e30e807fSPeter Brune 284e30e807fSPeter Brune /* 285e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 286e30e807fSPeter Brune 287285ae305SPeter Brune Right now this assumes one subdomain per processor. 288285ae305SPeter Brune 289e30e807fSPeter Brune */ 2900adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) 2910adebc6cSBarry Smith { 292e30e807fSPeter Brune PetscErrorCode ierr; 293285ae305SPeter Brune DMDALocalInfo info,subinfo; 294e30e807fSPeter Brune DM subdm; 295285ae305SPeter Brune MatStencil upper,lower; 296285ae305SPeter Brune IS idis,isis,odis,osis,gdis; 297285ae305SPeter Brune Vec svec,dvec,slvec; 29840234c92SPeter Brune PetscInt xm,ym,zm,xs,ys,zs; 29990c77898SPeter Brune PetscInt i; 300e30e807fSPeter Brune 301e30e807fSPeter Brune PetscFunctionBegin; 302285ae305SPeter Brune 303e30e807fSPeter Brune /* allocate the arrays of scatters */ 304785e854fSJed Brown if (iscat) {ierr = PetscMalloc1(nsubdms,iscat);CHKERRQ(ierr);} 305785e854fSJed Brown if (oscat) {ierr = PetscMalloc1(nsubdms,oscat);CHKERRQ(ierr);} 306785e854fSJed Brown if (lscat) {ierr = PetscMalloc1(nsubdms,lscat);CHKERRQ(ierr);} 307e30e807fSPeter Brune 308285ae305SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 30990c77898SPeter Brune for (i = 0; i < nsubdms; i++) { 31090c77898SPeter Brune subdm = subdms[i]; 311285ae305SPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 31290c77898SPeter Brune ierr = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 313e30e807fSPeter Brune 314285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 31540234c92SPeter Brune lower.i = xs; 31640234c92SPeter Brune lower.j = ys; 31740234c92SPeter Brune lower.k = zs; 31840234c92SPeter Brune upper.i = xs+xm; 31940234c92SPeter Brune upper.j = ys+ym; 32040234c92SPeter Brune upper.k = zs+zm; 321285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); 322285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); 323e30e807fSPeter Brune 324285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 325285ae305SPeter Brune lower.i = subinfo.xs; 326285ae305SPeter Brune lower.j = subinfo.ys; 327285ae305SPeter Brune lower.k = subinfo.zs; 328285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 329285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 330285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 331285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); 332285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); 333e30e807fSPeter Brune 334285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 335285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 336285ae305SPeter Brune lower.i = subinfo.gxs; 337285ae305SPeter Brune lower.j = subinfo.gys; 338285ae305SPeter Brune lower.k = subinfo.gzs; 339285ae305SPeter Brune upper.i = subinfo.gxs+subinfo.gxm; 340285ae305SPeter Brune upper.j = subinfo.gys+subinfo.gym; 341285ae305SPeter Brune upper.k = subinfo.gzs+subinfo.gzm; 342e30e807fSPeter Brune 343285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); 344e30e807fSPeter Brune 345e30e807fSPeter Brune /* form the scatter */ 346e30e807fSPeter Brune ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 347e30e807fSPeter Brune ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 348e30e807fSPeter Brune ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 349e30e807fSPeter Brune 3509448b7f1SJunchao Zhang if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[i]);CHKERRQ(ierr);} 3519448b7f1SJunchao Zhang if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[i]);CHKERRQ(ierr);} 3529448b7f1SJunchao Zhang if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[i]);CHKERRQ(ierr);} 353e30e807fSPeter Brune 354e30e807fSPeter Brune ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 355e30e807fSPeter Brune ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 356e30e807fSPeter Brune ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 357e30e807fSPeter Brune 358e30e807fSPeter Brune ierr = ISDestroy(&idis);CHKERRQ(ierr); 359e30e807fSPeter Brune ierr = ISDestroy(&isis);CHKERRQ(ierr); 360e30e807fSPeter Brune 361e30e807fSPeter Brune ierr = ISDestroy(&odis);CHKERRQ(ierr); 362e30e807fSPeter Brune ierr = ISDestroy(&osis);CHKERRQ(ierr); 363e30e807fSPeter Brune 364e30e807fSPeter Brune ierr = ISDestroy(&gdis);CHKERRQ(ierr); 36590c77898SPeter Brune } 366e30e807fSPeter Brune PetscFunctionReturn(0); 367e30e807fSPeter Brune } 368e30e807fSPeter Brune 36990c77898SPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,PetscInt n,DM *subdm,IS **iis,IS **ois) 3700adebc6cSBarry Smith { 371e30e807fSPeter Brune PetscErrorCode ierr; 37290c77898SPeter Brune PetscInt i; 373e30e807fSPeter Brune DMDALocalInfo info,subinfo; 374285ae305SPeter Brune MatStencil lower,upper; 375e30e807fSPeter Brune 376e30e807fSPeter Brune PetscFunctionBegin; 377e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 378785e854fSJed Brown if (iis) {ierr = PetscMalloc1(n,iis);CHKERRQ(ierr);} 379785e854fSJed Brown if (ois) {ierr = PetscMalloc1(n,ois);CHKERRQ(ierr);} 380e30e807fSPeter Brune 38190c77898SPeter Brune for (i = 0;i < n; i++) { 38290c77898SPeter Brune ierr = DMDAGetLocalInfo(subdm[i],&subinfo);CHKERRQ(ierr); 38390c77898SPeter Brune if (iis) { 384285ae305SPeter Brune /* create the inner IS */ 385285ae305SPeter Brune lower.i = info.xs; 386285ae305SPeter Brune lower.j = info.ys; 387285ae305SPeter Brune lower.k = info.zs; 388285ae305SPeter Brune upper.i = info.xs+info.xm; 389285ae305SPeter Brune upper.j = info.ys+info.ym; 390285ae305SPeter Brune upper.k = info.zs+info.zm; 39190c77898SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*iis)[i]);CHKERRQ(ierr); 39290c77898SPeter Brune } 393e30e807fSPeter Brune 39490c77898SPeter Brune if (ois) { 395285ae305SPeter Brune /* create the outer IS */ 396285ae305SPeter Brune lower.i = subinfo.xs; 397285ae305SPeter Brune lower.j = subinfo.ys; 398285ae305SPeter Brune lower.k = subinfo.zs; 399285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 400285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 401285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 40290c77898SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*ois)[i]);CHKERRQ(ierr); 40390c77898SPeter Brune } 40490c77898SPeter Brune } 405e30e807fSPeter Brune PetscFunctionReturn(0); 406e30e807fSPeter Brune } 407e30e807fSPeter Brune 4080adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) 4090adebc6cSBarry Smith { 410e30e807fSPeter Brune PetscErrorCode ierr; 41190c77898SPeter Brune DM *sdm; 41290c77898SPeter Brune PetscInt n,i; 4130a696f66SPeter Brune 4140adebc6cSBarry Smith PetscFunctionBegin; 41590c77898SPeter Brune ierr = DMDASubDomainDA_Private(dm,&n,&sdm);CHKERRQ(ierr); 41690c77898SPeter Brune if (names) { 417785e854fSJed Brown ierr = PetscMalloc1(n,names);CHKERRQ(ierr); 418*ea78f98cSLisandro Dalcin for (i=0;i<n;i++) (*names)[i] = NULL; 419e30e807fSPeter Brune } 42090c77898SPeter Brune ierr = DMDASubDomainIS_Private(dm,n,sdm,iis,ois);CHKERRQ(ierr); 42190c77898SPeter Brune if (subdm) *subdm = sdm; 422e2c616c8SPeter Brune else { 423e2c616c8SPeter Brune for (i=0;i<n;i++) { 424e2c616c8SPeter Brune ierr = DMDestroy(&sdm[i]);CHKERRQ(ierr); 425e2c616c8SPeter Brune } 426e2c616c8SPeter Brune } 42790c77898SPeter Brune if (len) *len = n; 428e30e807fSPeter Brune PetscFunctionReturn(0); 429e30e807fSPeter Brune } 430