diff --git a/src/dm/impls/da/dadd.c b/src/dm/impls/da/dadd.c index df324f2..b9fc079 100644 --- a/src/dm/impls/da/dadd.c +++ b/src/dm/impls/da/dadd.c @@ -299,10 +299,11 @@ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM DMDALocalInfo info,subinfo; DM subdm; MatStencil upper,lower; - IS idis,isis,odis,osis,gdis; + IS idis,isis,odis,osis,gdis,gsdis; Vec svec,dvec,slvec; PetscInt xm,ym,zm,xs,ys,zs; PetscInt i; + DM_DA *da = (DM_DA*)dm->data; PetscFunctionBegin; @@ -315,7 +316,7 @@ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM 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); + 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; @@ -340,13 +341,19 @@ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM /* 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; + if (da->bx != DM_BOUNDARY_PERIODIC) lower.i = PetscMax(lower.i,0); lower.j = subinfo.gys; + if (da->by != DM_BOUNDARY_PERIODIC) lower.j = PetscMax(lower.j,0); lower.k = subinfo.gzs; + if (da->bz != DM_BOUNDARY_PERIODIC) lower.k = PetscMax(lower.k,0); upper.i = subinfo.gxs+subinfo.gxm; + if (da->bx != DM_BOUNDARY_PERIODIC) upper.i = PetscMin(upper.i,da->M-1); upper.j = subinfo.gys+subinfo.gym; + if (da->by != DM_BOUNDARY_PERIODIC) upper.j = PetscMin(upper.j,da->N-1); upper.k = subinfo.gzs+subinfo.gzm; - + if (da->bz != DM_BOUNDARY_PERIODIC) upper.k = PetscMin(upper.k,da->P-1); ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); + ierr = DMDACreatePatchIS(subdm,&lower,&upper,&gsdis);CHKERRQ(ierr); /* form the scatter */ ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); @@ -355,7 +362,7 @@ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM 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);} + if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,gsdis,&(*lscat)[i]);CHKERRQ(ierr);} ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); @@ -368,6 +375,7 @@ PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM ierr = ISDestroy(&osis);CHKERRQ(ierr); ierr = ISDestroy(&gdis);CHKERRQ(ierr); + ierr = ISDestroy(&gsdis);CHKERRQ(ierr); } PetscFunctionReturn(0); }