1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 347c6ae99SBarry Smith /*@ 4*e43dc8daSBarry 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 10*e43dc8daSBarry Smith . M - the global X size 11*e43dc8daSBarry Smith . N - the global Y size 12*e43dc8daSBarry Smith - P - the global Z size 1347c6ae99SBarry Smith 1447c6ae99SBarry Smith Level: intermediate 1547c6ae99SBarry Smith 16*e43dc8daSBarry Smith Developer Notes: Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points 17*e43dc8daSBarry Smith 18aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership() 1947c6ae99SBarry Smith @*/ 207087cfbeSBarry Smith PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 2147c6ae99SBarry Smith { 2247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith PetscFunctionBegin; 2547c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 2647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,M,2); 2747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,N,3); 2847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,P,4); 29ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 30*e43dc8daSBarry Smith if (M < 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive"); 31*e43dc8daSBarry Smith if (N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive"); 32*e43dc8daSBarry Smith if (P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Z direction must be positive"); 3347c6ae99SBarry Smith 3447c6ae99SBarry Smith dd->M = M; 3547c6ae99SBarry Smith dd->N = N; 3647c6ae99SBarry Smith dd->P = P; 3747c6ae99SBarry Smith PetscFunctionReturn(0); 3847c6ae99SBarry Smith } 3947c6ae99SBarry Smith 4047c6ae99SBarry Smith /*@ 41aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 4247c6ae99SBarry Smith 43aa219208SBarry Smith Logically Collective on DMDA 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith Input Parameters: 46aa219208SBarry Smith + da - the DMDA 4747c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE) 4847c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE) 4947c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE) 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith Level: intermediate 5247c6ae99SBarry Smith 53aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership() 5447c6ae99SBarry Smith @*/ 557087cfbeSBarry Smith PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 5647c6ae99SBarry Smith { 5747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 58e3f3e4b6SBarry Smith PetscErrorCode ierr; 5947c6ae99SBarry Smith 6047c6ae99SBarry Smith PetscFunctionBegin; 6147c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 6247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,m,2); 6347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,n,3); 6447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,p,4); 65ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 6647c6ae99SBarry Smith dd->m = m; 6747c6ae99SBarry Smith dd->n = n; 6847c6ae99SBarry Smith dd->p = p; 69c73cfb54SMatthew G. Knepley if (da->dim == 2) { 70d3be247eSBarry Smith PetscMPIInt size; 71d3be247eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 72e3f3e4b6SBarry Smith if ((dd->m > 0) && (dd->n < 0)) { 73e3f3e4b6SBarry Smith dd->n = size/dd->m; 74ce94432eSBarry 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); 75e3f3e4b6SBarry Smith } 76e3f3e4b6SBarry Smith if ((dd->n > 0) && (dd->m < 0)) { 77e3f3e4b6SBarry Smith dd->m = size/dd->n; 78ce94432eSBarry 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); 79e3f3e4b6SBarry Smith } 80e3f3e4b6SBarry Smith } 8147c6ae99SBarry Smith PetscFunctionReturn(0); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith 8447c6ae99SBarry Smith /*@ 851321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith Not collective 8847c6ae99SBarry Smith 8947c6ae99SBarry Smith Input Parameter: 90aa219208SBarry Smith + da - The DMDA 91bff4a2f0SMatthew G. Knepley - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith Level: intermediate 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith .keywords: distributed array, periodicity 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; 10347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1045a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bx,2); 1055a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,by,3); 1065a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,bz,4); 107ce94432eSBarry Smith if (da->setupcalled) SETERRQ(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 12547c6ae99SBarry Smith .keywords: distributed array, degrees of freedom 126fb6725baSMatthew G. Knepley .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA 12747c6ae99SBarry Smith @*/ 12854cfb0beSLisandro Dalcin PetscErrorCode DMDASetDof(DM da, PetscInt dof) 12947c6ae99SBarry Smith { 13047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith PetscFunctionBegin; 13347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 13454cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da,dof,2); 135ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 13647c6ae99SBarry Smith dd->w = dof; 1371411c6eeSJed Brown da->bs = dof; 13847c6ae99SBarry Smith PetscFunctionReturn(0); 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith 141fb6725baSMatthew G. Knepley /*@ 142fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 143fb6725baSMatthew G. Knepley 144fb6725baSMatthew G. Knepley Not collective 145fb6725baSMatthew G. Knepley 146fb6725baSMatthew G. Knepley Input Parameter: 147fb6725baSMatthew G. Knepley . da - The DMDA 148fb6725baSMatthew G. Knepley 149fb6725baSMatthew G. Knepley Output Parameter: 150fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom 151fb6725baSMatthew G. Knepley 152fb6725baSMatthew G. Knepley Level: intermediate 153fb6725baSMatthew G. Knepley 154fb6725baSMatthew G. Knepley .keywords: distributed array, degrees of freedom 155fb6725baSMatthew G. Knepley .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA 156fb6725baSMatthew G. Knepley @*/ 157fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 158fb6725baSMatthew G. Knepley { 159fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 160fb6725baSMatthew G. Knepley 161fb6725baSMatthew G. Knepley PetscFunctionBegin; 162fb6725baSMatthew G. Knepley PetscValidHeaderSpecific(da,DM_CLASSID,1); 163fb6725baSMatthew G. Knepley PetscValidPointer(dof,2); 164fb6725baSMatthew G. Knepley *dof = dd->w; 165fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 166fb6725baSMatthew G. Knepley } 167fb6725baSMatthew G. Knepley 1687ddda789SPeter Brune /*@ 1697ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1707ddda789SPeter Brune 1717ddda789SPeter Brune Not collective 1727ddda789SPeter Brune 1737ddda789SPeter Brune Input Parameters: 1747ddda789SPeter Brune . da - The DMDA 1757ddda789SPeter Brune 1767ddda789SPeter Brune Output Parameters: 1777ddda789SPeter Brune + x - Overlap in the x direction 1787ddda789SPeter Brune . y - Overlap in the y direction 1797ddda789SPeter Brune - z - Overlap in the z direction 1807ddda789SPeter Brune 1817ddda789SPeter Brune Level: intermediate 1827ddda789SPeter Brune 1837ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 1847ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA 1857ddda789SPeter Brune @*/ 1867ddda789SPeter Brune PetscErrorCode DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z) 1877ddda789SPeter Brune { 1887ddda789SPeter Brune DM_DA *dd = (DM_DA*)da->data; 1897ddda789SPeter Brune 1907ddda789SPeter Brune PetscFunctionBegin; 1917ddda789SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 1927ddda789SPeter Brune if (x) *x = dd->xol; 1937ddda789SPeter Brune if (y) *y = dd->yol; 1947ddda789SPeter Brune if (z) *z = dd->zol; 1957ddda789SPeter Brune PetscFunctionReturn(0); 1967ddda789SPeter Brune } 1977ddda789SPeter Brune 19888661749SPeter Brune /*@ 19988661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 20088661749SPeter Brune 20188661749SPeter Brune Not collective 20288661749SPeter Brune 2037ddda789SPeter Brune Input Parameters: 20488661749SPeter Brune + da - The DMDA 2057ddda789SPeter Brune . x - Overlap in the x direction 2067ddda789SPeter Brune . y - Overlap in the y direction 2077ddda789SPeter Brune - z - Overlap in the z direction 20888661749SPeter Brune 20988661749SPeter Brune Level: intermediate 21088661749SPeter Brune 2117ddda789SPeter Brune .keywords: distributed array, overlap, domain decomposition 2127ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA 21388661749SPeter Brune @*/ 2147ddda789SPeter Brune PetscErrorCode DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z) 21588661749SPeter Brune { 21688661749SPeter Brune DM_DA *dd = (DM_DA*)da->data; 21788661749SPeter Brune 21888661749SPeter Brune PetscFunctionBegin; 21988661749SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2207ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,x,2); 2217ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,y,3); 2227ddda789SPeter Brune PetscValidLogicalCollectiveInt(da,z,4); 2237ddda789SPeter Brune dd->xol = x; 2247ddda789SPeter Brune dd->yol = y; 2257ddda789SPeter Brune dd->zol = z; 22688661749SPeter Brune PetscFunctionReturn(0); 22788661749SPeter Brune } 22888661749SPeter Brune 2293e7870d2SPeter Brune 2303e7870d2SPeter Brune /*@ 2313e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2323e7870d2SPeter Brune 2333e7870d2SPeter Brune Not collective 2343e7870d2SPeter Brune 2353e7870d2SPeter Brune Input Parameters: 2363e7870d2SPeter Brune . da - The DMDA 2373e7870d2SPeter Brune 2383e7870d2SPeter Brune Output Parameters: 2393e7870d2SPeter Brune + Nsub - Number of local subdomains created upon decomposition 2403e7870d2SPeter Brune 2413e7870d2SPeter Brune Level: intermediate 2423e7870d2SPeter Brune 2433e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2443e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA 2453e7870d2SPeter Brune @*/ 2463e7870d2SPeter Brune PetscErrorCode DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub) 2473e7870d2SPeter Brune { 2483e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2493e7870d2SPeter Brune 2503e7870d2SPeter Brune PetscFunctionBegin; 2513e7870d2SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2523e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2533e7870d2SPeter Brune PetscFunctionReturn(0); 2543e7870d2SPeter Brune } 2553e7870d2SPeter Brune 2563e7870d2SPeter Brune /*@ 2573e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2583e7870d2SPeter Brune 2593e7870d2SPeter Brune Not collective 2603e7870d2SPeter Brune 2613e7870d2SPeter Brune Input Parameters: 2623e7870d2SPeter Brune + da - The DMDA 2633e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2643e7870d2SPeter Brune 2653e7870d2SPeter Brune Level: intermediate 2663e7870d2SPeter Brune 2673e7870d2SPeter Brune .keywords: distributed array, domain decomposition 2683e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA 2693e7870d2SPeter Brune @*/ 2703e7870d2SPeter Brune PetscErrorCode DMDASetNumLocalSubDomains(DM da,PetscInt Nsub) 2713e7870d2SPeter Brune { 2723e7870d2SPeter Brune DM_DA *dd = (DM_DA*)da->data; 2733e7870d2SPeter Brune 2743e7870d2SPeter Brune PetscFunctionBegin; 2753e7870d2SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 2763e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da,Nsub,2); 2773e7870d2SPeter Brune dd->Nsub = Nsub; 2783e7870d2SPeter Brune PetscFunctionReturn(0); 2793e7870d2SPeter Brune } 2803e7870d2SPeter Brune 281d886c4f4SPeter Brune /*@ 282d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 283d886c4f4SPeter Brune 284d886c4f4SPeter Brune Collective on DA 285d886c4f4SPeter Brune 286d886c4f4SPeter Brune Input Parameter: 287d886c4f4SPeter Brune + da - The DMDA 288d886c4f4SPeter Brune . xo - The offset in the x direction 289d886c4f4SPeter Brune . yo - The offset in the y direction 290d886c4f4SPeter Brune - zo - The offset in the z direction 291d886c4f4SPeter Brune 292d886c4f4SPeter Brune Level: intermediate 293d886c4f4SPeter Brune 294d886c4f4SPeter Brune Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without 295d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 296d886c4f4SPeter Brune 297d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 298d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 299d886c4f4SPeter Brune @*/ 30095c13181SPeter Brune PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 301d886c4f4SPeter Brune { 30295c13181SPeter Brune PetscErrorCode ierr; 303d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 304d886c4f4SPeter Brune 305d886c4f4SPeter Brune PetscFunctionBegin; 306d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 307d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da,xo,2); 30895c13181SPeter Brune PetscValidLogicalCollectiveInt(da,yo,3); 30995c13181SPeter Brune PetscValidLogicalCollectiveInt(da,zo,4); 31095c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Mo,5); 31195c13181SPeter Brune PetscValidLogicalCollectiveInt(da,No,6); 31295c13181SPeter Brune PetscValidLogicalCollectiveInt(da,Po,7); 313d886c4f4SPeter Brune dd->xo = xo; 314d886c4f4SPeter Brune dd->yo = yo; 315d886c4f4SPeter Brune dd->zo = zo; 31695c13181SPeter Brune dd->Mo = Mo; 31795c13181SPeter Brune dd->No = No; 31895c13181SPeter Brune dd->Po = Po; 31995c13181SPeter Brune 32095c13181SPeter Brune if (da->coordinateDM) { 32195c13181SPeter Brune ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 32295c13181SPeter Brune } 323d886c4f4SPeter Brune PetscFunctionReturn(0); 324d886c4f4SPeter Brune } 325d886c4f4SPeter Brune 326d886c4f4SPeter Brune /*@ 327d886c4f4SPeter Brune DMDAGetOffset - Gets the index offset of the DA. 328d886c4f4SPeter Brune 329d886c4f4SPeter Brune Not collective 330d886c4f4SPeter Brune 331d886c4f4SPeter Brune Input Parameter: 332d886c4f4SPeter Brune . da - The DMDA 333d886c4f4SPeter Brune 334d886c4f4SPeter Brune Output Parameters: 335d886c4f4SPeter Brune + xo - The offset in the x direction 336d886c4f4SPeter Brune . yo - The offset in the y direction 33795c13181SPeter Brune . zo - The offset in the z direction 33895c13181SPeter Brune . Mo - The global size in the x direction 33995c13181SPeter Brune . No - The global size in the y direction 34095c13181SPeter Brune - Po - The global size in the z direction 341d886c4f4SPeter Brune 342d886c4f4SPeter Brune Level: intermediate 343d886c4f4SPeter Brune 344d886c4f4SPeter Brune .keywords: distributed array, degrees of freedom 345d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray() 346d886c4f4SPeter Brune @*/ 34795c13181SPeter Brune PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 348d886c4f4SPeter Brune { 349d886c4f4SPeter Brune DM_DA *dd = (DM_DA*)da->data; 350d886c4f4SPeter Brune 351d886c4f4SPeter Brune PetscFunctionBegin; 352d886c4f4SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 353d886c4f4SPeter Brune if (xo) *xo = dd->xo; 354d886c4f4SPeter Brune if (yo) *yo = dd->yo; 355d886c4f4SPeter Brune if (zo) *zo = dd->zo; 35695c13181SPeter Brune if (Mo) *Mo = dd->Mo; 35795c13181SPeter Brune if (No) *No = dd->No; 35895c13181SPeter Brune if (Po) *Po = dd->Po; 359d886c4f4SPeter Brune PetscFunctionReturn(0); 360d886c4f4SPeter Brune } 361d886c4f4SPeter Brune 36240234c92SPeter Brune /*@ 36340234c92SPeter Brune DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM. 36440234c92SPeter Brune 36540234c92SPeter Brune Not collective 36640234c92SPeter Brune 36740234c92SPeter Brune Input Parameter: 36840234c92SPeter Brune . da - The DMDA 36940234c92SPeter Brune 37040234c92SPeter Brune Output Parameters: 37140234c92SPeter Brune + xs - The start of the region in x 37240234c92SPeter Brune . ys - The start of the region in y 37340234c92SPeter Brune . zs - The start of the region in z 37440234c92SPeter Brune . xs - The size of the region in x 37540234c92SPeter Brune . ys - The size of the region in y 37640234c92SPeter Brune . zs - The size of the region in z 37740234c92SPeter Brune 37840234c92SPeter Brune Level: intermediate 37940234c92SPeter Brune 38040234c92SPeter Brune .keywords: distributed array, degrees of freedom 38140234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 38240234c92SPeter Brune @*/ 38340234c92SPeter Brune PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 38440234c92SPeter Brune { 38540234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 38640234c92SPeter Brune 38740234c92SPeter Brune PetscFunctionBegin; 38840234c92SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 38940234c92SPeter Brune if (xs) *xs = dd->nonxs; 39040234c92SPeter Brune if (ys) *ys = dd->nonys; 39140234c92SPeter Brune if (zs) *zs = dd->nonzs; 39240234c92SPeter Brune if (xm) *xm = dd->nonxm; 39340234c92SPeter Brune if (ym) *ym = dd->nonym; 39440234c92SPeter Brune if (zm) *zm = dd->nonzm; 39540234c92SPeter Brune PetscFunctionReturn(0); 39640234c92SPeter Brune } 39740234c92SPeter Brune 39840234c92SPeter Brune 39940234c92SPeter Brune /*@ 40040234c92SPeter Brune DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM. 40140234c92SPeter Brune 40240234c92SPeter Brune Collective on DA 40340234c92SPeter Brune 40440234c92SPeter Brune Input Parameter: 40540234c92SPeter Brune + da - The DMDA 40640234c92SPeter Brune . xs - The start of the region in x 40740234c92SPeter Brune . ys - The start of the region in y 40840234c92SPeter Brune . zs - The start of the region in z 40940234c92SPeter Brune . xs - The size of the region in x 41040234c92SPeter Brune . ys - The size of the region in y 41140234c92SPeter Brune . zs - The size of the region in z 41240234c92SPeter Brune 41340234c92SPeter Brune Level: intermediate 41440234c92SPeter Brune 41540234c92SPeter Brune .keywords: distributed array, degrees of freedom 41640234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray() 41740234c92SPeter Brune @*/ 41840234c92SPeter Brune PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 41940234c92SPeter Brune { 42040234c92SPeter Brune DM_DA *dd = (DM_DA*)da->data; 42140234c92SPeter Brune 42240234c92SPeter Brune PetscFunctionBegin; 42340234c92SPeter Brune PetscValidHeaderSpecific(da,DM_CLASSID,1); 42440234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xs,2); 42540234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ys,3); 42640234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zs,4); 42740234c92SPeter Brune PetscValidLogicalCollectiveInt(da,xm,5); 42840234c92SPeter Brune PetscValidLogicalCollectiveInt(da,ym,6); 42940234c92SPeter Brune PetscValidLogicalCollectiveInt(da,zm,7); 43040234c92SPeter Brune dd->nonxs = xs; 43140234c92SPeter Brune dd->nonys = ys; 43240234c92SPeter Brune dd->nonzs = zs; 43340234c92SPeter Brune dd->nonxm = xm; 43440234c92SPeter Brune dd->nonym = ym; 43540234c92SPeter Brune dd->nonzm = zm; 43640234c92SPeter Brune 43740234c92SPeter Brune PetscFunctionReturn(0); 43840234c92SPeter Brune } 43988661749SPeter Brune 44047c6ae99SBarry Smith /*@ 441aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 44247c6ae99SBarry Smith 443aa219208SBarry Smith Logically Collective on DMDA 44447c6ae99SBarry Smith 44547c6ae99SBarry Smith Input Parameter: 446aa219208SBarry Smith + da - The DMDA 447aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith Level: intermediate 45047c6ae99SBarry Smith 45147c6ae99SBarry Smith .keywords: distributed array, stencil 45296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 45347c6ae99SBarry Smith @*/ 4547087cfbeSBarry Smith PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 45547c6ae99SBarry Smith { 45647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 45747c6ae99SBarry Smith 45847c6ae99SBarry Smith PetscFunctionBegin; 45947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 46047c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,stype,2); 461ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 46247c6ae99SBarry Smith dd->stencil_type = stype; 46347c6ae99SBarry Smith PetscFunctionReturn(0); 46447c6ae99SBarry Smith } 46547c6ae99SBarry Smith 466fb6725baSMatthew G. Knepley /*@ 467fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 468fb6725baSMatthew G. Knepley 469fb6725baSMatthew G. Knepley Not collective 470fb6725baSMatthew G. Knepley 471fb6725baSMatthew G. Knepley Input Parameter: 472fb6725baSMatthew G. Knepley . da - The DMDA 473fb6725baSMatthew G. Knepley 474fb6725baSMatthew G. Knepley Output Parameter: 475fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 476fb6725baSMatthew G. Knepley 477fb6725baSMatthew G. Knepley Level: intermediate 478fb6725baSMatthew G. Knepley 479fb6725baSMatthew G. Knepley .keywords: distributed array, stencil 480fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 481fb6725baSMatthew G. Knepley @*/ 482fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 483fb6725baSMatthew G. Knepley { 484fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA*)da->data; 485fb6725baSMatthew G. Knepley 486fb6725baSMatthew G. Knepley PetscFunctionBegin; 487fb6725baSMatthew G. Knepley PetscValidHeaderSpecific(da,DM_CLASSID,1); 488fb6725baSMatthew G. Knepley PetscValidPointer(stype,2); 489fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 490fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 491fb6725baSMatthew G. Knepley } 492fb6725baSMatthew G. Knepley 49347c6ae99SBarry Smith /*@ 494aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 49547c6ae99SBarry Smith 496aa219208SBarry Smith Logically Collective on DMDA 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith Input Parameter: 499aa219208SBarry Smith + da - The DMDA 50047c6ae99SBarry Smith - width - The stencil width 50147c6ae99SBarry Smith 50247c6ae99SBarry Smith Level: intermediate 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith .keywords: distributed array, stencil 50596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 50647c6ae99SBarry Smith @*/ 5077087cfbeSBarry Smith PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 50847c6ae99SBarry Smith { 50947c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 51047c6ae99SBarry Smith 51147c6ae99SBarry Smith PetscFunctionBegin; 51247c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 51347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,width,2); 514ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 51547c6ae99SBarry Smith dd->s = width; 51647c6ae99SBarry Smith PetscFunctionReturn(0); 51747c6ae99SBarry Smith } 51847c6ae99SBarry Smith 519fb6725baSMatthew G. Knepley /*@ 520fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 521fb6725baSMatthew G. Knepley 522fb6725baSMatthew G. Knepley Not collective 523fb6725baSMatthew G. Knepley 524fb6725baSMatthew G. Knepley Input Parameter: 525fb6725baSMatthew G. Knepley . da - The DMDA 526fb6725baSMatthew G. Knepley 527fb6725baSMatthew G. Knepley Output Parameter: 528fb6725baSMatthew G. Knepley . width - The stencil width 529fb6725baSMatthew G. Knepley 530fb6725baSMatthew G. Knepley Level: intermediate 531fb6725baSMatthew G. Knepley 532fb6725baSMatthew G. Knepley .keywords: distributed array, stencil 533fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA 534fb6725baSMatthew G. Knepley @*/ 535fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 536fb6725baSMatthew G. Knepley { 537fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *) da->data; 538fb6725baSMatthew G. Knepley 539fb6725baSMatthew G. Knepley PetscFunctionBegin; 540fb6725baSMatthew G. Knepley PetscValidHeaderSpecific(da,DM_CLASSID,1); 541fb6725baSMatthew G. Knepley PetscValidPointer(width,2); 542fb6725baSMatthew G. Knepley *width = dd->s; 543fb6725baSMatthew G. Knepley PetscFunctionReturn(0); 544fb6725baSMatthew G. Knepley } 545fb6725baSMatthew G. Knepley 546aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 54747c6ae99SBarry Smith { 54847c6ae99SBarry Smith PetscInt i,sum; 54947c6ae99SBarry Smith 55047c6ae99SBarry Smith PetscFunctionBegin; 551ce94432eSBarry Smith if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 55247c6ae99SBarry Smith for (i=sum=0; i<m; i++) sum += lx[i]; 553ce94432eSBarry Smith if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 55447c6ae99SBarry Smith PetscFunctionReturn(0); 55547c6ae99SBarry Smith } 55647c6ae99SBarry Smith 55747c6ae99SBarry Smith /*@ 558aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 55947c6ae99SBarry Smith 560aa219208SBarry Smith Logically Collective on DMDA 56147c6ae99SBarry Smith 56247c6ae99SBarry Smith Input Parameter: 563aa219208SBarry Smith + da - The DMDA 5640298fd71SBarry 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 5650298fd71SBarry 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 5660298fd71SBarry 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. 56747c6ae99SBarry Smith 56847c6ae99SBarry Smith Level: intermediate 56947c6ae99SBarry Smith 570e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 571e3f3e4b6SBarry Smith 57247c6ae99SBarry Smith .keywords: distributed array 57396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA 57447c6ae99SBarry Smith @*/ 5757087cfbeSBarry Smith PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 57647c6ae99SBarry Smith { 57747c6ae99SBarry Smith PetscErrorCode ierr; 57847c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 57947c6ae99SBarry Smith 58047c6ae99SBarry Smith PetscFunctionBegin; 58147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 582ce94432eSBarry Smith if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 58347c6ae99SBarry Smith if (lx) { 584ce94432eSBarry Smith if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 585aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 58647c6ae99SBarry Smith if (!dd->lx) { 587785e854fSJed Brown ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr); 58847c6ae99SBarry Smith } 58947c6ae99SBarry Smith ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 59047c6ae99SBarry Smith } 59147c6ae99SBarry Smith if (ly) { 592ce94432eSBarry Smith if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 593aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 59447c6ae99SBarry Smith if (!dd->ly) { 595785e854fSJed Brown ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr); 59647c6ae99SBarry Smith } 59747c6ae99SBarry Smith ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 59847c6ae99SBarry Smith } 59947c6ae99SBarry Smith if (lz) { 600ce94432eSBarry Smith if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 601aa219208SBarry Smith ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 60247c6ae99SBarry Smith if (!dd->lz) { 603785e854fSJed Brown ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr); 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 60647c6ae99SBarry Smith } 60747c6ae99SBarry Smith PetscFunctionReturn(0); 60847c6ae99SBarry Smith } 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith /*@ 611aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 612e727c939SJed Brown returned by DMCreateInterpolation() 61347c6ae99SBarry Smith 614aa219208SBarry Smith Logically Collective on DMDA 61547c6ae99SBarry Smith 61647c6ae99SBarry Smith Input Parameter: 61747c6ae99SBarry Smith + da - initial distributed array 618aa219208SBarry Smith . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith Level: intermediate 62147c6ae99SBarry Smith 622e727c939SJed Brown Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 62347c6ae99SBarry Smith 62447c6ae99SBarry Smith .keywords: distributed array, interpolation 62547c6ae99SBarry Smith 62696e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 62747c6ae99SBarry Smith @*/ 6287087cfbeSBarry Smith PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 62947c6ae99SBarry Smith { 63047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 63147c6ae99SBarry Smith 63247c6ae99SBarry Smith PetscFunctionBegin; 63347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 63447c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da,ctype,2); 63547c6ae99SBarry Smith dd->interptype = ctype; 63647c6ae99SBarry Smith PetscFunctionReturn(0); 63747c6ae99SBarry Smith } 63847c6ae99SBarry Smith 6392dde6fd4SLisandro Dalcin /*@ 6402dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 641e727c939SJed Brown used by DMCreateInterpolation() 6422dde6fd4SLisandro Dalcin 6432dde6fd4SLisandro Dalcin Not Collective 6442dde6fd4SLisandro Dalcin 6452dde6fd4SLisandro Dalcin Input Parameter: 6462dde6fd4SLisandro Dalcin . da - distributed array 6472dde6fd4SLisandro Dalcin 6482dde6fd4SLisandro Dalcin Output Parameter: 6492dde6fd4SLisandro Dalcin . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 6502dde6fd4SLisandro Dalcin 6512dde6fd4SLisandro Dalcin Level: intermediate 6522dde6fd4SLisandro Dalcin 6532dde6fd4SLisandro Dalcin .keywords: distributed array, interpolation 6542dde6fd4SLisandro Dalcin 655e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 6562dde6fd4SLisandro Dalcin @*/ 6572dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 6582dde6fd4SLisandro Dalcin { 6592dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 6602dde6fd4SLisandro Dalcin 6612dde6fd4SLisandro Dalcin PetscFunctionBegin; 6622dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 6632dde6fd4SLisandro Dalcin PetscValidPointer(ctype,2); 6642dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6652dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 6662dde6fd4SLisandro Dalcin } 66747c6ae99SBarry Smith 6686a119db4SBarry Smith /*@C 669aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 67047c6ae99SBarry Smith processes neighbors. 67147c6ae99SBarry Smith 67247c6ae99SBarry Smith Not Collective 67347c6ae99SBarry Smith 67447c6ae99SBarry Smith Input Parameter: 675aa219208SBarry Smith . da - the DMDA object 67647c6ae99SBarry Smith 67747c6ae99SBarry Smith Output Parameters: 67847c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 67947c6ae99SBarry Smith this process itself is in the list 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith Notes: In 2d the array is of length 9, in 3d of length 27 68247c6ae99SBarry Smith Not supported in 1d 683aa219208SBarry Smith Do not free the array, it is freed when the DMDA is destroyed. 68447c6ae99SBarry Smith 68547c6ae99SBarry Smith Fortran Notes: In fortran you must pass in an array of the appropriate length. 68647c6ae99SBarry Smith 68747c6ae99SBarry Smith Level: intermediate 68847c6ae99SBarry Smith 68947c6ae99SBarry Smith @*/ 6907087cfbeSBarry Smith PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 69147c6ae99SBarry Smith { 69247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 6935fd66863SKarl Rupp 69447c6ae99SBarry Smith PetscFunctionBegin; 69547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 69647c6ae99SBarry Smith *ranks = dd->neighbors; 69747c6ae99SBarry Smith PetscFunctionReturn(0); 69847c6ae99SBarry Smith } 69947c6ae99SBarry Smith 70047c6ae99SBarry Smith /*@C 701aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 70247c6ae99SBarry Smith 70347c6ae99SBarry Smith Not Collective 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith Input Parameter: 706aa219208SBarry Smith . da - the DMDA object 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith Output Parameter: 70947c6ae99SBarry Smith + lx - ownership along x direction (optional) 71047c6ae99SBarry Smith . ly - ownership along y direction (optional) 71147c6ae99SBarry Smith - lz - ownership along z direction (optional) 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith Level: intermediate 71447c6ae99SBarry Smith 715aa219208SBarry Smith Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 71647c6ae99SBarry Smith 71747c6ae99SBarry Smith In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 718aa219208SBarry Smith eighth arguments from DMDAGetInfo() 71947c6ae99SBarry Smith 72047c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 721aa219208SBarry Smith DMDA they came from still exists (has not been destroyed). 72247c6ae99SBarry Smith 723e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 724e3f3e4b6SBarry Smith 7251fd49c25SBarry Smith Not available from Fortran 7261fd49c25SBarry Smith 727aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 72847c6ae99SBarry Smith @*/ 7297087cfbeSBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 73047c6ae99SBarry Smith { 73147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith PetscFunctionBegin; 73447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 73547c6ae99SBarry Smith if (lx) *lx = dd->lx; 73647c6ae99SBarry Smith if (ly) *ly = dd->ly; 73747c6ae99SBarry Smith if (lz) *lz = dd->lz; 73847c6ae99SBarry Smith PetscFunctionReturn(0); 73947c6ae99SBarry Smith } 74047c6ae99SBarry Smith 74147c6ae99SBarry Smith /*@ 742aa219208SBarry Smith DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 74347c6ae99SBarry Smith 744aa219208SBarry Smith Logically Collective on DMDA 74547c6ae99SBarry Smith 74647c6ae99SBarry Smith Input Parameters: 747aa219208SBarry Smith + da - the DMDA object 74847c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 74947c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 75047c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 75147c6ae99SBarry Smith 75247c6ae99SBarry Smith Options Database: 75347c6ae99SBarry Smith + -da_refine_x - refinement ratio in x direction 75447c6ae99SBarry Smith . -da_refine_y - refinement ratio in y direction 75547c6ae99SBarry Smith - -da_refine_z - refinement ratio in z direction 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith Level: intermediate 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith Notes: Pass PETSC_IGNORE to leave a value unchanged 76047c6ae99SBarry Smith 761aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor() 76247c6ae99SBarry Smith @*/ 7637087cfbeSBarry Smith PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 76447c6ae99SBarry Smith { 76547c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 76647c6ae99SBarry Smith 76747c6ae99SBarry Smith PetscFunctionBegin; 76847c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 76947c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_x,2); 77047c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_y,3); 77147c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da,refine_z,4); 77247c6ae99SBarry Smith 77347c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 77447c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 77547c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 77647c6ae99SBarry Smith PetscFunctionReturn(0); 77747c6ae99SBarry Smith } 77847c6ae99SBarry Smith 77947c6ae99SBarry Smith /*@C 780aa219208SBarry Smith DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith Not Collective 78347c6ae99SBarry Smith 78447c6ae99SBarry Smith Input Parameter: 785aa219208SBarry Smith . da - the DMDA object 78647c6ae99SBarry Smith 78747c6ae99SBarry Smith Output Parameters: 78847c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 78947c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 79047c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith Level: intermediate 79347c6ae99SBarry Smith 7940298fd71SBarry Smith Notes: Pass NULL for values you do not need 79547c6ae99SBarry Smith 796aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor() 79747c6ae99SBarry Smith @*/ 7987087cfbeSBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 79947c6ae99SBarry Smith { 80047c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith PetscFunctionBegin; 80347c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 80447c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 80547c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 80647c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 80747c6ae99SBarry Smith PetscFunctionReturn(0); 80847c6ae99SBarry Smith } 80947c6ae99SBarry Smith 81047c6ae99SBarry Smith /*@C 811aa219208SBarry Smith DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 81247c6ae99SBarry Smith 813aa219208SBarry Smith Logically Collective on DMDA 81447c6ae99SBarry Smith 81547c6ae99SBarry Smith Input Parameters: 816aa219208SBarry Smith + da - the DMDA object 817aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith Level: developer 82047c6ae99SBarry Smith 821aa219208SBarry Smith Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 82247c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 82347c6ae99SBarry Smith 8241fd49c25SBarry Smith Not supported from Fortran 8251fd49c25SBarry Smith 826950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills() 82747c6ae99SBarry Smith @*/ 828b412c318SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*)) 82947c6ae99SBarry Smith { 83047c6ae99SBarry Smith PetscFunctionBegin; 83147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 83225296bd5SBarry Smith da->ops->creatematrix = f; 83347c6ae99SBarry Smith PetscFunctionReturn(0); 83447c6ae99SBarry Smith } 83547c6ae99SBarry Smith 83647c6ae99SBarry Smith /* 83747c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 83847c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 83947c6ae99SBarry Smith 84047c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 84147c6ae99SBarry Smith */ 84294013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 84347c6ae99SBarry Smith { 84447c6ae99SBarry Smith PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 84547c6ae99SBarry Smith PetscErrorCode ierr; 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith PetscFunctionBegin; 848ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 84947c6ae99SBarry Smith if (ratio == 1) { 85047c6ae99SBarry Smith ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 85147c6ae99SBarry Smith PetscFunctionReturn(0); 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith for (i=0; i<m; i++) totalc += lc[i]; 85447c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 85547c6ae99SBarry Smith for (i=0; i<m; i++) { 85647c6ae99SBarry Smith PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 85747c6ae99SBarry Smith if (i == m-1) lf[i] = want; 85847c6ae99SBarry Smith else { 8597aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 8607aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 8617aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 8627aca7175SJed Brown while ((startf+want)/ratio < nextc - stencil_width) want++; 8637aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 8647aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 8657aca7175SJed Brown while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 8667aca7175SJed Brown /* Make sure all constraints are satisfied */ 86730729d88SBarry Smith if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 868ce94432eSBarry 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"); 86947c6ae99SBarry Smith } 87047c6ae99SBarry Smith lf[i] = want; 87147c6ae99SBarry Smith startc += lc[i]; 87247c6ae99SBarry Smith startf += lf[i]; 87347c6ae99SBarry Smith remaining -= lf[i]; 87447c6ae99SBarry Smith } 87547c6ae99SBarry Smith PetscFunctionReturn(0); 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith 8782be375d4SJed Brown /* 8792be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 8802be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 8812be375d4SJed Brown 8822be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 8832be375d4SJed Brown */ 8842be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 8852be375d4SJed Brown { 8862be375d4SJed Brown PetscInt i,totalf,remaining,startc,startf; 8872be375d4SJed Brown PetscErrorCode ierr; 8882be375d4SJed Brown 8892be375d4SJed Brown PetscFunctionBegin; 890ce94432eSBarry Smith if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 8912be375d4SJed Brown if (ratio == 1) { 8922be375d4SJed Brown ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 8932be375d4SJed Brown PetscFunctionReturn(0); 8942be375d4SJed Brown } 8952be375d4SJed Brown for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 8962be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 8972be375d4SJed Brown for (i=0,startc=0,startf=0; i<m; i++) { 8982be375d4SJed Brown PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 8992be375d4SJed Brown if (i == m-1) lc[i] = want; 9002be375d4SJed Brown else { 9012be375d4SJed Brown const PetscInt nextf = startf+lf[i]; 9022be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 9032be375d4SJed Brown * node is within one stencil width. */ 9042be375d4SJed Brown while (nextf/ratio < startc+want-stencil_width) want--; 9052be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 9062be375d4SJed Brown * fine node is within one stencil width. */ 9072be375d4SJed Brown while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 9082be375d4SJed Brown if (want < 0 || want > remaining 909ce94432eSBarry 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"); 9102be375d4SJed Brown } 9112be375d4SJed Brown lc[i] = want; 9122be375d4SJed Brown startc += lc[i]; 9132be375d4SJed Brown startf += lf[i]; 9142be375d4SJed Brown remaining -= lc[i]; 9152be375d4SJed Brown } 9162be375d4SJed Brown PetscFunctionReturn(0); 9172be375d4SJed Brown } 9182be375d4SJed Brown 9197087cfbeSBarry Smith PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 92047c6ae99SBarry Smith { 92147c6ae99SBarry Smith PetscErrorCode ierr; 922c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 9239a42bb27SBarry Smith DM da2; 92447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 92547c6ae99SBarry Smith 92647c6ae99SBarry Smith PetscFunctionBegin; 92747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 92847c6ae99SBarry Smith PetscValidPointer(daref,3); 92947c6ae99SBarry Smith 930c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 931bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 93247c6ae99SBarry Smith M = dd->refine_x*dd->M; 93347c6ae99SBarry Smith } else { 93447c6ae99SBarry Smith M = 1 + dd->refine_x*(dd->M - 1); 93547c6ae99SBarry Smith } 936bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 937c73cfb54SMatthew G. Knepley if (dim > 1) { 93847c6ae99SBarry Smith N = dd->refine_y*dd->N; 93947c6ae99SBarry Smith } else { 9401860e6e9SBarry Smith N = 1; 9411860e6e9SBarry Smith } 9421860e6e9SBarry Smith } else { 94347c6ae99SBarry Smith N = 1 + dd->refine_y*(dd->N - 1); 94447c6ae99SBarry Smith } 945bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 946c73cfb54SMatthew G. Knepley if (dim > 2) { 94747c6ae99SBarry Smith P = dd->refine_z*dd->P; 94847c6ae99SBarry Smith } else { 9491860e6e9SBarry Smith P = 1; 9501860e6e9SBarry Smith } 9511860e6e9SBarry Smith } else { 95247c6ae99SBarry Smith P = 1 + dd->refine_z*(dd->P - 1); 95347c6ae99SBarry Smith } 954ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 955397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 956c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 957397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 958397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 959397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 960397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 961397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 962397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 963c73cfb54SMatthew G. Knepley if (dim == 3) { 96447c6ae99SBarry Smith PetscInt *lx,*ly,*lz; 965dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);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); 967bff4a2f0SMatthew 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); 968bff4a2f0SMatthew 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); 969397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 97047c6ae99SBarry Smith ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 971c73cfb54SMatthew G. Knepley } else if (dim == 2) { 97247c6ae99SBarry Smith PetscInt *lx,*ly; 973dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 974bff4a2f0SMatthew 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); 975bff4a2f0SMatthew 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); 9760298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 97747c6ae99SBarry Smith ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 978c73cfb54SMatthew G. Knepley } else if (dim == 1) { 97947c6ae99SBarry Smith PetscInt *lx; 980785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 981bff4a2f0SMatthew 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); 9820298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 98347c6ae99SBarry Smith ierr = PetscFree(lx);CHKERRQ(ierr); 98447c6ae99SBarry Smith } 98547c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 98647c6ae99SBarry Smith 98747c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 98825296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 98925296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 99047c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 99147c6ae99SBarry Smith dd2->interptype = dd->interptype; 99247c6ae99SBarry Smith 99347c6ae99SBarry Smith /* copy fill information if given */ 99447c6ae99SBarry Smith if (dd->dfill) { 995854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 99647c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 99747c6ae99SBarry Smith } 99847c6ae99SBarry Smith if (dd->ofill) { 999854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 100047c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 100147c6ae99SBarry Smith } 100247c6ae99SBarry Smith /* copy the refine information */ 1003397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1004397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1005397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 100647c6ae99SBarry Smith 1007897f7067SBarry Smith if (dd->refine_z_hier) { 1008897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) { 1009897f7067SBarry Smith dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; 1010897f7067SBarry Smith } 1011897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) { 1012897f7067SBarry Smith dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; 1013897f7067SBarry Smith } 1014897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1015897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1016897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1017897f7067SBarry Smith } 1018897f7067SBarry Smith if (dd->refine_y_hier) { 1019897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) { 1020897f7067SBarry Smith dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; 1021897f7067SBarry Smith } 1022897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) { 1023897f7067SBarry Smith dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; 1024897f7067SBarry Smith } 1025897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1026897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1027897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1028897f7067SBarry Smith } 1029897f7067SBarry Smith if (dd->refine_x_hier) { 1030897f7067SBarry Smith if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) { 1031897f7067SBarry Smith dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; 1032897f7067SBarry Smith } 1033897f7067SBarry Smith if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) { 1034897f7067SBarry Smith dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; 1035897f7067SBarry Smith } 1036897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1037897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1038897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1039897f7067SBarry Smith } 1040897f7067SBarry Smith 1041897f7067SBarry Smith 104247c6ae99SBarry Smith /* copy vector type information */ 1043c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 1044ddcf8b74SDave May 1045efd51863SBarry Smith dd2->lf = dd->lf; 1046efd51863SBarry Smith dd2->lj = dd->lj; 1047efd51863SBarry Smith 10486e87535bSJed Brown da2->leveldown = da->leveldown; 10496e87535bSJed Brown da2->levelup = da->levelup + 1; 10508865f1eaSKarl Rupp 10516e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 10526e87535bSJed Brown 1053ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10546636e97aSMatthew G Knepley if (da->coordinates) { 1055ddcf8b74SDave May DM cdaf,cdac; 1056ddcf8b74SDave May Vec coordsc,coordsf; 1057ddcf8b74SDave May Mat II; 1058ddcf8b74SDave May 10596636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 10606636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 10616636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 1062b61d3410SDave May /* force creation of the coordinate vector */ 1063b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 10646636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 10650298fd71SBarry Smith ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 1066ddcf8b74SDave May ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 1067fcfd50ebSBarry Smith ierr = MatDestroy(&II);CHKERRQ(ierr); 1068ddcf8b74SDave May } 1069397b6216SJed Brown 1070f3141302SJed Brown for (i=0; i<da->bs; i++) { 1071f3141302SJed Brown const char *fieldname; 1072f3141302SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1073f3141302SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1074f3141302SJed Brown } 1075397b6216SJed Brown 107647c6ae99SBarry Smith *daref = da2; 107747c6ae99SBarry Smith PetscFunctionReturn(0); 107847c6ae99SBarry Smith } 107947c6ae99SBarry Smith 1080397b6216SJed Brown 10817087cfbeSBarry Smith PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 108247c6ae99SBarry Smith { 108347c6ae99SBarry Smith PetscErrorCode ierr; 1084c73cfb54SMatthew G. Knepley PetscInt M,N,P,i,dim; 10859a42bb27SBarry Smith DM da2; 108647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data,*dd2; 108747c6ae99SBarry Smith 108847c6ae99SBarry Smith PetscFunctionBegin; 108947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 109047c6ae99SBarry Smith PetscValidPointer(daref,3); 109147c6ae99SBarry Smith 1092c73cfb54SMatthew G. Knepley ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 1093bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1094397b6216SJed Brown M = dd->M / dd->coarsen_x; 109547c6ae99SBarry Smith } else { 1096397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 109747c6ae99SBarry Smith } 1098bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1099c73cfb54SMatthew G. Knepley if (dim > 1) { 1100397b6216SJed Brown N = dd->N / dd->coarsen_y; 110147c6ae99SBarry Smith } else { 11021860e6e9SBarry Smith N = 1; 11031860e6e9SBarry Smith } 11041860e6e9SBarry Smith } else { 1105397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 110647c6ae99SBarry Smith } 1107bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1108c73cfb54SMatthew G. Knepley if (dim > 2) { 1109397b6216SJed Brown P = dd->P / dd->coarsen_z; 111047c6ae99SBarry Smith } else { 11111860e6e9SBarry Smith P = 1; 11121860e6e9SBarry Smith } 11131860e6e9SBarry Smith } else { 1114397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 111547c6ae99SBarry Smith } 1116ce94432eSBarry Smith ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 1117397b6216SJed Brown ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 1118c73cfb54SMatthew G. Knepley ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 1119397b6216SJed Brown ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 1120397b6216SJed Brown ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 1121397b6216SJed Brown ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 1122397b6216SJed Brown ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 1123397b6216SJed Brown ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 1124397b6216SJed Brown ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 1125c73cfb54SMatthew G. Knepley if (dim == 3) { 11262be375d4SJed Brown PetscInt *lx,*ly,*lz; 1127dcca6d9dSJed Brown ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 1128bff4a2f0SMatthew 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); 1129bff4a2f0SMatthew 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); 1130bff4a2f0SMatthew 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); 1131397b6216SJed Brown ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 11322be375d4SJed Brown ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 1133c73cfb54SMatthew G. Knepley } else if (dim == 2) { 11342be375d4SJed Brown PetscInt *lx,*ly; 1135dcca6d9dSJed Brown ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 1136bff4a2f0SMatthew 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); 1137bff4a2f0SMatthew 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); 11380298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 11392be375d4SJed Brown ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 1140c73cfb54SMatthew G. Knepley } else if (dim == 1) { 11412be375d4SJed Brown PetscInt *lx; 1142785e854fSJed Brown ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 1143bff4a2f0SMatthew 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); 11440298fd71SBarry Smith ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 11452be375d4SJed Brown ierr = PetscFree(lx);CHKERRQ(ierr); 114647c6ae99SBarry Smith } 114747c6ae99SBarry Smith dd2 = (DM_DA*)da2->data; 114847c6ae99SBarry Smith 11494dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 115025296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 115125296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 115247c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 115347c6ae99SBarry Smith dd2->interptype = dd->interptype; 115447c6ae99SBarry Smith 115547c6ae99SBarry Smith /* copy fill information if given */ 115647c6ae99SBarry Smith if (dd->dfill) { 1157854ce69bSBarry Smith ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 115847c6ae99SBarry Smith ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 115947c6ae99SBarry Smith } 116047c6ae99SBarry Smith if (dd->ofill) { 1161854ce69bSBarry Smith ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 116247c6ae99SBarry Smith ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 116347c6ae99SBarry Smith } 116447c6ae99SBarry Smith /* copy the refine information */ 1165397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1166397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1167397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 116847c6ae99SBarry Smith 1169897f7067SBarry Smith if (dd->refine_z_hier) { 1170897f7067SBarry Smith if (da->levelup - da->leveldown -1 > -1 && da->levelup - da->leveldown - 1< dd->refine_z_hier_n) { 1171897f7067SBarry Smith dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown - 1]; 1172897f7067SBarry Smith } 1173897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_z_hier_n) { 1174897f7067SBarry Smith dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown - 2]; 1175897f7067SBarry Smith } 1176897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 1177897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1178897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1179897f7067SBarry Smith } 1180897f7067SBarry Smith if (dd->refine_y_hier) { 1181897f7067SBarry Smith if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1< dd->refine_y_hier_n) { 1182897f7067SBarry Smith dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown - 1]; 1183897f7067SBarry Smith } 1184897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_y_hier_n) { 1185897f7067SBarry Smith dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown - 2]; 1186897f7067SBarry Smith } 1187897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 1188897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1189897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1190897f7067SBarry Smith } 1191897f7067SBarry Smith if (dd->refine_x_hier) { 1192897f7067SBarry Smith if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1 < dd->refine_x_hier_n) { 1193897f7067SBarry Smith dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown - 1]; 1194897f7067SBarry Smith } 1195897f7067SBarry Smith if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_x_hier_n) { 1196897f7067SBarry Smith dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown - 2]; 1197897f7067SBarry Smith } 1198897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 1199897f7067SBarry Smith ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1200897f7067SBarry Smith ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1201897f7067SBarry Smith } 1202897f7067SBarry Smith 120347c6ae99SBarry Smith /* copy vector type information */ 1204c0dedaeaSBarry Smith ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 120547c6ae99SBarry Smith 1206644e2e5bSBarry Smith dd2->lf = dd->lf; 1207644e2e5bSBarry Smith dd2->lj = dd->lj; 1208644e2e5bSBarry Smith 12096e87535bSJed Brown da2->leveldown = da->leveldown + 1; 12106e87535bSJed Brown da2->levelup = da->levelup; 12118865f1eaSKarl Rupp 12126e87535bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 12136e87535bSJed Brown 1214ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 12156636e97aSMatthew G Knepley if (da->coordinates) { 1216ddcf8b74SDave May DM cdaf,cdac; 1217ddcf8b74SDave May Vec coordsc,coordsf; 12186dbf9973SLawrence Mitchell Mat inject; 12196dbf9973SLawrence Mitchell VecScatter vscat; 1220ddcf8b74SDave May 12216636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); 12226636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); 12236636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); 1224b61d3410SDave May /* force creation of the coordinate vector */ 1225b61d3410SDave May ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 12266636e97aSMatthew G Knepley ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 1227ddcf8b74SDave May 1228e727c939SJed Brown ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 12296dbf9973SLawrence Mitchell ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr); 12306dbf9973SLawrence Mitchell ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12316dbf9973SLawrence Mitchell ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12326dbf9973SLawrence Mitchell ierr = MatDestroy(&inject);CHKERRQ(ierr); 1233ddcf8b74SDave May } 1234f98405f7SJed Brown 1235f98405f7SJed Brown for (i=0; i<da->bs; i++) { 1236f98405f7SJed Brown const char *fieldname; 1237f98405f7SJed Brown ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1238f98405f7SJed Brown ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1239f98405f7SJed Brown } 1240f98405f7SJed Brown 124147c6ae99SBarry Smith *daref = da2; 124247c6ae99SBarry Smith PetscFunctionReturn(0); 124347c6ae99SBarry Smith } 124447c6ae99SBarry Smith 12457087cfbeSBarry Smith PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 124647c6ae99SBarry Smith { 124747c6ae99SBarry Smith PetscErrorCode ierr; 124847c6ae99SBarry Smith PetscInt i,n,*refx,*refy,*refz; 124947c6ae99SBarry Smith 125047c6ae99SBarry Smith PetscFunctionBegin; 125147c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1252ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 125347c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 125447c6ae99SBarry Smith PetscValidPointer(daf,3); 125547c6ae99SBarry Smith 1256aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 1257dcca6d9dSJed Brown ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr); 125847c6ae99SBarry Smith for (i=0; i<nlevels; i++) { 1259aa219208SBarry Smith ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith n = nlevels; 1262c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 126347c6ae99SBarry Smith n = nlevels; 1264c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 126547c6ae99SBarry Smith n = nlevels; 1266c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 126747c6ae99SBarry Smith 1268aa219208SBarry Smith ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1269ce94432eSBarry Smith ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 127047c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1271aa219208SBarry Smith ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1272ce94432eSBarry Smith ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 127347c6ae99SBarry Smith } 127447c6ae99SBarry Smith ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 127547c6ae99SBarry Smith PetscFunctionReturn(0); 127647c6ae99SBarry Smith } 127747c6ae99SBarry Smith 12787087cfbeSBarry Smith PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 127947c6ae99SBarry Smith { 128047c6ae99SBarry Smith PetscErrorCode ierr; 128147c6ae99SBarry Smith PetscInt i; 128247c6ae99SBarry Smith 128347c6ae99SBarry Smith PetscFunctionBegin; 128447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 1285ce94432eSBarry Smith if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 128647c6ae99SBarry Smith if (nlevels == 0) PetscFunctionReturn(0); 128747c6ae99SBarry Smith PetscValidPointer(dac,3); 1288ce94432eSBarry Smith ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 128947c6ae99SBarry Smith for (i=1; i<nlevels; i++) { 1290ce94432eSBarry Smith ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 129147c6ae99SBarry Smith } 129247c6ae99SBarry Smith PetscFunctionReturn(0); 129347c6ae99SBarry Smith } 129462197512SBarry Smith 129562197512SBarry Smith #include <petscgll.h> 129662197512SBarry Smith 129762197512SBarry Smith PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscGLL *gll) 129862197512SBarry Smith { 129962197512SBarry Smith PetscErrorCode ierr; 130062197512SBarry Smith PetscInt i,j,n = gll->n,xs,xn,q; 130162197512SBarry Smith PetscScalar *xx; 130262197512SBarry Smith PetscReal h; 130362197512SBarry Smith Vec x; 130462197512SBarry Smith DM_DA *da = (DM_DA*)dm->data; 130562197512SBarry Smith 130662197512SBarry Smith PetscFunctionBegin; 130762197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 130862197512SBarry Smith ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 130962197512SBarry Smith q = (q-1)/(n-1); /* number of spectral elements */ 131062197512SBarry Smith h = 2.0/q; 131162197512SBarry Smith ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 131262197512SBarry Smith xs = xs/(n-1); 131362197512SBarry Smith xn = xn/(n-1); 131462197512SBarry Smith ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr); 131562197512SBarry Smith ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 131662197512SBarry Smith ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr); 131762197512SBarry Smith 131862197512SBarry Smith /* loop over local spectral elements */ 131962197512SBarry Smith for (j=xs; j<xs+xn; j++) { 132062197512SBarry Smith /* 132162197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 132262197512SBarry Smith */ 132362197512SBarry Smith for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) { 132462197512SBarry Smith xx[j*(n-1) + i] = -1.0 + h*j + h*(gll->nodes[i]+1.0)/2.; 132562197512SBarry Smith } 132662197512SBarry Smith } 132762197512SBarry Smith ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr); 132862197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic"); 132962197512SBarry Smith PetscFunctionReturn(0); 133062197512SBarry Smith } 133162197512SBarry Smith 13321fd49c25SBarry Smith /*@ 133362197512SBarry Smith 133462197512SBarry 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 133562197512SBarry Smith 133662197512SBarry Smith Collective on DM 133762197512SBarry Smith 133862197512SBarry Smith Input Parameters: 133962197512SBarry Smith + da - the DMDA object 134062197512SBarry Smith - gll - the GLL object 134162197512SBarry Smith 134262197512SBarry Smith Notes: the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 134362197512SBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is 134462197512SBarry Smith periodic or not. 134562197512SBarry Smith 1346edc382c3SSatish Balay Level: advanced 1347edc382c3SSatish Balay 134862197512SBarry Smith .seealso: DMDACreate(), PetscGLLCreate(), DMGetCoordinates() 134962197512SBarry Smith @*/ 135062197512SBarry Smith PetscErrorCode DMDASetGLLCoordinates(DM da,PetscGLL *gll) 135162197512SBarry Smith { 135262197512SBarry Smith PetscErrorCode ierr; 135362197512SBarry Smith 135462197512SBarry Smith PetscFunctionBegin; 135562197512SBarry Smith if (da->dim == 1) { 135662197512SBarry Smith ierr = DMDASetGLLCoordinates_1d(da,gll);CHKERRQ(ierr); 135762197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d"); 135862197512SBarry Smith PetscFunctionReturn(0); 135962197512SBarry Smith } 1360