1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 347c6ae99SBarry Smith /*@ 4e43dc8daSBarry Smith DMDASetSizes - Sets the number of grid points in the three dimensional directions 547c6ae99SBarry Smith 6aa219208SBarry Smith Logically Collective on DMDA 747c6ae99SBarry Smith 847c6ae99SBarry Smith Input Parameters: 9aa219208SBarry Smith + da - the DMDA 10e43dc8daSBarry Smith . M - the global X size 11e43dc8daSBarry Smith . N - the global Y size 12e43dc8daSBarry Smith - P - the global Z size 1347c6ae99SBarry Smith 1447c6ae99SBarry Smith Level: intermediate 1547c6ae99SBarry Smith 16628e512eSPatrick Sanan Developer Notes: 17628e512eSPatrick Sanan Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points 18e43dc8daSBarry Smith 19628e512eSPatrick Sanan .seealso: PetscSplitOwnership() 2047c6ae99SBarry Smith @*/ 217087cfbeSBarry Smith PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 2247c6ae99SBarry Smith { 2347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith PetscFunctionBegin; 26a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 2747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,M,2); 2847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,N,3); 2947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,P,4); 30ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 31e43dc8daSBarry Smith if (M < 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive"); 32e43dc8daSBarry Smith if (N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive"); 33e43dc8daSBarry Smith if (P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Z direction must be positive"); 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith dd->M = M; 3647c6ae99SBarry Smith dd->N = N; 3747c6ae99SBarry Smith dd->P = P; 3847c6ae99SBarry Smith PetscFunctionReturn(0); 3947c6ae99SBarry Smith } 4047c6ae99SBarry Smith 4147c6ae99SBarry Smith /*@ 42aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 4347c6ae99SBarry Smith 44aa219208SBarry Smith Logically Collective on DMDA 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith Input Parameters: 47aa219208SBarry Smith + da - the DMDA 4847c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE) 4947c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE) 5047c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE) 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith Level: intermediate 5347c6ae99SBarry Smith 54628e512eSPatrick Sanan .seealso: DMDASetSizes(), PetscSplitOwnership() 5547c6ae99SBarry Smith @*/ 567087cfbeSBarry Smith PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 5747c6ae99SBarry Smith { 5847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 59e3f3e4b6SBarry Smith PetscErrorCode ierr; 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith PetscFunctionBegin; 62a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 6347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 6447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 6547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 66ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 6747c6ae99SBarry Smith dd->m = m; 6847c6ae99SBarry Smith dd->n = n; 6947c6ae99SBarry Smith dd->p = p; 70c73cfb54SMatthew G. Knepley if (da->dim == 2) { 71d3be247eSBarry Smith PetscMPIInt size; 72d3be247eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 73e3f3e4b6SBarry Smith if ((dd->m > 0) && (dd->n < 0)) { 74e3f3e4b6SBarry Smith dd->n = size/dd->m; 75ce94432eSBarry 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); 76e3f3e4b6SBarry Smith } 77e3f3e4b6SBarry Smith if ((dd->n > 0) && (dd->m < 0)) { 78e3f3e4b6SBarry Smith dd->m = size/dd->n; 79ce94432eSBarry 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); 80e3f3e4b6SBarry Smith } 81e3f3e4b6SBarry Smith } 8247c6ae99SBarry Smith PetscFunctionReturn(0); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith /*@ 861321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 8747c6ae99SBarry Smith 8847c6ae99SBarry Smith Not collective 8947c6ae99SBarry Smith 9047c6ae99SBarry Smith Input Parameter: 91aa219208SBarry Smith + da - The DMDA 92bff4a2f0SMatthew G. Knepley - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC 9347c6ae99SBarry Smith 9447c6ae99SBarry Smith Level: intermediate 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith .keywords: distributed array, periodicity 97bff4a2f0SMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType 9847c6ae99SBarry Smith @*/ 99bff4a2f0SMatthew G. Knepley PetscErrorCode DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz) 10047c6ae99SBarry Smith { 10147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith PetscFunctionBegin; 104a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 1055a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bx,2); 1065a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,by,3); 1075a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bz,4); 108ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 1091321219cSEthan Coon dd->bx = bx; 1101321219cSEthan Coon dd->by = by; 1111321219cSEthan Coon dd->bz = bz; 11247c6ae99SBarry Smith PetscFunctionReturn(0); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith /*@ 116aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith Not collective 11947c6ae99SBarry Smith 12059f3ab6dSMatthew G. Knepley Input Parameters: 121aa219208SBarry Smith + da - The DMDA 12247c6ae99SBarry Smith - dof - Number of degrees of freedom 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith Level: intermediate 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 127fb6725baSMatthew G. Knepley .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA 12847c6ae99SBarry Smith @*/ 12954cfb0beSLisandro Dalcin PetscErrorCode DMDASetDof(DM da, PetscInt dof) 13047c6ae99SBarry Smith { 13147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith PetscFunctionBegin; 134a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 13554cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da,dof,2); 136ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 13747c6ae99SBarry Smith dd->w = dof; 1381411c6eeSJed Brown da->bs = dof; 13947c6ae99SBarry Smith PetscFunctionReturn(0); 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith 142fb6725baSMatthew G. Knepley /*@ 143fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 144fb6725baSMatthew G. Knepley 145fb6725baSMatthew G. Knepley Not collective 146fb6725baSMatthew G. Knepley 147fb6725baSMatthew G. Knepley Input Parameter: 148fb6725baSMatthew G. Knepley . da - The DMDA 149fb6725baSMatthew G. Knepley 150fb6725baSMatthew G. Knepley Output Parameter: 151fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom 152fb6725baSMatthew G. Knepley 153fb6725baSMatthew G. Knepley Level: intermediate 154fb6725baSMatthew G. Knepley 155fb6725baSMatthew G. Knepley .keywords: distributed array, degrees of freedom 156fb6725baSMatthew G. Knepley .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA 157fb6725baSMatthew G. Knepley @*/ 158fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 159fb6725baSMatthew G. Knepley { 160fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 161fb6725baSMatthew G. Knepley 162fb6725baSMatthew G. Knepley PetscFunctionBegin; 163a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 164fb6725baSMatthew G. Knepley PetscValidPointer(dof,2); 165fb6725baSMatthew G. Knepley *dof = dd->w; 166fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 167fb6725baSMatthew G. Knepley } 168fb6725baSMatthew G. Knepley 1697ddda789SPeter Brune /*@ 1707ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1717ddda789SPeter Brune 1727ddda789SPeter Brune Not collective 1737ddda789SPeter Brune 1747ddda789SPeter Brune Input Parameters: 1757ddda789SPeter Brune . da - The DMDA 1767ddda789SPeter Brune 1777ddda789SPeter Brune Output Parameters: 1787ddda789SPeter Brune + x - Overlap in the x direction 1797ddda789SPeter Brune . y - Overlap in the y direction 1807ddda789SPeter Brune - z - Overlap in the z direction 1817ddda789SPeter Brune 1827ddda789SPeter Brune Level: intermediate 1837ddda789SPeter Brune 1847ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 1857ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA 1867ddda789SPeter Brune @*/ 1877ddda789SPeter Brune PetscErrorCode DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z) 1887ddda789SPeter Brune { 1897ddda789SPeter Brune DM_DA *dd = (DM_DA*)da->data; 1907ddda789SPeter Brune 1917ddda789SPeter Brune PetscFunctionBegin; 192a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 1937ddda789SPeter Brune if (x) *x = dd->xol; 1947ddda789SPeter Brune if (y) *y = dd->yol; 1957ddda789SPeter Brune if (z) *z = dd->zol; 1967ddda789SPeter Brune PetscFunctionReturn(0); 1977ddda789SPeter Brune } 1987ddda789SPeter Brune 19988661749SPeter Brune /*@ 20088661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 20188661749SPeter Brune 20288661749SPeter Brune Not collective 20388661749SPeter Brune 2047ddda789SPeter Brune Input Parameters: 20588661749SPeter Brune + da - The DMDA 2067ddda789SPeter Brune . x - Overlap in the x direction 2077ddda789SPeter Brune . y - Overlap in the y direction 2087ddda789SPeter Brune - z - Overlap in the z direction 20988661749SPeter Brune 21088661749SPeter Brune Level: intermediate 21188661749SPeter Brune 2127ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 2137ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA 21488661749SPeter Brune @*/ 2157ddda789SPeter Brune PetscErrorCode DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z) 21688661749SPeter Brune { 21788661749SPeter Brune DM_DA *dd = (DM_DA*)da->data; 21888661749SPeter Brune 21988661749SPeter Brune PetscFunctionBegin; 220a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2217ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,x,2); 2227ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,y,3); 2237ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,z,4); 2247ddda789SPeter Brune dd->xol = x; 2257ddda789SPeter Brune dd->yol = y; 2267ddda789SPeter Brune dd->zol = z; 22788661749SPeter Brune PetscFunctionReturn(0); 22888661749SPeter Brune } 22988661749SPeter Brune 2303e7870d2SPeter Brune 2313e7870d2SPeter Brune /*@ 2323e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2333e7870d2SPeter Brune 2343e7870d2SPeter Brune Not collective 2353e7870d2SPeter Brune 2363e7870d2SPeter Brune Input Parameters: 2373e7870d2SPeter Brune . da - The DMDA 2383e7870d2SPeter Brune 2393e7870d2SPeter Brune Output Parameters: 240*a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition 2413e7870d2SPeter Brune 2423e7870d2SPeter Brune Level: intermediate 2433e7870d2SPeter Brune 2443e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2453e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA 2463e7870d2SPeter Brune @*/ 2473e7870d2SPeter Brune PetscErrorCode DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub) 2483e7870d2SPeter Brune { 2493e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2503e7870d2SPeter Brune 2513e7870d2SPeter Brune PetscFunctionBegin; 252a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2533e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2543e7870d2SPeter Brune PetscFunctionReturn(0); 2553e7870d2SPeter Brune } 2563e7870d2SPeter Brune 2573e7870d2SPeter Brune /*@ 2583e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2593e7870d2SPeter Brune 2603e7870d2SPeter Brune Not collective 2613e7870d2SPeter Brune 2623e7870d2SPeter Brune Input Parameters: 2633e7870d2SPeter Brune + da - The DMDA 2643e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2653e7870d2SPeter Brune 2663e7870d2SPeter Brune Level: intermediate 2673e7870d2SPeter Brune 2683e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2693e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA 2703e7870d2SPeter Brune @*/ 2713e7870d2SPeter Brune PetscErrorCode DMDASetNumLocalSubDomains(DM da,PetscInt Nsub) 2723e7870d2SPeter Brune { 2733e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2743e7870d2SPeter Brune 2753e7870d2SPeter Brune PetscFunctionBegin; 276a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2773e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da,Nsub,2); 2783e7870d2SPeter Brune dd->Nsub = Nsub; 2793e7870d2SPeter Brune PetscFunctionReturn(0); 2803e7870d2SPeter Brune } 2813e7870d2SPeter Brune 282d886c4f4SPeter Brune /*@ 283d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 284d886c4f4SPeter Brune 285d886c4f4SPeter Brune Collective on DA 286d886c4f4SPeter Brune 287d886c4f4SPeter Brune Input Parameter: 288d886c4f4SPeter Brune + da - The DMDA 289d886c4f4SPeter Brune . xo - The offset in the x direction 290d886c4f4SPeter Brune . yo - The offset in the y direction 291d886c4f4SPeter Brune - zo - The offset in the z direction 292d886c4f4SPeter Brune 293d886c4f4SPeter Brune Level: intermediate 294d886c4f4SPeter Brune 29595452b02SPatrick Sanan Notes: 29695452b02SPatrick Sanan This is used primarily to overlap a computation on a local DA with that on a global DA without 297d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 298d886c4f4SPeter Brune 299d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 300d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 301d886c4f4SPeter Brune @*/ 30295c13181SPeter Brune PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 303d886c4f4SPeter Brune { 30495c13181SPeter Brune PetscErrorCode ierr; 305d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 306d886c4f4SPeter Brune 307d886c4f4SPeter Brune PetscFunctionBegin; 308a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 309d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da,xo,2); 31095c13181SPeter Brune PetscValidLogicalCollectiveInt(da,yo,3); 31195c13181SPeter Brune PetscValidLogicalCollectiveInt(da,zo,4); 31295c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Mo,5); 31395c13181SPeter Brune PetscValidLogicalCollectiveInt(da,No,6); 31495c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Po,7); 315d886c4f4SPeter Brune dd->xo = xo; 316d886c4f4SPeter Brune dd->yo = yo; 317d886c4f4SPeter Brune dd->zo = zo; 31895c13181SPeter Brune dd->Mo = Mo; 31995c13181SPeter Brune dd->No = No; 32095c13181SPeter Brune dd->Po = Po; 32195c13181SPeter Brune 32295c13181SPeter Brune if (da->coordinateDM) { 32395c13181SPeter Brune ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 32495c13181SPeter Brune } 325d886c4f4SPeter Brune PetscFunctionReturn(0); 326d886c4f4SPeter Brune } 327d886c4f4SPeter Brune 328d886c4f4SPeter Brune /*@ 329d886c4f4SPeter Brune DMDAGetOffset - Gets the index offset of the DA. 330d886c4f4SPeter Brune 331d886c4f4SPeter Brune Not collective 332d886c4f4SPeter Brune 333d886c4f4SPeter Brune Input Parameter: 334d886c4f4SPeter Brune . da - The DMDA 335d886c4f4SPeter Brune 336d886c4f4SPeter Brune Output Parameters: 337d886c4f4SPeter Brune + xo - The offset in the x direction 338d886c4f4SPeter Brune . yo - The offset in the y direction 33995c13181SPeter Brune . zo - The offset in the z direction 34095c13181SPeter Brune . Mo - The global size in the x direction 34195c13181SPeter Brune . No - The global size in the y direction 34295c13181SPeter Brune - Po - The global size in the z direction 343d886c4f4SPeter Brune 344d886c4f4SPeter Brune Level: intermediate 345d886c4f4SPeter Brune 346d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 347d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray() 348d886c4f4SPeter Brune @*/ 34995c13181SPeter Brune PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 350d886c4f4SPeter Brune { 351d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 352d886c4f4SPeter Brune 353d886c4f4SPeter Brune PetscFunctionBegin; 354a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 355d886c4f4SPeter Brune if (xo) *xo = dd->xo; 356d886c4f4SPeter Brune if (yo) *yo = dd->yo; 357d886c4f4SPeter Brune if (zo) *zo = dd->zo; 35895c13181SPeter Brune if (Mo) *Mo = dd->Mo; 35995c13181SPeter Brune if (No) *No = dd->No; 36095c13181SPeter Brune if (Po) *Po = dd->Po; 361d886c4f4SPeter Brune PetscFunctionReturn(0); 362d886c4f4SPeter Brune } 363d886c4f4SPeter Brune 36440234c92SPeter Brune /*@ 36540234c92SPeter Brune DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM. 36640234c92SPeter Brune 36740234c92SPeter Brune Not collective 36840234c92SPeter Brune 36940234c92SPeter Brune Input Parameter: 37040234c92SPeter Brune . da - The DMDA 37140234c92SPeter Brune 37240234c92SPeter Brune Output Parameters: 37340234c92SPeter Brune + xs - The start of the region in x 37440234c92SPeter Brune . ys - The start of the region in y 37540234c92SPeter Brune . zs - The start of the region in z 37640234c92SPeter Brune . xs - The size of the region in x 37740234c92SPeter Brune . ys - The size of the region in y 378*a2b725a8SWilliam Gropp - zs - The size of the region in z 37940234c92SPeter Brune 38040234c92SPeter Brune Level: intermediate 38140234c92SPeter Brune 38240234c92SPeter Brune .keywords: distributed array, degrees of freedom 38340234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 38440234c92SPeter Brune @*/ 38540234c92SPeter Brune PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 38640234c92SPeter Brune { 38740234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 38840234c92SPeter Brune 38940234c92SPeter Brune PetscFunctionBegin; 390a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 39140234c92SPeter Brune if (xs) *xs = dd->nonxs; 39240234c92SPeter Brune if (ys) *ys = dd->nonys; 39340234c92SPeter Brune if (zs) *zs = dd->nonzs; 39440234c92SPeter Brune if (xm) *xm = dd->nonxm; 39540234c92SPeter Brune if (ym) *ym = dd->nonym; 39640234c92SPeter Brune if (zm) *zm = dd->nonzm; 39740234c92SPeter Brune PetscFunctionReturn(0); 39840234c92SPeter Brune } 39940234c92SPeter Brune 40040234c92SPeter Brune 40140234c92SPeter Brune /*@ 40240234c92SPeter Brune DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM. 40340234c92SPeter Brune 40440234c92SPeter Brune Collective on DA 40540234c92SPeter Brune 40640234c92SPeter Brune Input Parameter: 40740234c92SPeter Brune + da - The DMDA 40840234c92SPeter Brune . xs - The start of the region in x 40940234c92SPeter Brune . ys - The start of the region in y 41040234c92SPeter Brune . zs - The start of the region in z 41140234c92SPeter Brune . xs - The size of the region in x 41240234c92SPeter Brune . ys - The size of the region in y 413*a2b725a8SWilliam Gropp - zs - The size of the region in z 41440234c92SPeter Brune 41540234c92SPeter Brune Level: intermediate 41640234c92SPeter Brune 41740234c92SPeter Brune .keywords: distributed array, degrees of freedom 41840234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 41940234c92SPeter Brune @*/ 42040234c92SPeter Brune PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 42140234c92SPeter Brune { 42240234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 42340234c92SPeter Brune 42440234c92SPeter Brune PetscFunctionBegin; 425a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 42640234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xs,2); 42740234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ys,3); 42840234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zs,4); 42940234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xm,5); 43040234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ym,6); 43140234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zm,7); 43240234c92SPeter Brune dd->nonxs = xs; 43340234c92SPeter Brune dd->nonys = ys; 43440234c92SPeter Brune dd->nonzs = zs; 43540234c92SPeter Brune dd->nonxm = xm; 43640234c92SPeter Brune dd->nonym = ym; 43740234c92SPeter Brune dd->nonzm = zm; 43840234c92SPeter Brune 43940234c92SPeter Brune PetscFunctionReturn(0); 44040234c92SPeter Brune } 44188661749SPeter Brune 44247c6ae99SBarry Smith /*@ 443aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 44447c6ae99SBarry Smith 445aa219208SBarry Smith Logically Collective on DMDA 44647c6ae99SBarry Smith 44747c6ae99SBarry Smith Input Parameter: 448aa219208SBarry Smith + da - The DMDA 449aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 45047c6ae99SBarry Smith 45147c6ae99SBarry Smith Level: intermediate 45247c6ae99SBarry Smith 45347c6ae99SBarry Smith .keywords: distributed array, stencil 45496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 45547c6ae99SBarry Smith @*/ 4567087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 45747c6ae99SBarry Smith { 45847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 45947c6ae99SBarry Smith 46047c6ae99SBarry Smith PetscFunctionBegin; 461a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 46247c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 463ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 46447c6ae99SBarry Smith dd->stencil_type = stype; 46547c6ae99SBarry Smith PetscFunctionReturn(0); 46647c6ae99SBarry Smith } 46747c6ae99SBarry Smith 468fb6725baSMatthew G. Knepley /*@ 469fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 470fb6725baSMatthew G. Knepley 471fb6725baSMatthew G. Knepley Not collective 472fb6725baSMatthew G. Knepley 473fb6725baSMatthew G. Knepley Input Parameter: 474fb6725baSMatthew G. Knepley . da - The DMDA 475fb6725baSMatthew G. Knepley 476fb6725baSMatthew G. Knepley Output Parameter: 477fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 478fb6725baSMatthew G. Knepley 479fb6725baSMatthew G. Knepley Level: intermediate 480fb6725baSMatthew G. Knepley 481fb6725baSMatthew G. Knepley .keywords: distributed array, stencil 482fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 483fb6725baSMatthew G. Knepley @*/ 484fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 485fb6725baSMatthew G. Knepley { 486fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA*)da->data; 487fb6725baSMatthew G. Knepley 488fb6725baSMatthew G. Knepley PetscFunctionBegin; 489a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 490fb6725baSMatthew G. Knepley PetscValidPointer(stype,2); 491fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 492fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 493fb6725baSMatthew G. Knepley } 494fb6725baSMatthew G. Knepley 49547c6ae99SBarry Smith /*@ 496aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 49747c6ae99SBarry Smith 498aa219208SBarry Smith Logically Collective on DMDA 49947c6ae99SBarry Smith 50047c6ae99SBarry Smith Input Parameter: 501aa219208SBarry Smith + da - The DMDA 50247c6ae99SBarry Smith - width - The stencil width 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith Level: intermediate 50547c6ae99SBarry Smith 50647c6ae99SBarry Smith .keywords: distributed array, stencil 50796e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 50847c6ae99SBarry Smith @*/ 5097087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 51047c6ae99SBarry Smith { 51147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 51247c6ae99SBarry Smith 51347c6ae99SBarry Smith PetscFunctionBegin; 514a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 51547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 516ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 51747c6ae99SBarry Smith dd->s = width; 51847c6ae99SBarry Smith PetscFunctionReturn(0); 51947c6ae99SBarry Smith } 52047c6ae99SBarry Smith 521fb6725baSMatthew G. Knepley /*@ 522fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 523fb6725baSMatthew G. Knepley 524fb6725baSMatthew G. Knepley Not collective 525fb6725baSMatthew G. Knepley 526fb6725baSMatthew G. Knepley Input Parameter: 527fb6725baSMatthew G. Knepley . da - The DMDA 528fb6725baSMatthew G. Knepley 529fb6725baSMatthew G. Knepley Output Parameter: 530fb6725baSMatthew G. Knepley . width - The stencil width 531fb6725baSMatthew G. Knepley 532fb6725baSMatthew G. Knepley Level: intermediate 533fb6725baSMatthew G. Knepley 534fb6725baSMatthew G. Knepley .keywords: distributed array, stencil 535fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 536fb6725baSMatthew G. Knepley @*/ 537fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 538fb6725baSMatthew G. Knepley { 539fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 540fb6725baSMatthew G. Knepley 541fb6725baSMatthew G. Knepley PetscFunctionBegin; 542a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 543fb6725baSMatthew G. Knepley PetscValidPointer(width,2); 544fb6725baSMatthew G. Knepley *width = dd->s; 545fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 546fb6725baSMatthew G. Knepley } 547fb6725baSMatthew G. Knepley 548aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 54947c6ae99SBarry Smith { 55047c6ae99SBarry Smith PetscInt i,sum; 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith PetscFunctionBegin; 553ce94432eSBarry Smith if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 55447c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 555ce94432eSBarry Smith if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 55647c6ae99SBarry Smith PetscFunctionReturn(0); 55747c6ae99SBarry Smith } 55847c6ae99SBarry Smith 55947c6ae99SBarry Smith /*@ 560aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 56147c6ae99SBarry Smith 562aa219208SBarry Smith Logically Collective on DMDA 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith Input Parameter: 565aa219208SBarry Smith + da - The DMDA 5660298fd71SBarry 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 5670298fd71SBarry 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 5680298fd71SBarry 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. 56947c6ae99SBarry Smith 57047c6ae99SBarry Smith Level: intermediate 57147c6ae99SBarry Smith 572e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 573e3f3e4b6SBarry Smith 57447c6ae99SBarry Smith .keywords: distributed array 57596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 57647c6ae99SBarry Smith @*/ 5777087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 57847c6ae99SBarry Smith { 57947c6ae99SBarry Smith PetscErrorCode ierr; 58047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 58147c6ae99SBarry Smith 58247c6ae99SBarry Smith PetscFunctionBegin; 583a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 584ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 58547c6ae99SBarry Smith if (lx) { 586ce94432eSBarry Smith if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 587aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 58847c6ae99SBarry Smith if (!dd->lx) { 589785e854fSJed Brown ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr); 59047c6ae99SBarry Smith } 59147c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith if (ly) { 594ce94432eSBarry Smith if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 595aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 59647c6ae99SBarry Smith if (!dd->ly) { 597785e854fSJed Brown ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr); 59847c6ae99SBarry Smith } 59947c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 60047c6ae99SBarry Smith } 60147c6ae99SBarry Smith if (lz) { 602ce94432eSBarry Smith if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 603aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 60447c6ae99SBarry Smith if (!dd->lz) { 605785e854fSJed Brown ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr); 60647c6ae99SBarry Smith } 60747c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 60847c6ae99SBarry Smith } 60947c6ae99SBarry Smith PetscFunctionReturn(0); 61047c6ae99SBarry Smith } 61147c6ae99SBarry Smith 61247c6ae99SBarry Smith /*@ 613aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 614e727c939SJed Brown returned by DMCreateInterpolation() 61547c6ae99SBarry Smith 616aa219208SBarry Smith Logically Collective on DMDA 61747c6ae99SBarry Smith 61847c6ae99SBarry Smith Input Parameter: 61947c6ae99SBarry Smith + da - initial distributed array 620*a2b725a8SWilliam Gropp - ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 62147c6ae99SBarry Smith 62247c6ae99SBarry Smith Level: intermediate 62347c6ae99SBarry Smith 62495452b02SPatrick Sanan Notes: 62595452b02SPatrick Sanan you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith .keywords: distributed array, interpolation 62847c6ae99SBarry Smith 62996e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 63047c6ae99SBarry Smith @*/ 6317087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 63247c6ae99SBarry Smith { 63347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 63447c6ae99SBarry Smith 63547c6ae99SBarry Smith PetscFunctionBegin; 636a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 63747c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 63847c6ae99SBarry Smith dd->interptype = ctype; 63947c6ae99SBarry Smith PetscFunctionReturn(0); 64047c6ae99SBarry Smith } 64147c6ae99SBarry Smith 6422dde6fd4SLisandro Dalcin /*@ 6432dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 644e727c939SJed Brown used by DMCreateInterpolation() 6452dde6fd4SLisandro Dalcin 6462dde6fd4SLisandro Dalcin Not Collective 6472dde6fd4SLisandro Dalcin 6482dde6fd4SLisandro Dalcin Input Parameter: 6492dde6fd4SLisandro Dalcin . da - distributed array 6502dde6fd4SLisandro Dalcin 6512dde6fd4SLisandro Dalcin Output Parameter: 6522dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 6532dde6fd4SLisandro Dalcin 6542dde6fd4SLisandro Dalcin Level: intermediate 6552dde6fd4SLisandro Dalcin 6562dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 6572dde6fd4SLisandro Dalcin 658e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 6592dde6fd4SLisandro Dalcin @*/ 6602dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 6612dde6fd4SLisandro Dalcin { 6622dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 6632dde6fd4SLisandro Dalcin 6642dde6fd4SLisandro Dalcin PetscFunctionBegin; 665a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 6662dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 6672dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6682dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 6692dde6fd4SLisandro Dalcin } 67047c6ae99SBarry Smith 6716a119db4SBarry Smith /*@C 672aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 67347c6ae99SBarry Smith processes neighbors. 67447c6ae99SBarry Smith 67547c6ae99SBarry Smith Not Collective 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith Input Parameter: 678aa219208SBarry Smith . da - the DMDA object 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith Output Parameters: 68147c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 68247c6ae99SBarry Smith this process itself is in the list 68347c6ae99SBarry Smith 68495452b02SPatrick Sanan Notes: 68595452b02SPatrick Sanan In 2d the array is of length 9, in 3d of length 27 68647c6ae99SBarry Smith Not supported in 1d 687aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 68847c6ae99SBarry Smith 68995452b02SPatrick Sanan Fortran Notes: 69095452b02SPatrick Sanan In fortran you must pass in an array of the appropriate length. 69147c6ae99SBarry Smith 69247c6ae99SBarry Smith Level: intermediate 69347c6ae99SBarry Smith 69447c6ae99SBarry Smith @*/ 6957087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 69647c6ae99SBarry Smith { 69747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 6985fd66863SKarl Rupp 69947c6ae99SBarry Smith PetscFunctionBegin; 700a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 70147c6ae99SBarry Smith *ranks = dd->neighbors; 70247c6ae99SBarry Smith PetscFunctionReturn(0); 70347c6ae99SBarry Smith } 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith /*@C 706aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith Not Collective 70947c6ae99SBarry Smith 71047c6ae99SBarry Smith Input Parameter: 711aa219208SBarry Smith . da - the DMDA object 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith Output Parameter: 71447c6ae99SBarry Smith + lx - ownership along x direction (optional) 71547c6ae99SBarry Smith . ly - ownership along y direction (optional) 71647c6ae99SBarry Smith - lz - ownership along z direction (optional) 71747c6ae99SBarry Smith 71847c6ae99SBarry Smith Level: intermediate 71947c6ae99SBarry Smith 720aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 72147c6ae99SBarry Smith 72247c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 723aa219208SBarry Smith eighth arguments from DMDAGetInfo() 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 726aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 72747c6ae99SBarry Smith 728e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 729e3f3e4b6SBarry Smith 730aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 73147c6ae99SBarry Smith @*/ 7327087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 73347c6ae99SBarry Smith { 73447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 73547c6ae99SBarry Smith 73647c6ae99SBarry Smith PetscFunctionBegin; 737a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 73847c6ae99SBarry Smith if (lx) *lx = dd->lx; 73947c6ae99SBarry Smith if (ly) *ly = dd->ly; 74047c6ae99SBarry Smith if (lz) *lz = dd->lz; 74147c6ae99SBarry Smith PetscFunctionReturn(0); 74247c6ae99SBarry Smith } 74347c6ae99SBarry Smith 74447c6ae99SBarry Smith /*@ 745aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 74647c6ae99SBarry Smith 747aa219208SBarry Smith Logically Collective on DMDA 74847c6ae99SBarry Smith 74947c6ae99SBarry Smith Input Parameters: 750aa219208SBarry Smith + da - the DMDA object 75147c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 75247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 75347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith Options Database: 75647c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 75747c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 75847c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 75947c6ae99SBarry Smith 76047c6ae99SBarry Smith Level: intermediate 76147c6ae99SBarry Smith 76295452b02SPatrick Sanan Notes: 76395452b02SPatrick Sanan Pass PETSC_IGNORE to leave a value unchanged 76447c6ae99SBarry Smith 765aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 76647c6ae99SBarry Smith @*/ 7677087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 76847c6ae99SBarry Smith { 76947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith PetscFunctionBegin; 772a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 77347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 77447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 77547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 77647c6ae99SBarry Smith 77747c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 77847c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 77947c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 78047c6ae99SBarry Smith PetscFunctionReturn(0); 78147c6ae99SBarry Smith } 78247c6ae99SBarry Smith 78347c6ae99SBarry Smith /*@C 784aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith Not Collective 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith Input Parameter: 789aa219208SBarry Smith . da - the DMDA object 79047c6ae99SBarry Smith 79147c6ae99SBarry Smith Output Parameters: 79247c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 79347c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 79447c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith Level: intermediate 79747c6ae99SBarry Smith 79895452b02SPatrick Sanan Notes: 79995452b02SPatrick Sanan Pass NULL for values you do not need 80047c6ae99SBarry Smith 801aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 80247c6ae99SBarry Smith @*/ 8037087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 80447c6ae99SBarry Smith { 80547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith PetscFunctionBegin; 808a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 80947c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 81047c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 81147c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 81247c6ae99SBarry Smith PetscFunctionReturn(0); 81347c6ae99SBarry Smith } 81447c6ae99SBarry Smith 81547c6ae99SBarry Smith /*@C 816aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 81747c6ae99SBarry Smith 818aa219208SBarry Smith Logically Collective on DMDA 81947c6ae99SBarry Smith 82047c6ae99SBarry Smith Input Parameters: 821aa219208SBarry Smith + da - the DMDA object 822aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith Level: developer 82547c6ae99SBarry Smith 82695452b02SPatrick Sanan Notes: 82795452b02SPatrick Sanan See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 82847c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 82947c6ae99SBarry Smith 8301fd49c25SBarry Smith Not supported from Fortran 8311fd49c25SBarry Smith 832950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 83347c6ae99SBarry Smith @*/ 834b412c318SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*)) 83547c6ae99SBarry Smith { 83647c6ae99SBarry Smith PetscFunctionBegin; 837a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 83825296bd5SBarry Smith da->ops->creatematrix = f; 83947c6ae99SBarry Smith PetscFunctionReturn(0); 84047c6ae99SBarry Smith } 84147c6ae99SBarry Smith 84247c6ae99SBarry Smith /* 84347c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 84447c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 84547c6ae99SBarry Smith 84647c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 84747c6ae99SBarry Smith */ 84894013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 84947c6ae99SBarry Smith { 85047c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 85147c6ae99SBarry Smith PetscErrorCode ierr; 85247c6ae99SBarry Smith 85347c6ae99SBarry Smith PetscFunctionBegin; 854ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 85547c6ae99SBarry Smith if (ratio == 1) { 85647c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 85747c6ae99SBarry Smith PetscFunctionReturn(0); 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 86047c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 86147c6ae99SBarry Smith for (i=0; i<m; i++) { 86247c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 86347c6ae99SBarry Smith if (i == m-1) lf[i] = want; 86447c6ae99SBarry Smith else { 8657aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 8667aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 8677aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 8687aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 8697aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 8707aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 8717aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 8727aca7175SJed Brown /* Make sure all constraints are satisfied */ 87330729d88SBarry Smith if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 874ce94432eSBarry 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"); 87547c6ae99SBarry Smith } 87647c6ae99SBarry Smith lf[i] = want; 87747c6ae99SBarry Smith startc += lc[i]; 87847c6ae99SBarry Smith startf += lf[i]; 87947c6ae99SBarry Smith remaining -= lf[i]; 88047c6ae99SBarry Smith } 88147c6ae99SBarry Smith PetscFunctionReturn(0); 88247c6ae99SBarry Smith } 88347c6ae99SBarry Smith 8842be375d4SJed Brown /* 8852be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 8862be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 8872be375d4SJed Brown 8882be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 8892be375d4SJed Brown */ 8902be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 8912be375d4SJed Brown { 8922be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 8932be375d4SJed Brown PetscErrorCode ierr; 8942be375d4SJed Brown 8952be375d4SJed Brown PetscFunctionBegin; 896ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 8972be375d4SJed Brown if (ratio == 1) { 8982be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 8992be375d4SJed Brown PetscFunctionReturn(0); 9002be375d4SJed Brown } 9012be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 9022be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 9032be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 9042be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 9052be375d4SJed Brown if (i == m-1) lc[i] = want; 9062be375d4SJed Brown else { 9072be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 9082be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 9092be375d4SJed Brown * node is within one stencil width. */ 9102be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 9112be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 9122be375d4SJed Brown * fine node is within one stencil width. */ 9132be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 9142be375d4SJed Brown if (want < 0 || want > remaining 915ce94432eSBarry 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"); 9162be375d4SJed Brown } 9172be375d4SJed Brown lc[i] = want; 9182be375d4SJed Brown startc += lc[i]; 9192be375d4SJed Brown startf += lf[i]; 9202be375d4SJed Brown remaining -= lc[i]; 9212be375d4SJed Brown } 9222be375d4SJed Brown PetscFunctionReturn(0); 9232be375d4SJed Brown } 9242be375d4SJed Brown 9257087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 92647c6ae99SBarry Smith { 92747c6ae99SBarry Smith PetscErrorCode ierr; 928c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 9299a42bb27SBarry Smith DM da2; 93047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith PetscFunctionBegin; 933a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 93447c6ae99SBarry Smith PetscValidPointer(daref,3); 93547c6ae99SBarry Smith 936c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 937bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 93847c6ae99SBarry Smith M = dd->refine_x*dd->M; 93947c6ae99SBarry Smith } else { 94047c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 94147c6ae99SBarry Smith } 942bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 943c73cfb54SMatthew G. Knepley if (dim > 1) { 94447c6ae99SBarry Smith N = dd->refine_y*dd->N; 94547c6ae99SBarry Smith } else { 9461860e6e9SBarry Smith N = 1; 9471860e6e9SBarry Smith } 9481860e6e9SBarry Smith } else { 94947c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 95047c6ae99SBarry Smith } 951bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 952c73cfb54SMatthew G. Knepley if (dim > 2) { 95347c6ae99SBarry Smith P = dd->refine_z*dd->P; 95447c6ae99SBarry Smith } else { 9551860e6e9SBarry Smith P = 1; 9561860e6e9SBarry Smith } 9571860e6e9SBarry Smith } else { 95847c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 95947c6ae99SBarry Smith } 960ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 961397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 962c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 963397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 964397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 965397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 966397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 967397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 968397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 969c73cfb54SMatthew G. Knepley if (dim == 3) { 97047c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 971dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 972bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 973bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 974bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 975397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 97647c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 977c73cfb54SMatthew G. Knepley } else if (dim == 2) { 97847c6ae99SBarry Smith PetscInt *lx,*ly; 979dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 980bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 981bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 9820298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 98347c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 984c73cfb54SMatthew G. Knepley } else if (dim == 1) { 98547c6ae99SBarry Smith PetscInt *lx; 986785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 987bff4a2f0SMatthew G. Knepley ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 9880298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 98947c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 99047c6ae99SBarry Smith } 99147c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 99247c6ae99SBarry Smith 99347c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 99425296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 99525296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 99647c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 99747c6ae99SBarry Smith dd2->interptype = dd->interptype; 99847c6ae99SBarry Smith 99947c6ae99SBarry Smith /* copy fill information if given */ 100047c6ae99SBarry Smith if (dd->dfill) { 1001854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 100247c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 100347c6ae99SBarry Smith } 100447c6ae99SBarry Smith if (dd->ofill) { 1005854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 100647c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 100747c6ae99SBarry Smith } 100847c6ae99SBarry Smith /* copy the refine information */ 1009397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1010397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1011397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 101247c6ae99SBarry Smith 1013897f7067SBarry Smith if (dd->refine_z_hier) { 1014897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) { 1015897f7067SBarry Smith dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; 1016897f7067SBarry Smith } 1017897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) { 1018897f7067SBarry Smith dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; 1019897f7067SBarry Smith } 1020897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1021897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1022897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1023897f7067SBarry Smith } 1024897f7067SBarry Smith if (dd->refine_y_hier) { 1025897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) { 1026897f7067SBarry Smith dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; 1027897f7067SBarry Smith } 1028897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) { 1029897f7067SBarry Smith dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; 1030897f7067SBarry Smith } 1031897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1032897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1033897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1034897f7067SBarry Smith } 1035897f7067SBarry Smith if (dd->refine_x_hier) { 1036897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) { 1037897f7067SBarry Smith dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; 1038897f7067SBarry Smith } 1039897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) { 1040897f7067SBarry Smith dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; 1041897f7067SBarry Smith } 1042897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1043897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1044897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1045897f7067SBarry Smith } 1046897f7067SBarry Smith 1047897f7067SBarry Smith 104847c6ae99SBarry Smith /* copy vector type information */ 1049c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 1050ddcf8b74SDave May 1051efd51863SBarry Smith dd2->lf = dd->lf; 1052efd51863SBarry Smith dd2->lj = dd->lj; 1053efd51863SBarry Smith 10546e87535bSJed Brown da2->leveldown = da->leveldown; 10556e87535bSJed Brown da2->levelup = da->levelup + 1; 10568865f1eaSKarl Rupp 10576e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 10586e87535bSJed Brown 1059ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10606636e97aSMatthew G Knepley if (da->coordinates) { 1061ddcf8b74SDave May DM cdaf,cdac; 1062ddcf8b74SDave May Vec coordsc,coordsf; 1063ddcf8b74SDave May Mat II; 1064ddcf8b74SDave May 10656636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 10666636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 10676636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 1068b61d3410SDave May /* force creation of the coordinate vector */ 1069b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 10706636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 10710298fd71SBarry Smith ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 1072ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 1073fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 1074ddcf8b74SDave May } 1075397b6216SJed Brown 1076f3141302SJed Brown for (i=0; i<da->bs; i++) { 1077f3141302SJed Brown const char *fieldname; 1078f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1079f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1080f3141302SJed Brown } 1081397b6216SJed Brown 108247c6ae99SBarry Smith *daref = da2; 108347c6ae99SBarry Smith PetscFunctionReturn(0); 108447c6ae99SBarry Smith } 108547c6ae99SBarry Smith 1086397b6216SJed Brown 10877087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 108847c6ae99SBarry Smith { 108947c6ae99SBarry Smith PetscErrorCode ierr; 1090c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 10919a42bb27SBarry Smith DM da2; 109247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 109347c6ae99SBarry Smith 109447c6ae99SBarry Smith PetscFunctionBegin; 1095a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 109647c6ae99SBarry Smith PetscValidPointer(daref,3); 109747c6ae99SBarry Smith 1098c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 1099bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1100397b6216SJed Brown M = dd->M / dd->coarsen_x; 110147c6ae99SBarry Smith } else { 1102397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 110347c6ae99SBarry Smith } 1104bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1105c73cfb54SMatthew G. Knepley if (dim > 1) { 1106397b6216SJed Brown N = dd->N / dd->coarsen_y; 110747c6ae99SBarry Smith } else { 11081860e6e9SBarry Smith N = 1; 11091860e6e9SBarry Smith } 11101860e6e9SBarry Smith } else { 1111397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 111247c6ae99SBarry Smith } 1113bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1114c73cfb54SMatthew G. Knepley if (dim > 2) { 1115397b6216SJed Brown P = dd->P / dd->coarsen_z; 111647c6ae99SBarry Smith } else { 11171860e6e9SBarry Smith P = 1; 11181860e6e9SBarry Smith } 11191860e6e9SBarry Smith } else { 1120397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 112147c6ae99SBarry Smith } 1122ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 1123397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 1124c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 1125397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 1126397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 1127397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 1128397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 1129397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 1130397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 1131c73cfb54SMatthew G. Knepley if (dim == 3) { 11322be375d4SJed Brown PetscInt *lx,*ly,*lz; 1133dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 1134bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1135bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 1136bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 1137397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 11382be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 1139c73cfb54SMatthew G. Knepley } else if (dim == 2) { 11402be375d4SJed Brown PetscInt *lx,*ly; 1141dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 1142bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1143bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 11440298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 11452be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 1146c73cfb54SMatthew G. Knepley } else if (dim == 1) { 11472be375d4SJed Brown PetscInt *lx; 1148785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 1149bff4a2f0SMatthew G. Knepley ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 11500298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 11512be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 115247c6ae99SBarry Smith } 115347c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 115447c6ae99SBarry Smith 11554dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 115625296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 115725296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 115847c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 115947c6ae99SBarry Smith dd2->interptype = dd->interptype; 116047c6ae99SBarry Smith 116147c6ae99SBarry Smith /* copy fill information if given */ 116247c6ae99SBarry Smith if (dd->dfill) { 1163854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 116447c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith if (dd->ofill) { 1167854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 116847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 116947c6ae99SBarry Smith } 117047c6ae99SBarry Smith /* copy the refine information */ 1171397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1172397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1173397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 117447c6ae99SBarry Smith 1175897f7067SBarry Smith if (dd->refine_z_hier) { 1176897f7067SBarry Smith if (da->levelup - da->leveldown -1 > -1 && da->levelup - da->leveldown - 1< dd->refine_z_hier_n) { 1177897f7067SBarry Smith dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown - 1]; 1178897f7067SBarry Smith } 1179897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_z_hier_n) { 1180897f7067SBarry Smith dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown - 2]; 1181897f7067SBarry Smith } 1182897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1183897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1184897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1185897f7067SBarry Smith } 1186897f7067SBarry Smith if (dd->refine_y_hier) { 1187897f7067SBarry Smith if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1< dd->refine_y_hier_n) { 1188897f7067SBarry Smith dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown - 1]; 1189897f7067SBarry Smith } 1190897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_y_hier_n) { 1191897f7067SBarry Smith dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown - 2]; 1192897f7067SBarry Smith } 1193897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1194897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1195897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1196897f7067SBarry Smith } 1197897f7067SBarry Smith if (dd->refine_x_hier) { 1198897f7067SBarry Smith if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1 < dd->refine_x_hier_n) { 1199897f7067SBarry Smith dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown - 1]; 1200897f7067SBarry Smith } 1201897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_x_hier_n) { 1202897f7067SBarry Smith dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown - 2]; 1203897f7067SBarry Smith } 1204897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1205897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1206897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1207897f7067SBarry Smith } 1208897f7067SBarry Smith 120947c6ae99SBarry Smith /* copy vector type information */ 1210c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 121147c6ae99SBarry Smith 1212644e2e5bSBarry Smith dd2->lf = dd->lf; 1213644e2e5bSBarry Smith dd2->lj = dd->lj; 1214644e2e5bSBarry Smith 12156e87535bSJed Brown da2->leveldown = da->leveldown + 1; 12166e87535bSJed Brown da2->levelup = da->levelup; 12178865f1eaSKarl Rupp 12186e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 12196e87535bSJed Brown 1220ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 12216636e97aSMatthew G Knepley if (da->coordinates) { 1222ddcf8b74SDave May DM cdaf,cdac; 1223ddcf8b74SDave May Vec coordsc,coordsf; 12246dbf9973SLawrence Mitchell Mat inject; 12256dbf9973SLawrence Mitchell VecScatter vscat; 1226ddcf8b74SDave May 12276636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); 12286636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); 12296636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); 1230b61d3410SDave May /* force creation of the coordinate vector */ 1231b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 12326636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 1233ddcf8b74SDave May 1234e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 12356dbf9973SLawrence Mitchell ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr); 12366dbf9973SLawrence Mitchell ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12376dbf9973SLawrence Mitchell ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12386dbf9973SLawrence Mitchell ierr = MatDestroy(&inject);CHKERRQ(ierr); 1239ddcf8b74SDave May } 1240f98405f7SJed Brown 1241f98405f7SJed Brown for (i=0; i<da->bs; i++) { 1242f98405f7SJed Brown const char *fieldname; 1243f98405f7SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1244f98405f7SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1245f98405f7SJed Brown } 1246f98405f7SJed Brown 124747c6ae99SBarry Smith *daref = da2; 124847c6ae99SBarry Smith PetscFunctionReturn(0); 124947c6ae99SBarry Smith } 125047c6ae99SBarry Smith 12517087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 125247c6ae99SBarry Smith { 125347c6ae99SBarry Smith PetscErrorCode ierr; 125447c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 125547c6ae99SBarry Smith 125647c6ae99SBarry Smith PetscFunctionBegin; 125747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1258ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 125947c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 126047c6ae99SBarry Smith PetscValidPointer(daf,3); 126147c6ae99SBarry Smith 1262aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 1263dcca6d9dSJed Brown ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr); 126447c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 1265aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 126647c6ae99SBarry Smith } 126747c6ae99SBarry Smith n = nlevels; 1268c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 126947c6ae99SBarry Smith n = nlevels; 1270c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 127147c6ae99SBarry Smith n = nlevels; 1272c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 127347c6ae99SBarry Smith 1274aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1275ce94432eSBarry Smith ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 127647c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1277aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1278ce94432eSBarry Smith ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 127947c6ae99SBarry Smith } 128047c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 128147c6ae99SBarry Smith PetscFunctionReturn(0); 128247c6ae99SBarry Smith } 128347c6ae99SBarry Smith 12847087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 128547c6ae99SBarry Smith { 128647c6ae99SBarry Smith PetscErrorCode ierr; 128747c6ae99SBarry Smith PetscInt i; 128847c6ae99SBarry Smith 128947c6ae99SBarry Smith PetscFunctionBegin; 129047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1291ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 129247c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 129347c6ae99SBarry Smith PetscValidPointer(dac,3); 1294ce94432eSBarry Smith ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 129547c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1296ce94432eSBarry Smith ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 129747c6ae99SBarry Smith } 129847c6ae99SBarry Smith PetscFunctionReturn(0); 129947c6ae99SBarry Smith } 130062197512SBarry Smith 1301f8f0706aSKarl Rupp #include <petscgll.h> 1302f8f0706aSKarl Rupp 1303f8f0706aSKarl Rupp PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscGLL *gll) 130462197512SBarry Smith { 130562197512SBarry Smith PetscErrorCode ierr; 1306f8f0706aSKarl Rupp PetscInt i,j,n = gll->n,xs,xn,q; 130762197512SBarry Smith PetscScalar *xx; 130862197512SBarry Smith PetscReal h; 130962197512SBarry Smith Vec x; 131062197512SBarry Smith DM_DA *da = (DM_DA*)dm->data; 131162197512SBarry Smith 131262197512SBarry Smith PetscFunctionBegin; 131362197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 131462197512SBarry Smith ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 131562197512SBarry Smith q = (q-1)/(n-1); /* number of spectral elements */ 131662197512SBarry Smith h = 2.0/q; 131762197512SBarry Smith ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 131862197512SBarry Smith xs = xs/(n-1); 131962197512SBarry Smith xn = xn/(n-1); 132062197512SBarry Smith ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr); 132162197512SBarry Smith ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 132262197512SBarry Smith ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr); 132362197512SBarry Smith 132462197512SBarry Smith /* loop over local spectral elements */ 132562197512SBarry Smith for (j=xs; j<xs+xn; j++) { 132662197512SBarry Smith /* 132762197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 132862197512SBarry Smith */ 132962197512SBarry Smith for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) { 1330f8f0706aSKarl Rupp xx[j*(n-1) + i] = -1.0 + h*j + h*(gll->nodes[i]+1.0)/2.; 133162197512SBarry Smith } 133262197512SBarry Smith } 133362197512SBarry Smith ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr); 133462197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic"); 133562197512SBarry Smith PetscFunctionReturn(0); 133662197512SBarry Smith } 133762197512SBarry Smith 13381fd49c25SBarry Smith /*@ 133962197512SBarry Smith 134062197512SBarry Smith DMDASetGLLCoordinates - Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points 134162197512SBarry Smith 134262197512SBarry Smith Collective on DM 134362197512SBarry Smith 134462197512SBarry Smith Input Parameters: 134562197512SBarry Smith + da - the DMDA object 1346f8f0706aSKarl Rupp - gll - the GLL object 134762197512SBarry Smith 134895452b02SPatrick Sanan Notes: 134995452b02SPatrick Sanan the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 135062197512SBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is 135162197512SBarry Smith periodic or not. 135262197512SBarry Smith 1353edc382c3SSatish Balay Level: advanced 1354edc382c3SSatish Balay 1355f8f0706aSKarl Rupp .seealso: DMDACreate(), PetscGLLCreate(), DMGetCoordinates() 135662197512SBarry Smith @*/ 1357f8f0706aSKarl Rupp PetscErrorCode DMDASetGLLCoordinates(DM da,PetscGLL *gll) 135862197512SBarry Smith { 135962197512SBarry Smith PetscErrorCode ierr; 136062197512SBarry Smith 136162197512SBarry Smith PetscFunctionBegin; 136262197512SBarry Smith if (da->dim == 1) { 1363f8f0706aSKarl Rupp ierr = DMDASetGLLCoordinates_1d(da,gll);CHKERRQ(ierr); 136462197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d"); 136562197512SBarry Smith PetscFunctionReturn(0); 136662197512SBarry Smith } 13677c3cd84eSPatrick Sanan 13687c3cd84eSPatrick Sanan PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1,DM dm2,PetscBool *compatible,PetscBool *set) 13697c3cd84eSPatrick Sanan { 13707c3cd84eSPatrick Sanan PetscErrorCode ierr; 13717c3cd84eSPatrick Sanan DM_DA *dd1 = (DM_DA*)da1->data,*dd2; 13727c3cd84eSPatrick Sanan DM da2; 13737c3cd84eSPatrick Sanan DMType dmtype2; 13747c3cd84eSPatrick Sanan PetscBool isda,compatibleLocal; 13757c3cd84eSPatrick Sanan PetscInt i; 13767c3cd84eSPatrick Sanan 13777c3cd84eSPatrick Sanan PetscFunctionBegin; 13787c3cd84eSPatrick Sanan if (!da1->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da1),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on first DM before DMGetCompatibility()"); 13797c3cd84eSPatrick Sanan ierr = DMGetType(dm2,&dmtype2);CHKERRQ(ierr); 13807c3cd84eSPatrick Sanan ierr = PetscStrcmp(dmtype2,DMDA,&isda);CHKERRQ(ierr); 13817c3cd84eSPatrick Sanan if (isda) { 13827c3cd84eSPatrick Sanan da2 = dm2; 13837c3cd84eSPatrick Sanan dd2 = (DM_DA*)da2->data; 13847c3cd84eSPatrick Sanan if (!da2->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da2),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on second DM before DMGetCompatibility()"); 1385110623a0SKarl Rupp compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1386c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1387110623a0SKarl Rupp /* Global size ranks Boundary type */ 1388c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1389c790a739SKarl Rupp if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1390c790a739SKarl Rupp if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 13917c3cd84eSPatrick Sanan if (compatibleLocal) { 13927c3cd84eSPatrick Sanan for (i=0; i<dd1->m; ++i) { 1393c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ 13947c3cd84eSPatrick Sanan } 13957c3cd84eSPatrick Sanan } 13967c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 1) { 13977c3cd84eSPatrick Sanan for (i=0; i<dd1->n; ++i) { 1398c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 13997c3cd84eSPatrick Sanan } 14007c3cd84eSPatrick Sanan } 14017c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 2) { 14027c3cd84eSPatrick Sanan for (i=0; i<dd1->p; ++i) { 1403c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 14047c3cd84eSPatrick Sanan } 14057c3cd84eSPatrick Sanan } 14067c3cd84eSPatrick Sanan *compatible = compatibleLocal; 14077c3cd84eSPatrick Sanan *set = PETSC_TRUE; 14087c3cd84eSPatrick Sanan } else { 14097c3cd84eSPatrick Sanan /* Decline to determine compatibility with other DM types */ 14107c3cd84eSPatrick Sanan *set = PETSC_FALSE; 14117c3cd84eSPatrick Sanan } 14127c3cd84eSPatrick Sanan PetscFunctionReturn(0); 14137c3cd84eSPatrick Sanan } 1414