14035e84dSBarry Smith #include <petsc-private/dmdaimpl.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); 25ce94432eSBarry 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); 56ce94432eSBarry 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); 91ce94432eSBarry 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; 97ce94432eSBarry 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; 100ce94432eSBarry 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; 104ce94432eSBarry 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); 135ce94432eSBarry 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); 165ce94432eSBarry 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 2363e7870d2SPeter Brune 2373e7870d2SPeter Brune #undef __FUNCT__ 2383e7870d2SPeter Brune #define __FUNCT__ "DMDAGetNumLocalSubDomains" 2393e7870d2SPeter Brune /*@ 2403e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2413e7870d2SPeter Brune 2423e7870d2SPeter Brune Not collective 2433e7870d2SPeter Brune 2443e7870d2SPeter Brune Input Parameters: 2453e7870d2SPeter Brune . da - The DMDA 2463e7870d2SPeter Brune 2473e7870d2SPeter Brune Output Parameters: 2483e7870d2SPeter Brune + Nsub - Number of local subdomains created upon decomposition 2493e7870d2SPeter Brune 2503e7870d2SPeter Brune Level: intermediate 2513e7870d2SPeter Brune 2523e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2533e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA 2543e7870d2SPeter Brune @*/ 2553e7870d2SPeter Brune PetscErrorCode DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub) 2563e7870d2SPeter Brune { 2573e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2583e7870d2SPeter Brune 2593e7870d2SPeter Brune PetscFunctionBegin; 2603e7870d2SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2613e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2623e7870d2SPeter Brune PetscFunctionReturn(0); 2633e7870d2SPeter Brune } 2643e7870d2SPeter Brune 2653e7870d2SPeter Brune #undef __FUNCT__ 2663e7870d2SPeter Brune #define __FUNCT__ "DMDASetNumLocalSubDomains" 2673e7870d2SPeter Brune /*@ 2683e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2693e7870d2SPeter Brune 2703e7870d2SPeter Brune Not collective 2713e7870d2SPeter Brune 2723e7870d2SPeter Brune Input Parameters: 2733e7870d2SPeter Brune + da - The DMDA 2743e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2753e7870d2SPeter Brune 2763e7870d2SPeter Brune Level: intermediate 2773e7870d2SPeter Brune 2783e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2793e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA 2803e7870d2SPeter Brune @*/ 2813e7870d2SPeter Brune PetscErrorCode DMDASetNumLocalSubDomains(DM da,PetscInt Nsub) 2823e7870d2SPeter Brune { 2833e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2843e7870d2SPeter Brune 2853e7870d2SPeter Brune PetscFunctionBegin; 2863e7870d2SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2873e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da,Nsub,2); 2883e7870d2SPeter Brune dd->Nsub = Nsub; 2893e7870d2SPeter Brune PetscFunctionReturn(0); 2903e7870d2SPeter Brune } 2913e7870d2SPeter Brune 292d886c4f4SPeter Brune #undef __FUNCT__ 293d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset" 294d886c4f4SPeter Brune /*@ 295d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 296d886c4f4SPeter Brune 297d886c4f4SPeter Brune Collective on DA 298d886c4f4SPeter Brune 299d886c4f4SPeter Brune Input Parameter: 300d886c4f4SPeter Brune + da - The DMDA 301d886c4f4SPeter Brune . xo - The offset in the x direction 302d886c4f4SPeter Brune . yo - The offset in the y direction 303d886c4f4SPeter Brune - zo - The offset in the z direction 304d886c4f4SPeter Brune 305d886c4f4SPeter Brune Level: intermediate 306d886c4f4SPeter Brune 307d886c4f4SPeter Brune Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without 308d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 309d886c4f4SPeter Brune 310d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 311d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 312d886c4f4SPeter Brune @*/ 31395c13181SPeter Brune PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 314d886c4f4SPeter Brune { 31595c13181SPeter Brune PetscErrorCode ierr; 316d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 317d886c4f4SPeter Brune 318d886c4f4SPeter Brune PetscFunctionBegin; 319d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 320d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da,xo,2); 32195c13181SPeter Brune PetscValidLogicalCollectiveInt(da,yo,3); 32295c13181SPeter Brune PetscValidLogicalCollectiveInt(da,zo,4); 32395c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Mo,5); 32495c13181SPeter Brune PetscValidLogicalCollectiveInt(da,No,6); 32595c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Po,7); 326d886c4f4SPeter Brune dd->xo = xo; 327d886c4f4SPeter Brune dd->yo = yo; 328d886c4f4SPeter Brune dd->zo = zo; 32995c13181SPeter Brune dd->Mo = Mo; 33095c13181SPeter Brune dd->No = No; 33195c13181SPeter Brune dd->Po = Po; 33295c13181SPeter Brune 33395c13181SPeter Brune if (da->coordinateDM) { 33495c13181SPeter Brune ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 33595c13181SPeter Brune } 336d886c4f4SPeter Brune PetscFunctionReturn(0); 337d886c4f4SPeter Brune } 338d886c4f4SPeter Brune 339d886c4f4SPeter Brune #undef __FUNCT__ 340d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset" 341d886c4f4SPeter Brune /*@ 342d886c4f4SPeter Brune DMDAGetOffset - Gets the index offset of the DA. 343d886c4f4SPeter Brune 344d886c4f4SPeter Brune Not collective 345d886c4f4SPeter Brune 346d886c4f4SPeter Brune Input Parameter: 347d886c4f4SPeter Brune . da - The DMDA 348d886c4f4SPeter Brune 349d886c4f4SPeter Brune Output Parameters: 350d886c4f4SPeter Brune + xo - The offset in the x direction 351d886c4f4SPeter Brune . yo - The offset in the y direction 35295c13181SPeter Brune . zo - The offset in the z direction 35395c13181SPeter Brune . Mo - The global size in the x direction 35495c13181SPeter Brune . No - The global size in the y direction 35595c13181SPeter Brune - Po - The global size in the z direction 356d886c4f4SPeter Brune 357d886c4f4SPeter Brune Level: intermediate 358d886c4f4SPeter Brune 359d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 360d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray() 361d886c4f4SPeter Brune @*/ 36295c13181SPeter Brune PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 363d886c4f4SPeter Brune { 364d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 365d886c4f4SPeter Brune 366d886c4f4SPeter Brune PetscFunctionBegin; 367d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 368d886c4f4SPeter Brune if (xo) *xo = dd->xo; 369d886c4f4SPeter Brune if (yo) *yo = dd->yo; 370d886c4f4SPeter Brune if (zo) *zo = dd->zo; 37195c13181SPeter Brune if (Mo) *Mo = dd->Mo; 37295c13181SPeter Brune if (No) *No = dd->No; 37395c13181SPeter Brune if (Po) *Po = dd->Po; 374d886c4f4SPeter Brune PetscFunctionReturn(0); 375d886c4f4SPeter Brune } 376d886c4f4SPeter Brune 37740234c92SPeter Brune #undef __FUNCT__ 37840234c92SPeter Brune #define __FUNCT__ "DMDAGetNonOverlappingRegion" 37940234c92SPeter Brune /*@ 38040234c92SPeter Brune DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM. 38140234c92SPeter Brune 38240234c92SPeter Brune Not collective 38340234c92SPeter Brune 38440234c92SPeter Brune Input Parameter: 38540234c92SPeter Brune . da - The DMDA 38640234c92SPeter Brune 38740234c92SPeter Brune Output Parameters: 38840234c92SPeter Brune + xs - The start of the region in x 38940234c92SPeter Brune . ys - The start of the region in y 39040234c92SPeter Brune . zs - The start of the region in z 39140234c92SPeter Brune . xs - The size of the region in x 39240234c92SPeter Brune . ys - The size of the region in y 39340234c92SPeter Brune . zs - The size of the region in z 39440234c92SPeter Brune 39540234c92SPeter Brune Level: intermediate 39640234c92SPeter Brune 39740234c92SPeter Brune .keywords: distributed array, degrees of freedom 39840234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 39940234c92SPeter Brune @*/ 40040234c92SPeter Brune PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 40140234c92SPeter Brune { 40240234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 40340234c92SPeter Brune 40440234c92SPeter Brune PetscFunctionBegin; 40540234c92SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 40640234c92SPeter Brune if (xs) *xs = dd->nonxs; 40740234c92SPeter Brune if (ys) *ys = dd->nonys; 40840234c92SPeter Brune if (zs) *zs = dd->nonzs; 40940234c92SPeter Brune if (xm) *xm = dd->nonxm; 41040234c92SPeter Brune if (ym) *ym = dd->nonym; 41140234c92SPeter Brune if (zm) *zm = dd->nonzm; 41240234c92SPeter Brune PetscFunctionReturn(0); 41340234c92SPeter Brune } 41440234c92SPeter Brune 41540234c92SPeter Brune 41640234c92SPeter Brune #undef __FUNCT__ 41740234c92SPeter Brune #define __FUNCT__ "DMDASetNonOverlappingRegion" 41840234c92SPeter Brune /*@ 41940234c92SPeter Brune DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM. 42040234c92SPeter Brune 42140234c92SPeter Brune Collective on DA 42240234c92SPeter Brune 42340234c92SPeter Brune Input Parameter: 42440234c92SPeter Brune + da - The DMDA 42540234c92SPeter Brune . xs - The start of the region in x 42640234c92SPeter Brune . ys - The start of the region in y 42740234c92SPeter Brune . zs - The start of the region in z 42840234c92SPeter Brune . xs - The size of the region in x 42940234c92SPeter Brune . ys - The size of the region in y 43040234c92SPeter Brune . zs - The size of the region in z 43140234c92SPeter Brune 43240234c92SPeter Brune Level: intermediate 43340234c92SPeter Brune 43440234c92SPeter Brune .keywords: distributed array, degrees of freedom 43540234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 43640234c92SPeter Brune @*/ 43740234c92SPeter Brune PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 43840234c92SPeter Brune { 43940234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 44040234c92SPeter Brune 44140234c92SPeter Brune PetscFunctionBegin; 44240234c92SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 44340234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xs,2); 44440234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ys,3); 44540234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zs,4); 44640234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xm,5); 44740234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ym,6); 44840234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zm,7); 44940234c92SPeter Brune dd->nonxs = xs; 45040234c92SPeter Brune dd->nonys = ys; 45140234c92SPeter Brune dd->nonzs = zs; 45240234c92SPeter Brune dd->nonxm = xm; 45340234c92SPeter Brune dd->nonym = ym; 45440234c92SPeter Brune dd->nonzm = zm; 45540234c92SPeter Brune 45640234c92SPeter Brune PetscFunctionReturn(0); 45740234c92SPeter Brune } 45888661749SPeter Brune 45988661749SPeter Brune #undef __FUNCT__ 460aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType" 46147c6ae99SBarry Smith /*@ 462aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 46347c6ae99SBarry Smith 464aa219208SBarry Smith Logically Collective on DMDA 46547c6ae99SBarry Smith 46647c6ae99SBarry Smith Input Parameter: 467aa219208SBarry Smith + da - The DMDA 468aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 46947c6ae99SBarry Smith 47047c6ae99SBarry Smith Level: intermediate 47147c6ae99SBarry Smith 47247c6ae99SBarry Smith .keywords: distributed array, stencil 47396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 47447c6ae99SBarry Smith @*/ 4757087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 47647c6ae99SBarry Smith { 47747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 47847c6ae99SBarry Smith 47947c6ae99SBarry Smith PetscFunctionBegin; 48047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 48147c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 482ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 48347c6ae99SBarry Smith dd->stencil_type = stype; 48447c6ae99SBarry Smith PetscFunctionReturn(0); 48547c6ae99SBarry Smith } 48647c6ae99SBarry Smith 48747c6ae99SBarry Smith #undef __FUNCT__ 488aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth" 48947c6ae99SBarry Smith /*@ 490aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 49147c6ae99SBarry Smith 492aa219208SBarry Smith Logically Collective on DMDA 49347c6ae99SBarry Smith 49447c6ae99SBarry Smith Input Parameter: 495aa219208SBarry Smith + da - The DMDA 49647c6ae99SBarry Smith - width - The stencil width 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith Level: intermediate 49947c6ae99SBarry Smith 50047c6ae99SBarry Smith .keywords: distributed array, stencil 50196e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 50247c6ae99SBarry Smith @*/ 5037087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 50447c6ae99SBarry Smith { 50547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 50647c6ae99SBarry Smith 50747c6ae99SBarry Smith PetscFunctionBegin; 50847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 50947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 510ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 51147c6ae99SBarry Smith dd->s = width; 51247c6ae99SBarry Smith PetscFunctionReturn(0); 51347c6ae99SBarry Smith } 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith #undef __FUNCT__ 516aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 517aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 51847c6ae99SBarry Smith { 51947c6ae99SBarry Smith PetscInt i,sum; 52047c6ae99SBarry Smith 52147c6ae99SBarry Smith PetscFunctionBegin; 522ce94432eSBarry Smith if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 52347c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 524ce94432eSBarry Smith if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 52547c6ae99SBarry Smith PetscFunctionReturn(0); 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith 52847c6ae99SBarry Smith #undef __FUNCT__ 529aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges" 53047c6ae99SBarry Smith /*@ 531aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 53247c6ae99SBarry Smith 533aa219208SBarry Smith Logically Collective on DMDA 53447c6ae99SBarry Smith 53547c6ae99SBarry Smith Input Parameter: 536aa219208SBarry Smith + da - The DMDA 5370298fd71SBarry 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 5380298fd71SBarry 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 5390298fd71SBarry 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. 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith Level: intermediate 54247c6ae99SBarry Smith 543e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 544e3f3e4b6SBarry Smith 54547c6ae99SBarry Smith .keywords: distributed array 54696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 54747c6ae99SBarry Smith @*/ 5487087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 54947c6ae99SBarry Smith { 55047c6ae99SBarry Smith PetscErrorCode ierr; 55147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith PetscFunctionBegin; 55447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 555ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 55647c6ae99SBarry Smith if (lx) { 557ce94432eSBarry Smith if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 558aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 55947c6ae99SBarry Smith if (!dd->lx) { 56047c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 56147c6ae99SBarry Smith } 56247c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 56347c6ae99SBarry Smith } 56447c6ae99SBarry Smith if (ly) { 565ce94432eSBarry Smith if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 566aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 56747c6ae99SBarry Smith if (!dd->ly) { 56847c6ae99SBarry Smith ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 57147c6ae99SBarry Smith } 57247c6ae99SBarry Smith if (lz) { 573ce94432eSBarry Smith if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 574aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 57547c6ae99SBarry Smith if (!dd->lz) { 57647c6ae99SBarry Smith ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 57947c6ae99SBarry Smith } 58047c6ae99SBarry Smith PetscFunctionReturn(0); 58147c6ae99SBarry Smith } 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith #undef __FUNCT__ 584aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType" 58547c6ae99SBarry Smith /*@ 586aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 587e727c939SJed Brown returned by DMCreateInterpolation() 58847c6ae99SBarry Smith 589aa219208SBarry Smith Logically Collective on DMDA 59047c6ae99SBarry Smith 59147c6ae99SBarry Smith Input Parameter: 59247c6ae99SBarry Smith + da - initial distributed array 593aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 59447c6ae99SBarry Smith 59547c6ae99SBarry Smith Level: intermediate 59647c6ae99SBarry Smith 597e727c939SJed Brown Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 59847c6ae99SBarry Smith 59947c6ae99SBarry Smith .keywords: distributed array, interpolation 60047c6ae99SBarry Smith 60196e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 60247c6ae99SBarry Smith @*/ 6037087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 60447c6ae99SBarry Smith { 60547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 60647c6ae99SBarry Smith 60747c6ae99SBarry Smith PetscFunctionBegin; 60847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 60947c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 61047c6ae99SBarry Smith dd->interptype = ctype; 61147c6ae99SBarry Smith PetscFunctionReturn(0); 61247c6ae99SBarry Smith } 61347c6ae99SBarry Smith 6142dde6fd4SLisandro Dalcin #undef __FUNCT__ 6152dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType" 6162dde6fd4SLisandro Dalcin /*@ 6172dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 618e727c939SJed Brown used by DMCreateInterpolation() 6192dde6fd4SLisandro Dalcin 6202dde6fd4SLisandro Dalcin Not Collective 6212dde6fd4SLisandro Dalcin 6222dde6fd4SLisandro Dalcin Input Parameter: 6232dde6fd4SLisandro Dalcin . da - distributed array 6242dde6fd4SLisandro Dalcin 6252dde6fd4SLisandro Dalcin Output Parameter: 6262dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 6272dde6fd4SLisandro Dalcin 6282dde6fd4SLisandro Dalcin Level: intermediate 6292dde6fd4SLisandro Dalcin 6302dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 6312dde6fd4SLisandro Dalcin 632e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 6332dde6fd4SLisandro Dalcin @*/ 6342dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 6352dde6fd4SLisandro Dalcin { 6362dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 6372dde6fd4SLisandro Dalcin 6382dde6fd4SLisandro Dalcin PetscFunctionBegin; 6392dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 6402dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 6412dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6422dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 6432dde6fd4SLisandro Dalcin } 64447c6ae99SBarry Smith 64547c6ae99SBarry Smith #undef __FUNCT__ 646aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors" 64747c6ae99SBarry Smith /*@C 648aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 64947c6ae99SBarry Smith processes neighbors. 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith Not Collective 65247c6ae99SBarry Smith 65347c6ae99SBarry Smith Input Parameter: 654aa219208SBarry Smith . da - the DMDA object 65547c6ae99SBarry Smith 65647c6ae99SBarry Smith Output Parameters: 65747c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 65847c6ae99SBarry Smith this process itself is in the list 65947c6ae99SBarry Smith 66047c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 66147c6ae99SBarry Smith Not supported in 1d 662aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 66347c6ae99SBarry Smith 66447c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 66547c6ae99SBarry Smith 66647c6ae99SBarry Smith Level: intermediate 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith @*/ 6697087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 67047c6ae99SBarry Smith { 67147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 6725fd66863SKarl Rupp 67347c6ae99SBarry Smith PetscFunctionBegin; 67447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 67547c6ae99SBarry Smith *ranks = dd->neighbors; 67647c6ae99SBarry Smith PetscFunctionReturn(0); 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith 67947c6ae99SBarry Smith #undef __FUNCT__ 680aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges" 68147c6ae99SBarry Smith /*@C 682aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith Not Collective 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith Input Parameter: 687aa219208SBarry Smith . da - the DMDA object 68847c6ae99SBarry Smith 68947c6ae99SBarry Smith Output Parameter: 69047c6ae99SBarry Smith + lx - ownership along x direction (optional) 69147c6ae99SBarry Smith . ly - ownership along y direction (optional) 69247c6ae99SBarry Smith - lz - ownership along z direction (optional) 69347c6ae99SBarry Smith 69447c6ae99SBarry Smith Level: intermediate 69547c6ae99SBarry Smith 696aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 69747c6ae99SBarry Smith 69847c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 699aa219208SBarry Smith eighth arguments from DMDAGetInfo() 70047c6ae99SBarry Smith 70147c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 702aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 70347c6ae99SBarry Smith 704e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 705e3f3e4b6SBarry Smith 706aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 70747c6ae99SBarry Smith @*/ 7087087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 70947c6ae99SBarry Smith { 71047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 71147c6ae99SBarry Smith 71247c6ae99SBarry Smith PetscFunctionBegin; 71347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 71447c6ae99SBarry Smith if (lx) *lx = dd->lx; 71547c6ae99SBarry Smith if (ly) *ly = dd->ly; 71647c6ae99SBarry Smith if (lz) *lz = dd->lz; 71747c6ae99SBarry Smith PetscFunctionReturn(0); 71847c6ae99SBarry Smith } 71947c6ae99SBarry Smith 72047c6ae99SBarry Smith #undef __FUNCT__ 721aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor" 72247c6ae99SBarry Smith /*@ 723aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 72447c6ae99SBarry Smith 725aa219208SBarry Smith Logically Collective on DMDA 72647c6ae99SBarry Smith 72747c6ae99SBarry Smith Input Parameters: 728aa219208SBarry Smith + da - the DMDA object 72947c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 73047c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 73147c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith Options Database: 73447c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 73547c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 73647c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 73747c6ae99SBarry Smith 73847c6ae99SBarry Smith Level: intermediate 73947c6ae99SBarry Smith 74047c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 74147c6ae99SBarry Smith 742aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 74347c6ae99SBarry Smith @*/ 7447087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 74547c6ae99SBarry Smith { 74647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 74747c6ae99SBarry Smith 74847c6ae99SBarry Smith PetscFunctionBegin; 74947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 75047c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 75147c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 75247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 75347c6ae99SBarry Smith 75447c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 75547c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 75647c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 75747c6ae99SBarry Smith PetscFunctionReturn(0); 75847c6ae99SBarry Smith } 75947c6ae99SBarry Smith 76047c6ae99SBarry Smith #undef __FUNCT__ 761aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor" 76247c6ae99SBarry Smith /*@C 763aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 76447c6ae99SBarry Smith 76547c6ae99SBarry Smith Not Collective 76647c6ae99SBarry Smith 76747c6ae99SBarry Smith Input Parameter: 768aa219208SBarry Smith . da - the DMDA object 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith Output Parameters: 77147c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 77247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 77347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith Level: intermediate 77647c6ae99SBarry Smith 7770298fd71SBarry Smith Notes: Pass NULL for values you do not need 77847c6ae99SBarry Smith 779aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 78047c6ae99SBarry Smith @*/ 7817087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 78247c6ae99SBarry Smith { 78347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith PetscFunctionBegin; 78647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 78747c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 78847c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 78947c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 79047c6ae99SBarry Smith PetscFunctionReturn(0); 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith 79347c6ae99SBarry Smith #undef __FUNCT__ 794aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix" 79547c6ae99SBarry Smith /*@C 796aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 79747c6ae99SBarry Smith 798aa219208SBarry Smith Logically Collective on DMDA 79947c6ae99SBarry Smith 80047c6ae99SBarry Smith Input Parameters: 801aa219208SBarry Smith + da - the DMDA object 802aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 80347c6ae99SBarry Smith 80447c6ae99SBarry Smith Level: developer 80547c6ae99SBarry Smith 806aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 80747c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 80847c6ae99SBarry Smith 809950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 81047c6ae99SBarry Smith @*/ 811b412c318SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*)) 81247c6ae99SBarry Smith { 81347c6ae99SBarry Smith PetscFunctionBegin; 81447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 81525296bd5SBarry Smith da->ops->creatematrix = f; 81647c6ae99SBarry Smith PetscFunctionReturn(0); 81747c6ae99SBarry Smith } 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith #undef __FUNCT__ 82094013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges" 82147c6ae99SBarry Smith /* 82247c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 82347c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 82447c6ae99SBarry Smith 82547c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 82647c6ae99SBarry Smith */ 82794013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 82847c6ae99SBarry Smith { 82947c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 83047c6ae99SBarry Smith PetscErrorCode ierr; 83147c6ae99SBarry Smith 83247c6ae99SBarry Smith PetscFunctionBegin; 833ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 83447c6ae99SBarry Smith if (ratio == 1) { 83547c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 83647c6ae99SBarry Smith PetscFunctionReturn(0); 83747c6ae99SBarry Smith } 83847c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 83947c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 84047c6ae99SBarry Smith for (i=0; i<m; i++) { 84147c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 84247c6ae99SBarry Smith if (i == m-1) lf[i] = want; 84347c6ae99SBarry Smith else { 8447aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 8457aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 8467aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 8477aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 8487aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 8497aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 8507aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 8517aca7175SJed Brown /* Make sure all constraints are satisfied */ 85230729d88SBarry Smith if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 853ce94432eSBarry 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"); 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith lf[i] = want; 85647c6ae99SBarry Smith startc += lc[i]; 85747c6ae99SBarry Smith startf += lf[i]; 85847c6ae99SBarry Smith remaining -= lf[i]; 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith PetscFunctionReturn(0); 86147c6ae99SBarry Smith } 86247c6ae99SBarry Smith 86347c6ae99SBarry Smith #undef __FUNCT__ 8642be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges" 8652be375d4SJed Brown /* 8662be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 8672be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 8682be375d4SJed Brown 8692be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 8702be375d4SJed Brown */ 8712be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 8722be375d4SJed Brown { 8732be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 8742be375d4SJed Brown PetscErrorCode ierr; 8752be375d4SJed Brown 8762be375d4SJed Brown PetscFunctionBegin; 877ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 8782be375d4SJed Brown if (ratio == 1) { 8792be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 8802be375d4SJed Brown PetscFunctionReturn(0); 8812be375d4SJed Brown } 8822be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 8832be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 8842be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 8852be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 8862be375d4SJed Brown if (i == m-1) lc[i] = want; 8872be375d4SJed Brown else { 8882be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 8892be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 8902be375d4SJed Brown * node is within one stencil width. */ 8912be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 8922be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 8932be375d4SJed Brown * fine node is within one stencil width. */ 8942be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 8952be375d4SJed Brown if (want < 0 || want > remaining 896ce94432eSBarry 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"); 8972be375d4SJed Brown } 8982be375d4SJed Brown lc[i] = want; 8992be375d4SJed Brown startc += lc[i]; 9002be375d4SJed Brown startf += lf[i]; 9012be375d4SJed Brown remaining -= lc[i]; 9022be375d4SJed Brown } 9032be375d4SJed Brown PetscFunctionReturn(0); 9042be375d4SJed Brown } 9052be375d4SJed Brown 9062be375d4SJed Brown #undef __FUNCT__ 90794013140SBarry Smith #define __FUNCT__ "DMRefine_DA" 9087087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 90947c6ae99SBarry Smith { 91047c6ae99SBarry Smith PetscErrorCode ierr; 911f3141302SJed Brown PetscInt M,N,P,i; 9129a42bb27SBarry Smith DM da2; 91347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 91447c6ae99SBarry Smith 91547c6ae99SBarry Smith PetscFunctionBegin; 91647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 91747c6ae99SBarry Smith PetscValidPointer(daref,3); 91847c6ae99SBarry Smith 9191321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 92047c6ae99SBarry Smith M = dd->refine_x*dd->M; 92147c6ae99SBarry Smith } else { 92247c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 92347c6ae99SBarry Smith } 9241321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 9251860e6e9SBarry Smith if (dd->dim > 1) { 92647c6ae99SBarry Smith N = dd->refine_y*dd->N; 92747c6ae99SBarry Smith } else { 9281860e6e9SBarry Smith N = 1; 9291860e6e9SBarry Smith } 9301860e6e9SBarry Smith } else { 93147c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 93247c6ae99SBarry Smith } 9331321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 9341860e6e9SBarry Smith if (dd->dim > 2) { 93547c6ae99SBarry Smith P = dd->refine_z*dd->P; 93647c6ae99SBarry Smith } else { 9371860e6e9SBarry Smith P = 1; 9381860e6e9SBarry Smith } 9391860e6e9SBarry Smith } else { 94047c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 94147c6ae99SBarry Smith } 942ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 943397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 944397b6216SJed Brown ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr); 945397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 946397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 947397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 948397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 949397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 950397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 95147c6ae99SBarry Smith if (dd->dim == 3) { 95247c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 953*dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 9541321219cSEthan 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); 9551321219cSEthan 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); 9561321219cSEthan 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); 957397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 95847c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 95947c6ae99SBarry Smith } else if (dd->dim == 2) { 96047c6ae99SBarry Smith PetscInt *lx,*ly; 961*dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 9621321219cSEthan 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); 9631321219cSEthan 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); 9640298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 96547c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 96647c6ae99SBarry Smith } else if (dd->dim == 1) { 96747c6ae99SBarry Smith PetscInt *lx; 96847c6ae99SBarry Smith ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 9691321219cSEthan 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); 9700298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 97147c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 97247c6ae99SBarry Smith } 97347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 97447c6ae99SBarry Smith 97547c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 97625296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 97725296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 97847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 97947c6ae99SBarry Smith dd2->interptype = dd->interptype; 98047c6ae99SBarry Smith 98147c6ae99SBarry Smith /* copy fill information if given */ 98247c6ae99SBarry Smith if (dd->dfill) { 98347c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 98447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 98547c6ae99SBarry Smith } 98647c6ae99SBarry Smith if (dd->ofill) { 98747c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 98847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 98947c6ae99SBarry Smith } 99047c6ae99SBarry Smith /* copy the refine information */ 991397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 992397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 993397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 99447c6ae99SBarry Smith 99547c6ae99SBarry Smith /* copy vector type information */ 996c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 997ddcf8b74SDave May 998efd51863SBarry Smith dd2->lf = dd->lf; 999efd51863SBarry Smith dd2->lj = dd->lj; 1000efd51863SBarry Smith 10016e87535bSJed Brown da2->leveldown = da->leveldown; 10026e87535bSJed Brown da2->levelup = da->levelup + 1; 10038865f1eaSKarl Rupp 10046e87535bSJed Brown ierr = DMSetFromOptions(da2);CHKERRQ(ierr); 10056e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 1006146574abSBarry Smith ierr = DMViewFromOptions(da2,NULL,"-dm_view");CHKERRQ(ierr); 10076e87535bSJed Brown 1008ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10096636e97aSMatthew G Knepley if (da->coordinates) { 1010ddcf8b74SDave May DM cdaf,cdac; 1011ddcf8b74SDave May Vec coordsc,coordsf; 1012ddcf8b74SDave May Mat II; 1013ddcf8b74SDave May 10146636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 10156636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 10166636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 1017b61d3410SDave May /* force creation of the coordinate vector */ 1018b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 10196636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 10200298fd71SBarry Smith ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 1021ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 1022fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 1023ddcf8b74SDave May } 1024397b6216SJed Brown 1025f3141302SJed Brown for (i=0; i<da->bs; i++) { 1026f3141302SJed Brown const char *fieldname; 1027f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1028f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1029f3141302SJed Brown } 1030397b6216SJed Brown 103147c6ae99SBarry Smith *daref = da2; 103247c6ae99SBarry Smith PetscFunctionReturn(0); 103347c6ae99SBarry Smith } 103447c6ae99SBarry Smith 1035397b6216SJed Brown 103647c6ae99SBarry Smith #undef __FUNCT__ 103794013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA" 10387087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 103947c6ae99SBarry Smith { 104047c6ae99SBarry Smith PetscErrorCode ierr; 1041397b6216SJed Brown PetscInt M,N,P,i; 10429a42bb27SBarry Smith DM da2; 104347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 104447c6ae99SBarry Smith 104547c6ae99SBarry Smith PetscFunctionBegin; 104647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 104747c6ae99SBarry Smith PetscValidPointer(daref,3); 104847c6ae99SBarry Smith 10491321219cSEthan Coon if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1050397b6216SJed Brown M = dd->M / dd->coarsen_x; 105147c6ae99SBarry Smith } else { 1052397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 105347c6ae99SBarry Smith } 10541321219cSEthan Coon if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 10551860e6e9SBarry Smith if (dd->dim > 1) { 1056397b6216SJed Brown N = dd->N / dd->coarsen_y; 105747c6ae99SBarry Smith } else { 10581860e6e9SBarry Smith N = 1; 10591860e6e9SBarry Smith } 10601860e6e9SBarry Smith } else { 1061397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 106247c6ae99SBarry Smith } 10631321219cSEthan Coon if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 10641860e6e9SBarry Smith if (dd->dim > 2) { 1065397b6216SJed Brown P = dd->P / dd->coarsen_z; 106647c6ae99SBarry Smith } else { 10671860e6e9SBarry Smith P = 1; 10681860e6e9SBarry Smith } 10691860e6e9SBarry Smith } else { 1070397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 107147c6ae99SBarry Smith } 1072ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 1073397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 1074397b6216SJed Brown ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr); 1075397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 1076397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 1077397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 1078397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 1079397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 1080397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 108147c6ae99SBarry Smith if (dd->dim == 3) { 10822be375d4SJed Brown PetscInt *lx,*ly,*lz; 1083*dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 1084397b6216SJed 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); 1085397b6216SJed 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); 1086397b6216SJed 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); 1087397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 10882be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 108947c6ae99SBarry Smith } else if (dd->dim == 2) { 10902be375d4SJed Brown PetscInt *lx,*ly; 1091*dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 1092397b6216SJed 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); 1093397b6216SJed 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); 10940298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 10952be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 109647c6ae99SBarry Smith } else if (dd->dim == 1) { 10972be375d4SJed Brown PetscInt *lx; 10982be375d4SJed Brown ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 1099397b6216SJed 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); 11000298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 11012be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 110247c6ae99SBarry Smith } 110347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 110447c6ae99SBarry Smith 11054dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 110625296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 110725296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 110847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 110947c6ae99SBarry Smith dd2->interptype = dd->interptype; 111047c6ae99SBarry Smith 111147c6ae99SBarry Smith /* copy fill information if given */ 111247c6ae99SBarry Smith if (dd->dfill) { 111347c6ae99SBarry Smith ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 111447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 111547c6ae99SBarry Smith } 111647c6ae99SBarry Smith if (dd->ofill) { 111747c6ae99SBarry Smith ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 111847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 111947c6ae99SBarry Smith } 112047c6ae99SBarry Smith /* copy the refine information */ 1121397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1122397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1123397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 112447c6ae99SBarry Smith 112547c6ae99SBarry Smith /* copy vector type information */ 1126c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 112747c6ae99SBarry Smith 1128644e2e5bSBarry Smith dd2->lf = dd->lf; 1129644e2e5bSBarry Smith dd2->lj = dd->lj; 1130644e2e5bSBarry Smith 11316e87535bSJed Brown da2->leveldown = da->leveldown + 1; 11326e87535bSJed Brown da2->levelup = da->levelup; 11338865f1eaSKarl Rupp 11346e87535bSJed Brown ierr = DMSetFromOptions(da2);CHKERRQ(ierr); 11356e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 1136146574abSBarry Smith ierr = DMViewFromOptions(da2,NULL,"-dm_view");CHKERRQ(ierr); 11376e87535bSJed Brown 1138ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 11396636e97aSMatthew G Knepley if (da->coordinates) { 1140ddcf8b74SDave May DM cdaf,cdac; 1141ddcf8b74SDave May Vec coordsc,coordsf; 1142ddcf8b74SDave May VecScatter inject; 1143ddcf8b74SDave May 11446636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); 11456636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); 11466636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); 1147b61d3410SDave May /* force creation of the coordinate vector */ 1148b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 11496636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 1150ddcf8b74SDave May 1151e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 1152ddcf8b74SDave May ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1153ddcf8b74SDave May ierr = VecScatterEnd(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1154fcfd50ebSBarry Smith ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 1155ddcf8b74SDave May } 1156f98405f7SJed Brown 1157f98405f7SJed Brown for (i=0; i<da->bs; i++) { 1158f98405f7SJed Brown const char *fieldname; 1159f98405f7SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1160f98405f7SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1161f98405f7SJed Brown } 1162f98405f7SJed Brown 116347c6ae99SBarry Smith *daref = da2; 116447c6ae99SBarry Smith PetscFunctionReturn(0); 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith 116747c6ae99SBarry Smith #undef __FUNCT__ 116894013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA" 11697087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 117047c6ae99SBarry Smith { 117147c6ae99SBarry Smith PetscErrorCode ierr; 117247c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 117347c6ae99SBarry Smith 117447c6ae99SBarry Smith PetscFunctionBegin; 117547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1176ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 117747c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 117847c6ae99SBarry Smith PetscValidPointer(daf,3); 117947c6ae99SBarry Smith 1180aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 1181*dcca6d9dSJed Brown ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr); 118247c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 1183aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 118447c6ae99SBarry Smith } 118547c6ae99SBarry Smith n = nlevels; 11860298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 118747c6ae99SBarry Smith n = nlevels; 11880298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 118947c6ae99SBarry Smith n = nlevels; 11900298fd71SBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 119147c6ae99SBarry Smith 1192aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1193ce94432eSBarry Smith ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 119447c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1195aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1196ce94432eSBarry Smith ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 119747c6ae99SBarry Smith } 119847c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 119947c6ae99SBarry Smith PetscFunctionReturn(0); 120047c6ae99SBarry Smith } 120147c6ae99SBarry Smith 120247c6ae99SBarry Smith #undef __FUNCT__ 120394013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA" 12047087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 120547c6ae99SBarry Smith { 120647c6ae99SBarry Smith PetscErrorCode ierr; 120747c6ae99SBarry Smith PetscInt i; 120847c6ae99SBarry Smith 120947c6ae99SBarry Smith PetscFunctionBegin; 121047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1211ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 121247c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 121347c6ae99SBarry Smith PetscValidPointer(dac,3); 1214ce94432eSBarry Smith ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 121547c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1216ce94432eSBarry Smith ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 121747c6ae99SBarry Smith } 121847c6ae99SBarry Smith PetscFunctionReturn(0); 121947c6ae99SBarry Smith } 1220