1e30e807fSPeter Brune #include <petsc-private/daimpl.h> 2e30e807fSPeter Brune 3e30e807fSPeter Brune #undef __FUNCT__ 4e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private" 5e30e807fSPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) { 6e30e807fSPeter Brune DM da; 7e30e807fSPeter Brune DM_DA *dd; 8e30e807fSPeter Brune PetscErrorCode ierr; 9e30e807fSPeter Brune DMDALocalInfo info; 10e30e807fSPeter Brune PetscReal lmin[3],lmax[3]; 11e30e807fSPeter Brune PetscInt xsize,ysize,zsize; 12e30e807fSPeter Brune PetscInt xo,yo,zo; 13e30e807fSPeter Brune 14e30e807fSPeter Brune PetscFunctionBegin; 15e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 16e30e807fSPeter Brune ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr); 17e30e807fSPeter Brune ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr); 18e30e807fSPeter Brune 19e30e807fSPeter Brune ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr); 20e30e807fSPeter Brune ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr); 21e30e807fSPeter Brune 22e30e807fSPeter Brune ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr); 23e30e807fSPeter Brune ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr); 24e30e807fSPeter Brune 25e30e807fSPeter Brune dd = (DM_DA *)dm->data; 26e30e807fSPeter Brune 27e30e807fSPeter Brune /* local with overlap */ 28e30e807fSPeter Brune xsize = info.xm; 29e30e807fSPeter Brune ysize = info.ym; 30e30e807fSPeter Brune zsize = info.zm; 31e30e807fSPeter Brune xo = info.xs; 32e30e807fSPeter Brune yo = info.ys; 33e30e807fSPeter Brune zo = info.zs; 34e30e807fSPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) { 35e30e807fSPeter Brune xsize += dd->overlap; 36e30e807fSPeter Brune xo -= dd->overlap; 37e30e807fSPeter Brune } 38e30e807fSPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) { 39e30e807fSPeter Brune ysize += dd->overlap; 40e30e807fSPeter Brune yo -= dd->overlap; 41e30e807fSPeter Brune } 42e30e807fSPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) { 43e30e807fSPeter Brune zsize += dd->overlap; 44e30e807fSPeter Brune zo -= dd->overlap; 45e30e807fSPeter Brune } 46e30e807fSPeter Brune 47e30e807fSPeter Brune if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) { 48e30e807fSPeter Brune xsize += dd->overlap; 49e30e807fSPeter Brune } 50e30e807fSPeter Brune if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) { 51e30e807fSPeter Brune ysize += dd->overlap; 52e30e807fSPeter Brune } 53e30e807fSPeter Brune if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) { 54e30e807fSPeter Brune zsize += dd->overlap; 55e30e807fSPeter Brune } 56e30e807fSPeter Brune 57e30e807fSPeter Brune ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr); 58e30e807fSPeter Brune ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr); 59e30e807fSPeter Brune ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 60e30e807fSPeter Brune 61e30e807fSPeter Brune /* set up as a block instead */ 62e30e807fSPeter Brune ierr = DMSetUp(da);CHKERRQ(ierr); 63e30e807fSPeter Brune ierr = DMDASetOffset(da,xo,yo,zo);CHKERRQ(ierr); 64e30e807fSPeter Brune 65e30e807fSPeter Brune 66e30e807fSPeter Brune /* todo - nonuniform coordinates */ 67e30e807fSPeter Brune ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 68e30e807fSPeter Brune ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr); 69e30e807fSPeter Brune 70e30e807fSPeter Brune dd = (DM_DA *)da->data; 71e30e807fSPeter Brune dd->Mo = info.mx; 72e30e807fSPeter Brune dd->No = info.my; 73e30e807fSPeter Brune dd->Po = info.mz; 74e30e807fSPeter Brune 75e30e807fSPeter Brune *dddm = da; 76e30e807fSPeter Brune 77e30e807fSPeter Brune PetscFunctionReturn(0); 78e30e807fSPeter Brune } 79e30e807fSPeter Brune 80e30e807fSPeter Brune #undef __FUNCT__ 81e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 82e30e807fSPeter Brune /* 83e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 84e30e807fSPeter Brune 85e30e807fSPeter Brune */ 86e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) { 87e30e807fSPeter Brune PetscErrorCode ierr; 88e30e807fSPeter Brune DMDALocalInfo dinfo,sinfo; 89e30e807fSPeter Brune IS isis,idis,osis,odis,gsis,gdis; 90e30e807fSPeter Brune PetscInt *ididx,*isidx,*odidx,*osidx,*gdidx,*gsidx,*idx_global,n_global,*idx_sub,n_sub; 91e30e807fSPeter Brune PetscInt l,i,j,k,d,n_i,n_o,n_g,sl,dl,di,dj,dk,si,sj,sk; 92e30e807fSPeter Brune Vec dvec,svec,slvec; 93e30e807fSPeter Brune DM subdm; 94e30e807fSPeter Brune 95e30e807fSPeter Brune PetscFunctionBegin; 96e30e807fSPeter Brune 97e30e807fSPeter Brune /* allocate the arrays of scatters */ 98*050c5e14SJed Brown if (iscat) {ierr = PetscMalloc(sizeof(VecScatter *),iscat);CHKERRQ(ierr);} 99*050c5e14SJed Brown if (oscat) {ierr = PetscMalloc(sizeof(VecScatter *),oscat);CHKERRQ(ierr);} 100*050c5e14SJed Brown if (lscat) {ierr = PetscMalloc(sizeof(VecScatter *),lscat);CHKERRQ(ierr);} 101e30e807fSPeter Brune 102e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&dinfo);CHKERRQ(ierr); 103e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(dm,&n_global,&idx_global);CHKERRQ(ierr); 104e30e807fSPeter Brune for (l = 0;l < nsubdms;l++) { 105e30e807fSPeter Brune n_i = 0; 106e30e807fSPeter Brune n_o = 0; 107e30e807fSPeter Brune n_g = 0; 108e30e807fSPeter Brune subdm = subdms[l]; 109e30e807fSPeter Brune ierr = DMDAGetLocalInfo(subdm,&sinfo);CHKERRQ(ierr); 110e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(subdm,&n_sub,&idx_sub);CHKERRQ(ierr); 111e30e807fSPeter Brune /* count the three region sizes */ 112e30e807fSPeter Brune for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) { 113e30e807fSPeter Brune for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) { 114e30e807fSPeter Brune for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) { 115e30e807fSPeter Brune for (d=0;d<sinfo.dof;d++) { 116e30e807fSPeter Brune if (k >= sinfo.zs && j >= sinfo.ys && i >= sinfo.xs && 117e30e807fSPeter Brune k < sinfo.zs+sinfo.zm && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) { 118e30e807fSPeter Brune 119e30e807fSPeter Brune /* interior - subinterior overlap */ 120e30e807fSPeter Brune if (k >= dinfo.zs && j >= dinfo.ys && i >= dinfo.xs && 121e30e807fSPeter Brune k < dinfo.zs+dinfo.zm && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) { 122e30e807fSPeter Brune n_i++; 123e30e807fSPeter Brune } 124e30e807fSPeter Brune /* ghost - subinterior overlap */ 125e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 126e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 127e30e807fSPeter Brune n_o++; 128e30e807fSPeter Brune } 129e30e807fSPeter Brune } 130e30e807fSPeter Brune 131e30e807fSPeter Brune /* ghost - subghost overlap */ 132e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 133e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 134e30e807fSPeter Brune n_g++; 135e30e807fSPeter Brune } 136e30e807fSPeter Brune } 137e30e807fSPeter Brune } 138e30e807fSPeter Brune } 139e30e807fSPeter Brune } 140e30e807fSPeter Brune 141e30e807fSPeter Brune if (n_g == 0) SETERRQ(((PetscObject)subdm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Processor-local domain and subdomain do not intersect!"); 142e30e807fSPeter Brune 143e30e807fSPeter Brune /* local and subdomain local index set indices */ 144e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_i,&ididx);CHKERRQ(ierr); 145e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_i,&isidx);CHKERRQ(ierr); 146e30e807fSPeter Brune 147e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_o,&odidx);CHKERRQ(ierr); 148e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_o,&osidx);CHKERRQ(ierr); 149e30e807fSPeter Brune 150e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gdidx);CHKERRQ(ierr); 151e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gsidx);CHKERRQ(ierr); 152e30e807fSPeter Brune 153e30e807fSPeter Brune n_i = 0; n_o = 0;n_g = 0; 154e30e807fSPeter Brune for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) { 155e30e807fSPeter Brune for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) { 156e30e807fSPeter Brune for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) { 157e30e807fSPeter Brune for (d=0;d<sinfo.dof;d++) { 158e30e807fSPeter Brune si = i - sinfo.gxs; 159e30e807fSPeter Brune sj = j - sinfo.gys; 160e30e807fSPeter Brune sk = k - sinfo.gzs; 161e30e807fSPeter Brune sl = d + sinfo.dof*(si + sinfo.gxm*(sj + sinfo.gym*sk)); 162e30e807fSPeter Brune di = i - dinfo.gxs; 163e30e807fSPeter Brune dj = j - dinfo.gys; 164e30e807fSPeter Brune dk = k - dinfo.gzs; 165e30e807fSPeter Brune dl = d + dinfo.dof*(di + dinfo.gxm*(dj + dinfo.gym*dk)); 166e30e807fSPeter Brune 167e30e807fSPeter Brune if (k >= sinfo.zs && j >= sinfo.ys && i >= sinfo.xs && 168e30e807fSPeter Brune k < sinfo.zs+sinfo.zm && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) { 169e30e807fSPeter Brune 170e30e807fSPeter Brune /* interior - subinterior overlap */ 171e30e807fSPeter Brune if (k >= dinfo.zs && j >= dinfo.ys && i >= dinfo.xs && 172e30e807fSPeter Brune k < dinfo.zs+dinfo.zm && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) { 173e30e807fSPeter Brune ididx[n_i] = idx_global[dl]; 174e30e807fSPeter Brune isidx[n_i] = idx_sub[sl]; 175e30e807fSPeter Brune n_i++; 176e30e807fSPeter Brune } 177e30e807fSPeter Brune /* ghost - subinterior overlap */ 178e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 179e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 180e30e807fSPeter Brune odidx[n_o] = idx_global[dl]; 181e30e807fSPeter Brune osidx[n_o] = idx_sub[sl]; 182e30e807fSPeter Brune n_o++; 183e30e807fSPeter Brune } 184e30e807fSPeter Brune } 185e30e807fSPeter Brune 186e30e807fSPeter Brune /* ghost - subghost overlap */ 187e30e807fSPeter Brune if (k >= dinfo.gzs && j >= dinfo.gys && i >= dinfo.gxs && 188e30e807fSPeter Brune k < dinfo.gzs+dinfo.gzm && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) { 189e30e807fSPeter Brune gdidx[n_g] = idx_global[dl]; 190e30e807fSPeter Brune gsidx[n_g] = sl; 191e30e807fSPeter Brune n_g++; 192e30e807fSPeter Brune } 193e30e807fSPeter Brune } 194e30e807fSPeter Brune } 195e30e807fSPeter Brune } 196e30e807fSPeter Brune } 197e30e807fSPeter Brune 198e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,ididx,PETSC_OWN_POINTER,&idis);CHKERRQ(ierr); 199e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,isidx,PETSC_OWN_POINTER,&isis);CHKERRQ(ierr); 200e30e807fSPeter Brune 201e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,odidx,PETSC_OWN_POINTER,&odis);CHKERRQ(ierr); 202e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,osidx,PETSC_OWN_POINTER,&osis);CHKERRQ(ierr); 203e30e807fSPeter Brune 204e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gdidx,PETSC_OWN_POINTER,&gdis);CHKERRQ(ierr); 205e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gsidx,PETSC_OWN_POINTER,&gsis);CHKERRQ(ierr); 206e30e807fSPeter Brune 207e30e807fSPeter Brune /* form the scatter */ 208e30e807fSPeter Brune ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 209e30e807fSPeter Brune ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 210e30e807fSPeter Brune ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 211e30e807fSPeter Brune 212e30e807fSPeter Brune if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[l]);CHKERRQ(ierr);} 213e30e807fSPeter Brune if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[l]);CHKERRQ(ierr);} 214e30e807fSPeter Brune if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,gsis,&(*lscat)[l]);CHKERRQ(ierr);} 215e30e807fSPeter Brune 216e30e807fSPeter Brune ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 217e30e807fSPeter Brune ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 218e30e807fSPeter Brune ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 219e30e807fSPeter Brune 220e30e807fSPeter Brune ierr = ISDestroy(&idis);CHKERRQ(ierr); 221e30e807fSPeter Brune ierr = ISDestroy(&isis);CHKERRQ(ierr); 222e30e807fSPeter Brune 223e30e807fSPeter Brune ierr = ISDestroy(&odis);CHKERRQ(ierr); 224e30e807fSPeter Brune ierr = ISDestroy(&osis);CHKERRQ(ierr); 225e30e807fSPeter Brune 226e30e807fSPeter Brune ierr = ISDestroy(&gdis);CHKERRQ(ierr); 227e30e807fSPeter Brune ierr = ISDestroy(&gsis);CHKERRQ(ierr); 228e30e807fSPeter Brune } 229e30e807fSPeter Brune PetscFunctionReturn(0); 230e30e807fSPeter Brune 231e30e807fSPeter Brune } 232e30e807fSPeter Brune 233e30e807fSPeter Brune #undef __FUNCT__ 234e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private" 235e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up 236e30e807fSPeter Brune 237e30e807fSPeter Brune ---------- 238e30e807fSPeter Brune ---------- 239e30e807fSPeter Brune --++++++o= 240e30e807fSPeter Brune --++++++o= 241e30e807fSPeter Brune --++++++o= 242e30e807fSPeter Brune --++++++o= 243e30e807fSPeter Brune --ooooooo= 244e30e807fSPeter Brune --======== 245e30e807fSPeter Brune 246e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's: 247e30e807fSPeter Brune 248e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up 249e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap 250e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter() 251e30e807fSPeter Brune 4. -: In the nonshared ghost region 252e30e807fSPeter Brune */ 253e30e807fSPeter Brune 254e30e807fSPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) { 255e30e807fSPeter Brune PetscErrorCode ierr; 256e30e807fSPeter Brune DMDALocalInfo info,subinfo; 257e30e807fSPeter Brune PetscInt *iiidx,*oiidx,*gidx,gindx; 258e30e807fSPeter Brune PetscInt i,j,k,d,n,nsub,nover,llindx,lindx,li,lj,lk,gi,gj,gk; 259e30e807fSPeter Brune 260e30e807fSPeter Brune PetscFunctionBegin; 261e30e807fSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 262e30e807fSPeter Brune ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 263e30e807fSPeter Brune ierr = DMDAGetGlobalIndices(dm,&n,&gidx);CHKERRQ(ierr); 264e30e807fSPeter Brune /* todo -- overlap */ 265e30e807fSPeter Brune nsub = info.xm*info.ym*info.zm*info.dof; 266e30e807fSPeter Brune nover = subinfo.xm*subinfo.ym*subinfo.zm*subinfo.dof; 267e30e807fSPeter Brune /* iis is going to have size of the local problem's global part but have a lot of fill-in */ 268e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*nsub,&iiidx);CHKERRQ(ierr); 269e30e807fSPeter Brune /* ois is going to have size of the local problem's global part */ 270e30e807fSPeter Brune ierr = PetscMalloc(sizeof(PetscInt)*nover,&oiidx);CHKERRQ(ierr); 271e30e807fSPeter Brune /* loop over the ghost region of the subdm and fill in the indices */ 272e30e807fSPeter Brune for (k=subinfo.gzs;k<subinfo.gzs+subinfo.gzm;k++) { 273e30e807fSPeter Brune for (j=subinfo.gys;j<subinfo.gys+subinfo.gym;j++) { 274e30e807fSPeter Brune for (i=subinfo.gxs;i<subinfo.gxs+subinfo.gxm;i++) { 275e30e807fSPeter Brune for (d=0;d<subinfo.dof;d++) { 276e30e807fSPeter Brune li = i - subinfo.xs; 277e30e807fSPeter Brune lj = j - subinfo.ys; 278e30e807fSPeter Brune lk = k - subinfo.zs; 279e30e807fSPeter Brune lindx = d + subinfo.dof*(li + subinfo.xm*(lj + subinfo.ym*lk)); 280e30e807fSPeter Brune li = i - info.xs; 281e30e807fSPeter Brune lj = j - info.ys; 282e30e807fSPeter Brune lk = k - info.zs; 283e30e807fSPeter Brune llindx = d + info.dof*(li + info.xm*(lj + info.ym*lk)); 284e30e807fSPeter Brune gi = i - info.gxs; 285e30e807fSPeter Brune gj = j - info.gys; 286e30e807fSPeter Brune gk = k - info.gzs; 287e30e807fSPeter Brune gindx = d + info.dof*(gi + info.gxm*(gj + info.gym*gk)); 288e30e807fSPeter Brune 289e30e807fSPeter Brune /* check if the current point is inside the interior region */ 290e30e807fSPeter Brune if (k >= info.zs && j >= info.ys && i >= info.xs && 291e30e807fSPeter Brune k < info.zs+info.zm && j < info.ys+info.ym && i < info.xs+info.xm) { 292e30e807fSPeter Brune iiidx[llindx] = gidx[gindx]; 293e30e807fSPeter Brune oiidx[lindx] = gidx[gindx]; 294e30e807fSPeter Brune /* overlap region */ 295e30e807fSPeter Brune } else if (k >= subinfo.zs && j >= subinfo.ys && i >= subinfo.xs && 296e30e807fSPeter Brune k < subinfo.zs+subinfo.zm && j < subinfo.ys+subinfo.ym && i < subinfo.xs+subinfo.xm) { 297e30e807fSPeter Brune oiidx[lindx] = gidx[gindx]; 298e30e807fSPeter Brune } 299e30e807fSPeter Brune } 300e30e807fSPeter Brune } 301e30e807fSPeter Brune } 302e30e807fSPeter Brune } 303e30e807fSPeter Brune 304e30e807fSPeter Brune /* create the index sets */ 305e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,nsub,iiidx,PETSC_OWN_POINTER,iis);CHKERRQ(ierr); 306e30e807fSPeter Brune ierr = ISCreateGeneral(PETSC_COMM_SELF,nover,oiidx,PETSC_OWN_POINTER,ois);CHKERRQ(ierr); 307e30e807fSPeter Brune PetscFunctionReturn(0); 308e30e807fSPeter Brune } 309e30e807fSPeter Brune 310e30e807fSPeter Brune #undef __FUNCT__ 311e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA" 312e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) { 313e30e807fSPeter Brune PetscErrorCode ierr; 314e30e807fSPeter Brune IS iis0,ois0; 315e30e807fSPeter Brune DM subdm0; 316e30e807fSPeter Brune PetscFunctionBegin; 317e30e807fSPeter Brune if (len)*len = 1; 318e30e807fSPeter Brune 319e30e807fSPeter Brune if (iis) {ierr = PetscMalloc(sizeof(IS *),iis);CHKERRQ(ierr);} 320e30e807fSPeter Brune if (ois) {ierr = PetscMalloc(sizeof(IS *),ois);CHKERRQ(ierr);} 321e30e807fSPeter Brune if (subdm) {ierr = PetscMalloc(sizeof(DM *),subdm);CHKERRQ(ierr);} 322e30e807fSPeter Brune if (names) {ierr = PetscMalloc(sizeof(char *),names);CHKERRQ(ierr);} 323e30e807fSPeter Brune ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr); 324e30e807fSPeter Brune ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr); 325e30e807fSPeter Brune if (iis) { 326e30e807fSPeter Brune (*iis)[0] = iis0; 327e30e807fSPeter Brune } else { 328e30e807fSPeter Brune ierr = ISDestroy(&iis0);CHKERRQ(ierr); 329e30e807fSPeter Brune } 330e30e807fSPeter Brune if (ois) { 331e30e807fSPeter Brune (*ois)[0] = ois0; 332e30e807fSPeter Brune } else { 333e30e807fSPeter Brune ierr = ISDestroy(&ois0);CHKERRQ(ierr); 334e30e807fSPeter Brune } 335e30e807fSPeter Brune if (subdm) (*subdm)[0] = subdm0; 336e30e807fSPeter Brune if (names) (*names)[0] = 0; 337e30e807fSPeter Brune PetscFunctionReturn(0); 338e30e807fSPeter Brune } 339