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 6d083f849SBarry Smith Logically Collective on da 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); 30*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 31*7a8be351SBarry Smith PetscCheck(M >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive"); 32*7a8be351SBarry Smith PetscCheck(N >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive"); 33*7a8be351SBarry Smith PetscCheck(P >= 0,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 44d083f849SBarry Smith Logically Collective on da 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); 66*7a8be351SBarry Smith PetscCheck(!da->setupcalled,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; 72ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRMPI(ierr); 73e3f3e4b6SBarry Smith if ((dd->m > 0) && (dd->n < 0)) { 74e3f3e4b6SBarry Smith dd->n = size/dd->m; 75*7a8be351SBarry Smith PetscCheck(dd->n*dd->m == size,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; 79*7a8be351SBarry Smith PetscCheck(dd->n*dd->m == size,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 90d8d19677SJose E. Roman Input Parameters: 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 96bff4a2f0SMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType 9747c6ae99SBarry Smith @*/ 98bff4a2f0SMatthew G. Knepley PetscErrorCode DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz) 9947c6ae99SBarry Smith { 10047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 10147c6ae99SBarry Smith 10247c6ae99SBarry Smith PetscFunctionBegin; 103a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 1045a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bx,2); 1055a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,by,3); 1065a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bz,4); 107*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 1081321219cSEthan Coon dd->bx = bx; 1091321219cSEthan Coon dd->by = by; 1101321219cSEthan Coon dd->bz = bz; 11147c6ae99SBarry Smith PetscFunctionReturn(0); 11247c6ae99SBarry Smith } 11347c6ae99SBarry Smith 11447c6ae99SBarry Smith /*@ 115aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith Not collective 11847c6ae99SBarry Smith 11959f3ab6dSMatthew G. Knepley Input Parameters: 120aa219208SBarry Smith + da - The DMDA 12147c6ae99SBarry Smith - dof - Number of degrees of freedom 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith Level: intermediate 12447c6ae99SBarry Smith 125fb6725baSMatthew G. Knepley .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA 12647c6ae99SBarry Smith @*/ 12754cfb0beSLisandro Dalcin PetscErrorCode DMDASetDof(DM da, PetscInt dof) 12847c6ae99SBarry Smith { 12947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith PetscFunctionBegin; 132a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 13354cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da,dof,2); 134*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 13547c6ae99SBarry Smith dd->w = dof; 1361411c6eeSJed Brown da->bs = dof; 13747c6ae99SBarry Smith PetscFunctionReturn(0); 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith 140fb6725baSMatthew G. Knepley /*@ 141fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 142fb6725baSMatthew G. Knepley 143fb6725baSMatthew G. Knepley Not collective 144fb6725baSMatthew G. Knepley 145fb6725baSMatthew G. Knepley Input Parameter: 146fb6725baSMatthew G. Knepley . da - The DMDA 147fb6725baSMatthew G. Knepley 148fb6725baSMatthew G. Knepley Output Parameter: 149fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom 150fb6725baSMatthew G. Knepley 151fb6725baSMatthew G. Knepley Level: intermediate 152fb6725baSMatthew G. Knepley 153fb6725baSMatthew G. Knepley .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA 154fb6725baSMatthew G. Knepley @*/ 155fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 156fb6725baSMatthew G. Knepley { 157fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 158fb6725baSMatthew G. Knepley 159fb6725baSMatthew G. Knepley PetscFunctionBegin; 160a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 161fb6725baSMatthew G. Knepley PetscValidPointer(dof,2); 162fb6725baSMatthew G. Knepley *dof = dd->w; 163fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 164fb6725baSMatthew G. Knepley } 165fb6725baSMatthew G. Knepley 1667ddda789SPeter Brune /*@ 1677ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1687ddda789SPeter Brune 1697ddda789SPeter Brune Not collective 1707ddda789SPeter Brune 171f899ff85SJose E. Roman Input Parameter: 1727ddda789SPeter Brune . da - The DMDA 1737ddda789SPeter Brune 1747ddda789SPeter Brune Output Parameters: 1757ddda789SPeter Brune + x - Overlap in the x direction 1767ddda789SPeter Brune . y - Overlap in the y direction 1777ddda789SPeter Brune - z - Overlap in the z direction 1787ddda789SPeter Brune 1797ddda789SPeter Brune Level: intermediate 1807ddda789SPeter Brune 1817ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA 1827ddda789SPeter Brune @*/ 1837ddda789SPeter Brune PetscErrorCode DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z) 1847ddda789SPeter Brune { 1857ddda789SPeter Brune DM_DA *dd = (DM_DA*)da->data; 1867ddda789SPeter Brune 1877ddda789SPeter Brune PetscFunctionBegin; 188a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 1897ddda789SPeter Brune if (x) *x = dd->xol; 1907ddda789SPeter Brune if (y) *y = dd->yol; 1917ddda789SPeter Brune if (z) *z = dd->zol; 1927ddda789SPeter Brune PetscFunctionReturn(0); 1937ddda789SPeter Brune } 1947ddda789SPeter Brune 19588661749SPeter Brune /*@ 19688661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 19788661749SPeter Brune 19888661749SPeter Brune Not collective 19988661749SPeter Brune 2007ddda789SPeter Brune Input Parameters: 20188661749SPeter Brune + da - The DMDA 2027ddda789SPeter Brune . x - Overlap in the x direction 2037ddda789SPeter Brune . y - Overlap in the y direction 2047ddda789SPeter Brune - z - Overlap in the z direction 20588661749SPeter Brune 20688661749SPeter Brune Level: intermediate 20788661749SPeter Brune 2087ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA 20988661749SPeter Brune @*/ 2107ddda789SPeter Brune PetscErrorCode DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z) 21188661749SPeter Brune { 21288661749SPeter Brune DM_DA *dd = (DM_DA*)da->data; 21388661749SPeter Brune 21488661749SPeter Brune PetscFunctionBegin; 215a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2167ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,x,2); 2177ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,y,3); 2187ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,z,4); 2197ddda789SPeter Brune dd->xol = x; 2207ddda789SPeter Brune dd->yol = y; 2217ddda789SPeter Brune dd->zol = z; 22288661749SPeter Brune PetscFunctionReturn(0); 22388661749SPeter Brune } 22488661749SPeter Brune 2253e7870d2SPeter Brune /*@ 2263e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2273e7870d2SPeter Brune 2283e7870d2SPeter Brune Not collective 2293e7870d2SPeter Brune 2303e7870d2SPeter Brune Input Parameters: 2313e7870d2SPeter Brune . da - The DMDA 2323e7870d2SPeter Brune 2333e7870d2SPeter Brune Output Parameters: 234a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition 2353e7870d2SPeter Brune 2363e7870d2SPeter Brune Level: intermediate 2373e7870d2SPeter Brune 2383e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA 2393e7870d2SPeter Brune @*/ 2403e7870d2SPeter Brune PetscErrorCode DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub) 2413e7870d2SPeter Brune { 2423e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2433e7870d2SPeter Brune 2443e7870d2SPeter Brune PetscFunctionBegin; 245a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2463e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2473e7870d2SPeter Brune PetscFunctionReturn(0); 2483e7870d2SPeter Brune } 2493e7870d2SPeter Brune 2503e7870d2SPeter Brune /*@ 2513e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2523e7870d2SPeter Brune 2533e7870d2SPeter Brune Not collective 2543e7870d2SPeter Brune 2553e7870d2SPeter Brune Input Parameters: 2563e7870d2SPeter Brune + da - The DMDA 2573e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2583e7870d2SPeter Brune 2593e7870d2SPeter Brune Level: intermediate 2603e7870d2SPeter Brune 2613e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA 2623e7870d2SPeter Brune @*/ 2633e7870d2SPeter Brune PetscErrorCode DMDASetNumLocalSubDomains(DM da,PetscInt Nsub) 2643e7870d2SPeter Brune { 2653e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2663e7870d2SPeter Brune 2673e7870d2SPeter Brune PetscFunctionBegin; 268a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2693e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da,Nsub,2); 2703e7870d2SPeter Brune dd->Nsub = Nsub; 2713e7870d2SPeter Brune PetscFunctionReturn(0); 2723e7870d2SPeter Brune } 2733e7870d2SPeter Brune 274d886c4f4SPeter Brune /*@ 275d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 276d886c4f4SPeter Brune 277d886c4f4SPeter Brune Collective on DA 278d886c4f4SPeter Brune 279d8d19677SJose E. Roman Input Parameters: 280d886c4f4SPeter Brune + da - The DMDA 281d886c4f4SPeter Brune . xo - The offset in the x direction 282d886c4f4SPeter Brune . yo - The offset in the y direction 283d886c4f4SPeter Brune - zo - The offset in the z direction 284d886c4f4SPeter Brune 285d886c4f4SPeter Brune Level: intermediate 286d886c4f4SPeter Brune 28795452b02SPatrick Sanan Notes: 28895452b02SPatrick Sanan This is used primarily to overlap a computation on a local DA with that on a global DA without 289d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 290d886c4f4SPeter Brune 291d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 292d886c4f4SPeter Brune @*/ 29395c13181SPeter Brune PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 294d886c4f4SPeter Brune { 29595c13181SPeter Brune PetscErrorCode ierr; 296d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 297d886c4f4SPeter Brune 298d886c4f4SPeter Brune PetscFunctionBegin; 299a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 300d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da,xo,2); 30195c13181SPeter Brune PetscValidLogicalCollectiveInt(da,yo,3); 30295c13181SPeter Brune PetscValidLogicalCollectiveInt(da,zo,4); 30395c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Mo,5); 30495c13181SPeter Brune PetscValidLogicalCollectiveInt(da,No,6); 30595c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Po,7); 306d886c4f4SPeter Brune dd->xo = xo; 307d886c4f4SPeter Brune dd->yo = yo; 308d886c4f4SPeter Brune dd->zo = zo; 30995c13181SPeter Brune dd->Mo = Mo; 31095c13181SPeter Brune dd->No = No; 31195c13181SPeter Brune dd->Po = Po; 31295c13181SPeter Brune 31395c13181SPeter Brune if (da->coordinateDM) { 31495c13181SPeter Brune ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 31595c13181SPeter Brune } 316d886c4f4SPeter Brune PetscFunctionReturn(0); 317d886c4f4SPeter Brune } 318d886c4f4SPeter Brune 319d886c4f4SPeter Brune /*@ 320d886c4f4SPeter Brune DMDAGetOffset - Gets the index offset of the DA. 321d886c4f4SPeter Brune 322d886c4f4SPeter Brune Not collective 323d886c4f4SPeter Brune 324d886c4f4SPeter Brune Input Parameter: 325d886c4f4SPeter Brune . da - The DMDA 326d886c4f4SPeter Brune 327d886c4f4SPeter Brune Output Parameters: 328d886c4f4SPeter Brune + xo - The offset in the x direction 329d886c4f4SPeter Brune . yo - The offset in the y direction 33095c13181SPeter Brune . zo - The offset in the z direction 33195c13181SPeter Brune . Mo - The global size in the x direction 33295c13181SPeter Brune . No - The global size in the y direction 33395c13181SPeter Brune - Po - The global size in the z direction 334d886c4f4SPeter Brune 335d886c4f4SPeter Brune Level: intermediate 336d886c4f4SPeter Brune 337d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray() 338d886c4f4SPeter Brune @*/ 33995c13181SPeter Brune PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 340d886c4f4SPeter Brune { 341d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 342d886c4f4SPeter Brune 343d886c4f4SPeter Brune PetscFunctionBegin; 344a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 345d886c4f4SPeter Brune if (xo) *xo = dd->xo; 346d886c4f4SPeter Brune if (yo) *yo = dd->yo; 347d886c4f4SPeter Brune if (zo) *zo = dd->zo; 34895c13181SPeter Brune if (Mo) *Mo = dd->Mo; 34995c13181SPeter Brune if (No) *No = dd->No; 35095c13181SPeter Brune if (Po) *Po = dd->Po; 351d886c4f4SPeter Brune PetscFunctionReturn(0); 352d886c4f4SPeter Brune } 353d886c4f4SPeter Brune 35440234c92SPeter Brune /*@ 35540234c92SPeter Brune DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM. 35640234c92SPeter Brune 35740234c92SPeter Brune Not collective 35840234c92SPeter Brune 35940234c92SPeter Brune Input Parameter: 36040234c92SPeter Brune . da - The DMDA 36140234c92SPeter Brune 36240234c92SPeter Brune Output Parameters: 36340234c92SPeter Brune + xs - The start of the region in x 36440234c92SPeter Brune . ys - The start of the region in y 36540234c92SPeter Brune . zs - The start of the region in z 36640234c92SPeter Brune . xs - The size of the region in x 36740234c92SPeter Brune . ys - The size of the region in y 368a2b725a8SWilliam Gropp - zs - The size of the region in z 36940234c92SPeter Brune 37040234c92SPeter Brune Level: intermediate 37140234c92SPeter Brune 37240234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 37340234c92SPeter Brune @*/ 37440234c92SPeter Brune PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 37540234c92SPeter Brune { 37640234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 37740234c92SPeter Brune 37840234c92SPeter Brune PetscFunctionBegin; 379a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 38040234c92SPeter Brune if (xs) *xs = dd->nonxs; 38140234c92SPeter Brune if (ys) *ys = dd->nonys; 38240234c92SPeter Brune if (zs) *zs = dd->nonzs; 38340234c92SPeter Brune if (xm) *xm = dd->nonxm; 38440234c92SPeter Brune if (ym) *ym = dd->nonym; 38540234c92SPeter Brune if (zm) *zm = dd->nonzm; 38640234c92SPeter Brune PetscFunctionReturn(0); 38740234c92SPeter Brune } 38840234c92SPeter Brune 38940234c92SPeter Brune /*@ 39040234c92SPeter Brune DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM. 39140234c92SPeter Brune 39240234c92SPeter Brune Collective on DA 39340234c92SPeter Brune 394d8d19677SJose E. Roman Input Parameters: 39540234c92SPeter Brune + da - The DMDA 39640234c92SPeter Brune . xs - The start of the region in x 39740234c92SPeter Brune . ys - The start of the region in y 39840234c92SPeter Brune . zs - The start of the region in z 39940234c92SPeter Brune . xs - The size of the region in x 40040234c92SPeter Brune . ys - The size of the region in y 401a2b725a8SWilliam Gropp - zs - The size of the region in z 40240234c92SPeter Brune 40340234c92SPeter Brune Level: intermediate 40440234c92SPeter Brune 40540234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 40640234c92SPeter Brune @*/ 40740234c92SPeter Brune PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 40840234c92SPeter Brune { 40940234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 41040234c92SPeter Brune 41140234c92SPeter Brune PetscFunctionBegin; 412a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 41340234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xs,2); 41440234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ys,3); 41540234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zs,4); 41640234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xm,5); 41740234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ym,6); 41840234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zm,7); 41940234c92SPeter Brune dd->nonxs = xs; 42040234c92SPeter Brune dd->nonys = ys; 42140234c92SPeter Brune dd->nonzs = zs; 42240234c92SPeter Brune dd->nonxm = xm; 42340234c92SPeter Brune dd->nonym = ym; 42440234c92SPeter Brune dd->nonzm = zm; 42540234c92SPeter Brune 42640234c92SPeter Brune PetscFunctionReturn(0); 42740234c92SPeter Brune } 42888661749SPeter Brune 42947c6ae99SBarry Smith /*@ 430aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 43147c6ae99SBarry Smith 432d083f849SBarry Smith Logically Collective on da 43347c6ae99SBarry Smith 434d8d19677SJose E. Roman Input Parameters: 435aa219208SBarry Smith + da - The DMDA 436aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 43747c6ae99SBarry Smith 43847c6ae99SBarry Smith Level: intermediate 43947c6ae99SBarry Smith 44096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 44147c6ae99SBarry Smith @*/ 4427087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 44347c6ae99SBarry Smith { 44447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 44547c6ae99SBarry Smith 44647c6ae99SBarry Smith PetscFunctionBegin; 447a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 44847c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 449*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 45047c6ae99SBarry Smith dd->stencil_type = stype; 45147c6ae99SBarry Smith PetscFunctionReturn(0); 45247c6ae99SBarry Smith } 45347c6ae99SBarry Smith 454fb6725baSMatthew G. Knepley /*@ 455fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 456fb6725baSMatthew G. Knepley 457fb6725baSMatthew G. Knepley Not collective 458fb6725baSMatthew G. Knepley 459fb6725baSMatthew G. Knepley Input Parameter: 460fb6725baSMatthew G. Knepley . da - The DMDA 461fb6725baSMatthew G. Knepley 462fb6725baSMatthew G. Knepley Output Parameter: 463fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 464fb6725baSMatthew G. Knepley 465fb6725baSMatthew G. Knepley Level: intermediate 466fb6725baSMatthew G. Knepley 467fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 468fb6725baSMatthew G. Knepley @*/ 469fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 470fb6725baSMatthew G. Knepley { 471fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA*)da->data; 472fb6725baSMatthew G. Knepley 473fb6725baSMatthew G. Knepley PetscFunctionBegin; 474a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 475fb6725baSMatthew G. Knepley PetscValidPointer(stype,2); 476fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 477fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 478fb6725baSMatthew G. Knepley } 479fb6725baSMatthew G. Knepley 48047c6ae99SBarry Smith /*@ 481aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 48247c6ae99SBarry Smith 483d083f849SBarry Smith Logically Collective on da 48447c6ae99SBarry Smith 485d8d19677SJose E. Roman Input Parameters: 486aa219208SBarry Smith + da - The DMDA 48747c6ae99SBarry Smith - width - The stencil width 48847c6ae99SBarry Smith 48947c6ae99SBarry Smith Level: intermediate 49047c6ae99SBarry Smith 49196e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 49247c6ae99SBarry Smith @*/ 4937087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 49447c6ae99SBarry Smith { 49547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 49647c6ae99SBarry Smith 49747c6ae99SBarry Smith PetscFunctionBegin; 498a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 49947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 500*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 50147c6ae99SBarry Smith dd->s = width; 50247c6ae99SBarry Smith PetscFunctionReturn(0); 50347c6ae99SBarry Smith } 50447c6ae99SBarry Smith 505fb6725baSMatthew G. Knepley /*@ 506fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 507fb6725baSMatthew G. Knepley 508fb6725baSMatthew G. Knepley Not collective 509fb6725baSMatthew G. Knepley 510fb6725baSMatthew G. Knepley Input Parameter: 511fb6725baSMatthew G. Knepley . da - The DMDA 512fb6725baSMatthew G. Knepley 513fb6725baSMatthew G. Knepley Output Parameter: 514fb6725baSMatthew G. Knepley . width - The stencil width 515fb6725baSMatthew G. Knepley 516fb6725baSMatthew G. Knepley Level: intermediate 517fb6725baSMatthew G. Knepley 518fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 519fb6725baSMatthew G. Knepley @*/ 520fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 521fb6725baSMatthew G. Knepley { 522fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 523fb6725baSMatthew G. Knepley 524fb6725baSMatthew G. Knepley PetscFunctionBegin; 525a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 526fb6725baSMatthew G. Knepley PetscValidPointer(width,2); 527fb6725baSMatthew G. Knepley *width = dd->s; 528fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 529fb6725baSMatthew G. Knepley } 530fb6725baSMatthew G. Knepley 531aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 53247c6ae99SBarry Smith { 53347c6ae99SBarry Smith PetscInt i,sum; 53447c6ae99SBarry Smith 53547c6ae99SBarry Smith PetscFunctionBegin; 536*7a8be351SBarry Smith PetscCheck(M >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 53747c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 538*7a8be351SBarry Smith PetscCheck(sum == M,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 53947c6ae99SBarry Smith PetscFunctionReturn(0); 54047c6ae99SBarry Smith } 54147c6ae99SBarry Smith 54247c6ae99SBarry Smith /*@ 543aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 54447c6ae99SBarry Smith 545d083f849SBarry Smith Logically Collective on da 54647c6ae99SBarry Smith 547d8d19677SJose E. Roman Input Parameters: 548aa219208SBarry Smith + da - The DMDA 5490298fd71SBarry 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 5500298fd71SBarry 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 5510298fd71SBarry 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. 55247c6ae99SBarry Smith 55347c6ae99SBarry Smith Level: intermediate 55447c6ae99SBarry Smith 555e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 556e3f3e4b6SBarry Smith 55796e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 55847c6ae99SBarry Smith @*/ 5597087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 56047c6ae99SBarry Smith { 56147c6ae99SBarry Smith PetscErrorCode ierr; 56247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith PetscFunctionBegin; 565a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 566*7a8be351SBarry Smith PetscCheck(!da->setupcalled,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 56747c6ae99SBarry Smith if (lx) { 568*7a8be351SBarry Smith PetscCheck(dd->m >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 569aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 57047c6ae99SBarry Smith if (!dd->lx) { 571785e854fSJed Brown ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr); 57247c6ae99SBarry Smith } 573580bdb30SBarry Smith ierr = PetscArraycpy(dd->lx, lx, dd->m);CHKERRQ(ierr); 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith if (ly) { 576*7a8be351SBarry Smith PetscCheck(dd->n >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 577aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 57847c6ae99SBarry Smith if (!dd->ly) { 579785e854fSJed Brown ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr); 58047c6ae99SBarry Smith } 581580bdb30SBarry Smith ierr = PetscArraycpy(dd->ly, ly, dd->n);CHKERRQ(ierr); 58247c6ae99SBarry Smith } 58347c6ae99SBarry Smith if (lz) { 584*7a8be351SBarry Smith PetscCheck(dd->p >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 585aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 58647c6ae99SBarry Smith if (!dd->lz) { 587785e854fSJed Brown ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr); 58847c6ae99SBarry Smith } 589580bdb30SBarry Smith ierr = PetscArraycpy(dd->lz, lz, dd->p);CHKERRQ(ierr); 59047c6ae99SBarry Smith } 59147c6ae99SBarry Smith PetscFunctionReturn(0); 59247c6ae99SBarry Smith } 59347c6ae99SBarry Smith 59447c6ae99SBarry Smith /*@ 595aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 596e727c939SJed Brown returned by DMCreateInterpolation() 59747c6ae99SBarry Smith 598d083f849SBarry Smith Logically Collective on da 59947c6ae99SBarry Smith 600d8d19677SJose E. Roman Input Parameters: 60147c6ae99SBarry Smith + da - initial distributed array 602a2b725a8SWilliam Gropp - ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 60347c6ae99SBarry Smith 60447c6ae99SBarry Smith Level: intermediate 60547c6ae99SBarry Smith 60695452b02SPatrick Sanan Notes: 60795452b02SPatrick Sanan you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 60847c6ae99SBarry Smith 60996e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 61047c6ae99SBarry Smith @*/ 6117087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 61247c6ae99SBarry Smith { 61347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 61447c6ae99SBarry Smith 61547c6ae99SBarry Smith PetscFunctionBegin; 616a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 61747c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 61847c6ae99SBarry Smith dd->interptype = ctype; 61947c6ae99SBarry Smith PetscFunctionReturn(0); 62047c6ae99SBarry Smith } 62147c6ae99SBarry Smith 6222dde6fd4SLisandro Dalcin /*@ 6232dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 624e727c939SJed Brown used by DMCreateInterpolation() 6252dde6fd4SLisandro Dalcin 6262dde6fd4SLisandro Dalcin Not Collective 6272dde6fd4SLisandro Dalcin 6282dde6fd4SLisandro Dalcin Input Parameter: 6292dde6fd4SLisandro Dalcin . da - distributed array 6302dde6fd4SLisandro Dalcin 6312dde6fd4SLisandro Dalcin Output Parameter: 6322dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 6332dde6fd4SLisandro Dalcin 6342dde6fd4SLisandro Dalcin Level: intermediate 6352dde6fd4SLisandro Dalcin 636e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 6372dde6fd4SLisandro Dalcin @*/ 6382dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 6392dde6fd4SLisandro Dalcin { 6402dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 6412dde6fd4SLisandro Dalcin 6422dde6fd4SLisandro Dalcin PetscFunctionBegin; 643a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 6442dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 6452dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6462dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 6472dde6fd4SLisandro Dalcin } 64847c6ae99SBarry Smith 6496a119db4SBarry Smith /*@C 650aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 65147c6ae99SBarry Smith processes neighbors. 65247c6ae99SBarry Smith 65347c6ae99SBarry Smith Not Collective 65447c6ae99SBarry Smith 65547c6ae99SBarry Smith Input Parameter: 656aa219208SBarry Smith . da - the DMDA object 65747c6ae99SBarry Smith 65847c6ae99SBarry Smith Output Parameters: 65947c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 66047c6ae99SBarry Smith this process itself is in the list 66147c6ae99SBarry Smith 66295452b02SPatrick Sanan Notes: 66395452b02SPatrick Sanan In 2d the array is of length 9, in 3d of length 27 66447c6ae99SBarry Smith Not supported in 1d 665aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 66647c6ae99SBarry Smith 66795452b02SPatrick Sanan Fortran Notes: 66895452b02SPatrick Sanan In fortran you must pass in an array of the appropriate length. 66947c6ae99SBarry Smith 67047c6ae99SBarry Smith Level: intermediate 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith @*/ 6737087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 67447c6ae99SBarry Smith { 67547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 6765fd66863SKarl Rupp 67747c6ae99SBarry Smith PetscFunctionBegin; 678a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 67947c6ae99SBarry Smith *ranks = dd->neighbors; 68047c6ae99SBarry Smith PetscFunctionReturn(0); 68147c6ae99SBarry Smith } 68247c6ae99SBarry Smith 68347c6ae99SBarry Smith /*@C 684aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith Not Collective 68747c6ae99SBarry Smith 68847c6ae99SBarry Smith Input Parameter: 689aa219208SBarry Smith . da - the DMDA object 69047c6ae99SBarry Smith 691d8d19677SJose E. Roman Output Parameters: 69247c6ae99SBarry Smith + lx - ownership along x direction (optional) 69347c6ae99SBarry Smith . ly - ownership along y direction (optional) 69447c6ae99SBarry Smith - lz - ownership along z direction (optional) 69547c6ae99SBarry Smith 69647c6ae99SBarry Smith Level: intermediate 69747c6ae99SBarry Smith 698aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 69947c6ae99SBarry Smith 70047c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 701aa219208SBarry Smith eighth arguments from DMDAGetInfo() 70247c6ae99SBarry Smith 70347c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 704aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 70547c6ae99SBarry Smith 706e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 707e3f3e4b6SBarry Smith 708aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 70947c6ae99SBarry Smith @*/ 7107087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 71147c6ae99SBarry Smith { 71247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith PetscFunctionBegin; 715a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 71647c6ae99SBarry Smith if (lx) *lx = dd->lx; 71747c6ae99SBarry Smith if (ly) *ly = dd->ly; 71847c6ae99SBarry Smith if (lz) *lz = dd->lz; 71947c6ae99SBarry Smith PetscFunctionReturn(0); 72047c6ae99SBarry Smith } 72147c6ae99SBarry Smith 72247c6ae99SBarry Smith /*@ 723aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 72447c6ae99SBarry Smith 725d083f849SBarry Smith Logically Collective on da 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: 73448eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction 73548eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction 73648eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction 73748eeb7c8SBarry Smith - -da_refine <n> - refine the DMDA object n times when it is created. 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith Level: intermediate 74047c6ae99SBarry Smith 74195452b02SPatrick Sanan Notes: 74295452b02SPatrick Sanan Pass PETSC_IGNORE to leave a value unchanged 74347c6ae99SBarry Smith 744aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 74547c6ae99SBarry Smith @*/ 7467087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 74747c6ae99SBarry Smith { 74847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith PetscFunctionBegin; 751a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 75247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 75347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 75447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 75747c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 75847c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 75947c6ae99SBarry Smith PetscFunctionReturn(0); 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith 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 77795452b02SPatrick Sanan Notes: 77895452b02SPatrick Sanan Pass NULL for values you do not need 77947c6ae99SBarry Smith 780aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 78147c6ae99SBarry Smith @*/ 7827087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 78347c6ae99SBarry Smith { 78447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith PetscFunctionBegin; 787a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 78847c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 78947c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 79047c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 79147c6ae99SBarry Smith PetscFunctionReturn(0); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 79447c6ae99SBarry Smith /*@C 795aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 79647c6ae99SBarry Smith 797d083f849SBarry Smith Logically Collective on da 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith Input Parameters: 800aa219208SBarry Smith + da - the DMDA object 801aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith Level: developer 80447c6ae99SBarry Smith 80595452b02SPatrick Sanan Notes: 80695452b02SPatrick Sanan 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 8091fd49c25SBarry Smith Not supported from Fortran 8101fd49c25SBarry Smith 811950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 81247c6ae99SBarry Smith @*/ 813b412c318SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*)) 81447c6ae99SBarry Smith { 81547c6ae99SBarry Smith PetscFunctionBegin; 816a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 81725296bd5SBarry Smith da->ops->creatematrix = f; 81847c6ae99SBarry Smith PetscFunctionReturn(0); 81947c6ae99SBarry Smith } 82047c6ae99SBarry Smith 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; 833*7a8be351SBarry Smith PetscCheck(ratio >= 1,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 83447c6ae99SBarry Smith if (ratio == 1) { 835580bdb30SBarry Smith ierr = PetscArraycpy(lf,lc,m);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 8632be375d4SJed Brown /* 8642be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 8652be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 8662be375d4SJed Brown 8672be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 8682be375d4SJed Brown */ 8692be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 8702be375d4SJed Brown { 8712be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 8722be375d4SJed Brown PetscErrorCode ierr; 8732be375d4SJed Brown 8742be375d4SJed Brown PetscFunctionBegin; 875*7a8be351SBarry Smith PetscCheck(ratio >= 1,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 8762be375d4SJed Brown if (ratio == 1) { 877580bdb30SBarry Smith ierr = PetscArraycpy(lc,lf,m);CHKERRQ(ierr); 8782be375d4SJed Brown PetscFunctionReturn(0); 8792be375d4SJed Brown } 8802be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 8812be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 8822be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 8832be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 8842be375d4SJed Brown if (i == m-1) lc[i] = want; 8852be375d4SJed Brown else { 8862be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 8872be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 8882be375d4SJed Brown * node is within one stencil width. */ 8892be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 8902be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 8912be375d4SJed Brown * fine node is within one stencil width. */ 8922be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 8932be375d4SJed Brown if (want < 0 || want > remaining 894ce94432eSBarry 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"); 8952be375d4SJed Brown } 8962be375d4SJed Brown lc[i] = want; 8972be375d4SJed Brown startc += lc[i]; 8982be375d4SJed Brown startf += lf[i]; 8992be375d4SJed Brown remaining -= lc[i]; 9002be375d4SJed Brown } 9012be375d4SJed Brown PetscFunctionReturn(0); 9022be375d4SJed Brown } 9032be375d4SJed Brown 9047087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 90547c6ae99SBarry Smith { 90647c6ae99SBarry Smith PetscErrorCode ierr; 907c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 9089a42bb27SBarry Smith DM da2; 90947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 91047c6ae99SBarry Smith 91147c6ae99SBarry Smith PetscFunctionBegin; 912a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 91347c6ae99SBarry Smith PetscValidPointer(daref,3); 91447c6ae99SBarry Smith 915c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 916bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 91747c6ae99SBarry Smith M = dd->refine_x*dd->M; 91847c6ae99SBarry Smith } else { 91947c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 92047c6ae99SBarry Smith } 921bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 922c73cfb54SMatthew G. Knepley if (dim > 1) { 92347c6ae99SBarry Smith N = dd->refine_y*dd->N; 92447c6ae99SBarry Smith } else { 9251860e6e9SBarry Smith N = 1; 9261860e6e9SBarry Smith } 9271860e6e9SBarry Smith } else { 92847c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 92947c6ae99SBarry Smith } 930bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 931c73cfb54SMatthew G. Knepley if (dim > 2) { 93247c6ae99SBarry Smith P = dd->refine_z*dd->P; 93347c6ae99SBarry Smith } else { 9341860e6e9SBarry Smith P = 1; 9351860e6e9SBarry Smith } 9361860e6e9SBarry Smith } else { 93747c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 93847c6ae99SBarry Smith } 939ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 940397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 941c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 942397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 943397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 944397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 945397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 946397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 947397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 948c73cfb54SMatthew G. Knepley if (dim == 3) { 94947c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 950dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 951bff4a2f0SMatthew 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); 952bff4a2f0SMatthew 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); 953bff4a2f0SMatthew 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); 954397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 95547c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 956c73cfb54SMatthew G. Knepley } else if (dim == 2) { 95747c6ae99SBarry Smith PetscInt *lx,*ly; 958dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 959bff4a2f0SMatthew 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); 960bff4a2f0SMatthew 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); 9610298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 96247c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 963c73cfb54SMatthew G. Knepley } else if (dim == 1) { 96447c6ae99SBarry Smith PetscInt *lx; 965785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 966bff4a2f0SMatthew 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); 9670298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 96847c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 96947c6ae99SBarry Smith } 97047c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 97147c6ae99SBarry Smith 97247c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 97325296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 97425296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 97547c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 97647c6ae99SBarry Smith dd2->interptype = dd->interptype; 97747c6ae99SBarry Smith 97847c6ae99SBarry Smith /* copy fill information if given */ 97947c6ae99SBarry Smith if (dd->dfill) { 980854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 981580bdb30SBarry Smith ierr = PetscArraycpy(dd2->dfill,dd->dfill,dd->dfill[dd->w]+dd->w+1);CHKERRQ(ierr); 98247c6ae99SBarry Smith } 98347c6ae99SBarry Smith if (dd->ofill) { 984854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 985580bdb30SBarry Smith ierr = PetscArraycpy(dd2->ofill,dd->ofill,dd->ofill[dd->w]+dd->w+1);CHKERRQ(ierr); 98647c6ae99SBarry Smith } 98747c6ae99SBarry Smith /* copy the refine information */ 988397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 989397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 990397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 99147c6ae99SBarry Smith 992897f7067SBarry Smith if (dd->refine_z_hier) { 993897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) { 994897f7067SBarry Smith dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; 995897f7067SBarry Smith } 996897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) { 997897f7067SBarry Smith dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; 998897f7067SBarry Smith } 999897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1000897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1001580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n);CHKERRQ(ierr); 1002897f7067SBarry Smith } 1003897f7067SBarry Smith if (dd->refine_y_hier) { 1004897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) { 1005897f7067SBarry Smith dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; 1006897f7067SBarry Smith } 1007897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) { 1008897f7067SBarry Smith dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; 1009897f7067SBarry Smith } 1010897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1011897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1012580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n);CHKERRQ(ierr); 1013897f7067SBarry Smith } 1014897f7067SBarry Smith if (dd->refine_x_hier) { 1015897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) { 1016897f7067SBarry Smith dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; 1017897f7067SBarry Smith } 1018897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) { 1019897f7067SBarry Smith dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; 1020897f7067SBarry Smith } 1021897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1022897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1023580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n);CHKERRQ(ierr); 1024897f7067SBarry Smith } 1025897f7067SBarry Smith 102647c6ae99SBarry Smith /* copy vector type information */ 1027c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 1028ddcf8b74SDave May 1029efd51863SBarry Smith dd2->lf = dd->lf; 1030efd51863SBarry Smith dd2->lj = dd->lj; 1031efd51863SBarry Smith 10326e87535bSJed Brown da2->leveldown = da->leveldown; 10336e87535bSJed Brown da2->levelup = da->levelup + 1; 10348865f1eaSKarl Rupp 10356e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 10366e87535bSJed Brown 1037ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10386636e97aSMatthew G Knepley if (da->coordinates) { 1039ddcf8b74SDave May DM cdaf,cdac; 1040ddcf8b74SDave May Vec coordsc,coordsf; 1041ddcf8b74SDave May Mat II; 1042ddcf8b74SDave May 10436636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 10446636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 10456636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 1046b61d3410SDave May /* force creation of the coordinate vector */ 1047b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 10486636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 10490298fd71SBarry Smith ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 1050ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 1051fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 1052ddcf8b74SDave May } 1053397b6216SJed Brown 1054f3141302SJed Brown for (i=0; i<da->bs; i++) { 1055f3141302SJed Brown const char *fieldname; 1056f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1057f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1058f3141302SJed Brown } 1059397b6216SJed Brown 106047c6ae99SBarry Smith *daref = da2; 106147c6ae99SBarry Smith PetscFunctionReturn(0); 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith 1064a5bc1bf3SBarry Smith PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm,DM *dmc) 106547c6ae99SBarry Smith { 106647c6ae99SBarry Smith PetscErrorCode ierr; 1067c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 1068a5bc1bf3SBarry Smith DM dmc2; 1069a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA*)dmf->data,*dd2; 107047c6ae99SBarry Smith 107147c6ae99SBarry Smith PetscFunctionBegin; 1072a5bc1bf3SBarry Smith PetscValidHeaderSpecificType(dmf,DM_CLASSID,1,DMDA); 1073a5bc1bf3SBarry Smith PetscValidPointer(dmc,3); 107447c6ae99SBarry Smith 1075a5bc1bf3SBarry Smith ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1076bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1077397b6216SJed Brown M = dd->M / dd->coarsen_x; 107847c6ae99SBarry Smith } else { 1079397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 108047c6ae99SBarry Smith } 1081bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1082c73cfb54SMatthew G. Knepley if (dim > 1) { 1083397b6216SJed Brown N = dd->N / dd->coarsen_y; 108447c6ae99SBarry Smith } else { 10851860e6e9SBarry Smith N = 1; 10861860e6e9SBarry Smith } 10871860e6e9SBarry Smith } else { 1088397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 108947c6ae99SBarry Smith } 1090bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1091c73cfb54SMatthew G. Knepley if (dim > 2) { 1092397b6216SJed Brown P = dd->P / dd->coarsen_z; 109347c6ae99SBarry Smith } else { 10941860e6e9SBarry Smith P = 1; 10951860e6e9SBarry Smith } 10961860e6e9SBarry Smith } else { 1097397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 109847c6ae99SBarry Smith } 1099a5bc1bf3SBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)dmf),&dmc2);CHKERRQ(ierr); 1100a5bc1bf3SBarry Smith ierr = DMSetOptionsPrefix(dmc2,((PetscObject)dmf)->prefix);CHKERRQ(ierr); 1101a5bc1bf3SBarry Smith ierr = DMSetDimension(dmc2,dim);CHKERRQ(ierr); 1102a5bc1bf3SBarry Smith ierr = DMDASetSizes(dmc2,M,N,P);CHKERRQ(ierr); 1103a5bc1bf3SBarry Smith ierr = DMDASetNumProcs(dmc2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 1104a5bc1bf3SBarry Smith ierr = DMDASetBoundaryType(dmc2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 1105a5bc1bf3SBarry Smith ierr = DMDASetDof(dmc2,dd->w);CHKERRQ(ierr); 1106a5bc1bf3SBarry Smith ierr = DMDASetStencilType(dmc2,dd->stencil_type);CHKERRQ(ierr); 1107a5bc1bf3SBarry Smith ierr = DMDASetStencilWidth(dmc2,dd->s);CHKERRQ(ierr); 1108c73cfb54SMatthew G. Knepley if (dim == 3) { 11092be375d4SJed Brown PetscInt *lx,*ly,*lz; 1110dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 1111a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1112a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 1113a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 1114a5bc1bf3SBarry Smith ierr = DMDASetOwnershipRanges(dmc2,lx,ly,lz);CHKERRQ(ierr); 11152be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 1116c73cfb54SMatthew G. Knepley } else if (dim == 2) { 11172be375d4SJed Brown PetscInt *lx,*ly; 1118dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 1119a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1120a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 1121a5bc1bf3SBarry Smith ierr = DMDASetOwnershipRanges(dmc2,lx,ly,NULL);CHKERRQ(ierr); 11222be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 1123c73cfb54SMatthew G. Knepley } else if (dim == 1) { 11242be375d4SJed Brown PetscInt *lx; 1125785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 1126a5bc1bf3SBarry Smith ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1127a5bc1bf3SBarry Smith ierr = DMDASetOwnershipRanges(dmc2,lx,NULL,NULL);CHKERRQ(ierr); 11282be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 112947c6ae99SBarry Smith } 1130a5bc1bf3SBarry Smith dd2 = (DM_DA*)dmc2->data; 113147c6ae99SBarry Smith 11324dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 1133a5bc1bf3SBarry Smith /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 1134a5bc1bf3SBarry Smith dmc2->ops->creatematrix = dmf->ops->creatematrix; 1135a5bc1bf3SBarry Smith dmc2->ops->getcoloring = dmf->ops->getcoloring; 113647c6ae99SBarry Smith dd2->interptype = dd->interptype; 113747c6ae99SBarry Smith 113847c6ae99SBarry Smith /* copy fill information if given */ 113947c6ae99SBarry Smith if (dd->dfill) { 1140854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 1141580bdb30SBarry Smith ierr = PetscArraycpy(dd2->dfill,dd->dfill,dd->dfill[dd->w]+dd->w+1);CHKERRQ(ierr); 114247c6ae99SBarry Smith } 114347c6ae99SBarry Smith if (dd->ofill) { 1144854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 1145580bdb30SBarry Smith ierr = PetscArraycpy(dd2->ofill,dd->ofill,dd->ofill[dd->w]+dd->w+1);CHKERRQ(ierr); 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith /* copy the refine information */ 1148397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1149397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1150397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 115147c6ae99SBarry Smith 1152897f7067SBarry Smith if (dd->refine_z_hier) { 1153a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown -1 > -1 && dmf->levelup - dmf->leveldown - 1< dd->refine_z_hier_n) { 1154a5bc1bf3SBarry Smith dd2->refine_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 1]; 1155897f7067SBarry Smith } 1156a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_z_hier_n) { 1157a5bc1bf3SBarry Smith dd2->coarsen_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 2]; 1158897f7067SBarry Smith } 1159897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1160897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1161580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n);CHKERRQ(ierr); 1162897f7067SBarry Smith } 1163897f7067SBarry Smith if (dd->refine_y_hier) { 1164a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1< dd->refine_y_hier_n) { 1165a5bc1bf3SBarry Smith dd2->refine_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 1]; 1166897f7067SBarry Smith } 1167a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_y_hier_n) { 1168a5bc1bf3SBarry Smith dd2->coarsen_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 2]; 1169897f7067SBarry Smith } 1170897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1171897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1172580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n);CHKERRQ(ierr); 1173897f7067SBarry Smith } 1174897f7067SBarry Smith if (dd->refine_x_hier) { 1175a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_x_hier_n) { 1176a5bc1bf3SBarry Smith dd2->refine_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 1]; 1177897f7067SBarry Smith } 1178a5bc1bf3SBarry Smith if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_x_hier_n) { 1179a5bc1bf3SBarry Smith dd2->coarsen_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 2]; 1180897f7067SBarry Smith } 1181897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1182897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1183580bdb30SBarry Smith ierr = PetscArraycpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n);CHKERRQ(ierr); 1184897f7067SBarry Smith } 1185897f7067SBarry Smith 118647c6ae99SBarry Smith /* copy vector type information */ 1187a5bc1bf3SBarry Smith ierr = DMSetVecType(dmc2,dmf->vectype);CHKERRQ(ierr); 118847c6ae99SBarry Smith 1189644e2e5bSBarry Smith dd2->lf = dd->lf; 1190644e2e5bSBarry Smith dd2->lj = dd->lj; 1191644e2e5bSBarry Smith 1192a5bc1bf3SBarry Smith dmc2->leveldown = dmf->leveldown + 1; 1193a5bc1bf3SBarry Smith dmc2->levelup = dmf->levelup; 11948865f1eaSKarl Rupp 1195a5bc1bf3SBarry Smith ierr = DMSetUp(dmc2);CHKERRQ(ierr); 11966e87535bSJed Brown 1197ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 1198a5bc1bf3SBarry Smith if (dmf->coordinates) { 1199ddcf8b74SDave May DM cdaf,cdac; 1200ddcf8b74SDave May Vec coordsc,coordsf; 12016dbf9973SLawrence Mitchell Mat inject; 12026dbf9973SLawrence Mitchell VecScatter vscat; 1203ddcf8b74SDave May 1204a5bc1bf3SBarry Smith ierr = DMGetCoordinateDM(dmf,&cdaf);CHKERRQ(ierr); 1205a5bc1bf3SBarry Smith ierr = DMGetCoordinates(dmf,&coordsf);CHKERRQ(ierr); 1206a5bc1bf3SBarry Smith ierr = DMGetCoordinateDM(dmc2,&cdac);CHKERRQ(ierr); 1207b61d3410SDave May /* force creation of the coordinate vector */ 1208a5bc1bf3SBarry Smith ierr = DMDASetUniformCoordinates(dmc2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1209a5bc1bf3SBarry Smith ierr = DMGetCoordinates(dmc2,&coordsc);CHKERRQ(ierr); 1210ddcf8b74SDave May 1211e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 12126dbf9973SLawrence Mitchell ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr); 12136dbf9973SLawrence Mitchell ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12146dbf9973SLawrence Mitchell ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12156dbf9973SLawrence Mitchell ierr = MatDestroy(&inject);CHKERRQ(ierr); 1216ddcf8b74SDave May } 1217f98405f7SJed Brown 1218a5bc1bf3SBarry Smith for (i=0; i<dmf->bs; i++) { 1219f98405f7SJed Brown const char *fieldname; 1220a5bc1bf3SBarry Smith ierr = DMDAGetFieldName(dmf,i,&fieldname);CHKERRQ(ierr); 1221a5bc1bf3SBarry Smith ierr = DMDASetFieldName(dmc2,i,fieldname);CHKERRQ(ierr); 1222f98405f7SJed Brown } 1223f98405f7SJed Brown 1224a5bc1bf3SBarry Smith *dmc = dmc2; 122547c6ae99SBarry Smith PetscFunctionReturn(0); 122647c6ae99SBarry Smith } 122747c6ae99SBarry Smith 12287087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 122947c6ae99SBarry Smith { 123047c6ae99SBarry Smith PetscErrorCode ierr; 123147c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 123247c6ae99SBarry Smith 123347c6ae99SBarry Smith PetscFunctionBegin; 123447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1235*7a8be351SBarry Smith PetscCheck(nlevels >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 123647c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 123747c6ae99SBarry Smith PetscValidPointer(daf,3); 123847c6ae99SBarry Smith 1239aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 1240dcca6d9dSJed Brown ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr); 124147c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 1242aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 124347c6ae99SBarry Smith } 124447c6ae99SBarry Smith n = nlevels; 1245c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 124647c6ae99SBarry Smith n = nlevels; 1247c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 124847c6ae99SBarry Smith n = nlevels; 1249c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 125047c6ae99SBarry Smith 1251aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1252ce94432eSBarry Smith ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 125347c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1254aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1255ce94432eSBarry Smith ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 125847c6ae99SBarry Smith PetscFunctionReturn(0); 125947c6ae99SBarry Smith } 126047c6ae99SBarry Smith 12617087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 126247c6ae99SBarry Smith { 126347c6ae99SBarry Smith PetscErrorCode ierr; 126447c6ae99SBarry Smith PetscInt i; 126547c6ae99SBarry Smith 126647c6ae99SBarry Smith PetscFunctionBegin; 126747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1268*7a8be351SBarry Smith PetscCheck(nlevels >= 0,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 126947c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 127047c6ae99SBarry Smith PetscValidPointer(dac,3); 1271ce94432eSBarry Smith ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 127247c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1273ce94432eSBarry Smith ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 127447c6ae99SBarry Smith } 127547c6ae99SBarry Smith PetscFunctionReturn(0); 127647c6ae99SBarry Smith } 127762197512SBarry Smith 12788272889dSSatish Balay PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscInt n,PetscReal *nodes) 127962197512SBarry Smith { 128062197512SBarry Smith PetscErrorCode ierr; 12818272889dSSatish Balay PetscInt i,j,xs,xn,q; 128262197512SBarry Smith PetscScalar *xx; 128362197512SBarry Smith PetscReal h; 128462197512SBarry Smith Vec x; 128562197512SBarry Smith DM_DA *da = (DM_DA*)dm->data; 128662197512SBarry Smith 128762197512SBarry Smith PetscFunctionBegin; 128862197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 128962197512SBarry Smith ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 129062197512SBarry Smith q = (q-1)/(n-1); /* number of spectral elements */ 129162197512SBarry Smith h = 2.0/q; 129262197512SBarry Smith ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 129362197512SBarry Smith xs = xs/(n-1); 129462197512SBarry Smith xn = xn/(n-1); 129562197512SBarry Smith ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr); 129662197512SBarry Smith ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 129762197512SBarry Smith ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr); 129862197512SBarry Smith 129962197512SBarry Smith /* loop over local spectral elements */ 130062197512SBarry Smith for (j=xs; j<xs+xn; j++) { 130162197512SBarry Smith /* 130262197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 130362197512SBarry Smith */ 130462197512SBarry Smith for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) { 13058272889dSSatish Balay xx[j*(n-1) + i] = -1.0 + h*j + h*(nodes[i]+1.0)/2.; 130662197512SBarry Smith } 130762197512SBarry Smith } 130862197512SBarry Smith ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr); 130962197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic"); 131062197512SBarry Smith PetscFunctionReturn(0); 131162197512SBarry Smith } 131262197512SBarry Smith 13131fd49c25SBarry Smith /*@ 131462197512SBarry Smith 131562197512SBarry 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 131662197512SBarry Smith 1317d083f849SBarry Smith Collective on da 131862197512SBarry Smith 131962197512SBarry Smith Input Parameters: 132062197512SBarry Smith + da - the DMDA object 13218272889dSSatish Balay - n - the number of GLL nodes 13228272889dSSatish Balay - nodes - the GLL nodes 132362197512SBarry Smith 132495452b02SPatrick Sanan Notes: 132595452b02SPatrick Sanan the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 132662197512SBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is 132762197512SBarry Smith periodic or not. 132862197512SBarry Smith 1329edc382c3SSatish Balay Level: advanced 1330edc382c3SSatish Balay 13318272889dSSatish Balay .seealso: DMDACreate(), PetscDTGaussLobattoLegendreQuadrature(), DMGetCoordinates() 133262197512SBarry Smith @*/ 13338272889dSSatish Balay PetscErrorCode DMDASetGLLCoordinates(DM da,PetscInt n,PetscReal *nodes) 133462197512SBarry Smith { 133562197512SBarry Smith PetscErrorCode ierr; 133662197512SBarry Smith 133762197512SBarry Smith PetscFunctionBegin; 133862197512SBarry Smith if (da->dim == 1) { 13398272889dSSatish Balay ierr = DMDASetGLLCoordinates_1d(da,n,nodes);CHKERRQ(ierr); 134062197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d"); 134162197512SBarry Smith PetscFunctionReturn(0); 134262197512SBarry Smith } 13437c3cd84eSPatrick Sanan 13447c3cd84eSPatrick Sanan PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1,DM dm2,PetscBool *compatible,PetscBool *set) 13457c3cd84eSPatrick Sanan { 13467c3cd84eSPatrick Sanan PetscErrorCode ierr; 13477c3cd84eSPatrick Sanan DM_DA *dd1 = (DM_DA*)da1->data,*dd2; 13487c3cd84eSPatrick Sanan DM da2; 13497c3cd84eSPatrick Sanan DMType dmtype2; 13507c3cd84eSPatrick Sanan PetscBool isda,compatibleLocal; 13517c3cd84eSPatrick Sanan PetscInt i; 13527c3cd84eSPatrick Sanan 13537c3cd84eSPatrick Sanan PetscFunctionBegin; 1354*7a8be351SBarry Smith PetscCheck(da1->setupcalled,PetscObjectComm((PetscObject)da1),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on first DM before DMGetCompatibility()"); 13557c3cd84eSPatrick Sanan ierr = DMGetType(dm2,&dmtype2);CHKERRQ(ierr); 13567c3cd84eSPatrick Sanan ierr = PetscStrcmp(dmtype2,DMDA,&isda);CHKERRQ(ierr); 13577c3cd84eSPatrick Sanan if (isda) { 13587c3cd84eSPatrick Sanan da2 = dm2; 13597c3cd84eSPatrick Sanan dd2 = (DM_DA*)da2->data; 1360*7a8be351SBarry Smith PetscCheck(da2->setupcalled,PetscObjectComm((PetscObject)da2),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on second DM before DMGetCompatibility()"); 1361110623a0SKarl Rupp compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1362c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1363110623a0SKarl Rupp /* Global size ranks Boundary type */ 1364c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1365c790a739SKarl Rupp if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1366c790a739SKarl Rupp if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 13677c3cd84eSPatrick Sanan if (compatibleLocal) { 13687c3cd84eSPatrick Sanan for (i=0; i<dd1->m; ++i) { 1369c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ 13707c3cd84eSPatrick Sanan } 13717c3cd84eSPatrick Sanan } 13727c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 1) { 13737c3cd84eSPatrick Sanan for (i=0; i<dd1->n; ++i) { 1374c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 13757c3cd84eSPatrick Sanan } 13767c3cd84eSPatrick Sanan } 13777c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 2) { 13787c3cd84eSPatrick Sanan for (i=0; i<dd1->p; ++i) { 1379c790a739SKarl Rupp compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 13807c3cd84eSPatrick Sanan } 13817c3cd84eSPatrick Sanan } 13827c3cd84eSPatrick Sanan *compatible = compatibleLocal; 13837c3cd84eSPatrick Sanan *set = PETSC_TRUE; 13847c3cd84eSPatrick Sanan } else { 13857c3cd84eSPatrick Sanan /* Decline to determine compatibility with other DM types */ 13867c3cd84eSPatrick Sanan *set = PETSC_FALSE; 13877c3cd84eSPatrick Sanan } 13887c3cd84eSPatrick Sanan PetscFunctionReturn(0); 13897c3cd84eSPatrick Sanan } 1390