1 #include <petsc-private/daimpl.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 PetscFunctionBegin; 23 24 comm = ((PetscObject)da)->comm; 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,PETSC_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);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,PETSC_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,PETSC_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 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "DMDASubDomainDA_Private" 117 PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) { 118 DM da; 119 DM_DA *dd; 120 PetscErrorCode ierr; 121 DMDALocalInfo info; 122 PetscReal lmin[3],lmax[3]; 123 PetscInt xsize,ysize,zsize; 124 PetscInt xo,yo,zo; 125 126 PetscFunctionBegin; 127 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 128 ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr); 129 ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr); 130 131 ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr); 132 ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr); 133 134 ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr); 135 ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr); 136 137 dd = (DM_DA *)dm->data; 138 139 /* local with overlap */ 140 xsize = info.xm; 141 ysize = info.ym; 142 zsize = info.zm; 143 xo = info.xs; 144 yo = info.ys; 145 zo = info.zs; 146 if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) { 147 xsize += dd->overlap; 148 xo -= dd->overlap; 149 } 150 if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) { 151 ysize += dd->overlap; 152 yo -= dd->overlap; 153 } 154 if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) { 155 zsize += dd->overlap; 156 zo -= dd->overlap; 157 } 158 159 if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) { 160 xsize += dd->overlap; 161 } 162 if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) { 163 ysize += dd->overlap; 164 } 165 if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) { 166 zsize += dd->overlap; 167 } 168 169 ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr); 170 ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr); 171 ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr); 172 173 /* set up as a block instead */ 174 ierr = DMSetUp(da);CHKERRQ(ierr); 175 ierr = DMDASetOffset(da,xo,yo,zo);CHKERRQ(ierr); 176 177 178 /* todo - nonuniform coordinates */ 179 ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 180 ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr); 181 182 dd = (DM_DA *)da->data; 183 dd->Mo = info.mx; 184 dd->No = info.my; 185 dd->Po = info.mz; 186 187 *dddm = da; 188 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA" 194 /* 195 Fills the local vector problem on the subdomain from the global problem. 196 197 Right now this assumes one subdomain per processor. 198 199 */ 200 PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) { 201 PetscErrorCode ierr; 202 DMDALocalInfo info,subinfo; 203 DM subdm; 204 MatStencil upper,lower; 205 IS idis,isis,odis,osis,gdis; 206 Vec svec,dvec,slvec; 207 208 PetscFunctionBegin; 209 210 if (nsubdms != 1) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot have more than one subdomain per processor (yet)"); 211 212 /* allocate the arrays of scatters */ 213 if (iscat) {ierr = PetscMalloc(sizeof(VecScatter *),iscat);CHKERRQ(ierr);} 214 if (oscat) {ierr = PetscMalloc(sizeof(VecScatter *),oscat);CHKERRQ(ierr);} 215 if (lscat) {ierr = PetscMalloc(sizeof(VecScatter *),lscat);CHKERRQ(ierr);} 216 217 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 218 subdm = subdms[0]; 219 ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 220 221 /* create the global and subdomain index sets for the inner domain */ 222 /* TODO - make this actually support multiple subdomains -- subdomain needs to provide where it's nonoverlapping portion belongs */ 223 lower.i = info.xs; 224 lower.j = info.ys; 225 lower.k = info.zs; 226 upper.i = info.xs+info.xm; 227 upper.j = info.ys+info.ym; 228 upper.k = info.zs+info.zm; 229 ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr); 230 ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr); 231 232 /* create the global and subdomain index sets for the outer subdomain */ 233 lower.i = subinfo.xs; 234 lower.j = subinfo.ys; 235 lower.k = subinfo.zs; 236 upper.i = subinfo.xs+subinfo.xm; 237 upper.j = subinfo.ys+subinfo.ym; 238 upper.k = subinfo.zs+subinfo.zm; 239 ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr); 240 ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr); 241 242 /* global and subdomain ISes for the local indices of the subdomain */ 243 /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 244 lower.i = subinfo.gxs; 245 lower.j = subinfo.gys; 246 lower.k = subinfo.gzs; 247 upper.i = subinfo.gxs+subinfo.gxm; 248 upper.j = subinfo.gys+subinfo.gym; 249 upper.k = subinfo.gzs+subinfo.gzm; 250 251 ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr); 252 253 /* form the scatter */ 254 ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr); 255 ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr); 256 ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr); 257 258 if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[0]);CHKERRQ(ierr);} 259 if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[0]);CHKERRQ(ierr);} 260 if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,PETSC_NULL,&(*lscat)[0]);CHKERRQ(ierr);} 261 262 ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr); 263 ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr); 264 ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr); 265 266 ierr = ISDestroy(&idis);CHKERRQ(ierr); 267 ierr = ISDestroy(&isis);CHKERRQ(ierr); 268 269 ierr = ISDestroy(&odis);CHKERRQ(ierr); 270 ierr = ISDestroy(&osis);CHKERRQ(ierr); 271 272 ierr = ISDestroy(&gdis);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "DMDASubDomainIS_Private" 278 /* We have that the interior regions are going to be the same, but the ghost regions might not match up 279 280 ---------- 281 ---------- 282 --++++++o= 283 --++++++o= 284 --++++++o= 285 --++++++o= 286 --ooooooo= 287 --======== 288 289 Therefore, for each point in the overall, we must check if it's: 290 291 1. +: In the interior of the global dm; it lines up 292 2. o: In the overlap region -- for now the same as 1; no overlap 293 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter() 294 4. -: In the nonshared ghost region 295 */ 296 297 PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) { 298 PetscErrorCode ierr; 299 DMDALocalInfo info,subinfo; 300 MatStencil lower,upper; 301 302 PetscFunctionBegin; 303 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 304 ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr); 305 306 /* create the inner IS */ 307 lower.i = info.xs; 308 lower.j = info.ys; 309 lower.k = info.zs; 310 upper.i = info.xs+info.xm; 311 upper.j = info.ys+info.ym; 312 upper.k = info.zs+info.zm; 313 314 ierr = DMDACreatePatchIS(dm,&lower,&upper,iis);CHKERRQ(ierr); 315 316 /* create the outer IS */ 317 lower.i = subinfo.xs; 318 lower.j = subinfo.ys; 319 lower.k = subinfo.zs; 320 upper.i = subinfo.xs+subinfo.xm; 321 upper.j = subinfo.ys+subinfo.ym; 322 upper.k = subinfo.zs+subinfo.zm; 323 ierr = DMDACreatePatchIS(dm,&lower,&upper,ois);CHKERRQ(ierr); 324 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "DMCreateDomainDecomposition_DA" 330 PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) { 331 PetscErrorCode ierr; 332 IS iis0,ois0; 333 DM subdm0; 334 PetscFunctionBegin; 335 if (len)*len = 1; 336 337 if (iis) {ierr = PetscMalloc(sizeof(IS *),iis);CHKERRQ(ierr);} 338 if (ois) {ierr = PetscMalloc(sizeof(IS *),ois);CHKERRQ(ierr);} 339 if (subdm) {ierr = PetscMalloc(sizeof(DM *),subdm);CHKERRQ(ierr);} 340 if (names) {ierr = PetscMalloc(sizeof(char *),names);CHKERRQ(ierr);} 341 ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr); 342 ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr); 343 if (iis) { 344 (*iis)[0] = iis0; 345 } else { 346 ierr = ISDestroy(&iis0);CHKERRQ(ierr); 347 } 348 if (ois) { 349 (*ois)[0] = ois0; 350 } else { 351 ierr = ISDestroy(&ois0);CHKERRQ(ierr); 352 } 353 if (subdm) (*subdm)[0] = subdm0; 354 if (names) (*names)[0] = 0; 355 PetscFunctionReturn(0); 356 } 357