1b45d2f2cSJed Brown #include <petsc-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); 25*ce94432eSBarry Smith if (dd->dim > 0 && dim != dd->dim) SETERRQ2(PetscObjectComm((PetscObject)da),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); 56*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 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 @*/ 817087cfbeSBarry Smith PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 8247c6ae99SBarry Smith { 8347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 84e3f3e4b6SBarry Smith PetscErrorCode ierr; 8547c6ae99SBarry Smith 8647c6ae99SBarry Smith PetscFunctionBegin; 8747c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 8847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 8947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 9047c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 91*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 9247c6ae99SBarry Smith dd->m = m; 9347c6ae99SBarry Smith dd->n = n; 9447c6ae99SBarry Smith dd->p = p; 95e3f3e4b6SBarry Smith if (dd->dim == 2) { 96e3f3e4b6SBarry Smith PetscMPIInt size; 97*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 98e3f3e4b6SBarry Smith if ((dd->m > 0) && (dd->n < 0)) { 99e3f3e4b6SBarry Smith dd->n = size/dd->m; 100*ce94432eSBarry Smith if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size); 101e3f3e4b6SBarry Smith } 102e3f3e4b6SBarry Smith if ((dd->n > 0) && (dd->m < 0)) { 103e3f3e4b6SBarry Smith dd->m = size/dd->n; 104*ce94432eSBarry Smith if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size); 105e3f3e4b6SBarry Smith } 106e3f3e4b6SBarry Smith } 10747c6ae99SBarry Smith PetscFunctionReturn(0); 10847c6ae99SBarry Smith } 10947c6ae99SBarry Smith 11047c6ae99SBarry Smith #undef __FUNCT__ 1113e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType" 11247c6ae99SBarry Smith /*@ 1131321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith Not collective 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith Input Parameter: 118aa219208SBarry Smith + da - The DMDA 1191321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC 12047c6ae99SBarry Smith 12147c6ae99SBarry Smith Level: intermediate 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith .keywords: distributed array, periodicity 12496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType 12547c6ae99SBarry Smith @*/ 1261321219cSEthan Coon PetscErrorCode DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz) 12747c6ae99SBarry Smith { 12847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith PetscFunctionBegin; 13147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1325a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bx,2); 1335a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,by,3); 1345a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bz,4); 135*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 1361321219cSEthan Coon dd->bx = bx; 1371321219cSEthan Coon dd->by = by; 1381321219cSEthan Coon dd->bz = bz; 13947c6ae99SBarry Smith PetscFunctionReturn(0); 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith #undef __FUNCT__ 143aa219208SBarry Smith #define __FUNCT__ "DMDASetDof" 14447c6ae99SBarry Smith /*@ 145aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 14647c6ae99SBarry Smith 14747c6ae99SBarry Smith Not collective 14847c6ae99SBarry Smith 14947c6ae99SBarry Smith Input Parameter: 150aa219208SBarry Smith + da - The DMDA 15147c6ae99SBarry Smith - dof - Number of degrees of freedom 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith Level: intermediate 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 15696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 15747c6ae99SBarry Smith @*/ 15854cfb0beSLisandro Dalcin PetscErrorCode DMDASetDof(DM da, PetscInt dof) 15947c6ae99SBarry Smith { 16047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 16147c6ae99SBarry Smith 16247c6ae99SBarry Smith PetscFunctionBegin; 16347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 16454cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da,dof,2); 165*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 16647c6ae99SBarry Smith dd->w = dof; 1671411c6eeSJed Brown da->bs = dof; 16847c6ae99SBarry Smith PetscFunctionReturn(0); 16947c6ae99SBarry Smith } 17047c6ae99SBarry Smith 17147c6ae99SBarry Smith #undef __FUNCT__ 1727ddda789SPeter Brune #define __FUNCT__ "DMDAGetOverlap" 1737ddda789SPeter Brune /*@ 1747ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1757ddda789SPeter Brune 1767ddda789SPeter Brune Not collective 1777ddda789SPeter Brune 1787ddda789SPeter Brune Input Parameters: 1797ddda789SPeter Brune . da - The DMDA 1807ddda789SPeter Brune 1817ddda789SPeter Brune Output Parameters: 1827ddda789SPeter Brune + x - Overlap in the x direction 1837ddda789SPeter Brune . y - Overlap in the y direction 1847ddda789SPeter Brune - z - Overlap in the z direction 1857ddda789SPeter Brune 1867ddda789SPeter Brune Level: intermediate 1877ddda789SPeter Brune 1887ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 1897ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA 1907ddda789SPeter Brune @*/ 1917ddda789SPeter Brune PetscErrorCode DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z) 1927ddda789SPeter Brune { 1937ddda789SPeter Brune DM_DA *dd = (DM_DA*)da->data; 1947ddda789SPeter Brune 1957ddda789SPeter Brune PetscFunctionBegin; 1967ddda789SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 1977ddda789SPeter Brune if (x) *x = dd->xol; 1987ddda789SPeter Brune if (y) *y = dd->yol; 1997ddda789SPeter Brune if (z) *z = dd->zol; 2007ddda789SPeter Brune PetscFunctionReturn(0); 2017ddda789SPeter Brune } 2027ddda789SPeter Brune 2037ddda789SPeter Brune #undef __FUNCT__ 20488661749SPeter Brune #define __FUNCT__ "DMDASetOverlap" 20588661749SPeter Brune /*@ 20688661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 20788661749SPeter Brune 20888661749SPeter Brune Not collective 20988661749SPeter Brune 2107ddda789SPeter Brune Input Parameters: 21188661749SPeter Brune + da - The DMDA 2127ddda789SPeter Brune . x - Overlap in the x direction 2137ddda789SPeter Brune . y - Overlap in the y direction 2147ddda789SPeter Brune - z - Overlap in the z direction 21588661749SPeter Brune 21688661749SPeter Brune Level: intermediate 21788661749SPeter Brune 2187ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 2197ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA 22088661749SPeter Brune @*/ 2217ddda789SPeter Brune PetscErrorCode DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z) 22288661749SPeter Brune { 22388661749SPeter Brune DM_DA *dd = (DM_DA*)da->data; 22488661749SPeter Brune 22588661749SPeter Brune PetscFunctionBegin; 22688661749SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2277ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,x,2); 2287ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,y,3); 2297ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,z,4); 2307ddda789SPeter Brune dd->xol = x; 2317ddda789SPeter Brune dd->yol = y; 2327ddda789SPeter Brune dd->zol = z; 23388661749SPeter Brune PetscFunctionReturn(0); 23488661749SPeter Brune } 23588661749SPeter Brune 236d886c4f4SPeter Brune #undef __FUNCT__ 237d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset" 238d886c4f4SPeter Brune /*@ 239d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 240d886c4f4SPeter Brune 241d886c4f4SPeter Brune Collective on DA 242d886c4f4SPeter Brune 243d886c4f4SPeter Brune Input Parameter: 244d886c4f4SPeter Brune + da - The DMDA 245d886c4f4SPeter Brune . xo - The offset in the x direction 246d886c4f4SPeter Brune . yo - The offset in the y direction 247d886c4f4SPeter Brune - zo - The offset in the z direction 248d886c4f4SPeter Brune 249d886c4f4SPeter Brune Level: intermediate 250d886c4f4SPeter Brune 251d886c4f4SPeter Brune Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without 252d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 253d886c4f4SPeter Brune 254d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 255d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 256d886c4f4SPeter Brune @*/ 25795c13181SPeter Brune PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 258d886c4f4SPeter Brune { 25995c13181SPeter Brune PetscErrorCode ierr; 260d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 261d886c4f4SPeter Brune 262d886c4f4SPeter Brune PetscFunctionBegin; 263d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 264d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da,xo,2); 26595c13181SPeter Brune PetscValidLogicalCollectiveInt(da,yo,3); 26695c13181SPeter Brune PetscValidLogicalCollectiveInt(da,zo,4); 26795c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Mo,5); 26895c13181SPeter Brune PetscValidLogicalCollectiveInt(da,No,6); 26995c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Po,7); 270d886c4f4SPeter Brune dd->xo = xo; 271d886c4f4SPeter Brune dd->yo = yo; 272d886c4f4SPeter Brune dd->zo = zo; 27395c13181SPeter Brune dd->Mo = Mo; 27495c13181SPeter Brune dd->No = No; 27595c13181SPeter Brune dd->Po = Po; 27695c13181SPeter Brune 27795c13181SPeter Brune if (da->coordinateDM) { 27895c13181SPeter Brune ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 27995c13181SPeter Brune } 280d886c4f4SPeter Brune PetscFunctionReturn(0); 281d886c4f4SPeter Brune } 282d886c4f4SPeter Brune 283d886c4f4SPeter Brune #undef __FUNCT__ 284d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset" 285d886c4f4SPeter Brune /*@ 286d886c4f4SPeter Brune DMDAGetOffset - Gets the index offset of the DA. 287d886c4f4SPeter Brune 288d886c4f4SPeter Brune Not collective 289d886c4f4SPeter Brune 290d886c4f4SPeter Brune Input Parameter: 291d886c4f4SPeter Brune . da - The DMDA 292d886c4f4SPeter Brune 293d886c4f4SPeter Brune Output Parameters: 294d886c4f4SPeter Brune + xo - The offset in the x direction 295d886c4f4SPeter Brune . yo - The offset in the y direction 29695c13181SPeter Brune . zo - The offset in the z direction 29795c13181SPeter Brune . Mo - The global size in the x direction 29895c13181SPeter Brune . No - The global size in the y direction 29995c13181SPeter Brune - Po - The global size in the z direction 300d886c4f4SPeter Brune 301d886c4f4SPeter Brune Level: intermediate 302d886c4f4SPeter Brune 303d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 304d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray() 305d886c4f4SPeter Brune @*/ 30695c13181SPeter Brune PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 307d886c4f4SPeter Brune { 308d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 309d886c4f4SPeter Brune 310d886c4f4SPeter Brune PetscFunctionBegin; 311d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 312d886c4f4SPeter Brune if (xo) *xo = dd->xo; 313d886c4f4SPeter Brune if (yo) *yo = dd->yo; 314d886c4f4SPeter Brune if (zo) *zo = dd->zo; 31595c13181SPeter Brune if (Mo) *Mo = dd->Mo; 31695c13181SPeter Brune if (No) *No = dd->No; 31795c13181SPeter Brune if (Po) *Po = dd->Po; 318d886c4f4SPeter Brune PetscFunctionReturn(0); 319d886c4f4SPeter Brune } 320d886c4f4SPeter Brune 32188661749SPeter Brune 32288661749SPeter Brune #undef __FUNCT__ 323aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType" 32447c6ae99SBarry Smith /*@ 325aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 32647c6ae99SBarry Smith 327aa219208SBarry Smith Logically Collective on DMDA 32847c6ae99SBarry Smith 32947c6ae99SBarry Smith Input Parameter: 330aa219208SBarry Smith + da - The DMDA 331aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 33247c6ae99SBarry Smith 33347c6ae99SBarry Smith Level: intermediate 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith .keywords: distributed array, stencil 33696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 33747c6ae99SBarry Smith @*/ 3387087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 33947c6ae99SBarry Smith { 34047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith PetscFunctionBegin; 34347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 34447c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 345*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 34647c6ae99SBarry Smith dd->stencil_type = stype; 34747c6ae99SBarry Smith PetscFunctionReturn(0); 34847c6ae99SBarry Smith } 34947c6ae99SBarry Smith 35047c6ae99SBarry Smith #undef __FUNCT__ 351aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth" 35247c6ae99SBarry Smith /*@ 353aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 35447c6ae99SBarry Smith 355aa219208SBarry Smith Logically Collective on DMDA 35647c6ae99SBarry Smith 35747c6ae99SBarry Smith Input Parameter: 358aa219208SBarry Smith + da - The DMDA 35947c6ae99SBarry Smith - width - The stencil width 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith Level: intermediate 36247c6ae99SBarry Smith 36347c6ae99SBarry Smith .keywords: distributed array, stencil 36496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 36547c6ae99SBarry Smith @*/ 3667087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 36747c6ae99SBarry Smith { 36847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 36947c6ae99SBarry Smith 37047c6ae99SBarry Smith PetscFunctionBegin; 37147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 37247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 373*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 37447c6ae99SBarry Smith dd->s = width; 37547c6ae99SBarry Smith PetscFunctionReturn(0); 37647c6ae99SBarry Smith } 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith #undef __FUNCT__ 379aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 380aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 38147c6ae99SBarry Smith { 38247c6ae99SBarry Smith PetscInt i,sum; 38347c6ae99SBarry Smith 38447c6ae99SBarry Smith PetscFunctionBegin; 385*ce94432eSBarry Smith if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 38647c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 387*ce94432eSBarry Smith if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 38847c6ae99SBarry Smith PetscFunctionReturn(0); 38947c6ae99SBarry Smith } 39047c6ae99SBarry Smith 39147c6ae99SBarry Smith #undef __FUNCT__ 392aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges" 39347c6ae99SBarry Smith /*@ 394aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 39547c6ae99SBarry Smith 396aa219208SBarry Smith Logically Collective on DMDA 39747c6ae99SBarry Smith 39847c6ae99SBarry Smith Input Parameter: 399aa219208SBarry Smith + da - The DMDA 4000298fd71SBarry Smith . lx - array containing number of nodes in the X direction on each process, or NULL. If non-null, must be of length da->m 4010298fd71SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or NULL. If non-null, must be of length da->n 4020298fd71SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or NULL. If non-null, must be of length da->p. 40347c6ae99SBarry Smith 40447c6ae99SBarry Smith Level: intermediate 40547c6ae99SBarry Smith 406e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 407e3f3e4b6SBarry Smith 40847c6ae99SBarry Smith .keywords: distributed array 40996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 41047c6ae99SBarry Smith @*/ 4117087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 41247c6ae99SBarry Smith { 41347c6ae99SBarry Smith PetscErrorCode ierr; 41447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 41547c6ae99SBarry Smith 41647c6ae99SBarry Smith PetscFunctionBegin; 41747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 418*ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 41947c6ae99SBarry Smith if (lx) { 420*ce94432eSBarry Smith if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 421aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 42247c6ae99SBarry Smith if (!dd->lx) { 42347c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 42447c6ae99SBarry Smith } 42547c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 42647c6ae99SBarry Smith } 42747c6ae99SBarry Smith if (ly) { 428*ce94432eSBarry Smith if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 429aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 43047c6ae99SBarry Smith if (!dd->ly) { 43147c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 43247c6ae99SBarry Smith } 43347c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 43447c6ae99SBarry Smith } 43547c6ae99SBarry Smith if (lz) { 436*ce94432eSBarry Smith if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 437aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 43847c6ae99SBarry Smith if (!dd->lz) { 43947c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 44047c6ae99SBarry Smith } 44147c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith PetscFunctionReturn(0); 44447c6ae99SBarry Smith } 44547c6ae99SBarry Smith 44647c6ae99SBarry Smith #undef __FUNCT__ 447aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType" 44847c6ae99SBarry Smith /*@ 449aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 450e727c939SJed Brown returned by DMCreateInterpolation() 45147c6ae99SBarry Smith 452aa219208SBarry Smith Logically Collective on DMDA 45347c6ae99SBarry Smith 45447c6ae99SBarry Smith Input Parameter: 45547c6ae99SBarry Smith + da - initial distributed array 456aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 45747c6ae99SBarry Smith 45847c6ae99SBarry Smith Level: intermediate 45947c6ae99SBarry Smith 460e727c939SJed Brown Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 46147c6ae99SBarry Smith 46247c6ae99SBarry Smith .keywords: distributed array, interpolation 46347c6ae99SBarry Smith 46496e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 46547c6ae99SBarry Smith @*/ 4667087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 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 PetscValidLogicalCollectiveEnum(da,ctype,2); 47347c6ae99SBarry Smith dd->interptype = ctype; 47447c6ae99SBarry Smith PetscFunctionReturn(0); 47547c6ae99SBarry Smith } 47647c6ae99SBarry Smith 4772dde6fd4SLisandro Dalcin #undef __FUNCT__ 4782dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType" 4792dde6fd4SLisandro Dalcin /*@ 4802dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 481e727c939SJed Brown used by DMCreateInterpolation() 4822dde6fd4SLisandro Dalcin 4832dde6fd4SLisandro Dalcin Not Collective 4842dde6fd4SLisandro Dalcin 4852dde6fd4SLisandro Dalcin Input Parameter: 4862dde6fd4SLisandro Dalcin . da - distributed array 4872dde6fd4SLisandro Dalcin 4882dde6fd4SLisandro Dalcin Output Parameter: 4892dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 4902dde6fd4SLisandro Dalcin 4912dde6fd4SLisandro Dalcin Level: intermediate 4922dde6fd4SLisandro Dalcin 4932dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 4942dde6fd4SLisandro Dalcin 495e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 4962dde6fd4SLisandro Dalcin @*/ 4972dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 4982dde6fd4SLisandro Dalcin { 4992dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 5002dde6fd4SLisandro Dalcin 5012dde6fd4SLisandro Dalcin PetscFunctionBegin; 5022dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 5032dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 5042dde6fd4SLisandro Dalcin *ctype = dd->interptype; 5052dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 5062dde6fd4SLisandro Dalcin } 50747c6ae99SBarry Smith 50847c6ae99SBarry Smith #undef __FUNCT__ 509aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors" 51047c6ae99SBarry Smith /*@C 511aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 51247c6ae99SBarry Smith processes neighbors. 51347c6ae99SBarry Smith 51447c6ae99SBarry Smith Not Collective 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith Input Parameter: 517aa219208SBarry Smith . da - the DMDA object 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith Output Parameters: 52047c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 52147c6ae99SBarry Smith this process itself is in the list 52247c6ae99SBarry Smith 52347c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 52447c6ae99SBarry Smith Not supported in 1d 525aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 52647c6ae99SBarry Smith 52747c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith Level: intermediate 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith @*/ 5327087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 53347c6ae99SBarry Smith { 53447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 5355fd66863SKarl Rupp 53647c6ae99SBarry Smith PetscFunctionBegin; 53747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 53847c6ae99SBarry Smith *ranks = dd->neighbors; 53947c6ae99SBarry Smith PetscFunctionReturn(0); 54047c6ae99SBarry Smith } 54147c6ae99SBarry Smith 54247c6ae99SBarry Smith #undef __FUNCT__ 543aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges" 54447c6ae99SBarry Smith /*@C 545aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith Not Collective 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith Input Parameter: 550aa219208SBarry Smith . da - the DMDA object 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith Output Parameter: 55347c6ae99SBarry Smith + lx - ownership along x direction (optional) 55447c6ae99SBarry Smith . ly - ownership along y direction (optional) 55547c6ae99SBarry Smith - lz - ownership along z direction (optional) 55647c6ae99SBarry Smith 55747c6ae99SBarry Smith Level: intermediate 55847c6ae99SBarry Smith 559aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 56047c6ae99SBarry Smith 56147c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 562aa219208SBarry Smith eighth arguments from DMDAGetInfo() 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 565aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 56647c6ae99SBarry Smith 567e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 568e3f3e4b6SBarry Smith 569aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 57047c6ae99SBarry Smith @*/ 5717087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 57247c6ae99SBarry Smith { 57347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 57447c6ae99SBarry Smith 57547c6ae99SBarry Smith PetscFunctionBegin; 57647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 57747c6ae99SBarry Smith if (lx) *lx = dd->lx; 57847c6ae99SBarry Smith if (ly) *ly = dd->ly; 57947c6ae99SBarry Smith if (lz) *lz = dd->lz; 58047c6ae99SBarry Smith PetscFunctionReturn(0); 58147c6ae99SBarry Smith } 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith #undef __FUNCT__ 584aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor" 58547c6ae99SBarry Smith /*@ 586aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 58747c6ae99SBarry Smith 588aa219208SBarry Smith Logically Collective on DMDA 58947c6ae99SBarry Smith 59047c6ae99SBarry Smith Input Parameters: 591aa219208SBarry Smith + da - the DMDA object 59247c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 59347c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 59447c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 59547c6ae99SBarry Smith 59647c6ae99SBarry Smith Options Database: 59747c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 59847c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 59947c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 60047c6ae99SBarry Smith 60147c6ae99SBarry Smith Level: intermediate 60247c6ae99SBarry Smith 60347c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 60447c6ae99SBarry Smith 605aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 60647c6ae99SBarry Smith @*/ 6077087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 60847c6ae99SBarry Smith { 60947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 61047c6ae99SBarry Smith 61147c6ae99SBarry Smith PetscFunctionBegin; 61247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 61347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 61447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 61547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 61647c6ae99SBarry Smith 61747c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 61847c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 61947c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 62047c6ae99SBarry Smith PetscFunctionReturn(0); 62147c6ae99SBarry Smith } 62247c6ae99SBarry Smith 62347c6ae99SBarry Smith #undef __FUNCT__ 624aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor" 62547c6ae99SBarry Smith /*@C 626aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 62747c6ae99SBarry Smith 62847c6ae99SBarry Smith Not Collective 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith Input Parameter: 631aa219208SBarry Smith . da - the DMDA object 63247c6ae99SBarry Smith 63347c6ae99SBarry Smith Output Parameters: 63447c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 63547c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 63647c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 63747c6ae99SBarry Smith 63847c6ae99SBarry Smith Level: intermediate 63947c6ae99SBarry Smith 6400298fd71SBarry Smith Notes: Pass NULL for values you do not need 64147c6ae99SBarry Smith 642aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 64347c6ae99SBarry Smith @*/ 6447087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 64547c6ae99SBarry Smith { 64647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith PetscFunctionBegin; 64947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 65047c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 65147c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 65247c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 65347c6ae99SBarry Smith PetscFunctionReturn(0); 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith 65647c6ae99SBarry Smith #undef __FUNCT__ 657aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix" 65847c6ae99SBarry Smith /*@C 659aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 66047c6ae99SBarry Smith 661aa219208SBarry Smith Logically Collective on DMDA 66247c6ae99SBarry Smith 66347c6ae99SBarry Smith Input Parameters: 664aa219208SBarry Smith + da - the DMDA object 665aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith Level: developer 66847c6ae99SBarry Smith 669aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 67047c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 67147c6ae99SBarry Smith 672950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 67347c6ae99SBarry Smith @*/ 67419fd82e9SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*)) 67547c6ae99SBarry Smith { 67647c6ae99SBarry Smith PetscFunctionBegin; 67747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 67825296bd5SBarry Smith da->ops->creatematrix = f; 67947c6ae99SBarry Smith PetscFunctionReturn(0); 68047c6ae99SBarry Smith } 68147c6ae99SBarry Smith 68247c6ae99SBarry Smith #undef __FUNCT__ 68394013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges" 68447c6ae99SBarry Smith /* 68547c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 68647c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 68747c6ae99SBarry Smith 68847c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 68947c6ae99SBarry Smith */ 69094013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 69147c6ae99SBarry Smith { 69247c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 69347c6ae99SBarry Smith PetscErrorCode ierr; 69447c6ae99SBarry Smith 69547c6ae99SBarry Smith PetscFunctionBegin; 696*ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 69747c6ae99SBarry Smith if (ratio == 1) { 69847c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 69947c6ae99SBarry Smith PetscFunctionReturn(0); 70047c6ae99SBarry Smith } 70147c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 70247c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 70347c6ae99SBarry Smith for (i=0; i<m; i++) { 70447c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 70547c6ae99SBarry Smith if (i == m-1) lf[i] = want; 70647c6ae99SBarry Smith else { 7077aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 7087aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 7097aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 7107aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 7117aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 7127aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 7137aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 7147aca7175SJed Brown /* Make sure all constraints are satisfied */ 71530729d88SBarry Smith if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 716*ce94432eSBarry Smith || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 71747c6ae99SBarry Smith } 71847c6ae99SBarry Smith lf[i] = want; 71947c6ae99SBarry Smith startc += lc[i]; 72047c6ae99SBarry Smith startf += lf[i]; 72147c6ae99SBarry Smith remaining -= lf[i]; 72247c6ae99SBarry Smith } 72347c6ae99SBarry Smith PetscFunctionReturn(0); 72447c6ae99SBarry Smith } 72547c6ae99SBarry Smith 72647c6ae99SBarry Smith #undef __FUNCT__ 7272be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges" 7282be375d4SJed Brown /* 7292be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 7302be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 7312be375d4SJed Brown 7322be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 7332be375d4SJed Brown */ 7342be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 7352be375d4SJed Brown { 7362be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 7372be375d4SJed Brown PetscErrorCode ierr; 7382be375d4SJed Brown 7392be375d4SJed Brown PetscFunctionBegin; 740*ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 7412be375d4SJed Brown if (ratio == 1) { 7422be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 7432be375d4SJed Brown PetscFunctionReturn(0); 7442be375d4SJed Brown } 7452be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 7462be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 7472be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 7482be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 7492be375d4SJed Brown if (i == m-1) lc[i] = want; 7502be375d4SJed Brown else { 7512be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 7522be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 7532be375d4SJed Brown * node is within one stencil width. */ 7542be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 7552be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 7562be375d4SJed Brown * fine node is within one stencil width. */ 7572be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 7582be375d4SJed Brown if (want < 0 || want > remaining 759*ce94432eSBarry Smith || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 7602be375d4SJed Brown } 7612be375d4SJed Brown lc[i] = want; 7622be375d4SJed Brown startc += lc[i]; 7632be375d4SJed Brown startf += lf[i]; 7642be375d4SJed Brown remaining -= lc[i]; 7652be375d4SJed Brown } 7662be375d4SJed Brown PetscFunctionReturn(0); 7672be375d4SJed Brown } 7682be375d4SJed Brown 7692be375d4SJed Brown #undef __FUNCT__ 77094013140SBarry Smith #define __FUNCT__ "DMRefine_DA" 7717087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 77247c6ae99SBarry Smith { 77347c6ae99SBarry Smith PetscErrorCode ierr; 774f3141302SJed Brown PetscInt M,N,P,i; 7759a42bb27SBarry Smith DM da2; 77647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 77747c6ae99SBarry Smith 77847c6ae99SBarry Smith PetscFunctionBegin; 77947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 78047c6ae99SBarry Smith PetscValidPointer(daref,3); 78147c6ae99SBarry Smith 7821321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 78347c6ae99SBarry Smith M = dd->refine_x*dd->M; 78447c6ae99SBarry Smith } else { 78547c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 78647c6ae99SBarry Smith } 7871321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 7881860e6e9SBarry Smith if (dd->dim > 1) { 78947c6ae99SBarry Smith N = dd->refine_y*dd->N; 79047c6ae99SBarry Smith } else { 7911860e6e9SBarry Smith N = 1; 7921860e6e9SBarry Smith } 7931860e6e9SBarry Smith } else { 79447c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 79547c6ae99SBarry Smith } 7961321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 7971860e6e9SBarry Smith if (dd->dim > 2) { 79847c6ae99SBarry Smith P = dd->refine_z*dd->P; 79947c6ae99SBarry Smith } else { 8001860e6e9SBarry Smith P = 1; 8011860e6e9SBarry Smith } 8021860e6e9SBarry Smith } else { 80347c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 80447c6ae99SBarry Smith } 805*ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 806397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 807397b6216SJed Brown ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr); 808397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 809397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 810397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 811397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 812397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 813397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 81447c6ae99SBarry Smith if (dd->dim == 3) { 81547c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 81647c6ae99SBarry Smith ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 8171321219cSEthan 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); 8181321219cSEthan 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); 8191321219cSEthan 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); 820397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 82147c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 82247c6ae99SBarry Smith } else if (dd->dim == 2) { 82347c6ae99SBarry Smith PetscInt *lx,*ly; 82447c6ae99SBarry Smith ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 8251321219cSEthan 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); 8261321219cSEthan 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); 8270298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 82847c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 82947c6ae99SBarry Smith } else if (dd->dim == 1) { 83047c6ae99SBarry Smith PetscInt *lx; 83147c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 8321321219cSEthan 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); 8330298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 83447c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 83747c6ae99SBarry Smith 83847c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 83925296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 84025296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 84147c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 84247c6ae99SBarry Smith dd2->interptype = dd->interptype; 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith /* copy fill information if given */ 84547c6ae99SBarry Smith if (dd->dfill) { 84647c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 84747c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 84847c6ae99SBarry Smith } 84947c6ae99SBarry Smith if (dd->ofill) { 85047c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 85147c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith /* copy the refine information */ 854397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 855397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 856397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 85747c6ae99SBarry Smith 85847c6ae99SBarry Smith /* copy vector type information */ 85947c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 86019fd82e9SBarry Smith ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr); 861ddcf8b74SDave May 862efd51863SBarry Smith dd2->lf = dd->lf; 863efd51863SBarry Smith dd2->lj = dd->lj; 864efd51863SBarry Smith 8656e87535bSJed Brown da2->leveldown = da->leveldown; 8666e87535bSJed Brown da2->levelup = da->levelup + 1; 8678865f1eaSKarl Rupp 8686e87535bSJed Brown ierr = DMSetFromOptions(da2);CHKERRQ(ierr); 8696e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 870ca266f36SBarry Smith ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr); 8716e87535bSJed Brown 872ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 8736636e97aSMatthew G Knepley if (da->coordinates) { 874ddcf8b74SDave May DM cdaf,cdac; 875ddcf8b74SDave May Vec coordsc,coordsf; 876ddcf8b74SDave May Mat II; 877ddcf8b74SDave May 8786636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 8796636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 8806636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 881b61d3410SDave May /* force creation of the coordinate vector */ 882b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 8836636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 8840298fd71SBarry Smith ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 885ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 886fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 887ddcf8b74SDave May } 888397b6216SJed Brown 889f3141302SJed Brown for (i=0; i<da->bs; i++) { 890f3141302SJed Brown const char *fieldname; 891f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 892f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 893f3141302SJed Brown } 894397b6216SJed Brown 89547c6ae99SBarry Smith *daref = da2; 89647c6ae99SBarry Smith PetscFunctionReturn(0); 89747c6ae99SBarry Smith } 89847c6ae99SBarry Smith 899397b6216SJed Brown 90047c6ae99SBarry Smith #undef __FUNCT__ 90194013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA" 9027087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 90347c6ae99SBarry Smith { 90447c6ae99SBarry Smith PetscErrorCode ierr; 905397b6216SJed Brown PetscInt M,N,P,i; 9069a42bb27SBarry Smith DM da2; 90747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 90847c6ae99SBarry Smith 90947c6ae99SBarry Smith PetscFunctionBegin; 91047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 91147c6ae99SBarry Smith PetscValidPointer(daref,3); 91247c6ae99SBarry Smith 9131321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 914397b6216SJed Brown M = dd->M / dd->coarsen_x; 91547c6ae99SBarry Smith } else { 916397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 91747c6ae99SBarry Smith } 9181321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 9191860e6e9SBarry Smith if (dd->dim > 1) { 920397b6216SJed Brown N = dd->N / dd->coarsen_y; 92147c6ae99SBarry Smith } else { 9221860e6e9SBarry Smith N = 1; 9231860e6e9SBarry Smith } 9241860e6e9SBarry Smith } else { 925397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 92647c6ae99SBarry Smith } 9271321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 9281860e6e9SBarry Smith if (dd->dim > 2) { 929397b6216SJed Brown P = dd->P / dd->coarsen_z; 93047c6ae99SBarry Smith } else { 9311860e6e9SBarry Smith P = 1; 9321860e6e9SBarry Smith } 9331860e6e9SBarry Smith } else { 934397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 93547c6ae99SBarry Smith } 936*ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 937397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 938397b6216SJed Brown ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr); 939397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 940397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 941397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 942397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 943397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 944397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 94547c6ae99SBarry Smith if (dd->dim == 3) { 9462be375d4SJed Brown PetscInt *lx,*ly,*lz; 9472be375d4SJed Brown ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 948397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 949397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 950397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 951397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 9522be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 95347c6ae99SBarry Smith } else if (dd->dim == 2) { 9542be375d4SJed Brown PetscInt *lx,*ly; 9552be375d4SJed Brown ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 956397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 957397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 9580298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 9592be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 96047c6ae99SBarry Smith } else if (dd->dim == 1) { 9612be375d4SJed Brown PetscInt *lx; 9622be375d4SJed Brown ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 963397b6216SJed Brown ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 9640298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 9652be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 96647c6ae99SBarry Smith } 96747c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 96847c6ae99SBarry Smith 9694dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 97025296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 97125296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 97247c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 97347c6ae99SBarry Smith dd2->interptype = dd->interptype; 97447c6ae99SBarry Smith 97547c6ae99SBarry Smith /* copy fill information if given */ 97647c6ae99SBarry Smith if (dd->dfill) { 97747c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 97847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 97947c6ae99SBarry Smith } 98047c6ae99SBarry Smith if (dd->ofill) { 98147c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 98247c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 98347c6ae99SBarry Smith } 98447c6ae99SBarry Smith /* copy the refine information */ 985397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 986397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 987397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 98847c6ae99SBarry Smith 98947c6ae99SBarry Smith /* copy vector type information */ 99047c6ae99SBarry Smith ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 99119fd82e9SBarry Smith ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr); 99247c6ae99SBarry Smith 993644e2e5bSBarry Smith dd2->lf = dd->lf; 994644e2e5bSBarry Smith dd2->lj = dd->lj; 995644e2e5bSBarry Smith 9966e87535bSJed Brown da2->leveldown = da->leveldown + 1; 9976e87535bSJed Brown da2->levelup = da->levelup; 9988865f1eaSKarl Rupp 9996e87535bSJed Brown ierr = DMSetFromOptions(da2);CHKERRQ(ierr); 10006e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 1001ca266f36SBarry Smith ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr); 10026e87535bSJed Brown 1003ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 10046636e97aSMatthew G Knepley if (da->coordinates) { 1005ddcf8b74SDave May DM cdaf,cdac; 1006ddcf8b74SDave May Vec coordsc,coordsf; 1007ddcf8b74SDave May VecScatter inject; 1008ddcf8b74SDave May 10096636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); 10106636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); 10116636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); 1012b61d3410SDave May /* force creation of the coordinate vector */ 1013b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 10146636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 1015ddcf8b74SDave May 1016e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 1017ddcf8b74SDave May ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1018ddcf8b74SDave May ierr = VecScatterEnd(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1019fcfd50ebSBarry Smith ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 1020ddcf8b74SDave May } 1021f98405f7SJed Brown 1022f98405f7SJed Brown for (i=0; i<da->bs; i++) { 1023f98405f7SJed Brown const char *fieldname; 1024f98405f7SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1025f98405f7SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1026f98405f7SJed Brown } 1027f98405f7SJed Brown 102847c6ae99SBarry Smith *daref = da2; 102947c6ae99SBarry Smith PetscFunctionReturn(0); 103047c6ae99SBarry Smith } 103147c6ae99SBarry Smith 103247c6ae99SBarry Smith #undef __FUNCT__ 103394013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA" 10347087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 103547c6ae99SBarry Smith { 103647c6ae99SBarry Smith PetscErrorCode ierr; 103747c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 103847c6ae99SBarry Smith 103947c6ae99SBarry Smith PetscFunctionBegin; 104047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1041*ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 104247c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 104347c6ae99SBarry Smith PetscValidPointer(daf,3); 104447c6ae99SBarry Smith 1045aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 104647c6ae99SBarry Smith ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr); 104747c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 1048aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 104947c6ae99SBarry Smith } 105047c6ae99SBarry Smith n = nlevels; 10510298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 105247c6ae99SBarry Smith n = nlevels; 10530298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 105447c6ae99SBarry Smith n = nlevels; 10550298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 105647c6ae99SBarry Smith 1057aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1058*ce94432eSBarry Smith ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 105947c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1060aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1061*ce94432eSBarry Smith ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 106447c6ae99SBarry Smith PetscFunctionReturn(0); 106547c6ae99SBarry Smith } 106647c6ae99SBarry Smith 106747c6ae99SBarry Smith #undef __FUNCT__ 106894013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA" 10697087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 107047c6ae99SBarry Smith { 107147c6ae99SBarry Smith PetscErrorCode ierr; 107247c6ae99SBarry Smith PetscInt i; 107347c6ae99SBarry Smith 107447c6ae99SBarry Smith PetscFunctionBegin; 107547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1076*ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 107747c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 107847c6ae99SBarry Smith PetscValidPointer(dac,3); 1079*ce94432eSBarry Smith ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 108047c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1081*ce94432eSBarry Smith ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 108247c6ae99SBarry Smith } 108347c6ae99SBarry Smith PetscFunctionReturn(0); 108447c6ae99SBarry Smith } 1085