1 #include <petsc-private/dmdaimpl.h> 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMDACreatePatchIS" 5 6 PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) 7 { 8 PetscErrorCode ierr; 9 PetscInt i,j,k,idx; 10 PetscInt ii,jj,kk; 11 Vec v; 12 PetscInt n,pn,bs; 13 PetscMPIInt rank; 14 PetscSF sf,psf; 15 PetscLayout map; 16 MPI_Comm comm; 17 PetscInt *natidx,*globidx,*leafidx; 18 PetscInt *pnatidx,*pleafidx; 19 PetscInt base; 20 PetscInt ox,oy,oz; 21 DM_DA *dd; 22 23 PetscFunctionBegin; 24 ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 25 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 26 dd = (DM_DA*)da->data; 27 28 /* construct the natural mapping */ 29 ierr = DMGetGlobalVector(da,&v);CHKERRQ(ierr); 30 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 31 ierr = VecGetOwnershipRange(v,&base,NULL);CHKERRQ(ierr); 32 ierr = VecGetBlockSize(v,&bs);CHKERRQ(ierr); 33 ierr = DMRestoreGlobalVector(da,&v);CHKERRQ(ierr); 34 35 /* construct the layout */ 36 ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr); 37 ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 38 ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr); 39 ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 40 41 /* construct the list of natural indices on this process when PETSc ordering is considered */ 42 ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr); 43 ierr = PetscMalloc(sizeof(PetscInt)*n,&natidx);CHKERRQ(ierr); 44 ierr = PetscMalloc(sizeof(PetscInt)*n,&globidx);CHKERRQ(ierr); 45 ierr = PetscMalloc(sizeof(PetscInt)*n,&leafidx);CHKERRQ(ierr); 46 idx = 0; 47 for (k=dd->zs; k<dd->ze; k++) { 48 for (j=dd->ys; j<dd->ye; j++) { 49 for (i=dd->xs; i<dd->xe; i++) { 50 natidx[idx] = i + dd->w*(j*dd->M + k*dd->M*dd->N); 51 globidx[idx] = base + idx; 52 leafidx[idx] = 0; 53 idx++; 54 } 55 } 56 } 57 58 if (idx != n) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong."); 59 60 /* construct the SF going from the natural indices to the local set of PETSc indices */ 61 ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 62 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 63 ierr = PetscSFSetGraphLayout(sf,map,n,NULL,PETSC_OWN_POINTER,natidx);CHKERRQ(ierr); 64 65 /* broadcast the global indices over to the corresponding natural indices */ 66 ierr = PetscSFGatherBegin(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 67 ierr = PetscSFGatherEnd(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr); 68 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 69 70 71 pn = dd->w*(upper->k - lower->k)*(upper->j - lower->j)*(upper->i - lower->i); 72 ierr = PetscMalloc(sizeof(PetscInt)*pn,&pnatidx);CHKERRQ(ierr); 73 ierr = PetscMalloc(sizeof(PetscInt)*pn,&pleafidx);CHKERRQ(ierr); 74 idx = 0; 75 for (k=lower->k-oz; k<upper->k-oz; k++) { 76 for (j=lower->j-oy; j<upper->j-oy; j++) { 77 for (i=dd->w*(lower->i-ox); i<dd->w*(upper->i-ox); i++) { 78 ii = i % (dd->w*dd->M); 79 jj = j % dd->N; 80 kk = k % dd->P; 81 if (ii < 0) ii = dd->w*dd->M + ii; 82 if (jj < 0) jj = dd->N + jj; 83 if (kk < 0) kk = dd->P + kk; 84 pnatidx[idx] = ii + dd->w*(jj*dd->M + kk*dd->M*dd->N); 85 idx++; 86 } 87 } 88 } 89 90 if (idx != pn) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong"); 91 92 ierr = PetscSFCreate(comm,&psf);CHKERRQ(ierr); 93 ierr = PetscSFSetFromOptions(psf);CHKERRQ(ierr); 94 ierr = PetscSFSetGraphLayout(psf,map,pn,NULL,PETSC_OWN_POINTER,pnatidx);CHKERRQ(ierr); 95 96 /* broadcast the global indices through to the patch */ 97 ierr = PetscSFBcastBegin(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 98 ierr = PetscSFBcastEnd(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr); 99 100 /* create the IS */ 101 ierr = ISCreateGeneral(comm,pn,pleafidx,PETSC_OWN_POINTER,is);CHKERRQ(ierr); 102 103 ierr = PetscSFDestroy(&psf);CHKERRQ(ierr); 104 105 ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 106 107 ierr = PetscFree(globidx);CHKERRQ(ierr); 108 ierr = PetscFree(leafidx);CHKERRQ(ierr); 109 ierr = PetscFree(natidx);CHKERRQ(ierr); 110 ierr = PetscFree(pnatidx);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 #undef __FUNCT__ 115 #define __FUNCT__ "DMDASubDomainDA_Private" 116 PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) 117 { 118 DM da; 119 PetscErrorCode ierr; 120 DMDALocalInfo info; 121 PetscReal lmin[3],lmax[3]; 122 PetscInt xsize,ysize,zsize; 123 PetscInt xo,yo,zo; 124 PetscInt xol,yol,zol; 125 126 PetscFunctionBegin; 127 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 128 ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr); 129 130 ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr); 131 ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr); 132 ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr); 133 ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr); 134 135 ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr); 136 ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr); 137 138 /* local with overlap */ 139 xsize = info.xm; 140 ysize = info.ym; 141 zsize = info.zm; 142 xo = info.xs; 143 yo = info.ys; 144 zo = info.zs; 145 if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) { 146 xsize += xol; 147 xo -= xol; 148 } 149 if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) { 150 ysize += yol; 151 yo -= yol; 152 } 153 if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) { 154 zsize += zol; 155 zo -= zol; 156 } 157 158 if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) xsize += xol; 159 if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) ysize += yol; 160 if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) zsize += zol; 161 162 ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr); 163 ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr); 164 ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 165 166 /* set up as a block instead */ 167 ierr = DMSetUp(da);CHKERRQ(ierr); 168 169 /* todo - nonuniform coordinates */ 170 ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 171 ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr); 172 173 /* nonoverlapping region */ 174 ierr = DMDASetNonOverlappingRegion(da,info.xs,info.ys,info.zs,info.xm,info.ym,info.zm);CHKERRQ(ierr); 175 176 /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 177 ierr = DMDASetOffset(da,xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr); 178 *dddm = da; 179 PetscFunctionReturn(0); 180 } 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 184 /* 185 Fills the local vector problem on the subdomain from the global problem. 186 187 Right now this assumes one subdomain per processor. 188 189 */ 190 PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) 191 { 192 PetscErrorCode ierr; 193 DMDALocalInfo info,subinfo; 194 DM subdm; 195 MatStencil upper,lower; 196 IS idis,isis,odis,osis,gdis; 197 Vec svec,dvec,slvec; 198 PetscInt xm,ym,zm,xs,ys,zs; 199 200 PetscFunctionBegin; 201 if (nsubdms != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have more than one subdomain per processor (yet)"); 202 203 /* allocate the arrays of scatters */ 204 if (iscat) {ierr = PetscMalloc(sizeof(VecScatter*),iscat);CHKERRQ(ierr);} 205 if (oscat) {ierr = PetscMalloc(sizeof(VecScatter*),oscat);CHKERRQ(ierr);} 206 if (lscat) {ierr = PetscMalloc(sizeof(VecScatter*),lscat);CHKERRQ(ierr);} 207 208 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 209 subdm = subdms[0]; 210 ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 211 212 /* create the global and subdomain index sets for the inner domain */ 213 /* TODO - make this actually support multiple subdomains -- subdomain needs to provide where it's nonoverlapping portion belongs */ 214 ierr = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 215 lower.i = xs; 216 lower.j = ys; 217 lower.k = zs; 218 upper.i = xs+xm; 219 upper.j = ys+ym; 220 upper.k = zs+zm; 221 ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); 222 ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); 223 224 /* create the global and subdomain index sets for the outer subdomain */ 225 lower.i = subinfo.xs; 226 lower.j = subinfo.ys; 227 lower.k = subinfo.zs; 228 upper.i = subinfo.xs+subinfo.xm; 229 upper.j = subinfo.ys+subinfo.ym; 230 upper.k = subinfo.zs+subinfo.zm; 231 ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); 232 ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); 233 234 /* global and subdomain ISes for the local indices of the subdomain */ 235 /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 236 lower.i = subinfo.gxs; 237 lower.j = subinfo.gys; 238 lower.k = subinfo.gzs; 239 upper.i = subinfo.gxs+subinfo.gxm; 240 upper.j = subinfo.gys+subinfo.gym; 241 upper.k = subinfo.gzs+subinfo.gzm; 242 243 ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); 244 245 /* form the scatter */ 246 ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 247 ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 248 ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 249 250 if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[0]);CHKERRQ(ierr);} 251 if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[0]);CHKERRQ(ierr);} 252 if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[0]);CHKERRQ(ierr);} 253 254 ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 255 ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 256 ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 257 258 ierr = ISDestroy(&idis);CHKERRQ(ierr); 259 ierr = ISDestroy(&isis);CHKERRQ(ierr); 260 261 ierr = ISDestroy(&odis);CHKERRQ(ierr); 262 ierr = ISDestroy(&osis);CHKERRQ(ierr); 263 264 ierr = ISDestroy(&gdis);CHKERRQ(ierr); 265 PetscFunctionReturn(0); 266 } 267 268 #undef __FUNCT__ 269 #define __FUNCT__ "DMDASubDomainIS_Private" 270 /* We have that the interior regions are going to be the same, but the ghost regions might not match up 271 272 ---------- 273 ---------- 274 --++++++o= 275 --++++++o= 276 --++++++o= 277 --++++++o= 278 --ooooooo= 279 --======== 280 281 Therefore, for each point in the overall, we must check if it's: 282 283 1. +: In the interior of the global dm; it lines up 284 2. o: In the overlap region -- for now the same as 1; no overlap 285 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter() 286 4. -: In the nonshared ghost region 287 */ 288 289 PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) 290 { 291 PetscErrorCode ierr; 292 DMDALocalInfo info,subinfo; 293 MatStencil lower,upper; 294 295 PetscFunctionBegin; 296 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 297 ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 298 299 /* create the inner IS */ 300 lower.i = info.xs; 301 lower.j = info.ys; 302 lower.k = info.zs; 303 upper.i = info.xs+info.xm; 304 upper.j = info.ys+info.ym; 305 upper.k = info.zs+info.zm; 306 307 ierr = DMDACreatePatchIS(dm,&lower,&upper,iis);CHKERRQ(ierr); 308 309 /* create the outer IS */ 310 lower.i = subinfo.xs; 311 lower.j = subinfo.ys; 312 lower.k = subinfo.zs; 313 upper.i = subinfo.xs+subinfo.xm; 314 upper.j = subinfo.ys+subinfo.ym; 315 upper.k = subinfo.zs+subinfo.zm; 316 ierr = DMDACreatePatchIS(dm,&lower,&upper,ois);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 #undef __FUNCT__ 321 #define __FUNCT__ "DMCreateDomainDecomposition_DA" 322 PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) 323 { 324 PetscErrorCode ierr; 325 IS iis0,ois0; 326 DM subdm0; 327 328 PetscFunctionBegin; 329 if (len) *len = 1; 330 331 if (iis) {ierr = PetscMalloc(sizeof(IS*),iis);CHKERRQ(ierr);} 332 if (ois) {ierr = PetscMalloc(sizeof(IS*),ois);CHKERRQ(ierr);} 333 if (subdm) {ierr = PetscMalloc(sizeof(DM*),subdm);CHKERRQ(ierr);} 334 if (names) {ierr = PetscMalloc(sizeof(char*),names);CHKERRQ(ierr);} 335 ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr); 336 ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr); 337 if (iis) (*iis)[0] = iis0; 338 else { 339 ierr = ISDestroy(&iis0);CHKERRQ(ierr); 340 } 341 if (ois) (*ois)[0] = ois0; 342 else { 343 ierr = ISDestroy(&ois0);CHKERRQ(ierr); 344 } 345 if (subdm) (*subdm)[0] = subdm0; 346 if (names) (*names)[0] = 0; 347 PetscFunctionReturn(0); 348 } 349