1c6db04a5SJed Brown #include <private/daimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 347c6ae99SBarry Smith 447c6ae99SBarry Smith #undef __FUNCT__ 5aa219208SBarry Smith #define __FUNCT__ "DMDASetDim" 647c6ae99SBarry Smith /*@ 7aa219208SBarry Smith DMDASetDim - Sets the dimension 847c6ae99SBarry Smith 9aa219208SBarry Smith Collective on DMDA 1047c6ae99SBarry Smith 1147c6ae99SBarry Smith Input Parameters: 12aa219208SBarry Smith + da - the DMDA 1347c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE) 1447c6ae99SBarry Smith 1547c6ae99SBarry Smith Level: intermediate 1647c6ae99SBarry Smith 17aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes() 1847c6ae99SBarry Smith @*/ 197087cfbeSBarry Smith PetscErrorCode DMDASetDim(DM da, PetscInt dim) 2047c6ae99SBarry Smith { 2147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2247c6ae99SBarry Smith 2347c6ae99SBarry Smith PetscFunctionBegin; 2447c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 25aa219208SBarry Smith if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim); 2647c6ae99SBarry Smith dd->dim = dim; 2747c6ae99SBarry Smith PetscFunctionReturn(0); 2847c6ae99SBarry Smith } 2947c6ae99SBarry Smith 3047c6ae99SBarry Smith #undef __FUNCT__ 31aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes" 3247c6ae99SBarry Smith /*@ 33aa219208SBarry Smith DMDASetSizes - Sets the global sizes 3447c6ae99SBarry Smith 35aa219208SBarry Smith Logically Collective on DMDA 3647c6ae99SBarry Smith 3747c6ae99SBarry Smith Input Parameters: 38aa219208SBarry Smith + da - the DMDA 3947c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE) 4047c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE) 4147c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE) 4247c6ae99SBarry Smith 4347c6ae99SBarry Smith Level: intermediate 4447c6ae99SBarry Smith 45aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership() 4647c6ae99SBarry Smith @*/ 477087cfbeSBarry Smith PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 4847c6ae99SBarry Smith { 4947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith PetscFunctionBegin; 5247c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 5347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,M,2); 5447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,N,3); 5547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,P,4); 5647c6ae99SBarry Smith 5747c6ae99SBarry Smith dd->M = M; 5847c6ae99SBarry Smith dd->N = N; 5947c6ae99SBarry Smith dd->P = P; 6047c6ae99SBarry Smith PetscFunctionReturn(0); 6147c6ae99SBarry Smith } 6247c6ae99SBarry Smith 6347c6ae99SBarry Smith #undef __FUNCT__ 64aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs" 6547c6ae99SBarry Smith /*@ 66aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 6747c6ae99SBarry Smith 68aa219208SBarry Smith Logically Collective on DMDA 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith Input Parameters: 71aa219208SBarry Smith + da - the DMDA 7247c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE) 7347c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE) 7447c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE) 7547c6ae99SBarry Smith 7647c6ae99SBarry Smith Level: intermediate 7747c6ae99SBarry Smith 78aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership() 7947c6ae99SBarry Smith @*/ 807087cfbeSBarry Smith PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 8147c6ae99SBarry Smith { 8247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 8347c6ae99SBarry Smith 8447c6ae99SBarry Smith PetscFunctionBegin; 8547c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 8647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 8747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 8847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 8947c6ae99SBarry Smith dd->m = m; 9047c6ae99SBarry Smith dd->n = n; 9147c6ae99SBarry Smith dd->p = p; 9247c6ae99SBarry Smith PetscFunctionReturn(0); 9347c6ae99SBarry Smith } 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith #undef __FUNCT__ 963e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType" 9747c6ae99SBarry Smith /*@ 981321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith Not collective 10147c6ae99SBarry Smith 10247c6ae99SBarry Smith Input Parameter: 103aa219208SBarry Smith + da - The DMDA 1041321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC 10547c6ae99SBarry Smith 10647c6ae99SBarry Smith Level: intermediate 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith .keywords: distributed array, periodicity 10996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType 11047c6ae99SBarry Smith @*/ 1111321219cSEthan Coon PetscErrorCode DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz) 11247c6ae99SBarry Smith { 11347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith PetscFunctionBegin; 11647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1171321219cSEthan Coon PetscValidLogicalCollectiveInt(da,bx,2); 1181321219cSEthan Coon PetscValidLogicalCollectiveInt(da,by,3); 1191321219cSEthan Coon PetscValidLogicalCollectiveInt(da,bz,4); 1201321219cSEthan Coon dd->bx = bx; 1211321219cSEthan Coon dd->by = by; 1221321219cSEthan Coon dd->bz = bz; 12347c6ae99SBarry Smith PetscFunctionReturn(0); 12447c6ae99SBarry Smith } 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith #undef __FUNCT__ 127aa219208SBarry Smith #define __FUNCT__ "DMDASetDof" 12847c6ae99SBarry Smith /*@ 129aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith Not collective 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith Input Parameter: 134aa219208SBarry Smith + da - The DMDA 13547c6ae99SBarry Smith - dof - Number of degrees of freedom 13647c6ae99SBarry Smith 13747c6ae99SBarry Smith Level: intermediate 13847c6ae99SBarry Smith 13947c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 14096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 14147c6ae99SBarry Smith @*/ 1427087cfbeSBarry Smith PetscErrorCode DMDASetDof(DM da, int dof) 14347c6ae99SBarry Smith { 14447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 14547c6ae99SBarry Smith 14647c6ae99SBarry Smith PetscFunctionBegin; 14747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 14847c6ae99SBarry Smith dd->w = dof; 1491411c6eeSJed Brown da->bs = dof; 15047c6ae99SBarry Smith PetscFunctionReturn(0); 15147c6ae99SBarry Smith } 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith #undef __FUNCT__ 154aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType" 15547c6ae99SBarry Smith /*@ 156aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 15747c6ae99SBarry Smith 158aa219208SBarry Smith Logically Collective on DMDA 15947c6ae99SBarry Smith 16047c6ae99SBarry Smith Input Parameter: 161aa219208SBarry Smith + da - The DMDA 162aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 16347c6ae99SBarry Smith 16447c6ae99SBarry Smith Level: intermediate 16547c6ae99SBarry Smith 16647c6ae99SBarry Smith .keywords: distributed array, stencil 16796e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 16847c6ae99SBarry Smith @*/ 1697087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 17047c6ae99SBarry Smith { 17147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 17247c6ae99SBarry Smith 17347c6ae99SBarry Smith PetscFunctionBegin; 17447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 17547c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 17647c6ae99SBarry Smith dd->stencil_type = stype; 17747c6ae99SBarry Smith PetscFunctionReturn(0); 17847c6ae99SBarry Smith } 17947c6ae99SBarry Smith 18047c6ae99SBarry Smith #undef __FUNCT__ 181aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth" 18247c6ae99SBarry Smith /*@ 183aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 18447c6ae99SBarry Smith 185aa219208SBarry Smith Logically Collective on DMDA 18647c6ae99SBarry Smith 18747c6ae99SBarry Smith Input Parameter: 188aa219208SBarry Smith + da - The DMDA 18947c6ae99SBarry Smith - width - The stencil width 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith Level: intermediate 19247c6ae99SBarry Smith 19347c6ae99SBarry Smith .keywords: distributed array, stencil 19496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 19547c6ae99SBarry Smith @*/ 1967087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 19747c6ae99SBarry Smith { 19847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19947c6ae99SBarry Smith 20047c6ae99SBarry Smith PetscFunctionBegin; 20147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 20247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 20347c6ae99SBarry Smith dd->s = width; 20447c6ae99SBarry Smith PetscFunctionReturn(0); 20547c6ae99SBarry Smith } 20647c6ae99SBarry Smith 20747c6ae99SBarry Smith #undef __FUNCT__ 208aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 209aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 21047c6ae99SBarry Smith { 21147c6ae99SBarry Smith PetscInt i,sum; 21247c6ae99SBarry Smith 21347c6ae99SBarry Smith PetscFunctionBegin; 21447c6ae99SBarry Smith if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 21547c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 21647c6ae99SBarry Smith if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 21747c6ae99SBarry Smith PetscFunctionReturn(0); 21847c6ae99SBarry Smith } 21947c6ae99SBarry Smith 22047c6ae99SBarry Smith #undef __FUNCT__ 221aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges" 22247c6ae99SBarry Smith /*@ 223aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 22447c6ae99SBarry Smith 225aa219208SBarry Smith Logically Collective on DMDA 22647c6ae99SBarry Smith 22747c6ae99SBarry Smith Input Parameter: 228aa219208SBarry Smith + da - The DMDA 22947c6ae99SBarry Smith . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m 23047c6ae99SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n 23147c6ae99SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p. 23247c6ae99SBarry Smith 23347c6ae99SBarry Smith Level: intermediate 23447c6ae99SBarry Smith 23547c6ae99SBarry Smith .keywords: distributed array 23696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 23747c6ae99SBarry Smith @*/ 2387087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 23947c6ae99SBarry Smith { 24047c6ae99SBarry Smith PetscErrorCode ierr; 24147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 24247c6ae99SBarry Smith 24347c6ae99SBarry Smith PetscFunctionBegin; 24447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 24547c6ae99SBarry Smith if (lx) { 24647c6ae99SBarry Smith if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 247aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 24847c6ae99SBarry Smith if (!dd->lx) { 24947c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 25047c6ae99SBarry Smith } 25147c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 25247c6ae99SBarry Smith } 25347c6ae99SBarry Smith if (ly) { 25447c6ae99SBarry Smith if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 255aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 25647c6ae99SBarry Smith if (!dd->ly) { 25747c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 25847c6ae99SBarry Smith } 25947c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 26047c6ae99SBarry Smith } 26147c6ae99SBarry Smith if (lz) { 26247c6ae99SBarry Smith if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 263aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 26447c6ae99SBarry Smith if (!dd->lz) { 26547c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 26647c6ae99SBarry Smith } 26747c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 26847c6ae99SBarry Smith } 26947c6ae99SBarry Smith PetscFunctionReturn(0); 27047c6ae99SBarry Smith } 27147c6ae99SBarry Smith 27247c6ae99SBarry Smith #undef __FUNCT__ 273aa219208SBarry Smith #define __FUNCT__ "DMDACreateOwnershipRanges" 27447c6ae99SBarry Smith /* 275aa219208SBarry Smith Ensure that da->lx, ly, and lz exist. Collective on DMDA. 27647c6ae99SBarry Smith */ 2777087cfbeSBarry Smith PetscErrorCode DMDACreateOwnershipRanges(DM da) 27847c6ae99SBarry Smith { 27947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 28047c6ae99SBarry Smith PetscErrorCode ierr; 28147c6ae99SBarry Smith PetscInt n; 28247c6ae99SBarry Smith MPI_Comm comm; 28347c6ae99SBarry Smith PetscMPIInt size; 28447c6ae99SBarry Smith 28547c6ae99SBarry Smith PetscFunctionBegin; 28647c6ae99SBarry Smith if (!dd->lx) { 28747c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr); 288aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr); 28947c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 29047c6ae99SBarry Smith n = (dd->xe-dd->xs)/dd->w; 29147c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr); 29247c6ae99SBarry Smith } 29347c6ae99SBarry Smith if (dd->dim > 1 && !dd->ly) { 29447c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr); 295aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr); 29647c6ae99SBarry Smith n = dd->ye-dd->ys; 29747c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr); 29847c6ae99SBarry Smith } 29947c6ae99SBarry Smith if (dd->dim > 2 && !dd->lz) { 30047c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr); 301aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr); 30247c6ae99SBarry Smith n = dd->ze-dd->zs; 30347c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);CHKERRQ(ierr); 30447c6ae99SBarry Smith } 30547c6ae99SBarry Smith PetscFunctionReturn(0); 30647c6ae99SBarry Smith } 30747c6ae99SBarry Smith 30847c6ae99SBarry Smith #undef __FUNCT__ 309aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType" 31047c6ae99SBarry Smith /*@ 311aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 31294013140SBarry Smith returned by DMGetInterpolation() 31347c6ae99SBarry Smith 314aa219208SBarry Smith Logically Collective on DMDA 31547c6ae99SBarry Smith 31647c6ae99SBarry Smith Input Parameter: 31747c6ae99SBarry Smith + da - initial distributed array 318aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 31947c6ae99SBarry Smith 32047c6ae99SBarry Smith Level: intermediate 32147c6ae99SBarry Smith 322aa219208SBarry Smith Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation() 32347c6ae99SBarry Smith 32447c6ae99SBarry Smith .keywords: distributed array, interpolation 32547c6ae99SBarry Smith 32696e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 32747c6ae99SBarry Smith @*/ 3287087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 32947c6ae99SBarry Smith { 33047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 33147c6ae99SBarry Smith 33247c6ae99SBarry Smith PetscFunctionBegin; 33347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 33447c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 33547c6ae99SBarry Smith dd->interptype = ctype; 33647c6ae99SBarry Smith PetscFunctionReturn(0); 33747c6ae99SBarry Smith } 33847c6ae99SBarry Smith 339*2dde6fd4SLisandro Dalcin #undef __FUNCT__ 340*2dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType" 341*2dde6fd4SLisandro Dalcin /*@ 342*2dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 343*2dde6fd4SLisandro Dalcin used by DMGetInterpolation() 344*2dde6fd4SLisandro Dalcin 345*2dde6fd4SLisandro Dalcin Not Collective 346*2dde6fd4SLisandro Dalcin 347*2dde6fd4SLisandro Dalcin Input Parameter: 348*2dde6fd4SLisandro Dalcin . da - distributed array 349*2dde6fd4SLisandro Dalcin 350*2dde6fd4SLisandro Dalcin Output Parameter: 351*2dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 352*2dde6fd4SLisandro Dalcin 353*2dde6fd4SLisandro Dalcin Level: intermediate 354*2dde6fd4SLisandro Dalcin 355*2dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 356*2dde6fd4SLisandro Dalcin 357*2dde6fd4SLisandro Dalcin .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMGetInterpolation() 358*2dde6fd4SLisandro Dalcin @*/ 359*2dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 360*2dde6fd4SLisandro Dalcin { 361*2dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 362*2dde6fd4SLisandro Dalcin 363*2dde6fd4SLisandro Dalcin PetscFunctionBegin; 364*2dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 365*2dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 366*2dde6fd4SLisandro Dalcin *ctype = dd->interptype; 367*2dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 368*2dde6fd4SLisandro Dalcin } 36947c6ae99SBarry Smith 37047c6ae99SBarry Smith #undef __FUNCT__ 371aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors" 37247c6ae99SBarry Smith /*@C 373aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 37447c6ae99SBarry Smith processes neighbors. 37547c6ae99SBarry Smith 37647c6ae99SBarry Smith Not Collective 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith Input Parameter: 379aa219208SBarry Smith . da - the DMDA object 38047c6ae99SBarry Smith 38147c6ae99SBarry Smith Output Parameters: 38247c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 38347c6ae99SBarry Smith this process itself is in the list 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 38647c6ae99SBarry Smith Not supported in 1d 387aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 38847c6ae99SBarry Smith 38947c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 39047c6ae99SBarry Smith 39147c6ae99SBarry Smith Level: intermediate 39247c6ae99SBarry Smith 39347c6ae99SBarry Smith @*/ 3947087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 39547c6ae99SBarry Smith { 39647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 39747c6ae99SBarry Smith PetscFunctionBegin; 39847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 39947c6ae99SBarry Smith *ranks = dd->neighbors; 40047c6ae99SBarry Smith PetscFunctionReturn(0); 40147c6ae99SBarry Smith } 40247c6ae99SBarry Smith 40347c6ae99SBarry Smith #undef __FUNCT__ 404aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges" 40547c6ae99SBarry Smith /*@C 406aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 40747c6ae99SBarry Smith 40847c6ae99SBarry Smith Not Collective 40947c6ae99SBarry Smith 41047c6ae99SBarry Smith Input Parameter: 411aa219208SBarry Smith . da - the DMDA object 41247c6ae99SBarry Smith 41347c6ae99SBarry Smith Output Parameter: 41447c6ae99SBarry Smith + lx - ownership along x direction (optional) 41547c6ae99SBarry Smith . ly - ownership along y direction (optional) 41647c6ae99SBarry Smith - lz - ownership along z direction (optional) 41747c6ae99SBarry Smith 41847c6ae99SBarry Smith Level: intermediate 41947c6ae99SBarry Smith 420aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 42147c6ae99SBarry Smith 42247c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 423aa219208SBarry Smith eighth arguments from DMDAGetInfo() 42447c6ae99SBarry Smith 42547c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 426aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 42747c6ae99SBarry Smith 428aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 42947c6ae99SBarry Smith @*/ 4307087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 43147c6ae99SBarry Smith { 43247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 43347c6ae99SBarry Smith 43447c6ae99SBarry Smith PetscFunctionBegin; 43547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 43647c6ae99SBarry Smith if (lx) *lx = dd->lx; 43747c6ae99SBarry Smith if (ly) *ly = dd->ly; 43847c6ae99SBarry Smith if (lz) *lz = dd->lz; 43947c6ae99SBarry Smith PetscFunctionReturn(0); 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith #undef __FUNCT__ 443aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor" 44447c6ae99SBarry Smith /*@ 445aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 44647c6ae99SBarry Smith 447aa219208SBarry Smith Logically Collective on DMDA 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith Input Parameters: 450aa219208SBarry Smith + da - the DMDA object 45147c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 45247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 45347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 45447c6ae99SBarry Smith 45547c6ae99SBarry Smith Options Database: 45647c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 45747c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 45847c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 45947c6ae99SBarry Smith 46047c6ae99SBarry Smith Level: intermediate 46147c6ae99SBarry Smith 46247c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 46347c6ae99SBarry Smith 464aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 46547c6ae99SBarry Smith @*/ 4667087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 46747c6ae99SBarry Smith { 46847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 46947c6ae99SBarry Smith 47047c6ae99SBarry Smith PetscFunctionBegin; 47147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 47247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 47347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 47447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 47547c6ae99SBarry Smith 47647c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 47747c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 47847c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 47947c6ae99SBarry Smith PetscFunctionReturn(0); 48047c6ae99SBarry Smith } 48147c6ae99SBarry Smith 48247c6ae99SBarry Smith #undef __FUNCT__ 483aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor" 48447c6ae99SBarry Smith /*@C 485aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 48647c6ae99SBarry Smith 48747c6ae99SBarry Smith Not Collective 48847c6ae99SBarry Smith 48947c6ae99SBarry Smith Input Parameter: 490aa219208SBarry Smith . da - the DMDA object 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith Output Parameters: 49347c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 49447c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 49547c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 49647c6ae99SBarry Smith 49747c6ae99SBarry Smith Level: intermediate 49847c6ae99SBarry Smith 49947c6ae99SBarry Smith Notes: Pass PETSC_NULL for values you do not need 50047c6ae99SBarry Smith 501aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 50247c6ae99SBarry Smith @*/ 5037087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 50447c6ae99SBarry Smith { 50547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 50647c6ae99SBarry Smith 50747c6ae99SBarry Smith PetscFunctionBegin; 50847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 50947c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 51047c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 51147c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 51247c6ae99SBarry Smith PetscFunctionReturn(0); 51347c6ae99SBarry Smith } 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith #undef __FUNCT__ 516aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix" 51747c6ae99SBarry Smith /*@C 518aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 51947c6ae99SBarry Smith 520aa219208SBarry Smith Logically Collective on DMDA 52147c6ae99SBarry Smith 52247c6ae99SBarry Smith Input Parameters: 523aa219208SBarry Smith + da - the DMDA object 524aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith Level: developer 52747c6ae99SBarry Smith 528aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 52947c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 53047c6ae99SBarry Smith 531aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills() 53247c6ae99SBarry Smith @*/ 5337087cfbeSBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*)) 53447c6ae99SBarry Smith { 53547c6ae99SBarry Smith PetscFunctionBegin; 53647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 53747c6ae99SBarry Smith da->ops->getmatrix = f; 53847c6ae99SBarry Smith PetscFunctionReturn(0); 53947c6ae99SBarry Smith } 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith #undef __FUNCT__ 54294013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges" 54347c6ae99SBarry Smith /* 54447c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 54547c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 54847c6ae99SBarry Smith */ 54994013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 55047c6ae99SBarry Smith { 55147c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 55247c6ae99SBarry Smith PetscErrorCode ierr; 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith PetscFunctionBegin; 55547c6ae99SBarry Smith if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 55647c6ae99SBarry Smith if (ratio == 1) { 55747c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 55847c6ae99SBarry Smith PetscFunctionReturn(0); 55947c6ae99SBarry Smith } 56047c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 56147c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 56247c6ae99SBarry Smith for (i=0; i<m; i++) { 56347c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 56447c6ae99SBarry Smith if (i == m-1) lf[i] = want; 56547c6ae99SBarry Smith else { 5667aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 5677aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 5687aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 5697aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 5707aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 5717aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 5727aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 5737aca7175SJed Brown /* Make sure all constraints are satisfied */ 5747aca7175SJed Brown if (want < 0 || want > remaining 5757aca7175SJed Brown || ((startf+want)/ratio < nextc - stencil_width) 5767aca7175SJed Brown || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) 5777aca7175SJed Brown SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 57847c6ae99SBarry Smith } 57947c6ae99SBarry Smith lf[i] = want; 58047c6ae99SBarry Smith startc += lc[i]; 58147c6ae99SBarry Smith startf += lf[i]; 58247c6ae99SBarry Smith remaining -= lf[i]; 58347c6ae99SBarry Smith } 58447c6ae99SBarry Smith PetscFunctionReturn(0); 58547c6ae99SBarry Smith } 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith #undef __FUNCT__ 5882be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges" 5892be375d4SJed Brown /* 5902be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 5912be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 5922be375d4SJed Brown 5932be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 5942be375d4SJed Brown */ 5952be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 5962be375d4SJed Brown { 5972be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 5982be375d4SJed Brown PetscErrorCode ierr; 5992be375d4SJed Brown 6002be375d4SJed Brown PetscFunctionBegin; 6012be375d4SJed Brown if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 6022be375d4SJed Brown if (ratio == 1) { 6032be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 6042be375d4SJed Brown PetscFunctionReturn(0); 6052be375d4SJed Brown } 6062be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 6072be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 6082be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 6092be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 6102be375d4SJed Brown if (i == m-1) lc[i] = want; 6112be375d4SJed Brown else { 6122be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 6132be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 6142be375d4SJed Brown * node is within one stencil width. */ 6152be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 6162be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 6172be375d4SJed Brown * fine node is within one stencil width. */ 6182be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 6192be375d4SJed Brown if (want < 0 || want > remaining 6202be375d4SJed Brown || (nextf/ratio < startc+want-stencil_width) 6212be375d4SJed Brown || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) 6222be375d4SJed Brown SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 6232be375d4SJed Brown } 6242be375d4SJed Brown lc[i] = want; 6252be375d4SJed Brown startc += lc[i]; 6262be375d4SJed Brown startf += lf[i]; 6272be375d4SJed Brown remaining -= lc[i]; 6282be375d4SJed Brown } 6292be375d4SJed Brown PetscFunctionReturn(0); 6302be375d4SJed Brown } 6312be375d4SJed Brown 6322be375d4SJed Brown #undef __FUNCT__ 63394013140SBarry Smith #define __FUNCT__ "DMRefine_DA" 6347087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 63547c6ae99SBarry Smith { 63647c6ae99SBarry Smith PetscErrorCode ierr; 63747c6ae99SBarry Smith PetscInt M,N,P; 6389a42bb27SBarry Smith DM da2; 63947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 64047c6ae99SBarry Smith 64147c6ae99SBarry Smith PetscFunctionBegin; 64247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 64347c6ae99SBarry Smith PetscValidPointer(daref,3); 64447c6ae99SBarry Smith 6451321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 64647c6ae99SBarry Smith M = dd->refine_x*dd->M; 64747c6ae99SBarry Smith } else { 64847c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 64947c6ae99SBarry Smith } 6501321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 65147c6ae99SBarry Smith N = dd->refine_y*dd->N; 65247c6ae99SBarry Smith } else { 65347c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 65447c6ae99SBarry Smith } 6551321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 65647c6ae99SBarry Smith P = dd->refine_z*dd->P; 65747c6ae99SBarry Smith } else { 65847c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 65947c6ae99SBarry Smith } 660aa219208SBarry Smith ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr); 66147c6ae99SBarry Smith if (dd->dim == 3) { 66247c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 66347c6ae99SBarry Smith ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 6641321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6651321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 6661321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 6671321219cSEthan Coon ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 66847c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 66947c6ae99SBarry Smith } else if (dd->dim == 2) { 67047c6ae99SBarry Smith PetscInt *lx,*ly; 67147c6ae99SBarry Smith ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 6721321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6731321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 6741321219cSEthan Coon ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 67547c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 67647c6ae99SBarry Smith } else if (dd->dim == 1) { 67747c6ae99SBarry Smith PetscInt *lx; 67847c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 6791321219cSEthan Coon ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 6801321219cSEthan Coon ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 68147c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 68247c6ae99SBarry Smith } 68347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 68447c6ae99SBarry Smith 68547c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 68647c6ae99SBarry Smith da2->ops->getmatrix = da->ops->getmatrix; 68747c6ae99SBarry Smith da2->ops->getinterpolation = da->ops->getinterpolation; 68847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 68947c6ae99SBarry Smith dd2->interptype = dd->interptype; 69047c6ae99SBarry Smith 69147c6ae99SBarry Smith /* copy fill information if given */ 69247c6ae99SBarry Smith if (dd->dfill) { 69347c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 69447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 69547c6ae99SBarry Smith } 69647c6ae99SBarry Smith if (dd->ofill) { 69747c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 69847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith /* copy the refine information */ 70147c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 70247c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 70347c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith /* copy vector type information */ 70647c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 70747c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 708ddcf8b74SDave May 709ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 710ddcf8b74SDave May if (dd->coordinates) { 711ddcf8b74SDave May DM cdaf,cdac; 712ddcf8b74SDave May Vec coordsc,coordsf; 713ddcf8b74SDave May Mat II; 714ddcf8b74SDave May 715ddcf8b74SDave May ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr); 716ddcf8b74SDave May ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr); 717b61d3410SDave May ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr); 718b61d3410SDave May /* force creation of the coordinate vector */ 719b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 720ddcf8b74SDave May ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 721ddcf8b74SDave May ierr = DMGetInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr); 722ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 723fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 724ddcf8b74SDave May } 72547c6ae99SBarry Smith *daref = da2; 72647c6ae99SBarry Smith PetscFunctionReturn(0); 72747c6ae99SBarry Smith } 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith #undef __FUNCT__ 73094013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA" 7317087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 73247c6ae99SBarry Smith { 73347c6ae99SBarry Smith PetscErrorCode ierr; 73447c6ae99SBarry Smith PetscInt M,N,P; 7359a42bb27SBarry Smith DM da2; 73647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 73747c6ae99SBarry Smith 73847c6ae99SBarry Smith PetscFunctionBegin; 73947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 74047c6ae99SBarry Smith PetscValidPointer(daref,3); 74147c6ae99SBarry Smith 7421321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 74347c6ae99SBarry Smith M = dd->M / dd->refine_x; 74447c6ae99SBarry Smith } else { 74547c6ae99SBarry Smith M = 1 + (dd->M - 1) / dd->refine_x; 74647c6ae99SBarry Smith } 7471321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 74847c6ae99SBarry Smith N = dd->N / dd->refine_y; 74947c6ae99SBarry Smith } else { 75047c6ae99SBarry Smith N = 1 + (dd->N - 1) / dd->refine_y; 75147c6ae99SBarry Smith } 7521321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 75347c6ae99SBarry Smith P = dd->P / dd->refine_z; 75447c6ae99SBarry Smith } else { 75547c6ae99SBarry Smith P = 1 + (dd->P - 1) / dd->refine_z; 75647c6ae99SBarry Smith } 7572be375d4SJed Brown ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr); 75847c6ae99SBarry Smith if (dd->dim == 3) { 7592be375d4SJed Brown PetscInt *lx,*ly,*lz; 7602be375d4SJed Brown ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 7611321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7621321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 7631321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 7641321219cSEthan Coon ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 7652be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 76647c6ae99SBarry Smith } else if (dd->dim == 2) { 7672be375d4SJed Brown PetscInt *lx,*ly; 7682be375d4SJed Brown ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 7691321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7701321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 7711321219cSEthan Coon ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 7722be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 77347c6ae99SBarry Smith } else if (dd->dim == 1) { 7742be375d4SJed Brown PetscInt *lx; 7752be375d4SJed Brown ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 7761321219cSEthan Coon ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 7771321219cSEthan Coon ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 7782be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 77947c6ae99SBarry Smith } 78047c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 78347c6ae99SBarry Smith da2->ops->getmatrix = da->ops->getmatrix; 78447c6ae99SBarry Smith da2->ops->getinterpolation = da->ops->getinterpolation; 78547c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 78647c6ae99SBarry Smith dd2->interptype = dd->interptype; 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith /* copy fill information if given */ 78947c6ae99SBarry Smith if (dd->dfill) { 79047c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 79147c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith if (dd->ofill) { 79447c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 79547c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 79647c6ae99SBarry Smith } 79747c6ae99SBarry Smith /* copy the refine information */ 79847c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 79947c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 80047c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith /* copy vector type information */ 80347c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 80447c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 80547c6ae99SBarry Smith 806ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 807ddcf8b74SDave May if (dd->coordinates) { 808ddcf8b74SDave May DM cdaf,cdac; 809ddcf8b74SDave May Vec coordsc,coordsf; 810ddcf8b74SDave May VecScatter inject; 811ddcf8b74SDave May 812ddcf8b74SDave May ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr); 813ddcf8b74SDave May ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr); 814b61d3410SDave May ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr); 815b61d3410SDave May /* force creation of the coordinate vector */ 816b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 817ddcf8b74SDave May ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 818ddcf8b74SDave May 819ddcf8b74SDave May ierr = DMGetInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 820ddcf8b74SDave May ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 821ddcf8b74SDave May ierr = VecScatterEnd(inject ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 822fcfd50ebSBarry Smith ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 823ddcf8b74SDave May } 82447c6ae99SBarry Smith *daref = da2; 82547c6ae99SBarry Smith PetscFunctionReturn(0); 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith 82847c6ae99SBarry Smith #undef __FUNCT__ 82994013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA" 8307087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 83147c6ae99SBarry Smith { 83247c6ae99SBarry Smith PetscErrorCode ierr; 83347c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 83447c6ae99SBarry Smith 83547c6ae99SBarry Smith PetscFunctionBegin; 83647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 83747c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 83847c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 83947c6ae99SBarry Smith PetscValidPointer(daf,3); 84047c6ae99SBarry Smith 841aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 84247c6ae99SBarry Smith ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr); 84347c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 844aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith n = nlevels; 84747c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr); 84847c6ae99SBarry Smith n = nlevels; 84947c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr); 85047c6ae99SBarry Smith n = nlevels; 85147c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr); 85247c6ae99SBarry Smith 853aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 85494013140SBarry Smith ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr); 85547c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 856aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 85794013140SBarry Smith ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr); 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 86047c6ae99SBarry Smith PetscFunctionReturn(0); 86147c6ae99SBarry Smith } 86247c6ae99SBarry Smith 86347c6ae99SBarry Smith #undef __FUNCT__ 86494013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA" 8657087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 86647c6ae99SBarry Smith { 86747c6ae99SBarry Smith PetscErrorCode ierr; 86847c6ae99SBarry Smith PetscInt i; 86947c6ae99SBarry Smith 87047c6ae99SBarry Smith PetscFunctionBegin; 87147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 87247c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 87347c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 87447c6ae99SBarry Smith PetscValidPointer(dac,3); 87594013140SBarry Smith ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr); 87647c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 87794013140SBarry Smith ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr); 87847c6ae99SBarry Smith } 87947c6ae99SBarry Smith PetscFunctionReturn(0); 88047c6ae99SBarry Smith } 881