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 620f4b53cSBarry 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 1660225df5SJacob Faibussowitsch Developer Notes: 17628e512eSPatrick Sanan Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points 18e43dc8daSBarry Smith 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 4420f4b53cSBarry 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 8720f4b53cSBarry Smith Not Collective 8847c6ae99SBarry Smith 89d8d19677SJose E. Roman Input Parameters: 90dce8aebaSBarry Smith + da - The `DMDA` 91*a4e35b19SJacob Faibussowitsch . bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 92*a4e35b19SJacob Faibussowitsch . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 93*a4e35b19SJacob Faibussowitsch - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith Level: intermediate 9647c6ae99SBarry Smith 97dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType` 9847c6ae99SBarry Smith @*/ 99d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz) 100d71ae5a4SJacob Faibussowitsch { 10147c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith PetscFunctionBegin; 104a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1055a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bx, 2); 1065a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, by, 3); 1075a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bz, 4); 1087a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 1091321219cSEthan Coon dd->bx = bx; 1101321219cSEthan Coon dd->by = by; 1111321219cSEthan Coon dd->bz = bz; 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith /*@ 116aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 11747c6ae99SBarry Smith 11820f4b53cSBarry Smith Not Collective 11947c6ae99SBarry Smith 12059f3ab6dSMatthew G. Knepley Input Parameters: 121dce8aebaSBarry Smith + da - The `DMDA` 12247c6ae99SBarry Smith - dof - Number of degrees of freedom 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith Level: intermediate 12547c6ae99SBarry Smith 126dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()` 12747c6ae99SBarry Smith @*/ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetDof(DM da, PetscInt dof) 129d71ae5a4SJacob Faibussowitsch { 13047c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith PetscFunctionBegin; 133a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 13454cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da, dof, 2); 1357a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 13647c6ae99SBarry Smith dd->w = dof; 1371411c6eeSJed Brown da->bs = dof; 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith 141fb6725baSMatthew G. Knepley /*@ 142fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 143fb6725baSMatthew G. Knepley 14420f4b53cSBarry Smith Not Collective 145fb6725baSMatthew G. Knepley 146fb6725baSMatthew G. Knepley Input Parameter: 147dce8aebaSBarry Smith . 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 154dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetDof()`, `DMDACreate()`, `DMDestroy()` 155fb6725baSMatthew G. Knepley @*/ 156d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 157d71ae5a4SJacob Faibussowitsch { 158fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 159fb6725baSMatthew G. Knepley 160fb6725baSMatthew G. Knepley PetscFunctionBegin; 161a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1624f572ea9SToby Isaac PetscAssertPointer(dof, 2); 163fb6725baSMatthew G. Knepley *dof = dd->w; 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165fb6725baSMatthew G. Knepley } 166fb6725baSMatthew G. Knepley 1677ddda789SPeter Brune /*@ 1687ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 1697ddda789SPeter Brune 17020f4b53cSBarry Smith Not Collective 1717ddda789SPeter Brune 172f899ff85SJose E. Roman Input Parameter: 173dce8aebaSBarry Smith . da - The `DMDA` 1747ddda789SPeter Brune 1757ddda789SPeter Brune Output Parameters: 1767ddda789SPeter Brune + x - Overlap in the x direction 1777ddda789SPeter Brune . y - Overlap in the y direction 1787ddda789SPeter Brune - z - Overlap in the z direction 1797ddda789SPeter Brune 1807ddda789SPeter Brune Level: intermediate 1817ddda789SPeter Brune 182dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetOverlap()` 1837ddda789SPeter Brune @*/ 184d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOverlap(DM da, PetscInt *x, PetscInt *y, PetscInt *z) 185d71ae5a4SJacob Faibussowitsch { 1867ddda789SPeter Brune DM_DA *dd = (DM_DA *)da->data; 1877ddda789SPeter Brune 1887ddda789SPeter Brune PetscFunctionBegin; 189a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1907ddda789SPeter Brune if (x) *x = dd->xol; 1917ddda789SPeter Brune if (y) *y = dd->yol; 1927ddda789SPeter Brune if (z) *z = dd->zol; 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1947ddda789SPeter Brune } 1957ddda789SPeter Brune 19688661749SPeter Brune /*@ 19788661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 19888661749SPeter Brune 19920f4b53cSBarry Smith Not Collective 20088661749SPeter Brune 2017ddda789SPeter Brune Input Parameters: 202dce8aebaSBarry Smith + da - The `DMDA` 2037ddda789SPeter Brune . x - Overlap in the x direction 2047ddda789SPeter Brune . y - Overlap in the y direction 2057ddda789SPeter Brune - z - Overlap in the z direction 20688661749SPeter Brune 20788661749SPeter Brune Level: intermediate 20888661749SPeter Brune 209dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetOverlap()` 21088661749SPeter Brune @*/ 211d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z) 212d71ae5a4SJacob Faibussowitsch { 21388661749SPeter Brune DM_DA *dd = (DM_DA *)da->data; 21488661749SPeter Brune 21588661749SPeter Brune PetscFunctionBegin; 216a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2177ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, x, 2); 2187ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, y, 3); 2197ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, z, 4); 2207ddda789SPeter Brune dd->xol = x; 2217ddda789SPeter Brune dd->yol = y; 2227ddda789SPeter Brune dd->zol = z; 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22488661749SPeter Brune } 22588661749SPeter Brune 2263e7870d2SPeter Brune /*@ 2273e7870d2SPeter Brune DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 2283e7870d2SPeter Brune 22920f4b53cSBarry Smith Not Collective 2303e7870d2SPeter Brune 2312fe279fdSBarry Smith Input Parameter: 232dce8aebaSBarry Smith . da - The `DMDA` 2333e7870d2SPeter Brune 2342fe279fdSBarry Smith Output Parameter: 235a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition 2363e7870d2SPeter Brune 2373e7870d2SPeter Brune Level: intermediate 2383e7870d2SPeter Brune 239dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()` 2403e7870d2SPeter Brune @*/ 241d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub) 242d71ae5a4SJacob Faibussowitsch { 2433e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2443e7870d2SPeter Brune 2453e7870d2SPeter Brune PetscFunctionBegin; 246a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2473e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2493e7870d2SPeter Brune } 2503e7870d2SPeter Brune 2513e7870d2SPeter Brune /*@ 2523e7870d2SPeter Brune DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 2533e7870d2SPeter Brune 25420f4b53cSBarry Smith Not Collective 2553e7870d2SPeter Brune 2563e7870d2SPeter Brune Input Parameters: 257dce8aebaSBarry Smith + da - The `DMDA` 2583e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2593e7870d2SPeter Brune 2603e7870d2SPeter Brune Level: intermediate 2613e7870d2SPeter Brune 262dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()` 2633e7870d2SPeter Brune @*/ 264d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub) 265d71ae5a4SJacob Faibussowitsch { 2663e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2673e7870d2SPeter Brune 2683e7870d2SPeter Brune PetscFunctionBegin; 269a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2703e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da, Nsub, 2); 2713e7870d2SPeter Brune dd->Nsub = Nsub; 2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2733e7870d2SPeter Brune } 2743e7870d2SPeter Brune 275d886c4f4SPeter Brune /*@ 276d886c4f4SPeter Brune DMDASetOffset - Sets the index offset of the DA. 277d886c4f4SPeter Brune 27820f4b53cSBarry Smith Collective 279d886c4f4SPeter Brune 280d8d19677SJose E. Roman Input Parameters: 281dce8aebaSBarry Smith + da - The `DMDA` 282d886c4f4SPeter Brune . xo - The offset in the x direction 283d886c4f4SPeter Brune . yo - The offset in the y direction 284*a4e35b19SJacob Faibussowitsch . zo - The offset in the z direction 285*a4e35b19SJacob Faibussowitsch . Mo - The problem offset in the x direction 286*a4e35b19SJacob Faibussowitsch . No - The problem offset in the y direction 287*a4e35b19SJacob Faibussowitsch - Po - The problem offset in the z direction 288d886c4f4SPeter Brune 289d886c4f4SPeter Brune Level: intermediate 290d886c4f4SPeter Brune 291dce8aebaSBarry Smith Note: 292dce8aebaSBarry Smith This is used primarily to overlap a computation on a local `DMDA` with that on a global `DMDA` without 293d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 294d886c4f4SPeter Brune 295dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 296d886c4f4SPeter Brune @*/ 297d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 298d71ae5a4SJacob Faibussowitsch { 299d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 300d886c4f4SPeter Brune 301d886c4f4SPeter Brune PetscFunctionBegin; 302a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 303d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da, xo, 2); 30495c13181SPeter Brune PetscValidLogicalCollectiveInt(da, yo, 3); 30595c13181SPeter Brune PetscValidLogicalCollectiveInt(da, zo, 4); 30695c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Mo, 5); 30795c13181SPeter Brune PetscValidLogicalCollectiveInt(da, No, 6); 30895c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Po, 7); 309d886c4f4SPeter Brune dd->xo = xo; 310d886c4f4SPeter Brune dd->yo = yo; 311d886c4f4SPeter Brune dd->zo = zo; 31295c13181SPeter Brune dd->Mo = Mo; 31395c13181SPeter Brune dd->No = No; 31495c13181SPeter Brune dd->Po = Po; 31595c13181SPeter Brune 3166858538eSMatthew G. Knepley if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po)); 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 318d886c4f4SPeter Brune } 319d886c4f4SPeter Brune 320d886c4f4SPeter Brune /*@ 321dce8aebaSBarry Smith DMDAGetOffset - Gets the index offset of the `DMDA`. 322d886c4f4SPeter Brune 32320f4b53cSBarry Smith Not Collective 324d886c4f4SPeter Brune 325d886c4f4SPeter Brune Input Parameter: 326dce8aebaSBarry Smith . da - The `DMDA` 327d886c4f4SPeter Brune 328d886c4f4SPeter Brune Output Parameters: 329d886c4f4SPeter Brune + xo - The offset in the x direction 330d886c4f4SPeter Brune . yo - The offset in the y direction 33195c13181SPeter Brune . zo - The offset in the z direction 33295c13181SPeter Brune . Mo - The global size in the x direction 33395c13181SPeter Brune . No - The global size in the y direction 33495c13181SPeter Brune - Po - The global size in the z direction 335d886c4f4SPeter Brune 336d886c4f4SPeter Brune Level: intermediate 337d886c4f4SPeter Brune 338dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetOffset()`, `DMDAVecGetArray()` 339d886c4f4SPeter Brune @*/ 340d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo, PetscInt *Mo, PetscInt *No, PetscInt *Po) 341d71ae5a4SJacob Faibussowitsch { 342d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 343d886c4f4SPeter Brune 344d886c4f4SPeter Brune PetscFunctionBegin; 345a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 346d886c4f4SPeter Brune if (xo) *xo = dd->xo; 347d886c4f4SPeter Brune if (yo) *yo = dd->yo; 348d886c4f4SPeter Brune if (zo) *zo = dd->zo; 34995c13181SPeter Brune if (Mo) *Mo = dd->Mo; 35095c13181SPeter Brune if (No) *No = dd->No; 35195c13181SPeter Brune if (Po) *Po = dd->Po; 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 353d886c4f4SPeter Brune } 354d886c4f4SPeter Brune 35540234c92SPeter Brune /*@ 356dce8aebaSBarry Smith DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain `DM`. 35740234c92SPeter Brune 35820f4b53cSBarry Smith Not Collective 35940234c92SPeter Brune 36040234c92SPeter Brune Input Parameter: 361dce8aebaSBarry Smith . da - The `DMDA` 36240234c92SPeter Brune 36340234c92SPeter Brune Output Parameters: 36440234c92SPeter Brune + xs - The start of the region in x 36540234c92SPeter Brune . ys - The start of the region in y 36640234c92SPeter Brune . zs - The start of the region in z 367*a4e35b19SJacob Faibussowitsch . xm - The size of the region in x 368*a4e35b19SJacob Faibussowitsch . ym - The size of the region in y 369*a4e35b19SJacob Faibussowitsch - zm - The size of the region in z 37040234c92SPeter Brune 37140234c92SPeter Brune Level: intermediate 37240234c92SPeter Brune 373dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 37440234c92SPeter Brune @*/ 375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 376d71ae5a4SJacob Faibussowitsch { 37740234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 37840234c92SPeter Brune 37940234c92SPeter Brune PetscFunctionBegin; 380a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 38140234c92SPeter Brune if (xs) *xs = dd->nonxs; 38240234c92SPeter Brune if (ys) *ys = dd->nonys; 38340234c92SPeter Brune if (zs) *zs = dd->nonzs; 38440234c92SPeter Brune if (xm) *xm = dd->nonxm; 38540234c92SPeter Brune if (ym) *ym = dd->nonym; 38640234c92SPeter Brune if (zm) *zm = dd->nonzm; 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38840234c92SPeter Brune } 38940234c92SPeter Brune 39040234c92SPeter Brune /*@ 391dce8aebaSBarry Smith DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain `DM`. 39240234c92SPeter Brune 39320f4b53cSBarry Smith Collective 39440234c92SPeter Brune 395d8d19677SJose E. Roman Input Parameters: 396dce8aebaSBarry Smith + da - The `DMDA` 39740234c92SPeter Brune . xs - The start of the region in x 39840234c92SPeter Brune . ys - The start of the region in y 39940234c92SPeter Brune . zs - The start of the region in z 400*a4e35b19SJacob Faibussowitsch . xm - The size of the region in x 401*a4e35b19SJacob Faibussowitsch . ym - The size of the region in y 402*a4e35b19SJacob Faibussowitsch - zm - The size of the region in z 40340234c92SPeter Brune 40440234c92SPeter Brune Level: intermediate 40540234c92SPeter Brune 406dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 40740234c92SPeter Brune @*/ 408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 409d71ae5a4SJacob Faibussowitsch { 41040234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 41140234c92SPeter Brune 41240234c92SPeter Brune PetscFunctionBegin; 413a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 41440234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xs, 2); 41540234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ys, 3); 41640234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zs, 4); 41740234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xm, 5); 41840234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ym, 6); 41940234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zm, 7); 42040234c92SPeter Brune dd->nonxs = xs; 42140234c92SPeter Brune dd->nonys = ys; 42240234c92SPeter Brune dd->nonzs = zs; 42340234c92SPeter Brune dd->nonxm = xm; 42440234c92SPeter Brune dd->nonym = ym; 42540234c92SPeter Brune dd->nonzm = zm; 42640234c92SPeter Brune 4273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42840234c92SPeter Brune } 42988661749SPeter Brune 43047c6ae99SBarry Smith /*@ 431aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 43247c6ae99SBarry Smith 43320f4b53cSBarry Smith Logically Collective 43447c6ae99SBarry Smith 435d8d19677SJose E. Roman Input Parameters: 436dce8aebaSBarry Smith + da - The `DMDA` 437dce8aebaSBarry Smith - stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 43847c6ae99SBarry Smith 43947c6ae99SBarry Smith Level: intermediate 44047c6ae99SBarry Smith 441dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 44247c6ae99SBarry Smith @*/ 443d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 444d71ae5a4SJacob Faibussowitsch { 44547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 44647c6ae99SBarry Smith 44747c6ae99SBarry Smith PetscFunctionBegin; 448a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 44947c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, stype, 2); 4507a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 45147c6ae99SBarry Smith dd->stencil_type = stype; 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45347c6ae99SBarry Smith } 45447c6ae99SBarry Smith 455fb6725baSMatthew G. Knepley /*@ 456fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 457fb6725baSMatthew G. Knepley 45820f4b53cSBarry Smith Not Collective 459fb6725baSMatthew G. Knepley 460fb6725baSMatthew G. Knepley Input Parameter: 461dce8aebaSBarry Smith . da - The `DMDA` 462fb6725baSMatthew G. Knepley 463fb6725baSMatthew G. Knepley Output Parameter: 464dce8aebaSBarry Smith . stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 465fb6725baSMatthew G. Knepley 466fb6725baSMatthew G. Knepley Level: intermediate 467fb6725baSMatthew G. Knepley 468dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 469fb6725baSMatthew G. Knepley @*/ 470d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 471d71ae5a4SJacob Faibussowitsch { 472fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 473fb6725baSMatthew G. Knepley 474fb6725baSMatthew G. Knepley PetscFunctionBegin; 475a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 4764f572ea9SToby Isaac PetscAssertPointer(stype, 2); 477fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 479fb6725baSMatthew G. Knepley } 480fb6725baSMatthew G. Knepley 48147c6ae99SBarry Smith /*@ 482aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 48347c6ae99SBarry Smith 48420f4b53cSBarry Smith Logically Collective 48547c6ae99SBarry Smith 486d8d19677SJose E. Roman Input Parameters: 487dce8aebaSBarry Smith + da - The `DMDA` 48847c6ae99SBarry Smith - width - The stencil width 48947c6ae99SBarry Smith 49047c6ae99SBarry Smith Level: intermediate 49147c6ae99SBarry Smith 492dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 49347c6ae99SBarry Smith @*/ 494d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 495d71ae5a4SJacob Faibussowitsch { 49647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 49747c6ae99SBarry Smith 49847c6ae99SBarry Smith PetscFunctionBegin; 499a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 50047c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, width, 2); 5017a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 50247c6ae99SBarry Smith dd->s = width; 5033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith 506fb6725baSMatthew G. Knepley /*@ 507fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 508fb6725baSMatthew G. Knepley 50920f4b53cSBarry Smith Not Collective 510fb6725baSMatthew G. Knepley 511fb6725baSMatthew G. Knepley Input Parameter: 512dce8aebaSBarry Smith . da - The `DMDA` 513fb6725baSMatthew G. Knepley 514fb6725baSMatthew G. Knepley Output Parameter: 515fb6725baSMatthew G. Knepley . width - The stencil width 516fb6725baSMatthew G. Knepley 517fb6725baSMatthew G. Knepley Level: intermediate 518fb6725baSMatthew G. Knepley 519dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 520fb6725baSMatthew G. Knepley @*/ 521d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 522d71ae5a4SJacob Faibussowitsch { 523fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 524fb6725baSMatthew G. Knepley 525fb6725baSMatthew G. Knepley PetscFunctionBegin; 526a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5274f572ea9SToby Isaac PetscAssertPointer(width, 2); 528fb6725baSMatthew G. Knepley *width = dd->s; 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 530fb6725baSMatthew G. Knepley } 531fb6725baSMatthew G. Knepley 532d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[]) 533d71ae5a4SJacob Faibussowitsch { 53447c6ae99SBarry Smith PetscInt i, sum; 53547c6ae99SBarry Smith 53647c6ae99SBarry Smith PetscFunctionBegin; 5377a8be351SBarry Smith PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set"); 53847c6ae99SBarry Smith for (i = sum = 0; i < m; i++) sum += lx[i]; 53963a3b9bcSJacob 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); 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54147c6ae99SBarry Smith } 54247c6ae99SBarry Smith 54347c6ae99SBarry Smith /*@ 544aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 54547c6ae99SBarry Smith 54620f4b53cSBarry Smith Logically Collective 54747c6ae99SBarry Smith 548d8d19677SJose E. Roman Input Parameters: 549dce8aebaSBarry Smith + da - The `DMDA` 5500298fd71SBarry 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 5510298fd71SBarry 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 5520298fd71SBarry 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. 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith Level: intermediate 55547c6ae99SBarry Smith 556*a4e35b19SJacob Faibussowitsch Note: 557*a4e35b19SJacob Faibussowitsch These numbers are NOT multiplied by the number of dof per node. 558e3f3e4b6SBarry Smith 559dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 56047c6ae99SBarry Smith @*/ 561d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 562d71ae5a4SJacob Faibussowitsch { 56347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 56447c6ae99SBarry Smith 56547c6ae99SBarry Smith PetscFunctionBegin; 566a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5677a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 56847c6ae99SBarry Smith if (lx) { 5697a8be351SBarry Smith PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5709566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx)); 57148a46eb9SPierre Jolivet if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx)); 5729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lx, lx, dd->m)); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith if (ly) { 5757a8be351SBarry Smith PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5769566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly)); 57748a46eb9SPierre Jolivet if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly)); 5789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->ly, ly, dd->n)); 57947c6ae99SBarry Smith } 58047c6ae99SBarry Smith if (lz) { 5817a8be351SBarry Smith PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs"); 5829566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz)); 58348a46eb9SPierre Jolivet if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz)); 5849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lz, lz, dd->p)); 58547c6ae99SBarry Smith } 5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58747c6ae99SBarry Smith } 58847c6ae99SBarry Smith 58947c6ae99SBarry Smith /*@ 590aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 591dce8aebaSBarry Smith returned by `DMCreateInterpolation()` 59247c6ae99SBarry Smith 59320f4b53cSBarry Smith Logically Collective 59447c6ae99SBarry Smith 595d8d19677SJose E. Roman Input Parameters: 59647c6ae99SBarry Smith + da - initial distributed array 597dce8aebaSBarry Smith - ctype - `DMDA_Q1` and `DMDA_Q0` are currently the only supported forms 59847c6ae99SBarry Smith 59947c6ae99SBarry Smith Level: intermediate 60047c6ae99SBarry Smith 601dce8aebaSBarry Smith Note: 602dce8aebaSBarry Smith You should call this on the coarser of the two `DMDA` you pass to `DMCreateInterpolation()` 60347c6ae99SBarry Smith 604dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDAInterpolationType` 60547c6ae99SBarry Smith @*/ 606d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype) 607d71ae5a4SJacob Faibussowitsch { 60847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith PetscFunctionBegin; 611a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 61247c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, ctype, 2); 61347c6ae99SBarry Smith dd->interptype = ctype; 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61547c6ae99SBarry Smith } 61647c6ae99SBarry Smith 6172dde6fd4SLisandro Dalcin /*@ 6182dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 619dce8aebaSBarry Smith used by `DMCreateInterpolation()` 6202dde6fd4SLisandro Dalcin 6212dde6fd4SLisandro Dalcin Not Collective 6222dde6fd4SLisandro Dalcin 6232dde6fd4SLisandro Dalcin Input Parameter: 6242dde6fd4SLisandro Dalcin . da - distributed array 6252dde6fd4SLisandro Dalcin 6262dde6fd4SLisandro Dalcin Output Parameter: 627dce8aebaSBarry Smith . ctype - interpolation type (`DMDA_Q1` and `DMDA_Q0` are currently the only supported forms) 6282dde6fd4SLisandro Dalcin 6292dde6fd4SLisandro Dalcin Level: intermediate 6302dde6fd4SLisandro Dalcin 631dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()` 6322dde6fd4SLisandro Dalcin @*/ 633d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype) 634d71ae5a4SJacob Faibussowitsch { 6352dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 6362dde6fd4SLisandro Dalcin 6372dde6fd4SLisandro Dalcin PetscFunctionBegin; 638a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 6394f572ea9SToby Isaac PetscAssertPointer(ctype, 2); 6402dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6422dde6fd4SLisandro Dalcin } 64347c6ae99SBarry Smith 6446a119db4SBarry Smith /*@C 645aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 64647c6ae99SBarry Smith processes neighbors. 64747c6ae99SBarry Smith 64847c6ae99SBarry Smith Not Collective 64947c6ae99SBarry Smith 65047c6ae99SBarry Smith Input Parameter: 651dce8aebaSBarry Smith . da - the `DMDA` object 65247c6ae99SBarry Smith 6532fe279fdSBarry Smith Output Parameter: 65447c6ae99SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 65547c6ae99SBarry Smith this process itself is in the list 65647c6ae99SBarry Smith 65747c6ae99SBarry Smith Level: intermediate 65847c6ae99SBarry Smith 659dce8aebaSBarry Smith Notes: 660dce8aebaSBarry Smith In 2d the array is of length 9, in 3d of length 27 661dce8aebaSBarry Smith 662dce8aebaSBarry Smith Not supported in 1d 663dce8aebaSBarry Smith 664dce8aebaSBarry Smith Do not free the array, it is freed when the `DMDA` is destroyed. 665dce8aebaSBarry Smith 66660225df5SJacob Faibussowitsch Fortran Notes: 667dce8aebaSBarry Smith In fortran you must pass in an array of the appropriate length. 668dce8aebaSBarry Smith 669dce8aebaSBarry Smith .seealso: `DMDA`, `DM` 67047c6ae99SBarry Smith @*/ 671d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[]) 672d71ae5a4SJacob Faibussowitsch { 67347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 6745fd66863SKarl Rupp 67547c6ae99SBarry Smith PetscFunctionBegin; 676a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 67747c6ae99SBarry Smith *ranks = dd->neighbors; 6783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67947c6ae99SBarry Smith } 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith /*@C 682aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith Not Collective 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith Input Parameter: 687dce8aebaSBarry Smith . da - the `DMDA` object 68847c6ae99SBarry Smith 689d8d19677SJose E. Roman Output Parameters: 69047c6ae99SBarry Smith + lx - ownership along x direction (optional) 69147c6ae99SBarry Smith . ly - ownership along y direction (optional) 69247c6ae99SBarry Smith - lz - ownership along z direction (optional) 69347c6ae99SBarry Smith 69447c6ae99SBarry Smith Level: intermediate 69547c6ae99SBarry Smith 696dce8aebaSBarry Smith Note: 697dce8aebaSBarry Smith These correspond to the optional final arguments passed to `DMDACreate()`, `DMDACreate2d()`, `DMDACreate3d()` 69847c6ae99SBarry Smith 69947c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 700dce8aebaSBarry Smith `DMDA` they came from still exists (has not been destroyed). 70147c6ae99SBarry Smith 702e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 703e3f3e4b6SBarry Smith 70460225df5SJacob Faibussowitsch Fortran Notes: 70520f4b53cSBarry Smith In Fortran one must pass in arrays `lx`, `ly`, and `lz` that are long enough to hold the values; the sixth, seventh and 706dce8aebaSBarry Smith eighth arguments from `DMDAGetInfo()` 707dce8aebaSBarry Smith 708dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()` 70947c6ae99SBarry Smith @*/ 710d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOwnershipRanges(DM da, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[]) 711d71ae5a4SJacob Faibussowitsch { 71247c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith PetscFunctionBegin; 715a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 71647c6ae99SBarry Smith if (lx) *lx = dd->lx; 71747c6ae99SBarry Smith if (ly) *ly = dd->ly; 71847c6ae99SBarry Smith if (lz) *lz = dd->lz; 7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72047c6ae99SBarry Smith } 72147c6ae99SBarry Smith 72247c6ae99SBarry Smith /*@ 723dce8aebaSBarry Smith DMDASetRefinementFactor - Set the ratios that the `DMDA` grid is refined 72447c6ae99SBarry Smith 72520f4b53cSBarry Smith Logically Collective 72647c6ae99SBarry Smith 72747c6ae99SBarry Smith Input Parameters: 728dce8aebaSBarry Smith + da - the `DMDA` object 72947c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 73047c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 73147c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 73247c6ae99SBarry Smith 733dce8aebaSBarry Smith Options Database Keys: 73448eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction 73548eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction 73648eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction 73748eeb7c8SBarry Smith - -da_refine <n> - refine the DMDA object n times when it is created. 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith Level: intermediate 74047c6ae99SBarry Smith 741dce8aebaSBarry Smith Note: 742dce8aebaSBarry Smith Pass `PETSC_IGNORE` to leave a value unchanged 74347c6ae99SBarry Smith 744dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMRefine()`, `DMDAGetRefinementFactor()` 74547c6ae99SBarry Smith @*/ 746d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z) 747d71ae5a4SJacob Faibussowitsch { 74847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 74947c6ae99SBarry Smith 75047c6ae99SBarry Smith PetscFunctionBegin; 751a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 75247c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_x, 2); 75347c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_y, 3); 75447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_z, 4); 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 75747c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 75847c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 7593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76047c6ae99SBarry Smith } 76147c6ae99SBarry Smith 76247c6ae99SBarry Smith /*@C 763dce8aebaSBarry Smith DMDAGetRefinementFactor - Gets the ratios that the `DMDA` grid is refined 76447c6ae99SBarry Smith 76547c6ae99SBarry Smith Not Collective 76647c6ae99SBarry Smith 76747c6ae99SBarry Smith Input Parameter: 768dce8aebaSBarry Smith . da - the `DMDA` object 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith Output Parameters: 77147c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 77247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 77347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 77447c6ae99SBarry Smith 77547c6ae99SBarry Smith Level: intermediate 77647c6ae99SBarry Smith 777dce8aebaSBarry Smith Note: 77820f4b53cSBarry Smith Pass `NULL` for values you do not need 77947c6ae99SBarry Smith 780dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMRefine()`, `DMDASetRefinementFactor()` 78147c6ae99SBarry Smith @*/ 782d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y, PetscInt *refine_z) 783d71ae5a4SJacob Faibussowitsch { 78447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith PetscFunctionBegin; 787a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 78847c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 78947c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 79047c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 7913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 79447c6ae99SBarry Smith /*@C 795dce8aebaSBarry Smith DMDASetGetMatrix - Sets the routine used by the `DMDA` to allocate a matrix. 79647c6ae99SBarry Smith 79720f4b53cSBarry Smith Logically Collective; No Fortran Support 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith Input Parameters: 800dce8aebaSBarry Smith + da - the `DMDA` object 801aa219208SBarry Smith - f - the function that allocates the matrix for that specific DMDA 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith Level: developer 80447c6ae99SBarry Smith 805dce8aebaSBarry Smith Note: 806dce8aebaSBarry Smith See `DMDASetBlockFills()` that provides a simple way to provide the nonzero structure for 80747c6ae99SBarry Smith the diagonal and off-diagonal blocks of the matrix 80847c6ae99SBarry Smith 809dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateMatrix()`, `DMDASetBlockFills()` 81047c6ae99SBarry Smith @*/ 811d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM, Mat *)) 812d71ae5a4SJacob Faibussowitsch { 81347c6ae99SBarry Smith PetscFunctionBegin; 814a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 81525296bd5SBarry Smith da->ops->creatematrix = f; 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81747c6ae99SBarry Smith } 81847c6ae99SBarry Smith 81938fb4e8eSJunchao Zhang /*@ 820dce8aebaSBarry Smith DMDAMapMatStencilToGlobal - Map a list of `MatStencil` on a grid to global indices. 82138fb4e8eSJunchao Zhang 82238fb4e8eSJunchao Zhang Not Collective 82338fb4e8eSJunchao Zhang 82438fb4e8eSJunchao Zhang Input Parameters: 825dce8aebaSBarry Smith + da - the `DMDA` object 82638fb4e8eSJunchao Zhang . m - number of MatStencils 82738fb4e8eSJunchao Zhang - idxm - grid points (and component number when dof > 1) 82838fb4e8eSJunchao Zhang 8291179163eSBarry Smith Output Parameter: 8301179163eSBarry Smith . gidxm - global row indices 8311179163eSBarry Smith 8321179163eSBarry Smith Level: intermediate 83338fb4e8eSJunchao Zhang 834dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `MatStencil` 83538fb4e8eSJunchao Zhang @*/ 836d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[]) 837d71ae5a4SJacob Faibussowitsch { 83838fb4e8eSJunchao Zhang const DM_DA *dd = (const DM_DA *)da->data; 83938fb4e8eSJunchao Zhang const PetscInt *dxm = (const PetscInt *)idxm; 84038fb4e8eSJunchao Zhang PetscInt i, j, sdim, tmp, dim; 84138fb4e8eSJunchao Zhang PetscInt dims[4], starts[4], dims2[3], starts2[3], dof = dd->w; 84238fb4e8eSJunchao Zhang ISLocalToGlobalMapping ltog; 84338fb4e8eSJunchao Zhang 84438fb4e8eSJunchao Zhang PetscFunctionBegin; 8453ba16761SJacob Faibussowitsch if (m <= 0) PetscFunctionReturn(PETSC_SUCCESS); 84638fb4e8eSJunchao Zhang 8472f27f4d1SJunchao Zhang /* Code adapted from DMDAGetGhostCorners() */ 84838fb4e8eSJunchao Zhang starts2[0] = dd->Xs / dof + dd->xo; 84938fb4e8eSJunchao Zhang starts2[1] = dd->Ys + dd->yo; 85038fb4e8eSJunchao Zhang starts2[2] = dd->Zs + dd->zo; 85138fb4e8eSJunchao Zhang dims2[0] = (dd->Xe - dd->Xs) / dof; 85238fb4e8eSJunchao Zhang dims2[1] = (dd->Ye - dd->Ys); 85338fb4e8eSJunchao Zhang dims2[2] = (dd->Ze - dd->Zs); 85438fb4e8eSJunchao Zhang 8552f27f4d1SJunchao Zhang /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */ 8562f27f4d1SJunchao Zhang dim = da->dim; /* DA dim: 1 to 3 */ 8572f27f4d1SJunchao Zhang sdim = dim + (dof > 1 ? 1 : 0); /* Dimensions in MatStencil's (k,j,i,c) view */ 8582f27f4d1SJunchao Zhang for (i = 0; i < dim; i++) { /* Reverse the order and also skip the unused dimensions */ 8592f27f4d1SJunchao 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 */ 86038fb4e8eSJunchao Zhang starts[i] = starts2[dim - i - 1]; 86138fb4e8eSJunchao Zhang } 8622f27f4d1SJunchao Zhang starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */ 86338fb4e8eSJunchao Zhang dims[dim] = dof; 86438fb4e8eSJunchao Zhang 86538fb4e8eSJunchao Zhang /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */ 86638fb4e8eSJunchao Zhang for (i = 0; i < m; i++) { 8672f27f4d1SJunchao Zhang dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */ 8682f27f4d1SJunchao Zhang tmp = 0; 8692f27f4d1SJunchao Zhang for (j = 0; j < sdim; j++) { /* Iter over, ex. j,i or j,i,c in 2D */ 8702f27f4d1SJunchao 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 */ 8712f27f4d1SJunchao Zhang else tmp = tmp * dims[j] + (dxm[j] - starts[j]); 87238fb4e8eSJunchao Zhang } 87338fb4e8eSJunchao Zhang gidxm[i] = tmp; 8742f27f4d1SJunchao Zhang /* Move to the next MatStencil point */ 8752f27f4d1SJunchao Zhang if (dof > 1) dxm += sdim; /* c is already counted in sdim */ 8762f27f4d1SJunchao Zhang else dxm += sdim + 1; /* skip the unused c */ 87738fb4e8eSJunchao Zhang } 87838fb4e8eSJunchao Zhang 87938fb4e8eSJunchao Zhang /* Map local indices to global indices */ 88038fb4e8eSJunchao Zhang PetscCall(DMGetLocalToGlobalMapping(da, <og)); 88138fb4e8eSJunchao Zhang PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm)); 8823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 88338fb4e8eSJunchao Zhang } 88438fb4e8eSJunchao Zhang 88547c6ae99SBarry Smith /* 88647c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 88747c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 89047c6ae99SBarry Smith */ 891d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[]) 892d71ae5a4SJacob Faibussowitsch { 89347c6ae99SBarry Smith PetscInt i, totalc = 0, remaining, startc = 0, startf = 0; 89447c6ae99SBarry Smith 89547c6ae99SBarry Smith PetscFunctionBegin; 89663a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 89747c6ae99SBarry Smith if (ratio == 1) { 8989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lf, lc, m)); 8993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90047c6ae99SBarry Smith } 90147c6ae99SBarry Smith for (i = 0; i < m; i++) totalc += lc[i]; 90247c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 90347c6ae99SBarry Smith for (i = 0; i < m; i++) { 90447c6ae99SBarry Smith PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 90547c6ae99SBarry Smith if (i == m - 1) lf[i] = want; 90647c6ae99SBarry Smith else { 9077aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 9087aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 9097aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 9107aca7175SJed Brown while ((startf + want) / ratio < nextc - stencil_width) want++; 9117aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 9127aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 9137aca7175SJed Brown while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--; 9147aca7175SJed Brown /* Make sure all constraints are satisfied */ 9159371c9d4SSatish Balay if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width)) 9169371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range"); 91747c6ae99SBarry Smith } 91847c6ae99SBarry Smith lf[i] = want; 91947c6ae99SBarry Smith startc += lc[i]; 92047c6ae99SBarry Smith startf += lf[i]; 92147c6ae99SBarry Smith remaining -= lf[i]; 92247c6ae99SBarry Smith } 9233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92447c6ae99SBarry Smith } 92547c6ae99SBarry Smith 9262be375d4SJed Brown /* 9272be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 9282be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 9292be375d4SJed Brown 9302be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 9312be375d4SJed Brown */ 932d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[]) 933d71ae5a4SJacob Faibussowitsch { 9342be375d4SJed Brown PetscInt i, totalf, remaining, startc, startf; 9352be375d4SJed Brown 9362be375d4SJed Brown PetscFunctionBegin; 93763a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 9382be375d4SJed Brown if (ratio == 1) { 9399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lc, lf, m)); 9403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9412be375d4SJed Brown } 9422be375d4SJed Brown for (i = 0, totalf = 0; i < m; i++) totalf += lf[i]; 9432be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 9442be375d4SJed Brown for (i = 0, startc = 0, startf = 0; i < m; i++) { 9452be375d4SJed Brown PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 9462be375d4SJed Brown if (i == m - 1) lc[i] = want; 9472be375d4SJed Brown else { 9482be375d4SJed Brown const PetscInt nextf = startf + lf[i]; 9492be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 9502be375d4SJed Brown * node is within one stencil width. */ 9512be375d4SJed Brown while (nextf / ratio < startc + want - stencil_width) want--; 9522be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 9532be375d4SJed Brown * fine node is within one stencil width. */ 9542be375d4SJed Brown while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++; 9559371c9d4SSatish Balay if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width)) 9569371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range"); 9572be375d4SJed Brown } 9582be375d4SJed Brown lc[i] = want; 9592be375d4SJed Brown startc += lc[i]; 9602be375d4SJed Brown startf += lf[i]; 9612be375d4SJed Brown remaining -= lc[i]; 9622be375d4SJed Brown } 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9642be375d4SJed Brown } 9652be375d4SJed Brown 966d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref) 967d71ae5a4SJacob Faibussowitsch { 968c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 9696858538eSMatthew G. Knepley Vec coordsc, coordsf; 9709a42bb27SBarry Smith DM da2; 97147c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data, *dd2; 97247c6ae99SBarry Smith 97347c6ae99SBarry Smith PetscFunctionBegin; 974a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 9754f572ea9SToby Isaac PetscAssertPointer(daref, 3); 97647c6ae99SBarry Smith 9779566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 978bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 97947c6ae99SBarry Smith M = dd->refine_x * dd->M; 98047c6ae99SBarry Smith } else { 98147c6ae99SBarry Smith M = 1 + dd->refine_x * (dd->M - 1); 98247c6ae99SBarry Smith } 983bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 984c73cfb54SMatthew G. Knepley if (dim > 1) { 98547c6ae99SBarry Smith N = dd->refine_y * dd->N; 98647c6ae99SBarry Smith } else { 9871860e6e9SBarry Smith N = 1; 9881860e6e9SBarry Smith } 9891860e6e9SBarry Smith } else { 99047c6ae99SBarry Smith N = 1 + dd->refine_y * (dd->N - 1); 99147c6ae99SBarry Smith } 992bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 993c73cfb54SMatthew G. Knepley if (dim > 2) { 99447c6ae99SBarry Smith P = dd->refine_z * dd->P; 99547c6ae99SBarry Smith } else { 9961860e6e9SBarry Smith P = 1; 9971860e6e9SBarry Smith } 9981860e6e9SBarry Smith } else { 99947c6ae99SBarry Smith P = 1 + dd->refine_z * (dd->P - 1); 100047c6ae99SBarry Smith } 10019566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2)); 10029566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix)); 10039566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da2, dim)); 10049566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da2, M, N, P)); 10059566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p)); 10069566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz)); 10079566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da2, dd->w)); 10089566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da2, dd->stencil_type)); 10099566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da2, dd->s)); 1010c73cfb54SMatthew G. Knepley if (dim == 3) { 101147c6ae99SBarry Smith PetscInt *lx, *ly, *lz; 10129566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 10139566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10149566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10159566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz)); 10169566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz)); 10179566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1018c73cfb54SMatthew G. Knepley } else if (dim == 2) { 101947c6ae99SBarry Smith PetscInt *lx, *ly; 10209566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 10219566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10229566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10239566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL)); 10249566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1025c73cfb54SMatthew G. Knepley } else if (dim == 1) { 102647c6ae99SBarry Smith PetscInt *lx; 10279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 10289566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10299566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL)); 10309566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 103147c6ae99SBarry Smith } 103247c6ae99SBarry Smith dd2 = (DM_DA *)da2->data; 103347c6ae99SBarry Smith 103447c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 103525296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 103625296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 103747c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 103847c6ae99SBarry Smith dd2->interptype = dd->interptype; 103947c6ae99SBarry Smith 104047c6ae99SBarry Smith /* copy fill information if given */ 104147c6ae99SBarry Smith if (dd->dfill) { 10429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 10439566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 104447c6ae99SBarry Smith } 104547c6ae99SBarry Smith if (dd->ofill) { 10469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 10479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 104847c6ae99SBarry Smith } 104947c6ae99SBarry Smith /* copy the refine information */ 1050397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1051397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1052397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 105347c6ae99SBarry Smith 1054897f7067SBarry Smith if (dd->refine_z_hier) { 1055ad540459SPierre 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]; 1056ad540459SPierre 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]; 1057897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 10589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 10599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1060897f7067SBarry Smith } 1061897f7067SBarry Smith if (dd->refine_y_hier) { 1062ad540459SPierre 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]; 1063ad540459SPierre 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]; 1064897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 10659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 10669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1067897f7067SBarry Smith } 1068897f7067SBarry Smith if (dd->refine_x_hier) { 1069ad540459SPierre 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]; 1070ad540459SPierre 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]; 1071897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 10729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 10739566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1074897f7067SBarry Smith } 1075897f7067SBarry Smith 107647c6ae99SBarry Smith /* copy vector type information */ 10779566063dSJacob Faibussowitsch PetscCall(DMSetVecType(da2, da->vectype)); 1078ddcf8b74SDave May 1079efd51863SBarry Smith dd2->lf = dd->lf; 1080efd51863SBarry Smith dd2->lj = dd->lj; 1081efd51863SBarry Smith 10826e87535bSJed Brown da2->leveldown = da->leveldown; 10836e87535bSJed Brown da2->levelup = da->levelup + 1; 10848865f1eaSKarl Rupp 10859566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2)); 10866e87535bSJed Brown 1087ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 10886858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(da, &coordsc)); 10896858538eSMatthew G. Knepley if (coordsc) { 1090ddcf8b74SDave May DM cdaf, cdac; 1091ddcf8b74SDave May Mat II; 1092ddcf8b74SDave May 10939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cdac)); 10949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da2, &cdaf)); 1095b61d3410SDave May /* force creation of the coordinate vector */ 10969566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 10979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da2, &coordsf)); 10989566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL)); 10999566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 11009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 1101ddcf8b74SDave May } 1102397b6216SJed Brown 1103f3141302SJed Brown for (i = 0; i < da->bs; i++) { 1104f3141302SJed Brown const char *fieldname; 11059566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, i, &fieldname)); 11069566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2, i, fieldname)); 1107f3141302SJed Brown } 1108397b6216SJed Brown 110947c6ae99SBarry Smith *daref = da2; 11103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 111147c6ae99SBarry Smith } 111247c6ae99SBarry Smith 1113d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc) 1114d71ae5a4SJacob Faibussowitsch { 1115c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 11166858538eSMatthew G. Knepley Vec coordsc, coordsf; 1117a5bc1bf3SBarry Smith DM dmc2; 1118a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA *)dmf->data, *dd2; 111947c6ae99SBarry Smith 112047c6ae99SBarry Smith PetscFunctionBegin; 1121a5bc1bf3SBarry Smith PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA); 11224f572ea9SToby Isaac PetscAssertPointer(dmc, 3); 112347c6ae99SBarry Smith 11249566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmf, &dim)); 1125bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1126397b6216SJed Brown M = dd->M / dd->coarsen_x; 112747c6ae99SBarry Smith } else { 1128397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 112947c6ae99SBarry Smith } 1130bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1131c73cfb54SMatthew G. Knepley if (dim > 1) { 1132397b6216SJed Brown N = dd->N / dd->coarsen_y; 113347c6ae99SBarry Smith } else { 11341860e6e9SBarry Smith N = 1; 11351860e6e9SBarry Smith } 11361860e6e9SBarry Smith } else { 1137397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 113847c6ae99SBarry Smith } 1139bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1140c73cfb54SMatthew G. Knepley if (dim > 2) { 1141397b6216SJed Brown P = dd->P / dd->coarsen_z; 114247c6ae99SBarry Smith } else { 11431860e6e9SBarry Smith P = 1; 11441860e6e9SBarry Smith } 11451860e6e9SBarry Smith } else { 1146397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 114747c6ae99SBarry Smith } 11489566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2)); 11499566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix)); 11509566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dmc2, dim)); 11519566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(dmc2, M, N, P)); 11529566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p)); 11539566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz)); 11549566063dSJacob Faibussowitsch PetscCall(DMDASetDof(dmc2, dd->w)); 11559566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(dmc2, dd->stencil_type)); 11569566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(dmc2, dd->s)); 1157c73cfb54SMatthew G. Knepley if (dim == 3) { 11582be375d4SJed Brown PetscInt *lx, *ly, *lz; 11599566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 11609566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11619566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 11629566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz)); 11639566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz)); 11649566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1165c73cfb54SMatthew G. Knepley } else if (dim == 2) { 11662be375d4SJed Brown PetscInt *lx, *ly; 11679566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 11689566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11699566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 11709566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL)); 11719566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1172c73cfb54SMatthew G. Knepley } else if (dim == 1) { 11732be375d4SJed Brown PetscInt *lx; 11749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 11759566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 11769566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL)); 11779566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 117847c6ae99SBarry Smith } 1179a5bc1bf3SBarry Smith dd2 = (DM_DA *)dmc2->data; 118047c6ae99SBarry Smith 11814dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 1182a5bc1bf3SBarry Smith /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 1183a5bc1bf3SBarry Smith dmc2->ops->creatematrix = dmf->ops->creatematrix; 1184a5bc1bf3SBarry Smith dmc2->ops->getcoloring = dmf->ops->getcoloring; 118547c6ae99SBarry Smith dd2->interptype = dd->interptype; 118647c6ae99SBarry Smith 118747c6ae99SBarry Smith /* copy fill information if given */ 118847c6ae99SBarry Smith if (dd->dfill) { 11899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 11909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 119147c6ae99SBarry Smith } 119247c6ae99SBarry Smith if (dd->ofill) { 11939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 11949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 119547c6ae99SBarry Smith } 119647c6ae99SBarry Smith /* copy the refine information */ 1197397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1198397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1199397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 120047c6ae99SBarry Smith 1201897f7067SBarry Smith if (dd->refine_z_hier) { 1202ad540459SPierre 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]; 1203ad540459SPierre 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]; 1204897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 12059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 12069566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1207897f7067SBarry Smith } 1208897f7067SBarry Smith if (dd->refine_y_hier) { 1209ad540459SPierre 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]; 1210ad540459SPierre 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]; 1211897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 12129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 12139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1214897f7067SBarry Smith } 1215897f7067SBarry Smith if (dd->refine_x_hier) { 1216ad540459SPierre 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]; 1217ad540459SPierre 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]; 1218897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 12199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 12209566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1221897f7067SBarry Smith } 1222897f7067SBarry Smith 122347c6ae99SBarry Smith /* copy vector type information */ 12249566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dmc2, dmf->vectype)); 122547c6ae99SBarry Smith 1226644e2e5bSBarry Smith dd2->lf = dd->lf; 1227644e2e5bSBarry Smith dd2->lj = dd->lj; 1228644e2e5bSBarry Smith 1229a5bc1bf3SBarry Smith dmc2->leveldown = dmf->leveldown + 1; 1230a5bc1bf3SBarry Smith dmc2->levelup = dmf->levelup; 12318865f1eaSKarl Rupp 12329566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmc2)); 12336e87535bSJed Brown 1234ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 12356858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dmf, &coordsf)); 12366858538eSMatthew G. Knepley if (coordsf) { 1237ddcf8b74SDave May DM cdaf, cdac; 12386dbf9973SLawrence Mitchell Mat inject; 12396dbf9973SLawrence Mitchell VecScatter vscat; 1240ddcf8b74SDave May 12419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmf, &cdaf)); 12429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmc2, &cdac)); 1243b61d3410SDave May /* force creation of the coordinate vector */ 12449566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 12459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmc2, &coordsc)); 1246ddcf8b74SDave May 12479566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(cdac, cdaf, &inject)); 12489566063dSJacob Faibussowitsch PetscCall(MatScatterGetVecScatter(inject, &vscat)); 12499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inject)); 1252ddcf8b74SDave May } 1253f98405f7SJed Brown 1254a5bc1bf3SBarry Smith for (i = 0; i < dmf->bs; i++) { 1255f98405f7SJed Brown const char *fieldname; 12569566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(dmf, i, &fieldname)); 12579566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(dmc2, i, fieldname)); 1258f98405f7SJed Brown } 1259f98405f7SJed Brown 1260a5bc1bf3SBarry Smith *dmc = dmc2; 12613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126247c6ae99SBarry Smith } 126347c6ae99SBarry Smith 1264d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[]) 1265d71ae5a4SJacob Faibussowitsch { 126647c6ae99SBarry Smith PetscInt i, n, *refx, *refy, *refz; 126747c6ae99SBarry Smith 126847c6ae99SBarry Smith PetscFunctionBegin; 126947c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 12707a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 12713ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 12724f572ea9SToby Isaac PetscAssertPointer(daf, 3); 127347c6ae99SBarry Smith 1274aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 12759566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz)); 127648a46eb9SPierre Jolivet for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i])); 127747c6ae99SBarry Smith n = nlevels; 12789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL)); 127947c6ae99SBarry Smith n = nlevels; 12809566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL)); 128147c6ae99SBarry Smith n = nlevels; 12829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL)); 128347c6ae99SBarry Smith 12849566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0])); 12859566063dSJacob Faibussowitsch PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0])); 128647c6ae99SBarry Smith for (i = 1; i < nlevels; i++) { 12879566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i])); 12889566063dSJacob Faibussowitsch PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i])); 128947c6ae99SBarry Smith } 12909566063dSJacob Faibussowitsch PetscCall(PetscFree3(refx, refy, refz)); 12913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129247c6ae99SBarry Smith } 129347c6ae99SBarry Smith 1294d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[]) 1295d71ae5a4SJacob Faibussowitsch { 129647c6ae99SBarry Smith PetscInt i; 129747c6ae99SBarry Smith 129847c6ae99SBarry Smith PetscFunctionBegin; 129947c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 13007a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 13013ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 13024f572ea9SToby Isaac PetscAssertPointer(dac, 3); 13039566063dSJacob Faibussowitsch PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0])); 130448a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i])); 13053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130647c6ae99SBarry Smith } 130762197512SBarry Smith 1308d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes) 1309d71ae5a4SJacob Faibussowitsch { 13108272889dSSatish Balay PetscInt i, j, xs, xn, q; 131162197512SBarry Smith PetscScalar *xx; 131262197512SBarry Smith PetscReal h; 131362197512SBarry Smith Vec x; 131462197512SBarry Smith DM_DA *da = (DM_DA *)dm->data; 131562197512SBarry Smith 131662197512SBarry Smith PetscFunctionBegin; 131762197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 13189566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 131962197512SBarry Smith q = (q - 1) / (n - 1); /* number of spectral elements */ 132062197512SBarry Smith h = 2.0 / q; 13219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL)); 132262197512SBarry Smith xs = xs / (n - 1); 132362197512SBarry Smith xn = xn / (n - 1); 13249566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.)); 13259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &x)); 13269566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, x, &xx)); 132762197512SBarry Smith 132862197512SBarry Smith /* loop over local spectral elements */ 132962197512SBarry Smith for (j = xs; j < xs + xn; j++) { 133062197512SBarry Smith /* 133162197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 133262197512SBarry Smith */ 1333ad540459SPierre 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.; 133462197512SBarry Smith } 13359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, x, &xx)); 133662197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic"); 13373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133862197512SBarry Smith } 133962197512SBarry Smith 13401fd49c25SBarry Smith /*@ 134162197512SBarry Smith 134262197512SBarry 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 134362197512SBarry Smith 134420f4b53cSBarry Smith Collective 134562197512SBarry Smith 134662197512SBarry Smith Input Parameters: 1347dce8aebaSBarry Smith + da - the `DMDA` object 13482fe279fdSBarry Smith . n - the number of GLL nodes 13498272889dSSatish Balay - nodes - the GLL nodes 135062197512SBarry Smith 1351edc382c3SSatish Balay Level: advanced 1352edc382c3SSatish Balay 1353dce8aebaSBarry Smith Note: 1354dce8aebaSBarry Smith The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 1355dce8aebaSBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the `DM` is 1356dce8aebaSBarry Smith periodic or not. 1357dce8aebaSBarry Smith 1358dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()` 135962197512SBarry Smith @*/ 1360d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes) 1361d71ae5a4SJacob Faibussowitsch { 136262197512SBarry Smith PetscFunctionBegin; 136362197512SBarry Smith if (da->dim == 1) { 13649566063dSJacob Faibussowitsch PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes)); 136562197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d"); 13663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136762197512SBarry Smith } 13687c3cd84eSPatrick Sanan 1369d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set) 1370d71ae5a4SJacob Faibussowitsch { 13717c3cd84eSPatrick Sanan DM_DA *dd1 = (DM_DA *)da1->data, *dd2; 13727c3cd84eSPatrick Sanan DM da2; 13737c3cd84eSPatrick Sanan DMType dmtype2; 13747c3cd84eSPatrick Sanan PetscBool isda, compatibleLocal; 13757c3cd84eSPatrick Sanan PetscInt i; 13767c3cd84eSPatrick Sanan 13777c3cd84eSPatrick Sanan PetscFunctionBegin; 13787a8be351SBarry Smith PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()"); 13799566063dSJacob Faibussowitsch PetscCall(DMGetType(dm2, &dmtype2)); 13809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(dmtype2, DMDA, &isda)); 13817c3cd84eSPatrick Sanan if (isda) { 13827c3cd84eSPatrick Sanan da2 = dm2; 13837c3cd84eSPatrick Sanan dd2 = (DM_DA *)da2->data; 13847a8be351SBarry Smith PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()"); 1385110623a0SKarl Rupp compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1386c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1387110623a0SKarl Rupp /* Global size ranks Boundary type */ 1388c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1389c790a739SKarl Rupp if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1390c790a739SKarl Rupp if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 13917c3cd84eSPatrick Sanan if (compatibleLocal) { 13929371c9d4SSatish Balay for (i = 0; i < dd1->m; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ } 13937c3cd84eSPatrick Sanan } 13947c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 1) { 1395ad540459SPierre Jolivet for (i = 0; i < dd1->n; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 13967c3cd84eSPatrick Sanan } 13977c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 2) { 1398ad540459SPierre Jolivet for (i = 0; i < dd1->p; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 13997c3cd84eSPatrick Sanan } 14007c3cd84eSPatrick Sanan *compatible = compatibleLocal; 14017c3cd84eSPatrick Sanan *set = PETSC_TRUE; 14027c3cd84eSPatrick Sanan } else { 14037c3cd84eSPatrick Sanan /* Decline to determine compatibility with other DM types */ 14047c3cd84eSPatrick Sanan *set = PETSC_FALSE; 14057c3cd84eSPatrick Sanan } 14063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14077c3cd84eSPatrick Sanan } 1408