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 6*20f4b53cSBarry Smith Logically Collective 747c6ae99SBarry Smith 847c6ae99SBarry Smith Input Parameters: 9dce8aebaSBarry 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 16dce8aebaSBarry Smith Developer Note: 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 19dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `PetscSplitOwnership()` 2047c6ae99SBarry Smith @*/ 21d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 22d71ae5a4SJacob Faibussowitsch { 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); 307a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 317a8be351SBarry Smith PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in X direction must be positive"); 327a8be351SBarry Smith PetscCheck(N >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in Y direction must be positive"); 337a8be351SBarry 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; 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3947c6ae99SBarry Smith } 4047c6ae99SBarry Smith 4147c6ae99SBarry Smith /*@ 42aa219208SBarry Smith DMDASetNumProcs - Sets the number of processes in each dimension 4347c6ae99SBarry Smith 44*20f4b53cSBarry Smith Logically Collective 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith Input Parameters: 47dce8aebaSBarry Smith + da - the `DMDA` 48dce8aebaSBarry Smith . m - the number of X procs (or `PETSC_DECIDE`) 49dce8aebaSBarry Smith . n - the number of Y procs (or `PETSC_DECIDE`) 50dce8aebaSBarry Smith - p - the number of Z procs (or `PETSC_DECIDE`) 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith Level: intermediate 5347c6ae99SBarry Smith 54dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetSizes()`, `PetscSplitOwnership()` 5547c6ae99SBarry Smith @*/ 56d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 57d71ae5a4SJacob Faibussowitsch { 5847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 5947c6ae99SBarry Smith 6047c6ae99SBarry Smith PetscFunctionBegin; 61a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 6247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, m, 2); 6347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, n, 3); 6447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, p, 4); 657a8be351SBarry Smith PetscCheck(!da->setupcalled, 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; 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size)); 72e3f3e4b6SBarry Smith if ((dd->m > 0) && (dd->n < 0)) { 73e3f3e4b6SBarry Smith dd->n = size / dd->m; 7463a3b9bcSJacob Faibussowitsch PetscCheck(dd->n * dd->m == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " 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; 7863a3b9bcSJacob Faibussowitsch PetscCheck(dd->n * dd->m == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " processes in Y direction not divisible into comm size %d", n, size); 79e3f3e4b6SBarry Smith } 80e3f3e4b6SBarry Smith } 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith 8447c6ae99SBarry Smith /*@ 851321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 8647c6ae99SBarry Smith 87*20f4b53cSBarry Smith Not Collective 8847c6ae99SBarry Smith 89d8d19677SJose E. Roman Input Parameters: 90dce8aebaSBarry Smith + da - The `DMDA` 91dce8aebaSBarry Smith - bx,by,bz - One of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith Level: intermediate 9447c6ae99SBarry Smith 95dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType` 9647c6ae99SBarry Smith @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz) 98d71ae5a4SJacob Faibussowitsch { 9947c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith PetscFunctionBegin; 102a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1035a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bx, 2); 1045a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, by, 3); 1055a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bz, 4); 1067a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 1071321219cSEthan Coon dd->bx = bx; 1081321219cSEthan Coon dd->by = by; 1091321219cSEthan Coon dd->bz = bz; 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11147c6ae99SBarry Smith } 11247c6ae99SBarry Smith 11347c6ae99SBarry Smith /*@ 114aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 11547c6ae99SBarry Smith 116*20f4b53cSBarry Smith Not Collective 11747c6ae99SBarry Smith 11859f3ab6dSMatthew G. Knepley Input Parameters: 119dce8aebaSBarry Smith + da - The `DMDA` 12047c6ae99SBarry Smith - dof - Number of degrees of freedom 12147c6ae99SBarry Smith 12247c6ae99SBarry Smith Level: intermediate 12347c6ae99SBarry Smith 124dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()` 12547c6ae99SBarry Smith @*/ 126d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetDof(DM da, PetscInt dof) 127d71ae5a4SJacob Faibussowitsch { 12847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith PetscFunctionBegin; 131a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 13254cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da, dof, 2); 1337a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 13447c6ae99SBarry Smith dd->w = dof; 1351411c6eeSJed Brown da->bs = dof; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13747c6ae99SBarry Smith } 13847c6ae99SBarry Smith 139fb6725baSMatthew G. Knepley /*@ 140fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 141fb6725baSMatthew G. Knepley 142*20f4b53cSBarry Smith Not Collective 143fb6725baSMatthew G. Knepley 144fb6725baSMatthew G. Knepley Input Parameter: 145dce8aebaSBarry Smith . da - The `DMDA` 146fb6725baSMatthew G. Knepley 147fb6725baSMatthew G. Knepley Output Parameter: 148fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom 149fb6725baSMatthew G. Knepley 150fb6725baSMatthew G. Knepley Level: intermediate 151fb6725baSMatthew G. Knepley 152dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetDof()`, `DMDACreate()`, `DMDestroy()` 153fb6725baSMatthew G. Knepley @*/ 154d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 155d71ae5a4SJacob Faibussowitsch { 156fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 157fb6725baSMatthew G. Knepley 158fb6725baSMatthew G. Knepley PetscFunctionBegin; 159a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 160dadcf809SJacob Faibussowitsch PetscValidIntPointer(dof, 2); 161fb6725baSMatthew G. Knepley *dof = dd->w; 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163fb6725baSMatthew G. Knepley } 164fb6725baSMatthew G. Knepley 1657ddda789SPeter Brune /*@ 1667ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1677ddda789SPeter Brune 168*20f4b53cSBarry Smith Not Collective 1697ddda789SPeter Brune 170f899ff85SJose E. Roman Input Parameter: 171dce8aebaSBarry Smith . da - The `DMDA` 1727ddda789SPeter Brune 1737ddda789SPeter Brune Output Parameters: 1747ddda789SPeter Brune + x - Overlap in the x direction 1757ddda789SPeter Brune . y - Overlap in the y direction 1767ddda789SPeter Brune - z - Overlap in the z direction 1777ddda789SPeter Brune 1787ddda789SPeter Brune Level: intermediate 1797ddda789SPeter Brune 180dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetOverlap()` 1817ddda789SPeter Brune @*/ 182d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOverlap(DM da, PetscInt *x, PetscInt *y, PetscInt *z) 183d71ae5a4SJacob Faibussowitsch { 1847ddda789SPeter Brune DM_DA *dd = (DM_DA *)da->data; 1857ddda789SPeter Brune 1867ddda789SPeter Brune PetscFunctionBegin; 187a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1887ddda789SPeter Brune if (x) *x = dd->xol; 1897ddda789SPeter Brune if (y) *y = dd->yol; 1907ddda789SPeter Brune if (z) *z = dd->zol; 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1927ddda789SPeter Brune } 1937ddda789SPeter Brune 19488661749SPeter Brune /*@ 19588661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 19688661749SPeter Brune 197*20f4b53cSBarry Smith Not Collective 19888661749SPeter Brune 1997ddda789SPeter Brune Input Parameters: 200dce8aebaSBarry Smith + da - The `DMDA` 2017ddda789SPeter Brune . x - Overlap in the x direction 2027ddda789SPeter Brune . y - Overlap in the y direction 2037ddda789SPeter Brune - z - Overlap in the z direction 20488661749SPeter Brune 20588661749SPeter Brune Level: intermediate 20688661749SPeter Brune 207dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetOverlap()` 20888661749SPeter Brune @*/ 209d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z) 210d71ae5a4SJacob Faibussowitsch { 21188661749SPeter Brune DM_DA *dd = (DM_DA *)da->data; 21288661749SPeter Brune 21388661749SPeter Brune PetscFunctionBegin; 214a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2157ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, x, 2); 2167ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, y, 3); 2177ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, z, 4); 2187ddda789SPeter Brune dd->xol = x; 2197ddda789SPeter Brune dd->yol = y; 2207ddda789SPeter Brune dd->zol = z; 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22288661749SPeter Brune } 22388661749SPeter Brune 2243e7870d2SPeter Brune /*@ 2253e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2263e7870d2SPeter Brune 227*20f4b53cSBarry Smith Not Collective 2283e7870d2SPeter Brune 2293e7870d2SPeter Brune Input Parameters: 230dce8aebaSBarry Smith . da - The `DMDA` 2313e7870d2SPeter Brune 2323e7870d2SPeter Brune Output Parameters: 233a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition 2343e7870d2SPeter Brune 2353e7870d2SPeter Brune Level: intermediate 2363e7870d2SPeter Brune 237dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()` 2383e7870d2SPeter Brune @*/ 239d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub) 240d71ae5a4SJacob Faibussowitsch { 2413e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2423e7870d2SPeter Brune 2433e7870d2SPeter Brune PetscFunctionBegin; 244a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2453e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2473e7870d2SPeter Brune } 2483e7870d2SPeter Brune 2493e7870d2SPeter Brune /*@ 2503e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2513e7870d2SPeter Brune 252*20f4b53cSBarry Smith Not Collective 2533e7870d2SPeter Brune 2543e7870d2SPeter Brune Input Parameters: 255dce8aebaSBarry Smith + da - The `DMDA` 2563e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2573e7870d2SPeter Brune 2583e7870d2SPeter Brune Level: intermediate 2593e7870d2SPeter Brune 260dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()` 2613e7870d2SPeter Brune @*/ 262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub) 263d71ae5a4SJacob Faibussowitsch { 2643e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2653e7870d2SPeter Brune 2663e7870d2SPeter Brune PetscFunctionBegin; 267a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2683e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da, Nsub, 2); 2693e7870d2SPeter Brune dd->Nsub = Nsub; 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2713e7870d2SPeter Brune } 2723e7870d2SPeter Brune 273d886c4f4SPeter Brune /*@ 274d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 275d886c4f4SPeter Brune 276*20f4b53cSBarry Smith Collective 277d886c4f4SPeter Brune 278d8d19677SJose E. Roman Input Parameters: 279dce8aebaSBarry Smith + da - The `DMDA` 280d886c4f4SPeter Brune . xo - The offset in the x direction 281d886c4f4SPeter Brune . yo - The offset in the y direction 282d886c4f4SPeter Brune - zo - The offset in the z direction 283d886c4f4SPeter Brune 284d886c4f4SPeter Brune Level: intermediate 285d886c4f4SPeter Brune 286dce8aebaSBarry Smith Note: 287dce8aebaSBarry Smith This is used primarily to overlap a computation on a local `DMDA` with that on a global `DMDA` without 288d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 289d886c4f4SPeter Brune 290dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 291d886c4f4SPeter Brune @*/ 292d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 293d71ae5a4SJacob Faibussowitsch { 294d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 295d886c4f4SPeter Brune 296d886c4f4SPeter Brune PetscFunctionBegin; 297a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 298d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da, xo, 2); 29995c13181SPeter Brune PetscValidLogicalCollectiveInt(da, yo, 3); 30095c13181SPeter Brune PetscValidLogicalCollectiveInt(da, zo, 4); 30195c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Mo, 5); 30295c13181SPeter Brune PetscValidLogicalCollectiveInt(da, No, 6); 30395c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Po, 7); 304d886c4f4SPeter Brune dd->xo = xo; 305d886c4f4SPeter Brune dd->yo = yo; 306d886c4f4SPeter Brune dd->zo = zo; 30795c13181SPeter Brune dd->Mo = Mo; 30895c13181SPeter Brune dd->No = No; 30995c13181SPeter Brune dd->Po = Po; 31095c13181SPeter Brune 3116858538eSMatthew G. Knepley if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po)); 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313d886c4f4SPeter Brune } 314d886c4f4SPeter Brune 315d886c4f4SPeter Brune /*@ 316dce8aebaSBarry Smith DMDAGetOffset - Gets the index offset of the `DMDA`. 317d886c4f4SPeter Brune 318*20f4b53cSBarry Smith Not Collective 319d886c4f4SPeter Brune 320d886c4f4SPeter Brune Input Parameter: 321dce8aebaSBarry Smith . da - The `DMDA` 322d886c4f4SPeter Brune 323d886c4f4SPeter Brune Output Parameters: 324d886c4f4SPeter Brune + xo - The offset in the x direction 325d886c4f4SPeter Brune . yo - The offset in the y direction 32695c13181SPeter Brune . zo - The offset in the z direction 32795c13181SPeter Brune . Mo - The global size in the x direction 32895c13181SPeter Brune . No - The global size in the y direction 32995c13181SPeter Brune - Po - The global size in the z direction 330d886c4f4SPeter Brune 331d886c4f4SPeter Brune Level: intermediate 332d886c4f4SPeter Brune 333dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetOffset()`, `DMDAVecGetArray()` 334d886c4f4SPeter Brune @*/ 335d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo, PetscInt *Mo, PetscInt *No, PetscInt *Po) 336d71ae5a4SJacob Faibussowitsch { 337d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 338d886c4f4SPeter Brune 339d886c4f4SPeter Brune PetscFunctionBegin; 340a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 341d886c4f4SPeter Brune if (xo) *xo = dd->xo; 342d886c4f4SPeter Brune if (yo) *yo = dd->yo; 343d886c4f4SPeter Brune if (zo) *zo = dd->zo; 34495c13181SPeter Brune if (Mo) *Mo = dd->Mo; 34595c13181SPeter Brune if (No) *No = dd->No; 34695c13181SPeter Brune if (Po) *Po = dd->Po; 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348d886c4f4SPeter Brune } 349d886c4f4SPeter Brune 35040234c92SPeter Brune /*@ 351dce8aebaSBarry Smith DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain `DM`. 35240234c92SPeter Brune 353*20f4b53cSBarry Smith Not Collective 35440234c92SPeter Brune 35540234c92SPeter Brune Input Parameter: 356dce8aebaSBarry Smith . da - The `DMDA` 35740234c92SPeter Brune 35840234c92SPeter Brune Output Parameters: 35940234c92SPeter Brune + xs - The start of the region in x 36040234c92SPeter Brune . ys - The start of the region in y 36140234c92SPeter Brune . zs - The start of the region in z 36240234c92SPeter Brune . xs - The size of the region in x 36340234c92SPeter Brune . ys - The size of the region in y 364a2b725a8SWilliam Gropp - zs - The size of the region in z 36540234c92SPeter Brune 36640234c92SPeter Brune Level: intermediate 36740234c92SPeter Brune 368dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 36940234c92SPeter Brune @*/ 370d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 371d71ae5a4SJacob Faibussowitsch { 37240234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 37340234c92SPeter Brune 37440234c92SPeter Brune PetscFunctionBegin; 375a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 37640234c92SPeter Brune if (xs) *xs = dd->nonxs; 37740234c92SPeter Brune if (ys) *ys = dd->nonys; 37840234c92SPeter Brune if (zs) *zs = dd->nonzs; 37940234c92SPeter Brune if (xm) *xm = dd->nonxm; 38040234c92SPeter Brune if (ym) *ym = dd->nonym; 38140234c92SPeter Brune if (zm) *zm = dd->nonzm; 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38340234c92SPeter Brune } 38440234c92SPeter Brune 38540234c92SPeter Brune /*@ 386dce8aebaSBarry Smith DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain `DM`. 38740234c92SPeter Brune 388*20f4b53cSBarry Smith Collective 38940234c92SPeter Brune 390d8d19677SJose E. Roman Input Parameters: 391dce8aebaSBarry Smith + da - The `DMDA` 39240234c92SPeter Brune . xs - The start of the region in x 39340234c92SPeter Brune . ys - The start of the region in y 39440234c92SPeter Brune . zs - The start of the region in z 39540234c92SPeter Brune . xs - The size of the region in x 39640234c92SPeter Brune . ys - The size of the region in y 397a2b725a8SWilliam Gropp - zs - The size of the region in z 39840234c92SPeter Brune 39940234c92SPeter Brune Level: intermediate 40040234c92SPeter Brune 401dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 40240234c92SPeter Brune @*/ 403d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 404d71ae5a4SJacob Faibussowitsch { 40540234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 40640234c92SPeter Brune 40740234c92SPeter Brune PetscFunctionBegin; 408a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 40940234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xs, 2); 41040234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ys, 3); 41140234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zs, 4); 41240234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xm, 5); 41340234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ym, 6); 41440234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zm, 7); 41540234c92SPeter Brune dd->nonxs = xs; 41640234c92SPeter Brune dd->nonys = ys; 41740234c92SPeter Brune dd->nonzs = zs; 41840234c92SPeter Brune dd->nonxm = xm; 41940234c92SPeter Brune dd->nonym = ym; 42040234c92SPeter Brune dd->nonzm = zm; 42140234c92SPeter Brune 4223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42340234c92SPeter Brune } 42488661749SPeter Brune 42547c6ae99SBarry Smith /*@ 426aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 42747c6ae99SBarry Smith 428*20f4b53cSBarry Smith Logically Collective 42947c6ae99SBarry Smith 430d8d19677SJose E. Roman Input Parameters: 431dce8aebaSBarry Smith + da - The `DMDA` 432dce8aebaSBarry Smith - stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 43347c6ae99SBarry Smith 43447c6ae99SBarry Smith Level: intermediate 43547c6ae99SBarry Smith 436dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 43747c6ae99SBarry Smith @*/ 438d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 439d71ae5a4SJacob Faibussowitsch { 44047c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 44147c6ae99SBarry Smith 44247c6ae99SBarry Smith PetscFunctionBegin; 443a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 44447c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, stype, 2); 4457a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 44647c6ae99SBarry Smith dd->stencil_type = stype; 4473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44847c6ae99SBarry Smith } 44947c6ae99SBarry Smith 450fb6725baSMatthew G. Knepley /*@ 451fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 452fb6725baSMatthew G. Knepley 453*20f4b53cSBarry Smith Not Collective 454fb6725baSMatthew G. Knepley 455fb6725baSMatthew G. Knepley Input Parameter: 456dce8aebaSBarry Smith . da - The `DMDA` 457fb6725baSMatthew G. Knepley 458fb6725baSMatthew G. Knepley Output Parameter: 459dce8aebaSBarry Smith . stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 460fb6725baSMatthew G. Knepley 461fb6725baSMatthew G. Knepley Level: intermediate 462fb6725baSMatthew G. Knepley 463dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 464fb6725baSMatthew G. Knepley @*/ 465d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 466d71ae5a4SJacob Faibussowitsch { 467fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 468fb6725baSMatthew G. Knepley 469fb6725baSMatthew G. Knepley PetscFunctionBegin; 470a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 471fb6725baSMatthew G. Knepley PetscValidPointer(stype, 2); 472fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474fb6725baSMatthew G. Knepley } 475fb6725baSMatthew G. Knepley 47647c6ae99SBarry Smith /*@ 477aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 47847c6ae99SBarry Smith 479*20f4b53cSBarry Smith Logically Collective 48047c6ae99SBarry Smith 481d8d19677SJose E. Roman Input Parameters: 482dce8aebaSBarry Smith + da - The `DMDA` 48347c6ae99SBarry Smith - width - The stencil width 48447c6ae99SBarry Smith 48547c6ae99SBarry Smith Level: intermediate 48647c6ae99SBarry Smith 487dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 48847c6ae99SBarry Smith @*/ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 490d71ae5a4SJacob Faibussowitsch { 49147c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith PetscFunctionBegin; 494a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 49547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, width, 2); 4967a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 49747c6ae99SBarry Smith dd->s = width; 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith 501fb6725baSMatthew G. Knepley /*@ 502fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 503fb6725baSMatthew G. Knepley 504*20f4b53cSBarry Smith Not Collective 505fb6725baSMatthew G. Knepley 506fb6725baSMatthew G. Knepley Input Parameter: 507dce8aebaSBarry Smith . da - The `DMDA` 508fb6725baSMatthew G. Knepley 509fb6725baSMatthew G. Knepley Output Parameter: 510fb6725baSMatthew G. Knepley . width - The stencil width 511fb6725baSMatthew G. Knepley 512fb6725baSMatthew G. Knepley Level: intermediate 513fb6725baSMatthew G. Knepley 514dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 515fb6725baSMatthew G. Knepley @*/ 516d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 517d71ae5a4SJacob Faibussowitsch { 518fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 519fb6725baSMatthew G. Knepley 520fb6725baSMatthew G. Knepley PetscFunctionBegin; 521a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 522dadcf809SJacob Faibussowitsch PetscValidIntPointer(width, 2); 523fb6725baSMatthew G. Knepley *width = dd->s; 5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 525fb6725baSMatthew G. Knepley } 526fb6725baSMatthew G. Knepley 527d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[]) 528d71ae5a4SJacob Faibussowitsch { 52947c6ae99SBarry Smith PetscInt i, sum; 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith PetscFunctionBegin; 5327a8be351SBarry Smith PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set"); 53347c6ae99SBarry Smith for (i = sum = 0; i < m; i++) sum += lx[i]; 53463a3b9bcSJacob Faibussowitsch PetscCheck(sum == M, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "Ownership ranges sum to %" PetscInt_FMT " but global dimension is %" PetscInt_FMT, sum, M); 5353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53647c6ae99SBarry Smith } 53747c6ae99SBarry Smith 53847c6ae99SBarry Smith /*@ 539aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 54047c6ae99SBarry Smith 541*20f4b53cSBarry Smith Logically Collective 54247c6ae99SBarry Smith 543d8d19677SJose E. Roman Input Parameters: 544dce8aebaSBarry Smith + da - The `DMDA` 5450298fd71SBarry 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 5460298fd71SBarry 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 5470298fd71SBarry 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. 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith Level: intermediate 55047c6ae99SBarry Smith 551e3f3e4b6SBarry Smith Note: these numbers are NOT multiplied by the number of dof per node. 552e3f3e4b6SBarry Smith 553dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 55447c6ae99SBarry Smith @*/ 555d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 556d71ae5a4SJacob Faibussowitsch { 55747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 55847c6ae99SBarry Smith 55947c6ae99SBarry Smith PetscFunctionBegin; 560a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5617a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 56247c6ae99SBarry Smith if (lx) { 5637a8be351SBarry Smith PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5649566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx)); 56548a46eb9SPierre Jolivet if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx)); 5669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lx, lx, dd->m)); 56747c6ae99SBarry Smith } 56847c6ae99SBarry Smith if (ly) { 5697a8be351SBarry Smith PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5709566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly)); 57148a46eb9SPierre Jolivet if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly)); 5729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->ly, ly, dd->n)); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith if (lz) { 5757a8be351SBarry Smith PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5769566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz)); 57748a46eb9SPierre Jolivet if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz)); 5789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lz, lz, dd->p)); 57947c6ae99SBarry Smith } 5803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58147c6ae99SBarry Smith } 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith /*@ 584aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 585dce8aebaSBarry Smith returned by `DMCreateInterpolation()` 58647c6ae99SBarry Smith 587*20f4b53cSBarry Smith Logically Collective 58847c6ae99SBarry Smith 589d8d19677SJose E. Roman Input Parameters: 59047c6ae99SBarry Smith + da - initial distributed array 591dce8aebaSBarry Smith - ctype - `DMDA_Q1` and `DMDA_Q0` are currently the only supported forms 59247c6ae99SBarry Smith 59347c6ae99SBarry Smith Level: intermediate 59447c6ae99SBarry Smith 595dce8aebaSBarry Smith Note: 596dce8aebaSBarry Smith You should call this on the coarser of the two `DMDA` you pass to `DMCreateInterpolation()` 59747c6ae99SBarry Smith 598dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDAInterpolationType` 59947c6ae99SBarry Smith @*/ 600d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype) 601d71ae5a4SJacob Faibussowitsch { 60247c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 60347c6ae99SBarry Smith 60447c6ae99SBarry Smith PetscFunctionBegin; 605a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 60647c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, ctype, 2); 60747c6ae99SBarry Smith dd->interptype = ctype; 6083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60947c6ae99SBarry Smith } 61047c6ae99SBarry Smith 6112dde6fd4SLisandro Dalcin /*@ 6122dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 613dce8aebaSBarry Smith used by `DMCreateInterpolation()` 6142dde6fd4SLisandro Dalcin 6152dde6fd4SLisandro Dalcin Not Collective 6162dde6fd4SLisandro Dalcin 6172dde6fd4SLisandro Dalcin Input Parameter: 6182dde6fd4SLisandro Dalcin . da - distributed array 6192dde6fd4SLisandro Dalcin 6202dde6fd4SLisandro Dalcin Output Parameter: 621dce8aebaSBarry Smith . ctype - interpolation type (`DMDA_Q1` and `DMDA_Q0` are currently the only supported forms) 6222dde6fd4SLisandro Dalcin 6232dde6fd4SLisandro Dalcin Level: intermediate 6242dde6fd4SLisandro Dalcin 625dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()` 6262dde6fd4SLisandro Dalcin @*/ 627d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype) 628d71ae5a4SJacob Faibussowitsch { 6292dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 6302dde6fd4SLisandro Dalcin 6312dde6fd4SLisandro Dalcin PetscFunctionBegin; 632a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 6332dde6fd4SLisandro Dalcin PetscValidPointer(ctype, 2); 6342dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6362dde6fd4SLisandro Dalcin } 63747c6ae99SBarry Smith 6386a119db4SBarry Smith /*@C 639aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 64047c6ae99SBarry Smith processes neighbors. 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith Not Collective 64347c6ae99SBarry Smith 64447c6ae99SBarry Smith Input Parameter: 645dce8aebaSBarry Smith . da - the `DMDA` object 64647c6ae99SBarry Smith 64747c6ae99SBarry Smith Output Parameters: 64847c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 64947c6ae99SBarry Smith this process itself is in the list 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith Level: intermediate 65247c6ae99SBarry Smith 653dce8aebaSBarry Smith Notes: 654dce8aebaSBarry Smith In 2d the array is of length 9, in 3d of length 27 655dce8aebaSBarry Smith 656dce8aebaSBarry Smith Not supported in 1d 657dce8aebaSBarry Smith 658dce8aebaSBarry Smith Do not free the array, it is freed when the `DMDA` is destroyed. 659dce8aebaSBarry Smith 660dce8aebaSBarry Smith Fortran Note: 661dce8aebaSBarry Smith In fortran you must pass in an array of the appropriate length. 662dce8aebaSBarry Smith 663dce8aebaSBarry Smith .seealso: `DMDA`, `DM` 66447c6ae99SBarry Smith @*/ 665d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[]) 666d71ae5a4SJacob Faibussowitsch { 66747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 6685fd66863SKarl Rupp 66947c6ae99SBarry Smith PetscFunctionBegin; 670a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 67147c6ae99SBarry Smith *ranks = dd->neighbors; 6723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67347c6ae99SBarry Smith } 67447c6ae99SBarry Smith 67547c6ae99SBarry Smith /*@C 676aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith Not Collective 67947c6ae99SBarry Smith 68047c6ae99SBarry Smith Input Parameter: 681dce8aebaSBarry Smith . da - the `DMDA` object 68247c6ae99SBarry Smith 683d8d19677SJose E. Roman Output Parameters: 68447c6ae99SBarry Smith + lx - ownership along x direction (optional) 68547c6ae99SBarry Smith . ly - ownership along y direction (optional) 68647c6ae99SBarry Smith - lz - ownership along z direction (optional) 68747c6ae99SBarry Smith 68847c6ae99SBarry Smith Level: intermediate 68947c6ae99SBarry Smith 690dce8aebaSBarry Smith Note: 691dce8aebaSBarry Smith These correspond to the optional final arguments passed to `DMDACreate()`, `DMDACreate2d()`, `DMDACreate3d()` 69247c6ae99SBarry Smith 69347c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 694dce8aebaSBarry Smith `DMDA` they came from still exists (has not been destroyed). 69547c6ae99SBarry Smith 696e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 697e3f3e4b6SBarry Smith 698dce8aebaSBarry Smith Fortran Note: 699*20f4b53cSBarry Smith In Fortran one must pass in arrays `lx`, `ly`, and `lz` that are long enough to hold the values; the sixth, seventh and 700dce8aebaSBarry Smith eighth arguments from `DMDAGetInfo()` 701dce8aebaSBarry Smith 702dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()` 70347c6ae99SBarry Smith @*/ 704d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOwnershipRanges(DM da, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[]) 705d71ae5a4SJacob Faibussowitsch { 70647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 70747c6ae99SBarry Smith 70847c6ae99SBarry Smith PetscFunctionBegin; 709a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 71047c6ae99SBarry Smith if (lx) *lx = dd->lx; 71147c6ae99SBarry Smith if (ly) *ly = dd->ly; 71247c6ae99SBarry Smith if (lz) *lz = dd->lz; 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71447c6ae99SBarry Smith } 71547c6ae99SBarry Smith 71647c6ae99SBarry Smith /*@ 717dce8aebaSBarry Smith DMDASetRefinementFactor - Set the ratios that the `DMDA` grid is refined 71847c6ae99SBarry Smith 719*20f4b53cSBarry Smith Logically Collective 72047c6ae99SBarry Smith 72147c6ae99SBarry Smith Input Parameters: 722dce8aebaSBarry Smith + da - the `DMDA` object 72347c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 72447c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 72547c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 72647c6ae99SBarry Smith 727dce8aebaSBarry Smith Options Database Keys: 72848eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction 72948eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction 73048eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction 73148eeb7c8SBarry Smith - -da_refine <n> - refine the DMDA object n times when it is created. 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith Level: intermediate 73447c6ae99SBarry Smith 735dce8aebaSBarry Smith Note: 736dce8aebaSBarry Smith Pass `PETSC_IGNORE` to leave a value unchanged 73747c6ae99SBarry Smith 738dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMRefine()`, `DMDAGetRefinementFactor()` 73947c6ae99SBarry Smith @*/ 740d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z) 741d71ae5a4SJacob Faibussowitsch { 74247c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 74347c6ae99SBarry Smith 74447c6ae99SBarry Smith PetscFunctionBegin; 745a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 74647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_x, 2); 74747c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_y, 3); 74847c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_z, 4); 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 75147c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 75247c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 7533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75447c6ae99SBarry Smith } 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith /*@C 757dce8aebaSBarry Smith DMDAGetRefinementFactor - Gets the ratios that the `DMDA` grid is refined 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith Not Collective 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith Input Parameter: 762dce8aebaSBarry Smith . da - the `DMDA` object 76347c6ae99SBarry Smith 76447c6ae99SBarry Smith Output Parameters: 76547c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 76647c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 76747c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 76847c6ae99SBarry Smith 76947c6ae99SBarry Smith Level: intermediate 77047c6ae99SBarry Smith 771dce8aebaSBarry Smith Note: 772*20f4b53cSBarry Smith Pass `NULL` for values you do not need 77347c6ae99SBarry Smith 774dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMRefine()`, `DMDASetRefinementFactor()` 77547c6ae99SBarry Smith @*/ 776d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y, PetscInt *refine_z) 777d71ae5a4SJacob Faibussowitsch { 77847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 77947c6ae99SBarry Smith 78047c6ae99SBarry Smith PetscFunctionBegin; 781a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 78247c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 78347c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 78447c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 7853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78647c6ae99SBarry Smith } 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith /*@C 789dce8aebaSBarry Smith DMDASetGetMatrix - Sets the routine used by the `DMDA` to allocate a matrix. 79047c6ae99SBarry Smith 791*20f4b53cSBarry Smith Logically Collective; No Fortran Support 79247c6ae99SBarry Smith 79347c6ae99SBarry Smith Input Parameters: 794dce8aebaSBarry Smith + da - the `DMDA` object 795aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 79647c6ae99SBarry Smith 79747c6ae99SBarry Smith Level: developer 79847c6ae99SBarry Smith 799dce8aebaSBarry Smith Note: 800dce8aebaSBarry Smith See `DMDASetBlockFills()` that provides a simple way to provide the nonzero structure for 80147c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 80247c6ae99SBarry Smith 803dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateMatrix()`, `DMDASetBlockFills()` 80447c6ae99SBarry Smith @*/ 805d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM, Mat *)) 806d71ae5a4SJacob Faibussowitsch { 80747c6ae99SBarry Smith PetscFunctionBegin; 808a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 80925296bd5SBarry Smith da->ops->creatematrix = f; 8103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 81338fb4e8eSJunchao Zhang /*@ 814dce8aebaSBarry Smith DMDAMapMatStencilToGlobal - Map a list of `MatStencil` on a grid to global indices. 81538fb4e8eSJunchao Zhang 81638fb4e8eSJunchao Zhang Not Collective 81738fb4e8eSJunchao Zhang 81838fb4e8eSJunchao Zhang Input Parameters: 819dce8aebaSBarry Smith + da - the `DMDA` object 82038fb4e8eSJunchao Zhang . m - number of MatStencils 82138fb4e8eSJunchao Zhang - idxm - grid points (and component number when dof > 1) 82238fb4e8eSJunchao Zhang 8231179163eSBarry Smith Output Parameter: 8241179163eSBarry Smith . gidxm - global row indices 8251179163eSBarry Smith 8261179163eSBarry Smith Level: intermediate 82738fb4e8eSJunchao Zhang 828dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `MatStencil` 82938fb4e8eSJunchao Zhang @*/ 830d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[]) 831d71ae5a4SJacob Faibussowitsch { 83238fb4e8eSJunchao Zhang const DM_DA *dd = (const DM_DA *)da->data; 83338fb4e8eSJunchao Zhang const PetscInt *dxm = (const PetscInt *)idxm; 83438fb4e8eSJunchao Zhang PetscInt i, j, sdim, tmp, dim; 83538fb4e8eSJunchao Zhang PetscInt dims[4], starts[4], dims2[3], starts2[3], dof = dd->w; 83638fb4e8eSJunchao Zhang ISLocalToGlobalMapping ltog; 83738fb4e8eSJunchao Zhang 83838fb4e8eSJunchao Zhang PetscFunctionBegin; 8393ba16761SJacob Faibussowitsch if (m <= 0) PetscFunctionReturn(PETSC_SUCCESS); 84038fb4e8eSJunchao Zhang 8412f27f4d1SJunchao Zhang /* Code adapted from DMDAGetGhostCorners() */ 84238fb4e8eSJunchao Zhang starts2[0] = dd->Xs / dof + dd->xo; 84338fb4e8eSJunchao Zhang starts2[1] = dd->Ys + dd->yo; 84438fb4e8eSJunchao Zhang starts2[2] = dd->Zs + dd->zo; 84538fb4e8eSJunchao Zhang dims2[0] = (dd->Xe - dd->Xs) / dof; 84638fb4e8eSJunchao Zhang dims2[1] = (dd->Ye - dd->Ys); 84738fb4e8eSJunchao Zhang dims2[2] = (dd->Ze - dd->Zs); 84838fb4e8eSJunchao Zhang 8492f27f4d1SJunchao Zhang /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */ 8502f27f4d1SJunchao Zhang dim = da->dim; /* DA dim: 1 to 3 */ 8512f27f4d1SJunchao Zhang sdim = dim + (dof > 1 ? 1 : 0); /* Dimensions in MatStencil's (k,j,i,c) view */ 8522f27f4d1SJunchao Zhang for (i = 0; i < dim; i++) { /* Reverse the order and also skip the unused dimensions */ 8532f27f4d1SJunchao Zhang dims[i] = dims2[dim - i - 1]; /* ex. dims/starts[] are in order of {i} for 1D, {j,i} for 2D and {k,j,i} for 3D */ 85438fb4e8eSJunchao Zhang starts[i] = starts2[dim - i - 1]; 85538fb4e8eSJunchao Zhang } 8562f27f4d1SJunchao Zhang starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */ 85738fb4e8eSJunchao Zhang dims[dim] = dof; 85838fb4e8eSJunchao Zhang 85938fb4e8eSJunchao Zhang /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */ 86038fb4e8eSJunchao Zhang for (i = 0; i < m; i++) { 8612f27f4d1SJunchao Zhang dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */ 8622f27f4d1SJunchao Zhang tmp = 0; 8632f27f4d1SJunchao Zhang for (j = 0; j < sdim; j++) { /* Iter over, ex. j,i or j,i,c in 2D */ 8642f27f4d1SJunchao Zhang if (tmp < 0 || dxm[j] < starts[j] || dxm[j] >= (starts[j] + dims[j])) tmp = -1; /* Beyond the ghost region, therefore ignored with negative indices */ 8652f27f4d1SJunchao Zhang else tmp = tmp * dims[j] + (dxm[j] - starts[j]); 86638fb4e8eSJunchao Zhang } 86738fb4e8eSJunchao Zhang gidxm[i] = tmp; 8682f27f4d1SJunchao Zhang /* Move to the next MatStencil point */ 8692f27f4d1SJunchao Zhang if (dof > 1) dxm += sdim; /* c is already counted in sdim */ 8702f27f4d1SJunchao Zhang else dxm += sdim + 1; /* skip the unused c */ 87138fb4e8eSJunchao Zhang } 87238fb4e8eSJunchao Zhang 87338fb4e8eSJunchao Zhang /* Map local indices to global indices */ 87438fb4e8eSJunchao Zhang PetscCall(DMGetLocalToGlobalMapping(da, <og)); 87538fb4e8eSJunchao Zhang PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm)); 8763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87738fb4e8eSJunchao Zhang } 87838fb4e8eSJunchao Zhang 87947c6ae99SBarry Smith /* 88047c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 88147c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 88247c6ae99SBarry Smith 88347c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 88447c6ae99SBarry Smith */ 885d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[]) 886d71ae5a4SJacob Faibussowitsch { 88747c6ae99SBarry Smith PetscInt i, totalc = 0, remaining, startc = 0, startf = 0; 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith PetscFunctionBegin; 89063a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 89147c6ae99SBarry Smith if (ratio == 1) { 8929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lf, lc, m)); 8933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89447c6ae99SBarry Smith } 89547c6ae99SBarry Smith for (i = 0; i < m; i++) totalc += lc[i]; 89647c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 89747c6ae99SBarry Smith for (i = 0; i < m; i++) { 89847c6ae99SBarry Smith PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 89947c6ae99SBarry Smith if (i == m - 1) lf[i] = want; 90047c6ae99SBarry Smith else { 9017aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 9027aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 9037aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 9047aca7175SJed Brown while ((startf + want) / ratio < nextc - stencil_width) want++; 9057aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 9067aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 9077aca7175SJed Brown while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--; 9087aca7175SJed Brown /* Make sure all constraints are satisfied */ 9099371c9d4SSatish Balay if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width)) 9109371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range"); 91147c6ae99SBarry Smith } 91247c6ae99SBarry Smith lf[i] = want; 91347c6ae99SBarry Smith startc += lc[i]; 91447c6ae99SBarry Smith startf += lf[i]; 91547c6ae99SBarry Smith remaining -= lf[i]; 91647c6ae99SBarry Smith } 9173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91847c6ae99SBarry Smith } 91947c6ae99SBarry Smith 9202be375d4SJed Brown /* 9212be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 9222be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 9232be375d4SJed Brown 9242be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 9252be375d4SJed Brown */ 926d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[]) 927d71ae5a4SJacob Faibussowitsch { 9282be375d4SJed Brown PetscInt i, totalf, remaining, startc, startf; 9292be375d4SJed Brown 9302be375d4SJed Brown PetscFunctionBegin; 93163a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 9322be375d4SJed Brown if (ratio == 1) { 9339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lc, lf, m)); 9343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9352be375d4SJed Brown } 9362be375d4SJed Brown for (i = 0, totalf = 0; i < m; i++) totalf += lf[i]; 9372be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 9382be375d4SJed Brown for (i = 0, startc = 0, startf = 0; i < m; i++) { 9392be375d4SJed Brown PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 9402be375d4SJed Brown if (i == m - 1) lc[i] = want; 9412be375d4SJed Brown else { 9422be375d4SJed Brown const PetscInt nextf = startf + lf[i]; 9432be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 9442be375d4SJed Brown * node is within one stencil width. */ 9452be375d4SJed Brown while (nextf / ratio < startc + want - stencil_width) want--; 9462be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 9472be375d4SJed Brown * fine node is within one stencil width. */ 9482be375d4SJed Brown while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++; 9499371c9d4SSatish Balay if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width)) 9509371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range"); 9512be375d4SJed Brown } 9522be375d4SJed Brown lc[i] = want; 9532be375d4SJed Brown startc += lc[i]; 9542be375d4SJed Brown startf += lf[i]; 9552be375d4SJed Brown remaining -= lc[i]; 9562be375d4SJed Brown } 9573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9582be375d4SJed Brown } 9592be375d4SJed Brown 960d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref) 961d71ae5a4SJacob Faibussowitsch { 962c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 9636858538eSMatthew G. Knepley Vec coordsc, coordsf; 9649a42bb27SBarry Smith DM da2; 96547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data, *dd2; 96647c6ae99SBarry Smith 96747c6ae99SBarry Smith PetscFunctionBegin; 968a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 96947c6ae99SBarry Smith PetscValidPointer(daref, 3); 97047c6ae99SBarry Smith 9719566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 972bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 97347c6ae99SBarry Smith M = dd->refine_x * dd->M; 97447c6ae99SBarry Smith } else { 97547c6ae99SBarry Smith M = 1 + dd->refine_x * (dd->M - 1); 97647c6ae99SBarry Smith } 977bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 978c73cfb54SMatthew G. Knepley if (dim > 1) { 97947c6ae99SBarry Smith N = dd->refine_y * dd->N; 98047c6ae99SBarry Smith } else { 9811860e6e9SBarry Smith N = 1; 9821860e6e9SBarry Smith } 9831860e6e9SBarry Smith } else { 98447c6ae99SBarry Smith N = 1 + dd->refine_y * (dd->N - 1); 98547c6ae99SBarry Smith } 986bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 987c73cfb54SMatthew G. Knepley if (dim > 2) { 98847c6ae99SBarry Smith P = dd->refine_z * dd->P; 98947c6ae99SBarry Smith } else { 9901860e6e9SBarry Smith P = 1; 9911860e6e9SBarry Smith } 9921860e6e9SBarry Smith } else { 99347c6ae99SBarry Smith P = 1 + dd->refine_z * (dd->P - 1); 99447c6ae99SBarry Smith } 9959566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2)); 9969566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix)); 9979566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da2, dim)); 9989566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da2, M, N, P)); 9999566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p)); 10009566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz)); 10019566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da2, dd->w)); 10029566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da2, dd->stencil_type)); 10039566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da2, dd->s)); 1004c73cfb54SMatthew G. Knepley if (dim == 3) { 100547c6ae99SBarry Smith PetscInt *lx, *ly, *lz; 10069566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 10079566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10089566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10099566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz)); 10109566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz)); 10119566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1012c73cfb54SMatthew G. Knepley } else if (dim == 2) { 101347c6ae99SBarry Smith PetscInt *lx, *ly; 10149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 10159566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10169566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10179566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL)); 10189566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1019c73cfb54SMatthew G. Knepley } else if (dim == 1) { 102047c6ae99SBarry Smith PetscInt *lx; 10219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 10229566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10239566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL)); 10249566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith dd2 = (DM_DA *)da2->data; 102747c6ae99SBarry Smith 102847c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 102925296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 103025296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 103147c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 103247c6ae99SBarry Smith dd2->interptype = dd->interptype; 103347c6ae99SBarry Smith 103447c6ae99SBarry Smith /* copy fill information if given */ 103547c6ae99SBarry Smith if (dd->dfill) { 10369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 10379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 103847c6ae99SBarry Smith } 103947c6ae99SBarry Smith if (dd->ofill) { 10409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 10419566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 104247c6ae99SBarry Smith } 104347c6ae99SBarry Smith /* copy the refine information */ 1044397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1045397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1046397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 104747c6ae99SBarry Smith 1048897f7067SBarry Smith if (dd->refine_z_hier) { 1049ad540459SPierre Jolivet if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; 1050ad540459SPierre Jolivet if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; 1051897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 10529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 10539566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1054897f7067SBarry Smith } 1055897f7067SBarry Smith if (dd->refine_y_hier) { 1056ad540459SPierre Jolivet if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; 1057ad540459SPierre Jolivet if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; 1058897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 10599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 10609566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1061897f7067SBarry Smith } 1062897f7067SBarry Smith if (dd->refine_x_hier) { 1063ad540459SPierre Jolivet if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; 1064ad540459SPierre Jolivet if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; 1065897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 10669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 10679566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1068897f7067SBarry Smith } 1069897f7067SBarry Smith 107047c6ae99SBarry Smith /* copy vector type information */ 10719566063dSJacob Faibussowitsch PetscCall(DMSetVecType(da2, da->vectype)); 1072ddcf8b74SDave May 1073efd51863SBarry Smith dd2->lf = dd->lf; 1074efd51863SBarry Smith dd2->lj = dd->lj; 1075efd51863SBarry Smith 10766e87535bSJed Brown da2->leveldown = da->leveldown; 10776e87535bSJed Brown da2->levelup = da->levelup + 1; 10788865f1eaSKarl Rupp 10799566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2)); 10806e87535bSJed Brown 1081ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10826858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(da, &coordsc)); 10836858538eSMatthew G. Knepley if (coordsc) { 1084ddcf8b74SDave May DM cdaf, cdac; 1085ddcf8b74SDave May Mat II; 1086ddcf8b74SDave May 10879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cdac)); 10889566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da2, &cdaf)); 1089b61d3410SDave May /* force creation of the coordinate vector */ 10909566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 10919566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da2, &coordsf)); 10929566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL)); 10939566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 10949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 1095ddcf8b74SDave May } 1096397b6216SJed Brown 1097f3141302SJed Brown for (i = 0; i < da->bs; i++) { 1098f3141302SJed Brown const char *fieldname; 10999566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, i, &fieldname)); 11009566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2, i, fieldname)); 1101f3141302SJed Brown } 1102397b6216SJed Brown 110347c6ae99SBarry Smith *daref = da2; 11043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110547c6ae99SBarry Smith } 110647c6ae99SBarry Smith 1107d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc) 1108d71ae5a4SJacob Faibussowitsch { 1109c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 11106858538eSMatthew G. Knepley Vec coordsc, coordsf; 1111a5bc1bf3SBarry Smith DM dmc2; 1112a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA *)dmf->data, *dd2; 111347c6ae99SBarry Smith 111447c6ae99SBarry Smith PetscFunctionBegin; 1115a5bc1bf3SBarry Smith PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA); 1116a5bc1bf3SBarry Smith PetscValidPointer(dmc, 3); 111747c6ae99SBarry Smith 11189566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmf, &dim)); 1119bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1120397b6216SJed Brown M = dd->M / dd->coarsen_x; 112147c6ae99SBarry Smith } else { 1122397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 112347c6ae99SBarry Smith } 1124bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1125c73cfb54SMatthew G. Knepley if (dim > 1) { 1126397b6216SJed Brown N = dd->N / dd->coarsen_y; 112747c6ae99SBarry Smith } else { 11281860e6e9SBarry Smith N = 1; 11291860e6e9SBarry Smith } 11301860e6e9SBarry Smith } else { 1131397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 113247c6ae99SBarry Smith } 1133bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1134c73cfb54SMatthew G. Knepley if (dim > 2) { 1135397b6216SJed Brown P = dd->P / dd->coarsen_z; 113647c6ae99SBarry Smith } else { 11371860e6e9SBarry Smith P = 1; 11381860e6e9SBarry Smith } 11391860e6e9SBarry Smith } else { 1140397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 114147c6ae99SBarry Smith } 11429566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2)); 11439566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix)); 11449566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dmc2, dim)); 11459566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(dmc2, M, N, P)); 11469566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p)); 11479566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz)); 11489566063dSJacob Faibussowitsch PetscCall(DMDASetDof(dmc2, dd->w)); 11499566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(dmc2, dd->stencil_type)); 11509566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(dmc2, dd->s)); 1151c73cfb54SMatthew G. Knepley if (dim == 3) { 11522be375d4SJed Brown PetscInt *lx, *ly, *lz; 11539566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 11549566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11559566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 11569566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz)); 11579566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz)); 11589566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1159c73cfb54SMatthew G. Knepley } else if (dim == 2) { 11602be375d4SJed Brown PetscInt *lx, *ly; 11619566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 11629566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11639566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 11649566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL)); 11659566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1166c73cfb54SMatthew G. Knepley } else if (dim == 1) { 11672be375d4SJed Brown PetscInt *lx; 11689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 11699566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11709566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL)); 11719566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 117247c6ae99SBarry Smith } 1173a5bc1bf3SBarry Smith dd2 = (DM_DA *)dmc2->data; 117447c6ae99SBarry Smith 11754dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 1176a5bc1bf3SBarry Smith /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 1177a5bc1bf3SBarry Smith dmc2->ops->creatematrix = dmf->ops->creatematrix; 1178a5bc1bf3SBarry Smith dmc2->ops->getcoloring = dmf->ops->getcoloring; 117947c6ae99SBarry Smith dd2->interptype = dd->interptype; 118047c6ae99SBarry Smith 118147c6ae99SBarry Smith /* copy fill information if given */ 118247c6ae99SBarry Smith if (dd->dfill) { 11839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 11849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 118547c6ae99SBarry Smith } 118647c6ae99SBarry Smith if (dd->ofill) { 11879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 11889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 118947c6ae99SBarry Smith } 119047c6ae99SBarry Smith /* copy the refine information */ 1191397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1192397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1193397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 119447c6ae99SBarry Smith 1195897f7067SBarry Smith if (dd->refine_z_hier) { 1196ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_z_hier_n) dd2->refine_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 1]; 1197ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_z_hier_n) dd2->coarsen_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 2]; 1198897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 11999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 12009566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1201897f7067SBarry Smith } 1202897f7067SBarry Smith if (dd->refine_y_hier) { 1203ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_y_hier_n) dd2->refine_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 1]; 1204ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_y_hier_n) dd2->coarsen_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 2]; 1205897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 12069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 12079566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1208897f7067SBarry Smith } 1209897f7067SBarry Smith if (dd->refine_x_hier) { 1210ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_x_hier_n) dd2->refine_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 1]; 1211ad540459SPierre Jolivet if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_x_hier_n) dd2->coarsen_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 2]; 1212897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 12139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 12149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1215897f7067SBarry Smith } 1216897f7067SBarry Smith 121747c6ae99SBarry Smith /* copy vector type information */ 12189566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dmc2, dmf->vectype)); 121947c6ae99SBarry Smith 1220644e2e5bSBarry Smith dd2->lf = dd->lf; 1221644e2e5bSBarry Smith dd2->lj = dd->lj; 1222644e2e5bSBarry Smith 1223a5bc1bf3SBarry Smith dmc2->leveldown = dmf->leveldown + 1; 1224a5bc1bf3SBarry Smith dmc2->levelup = dmf->levelup; 12258865f1eaSKarl Rupp 12269566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmc2)); 12276e87535bSJed Brown 1228ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 12296858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dmf, &coordsf)); 12306858538eSMatthew G. Knepley if (coordsf) { 1231ddcf8b74SDave May DM cdaf, cdac; 12326dbf9973SLawrence Mitchell Mat inject; 12336dbf9973SLawrence Mitchell VecScatter vscat; 1234ddcf8b74SDave May 12359566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmf, &cdaf)); 12369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmc2, &cdac)); 1237b61d3410SDave May /* force creation of the coordinate vector */ 12389566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 12399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmc2, &coordsc)); 1240ddcf8b74SDave May 12419566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(cdac, cdaf, &inject)); 12429566063dSJacob Faibussowitsch PetscCall(MatScatterGetVecScatter(inject, &vscat)); 12439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inject)); 1246ddcf8b74SDave May } 1247f98405f7SJed Brown 1248a5bc1bf3SBarry Smith for (i = 0; i < dmf->bs; i++) { 1249f98405f7SJed Brown const char *fieldname; 12509566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(dmf, i, &fieldname)); 12519566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(dmc2, i, fieldname)); 1252f98405f7SJed Brown } 1253f98405f7SJed Brown 1254a5bc1bf3SBarry Smith *dmc = dmc2; 12553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125647c6ae99SBarry Smith } 125747c6ae99SBarry Smith 1258d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[]) 1259d71ae5a4SJacob Faibussowitsch { 126047c6ae99SBarry Smith PetscInt i, n, *refx, *refy, *refz; 126147c6ae99SBarry Smith 126247c6ae99SBarry Smith PetscFunctionBegin; 126347c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 12647a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 12653ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 126647c6ae99SBarry Smith PetscValidPointer(daf, 3); 126747c6ae99SBarry Smith 1268aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 12699566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz)); 127048a46eb9SPierre Jolivet for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i])); 127147c6ae99SBarry Smith n = nlevels; 12729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL)); 127347c6ae99SBarry Smith n = nlevels; 12749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL)); 127547c6ae99SBarry Smith n = nlevels; 12769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL)); 127747c6ae99SBarry Smith 12789566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0])); 12799566063dSJacob Faibussowitsch PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0])); 128047c6ae99SBarry Smith for (i = 1; i < nlevels; i++) { 12819566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i])); 12829566063dSJacob Faibussowitsch PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i])); 128347c6ae99SBarry Smith } 12849566063dSJacob Faibussowitsch PetscCall(PetscFree3(refx, refy, refz)); 12853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 128647c6ae99SBarry Smith } 128747c6ae99SBarry Smith 1288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[]) 1289d71ae5a4SJacob Faibussowitsch { 129047c6ae99SBarry Smith PetscInt i; 129147c6ae99SBarry Smith 129247c6ae99SBarry Smith PetscFunctionBegin; 129347c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 12947a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 12953ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 129647c6ae99SBarry Smith PetscValidPointer(dac, 3); 12979566063dSJacob Faibussowitsch PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0])); 129848a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i])); 12993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130047c6ae99SBarry Smith } 130162197512SBarry Smith 1302d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes) 1303d71ae5a4SJacob Faibussowitsch { 13048272889dSSatish Balay PetscInt i, j, xs, xn, q; 130562197512SBarry Smith PetscScalar *xx; 130662197512SBarry Smith PetscReal h; 130762197512SBarry Smith Vec x; 130862197512SBarry Smith DM_DA *da = (DM_DA *)dm->data; 130962197512SBarry Smith 131062197512SBarry Smith PetscFunctionBegin; 131162197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 13129566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 131362197512SBarry Smith q = (q - 1) / (n - 1); /* number of spectral elements */ 131462197512SBarry Smith h = 2.0 / q; 13159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL)); 131662197512SBarry Smith xs = xs / (n - 1); 131762197512SBarry Smith xn = xn / (n - 1); 13189566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.)); 13199566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &x)); 13209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, x, &xx)); 132162197512SBarry Smith 132262197512SBarry Smith /* loop over local spectral elements */ 132362197512SBarry Smith for (j = xs; j < xs + xn; j++) { 132462197512SBarry Smith /* 132562197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 132662197512SBarry Smith */ 1327ad540459SPierre Jolivet for (i = (j == xs && xs > 0) ? 1 : 0; i < n; i++) xx[j * (n - 1) + i] = -1.0 + h * j + h * (nodes[i] + 1.0) / 2.; 132862197512SBarry Smith } 13299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, x, &xx)); 133062197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic"); 13313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133262197512SBarry Smith } 133362197512SBarry Smith 13341fd49c25SBarry Smith /*@ 133562197512SBarry Smith 133662197512SBarry 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 133762197512SBarry Smith 1338*20f4b53cSBarry Smith Collective 133962197512SBarry Smith 134062197512SBarry Smith Input Parameters: 1341dce8aebaSBarry Smith + da - the `DMDA` object 13428272889dSSatish Balay - n - the number of GLL nodes 13438272889dSSatish Balay - nodes - the GLL nodes 134462197512SBarry Smith 1345edc382c3SSatish Balay Level: advanced 1346edc382c3SSatish Balay 1347dce8aebaSBarry Smith Note: 1348dce8aebaSBarry Smith The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 1349dce8aebaSBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the `DM` is 1350dce8aebaSBarry Smith periodic or not. 1351dce8aebaSBarry Smith 1352dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()` 135362197512SBarry Smith @*/ 1354d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes) 1355d71ae5a4SJacob Faibussowitsch { 135662197512SBarry Smith PetscFunctionBegin; 135762197512SBarry Smith if (da->dim == 1) { 13589566063dSJacob Faibussowitsch PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes)); 135962197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d"); 13603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136162197512SBarry Smith } 13627c3cd84eSPatrick Sanan 1363d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set) 1364d71ae5a4SJacob Faibussowitsch { 13657c3cd84eSPatrick Sanan DM_DA *dd1 = (DM_DA *)da1->data, *dd2; 13667c3cd84eSPatrick Sanan DM da2; 13677c3cd84eSPatrick Sanan DMType dmtype2; 13687c3cd84eSPatrick Sanan PetscBool isda, compatibleLocal; 13697c3cd84eSPatrick Sanan PetscInt i; 13707c3cd84eSPatrick Sanan 13717c3cd84eSPatrick Sanan PetscFunctionBegin; 13727a8be351SBarry Smith PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()"); 13739566063dSJacob Faibussowitsch PetscCall(DMGetType(dm2, &dmtype2)); 13749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(dmtype2, DMDA, &isda)); 13757c3cd84eSPatrick Sanan if (isda) { 13767c3cd84eSPatrick Sanan da2 = dm2; 13777c3cd84eSPatrick Sanan dd2 = (DM_DA *)da2->data; 13787a8be351SBarry Smith PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()"); 1379110623a0SKarl Rupp compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1380c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1381110623a0SKarl Rupp /* Global size ranks Boundary type */ 1382c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1383c790a739SKarl Rupp if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1384c790a739SKarl Rupp if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 13857c3cd84eSPatrick Sanan if (compatibleLocal) { 13869371c9d4SSatish Balay for (i = 0; i < dd1->m; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ } 13877c3cd84eSPatrick Sanan } 13887c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 1) { 1389ad540459SPierre Jolivet for (i = 0; i < dd1->n; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 13907c3cd84eSPatrick Sanan } 13917c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 2) { 1392ad540459SPierre Jolivet for (i = 0; i < dd1->p; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 13937c3cd84eSPatrick Sanan } 13947c3cd84eSPatrick Sanan *compatible = compatibleLocal; 13957c3cd84eSPatrick Sanan *set = PETSC_TRUE; 13967c3cd84eSPatrick Sanan } else { 13977c3cd84eSPatrick Sanan /* Decline to determine compatibility with other DM types */ 13987c3cd84eSPatrick Sanan *set = PETSC_FALSE; 13997c3cd84eSPatrick Sanan } 14003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14017c3cd84eSPatrick Sanan } 1402