104dc8af9SSatish Balay #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3e30e807fSPeter Brune #undef __FUNCT__ 4ff9846d9SPeter Brune #define __FUNCT__ "DMDACreatePatchIS" 5ff9846d9SPeter Brune 6ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) 7ff9846d9SPeter Brune { 8ff9846d9SPeter Brune PetscErrorCode ierr; 9ff9846d9SPeter Brune PetscInt i,j,k,idx; 10ff9846d9SPeter Brune PetscInt ii,jj,kk; 11ff9846d9SPeter Brune Vec v; 12ff9846d9SPeter Brune PetscInt n,pn,bs; 13ff9846d9SPeter Brune PetscMPIInt rank; 14ff9846d9SPeter Brune PetscSF sf,psf; 15ff9846d9SPeter Brune PetscLayout map; 16ff9846d9SPeter Brune MPI_Comm comm; 17ff9846d9SPeter Brune PetscInt *natidx,*globidx,*leafidx; 18ff9846d9SPeter Brune PetscInt *pnatidx,*pleafidx; 19ff9846d9SPeter Brune PetscInt base; 20ff9846d9SPeter Brune PetscInt ox,oy,oz; 21ff9846d9SPeter Brune DM_DA *dd; 22ff9846d9SPeter Brune 230adebc6cSBarry Smith PetscFunctionBegin; 24*ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 25ff9846d9SPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 26ff9846d9SPeter Brune dd = (DM_DA*)da->data; 27ff9846d9SPeter Brune 28ff9846d9SPeter Brune /* construct the natural mapping */ 29ff9846d9SPeter Brune ierr = DMGetGlobalVector(da,&v);CHKERRQ(ierr); 30ff9846d9SPeter Brune ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 310298fd71SBarry Smith ierr = VecGetOwnershipRange(v,&base,NULL);CHKERRQ(ierr); 32ff9846d9SPeter Brune ierr = VecGetBlockSize(v,&bs);CHKERRQ(ierr); 33ff9846d9SPeter Brune ierr = DMRestoreGlobalVector(da,&v);CHKERRQ(ierr); 34ff9846d9SPeter Brune 35ff9846d9SPeter Brune /* construct the layout */ 36ff9846d9SPeter Brune ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr); 37ff9846d9SPeter Brune ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 38ff9846d9SPeter Brune ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr); 39ff9846d9SPeter Brune ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 40ff9846d9SPeter Brune 41ff9846d9SPeter Brune /* construct the list of natural indices on this process when PETSc ordering is considered */ 420298fd71SBarry Smith ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr); 43ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&natidx);CHKERRQ(ierr); 44ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&globidx);CHKERRQ(ierr); 45ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&leafidx);CHKERRQ(ierr); 46ff9846d9SPeter Brune idx = 0; 47ff9846d9SPeter Brune for (k=dd->zs; k<dd->ze; k++) { 48ff9846d9SPeter Brune for (j=dd->ys; j<dd->ye; j++) { 49ff9846d9SPeter Brune for (i=dd->xs; i<dd->xe; i++) { 50ff9846d9SPeter Brune natidx[idx] = i + dd->w*(j*dd->M + k*dd->M*dd->N); 51ff9846d9SPeter Brune globidx[idx] = base + idx; 52ff9846d9SPeter Brune leafidx[idx] = 0; 53ff9846d9SPeter Brune idx++; 54ff9846d9SPeter Brune } 55ff9846d9SPeter Brune } 56ff9846d9SPeter Brune } 57ff9846d9SPeter Brune 58ff9846d9SPeter Brune if (idx != n) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong."); 59ff9846d9SPeter Brune 60ff9846d9SPeter Brune /* construct the SF going from the natural indices to the local set of PETSc indices */ 61ff9846d9SPeter Brune ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 62ff9846d9SPeter Brune ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 630298fd71SBarry Smith ierr = PetscSFSetGraphLayout(sf,map,n,NULL,PETSC_OWN_POINTER,natidx);CHKERRQ(ierr); 64ff9846d9SPeter Brune 65ff9846d9SPeter Brune /* broadcast the global indices over to the corresponding natural indices */ 66ff9846d9SPeter Brune ierr = PetscSFGatherBegin(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 67ff9846d9SPeter Brune ierr = PetscSFGatherEnd(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 68ff9846d9SPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 69ff9846d9SPeter Brune 70ff9846d9SPeter Brune 71ff9846d9SPeter Brune pn = dd->w*(upper->k - lower->k)*(upper->j - lower->j)*(upper->i - lower->i); 72ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*pn,&pnatidx);CHKERRQ(ierr); 73ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*pn,&pleafidx);CHKERRQ(ierr); 74ff9846d9SPeter Brune idx = 0; 75ff9846d9SPeter Brune for (k=lower->k-oz; k<upper->k-oz; k++) { 76ff9846d9SPeter Brune for (j=lower->j-oy; j<upper->j-oy; j++) { 77ff9846d9SPeter Brune for (i=dd->w*(lower->i-ox); i<dd->w*(upper->i-ox); i++) { 78ff9846d9SPeter Brune ii = i % (dd->w*dd->M); 79ff9846d9SPeter Brune jj = j % dd->N; 80ff9846d9SPeter Brune kk = k % dd->P; 81ff9846d9SPeter Brune if (ii < 0) ii = dd->w*dd->M + ii; 82ff9846d9SPeter Brune if (jj < 0) jj = dd->N + jj; 83ff9846d9SPeter Brune if (kk < 0) kk = dd->P + kk; 84ff9846d9SPeter Brune pnatidx[idx] = ii + dd->w*(jj*dd->M + kk*dd->M*dd->N); 85ff9846d9SPeter Brune idx++; 86ff9846d9SPeter Brune } 87ff9846d9SPeter Brune } 88ff9846d9SPeter Brune } 89ff9846d9SPeter Brune 90ff9846d9SPeter Brune if (idx != pn) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong"); 91ff9846d9SPeter Brune 92ff9846d9SPeter Brune ierr = PetscSFCreate(comm,&psf);CHKERRQ(ierr); 93ff9846d9SPeter Brune ierr = PetscSFSetFromOptions(psf);CHKERRQ(ierr); 940298fd71SBarry Smith ierr = PetscSFSetGraphLayout(psf,map,pn,NULL,PETSC_OWN_POINTER,pnatidx);CHKERRQ(ierr); 95ff9846d9SPeter Brune 96ff9846d9SPeter Brune /* broadcast the global indices through to the patch */ 97ff9846d9SPeter Brune ierr = PetscSFBcastBegin(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 98ff9846d9SPeter Brune ierr = PetscSFBcastEnd(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 99ff9846d9SPeter Brune 100ff9846d9SPeter Brune /* create the IS */ 101ff9846d9SPeter Brune ierr = ISCreateGeneral(comm,pn,pleafidx,PETSC_OWN_POINTER,is);CHKERRQ(ierr); 102ff9846d9SPeter Brune 103ff9846d9SPeter Brune ierr = PetscSFDestroy(&psf);CHKERRQ(ierr); 104ff9846d9SPeter Brune 105ff9846d9SPeter Brune ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 106ff9846d9SPeter Brune 107ff9846d9SPeter Brune ierr = PetscFree(globidx);CHKERRQ(ierr); 108ff9846d9SPeter Brune ierr = PetscFree(leafidx);CHKERRQ(ierr); 109ff9846d9SPeter Brune ierr = PetscFree(natidx);CHKERRQ(ierr); 110ff9846d9SPeter Brune ierr = PetscFree(pnatidx);CHKERRQ(ierr); 111ff9846d9SPeter Brune PetscFunctionReturn(0); 112ff9846d9SPeter Brune } 113ff9846d9SPeter Brune 114ff9846d9SPeter Brune #undef __FUNCT__ 115e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private" 1160adebc6cSBarry Smith PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) 1170adebc6cSBarry Smith { 118e30e807fSPeter Brune DM da; 119e30e807fSPeter Brune PetscErrorCode ierr; 120e30e807fSPeter Brune DMDALocalInfo info; 121e30e807fSPeter Brune PetscReal lmin[3],lmax[3]; 122ff9846d9SPeter Brune PetscInt xsize,ysize,zsize; 123e30e807fSPeter Brune PetscInt xo,yo,zo; 1247ddda789SPeter Brune PetscInt xol,yol,zol; 125e30e807fSPeter Brune 126e30e807fSPeter Brune PetscFunctionBegin; 127e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1287ddda789SPeter Brune ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr); 1297ddda789SPeter Brune 130e30e807fSPeter Brune ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr); 131e30e807fSPeter Brune ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr); 132e30e807fSPeter Brune ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr); 133e30e807fSPeter Brune ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr); 134e30e807fSPeter Brune 135e30e807fSPeter Brune ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr); 136e30e807fSPeter Brune ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr); 137e30e807fSPeter Brune 138e30e807fSPeter Brune /* local with overlap */ 139e30e807fSPeter Brune xsize = info.xm; 140e30e807fSPeter Brune ysize = info.ym; 141e30e807fSPeter Brune zsize = info.zm; 142e30e807fSPeter Brune xo = info.xs; 143e30e807fSPeter Brune yo = info.ys; 144e30e807fSPeter Brune zo = info.zs; 145e30e807fSPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) { 1467ddda789SPeter Brune xsize += xol; 1477ddda789SPeter Brune xo -= xol; 148e30e807fSPeter Brune } 149e30e807fSPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) { 1507ddda789SPeter Brune ysize += yol; 1517ddda789SPeter Brune yo -= yol; 152e30e807fSPeter Brune } 153e30e807fSPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) { 1547ddda789SPeter Brune zsize += zol; 1557ddda789SPeter Brune zo -= zol; 156e30e807fSPeter Brune } 157e30e807fSPeter Brune 1587ddda789SPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) xsize += xol; 1597ddda789SPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) ysize += yol; 1607ddda789SPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) zsize += zol; 1618865f1eaSKarl Rupp 162e30e807fSPeter Brune ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr); 163e30e807fSPeter Brune ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr); 164e30e807fSPeter Brune ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 165e30e807fSPeter Brune 166e30e807fSPeter Brune /* set up as a block instead */ 167e30e807fSPeter Brune ierr = DMSetUp(da);CHKERRQ(ierr); 168e30e807fSPeter Brune 169e30e807fSPeter Brune /* todo - nonuniform coordinates */ 170e30e807fSPeter Brune ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 171e30e807fSPeter Brune ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr); 172e30e807fSPeter Brune 17395c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 17495c13181SPeter Brune ierr = DMDASetOffset(da,xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr); 175e30e807fSPeter Brune 176e30e807fSPeter Brune *dddm = da; 177e30e807fSPeter Brune PetscFunctionReturn(0); 178e30e807fSPeter Brune } 179e30e807fSPeter Brune 180e30e807fSPeter Brune #undef __FUNCT__ 181e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 182e30e807fSPeter Brune /* 183e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 184e30e807fSPeter Brune 185285ae305SPeter Brune Right now this assumes one subdomain per processor. 186285ae305SPeter Brune 187e30e807fSPeter Brune */ 1880adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) 1890adebc6cSBarry Smith { 190e30e807fSPeter Brune PetscErrorCode ierr; 191285ae305SPeter Brune DMDALocalInfo info,subinfo; 192e30e807fSPeter Brune DM subdm; 193285ae305SPeter Brune MatStencil upper,lower; 194285ae305SPeter Brune IS idis,isis,odis,osis,gdis; 195285ae305SPeter Brune Vec svec,dvec,slvec; 196e30e807fSPeter Brune 197e30e807fSPeter Brune PetscFunctionBegin; 198*ce94432eSBarry Smith if (nsubdms != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have more than one subdomain per processor (yet)"); 199285ae305SPeter Brune 200e30e807fSPeter Brune /* allocate the arrays of scatters */ 201050c5e14SJed Brown if (iscat) {ierr = PetscMalloc(sizeof(VecScatter*),iscat);CHKERRQ(ierr);} 202050c5e14SJed Brown if (oscat) {ierr = PetscMalloc(sizeof(VecScatter*),oscat);CHKERRQ(ierr);} 203050c5e14SJed Brown if (lscat) {ierr = PetscMalloc(sizeof(VecScatter*),lscat);CHKERRQ(ierr);} 204e30e807fSPeter Brune 205285ae305SPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 206285ae305SPeter Brune subdm = subdms[0]; 207285ae305SPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 208e30e807fSPeter Brune 209285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 210285ae305SPeter Brune /* TODO - make this actually support multiple subdomains -- subdomain needs to provide where it's nonoverlapping portion belongs */ 211285ae305SPeter Brune lower.i = info.xs; 212285ae305SPeter Brune lower.j = info.ys; 213285ae305SPeter Brune lower.k = info.zs; 214285ae305SPeter Brune upper.i = info.xs+info.xm; 215285ae305SPeter Brune upper.j = info.ys+info.ym; 216285ae305SPeter Brune upper.k = info.zs+info.zm; 217285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); 218285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); 219e30e807fSPeter Brune 220285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 221285ae305SPeter Brune lower.i = subinfo.xs; 222285ae305SPeter Brune lower.j = subinfo.ys; 223285ae305SPeter Brune lower.k = subinfo.zs; 224285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 225285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 226285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 227285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); 228285ae305SPeter Brune ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); 229e30e807fSPeter Brune 230285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 231285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 232285ae305SPeter Brune lower.i = subinfo.gxs; 233285ae305SPeter Brune lower.j = subinfo.gys; 234285ae305SPeter Brune lower.k = subinfo.gzs; 235285ae305SPeter Brune upper.i = subinfo.gxs+subinfo.gxm; 236285ae305SPeter Brune upper.j = subinfo.gys+subinfo.gym; 237285ae305SPeter Brune upper.k = subinfo.gzs+subinfo.gzm; 238e30e807fSPeter Brune 239285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); 240e30e807fSPeter Brune 241e30e807fSPeter Brune /* form the scatter */ 242e30e807fSPeter Brune ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 243e30e807fSPeter Brune ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 244e30e807fSPeter Brune ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 245e30e807fSPeter Brune 246285ae305SPeter Brune if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[0]);CHKERRQ(ierr);} 247285ae305SPeter Brune if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[0]);CHKERRQ(ierr);} 2480298fd71SBarry Smith if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[0]);CHKERRQ(ierr);} 249e30e807fSPeter Brune 250e30e807fSPeter Brune ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 251e30e807fSPeter Brune ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 252e30e807fSPeter Brune ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 253e30e807fSPeter Brune 254e30e807fSPeter Brune ierr = ISDestroy(&idis);CHKERRQ(ierr); 255e30e807fSPeter Brune ierr = ISDestroy(&isis);CHKERRQ(ierr); 256e30e807fSPeter Brune 257e30e807fSPeter Brune ierr = ISDestroy(&odis);CHKERRQ(ierr); 258e30e807fSPeter Brune ierr = ISDestroy(&osis);CHKERRQ(ierr); 259e30e807fSPeter Brune 260e30e807fSPeter Brune ierr = ISDestroy(&gdis);CHKERRQ(ierr); 261e30e807fSPeter Brune PetscFunctionReturn(0); 262e30e807fSPeter Brune } 263e30e807fSPeter Brune 264e30e807fSPeter Brune #undef __FUNCT__ 265e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private" 266e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up 267e30e807fSPeter Brune 268e30e807fSPeter Brune ---------- 269e30e807fSPeter Brune ---------- 270e30e807fSPeter Brune --++++++o= 271e30e807fSPeter Brune --++++++o= 272e30e807fSPeter Brune --++++++o= 273e30e807fSPeter Brune --++++++o= 274e30e807fSPeter Brune --ooooooo= 275e30e807fSPeter Brune --======== 276e30e807fSPeter Brune 277e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's: 278e30e807fSPeter Brune 279e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up 280e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap 281e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter() 282e30e807fSPeter Brune 4. -: In the nonshared ghost region 283e30e807fSPeter Brune */ 284e30e807fSPeter Brune 2850adebc6cSBarry Smith PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) 2860adebc6cSBarry Smith { 287e30e807fSPeter Brune PetscErrorCode ierr; 288e30e807fSPeter Brune DMDALocalInfo info,subinfo; 289285ae305SPeter Brune MatStencil lower,upper; 290e30e807fSPeter Brune 291e30e807fSPeter Brune PetscFunctionBegin; 292e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 293e30e807fSPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 294e30e807fSPeter Brune 295285ae305SPeter Brune /* create the inner IS */ 296285ae305SPeter Brune lower.i = info.xs; 297285ae305SPeter Brune lower.j = info.ys; 298285ae305SPeter Brune lower.k = info.zs; 299285ae305SPeter Brune upper.i = info.xs+info.xm; 300285ae305SPeter Brune upper.j = info.ys+info.ym; 301285ae305SPeter Brune upper.k = info.zs+info.zm; 302e30e807fSPeter Brune 303285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,iis);CHKERRQ(ierr); 304285ae305SPeter Brune 305285ae305SPeter Brune /* create the outer IS */ 306285ae305SPeter Brune lower.i = subinfo.xs; 307285ae305SPeter Brune lower.j = subinfo.ys; 308285ae305SPeter Brune lower.k = subinfo.zs; 309285ae305SPeter Brune upper.i = subinfo.xs+subinfo.xm; 310285ae305SPeter Brune upper.j = subinfo.ys+subinfo.ym; 311285ae305SPeter Brune upper.k = subinfo.zs+subinfo.zm; 312285ae305SPeter Brune ierr = DMDACreatePatchIS(dm,&lower,&upper,ois);CHKERRQ(ierr); 313e30e807fSPeter Brune PetscFunctionReturn(0); 314e30e807fSPeter Brune } 315e30e807fSPeter Brune 31665db9045SDmitry Karpeev 31765db9045SDmitry Karpeev #undef __FUNCT__ 31865db9045SDmitry Karpeev #define __FUNCT__ "DMDASetUseDomainDecomposition" 31965db9045SDmitry Karpeev /*@ 32065db9045SDmitry Karpeev DMDASetUseDomainDecomposition - sets option controlling whether DMCreateDomainDecomposition() is used. 32165db9045SDmitry Karpeev 32265db9045SDmitry Karpeev Logically collective. 32365db9045SDmitry Karpeev 32465db9045SDmitry Karpeev Input Parameters: 32565db9045SDmitry Karpeev - dm - the DM object of type DA 3268c3ce3fcSDmitry Karpeev + flag - PETSC_TRUE or PETSC_FALSE according to whether DMCreateDomainDecomposition() returns geometrically-defined or null subdomains. 32765db9045SDmitry Karpeev 32865db9045SDmitry Karpeev .seealso DMDAGetUseDomainDecomposition() 329129686c5SSatish Balay @*/ 33065db9045SDmitry Karpeev PetscErrorCode DMDASetUseDomainDecomposition(DM dm, PetscBool flag) 33165db9045SDmitry Karpeev { 33285cbd38bSJed Brown DM_DA *dd = (DM_DA*)dm->data; 33365db9045SDmitry Karpeev PetscBool isda; 33465db9045SDmitry Karpeev PetscErrorCode ierr; 33585cbd38bSJed Brown 33665db9045SDmitry Karpeev PetscFunctionBegin; 33765db9045SDmitry Karpeev PetscValidHeaderSpecific(dm,DM_CLASSID,1); 33885cbd38bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr); 33965db9045SDmitry Karpeev if (!isda) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"DM must be of type DMDA"); 34065db9045SDmitry Karpeev dd->decompositiondm = flag; 34165db9045SDmitry Karpeev PetscFunctionReturn(0); 34265db9045SDmitry Karpeev } 34365db9045SDmitry Karpeev 34465db9045SDmitry Karpeev #undef __FUNCT__ 34565db9045SDmitry Karpeev #define __FUNCT__ "DMDAGetUseDomainDecomposition" 34665db9045SDmitry Karpeev /*@ 34765db9045SDmitry Karpeev DMDAGetUseDomainDecomposition - returns option controlling whether DMCreateDomainDecomposition() is used. 34865db9045SDmitry Karpeev 34965db9045SDmitry Karpeev Logically collective. 35065db9045SDmitry Karpeev 35165db9045SDmitry Karpeev Input Parameters: 35265db9045SDmitry Karpeev . dm - the DM object of type DA 35365db9045SDmitry Karpeev 35465db9045SDmitry Karpeev Output Parameters: 3558c3ce3fcSDmitry Karpeev . flag - PETSC_TRUE or PETSC_FALSE according to whether DMCreateDomainDecomposition returns geometrically-defined or null subdomains. 35665db9045SDmitry Karpeev 35765db9045SDmitry Karpeev .seealso DMDASetUseDomainDecomposition() 358129686c5SSatish Balay @*/ 35965db9045SDmitry Karpeev PetscErrorCode DMDAGetUseDomainDecomposition(DM dm,PetscBool *flag) 36065db9045SDmitry Karpeev { 36185cbd38bSJed Brown DM_DA *dd = (DM_DA*)dm->data; 36265db9045SDmitry Karpeev PetscBool isda; 36365db9045SDmitry Karpeev PetscErrorCode ierr; 36485cbd38bSJed Brown 36565db9045SDmitry Karpeev PetscFunctionBegin; 36665db9045SDmitry Karpeev PetscValidHeaderSpecific(dm,DM_CLASSID,1); 36785cbd38bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr); 36865db9045SDmitry Karpeev if (!isda) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"DM must be of type DMDA"); 36965db9045SDmitry Karpeev if (flag) *flag = dd->decompositiondm; 37065db9045SDmitry Karpeev PetscFunctionReturn(0); 37165db9045SDmitry Karpeev } 37265db9045SDmitry Karpeev 373e30e807fSPeter Brune #undef __FUNCT__ 374e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA" 3750adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) 3760adebc6cSBarry Smith { 377e30e807fSPeter Brune PetscErrorCode ierr; 378e30e807fSPeter Brune IS iis0,ois0; 379e30e807fSPeter Brune DM subdm0; 38006cef940SPeter Brune DM_DA *dd = (DM_DA*)dm->data; 3810a696f66SPeter Brune 3820adebc6cSBarry Smith PetscFunctionBegin; 3830a696f66SPeter Brune /* fix to enable PCASM default behavior as taking overlap from the matrix */ 3840a696f66SPeter Brune if (!dd->decompositiondm) { 3850a696f66SPeter Brune if (len) *len=0; 3860a696f66SPeter Brune if (names) *names=0; 3870a696f66SPeter Brune if (iis) *iis=0; 3880a696f66SPeter Brune if (ois) *ois=0; 3890a696f66SPeter Brune if (subdm) *subdm=0; 3900a696f66SPeter Brune PetscFunctionReturn(0); 3910a696f66SPeter Brune } 3920a696f66SPeter Brune 393e30e807fSPeter Brune if (len) *len = 1; 394e30e807fSPeter Brune 395e30e807fSPeter Brune if (iis) {ierr = PetscMalloc(sizeof(IS*),iis);CHKERRQ(ierr);} 396e30e807fSPeter Brune if (ois) {ierr = PetscMalloc(sizeof(IS*),ois);CHKERRQ(ierr);} 397e30e807fSPeter Brune if (subdm) {ierr = PetscMalloc(sizeof(DM*),subdm);CHKERRQ(ierr);} 398e30e807fSPeter Brune if (names) {ierr = PetscMalloc(sizeof(char*),names);CHKERRQ(ierr);} 399e30e807fSPeter Brune ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr); 400e30e807fSPeter Brune ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr); 4018865f1eaSKarl Rupp if (iis) (*iis)[0] = iis0; 4028865f1eaSKarl Rupp else { 403e30e807fSPeter Brune ierr = ISDestroy(&iis0);CHKERRQ(ierr); 404e30e807fSPeter Brune } 4058865f1eaSKarl Rupp if (ois) (*ois)[0] = ois0; 4068865f1eaSKarl Rupp else { 407e30e807fSPeter Brune ierr = ISDestroy(&ois0);CHKERRQ(ierr); 408e30e807fSPeter Brune } 409e30e807fSPeter Brune if (subdm) (*subdm)[0] = subdm0; 410e30e807fSPeter Brune if (names) (*names)[0] = 0; 411e30e807fSPeter Brune PetscFunctionReturn(0); 412e30e807fSPeter Brune } 413