/*@ DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields Collective on DMDA Input Parameters: + da - the distributed array - nfields - number of fields in new DMDA Output Parameter: . nda - the new DMDA Level: intermediate .keywords: distributed array, get, corners, nodes, local indices, coordinates .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates() @*/ PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) { PetscErrorCode ierr; DM_DA *dd = (DM_DA*)da->data; PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; const PetscInt *lx,*ly,*lz; DMBoundaryType bx,by,bz; DMDAStencilType stencil_type; PetscInt ox,oy,oz; PetscInt cl,rl; PetscFunctionBegin; dim = da->dim; M = dd->M; N = dd->N; P = dd->P; m = dd->m; n = dd->n; p = dd->p; s = dd->s; bx = dd->bx; by = dd->by; bz = dd->bz; stencil_type = dd->stencil_type; ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); if (dim == 1) { ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr); } else if (dim == 2) { ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr); } else if (dim == 3) { ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr); } ierr = DMSetUp(*nda);CHKERRQ(ierr); if (da->coordinates) { ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr); (*nda)->coordinates = da->coordinates; } /* allow for getting a reduced DA corresponding to a domain decomposition */ ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr); ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr); /* allow for getting a reduced DA corresponding to a coarsened DA */ ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr); ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr); (*nda)->levelup = rl; (*nda)->leveldown = cl; PetscFunctionReturn(0); }
PETSC_EXTERN void PETSC_STDCALL dmdagetoffset_(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po, int *__ierr ){ *__ierr = DMDAGetOffset( (DM)PetscToPointer((da) ),xo,yo,zo,Mo,No,Po); }
/*@ DMDACreatePatchIS - Creates an index set corresponding to a patch of the DA. Not Collective Input Parameters: + da - the DMDA . lower - a matstencil with i, j and k corresponding to the lower corner of the patch - upper - a matstencil with i, j and k corresponding to the upper corner of the patch Output Parameters: . is - the IS corresponding to the patch Level: developer .seealso: DMDACreateDomainDecomposition(), DMDACreateDomainDecompositionScatters() @*/ PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is) { PetscInt ms=0,ns=0,ps=0; PetscInt me=1,ne=1,pe=1; PetscInt mr=0,nr=0,pr=0; PetscInt ii,jj,kk; PetscInt si,sj,sk; PetscInt i,j,k,l,idx; PetscInt base; PetscInt xm=1,ym=1,zm=1; const PetscInt *lx,*ly,*lz; PetscInt ox,oy,oz; PetscInt m,n,p,M,N,P,dof; PetscInt nindices; PetscInt *indices; DM_DA *dd = (DM_DA*)da->data; PetscErrorCode ierr; PetscFunctionBegin; /* need to get the sizes of the actual DM rather than the "global" space of a subdomain DM */ M = dd->M;N = dd->N;P=dd->P; m = dd->m;n = dd->n;p=dd->p; dof = dd->w; ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr); ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); nindices = (upper->i - lower->i)*(upper->j - lower->j)*(upper->k - lower->k)*dof; ierr = PetscMalloc(sizeof(PetscInt)*nindices,&indices);CHKERRQ(ierr); /* start at index 0 on processor 0 */ mr = 0; nr = 0; pr = 0; ms = 0; ns = 0; ps = 0; if (lx) me = lx[0]; if (ly) ne = ly[0]; if (lz) pe = lz[0]; idx = 0; for (k=lower->k-oz;k<upper->k-oz;k++) { for (j=lower->j-oy;j < upper->j-oy;j++) { for (i=lower->i-ox;i < upper->i-ox;i++) { /* "actual" indices rather than ones outside of the domain */ ii = i; jj = j; kk = k; if (ii < 0) ii = M + ii; if (jj < 0) jj = N + jj; if (kk < 0) kk = P + kk; if (ii > M-1) ii = ii - M; if (jj > N-1) jj = jj - N; if (kk > P-1) kk = kk - P; /* gone out of processor range on x axis */ while(ii > me-1 || ii < ms) { if (mr == m-1) { ms = 0; me = lx[0]; mr = 0; } else { mr++; ms = me; me += lx[mr]; } } /* gone out of processor range on y axis */ while(jj > ne-1 || jj < ns) { if (nr == n-1) { ns = 0; ne = ly[0]; nr = 0; } else { nr++; ns = ne; ne += ly[nr]; } } /* gone out of processor range on z axis */ while(kk > pe-1 || kk < ps) { if (pr == p-1) { ps = 0; pe = lz[0]; pr = 0; } else { pr++; ps = pe; pe += lz[pr]; } } /* compute the vector base on owning processor */ xm = me - ms; ym = ne - ns; zm = pe - ps; base = ms*ym*zm + ns*M*zm + ps*M*N; /* compute the local coordinates on owning processor */ si = ii - ms; sj = jj - ns; sk = kk - ps; for (l=0;l<dof;l++) { indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk); idx++; } } } } ISCreateGeneral(PETSC_COMM_SELF,idx,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr); PetscFunctionReturn(0); }