147c6ae99SBarry Smith #define PETSCDM_DLL 2e1589f56SBarry Smith #include "private/daimpl.h" /*I "petscdm.h" I*/ 347c6ae99SBarry Smith 447c6ae99SBarry Smith 547c6ae99SBarry Smith #undef __FUNCT__ 6aa219208SBarry Smith #define __FUNCT__ "DMDASetDim" 747c6ae99SBarry Smith /*@ 8aa219208SBarry Smith DMDASetDim - Sets the dimension 947c6ae99SBarry Smith 10aa219208SBarry Smith Collective on DMDA 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith Input Parameters: 13aa219208SBarry Smith + da - the DMDA 1447c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE) 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith Level: intermediate 1747c6ae99SBarry Smith 18aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes() 1947c6ae99SBarry Smith @*/ 20aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetDim(DM da, PetscInt dim) 2147c6ae99SBarry Smith { 2247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith PetscFunctionBegin; 2547c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 26aa219208SBarry 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); 2747c6ae99SBarry Smith dd->dim = dim; 2847c6ae99SBarry Smith PetscFunctionReturn(0); 2947c6ae99SBarry Smith } 3047c6ae99SBarry Smith 3147c6ae99SBarry Smith #undef __FUNCT__ 32aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes" 3347c6ae99SBarry Smith /*@ 34aa219208SBarry Smith DMDASetSizes - Sets the global sizes 3547c6ae99SBarry Smith 36aa219208SBarry Smith Logically Collective on DMDA 3747c6ae99SBarry Smith 3847c6ae99SBarry Smith Input Parameters: 39aa219208SBarry Smith + da - the DMDA 4047c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE) 4147c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE) 4247c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE) 4347c6ae99SBarry Smith 4447c6ae99SBarry Smith Level: intermediate 4547c6ae99SBarry Smith 46aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership() 4747c6ae99SBarry Smith @*/ 48aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 4947c6ae99SBarry Smith { 5047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith PetscFunctionBegin; 5347c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 5447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,M,2); 5547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,N,3); 5647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,P,4); 5747c6ae99SBarry Smith 5847c6ae99SBarry Smith dd->M = M; 5947c6ae99SBarry Smith dd->N = N; 6047c6ae99SBarry Smith dd->P = P; 6147c6ae99SBarry Smith PetscFunctionReturn(0); 6247c6ae99SBarry Smith } 6347c6ae99SBarry Smith 6447c6ae99SBarry Smith #undef __FUNCT__ 65aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs" 6647c6ae99SBarry Smith /*@ 67aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 6847c6ae99SBarry Smith 69aa219208SBarry Smith Logically Collective on DMDA 7047c6ae99SBarry Smith 7147c6ae99SBarry Smith Input Parameters: 72aa219208SBarry Smith + da - the DMDA 7347c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE) 7447c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE) 7547c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE) 7647c6ae99SBarry Smith 7747c6ae99SBarry Smith Level: intermediate 7847c6ae99SBarry Smith 79aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership() 8047c6ae99SBarry Smith @*/ 81aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 8247c6ae99SBarry Smith { 8347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith PetscFunctionBegin; 8647c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 8747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 8847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 8947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 9047c6ae99SBarry Smith dd->m = m; 9147c6ae99SBarry Smith dd->n = n; 9247c6ae99SBarry Smith dd->p = p; 9347c6ae99SBarry Smith PetscFunctionReturn(0); 9447c6ae99SBarry Smith } 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith #undef __FUNCT__ 97aa219208SBarry Smith #define __FUNCT__ "DMDASetPeriodicity" 9847c6ae99SBarry Smith /*@ 99aa219208SBarry Smith DMDASetPeriodicity - Sets the type of periodicity 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith Not collective 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith Input Parameter: 104aa219208SBarry Smith + da - The DMDA 105aa219208SBarry Smith - ptype - One of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, DMDA_ZPERIODIC, DMDA_XYPERIODIC, DMDA_XZPERIODIC, DMDA_YZPERIODIC, or DMDA_XYZPERIODIC 10647c6ae99SBarry Smith 10747c6ae99SBarry Smith Level: intermediate 10847c6ae99SBarry Smith 10947c6ae99SBarry Smith .keywords: distributed array, periodicity 110aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDAPeriodicType 11147c6ae99SBarry Smith @*/ 112aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetPeriodicity(DM da, DMDAPeriodicType ptype) 11347c6ae99SBarry Smith { 11447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 11547c6ae99SBarry Smith 11647c6ae99SBarry Smith PetscFunctionBegin; 11747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 11847c6ae99SBarry Smith dd->wrap = ptype; 11947c6ae99SBarry Smith PetscFunctionReturn(0); 12047c6ae99SBarry Smith } 12147c6ae99SBarry Smith 12247c6ae99SBarry Smith #undef __FUNCT__ 123aa219208SBarry Smith #define __FUNCT__ "DMDASetDof" 12447c6ae99SBarry Smith /*@ 125aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith Not collective 12847c6ae99SBarry Smith 12947c6ae99SBarry Smith Input Parameter: 130aa219208SBarry Smith + da - The DMDA 13147c6ae99SBarry Smith - dof - Number of degrees of freedom 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith Level: intermediate 13447c6ae99SBarry Smith 13547c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 136aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 13747c6ae99SBarry Smith @*/ 138aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetDof(DM da, int dof) 13947c6ae99SBarry Smith { 14047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith PetscFunctionBegin; 14347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 14447c6ae99SBarry Smith dd->w = dof; 14547c6ae99SBarry Smith PetscFunctionReturn(0); 14647c6ae99SBarry Smith } 14747c6ae99SBarry Smith 14847c6ae99SBarry Smith #undef __FUNCT__ 149aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType" 15047c6ae99SBarry Smith /*@ 151aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 15247c6ae99SBarry Smith 153aa219208SBarry Smith Logically Collective on DMDA 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith Input Parameter: 156aa219208SBarry Smith + da - The DMDA 157aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 15847c6ae99SBarry Smith 15947c6ae99SBarry Smith Level: intermediate 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith .keywords: distributed array, stencil 162aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 16347c6ae99SBarry Smith @*/ 164aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilType(DM da, DMDAStencilType stype) 16547c6ae99SBarry Smith { 16647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 16747c6ae99SBarry Smith 16847c6ae99SBarry Smith PetscFunctionBegin; 16947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 17047c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 17147c6ae99SBarry Smith dd->stencil_type = stype; 17247c6ae99SBarry Smith PetscFunctionReturn(0); 17347c6ae99SBarry Smith } 17447c6ae99SBarry Smith 17547c6ae99SBarry Smith #undef __FUNCT__ 176aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth" 17747c6ae99SBarry Smith /*@ 178aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 17947c6ae99SBarry Smith 180aa219208SBarry Smith Logically Collective on DMDA 18147c6ae99SBarry Smith 18247c6ae99SBarry Smith Input Parameter: 183aa219208SBarry Smith + da - The DMDA 18447c6ae99SBarry Smith - width - The stencil width 18547c6ae99SBarry Smith 18647c6ae99SBarry Smith Level: intermediate 18747c6ae99SBarry Smith 18847c6ae99SBarry Smith .keywords: distributed array, stencil 189aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 19047c6ae99SBarry Smith @*/ 191aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilWidth(DM da, PetscInt width) 19247c6ae99SBarry Smith { 19347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 19447c6ae99SBarry Smith 19547c6ae99SBarry Smith PetscFunctionBegin; 19647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 19747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 19847c6ae99SBarry Smith dd->s = width; 19947c6ae99SBarry Smith PetscFunctionReturn(0); 20047c6ae99SBarry Smith } 20147c6ae99SBarry Smith 20247c6ae99SBarry Smith #undef __FUNCT__ 203aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 204aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 20547c6ae99SBarry Smith { 20647c6ae99SBarry Smith PetscInt i,sum; 20747c6ae99SBarry Smith 20847c6ae99SBarry Smith PetscFunctionBegin; 20947c6ae99SBarry Smith if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 21047c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 21147c6ae99SBarry Smith if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 21247c6ae99SBarry Smith PetscFunctionReturn(0); 21347c6ae99SBarry Smith } 21447c6ae99SBarry Smith 21547c6ae99SBarry Smith #undef __FUNCT__ 216aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges" 21747c6ae99SBarry Smith /*@ 218aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 21947c6ae99SBarry Smith 220aa219208SBarry Smith Logically Collective on DMDA 22147c6ae99SBarry Smith 22247c6ae99SBarry Smith Input Parameter: 223aa219208SBarry Smith + da - The DMDA 22447c6ae99SBarry 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 22547c6ae99SBarry 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 22647c6ae99SBarry 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. 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith Level: intermediate 22947c6ae99SBarry Smith 23047c6ae99SBarry Smith .keywords: distributed array 231aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 23247c6ae99SBarry Smith @*/ 233aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 23447c6ae99SBarry Smith { 23547c6ae99SBarry Smith PetscErrorCode ierr; 23647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 23747c6ae99SBarry Smith 23847c6ae99SBarry Smith PetscFunctionBegin; 23947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 24047c6ae99SBarry Smith if (lx) { 24147c6ae99SBarry Smith if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 242aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 24347c6ae99SBarry Smith if (!dd->lx) { 24447c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 24547c6ae99SBarry Smith } 24647c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 24747c6ae99SBarry Smith } 24847c6ae99SBarry Smith if (ly) { 24947c6ae99SBarry Smith if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 250aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 25147c6ae99SBarry Smith if (!dd->ly) { 25247c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 25347c6ae99SBarry Smith } 25447c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 25547c6ae99SBarry Smith } 25647c6ae99SBarry Smith if (lz) { 25747c6ae99SBarry Smith if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 258aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 25947c6ae99SBarry Smith if (!dd->lz) { 26047c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 26147c6ae99SBarry Smith } 26247c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 26347c6ae99SBarry Smith } 26447c6ae99SBarry Smith PetscFunctionReturn(0); 26547c6ae99SBarry Smith } 26647c6ae99SBarry Smith 26747c6ae99SBarry Smith #undef __FUNCT__ 268aa219208SBarry Smith #define __FUNCT__ "DMDACreateOwnershipRanges" 26947c6ae99SBarry Smith /* 270aa219208SBarry Smith Ensure that da->lx, ly, and lz exist. Collective on DMDA. 27147c6ae99SBarry Smith */ 272aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDACreateOwnershipRanges(DM da) 27347c6ae99SBarry Smith { 27447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 27547c6ae99SBarry Smith PetscErrorCode ierr; 27647c6ae99SBarry Smith PetscInt n; 27747c6ae99SBarry Smith MPI_Comm comm; 27847c6ae99SBarry Smith PetscMPIInt size; 27947c6ae99SBarry Smith 28047c6ae99SBarry Smith PetscFunctionBegin; 28147c6ae99SBarry Smith if (!dd->lx) { 28247c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr); 283aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr); 28447c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 28547c6ae99SBarry Smith n = (dd->xe-dd->xs)/dd->w; 28647c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr); 28747c6ae99SBarry Smith } 28847c6ae99SBarry Smith if (dd->dim > 1 && !dd->ly) { 28947c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr); 290aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr); 29147c6ae99SBarry Smith n = dd->ye-dd->ys; 29247c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr); 29347c6ae99SBarry Smith } 29447c6ae99SBarry Smith if (dd->dim > 2 && !dd->lz) { 29547c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr); 296aa219208SBarry Smith ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr); 29747c6ae99SBarry Smith n = dd->ze-dd->zs; 29847c6ae99SBarry Smith ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);CHKERRQ(ierr); 29947c6ae99SBarry Smith } 30047c6ae99SBarry Smith PetscFunctionReturn(0); 30147c6ae99SBarry Smith } 30247c6ae99SBarry Smith 30347c6ae99SBarry Smith #undef __FUNCT__ 304aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType" 30547c6ae99SBarry Smith /*@ 306aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 30794013140SBarry Smith returned by DMGetInterpolation() 30847c6ae99SBarry Smith 309aa219208SBarry Smith Logically Collective on DMDA 31047c6ae99SBarry Smith 31147c6ae99SBarry Smith Input Parameter: 31247c6ae99SBarry Smith + da - initial distributed array 313aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith Level: intermediate 31647c6ae99SBarry Smith 317aa219208SBarry Smith Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation() 31847c6ae99SBarry Smith 31947c6ae99SBarry Smith .keywords: distributed array, interpolation 32047c6ae99SBarry Smith 321aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 32247c6ae99SBarry Smith @*/ 323aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 32447c6ae99SBarry Smith { 32547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 32647c6ae99SBarry Smith 32747c6ae99SBarry Smith PetscFunctionBegin; 32847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 32947c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 33047c6ae99SBarry Smith dd->interptype = ctype; 33147c6ae99SBarry Smith PetscFunctionReturn(0); 33247c6ae99SBarry Smith } 33347c6ae99SBarry Smith 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith #undef __FUNCT__ 336aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors" 33747c6ae99SBarry Smith /*@C 338aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 33947c6ae99SBarry Smith processes neighbors. 34047c6ae99SBarry Smith 34147c6ae99SBarry Smith Not Collective 34247c6ae99SBarry Smith 34347c6ae99SBarry Smith Input Parameter: 344aa219208SBarry Smith . da - the DMDA object 34547c6ae99SBarry Smith 34647c6ae99SBarry Smith Output Parameters: 34747c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 34847c6ae99SBarry Smith this process itself is in the list 34947c6ae99SBarry Smith 35047c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 35147c6ae99SBarry Smith Not supported in 1d 352aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 35347c6ae99SBarry Smith 35447c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 35547c6ae99SBarry Smith 35647c6ae99SBarry Smith Level: intermediate 35747c6ae99SBarry Smith 35847c6ae99SBarry Smith @*/ 359aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 36047c6ae99SBarry Smith { 36147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 36247c6ae99SBarry Smith PetscFunctionBegin; 36347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 36447c6ae99SBarry Smith *ranks = dd->neighbors; 36547c6ae99SBarry Smith PetscFunctionReturn(0); 36647c6ae99SBarry Smith } 36747c6ae99SBarry Smith 36847c6ae99SBarry Smith /*@C 369aa219208SBarry Smith DMDASetElementType - Sets the element type to be returned by DMGetElements() 37047c6ae99SBarry Smith 37147c6ae99SBarry Smith Not Collective 37247c6ae99SBarry Smith 37347c6ae99SBarry Smith Input Parameter: 374aa219208SBarry Smith . da - the DMDA object 37547c6ae99SBarry Smith 37647c6ae99SBarry Smith Output Parameters: 377aa219208SBarry Smith . etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1 37847c6ae99SBarry Smith 37947c6ae99SBarry Smith Level: intermediate 38047c6ae99SBarry Smith 381aa219208SBarry Smith .seealso: DMDAElementType, DMDAGetElementType(), DMGetElements(), DMRestoreElements() 38247c6ae99SBarry Smith @*/ 38347c6ae99SBarry Smith #undef __FUNCT__ 384aa219208SBarry Smith #define __FUNCT__ "DMDASetElementType" 385aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetElementType(DM da, DMDAElementType etype) 38647c6ae99SBarry Smith { 38747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 38847c6ae99SBarry Smith PetscErrorCode ierr; 38947c6ae99SBarry Smith 39047c6ae99SBarry Smith PetscFunctionBegin; 39147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 39247c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,etype,2); 39347c6ae99SBarry Smith if (dd->elementtype != etype) { 39447c6ae99SBarry Smith ierr = PetscFree(dd->e);CHKERRQ(ierr); 39547c6ae99SBarry Smith dd->elementtype = etype; 39647c6ae99SBarry Smith dd->ne = 0; 39747c6ae99SBarry Smith dd->e = PETSC_NULL; 39847c6ae99SBarry Smith } 39947c6ae99SBarry Smith PetscFunctionReturn(0); 40047c6ae99SBarry Smith } 40147c6ae99SBarry Smith 40247c6ae99SBarry Smith /*@C 403aa219208SBarry Smith DMDAGetElementType - Gets the element type to be returned by DMGetElements() 40447c6ae99SBarry Smith 40547c6ae99SBarry Smith Not Collective 40647c6ae99SBarry Smith 40747c6ae99SBarry Smith Input Parameter: 408aa219208SBarry Smith . da - the DMDA object 40947c6ae99SBarry Smith 41047c6ae99SBarry Smith Output Parameters: 411aa219208SBarry Smith . etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1 41247c6ae99SBarry Smith 41347c6ae99SBarry Smith Level: intermediate 41447c6ae99SBarry Smith 415aa219208SBarry Smith .seealso: DMDAElementType, DMDASetElementType(), DMGetElements(), DMRestoreElements() 41647c6ae99SBarry Smith @*/ 41747c6ae99SBarry Smith #undef __FUNCT__ 418aa219208SBarry Smith #define __FUNCT__ "DMDAGetElementType" 419aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetElementType(DM da, DMDAElementType *etype) 42047c6ae99SBarry Smith { 42147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 42247c6ae99SBarry Smith PetscFunctionBegin; 42347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 42447c6ae99SBarry Smith PetscValidPointer(etype,2); 42547c6ae99SBarry Smith *etype = dd->elementtype; 42647c6ae99SBarry Smith PetscFunctionReturn(0); 42747c6ae99SBarry Smith } 42847c6ae99SBarry Smith 42947c6ae99SBarry Smith #undef __FUNCT__ 43047c6ae99SBarry Smith #define __FUNCT__ "DMGetElements" 43147c6ae99SBarry Smith /*@C 43247c6ae99SBarry Smith DMGetElements - Gets an array containing the indices (in local coordinates) 43347c6ae99SBarry Smith of all the local elements 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith Not Collective 43647c6ae99SBarry Smith 43747c6ae99SBarry Smith Input Parameter: 43847c6ae99SBarry Smith . dm - the DM object 43947c6ae99SBarry Smith 44047c6ae99SBarry Smith Output Parameters: 441*454e267fSLisandro Dalcin + nel - number of local elements 442*454e267fSLisandro Dalcin . nen - number of element nodes 44347c6ae99SBarry Smith - e - the indices of the elements vertices 44447c6ae99SBarry Smith 44547c6ae99SBarry Smith Level: intermediate 44647c6ae99SBarry Smith 44747c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMRestoreElements() 44847c6ae99SBarry Smith @*/ 449*454e267fSLisandro Dalcin PetscErrorCode PETSCDM_DLLEXPORT DMGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 45047c6ae99SBarry Smith { 45147c6ae99SBarry Smith PetscErrorCode ierr; 45247c6ae99SBarry Smith PetscFunctionBegin; 45347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 454*454e267fSLisandro Dalcin PetscValidIntPointer(nel,2); 455*454e267fSLisandro Dalcin PetscValidIntPointer(nen,3); 456*454e267fSLisandro Dalcin PetscValidPointer(e,4); 457*454e267fSLisandro Dalcin ierr = (dm->ops->getelements)(dm,nel,nen,e);CHKERRQ(ierr); 45847c6ae99SBarry Smith PetscFunctionReturn(0); 45947c6ae99SBarry Smith } 46047c6ae99SBarry Smith 46147c6ae99SBarry Smith #undef __FUNCT__ 46247c6ae99SBarry Smith #define __FUNCT__ "DMRestoreElements" 46347c6ae99SBarry Smith /*@C 46447c6ae99SBarry Smith DMRestoreElements - Returns an array containing the indices (in local coordinates) 46547c6ae99SBarry Smith of all the local elements obtained with DMGetElements() 46647c6ae99SBarry Smith 46747c6ae99SBarry Smith Not Collective 46847c6ae99SBarry Smith 46947c6ae99SBarry Smith Input Parameter: 47047c6ae99SBarry Smith + dm - the DM object 471*454e267fSLisandro Dalcin . nel - number of local elements 472*454e267fSLisandro Dalcin . nen - number of element nodes 47347c6ae99SBarry Smith - e - the indices of the elements vertices 47447c6ae99SBarry Smith 47547c6ae99SBarry Smith Level: intermediate 47647c6ae99SBarry Smith 47747c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMGetElements() 47847c6ae99SBarry Smith @*/ 479*454e267fSLisandro Dalcin PetscErrorCode PETSCDM_DLLEXPORT DMRestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 48047c6ae99SBarry Smith { 48147c6ae99SBarry Smith PetscErrorCode ierr; 48247c6ae99SBarry Smith PetscFunctionBegin; 48347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 484*454e267fSLisandro Dalcin PetscValidIntPointer(nel,2); 485*454e267fSLisandro Dalcin PetscValidIntPointer(nen,3); 486*454e267fSLisandro Dalcin PetscValidPointer(e,4); 48747c6ae99SBarry Smith if (dm->ops->restoreelements) { 488*454e267fSLisandro Dalcin ierr = (dm->ops->restoreelements)(dm,nel,nen,e);CHKERRQ(ierr); 48947c6ae99SBarry Smith } 49047c6ae99SBarry Smith PetscFunctionReturn(0); 49147c6ae99SBarry Smith } 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith #undef __FUNCT__ 494aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges" 49547c6ae99SBarry Smith /*@C 496aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith Not Collective 49947c6ae99SBarry Smith 50047c6ae99SBarry Smith Input Parameter: 501aa219208SBarry Smith . da - the DMDA object 50247c6ae99SBarry Smith 50347c6ae99SBarry Smith Output Parameter: 50447c6ae99SBarry Smith + lx - ownership along x direction (optional) 50547c6ae99SBarry Smith . ly - ownership along y direction (optional) 50647c6ae99SBarry Smith - lz - ownership along z direction (optional) 50747c6ae99SBarry Smith 50847c6ae99SBarry Smith Level: intermediate 50947c6ae99SBarry Smith 510aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 513aa219208SBarry Smith eighth arguments from DMDAGetInfo() 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 516aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 51747c6ae99SBarry Smith 518aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 51947c6ae99SBarry Smith @*/ 520aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 52147c6ae99SBarry Smith { 52247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 52347c6ae99SBarry Smith 52447c6ae99SBarry Smith PetscFunctionBegin; 52547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 52647c6ae99SBarry Smith if (lx) *lx = dd->lx; 52747c6ae99SBarry Smith if (ly) *ly = dd->ly; 52847c6ae99SBarry Smith if (lz) *lz = dd->lz; 52947c6ae99SBarry Smith PetscFunctionReturn(0); 53047c6ae99SBarry Smith } 53147c6ae99SBarry Smith 53247c6ae99SBarry Smith #undef __FUNCT__ 533aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor" 53447c6ae99SBarry Smith /*@ 535aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 53647c6ae99SBarry Smith 537aa219208SBarry Smith Logically Collective on DMDA 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith Input Parameters: 540aa219208SBarry Smith + da - the DMDA object 54147c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 54247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 54347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 54447c6ae99SBarry Smith 54547c6ae99SBarry Smith Options Database: 54647c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 54747c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 54847c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith Level: intermediate 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 55347c6ae99SBarry Smith 554aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 55547c6ae99SBarry Smith @*/ 556aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 55747c6ae99SBarry Smith { 55847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 55947c6ae99SBarry Smith 56047c6ae99SBarry Smith PetscFunctionBegin; 56147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 56247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 56347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 56447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 56547c6ae99SBarry Smith 56647c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 56747c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 56847c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 56947c6ae99SBarry Smith PetscFunctionReturn(0); 57047c6ae99SBarry Smith } 57147c6ae99SBarry Smith 57247c6ae99SBarry Smith #undef __FUNCT__ 573aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor" 57447c6ae99SBarry Smith /*@C 575aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 57647c6ae99SBarry Smith 57747c6ae99SBarry Smith Not Collective 57847c6ae99SBarry Smith 57947c6ae99SBarry Smith Input Parameter: 580aa219208SBarry Smith . da - the DMDA object 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith Output Parameters: 58347c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 58447c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 58547c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith Level: intermediate 58847c6ae99SBarry Smith 58947c6ae99SBarry Smith Notes: Pass PETSC_NULL for values you do not need 59047c6ae99SBarry Smith 591aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 59247c6ae99SBarry Smith @*/ 593aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 59447c6ae99SBarry Smith { 59547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 59647c6ae99SBarry Smith 59747c6ae99SBarry Smith PetscFunctionBegin; 59847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 59947c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 60047c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 60147c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 60247c6ae99SBarry Smith PetscFunctionReturn(0); 60347c6ae99SBarry Smith } 60447c6ae99SBarry Smith 60547c6ae99SBarry Smith #undef __FUNCT__ 606aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix" 60747c6ae99SBarry Smith /*@C 608aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 60947c6ae99SBarry Smith 610aa219208SBarry Smith Logically Collective on DMDA 61147c6ae99SBarry Smith 61247c6ae99SBarry Smith Input Parameters: 613aa219208SBarry Smith + da - the DMDA object 614aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 61547c6ae99SBarry Smith 61647c6ae99SBarry Smith Level: developer 61747c6ae99SBarry Smith 618aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 61947c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 62047c6ae99SBarry Smith 621aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills() 62247c6ae99SBarry Smith @*/ 623aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*)) 62447c6ae99SBarry Smith { 62547c6ae99SBarry Smith PetscFunctionBegin; 62647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 62747c6ae99SBarry Smith da->ops->getmatrix = f; 62847c6ae99SBarry Smith PetscFunctionReturn(0); 62947c6ae99SBarry Smith } 63047c6ae99SBarry Smith 63147c6ae99SBarry Smith #undef __FUNCT__ 63294013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges" 63347c6ae99SBarry Smith /* 63447c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 63547c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 63647c6ae99SBarry Smith 63747c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 63847c6ae99SBarry Smith */ 63994013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 64047c6ae99SBarry Smith { 64147c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 64247c6ae99SBarry Smith PetscErrorCode ierr; 64347c6ae99SBarry Smith 64447c6ae99SBarry Smith PetscFunctionBegin; 64547c6ae99SBarry Smith if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 64647c6ae99SBarry Smith if (ratio == 1) { 64747c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 64847c6ae99SBarry Smith PetscFunctionReturn(0); 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 65147c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 65247c6ae99SBarry Smith for (i=0; i<m; i++) { 65347c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 65447c6ae99SBarry Smith if (i == m-1) lf[i] = want; 65547c6ae99SBarry Smith else { 6567aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 6577aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 6587aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 6597aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 6607aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 6617aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 6627aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 6637aca7175SJed Brown /* Make sure all constraints are satisfied */ 6647aca7175SJed Brown if (want < 0 || want > remaining 6657aca7175SJed Brown || ((startf+want)/ratio < nextc - stencil_width) 6667aca7175SJed Brown || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) 6677aca7175SJed Brown SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 66847c6ae99SBarry Smith } 66947c6ae99SBarry Smith lf[i] = want; 67047c6ae99SBarry Smith startc += lc[i]; 67147c6ae99SBarry Smith startf += lf[i]; 67247c6ae99SBarry Smith remaining -= lf[i]; 67347c6ae99SBarry Smith } 67447c6ae99SBarry Smith PetscFunctionReturn(0); 67547c6ae99SBarry Smith } 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith #undef __FUNCT__ 6782be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges" 6792be375d4SJed Brown /* 6802be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 6812be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 6822be375d4SJed Brown 6832be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 6842be375d4SJed Brown */ 6852be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 6862be375d4SJed Brown { 6872be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 6882be375d4SJed Brown PetscErrorCode ierr; 6892be375d4SJed Brown 6902be375d4SJed Brown PetscFunctionBegin; 6912be375d4SJed Brown if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 6922be375d4SJed Brown if (ratio == 1) { 6932be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 6942be375d4SJed Brown PetscFunctionReturn(0); 6952be375d4SJed Brown } 6962be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 6972be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 6982be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 6992be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 7002be375d4SJed Brown if (i == m-1) lc[i] = want; 7012be375d4SJed Brown else { 7022be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 7032be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 7042be375d4SJed Brown * node is within one stencil width. */ 7052be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 7062be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 7072be375d4SJed Brown * fine node is within one stencil width. */ 7082be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 7092be375d4SJed Brown if (want < 0 || want > remaining 7102be375d4SJed Brown || (nextf/ratio < startc+want-stencil_width) 7112be375d4SJed Brown || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) 7122be375d4SJed Brown SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 7132be375d4SJed Brown } 7142be375d4SJed Brown lc[i] = want; 7152be375d4SJed Brown startc += lc[i]; 7162be375d4SJed Brown startf += lf[i]; 7172be375d4SJed Brown remaining -= lc[i]; 7182be375d4SJed Brown } 7192be375d4SJed Brown PetscFunctionReturn(0); 7202be375d4SJed Brown } 7212be375d4SJed Brown 7222be375d4SJed Brown #undef __FUNCT__ 72394013140SBarry Smith #define __FUNCT__ "DMRefine_DA" 72494013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 72547c6ae99SBarry Smith { 72647c6ae99SBarry Smith PetscErrorCode ierr; 72747c6ae99SBarry Smith PetscInt M,N,P; 7289a42bb27SBarry Smith DM da2; 72947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith PetscFunctionBegin; 73247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 73347c6ae99SBarry Smith PetscValidPointer(daref,3); 73447c6ae99SBarry Smith 735aa219208SBarry Smith if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 73647c6ae99SBarry Smith M = dd->refine_x*dd->M; 73747c6ae99SBarry Smith } else { 73847c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 73947c6ae99SBarry Smith } 740aa219208SBarry Smith if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 74147c6ae99SBarry Smith N = dd->refine_y*dd->N; 74247c6ae99SBarry Smith } else { 74347c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 74447c6ae99SBarry Smith } 745aa219208SBarry Smith if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 74647c6ae99SBarry Smith P = dd->refine_z*dd->P; 74747c6ae99SBarry Smith } else { 74847c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 74947c6ae99SBarry Smith } 750aa219208SBarry Smith ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr); 75147c6ae99SBarry Smith if (dd->dim == 3) { 75247c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 75347c6ae99SBarry Smith ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 754aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 755aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 756aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 757aa219208SBarry Smith ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 75847c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 75947c6ae99SBarry Smith } else if (dd->dim == 2) { 76047c6ae99SBarry Smith PetscInt *lx,*ly; 76147c6ae99SBarry Smith ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 762aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 763aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 764aa219208SBarry Smith ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 76547c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 76647c6ae99SBarry Smith } else if (dd->dim == 1) { 76747c6ae99SBarry Smith PetscInt *lx; 76847c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 769aa219208SBarry Smith ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 770aa219208SBarry Smith ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 77147c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 77247c6ae99SBarry Smith } 77347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 77647c6ae99SBarry Smith da2->ops->getmatrix = da->ops->getmatrix; 77747c6ae99SBarry Smith da2->ops->getinterpolation = da->ops->getinterpolation; 77847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 77947c6ae99SBarry Smith dd2->interptype = dd->interptype; 78047c6ae99SBarry Smith 78147c6ae99SBarry Smith /* copy fill information if given */ 78247c6ae99SBarry Smith if (dd->dfill) { 78347c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 78447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith if (dd->ofill) { 78747c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 78847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 78947c6ae99SBarry Smith } 79047c6ae99SBarry Smith /* copy the refine information */ 79147c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 79247c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 79347c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 79447c6ae99SBarry Smith 79547c6ae99SBarry Smith /* copy vector type information */ 79647c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 79747c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 79847c6ae99SBarry Smith *daref = da2; 79947c6ae99SBarry Smith PetscFunctionReturn(0); 80047c6ae99SBarry Smith } 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith #undef __FUNCT__ 80394013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA" 80494013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 80547c6ae99SBarry Smith { 80647c6ae99SBarry Smith PetscErrorCode ierr; 80747c6ae99SBarry Smith PetscInt M,N,P; 8089a42bb27SBarry Smith DM da2; 80947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 81047c6ae99SBarry Smith 81147c6ae99SBarry Smith PetscFunctionBegin; 81247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 81347c6ae99SBarry Smith PetscValidPointer(daref,3); 81447c6ae99SBarry Smith 815aa219208SBarry Smith if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 81647c6ae99SBarry Smith M = dd->M / dd->refine_x; 81747c6ae99SBarry Smith } else { 81847c6ae99SBarry Smith M = 1 + (dd->M - 1) / dd->refine_x; 81947c6ae99SBarry Smith } 820aa219208SBarry Smith if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 82147c6ae99SBarry Smith N = dd->N / dd->refine_y; 82247c6ae99SBarry Smith } else { 82347c6ae99SBarry Smith N = 1 + (dd->N - 1) / dd->refine_y; 82447c6ae99SBarry Smith } 825aa219208SBarry Smith if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){ 82647c6ae99SBarry Smith P = dd->P / dd->refine_z; 82747c6ae99SBarry Smith } else { 82847c6ae99SBarry Smith P = 1 + (dd->P - 1) / dd->refine_z; 82947c6ae99SBarry Smith } 8302be375d4SJed Brown ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr); 83147c6ae99SBarry Smith if (dd->dim == 3) { 8322be375d4SJed Brown PetscInt *lx,*ly,*lz; 8332be375d4SJed Brown ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 8342be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 8352be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 8362be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 8372be375d4SJed Brown ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 8382be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 83947c6ae99SBarry Smith } else if (dd->dim == 2) { 8402be375d4SJed Brown PetscInt *lx,*ly; 8412be375d4SJed Brown ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 8422be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 8432be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 8442be375d4SJed Brown ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 8452be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 84647c6ae99SBarry Smith } else if (dd->dim == 1) { 8472be375d4SJed Brown PetscInt *lx; 8482be375d4SJed Brown ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 8492be375d4SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 8502be375d4SJed Brown ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 8512be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 85647c6ae99SBarry Smith da2->ops->getmatrix = da->ops->getmatrix; 85747c6ae99SBarry Smith da2->ops->getinterpolation = da->ops->getinterpolation; 85847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 85947c6ae99SBarry Smith dd2->interptype = dd->interptype; 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith /* copy fill information if given */ 86247c6ae99SBarry Smith if (dd->dfill) { 86347c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 86447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 86547c6ae99SBarry Smith } 86647c6ae99SBarry Smith if (dd->ofill) { 86747c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 86847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith /* copy the refine information */ 87147c6ae99SBarry Smith dd2->refine_x = dd->refine_x; 87247c6ae99SBarry Smith dd2->refine_y = dd->refine_y; 87347c6ae99SBarry Smith dd2->refine_z = dd->refine_z; 87447c6ae99SBarry Smith 87547c6ae99SBarry Smith /* copy vector type information */ 87647c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 87747c6ae99SBarry Smith ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith *daref = da2; 88047c6ae99SBarry Smith PetscFunctionReturn(0); 88147c6ae99SBarry Smith } 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith #undef __FUNCT__ 88494013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA" 88594013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 88647c6ae99SBarry Smith { 88747c6ae99SBarry Smith PetscErrorCode ierr; 88847c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith PetscFunctionBegin; 89147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 89247c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 89347c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 89447c6ae99SBarry Smith PetscValidPointer(daf,3); 89547c6ae99SBarry Smith 896aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 89747c6ae99SBarry Smith ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr); 89847c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 899aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 90047c6ae99SBarry Smith } 90147c6ae99SBarry Smith n = nlevels; 90247c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr); 90347c6ae99SBarry Smith n = nlevels; 90447c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr); 90547c6ae99SBarry Smith n = nlevels; 90647c6ae99SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr); 90747c6ae99SBarry Smith 908aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 90994013140SBarry Smith ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr); 91047c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 911aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 91294013140SBarry Smith ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr); 91347c6ae99SBarry Smith } 91447c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 91547c6ae99SBarry Smith PetscFunctionReturn(0); 91647c6ae99SBarry Smith } 91747c6ae99SBarry Smith 91847c6ae99SBarry Smith #undef __FUNCT__ 91994013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA" 92094013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 92147c6ae99SBarry Smith { 92247c6ae99SBarry Smith PetscErrorCode ierr; 92347c6ae99SBarry Smith PetscInt i; 92447c6ae99SBarry Smith 92547c6ae99SBarry Smith PetscFunctionBegin; 92647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 92747c6ae99SBarry Smith if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 92847c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 92947c6ae99SBarry Smith PetscValidPointer(dac,3); 93094013140SBarry Smith ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr); 93147c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 93294013140SBarry Smith ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr); 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith PetscFunctionReturn(0); 93547c6ae99SBarry Smith } 936