1e30e807fSPeter Brune #include <petsc-private/daimpl.h> 2e30e807fSPeter Brune 3e30e807fSPeter Brune #undef __FUNCT__ 4*ff9846d9SPeter Brune #define __FUNCT__ "DMDACreatePatchIS" 5*ff9846d9SPeter Brune 6*ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) 7*ff9846d9SPeter Brune { 8*ff9846d9SPeter Brune PetscErrorCode ierr; 9*ff9846d9SPeter Brune PetscInt i,j,k,idx; 10*ff9846d9SPeter Brune PetscInt ii,jj,kk; 11*ff9846d9SPeter Brune Vec v; 12*ff9846d9SPeter Brune PetscInt n,pn,bs; 13*ff9846d9SPeter Brune PetscMPIInt rank; 14*ff9846d9SPeter Brune PetscSF sf,psf; 15*ff9846d9SPeter Brune PetscLayout map; 16*ff9846d9SPeter Brune MPI_Comm comm; 17*ff9846d9SPeter Brune PetscInt *natidx,*globidx,*leafidx; 18*ff9846d9SPeter Brune PetscInt *pnatidx,*pleafidx; 19*ff9846d9SPeter Brune PetscInt base; 20*ff9846d9SPeter Brune PetscInt ox,oy,oz; 21*ff9846d9SPeter Brune DM_DA *dd; 22*ff9846d9SPeter Brune PetscFunctionBegin; 23*ff9846d9SPeter Brune 24*ff9846d9SPeter Brune comm = ((PetscObject)da)->comm; 25*ff9846d9SPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 26*ff9846d9SPeter Brune dd = (DM_DA *)da->data; 27*ff9846d9SPeter Brune 28*ff9846d9SPeter Brune /* construct the natural mapping */ 29*ff9846d9SPeter Brune ierr = DMGetGlobalVector(da,&v);CHKERRQ(ierr); 30*ff9846d9SPeter Brune ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 31*ff9846d9SPeter Brune ierr = VecGetOwnershipRange(v,&base,PETSC_NULL);CHKERRQ(ierr); 32*ff9846d9SPeter Brune ierr = VecGetBlockSize(v,&bs);CHKERRQ(ierr); 33*ff9846d9SPeter Brune ierr = DMRestoreGlobalVector(da,&v);CHKERRQ(ierr); 34*ff9846d9SPeter Brune 35*ff9846d9SPeter Brune /* construct the layout */ 36*ff9846d9SPeter Brune ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr); 37*ff9846d9SPeter Brune ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 38*ff9846d9SPeter Brune ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr); 39*ff9846d9SPeter Brune ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 40*ff9846d9SPeter Brune 41*ff9846d9SPeter Brune /* construct the list of natural indices on this process when PETSc ordering is considered */ 42*ff9846d9SPeter Brune ierr = DMDAGetOffset(da,&ox,&oy,&oz);CHKERRQ(ierr); 43*ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&natidx);CHKERRQ(ierr); 44*ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&globidx);CHKERRQ(ierr); 45*ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n,&leafidx);CHKERRQ(ierr); 46*ff9846d9SPeter Brune idx = 0; 47*ff9846d9SPeter Brune for (k=dd->zs;k<dd->ze;k++) { 48*ff9846d9SPeter Brune for (j=dd->ys;j<dd->ye;j++) { 49*ff9846d9SPeter Brune for (i=dd->xs;i<dd->xe;i++) { 50*ff9846d9SPeter Brune natidx[idx] = i + dd->w*(j*dd->M + k*dd->M*dd->N); 51*ff9846d9SPeter Brune globidx[idx] = base + idx; 52*ff9846d9SPeter Brune leafidx[idx] = 0; 53*ff9846d9SPeter Brune idx++; 54*ff9846d9SPeter Brune } 55*ff9846d9SPeter Brune } 56*ff9846d9SPeter Brune } 57*ff9846d9SPeter Brune 58*ff9846d9SPeter Brune if (idx != n) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong."); 59*ff9846d9SPeter Brune 60*ff9846d9SPeter Brune /* construct the SF going from the natural indices to the local set of PETSc indices */ 61*ff9846d9SPeter Brune ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 62*ff9846d9SPeter Brune ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 63*ff9846d9SPeter Brune ierr = PetscSFSetGraphLayout(sf,map,n,PETSC_NULL,PETSC_OWN_POINTER,natidx);CHKERRQ(ierr); 64*ff9846d9SPeter Brune 65*ff9846d9SPeter Brune /* broadcast the global indices over to the corresponding natural indices */ 66*ff9846d9SPeter Brune ierr = PetscSFGatherBegin(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 67*ff9846d9SPeter Brune ierr = PetscSFGatherEnd(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 68*ff9846d9SPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 69*ff9846d9SPeter Brune 70*ff9846d9SPeter Brune 71*ff9846d9SPeter Brune pn = dd->w*(upper->k - lower->k)*(upper->j - lower->j)*(upper->i - lower->i); 72*ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*pn,&pnatidx);CHKERRQ(ierr); 73*ff9846d9SPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*pn,&pleafidx);CHKERRQ(ierr); 74*ff9846d9SPeter Brune idx = 0; 75*ff9846d9SPeter Brune for (k=lower->k-oz;k<upper->k-oz;k++) { 76*ff9846d9SPeter Brune for (j=lower->j-oy;j<upper->j-oy;j++) { 77*ff9846d9SPeter Brune for (i=dd->w*(lower->i-ox);i<dd->w*(upper->i-ox);i++) { 78*ff9846d9SPeter Brune ii = i % (dd->w*dd->M); 79*ff9846d9SPeter Brune jj = j % dd->N; 80*ff9846d9SPeter Brune kk = k % dd->P; 81*ff9846d9SPeter Brune if (ii < 0) ii = dd->w*dd->M + ii; 82*ff9846d9SPeter Brune if (jj < 0) jj = dd->N + jj; 83*ff9846d9SPeter Brune if (kk < 0) kk = dd->P + kk; 84*ff9846d9SPeter Brune pnatidx[idx] = ii + dd->w*(jj*dd->M + kk*dd->M*dd->N); 85*ff9846d9SPeter Brune idx++; 86*ff9846d9SPeter Brune } 87*ff9846d9SPeter Brune } 88*ff9846d9SPeter Brune } 89*ff9846d9SPeter Brune 90*ff9846d9SPeter Brune if (idx != pn) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong"); 91*ff9846d9SPeter Brune 92*ff9846d9SPeter Brune ierr = PetscSFCreate(comm,&psf);CHKERRQ(ierr); 93*ff9846d9SPeter Brune ierr = PetscSFSetFromOptions(psf);CHKERRQ(ierr); 94*ff9846d9SPeter Brune ierr = PetscSFSetGraphLayout(psf,map,pn,PETSC_NULL,PETSC_OWN_POINTER,pnatidx);CHKERRQ(ierr); 95*ff9846d9SPeter Brune 96*ff9846d9SPeter Brune /* broadcast the global indices through to the patch */ 97*ff9846d9SPeter Brune ierr = PetscSFBcastBegin(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 98*ff9846d9SPeter Brune ierr = PetscSFBcastEnd(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 99*ff9846d9SPeter Brune 100*ff9846d9SPeter Brune /* create the IS */ 101*ff9846d9SPeter Brune ierr = ISCreateGeneral(comm,pn,pleafidx,PETSC_OWN_POINTER,is);CHKERRQ(ierr); 102*ff9846d9SPeter Brune 103*ff9846d9SPeter Brune ierr = PetscSFDestroy(&psf);CHKERRQ(ierr); 104*ff9846d9SPeter Brune 105*ff9846d9SPeter Brune ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 106*ff9846d9SPeter Brune 107*ff9846d9SPeter Brune ierr = PetscFree(globidx);CHKERRQ(ierr); 108*ff9846d9SPeter Brune ierr = PetscFree(leafidx);CHKERRQ(ierr); 109*ff9846d9SPeter Brune ierr = PetscFree(natidx);CHKERRQ(ierr); 110*ff9846d9SPeter Brune ierr = PetscFree(pnatidx);CHKERRQ(ierr); 111*ff9846d9SPeter Brune 112*ff9846d9SPeter Brune PetscFunctionReturn(0); 113*ff9846d9SPeter Brune } 114*ff9846d9SPeter Brune 115*ff9846d9SPeter Brune #undef __FUNCT__ 116e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private" 117e30e807fSPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) { 118e30e807fSPeter Brune DM da; 119e30e807fSPeter Brune DM_DA *dd; 120e30e807fSPeter Brune PetscErrorCode ierr; 121e30e807fSPeter Brune DMDALocalInfo info; 122e30e807fSPeter Brune PetscReal lmin[3],lmax[3]; 123*ff9846d9SPeter Brune PetscInt xsize,ysize,zsize; 124e30e807fSPeter Brune PetscInt xo,yo,zo; 125e30e807fSPeter Brune 126e30e807fSPeter Brune PetscFunctionBegin; 127e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 128e30e807fSPeter Brune ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr); 129e30e807fSPeter Brune ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr); 130e30e807fSPeter Brune 131e30e807fSPeter Brune ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr); 132e30e807fSPeter Brune ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr); 133e30e807fSPeter Brune 134e30e807fSPeter Brune ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr); 135e30e807fSPeter Brune ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr); 136e30e807fSPeter Brune 137e30e807fSPeter Brune dd = (DM_DA *)dm->data; 138e30e807fSPeter Brune 139e30e807fSPeter Brune /* local with overlap */ 140e30e807fSPeter Brune xsize = info.xm; 141e30e807fSPeter Brune ysize = info.ym; 142e30e807fSPeter Brune zsize = info.zm; 143e30e807fSPeter Brune xo = info.xs; 144e30e807fSPeter Brune yo = info.ys; 145e30e807fSPeter Brune zo = info.zs; 146e30e807fSPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) { 147e30e807fSPeter Brune xsize += dd->overlap; 148e30e807fSPeter Brune xo -= dd->overlap; 149e30e807fSPeter Brune } 150e30e807fSPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) { 151e30e807fSPeter Brune ysize += dd->overlap; 152e30e807fSPeter Brune yo -= dd->overlap; 153e30e807fSPeter Brune } 154e30e807fSPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) { 155e30e807fSPeter Brune zsize += dd->overlap; 156e30e807fSPeter Brune zo -= dd->overlap; 157e30e807fSPeter Brune } 158e30e807fSPeter Brune 159e30e807fSPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) { 160e30e807fSPeter Brune xsize += dd->overlap; 161e30e807fSPeter Brune } 162e30e807fSPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) { 163e30e807fSPeter Brune ysize += dd->overlap; 164e30e807fSPeter Brune } 165e30e807fSPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) { 166e30e807fSPeter Brune zsize += dd->overlap; 167e30e807fSPeter Brune } 168e30e807fSPeter Brune 169e30e807fSPeter Brune ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr); 170e30e807fSPeter Brune ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr); 171e30e807fSPeter Brune ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 172e30e807fSPeter Brune 173e30e807fSPeter Brune /* set up as a block instead */ 174e30e807fSPeter Brune ierr = DMSetUp(da);CHKERRQ(ierr); 175e30e807fSPeter Brune ierr = DMDASetOffset(da,xo,yo,zo);CHKERRQ(ierr); 176e30e807fSPeter Brune 177e30e807fSPeter Brune 178e30e807fSPeter Brune /* todo - nonuniform coordinates */ 179e30e807fSPeter Brune ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 180e30e807fSPeter Brune ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr); 181e30e807fSPeter Brune 182e30e807fSPeter Brune dd = (DM_DA *)da->data; 183e30e807fSPeter Brune dd->Mo = info.mx; 184e30e807fSPeter Brune dd->No = info.my; 185e30e807fSPeter Brune dd->Po = info.mz; 186e30e807fSPeter Brune 187e30e807fSPeter Brune *dddm = da; 188e30e807fSPeter Brune 189e30e807fSPeter Brune PetscFunctionReturn(0); 190e30e807fSPeter Brune } 191e30e807fSPeter Brune 192e30e807fSPeter Brune #undef __FUNCT__ 193e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 194e30e807fSPeter Brune /* 195e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 196e30e807fSPeter Brune 197e30e807fSPeter Brune */ 198e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) { 199e30e807fSPeter Brune PetscErrorCode ierr; 200e30e807fSPeter Brune DMDALocalInfo dinfo,sinfo; 201e30e807fSPeter Brune IS isis,idis,osis,odis,gsis,gdis; 202e30e807fSPeter Brune PetscInt *ididx,*isidx,*odidx,*osidx,*gdidx,*gsidx,*idx_global,n_global,*idx_sub,n_sub; 203e30e807fSPeter Brune PetscInt l,i,j,k,d,n_i,n_o,n_g,sl,dl,di,dj,dk,si,sj,sk; 204e30e807fSPeter Brune Vec dvec,svec,slvec; 205e30e807fSPeter Brune DM subdm; 206e30e807fSPeter Brune 207e30e807fSPeter Brune PetscFunctionBegin; 208e30e807fSPeter Brune 209e30e807fSPeter Brune /* allocate the arrays of scatters */ 210050c5e14SJed Brown if (iscat) {ierr = PetscMalloc(sizeof(VecScatter *),iscat);CHKERRQ(ierr);} 211050c5e14SJed Brown if (oscat) {ierr = PetscMalloc(sizeof(VecScatter *),oscat);CHKERRQ(ierr);} 212050c5e14SJed Brown if (lscat) {ierr = PetscMalloc(sizeof(VecScatter *),lscat);CHKERRQ(ierr);} 213e30e807fSPeter Brune 214e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&dinfo);CHKERRQ(ierr); 215e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(dm,&n_global,&idx_global);CHKERRQ(ierr); 216e30e807fSPeter Brune for (l = 0;l < nsubdms;l++) { 217e30e807fSPeter Brune n_i = 0; 218e30e807fSPeter Brune n_o = 0; 219e30e807fSPeter Brune n_g = 0; 220e30e807fSPeter Brune subdm = subdms[l]; 221e30e807fSPeter Brune ierr = DMDAGetLocalInfo(subdm,&sinfo);CHKERRQ(ierr); 222e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(subdm,&n_sub,&idx_sub);CHKERRQ(ierr); 223e30e807fSPeter Brune /* count the three region sizes */ 224e30e807fSPeter Brune for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) { 225e30e807fSPeter Brune for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) { 226e30e807fSPeter Brune for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) { 227e30e807fSPeter Brune for (d=0;d<sinfo.dof;d++) { 228e30e807fSPeter Brune if (k >= sinfo.zs && j >= sinfo.ys && i >= sinfo.xs && 229e30e807fSPeter Brune k < sinfo.zs+sinfo.zm && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) { 230e30e807fSPeter Brune 231e30e807fSPeter Brune /* interior - subinterior overlap */ 232e30e807fSPeter Brune if (k >= dinfo.zs && j >= dinfo.ys && i >= dinfo.xs && 233e30e807fSPeter Brune k < dinfo.zs+dinfo.zm && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) { 234e30e807fSPeter Brune n_i++; 235e30e807fSPeter Brune } 236e30e807fSPeter Brune /* ghost - subinterior overlap */ 237e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 238e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 239e30e807fSPeter Brune n_o++; 240e30e807fSPeter Brune } 241e30e807fSPeter Brune } 242e30e807fSPeter Brune 243e30e807fSPeter Brune /* ghost - subghost overlap */ 244e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 245e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 246e30e807fSPeter Brune n_g++; 247e30e807fSPeter Brune } 248e30e807fSPeter Brune } 249e30e807fSPeter Brune } 250e30e807fSPeter Brune } 251e30e807fSPeter Brune } 252e30e807fSPeter Brune 253e30e807fSPeter Brune if (n_g == 0) SETERRQ(((PetscObject)subdm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Processor-local domain and subdomain do not intersect!"); 254e30e807fSPeter Brune 255e30e807fSPeter Brune /* local and subdomain local index set indices */ 256e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_i,&ididx);CHKERRQ(ierr); 257e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_i,&isidx);CHKERRQ(ierr); 258e30e807fSPeter Brune 259e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_o,&odidx);CHKERRQ(ierr); 260e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_o,&osidx);CHKERRQ(ierr); 261e30e807fSPeter Brune 262e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gdidx);CHKERRQ(ierr); 263e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gsidx);CHKERRQ(ierr); 264e30e807fSPeter Brune 265e30e807fSPeter Brune n_i = 0; n_o = 0;n_g = 0; 266e30e807fSPeter Brune for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) { 267e30e807fSPeter Brune for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) { 268e30e807fSPeter Brune for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) { 269e30e807fSPeter Brune for (d=0;d<sinfo.dof;d++) { 270e30e807fSPeter Brune si = i - sinfo.gxs; 271e30e807fSPeter Brune sj = j - sinfo.gys; 272e30e807fSPeter Brune sk = k - sinfo.gzs; 273e30e807fSPeter Brune sl = d + sinfo.dof*(si + sinfo.gxm*(sj + sinfo.gym*sk)); 274e30e807fSPeter Brune di = i - dinfo.gxs; 275e30e807fSPeter Brune dj = j - dinfo.gys; 276e30e807fSPeter Brune dk = k - dinfo.gzs; 277e30e807fSPeter Brune dl = d + dinfo.dof*(di + dinfo.gxm*(dj + dinfo.gym*dk)); 278e30e807fSPeter Brune 279e30e807fSPeter Brune if (k >= sinfo.zs && j >= sinfo.ys && i >= sinfo.xs && 280e30e807fSPeter Brune k < sinfo.zs+sinfo.zm && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) { 281e30e807fSPeter Brune 282e30e807fSPeter Brune /* interior - subinterior overlap */ 283e30e807fSPeter Brune if (k >= dinfo.zs && j >= dinfo.ys && i >= dinfo.xs && 284e30e807fSPeter Brune k < dinfo.zs+dinfo.zm && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) { 285e30e807fSPeter Brune ididx[n_i] = idx_global[dl]; 286e30e807fSPeter Brune isidx[n_i] = idx_sub[sl]; 287e30e807fSPeter Brune n_i++; 288e30e807fSPeter Brune } 289e30e807fSPeter Brune /* ghost - subinterior overlap */ 290e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 291e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 292e30e807fSPeter Brune odidx[n_o] = idx_global[dl]; 293e30e807fSPeter Brune osidx[n_o] = idx_sub[sl]; 294e30e807fSPeter Brune n_o++; 295e30e807fSPeter Brune } 296e30e807fSPeter Brune } 297e30e807fSPeter Brune 298e30e807fSPeter Brune /* ghost - subghost overlap */ 299e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 300e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 301e30e807fSPeter Brune gdidx[n_g] = idx_global[dl]; 302e30e807fSPeter Brune gsidx[n_g] = sl; 303e30e807fSPeter Brune n_g++; 304e30e807fSPeter Brune } 305e30e807fSPeter Brune } 306e30e807fSPeter Brune } 307e30e807fSPeter Brune } 308e30e807fSPeter Brune } 309e30e807fSPeter Brune 310e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,ididx,PETSC_OWN_POINTER,&idis);CHKERRQ(ierr); 311e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,isidx,PETSC_OWN_POINTER,&isis);CHKERRQ(ierr); 312e30e807fSPeter Brune 313e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,odidx,PETSC_OWN_POINTER,&odis);CHKERRQ(ierr); 314e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,osidx,PETSC_OWN_POINTER,&osis);CHKERRQ(ierr); 315e30e807fSPeter Brune 316e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gdidx,PETSC_OWN_POINTER,&gdis);CHKERRQ(ierr); 317e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gsidx,PETSC_OWN_POINTER,&gsis);CHKERRQ(ierr); 318e30e807fSPeter Brune 319e30e807fSPeter Brune /* form the scatter */ 320e30e807fSPeter Brune ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 321e30e807fSPeter Brune ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 322e30e807fSPeter Brune ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 323e30e807fSPeter Brune 324e30e807fSPeter Brune if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[l]);CHKERRQ(ierr);} 325e30e807fSPeter Brune if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[l]);CHKERRQ(ierr);} 326e30e807fSPeter Brune if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,gsis,&(*lscat)[l]);CHKERRQ(ierr);} 327e30e807fSPeter Brune 328e30e807fSPeter Brune ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 329e30e807fSPeter Brune ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 330e30e807fSPeter Brune ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 331e30e807fSPeter Brune 332e30e807fSPeter Brune ierr = ISDestroy(&idis);CHKERRQ(ierr); 333e30e807fSPeter Brune ierr = ISDestroy(&isis);CHKERRQ(ierr); 334e30e807fSPeter Brune 335e30e807fSPeter Brune ierr = ISDestroy(&odis);CHKERRQ(ierr); 336e30e807fSPeter Brune ierr = ISDestroy(&osis);CHKERRQ(ierr); 337e30e807fSPeter Brune 338e30e807fSPeter Brune ierr = ISDestroy(&gdis);CHKERRQ(ierr); 339e30e807fSPeter Brune ierr = ISDestroy(&gsis);CHKERRQ(ierr); 340e30e807fSPeter Brune } 341e30e807fSPeter Brune PetscFunctionReturn(0); 342e30e807fSPeter Brune 343e30e807fSPeter Brune } 344e30e807fSPeter Brune 345e30e807fSPeter Brune #undef __FUNCT__ 346e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private" 347e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up 348e30e807fSPeter Brune 349e30e807fSPeter Brune ---------- 350e30e807fSPeter Brune ---------- 351e30e807fSPeter Brune --++++++o= 352e30e807fSPeter Brune --++++++o= 353e30e807fSPeter Brune --++++++o= 354e30e807fSPeter Brune --++++++o= 355e30e807fSPeter Brune --ooooooo= 356e30e807fSPeter Brune --======== 357e30e807fSPeter Brune 358e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's: 359e30e807fSPeter Brune 360e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up 361e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap 362e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter() 363e30e807fSPeter Brune 4. -: In the nonshared ghost region 364e30e807fSPeter Brune */ 365e30e807fSPeter Brune 366e30e807fSPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) { 367e30e807fSPeter Brune PetscErrorCode ierr; 368e30e807fSPeter Brune DMDALocalInfo info,subinfo; 369e30e807fSPeter Brune PetscInt *iiidx,*oiidx,*gidx,gindx; 370e30e807fSPeter Brune PetscInt i,j,k,d,n,nsub,nover,llindx,lindx,li,lj,lk,gi,gj,gk; 371e30e807fSPeter Brune 372e30e807fSPeter Brune PetscFunctionBegin; 373e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 374e30e807fSPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 375e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(dm,&n,&gidx);CHKERRQ(ierr); 376e30e807fSPeter Brune /* todo -- overlap */ 377e30e807fSPeter Brune nsub = info.xm*info.ym*info.zm*info.dof; 378e30e807fSPeter Brune nover = subinfo.xm*subinfo.ym*subinfo.zm*subinfo.dof; 379e30e807fSPeter Brune /* iis is going to have size of the local problem's global part but have a lot of fill-in */ 380e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*nsub,&iiidx);CHKERRQ(ierr); 381e30e807fSPeter Brune /* ois is going to have size of the local problem's global part */ 382e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*nover,&oiidx);CHKERRQ(ierr); 383e30e807fSPeter Brune /* loop over the ghost region of the subdm and fill in the indices */ 384e30e807fSPeter Brune for (k=subinfo.gzs;k<subinfo.gzs+subinfo.gzm;k++) { 385e30e807fSPeter Brune for (j=subinfo.gys;j<subinfo.gys+subinfo.gym;j++) { 386e30e807fSPeter Brune for (i=subinfo.gxs;i<subinfo.gxs+subinfo.gxm;i++) { 387e30e807fSPeter Brune for (d=0;d<subinfo.dof;d++) { 388e30e807fSPeter Brune li = i - subinfo.xs; 389e30e807fSPeter Brune lj = j - subinfo.ys; 390e30e807fSPeter Brune lk = k - subinfo.zs; 391e30e807fSPeter Brune lindx = d + subinfo.dof*(li + subinfo.xm*(lj + subinfo.ym*lk)); 392e30e807fSPeter Brune li = i - info.xs; 393e30e807fSPeter Brune lj = j - info.ys; 394e30e807fSPeter Brune lk = k - info.zs; 395e30e807fSPeter Brune llindx = d + info.dof*(li + info.xm*(lj + info.ym*lk)); 396e30e807fSPeter Brune gi = i - info.gxs; 397e30e807fSPeter Brune gj = j - info.gys; 398e30e807fSPeter Brune gk = k - info.gzs; 399e30e807fSPeter Brune gindx = d + info.dof*(gi + info.gxm*(gj + info.gym*gk)); 400e30e807fSPeter Brune 401e30e807fSPeter Brune /* check if the current point is inside the interior region */ 402e30e807fSPeter Brune if (k >= info.zs && j >= info.ys && i >= info.xs && 403e30e807fSPeter Brune k < info.zs+info.zm && j < info.ys+info.ym && i < info.xs+info.xm) { 404e30e807fSPeter Brune iiidx[llindx] = gidx[gindx]; 405e30e807fSPeter Brune oiidx[lindx] = gidx[gindx]; 406e30e807fSPeter Brune /* overlap region */ 407e30e807fSPeter Brune } else if (k >= subinfo.zs && j >= subinfo.ys && i >= subinfo.xs && 408e30e807fSPeter Brune k < subinfo.zs+subinfo.zm && j < subinfo.ys+subinfo.ym && i < subinfo.xs+subinfo.xm) { 409e30e807fSPeter Brune oiidx[lindx] = gidx[gindx]; 410e30e807fSPeter Brune } 411e30e807fSPeter Brune } 412e30e807fSPeter Brune } 413e30e807fSPeter Brune } 414e30e807fSPeter Brune } 415e30e807fSPeter Brune 416e30e807fSPeter Brune /* create the index sets */ 417e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,nsub,iiidx,PETSC_OWN_POINTER,iis);CHKERRQ(ierr); 418e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,nover,oiidx,PETSC_OWN_POINTER,ois);CHKERRQ(ierr); 419e30e807fSPeter Brune PetscFunctionReturn(0); 420e30e807fSPeter Brune } 421e30e807fSPeter Brune 422e30e807fSPeter Brune #undef __FUNCT__ 423e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA" 424e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) { 425e30e807fSPeter Brune PetscErrorCode ierr; 426e30e807fSPeter Brune IS iis0,ois0; 427e30e807fSPeter Brune DM subdm0; 428e30e807fSPeter Brune PetscFunctionBegin; 429e30e807fSPeter Brune if (len)*len = 1; 430e30e807fSPeter Brune 431e30e807fSPeter Brune if (iis) {ierr = PetscMalloc(sizeof(IS *),iis);CHKERRQ(ierr);} 432e30e807fSPeter Brune if (ois) {ierr = PetscMalloc(sizeof(IS *),ois);CHKERRQ(ierr);} 433e30e807fSPeter Brune if (subdm) {ierr = PetscMalloc(sizeof(DM *),subdm);CHKERRQ(ierr);} 434e30e807fSPeter Brune if (names) {ierr = PetscMalloc(sizeof(char *),names);CHKERRQ(ierr);} 435e30e807fSPeter Brune ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr); 436e30e807fSPeter Brune ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr); 437e30e807fSPeter Brune if (iis) { 438e30e807fSPeter Brune (*iis)[0] = iis0; 439e30e807fSPeter Brune } else { 440e30e807fSPeter Brune ierr = ISDestroy(&iis0);CHKERRQ(ierr); 441e30e807fSPeter Brune } 442e30e807fSPeter Brune if (ois) { 443e30e807fSPeter Brune (*ois)[0] = ois0; 444e30e807fSPeter Brune } else { 445e30e807fSPeter Brune ierr = ISDestroy(&ois0);CHKERRQ(ierr); 446e30e807fSPeter Brune } 447e30e807fSPeter Brune if (subdm) (*subdm)[0] = subdm0; 448e30e807fSPeter Brune if (names) (*names)[0] = 0; 449e30e807fSPeter Brune PetscFunctionReturn(0); 450e30e807fSPeter Brune } 451