PETSC_EXTERN void PETSC_STDCALL dmdagetoverlap_(DM da,PetscInt *x,PetscInt *y,PetscInt *z, int *__ierr ){ *__ierr = DMDAGetOverlap( (DM)PetscToPointer((da) ),x,y,z); }
PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) { DM *da; PetscInt dim,size,i,j,k,idx; PetscErrorCode ierr; DMDALocalInfo info; PetscInt xsize,ysize,zsize; PetscInt xo,yo,zo; PetscInt xs,ys,zs; PetscInt xm=1,ym=1,zm=1; PetscInt xol,yol,zol; PetscInt m=1,n=1,p=1; PetscInt M,N,P; PetscInt pm,mtmp; PetscFunctionBegin; ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr); ierr = DMDAGetNumLocalSubDomains(dm,&size);CHKERRQ(ierr); ierr = PetscMalloc1(size,&da);CHKERRQ(ierr); if (nlocal) *nlocal = size; dim = info.dim; M = info.xm; N = info.ym; P = info.zm; if (dim == 1) { m = size; } else if (dim == 2) { m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N))); while (m > 0) { n = size/m; if (m*n*p == size) break; m--; } } else if (dim == 3) { n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); if (!n) n = 1; while (n > 0) { pm = size/n; if (n*pm == size) break; n--; } if (!n) n = 1; m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n))); if (!m) m = 1; while (m > 0) { p = size/(m*n); if (m*n*p == size) break; m--; } if (M > P && m < p) {mtmp = m; m = p; p = mtmp;} } zs = info.zs; idx = 0; for (k = 0; k < p; k++) { ys = info.ys; for (j = 0; j < n; j++) { xs = info.xs; for (i = 0; i < m; i++) { if (dim == 1) { xm = M/m + ((M % m) > i); } else if (dim == 2) { xm = M/m + ((M % m) > i); ym = N/n + ((N % n) > j); } else if (dim == 3) { xm = M/m + ((M % m) > i); ym = N/n + ((N % n) > j); zm = P/p + ((P % p) > k); } xsize = xm; ysize = ym; zsize = zm; xo = xs; yo = ys; zo = zs; ierr = DMDACreate(PETSC_COMM_SELF,&(da[idx]));CHKERRQ(ierr); ierr = DMSetOptionsPrefix(da[idx],"sub_");CHKERRQ(ierr); ierr = DMSetDimension(da[idx], info.dim);CHKERRQ(ierr); ierr = DMDASetDof(da[idx], info.dof);CHKERRQ(ierr); ierr = DMDASetStencilType(da[idx],info.st);CHKERRQ(ierr); ierr = DMDASetStencilWidth(da[idx],info.sw);CHKERRQ(ierr); if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) { xsize += xol; xo -= xol; } if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) { ysize += yol; yo -= yol; } if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) { zsize += zol; zo -= zol; } if (info.bx == DM_BOUNDARY_PERIODIC || (xs+xm != info.mx)) xsize += xol; if (info.by == DM_BOUNDARY_PERIODIC || (ys+ym != info.my)) ysize += yol; if (info.bz == DM_BOUNDARY_PERIODIC || (zs+zm != info.mz)) zsize += zol; if (info.bx != DM_BOUNDARY_PERIODIC) { if (xo < 0) { xsize += xo; xo = 0; } if (xo+xsize > info.mx-1) { xsize -= xo+xsize - info.mx; } } if (info.by != DM_BOUNDARY_PERIODIC) { if (yo < 0) { ysize += yo; yo = 0; } if (yo+ysize > info.my-1) { ysize -= yo+ysize - info.my; } } if (info.bz != DM_BOUNDARY_PERIODIC) { if (zo < 0) { zsize += zo; zo = 0; } if (zo+zsize > info.mz-1) { zsize -= zo+zsize - info.mz; } } ierr = DMDASetSizes(da[idx], xsize, ysize, zsize);CHKERRQ(ierr); ierr = DMDASetNumProcs(da[idx], 1, 1, 1);CHKERRQ(ierr); ierr = DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED);CHKERRQ(ierr); /* set up as a block instead */ ierr = DMSetUp(da[idx]);CHKERRQ(ierr); /* nonoverlapping region */ ierr = DMDASetNonOverlappingRegion(da[idx],xs,ys,zs,xm,ym,zm);CHKERRQ(ierr); /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ ierr = DMDASetOffset(da[idx],xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr); xs += xm; idx++; } ys += ym; } zs += zm; } *sdm = da; PetscFunctionReturn(0); }