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 1612b4a537SBarry 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 1912b4a537SBarry Smith .seealso: [](sec_struct), `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` 4812b4a537SBarry Smith . m - the number of X processes (or `PETSC_DECIDE`) 4912b4a537SBarry Smith . n - the number of Y processes (or `PETSC_DECIDE`) 5012b4a537SBarry Smith - p - the number of Z processes (or `PETSC_DECIDE`) 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith Level: intermediate 5347c6ae99SBarry Smith 5412b4a537SBarry Smith .seealso: [](sec_struct), `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 /*@ 85fe58071bSMatthew G. Knepley DMDAGetBoundaryType - Gets the type of ghost nodes on domain boundaries. 86fe58071bSMatthew G. Knepley 87fe58071bSMatthew G. Knepley Not Collective 88fe58071bSMatthew G. Knepley 89fe58071bSMatthew G. Knepley Input Parameter: 90fe58071bSMatthew G. Knepley . da - The `DMDA` 91fe58071bSMatthew G. Knepley 92fe58071bSMatthew G. Knepley Output Parameters: 93fe58071bSMatthew G. Knepley + bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 94fe58071bSMatthew G. Knepley . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 95fe58071bSMatthew G. Knepley - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 96fe58071bSMatthew G. Knepley 97fe58071bSMatthew G. Knepley Level: intermediate 98fe58071bSMatthew G. Knepley 99fe58071bSMatthew G. Knepley .seealso: [](sec_struct), `DMDASetBoundaryType()`, `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType`, `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 100fe58071bSMatthew G. Knepley @*/ 101fe58071bSMatthew G. Knepley PetscErrorCode DMDAGetBoundaryType(DM da, DMBoundaryType *bx, DMBoundaryType *by, DMBoundaryType *bz) 102fe58071bSMatthew G. Knepley { 103fe58071bSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 104fe58071bSMatthew G. Knepley 105fe58071bSMatthew G. Knepley PetscFunctionBegin; 106fe58071bSMatthew G. Knepley PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 107fe58071bSMatthew G. Knepley if (bx) PetscAssertPointer(bx, 2); 108fe58071bSMatthew G. Knepley if (by) PetscAssertPointer(by, 3); 109fe58071bSMatthew G. Knepley if (bz) PetscAssertPointer(bz, 4); 110fe58071bSMatthew G. Knepley if (bx) *bx = dd->bx; 111fe58071bSMatthew G. Knepley if (by) *by = dd->by; 112fe58071bSMatthew G. Knepley if (bz) *bz = dd->bz; 113fe58071bSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 114fe58071bSMatthew G. Knepley } 115fe58071bSMatthew G. Knepley 116fe58071bSMatthew G. Knepley /*@ 1171321219cSEthan Coon DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 11847c6ae99SBarry Smith 11920f4b53cSBarry Smith Not Collective 12047c6ae99SBarry Smith 121d8d19677SJose E. Roman Input Parameters: 122dce8aebaSBarry Smith + da - The `DMDA` 123a4e35b19SJacob Faibussowitsch . bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 124a4e35b19SJacob Faibussowitsch . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 125a4e35b19SJacob Faibussowitsch - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith Level: intermediate 12847c6ae99SBarry Smith 129fe58071bSMatthew G. Knepley .seealso: [](sec_struct), `DMDAGetBoundaryType()`, `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType`, `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC` 13047c6ae99SBarry Smith @*/ 131d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz) 132d71ae5a4SJacob Faibussowitsch { 13347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 13447c6ae99SBarry Smith 13547c6ae99SBarry Smith PetscFunctionBegin; 136a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1375a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bx, 2); 1385a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, by, 3); 1395a43f728SLisandro Dalcin PetscValidLogicalCollectiveEnum(da, bz, 4); 1407a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 1411321219cSEthan Coon dd->bx = bx; 1421321219cSEthan Coon dd->by = by; 1431321219cSEthan Coon dd->bz = bz; 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14547c6ae99SBarry Smith } 14647c6ae99SBarry Smith 14747c6ae99SBarry Smith /*@ 148aa219208SBarry Smith DMDASetDof - Sets the number of degrees of freedom per vertex 14947c6ae99SBarry Smith 15020f4b53cSBarry Smith Not Collective 15147c6ae99SBarry Smith 15259f3ab6dSMatthew G. Knepley Input Parameters: 153dce8aebaSBarry Smith + da - The `DMDA` 15412b4a537SBarry Smith - dof - Number of degrees of freedom per vertex 15547c6ae99SBarry Smith 15647c6ae99SBarry Smith Level: intermediate 15747c6ae99SBarry Smith 15812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()` 15947c6ae99SBarry Smith @*/ 160d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetDof(DM da, PetscInt dof) 161d71ae5a4SJacob Faibussowitsch { 16247c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 16347c6ae99SBarry Smith 16447c6ae99SBarry Smith PetscFunctionBegin; 165a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 16654cfb0beSLisandro Dalcin PetscValidLogicalCollectiveInt(da, dof, 2); 1677a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 16847c6ae99SBarry Smith dd->w = dof; 1691411c6eeSJed Brown da->bs = dof; 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17147c6ae99SBarry Smith } 17247c6ae99SBarry Smith 173fb6725baSMatthew G. Knepley /*@ 174fb6725baSMatthew G. Knepley DMDAGetDof - Gets the number of degrees of freedom per vertex 175fb6725baSMatthew G. Knepley 17620f4b53cSBarry Smith Not Collective 177fb6725baSMatthew G. Knepley 178fb6725baSMatthew G. Knepley Input Parameter: 179dce8aebaSBarry Smith . da - The `DMDA` 180fb6725baSMatthew G. Knepley 181fb6725baSMatthew G. Knepley Output Parameter: 18212b4a537SBarry Smith . dof - Number of degrees of freedom per vertex 183fb6725baSMatthew G. Knepley 184fb6725baSMatthew G. Knepley Level: intermediate 185fb6725baSMatthew G. Knepley 18612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetDof()`, `DMDACreate()`, `DMDestroy()` 187fb6725baSMatthew G. Knepley @*/ 188d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 189d71ae5a4SJacob Faibussowitsch { 190fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 191fb6725baSMatthew G. Knepley 192fb6725baSMatthew G. Knepley PetscFunctionBegin; 193a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1944f572ea9SToby Isaac PetscAssertPointer(dof, 2); 195fb6725baSMatthew G. Knepley *dof = dd->w; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197fb6725baSMatthew G. Knepley } 198fb6725baSMatthew G. Knepley 1997ddda789SPeter Brune /*@ 2007ddda789SPeter Brune DMDAGetOverlap - Gets the size of the per-processor overlap. 2017ddda789SPeter Brune 20220f4b53cSBarry Smith Not Collective 2037ddda789SPeter Brune 204f899ff85SJose E. Roman Input Parameter: 205dce8aebaSBarry Smith . da - The `DMDA` 2067ddda789SPeter Brune 2077ddda789SPeter Brune Output Parameters: 2087ddda789SPeter Brune + x - Overlap in the x direction 2097ddda789SPeter Brune . y - Overlap in the y direction 2107ddda789SPeter Brune - z - Overlap in the z direction 2117ddda789SPeter Brune 2127ddda789SPeter Brune Level: intermediate 2137ddda789SPeter Brune 21412b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetOverlap()` 2157ddda789SPeter Brune @*/ 216d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOverlap(DM da, PetscInt *x, PetscInt *y, PetscInt *z) 217d71ae5a4SJacob Faibussowitsch { 2187ddda789SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2197ddda789SPeter Brune 2207ddda789SPeter Brune PetscFunctionBegin; 221a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2227ddda789SPeter Brune if (x) *x = dd->xol; 2237ddda789SPeter Brune if (y) *y = dd->yol; 2247ddda789SPeter Brune if (z) *z = dd->zol; 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2267ddda789SPeter Brune } 2277ddda789SPeter Brune 22888661749SPeter Brune /*@ 22988661749SPeter Brune DMDASetOverlap - Sets the size of the per-processor overlap. 23088661749SPeter Brune 23120f4b53cSBarry Smith Not Collective 23288661749SPeter Brune 2337ddda789SPeter Brune Input Parameters: 234dce8aebaSBarry Smith + da - The `DMDA` 2357ddda789SPeter Brune . x - Overlap in the x direction 2367ddda789SPeter Brune . y - Overlap in the y direction 2377ddda789SPeter Brune - z - Overlap in the z direction 23888661749SPeter Brune 23988661749SPeter Brune Level: intermediate 24088661749SPeter Brune 24112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetOverlap()` 24288661749SPeter Brune @*/ 243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z) 244d71ae5a4SJacob Faibussowitsch { 24588661749SPeter Brune DM_DA *dd = (DM_DA *)da->data; 24688661749SPeter Brune 24788661749SPeter Brune PetscFunctionBegin; 248a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2497ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, x, 2); 2507ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, y, 3); 2517ddda789SPeter Brune PetscValidLogicalCollectiveInt(da, z, 4); 2527ddda789SPeter Brune dd->xol = x; 2537ddda789SPeter Brune dd->yol = y; 2547ddda789SPeter Brune dd->zol = z; 2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25688661749SPeter Brune } 25788661749SPeter Brune 2583e7870d2SPeter Brune /*@ 25912b4a537SBarry Smith DMDAGetNumLocalSubDomains - Gets the number of local subdomains that would be created upon decomposition. 2603e7870d2SPeter Brune 26120f4b53cSBarry Smith Not Collective 2623e7870d2SPeter Brune 2632fe279fdSBarry Smith Input Parameter: 264dce8aebaSBarry Smith . da - The `DMDA` 2653e7870d2SPeter Brune 2662fe279fdSBarry Smith Output Parameter: 267a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition 2683e7870d2SPeter Brune 2693e7870d2SPeter Brune Level: intermediate 2703e7870d2SPeter Brune 27112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()` 2723e7870d2SPeter Brune @*/ 273d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub) 274d71ae5a4SJacob Faibussowitsch { 2753e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2763e7870d2SPeter Brune 2773e7870d2SPeter Brune PetscFunctionBegin; 278a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2793e7870d2SPeter Brune if (Nsub) *Nsub = dd->Nsub; 2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2813e7870d2SPeter Brune } 2823e7870d2SPeter Brune 2833e7870d2SPeter Brune /*@ 28412b4a537SBarry Smith DMDASetNumLocalSubDomains - Sets the number of local subdomains to create when decomposing with `DMCreateDomainDecomposition()` 2853e7870d2SPeter Brune 28620f4b53cSBarry Smith Not Collective 2873e7870d2SPeter Brune 2883e7870d2SPeter Brune Input Parameters: 289dce8aebaSBarry Smith + da - The `DMDA` 2903e7870d2SPeter Brune - Nsub - The number of local subdomains requested 2913e7870d2SPeter Brune 2923e7870d2SPeter Brune Level: intermediate 2933e7870d2SPeter Brune 29412b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()` 2953e7870d2SPeter Brune @*/ 296d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub) 297d71ae5a4SJacob Faibussowitsch { 2983e7870d2SPeter Brune DM_DA *dd = (DM_DA *)da->data; 2993e7870d2SPeter Brune 3003e7870d2SPeter Brune PetscFunctionBegin; 301a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 3023e7870d2SPeter Brune PetscValidLogicalCollectiveInt(da, Nsub, 2); 3033e7870d2SPeter Brune dd->Nsub = Nsub; 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3053e7870d2SPeter Brune } 3063e7870d2SPeter Brune 307d886c4f4SPeter Brune /*@ 30812b4a537SBarry Smith DMDASetOffset - Sets the index offset of the `DMDA`. 309d886c4f4SPeter Brune 31020f4b53cSBarry Smith Collective 311d886c4f4SPeter Brune 312d8d19677SJose E. Roman Input Parameters: 313dce8aebaSBarry Smith + da - The `DMDA` 314d886c4f4SPeter Brune . xo - The offset in the x direction 315d886c4f4SPeter Brune . yo - The offset in the y direction 316a4e35b19SJacob Faibussowitsch . zo - The offset in the z direction 317a4e35b19SJacob Faibussowitsch . Mo - The problem offset in the x direction 318a4e35b19SJacob Faibussowitsch . No - The problem offset in the y direction 319a4e35b19SJacob Faibussowitsch - Po - The problem offset in the z direction 320d886c4f4SPeter Brune 32112b4a537SBarry Smith Level: developer 322d886c4f4SPeter Brune 323dce8aebaSBarry Smith Note: 324dce8aebaSBarry Smith This is used primarily to overlap a computation on a local `DMDA` with that on a global `DMDA` without 325d886c4f4SPeter Brune changing boundary conditions or subdomain features that depend upon the global offsets. 326d886c4f4SPeter Brune 32712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 328d886c4f4SPeter Brune @*/ 329d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 330d71ae5a4SJacob Faibussowitsch { 331d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 332d886c4f4SPeter Brune 333d886c4f4SPeter Brune PetscFunctionBegin; 334a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 335d886c4f4SPeter Brune PetscValidLogicalCollectiveInt(da, xo, 2); 33695c13181SPeter Brune PetscValidLogicalCollectiveInt(da, yo, 3); 33795c13181SPeter Brune PetscValidLogicalCollectiveInt(da, zo, 4); 33895c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Mo, 5); 33995c13181SPeter Brune PetscValidLogicalCollectiveInt(da, No, 6); 34095c13181SPeter Brune PetscValidLogicalCollectiveInt(da, Po, 7); 341d886c4f4SPeter Brune dd->xo = xo; 342d886c4f4SPeter Brune dd->yo = yo; 343d886c4f4SPeter Brune dd->zo = zo; 34495c13181SPeter Brune dd->Mo = Mo; 34595c13181SPeter Brune dd->No = No; 34695c13181SPeter Brune dd->Po = Po; 34795c13181SPeter Brune 3486858538eSMatthew G. Knepley if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po)); 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 350d886c4f4SPeter Brune } 351d886c4f4SPeter Brune 352d886c4f4SPeter Brune /*@ 353dce8aebaSBarry Smith DMDAGetOffset - Gets the index offset of the `DMDA`. 354d886c4f4SPeter Brune 35520f4b53cSBarry Smith Not Collective 356d886c4f4SPeter Brune 357d886c4f4SPeter Brune Input Parameter: 358dce8aebaSBarry Smith . da - The `DMDA` 359d886c4f4SPeter Brune 360d886c4f4SPeter Brune Output Parameters: 361d886c4f4SPeter Brune + xo - The offset in the x direction 362d886c4f4SPeter Brune . yo - The offset in the y direction 36395c13181SPeter Brune . zo - The offset in the z direction 36495c13181SPeter Brune . Mo - The global size in the x direction 36595c13181SPeter Brune . No - The global size in the y direction 36695c13181SPeter Brune - Po - The global size in the z direction 367d886c4f4SPeter Brune 36812b4a537SBarry Smith Level: developer 369d886c4f4SPeter Brune 37012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetOffset()`, `DMDAVecGetArray()` 371d886c4f4SPeter Brune @*/ 372d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo, PetscInt *Mo, PetscInt *No, PetscInt *Po) 373d71ae5a4SJacob Faibussowitsch { 374d886c4f4SPeter Brune DM_DA *dd = (DM_DA *)da->data; 375d886c4f4SPeter Brune 376d886c4f4SPeter Brune PetscFunctionBegin; 377a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 378d886c4f4SPeter Brune if (xo) *xo = dd->xo; 379d886c4f4SPeter Brune if (yo) *yo = dd->yo; 380d886c4f4SPeter Brune if (zo) *zo = dd->zo; 38195c13181SPeter Brune if (Mo) *Mo = dd->Mo; 38295c13181SPeter Brune if (No) *No = dd->No; 38395c13181SPeter Brune if (Po) *Po = dd->Po; 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 385d886c4f4SPeter Brune } 386d886c4f4SPeter Brune 38740234c92SPeter Brune /*@ 38812b4a537SBarry Smith DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain `DMDA`. 38940234c92SPeter Brune 39020f4b53cSBarry Smith Not Collective 39140234c92SPeter Brune 39240234c92SPeter Brune Input Parameter: 393dce8aebaSBarry Smith . da - The `DMDA` 39440234c92SPeter Brune 39540234c92SPeter Brune Output Parameters: 39640234c92SPeter Brune + xs - The start of the region in x 39740234c92SPeter Brune . ys - The start of the region in y 39840234c92SPeter Brune . zs - The start of the region in z 399a4e35b19SJacob Faibussowitsch . xm - The size of the region in x 400a4e35b19SJacob Faibussowitsch . ym - The size of the region in y 401a4e35b19SJacob Faibussowitsch - zm - The size of the region in z 40240234c92SPeter Brune 40340234c92SPeter Brune Level: intermediate 40440234c92SPeter Brune 40512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 40640234c92SPeter Brune @*/ 407d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 408d71ae5a4SJacob Faibussowitsch { 40940234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 41040234c92SPeter Brune 41140234c92SPeter Brune PetscFunctionBegin; 412a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 41340234c92SPeter Brune if (xs) *xs = dd->nonxs; 41440234c92SPeter Brune if (ys) *ys = dd->nonys; 41540234c92SPeter Brune if (zs) *zs = dd->nonzs; 41640234c92SPeter Brune if (xm) *xm = dd->nonxm; 41740234c92SPeter Brune if (ym) *ym = dd->nonym; 41840234c92SPeter Brune if (zm) *zm = dd->nonzm; 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42040234c92SPeter Brune } 42140234c92SPeter Brune 42240234c92SPeter Brune /*@ 42312b4a537SBarry Smith DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain `DMDA`. 42440234c92SPeter Brune 42520f4b53cSBarry Smith Collective 42640234c92SPeter Brune 427d8d19677SJose E. Roman Input Parameters: 428dce8aebaSBarry Smith + da - The `DMDA` 42940234c92SPeter Brune . xs - The start of the region in x 43040234c92SPeter Brune . ys - The start of the region in y 43140234c92SPeter Brune . zs - The start of the region in z 432a4e35b19SJacob Faibussowitsch . xm - The size of the region in x 433a4e35b19SJacob Faibussowitsch . ym - The size of the region in y 434a4e35b19SJacob Faibussowitsch - zm - The size of the region in z 43540234c92SPeter Brune 43640234c92SPeter Brune Level: intermediate 43740234c92SPeter Brune 43812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()` 43940234c92SPeter Brune @*/ 440d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 441d71ae5a4SJacob Faibussowitsch { 44240234c92SPeter Brune DM_DA *dd = (DM_DA *)da->data; 44340234c92SPeter Brune 44440234c92SPeter Brune PetscFunctionBegin; 445a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 44640234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xs, 2); 44740234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ys, 3); 44840234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zs, 4); 44940234c92SPeter Brune PetscValidLogicalCollectiveInt(da, xm, 5); 45040234c92SPeter Brune PetscValidLogicalCollectiveInt(da, ym, 6); 45140234c92SPeter Brune PetscValidLogicalCollectiveInt(da, zm, 7); 45240234c92SPeter Brune dd->nonxs = xs; 45340234c92SPeter Brune dd->nonys = ys; 45440234c92SPeter Brune dd->nonzs = zs; 45540234c92SPeter Brune dd->nonxm = xm; 45640234c92SPeter Brune dd->nonym = ym; 45740234c92SPeter Brune dd->nonzm = zm; 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45940234c92SPeter Brune } 46088661749SPeter Brune 46147c6ae99SBarry Smith /*@ 462aa219208SBarry Smith DMDASetStencilType - Sets the type of the communication stencil 46347c6ae99SBarry Smith 46420f4b53cSBarry Smith Logically Collective 46547c6ae99SBarry Smith 466d8d19677SJose E. Roman Input Parameters: 467dce8aebaSBarry Smith + da - The `DMDA` 468dce8aebaSBarry Smith - stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 46947c6ae99SBarry Smith 47047c6ae99SBarry Smith Level: intermediate 47147c6ae99SBarry Smith 47212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.` 47347c6ae99SBarry Smith @*/ 474d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 475d71ae5a4SJacob Faibussowitsch { 47647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 47747c6ae99SBarry Smith 47847c6ae99SBarry Smith PetscFunctionBegin; 479a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 48047c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, stype, 2); 4817a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 48247c6ae99SBarry Smith dd->stencil_type = stype; 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48447c6ae99SBarry Smith } 48547c6ae99SBarry Smith 486fb6725baSMatthew G. Knepley /*@ 487fb6725baSMatthew G. Knepley DMDAGetStencilType - Gets the type of the communication stencil 488fb6725baSMatthew G. Knepley 48920f4b53cSBarry Smith Not Collective 490fb6725baSMatthew G. Knepley 491fb6725baSMatthew G. Knepley Input Parameter: 492dce8aebaSBarry Smith . da - The `DMDA` 493fb6725baSMatthew G. Knepley 494fb6725baSMatthew G. Knepley Output Parameter: 495dce8aebaSBarry Smith . stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 496fb6725baSMatthew G. Knepley 497fb6725baSMatthew G. Knepley Level: intermediate 498fb6725baSMatthew G. Knepley 49912b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.` 500fb6725baSMatthew G. Knepley @*/ 501d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 502d71ae5a4SJacob Faibussowitsch { 503fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 504fb6725baSMatthew G. Knepley 505fb6725baSMatthew G. Knepley PetscFunctionBegin; 506a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5074f572ea9SToby Isaac PetscAssertPointer(stype, 2); 508fb6725baSMatthew G. Knepley *stype = dd->stencil_type; 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 510fb6725baSMatthew G. Knepley } 511fb6725baSMatthew G. Knepley 51247c6ae99SBarry Smith /*@ 513aa219208SBarry Smith DMDASetStencilWidth - Sets the width of the communication stencil 51447c6ae99SBarry Smith 51520f4b53cSBarry Smith Logically Collective 51647c6ae99SBarry Smith 517d8d19677SJose E. Roman Input Parameters: 518dce8aebaSBarry Smith + da - The `DMDA` 51947c6ae99SBarry Smith - width - The stencil width 52047c6ae99SBarry Smith 52147c6ae99SBarry Smith Level: intermediate 52247c6ae99SBarry Smith 52312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.` 52447c6ae99SBarry Smith @*/ 525d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 526d71ae5a4SJacob Faibussowitsch { 52747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith PetscFunctionBegin; 530a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 53147c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, width, 2); 5327a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 53347c6ae99SBarry Smith dd->s = width; 5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53547c6ae99SBarry Smith } 53647c6ae99SBarry Smith 537fb6725baSMatthew G. Knepley /*@ 538fb6725baSMatthew G. Knepley DMDAGetStencilWidth - Gets the width of the communication stencil 539fb6725baSMatthew G. Knepley 54020f4b53cSBarry Smith Not Collective 541fb6725baSMatthew G. Knepley 542fb6725baSMatthew G. Knepley Input Parameter: 543dce8aebaSBarry Smith . da - The `DMDA` 544fb6725baSMatthew G. Knepley 545fb6725baSMatthew G. Knepley Output Parameter: 546fb6725baSMatthew G. Knepley . width - The stencil width 547fb6725baSMatthew G. Knepley 548fb6725baSMatthew G. Knepley Level: intermediate 549fb6725baSMatthew G. Knepley 55012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.` 551fb6725baSMatthew G. Knepley @*/ 552d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 553d71ae5a4SJacob Faibussowitsch { 554fb6725baSMatthew G. Knepley DM_DA *dd = (DM_DA *)da->data; 555fb6725baSMatthew G. Knepley 556fb6725baSMatthew G. Knepley PetscFunctionBegin; 557a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5584f572ea9SToby Isaac PetscAssertPointer(width, 2); 559fb6725baSMatthew G. Knepley *width = dd->s; 5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 561fb6725baSMatthew G. Knepley } 562fb6725baSMatthew G. Knepley 563d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[]) 564d71ae5a4SJacob Faibussowitsch { 56547c6ae99SBarry Smith PetscInt i, sum; 56647c6ae99SBarry Smith 56747c6ae99SBarry Smith PetscFunctionBegin; 5687a8be351SBarry Smith PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set"); 56947c6ae99SBarry Smith for (i = sum = 0; i < m; i++) sum += lx[i]; 57063a3b9bcSJacob 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); 5713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57247c6ae99SBarry Smith } 57347c6ae99SBarry Smith 57447c6ae99SBarry Smith /*@ 575aa219208SBarry Smith DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 57647c6ae99SBarry Smith 57720f4b53cSBarry Smith Logically Collective 57847c6ae99SBarry Smith 579d8d19677SJose E. Roman Input Parameters: 580dce8aebaSBarry Smith + da - The `DMDA` 58112b4a537SBarry 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 58212b4a537SBarry 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 58312b4a537SBarry 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. 58447c6ae99SBarry Smith 58547c6ae99SBarry Smith Level: intermediate 58647c6ae99SBarry Smith 587a4e35b19SJacob Faibussowitsch Note: 588a4e35b19SJacob Faibussowitsch These numbers are NOT multiplied by the number of dof per node. 589e3f3e4b6SBarry Smith 59012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()` 59147c6ae99SBarry Smith @*/ 592d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 593d71ae5a4SJacob Faibussowitsch { 59447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 59547c6ae99SBarry Smith 59647c6ae99SBarry Smith PetscFunctionBegin; 597a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5987a8be351SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()"); 59947c6ae99SBarry Smith if (lx) { 60012b4a537SBarry Smith PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes"); 6019566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx)); 60248a46eb9SPierre Jolivet if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx)); 6039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lx, lx, dd->m)); 60447c6ae99SBarry Smith } 60547c6ae99SBarry Smith if (ly) { 60612b4a537SBarry Smith PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes"); 6079566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly)); 60848a46eb9SPierre Jolivet if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly)); 6099566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->ly, ly, dd->n)); 61047c6ae99SBarry Smith } 61147c6ae99SBarry Smith if (lz) { 61212b4a537SBarry Smith PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes"); 6139566063dSJacob Faibussowitsch PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz)); 61448a46eb9SPierre Jolivet if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz)); 6159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->lz, lz, dd->p)); 61647c6ae99SBarry Smith } 6173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61847c6ae99SBarry Smith } 61947c6ae99SBarry Smith 62047c6ae99SBarry Smith /*@ 621aa219208SBarry Smith DMDASetInterpolationType - Sets the type of interpolation that will be 622dce8aebaSBarry Smith returned by `DMCreateInterpolation()` 62347c6ae99SBarry Smith 62420f4b53cSBarry Smith Logically Collective 62547c6ae99SBarry Smith 626d8d19677SJose E. Roman Input Parameters: 62747c6ae99SBarry Smith + da - initial distributed array 628dce8aebaSBarry Smith - ctype - `DMDA_Q1` and `DMDA_Q0` are currently the only supported forms 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith Level: intermediate 63147c6ae99SBarry Smith 632dce8aebaSBarry Smith Note: 633dce8aebaSBarry Smith You should call this on the coarser of the two `DMDA` you pass to `DMCreateInterpolation()` 63447c6ae99SBarry Smith 63512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDAInterpolationType`, 63612b4a537SBarry Smith `DMDA_Q1`, `DMDA_Q0` 63747c6ae99SBarry Smith @*/ 638d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype) 639d71ae5a4SJacob Faibussowitsch { 64047c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith PetscFunctionBegin; 643a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 64447c6ae99SBarry Smith PetscValidLogicalCollectiveEnum(da, ctype, 2); 64547c6ae99SBarry Smith dd->interptype = ctype; 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64747c6ae99SBarry Smith } 64847c6ae99SBarry Smith 6492dde6fd4SLisandro Dalcin /*@ 6502dde6fd4SLisandro Dalcin DMDAGetInterpolationType - Gets the type of interpolation that will be 651dce8aebaSBarry Smith used by `DMCreateInterpolation()` 6522dde6fd4SLisandro Dalcin 6532dde6fd4SLisandro Dalcin Not Collective 6542dde6fd4SLisandro Dalcin 6552dde6fd4SLisandro Dalcin Input Parameter: 6562dde6fd4SLisandro Dalcin . da - distributed array 6572dde6fd4SLisandro Dalcin 6582dde6fd4SLisandro Dalcin Output Parameter: 659dce8aebaSBarry Smith . ctype - interpolation type (`DMDA_Q1` and `DMDA_Q0` are currently the only supported forms) 6602dde6fd4SLisandro Dalcin 6612dde6fd4SLisandro Dalcin Level: intermediate 6622dde6fd4SLisandro Dalcin 66312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()`, 66412b4a537SBarry Smith `DMDA_Q1`, `DMDA_Q0` 6652dde6fd4SLisandro Dalcin @*/ 666d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype) 667d71ae5a4SJacob Faibussowitsch { 6682dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA *)da->data; 6692dde6fd4SLisandro Dalcin 6702dde6fd4SLisandro Dalcin PetscFunctionBegin; 671a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 6724f572ea9SToby Isaac PetscAssertPointer(ctype, 2); 6732dde6fd4SLisandro Dalcin *ctype = dd->interptype; 6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6752dde6fd4SLisandro Dalcin } 67647c6ae99SBarry Smith 6776a119db4SBarry Smith /*@C 678aa219208SBarry Smith DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 67947c6ae99SBarry Smith processes neighbors. 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith Not Collective 68247c6ae99SBarry Smith 68347c6ae99SBarry Smith Input Parameter: 684dce8aebaSBarry Smith . da - the `DMDA` object 68547c6ae99SBarry Smith 6862fe279fdSBarry Smith Output Parameter: 68712b4a537SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. The process itself is in the list 68847c6ae99SBarry Smith 68947c6ae99SBarry Smith Level: intermediate 69047c6ae99SBarry Smith 691dce8aebaSBarry Smith Notes: 692*f13dfd9eSBarry Smith In 2d the `ranks` is of length 9, in 3d of length 27 693dce8aebaSBarry Smith 694dce8aebaSBarry Smith Not supported in 1d 695dce8aebaSBarry Smith 696dce8aebaSBarry Smith Do not free the array, it is freed when the `DMDA` is destroyed. 697dce8aebaSBarry Smith 69812b4a537SBarry Smith Fortran Note: 69912b4a537SBarry Smith Pass in an array of the appropriate length to contain the values 700dce8aebaSBarry Smith 70112b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM` 70247c6ae99SBarry Smith @*/ 703d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[]) 704d71ae5a4SJacob Faibussowitsch { 70547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 7065fd66863SKarl Rupp 70747c6ae99SBarry Smith PetscFunctionBegin; 708a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 70947c6ae99SBarry Smith *ranks = dd->neighbors; 7103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71147c6ae99SBarry Smith } 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith /*@C 714aa219208SBarry Smith DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 71547c6ae99SBarry Smith 71647c6ae99SBarry Smith Not Collective 71747c6ae99SBarry Smith 71847c6ae99SBarry Smith Input Parameter: 719dce8aebaSBarry Smith . da - the `DMDA` object 72047c6ae99SBarry Smith 721d8d19677SJose E. Roman Output Parameters: 72247c6ae99SBarry Smith + lx - ownership along x direction (optional) 72347c6ae99SBarry Smith . ly - ownership along y direction (optional) 72447c6ae99SBarry Smith - lz - ownership along z direction (optional) 72547c6ae99SBarry Smith 72647c6ae99SBarry Smith Level: intermediate 72747c6ae99SBarry Smith 728dce8aebaSBarry Smith Note: 729dce8aebaSBarry Smith These correspond to the optional final arguments passed to `DMDACreate()`, `DMDACreate2d()`, `DMDACreate3d()` 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith In C you should not free these arrays, nor change the values in them. They will only have valid values while the 732dce8aebaSBarry Smith `DMDA` they came from still exists (has not been destroyed). 73347c6ae99SBarry Smith 734e3f3e4b6SBarry Smith These numbers are NOT multiplied by the number of dof per node. 735e3f3e4b6SBarry Smith 73612b4a537SBarry Smith Fortran Note: 73712b4a537SBarry Smith Pass in arrays `lx`, `ly`, and `lz` of the appropriate length to hold the values; the sixth, seventh and 738dce8aebaSBarry Smith eighth arguments from `DMDAGetInfo()` 739dce8aebaSBarry Smith 74012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()` 74147c6ae99SBarry Smith @*/ 742d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOwnershipRanges(DM da, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[]) 743d71ae5a4SJacob Faibussowitsch { 74447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 74547c6ae99SBarry Smith 74647c6ae99SBarry Smith PetscFunctionBegin; 747a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 74847c6ae99SBarry Smith if (lx) *lx = dd->lx; 74947c6ae99SBarry Smith if (ly) *ly = dd->ly; 75047c6ae99SBarry Smith if (lz) *lz = dd->lz; 7513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75247c6ae99SBarry Smith } 75347c6ae99SBarry Smith 75447c6ae99SBarry Smith /*@ 755dce8aebaSBarry Smith DMDASetRefinementFactor - Set the ratios that the `DMDA` grid is refined 75647c6ae99SBarry Smith 75720f4b53cSBarry Smith Logically Collective 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith Input Parameters: 760dce8aebaSBarry Smith + da - the `DMDA` object 76147c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default) 76247c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 76347c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 76447c6ae99SBarry Smith 765dce8aebaSBarry Smith Options Database Keys: 76648eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction 76748eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction 76848eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction 76912b4a537SBarry Smith - -da_refine <n> - refine the `DMDA` object n times when it is created. 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith Level: intermediate 77247c6ae99SBarry Smith 773dce8aebaSBarry Smith Note: 774dce8aebaSBarry Smith Pass `PETSC_IGNORE` to leave a value unchanged 77547c6ae99SBarry Smith 77612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDAGetRefinementFactor()` 77747c6ae99SBarry Smith @*/ 778d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z) 779d71ae5a4SJacob Faibussowitsch { 78047c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 78147c6ae99SBarry Smith 78247c6ae99SBarry Smith PetscFunctionBegin; 783a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 78447c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_x, 2); 78547c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_y, 3); 78647c6ae99SBarry Smith PetscValidLogicalCollectiveInt(da, refine_z, 4); 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith if (refine_x > 0) dd->refine_x = refine_x; 78947c6ae99SBarry Smith if (refine_y > 0) dd->refine_y = refine_y; 79047c6ae99SBarry Smith if (refine_z > 0) dd->refine_z = refine_z; 7913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 79447c6ae99SBarry Smith /*@C 795dce8aebaSBarry Smith DMDAGetRefinementFactor - Gets the ratios that the `DMDA` grid is refined 79647c6ae99SBarry Smith 79747c6ae99SBarry Smith Not Collective 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith Input Parameter: 800dce8aebaSBarry Smith . da - the `DMDA` object 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith Output Parameters: 80347c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default) 80447c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default) 80547c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default) 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith Level: intermediate 80847c6ae99SBarry Smith 809dce8aebaSBarry Smith Note: 81020f4b53cSBarry Smith Pass `NULL` for values you do not need 81147c6ae99SBarry Smith 81212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDASetRefinementFactor()` 81347c6ae99SBarry Smith @*/ 814d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y, PetscInt *refine_z) 815d71ae5a4SJacob Faibussowitsch { 81647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 81747c6ae99SBarry Smith 81847c6ae99SBarry Smith PetscFunctionBegin; 819a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 82047c6ae99SBarry Smith if (refine_x) *refine_x = dd->refine_x; 82147c6ae99SBarry Smith if (refine_y) *refine_y = dd->refine_y; 82247c6ae99SBarry Smith if (refine_z) *refine_z = dd->refine_z; 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82447c6ae99SBarry Smith } 82547c6ae99SBarry Smith 82647c6ae99SBarry Smith /*@C 827dce8aebaSBarry Smith DMDASetGetMatrix - Sets the routine used by the `DMDA` to allocate a matrix. 82847c6ae99SBarry Smith 82920f4b53cSBarry Smith Logically Collective; No Fortran Support 83047c6ae99SBarry Smith 83147c6ae99SBarry Smith Input Parameters: 832dce8aebaSBarry Smith + da - the `DMDA` object 83312b4a537SBarry Smith - f - the function that allocates the matrix for that specific `DMDA` 83412b4a537SBarry Smith 83512b4a537SBarry Smith Calling sequence of `f`: 83612b4a537SBarry Smith + da - the `DMDA` object 83712b4a537SBarry Smith - A - the created matrix 83847c6ae99SBarry Smith 83947c6ae99SBarry Smith Level: developer 84047c6ae99SBarry Smith 84112b4a537SBarry Smith Notes: 84212b4a537SBarry Smith If the function is not provided a default function is used that uses the `DMDAStencilType`, `DMBoundaryType`, and value of `DMDASetStencilWidth()` 84312b4a537SBarry Smith to construct the matrix. 84447c6ae99SBarry Smith 84512b4a537SBarry Smith See `DMDASetBlockFills()` that provides a simple way to provide the nonzero structure for 84612b4a537SBarry Smith the diagonal and off-diagonal blocks of the matrix without providing a custom function 84712b4a537SBarry Smith 84812b4a537SBarry Smith Developer Note: 84912b4a537SBarry Smith This should be called `DMDASetCreateMatrix()` 85012b4a537SBarry Smith 85112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateMatrix()`, `DMDASetBlockFills()` 85247c6ae99SBarry Smith @*/ 85312b4a537SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM da, Mat *A)) 854d71ae5a4SJacob Faibussowitsch { 85547c6ae99SBarry Smith PetscFunctionBegin; 856a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 85725296bd5SBarry Smith da->ops->creatematrix = f; 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85947c6ae99SBarry Smith } 86047c6ae99SBarry Smith 86138fb4e8eSJunchao Zhang /*@ 862dce8aebaSBarry Smith DMDAMapMatStencilToGlobal - Map a list of `MatStencil` on a grid to global indices. 86338fb4e8eSJunchao Zhang 86438fb4e8eSJunchao Zhang Not Collective 86538fb4e8eSJunchao Zhang 86638fb4e8eSJunchao Zhang Input Parameters: 867dce8aebaSBarry Smith + da - the `DMDA` object 86812b4a537SBarry Smith . m - number of `MatStencil` to map 86938fb4e8eSJunchao Zhang - idxm - grid points (and component number when dof > 1) 87038fb4e8eSJunchao Zhang 8711179163eSBarry Smith Output Parameter: 8721179163eSBarry Smith . gidxm - global row indices 8731179163eSBarry Smith 8741179163eSBarry Smith Level: intermediate 87538fb4e8eSJunchao Zhang 87612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `MatStencil` 87738fb4e8eSJunchao Zhang @*/ 878d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[]) 879d71ae5a4SJacob Faibussowitsch { 88038fb4e8eSJunchao Zhang const DM_DA *dd = (const DM_DA *)da->data; 88138fb4e8eSJunchao Zhang const PetscInt *dxm = (const PetscInt *)idxm; 88238fb4e8eSJunchao Zhang PetscInt i, j, sdim, tmp, dim; 88338fb4e8eSJunchao Zhang PetscInt dims[4], starts[4], dims2[3], starts2[3], dof = dd->w; 88438fb4e8eSJunchao Zhang ISLocalToGlobalMapping ltog; 88538fb4e8eSJunchao Zhang 88638fb4e8eSJunchao Zhang PetscFunctionBegin; 8873ba16761SJacob Faibussowitsch if (m <= 0) PetscFunctionReturn(PETSC_SUCCESS); 88838fb4e8eSJunchao Zhang 8892f27f4d1SJunchao Zhang /* Code adapted from DMDAGetGhostCorners() */ 89038fb4e8eSJunchao Zhang starts2[0] = dd->Xs / dof + dd->xo; 89138fb4e8eSJunchao Zhang starts2[1] = dd->Ys + dd->yo; 89238fb4e8eSJunchao Zhang starts2[2] = dd->Zs + dd->zo; 89338fb4e8eSJunchao Zhang dims2[0] = (dd->Xe - dd->Xs) / dof; 89438fb4e8eSJunchao Zhang dims2[1] = (dd->Ye - dd->Ys); 89538fb4e8eSJunchao Zhang dims2[2] = (dd->Ze - dd->Zs); 89638fb4e8eSJunchao Zhang 8972f27f4d1SJunchao Zhang /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */ 8982f27f4d1SJunchao Zhang dim = da->dim; /* DA dim: 1 to 3 */ 8992f27f4d1SJunchao Zhang sdim = dim + (dof > 1 ? 1 : 0); /* Dimensions in MatStencil's (k,j,i,c) view */ 9002f27f4d1SJunchao Zhang for (i = 0; i < dim; i++) { /* Reverse the order and also skip the unused dimensions */ 9012f27f4d1SJunchao 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 */ 90238fb4e8eSJunchao Zhang starts[i] = starts2[dim - i - 1]; 90338fb4e8eSJunchao Zhang } 9042f27f4d1SJunchao Zhang starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */ 90538fb4e8eSJunchao Zhang dims[dim] = dof; 90638fb4e8eSJunchao Zhang 90738fb4e8eSJunchao Zhang /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */ 90838fb4e8eSJunchao Zhang for (i = 0; i < m; i++) { 9092f27f4d1SJunchao Zhang dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */ 9102f27f4d1SJunchao Zhang tmp = 0; 9112f27f4d1SJunchao Zhang for (j = 0; j < sdim; j++) { /* Iter over, ex. j,i or j,i,c in 2D */ 9122f27f4d1SJunchao 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 */ 9132f27f4d1SJunchao Zhang else tmp = tmp * dims[j] + (dxm[j] - starts[j]); 91438fb4e8eSJunchao Zhang } 91538fb4e8eSJunchao Zhang gidxm[i] = tmp; 9162f27f4d1SJunchao Zhang /* Move to the next MatStencil point */ 9172f27f4d1SJunchao Zhang if (dof > 1) dxm += sdim; /* c is already counted in sdim */ 9182f27f4d1SJunchao Zhang else dxm += sdim + 1; /* skip the unused c */ 91938fb4e8eSJunchao Zhang } 92038fb4e8eSJunchao Zhang 92138fb4e8eSJunchao Zhang /* Map local indices to global indices */ 92238fb4e8eSJunchao Zhang PetscCall(DMGetLocalToGlobalMapping(da, <og)); 92338fb4e8eSJunchao Zhang PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm)); 9243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92538fb4e8eSJunchao Zhang } 92638fb4e8eSJunchao Zhang 92747c6ae99SBarry Smith /* 92847c6ae99SBarry Smith Creates "balanced" ownership ranges after refinement, constrained by the need for the 92947c6ae99SBarry Smith fine grid boundaries to fall within one stencil width of the coarse partition. 93047c6ae99SBarry Smith 93147c6ae99SBarry Smith Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 93247c6ae99SBarry Smith */ 933d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[]) 934d71ae5a4SJacob Faibussowitsch { 93547c6ae99SBarry Smith PetscInt i, totalc = 0, remaining, startc = 0, startf = 0; 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith PetscFunctionBegin; 93863a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 93947c6ae99SBarry Smith if (ratio == 1) { 9409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lf, lc, m)); 9413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 94247c6ae99SBarry Smith } 94347c6ae99SBarry Smith for (i = 0; i < m; i++) totalc += lc[i]; 94447c6ae99SBarry Smith remaining = (!periodic) + ratio * (totalc - (!periodic)); 94547c6ae99SBarry Smith for (i = 0; i < m; i++) { 94647c6ae99SBarry Smith PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 94747c6ae99SBarry Smith if (i == m - 1) lf[i] = want; 94847c6ae99SBarry Smith else { 9497aca7175SJed Brown const PetscInt nextc = startc + lc[i]; 9507aca7175SJed Brown /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 9517aca7175SJed Brown * coarse stencil width of the first coarse node in the next subdomain. */ 9527aca7175SJed Brown while ((startf + want) / ratio < nextc - stencil_width) want++; 9537aca7175SJed Brown /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 9547aca7175SJed Brown * coarse stencil width of the last coarse node in the current subdomain. */ 9557aca7175SJed Brown while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--; 9567aca7175SJed Brown /* Make sure all constraints are satisfied */ 9579371c9d4SSatish Balay if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width)) 9589371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range"); 95947c6ae99SBarry Smith } 96047c6ae99SBarry Smith lf[i] = want; 96147c6ae99SBarry Smith startc += lc[i]; 96247c6ae99SBarry Smith startf += lf[i]; 96347c6ae99SBarry Smith remaining -= lf[i]; 96447c6ae99SBarry Smith } 9653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96647c6ae99SBarry Smith } 96747c6ae99SBarry Smith 9682be375d4SJed Brown /* 9692be375d4SJed Brown Creates "balanced" ownership ranges after coarsening, constrained by the need for the 9702be375d4SJed Brown fine grid boundaries to fall within one stencil width of the coarse partition. 9712be375d4SJed Brown 9722be375d4SJed Brown Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 9732be375d4SJed Brown */ 974d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[]) 975d71ae5a4SJacob Faibussowitsch { 9762be375d4SJed Brown PetscInt i, totalf, remaining, startc, startf; 9772be375d4SJed Brown 9782be375d4SJed Brown PetscFunctionBegin; 97963a3b9bcSJacob Faibussowitsch PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio); 9802be375d4SJed Brown if (ratio == 1) { 9819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lc, lf, m)); 9823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9832be375d4SJed Brown } 9842be375d4SJed Brown for (i = 0, totalf = 0; i < m; i++) totalf += lf[i]; 9852be375d4SJed Brown remaining = (!periodic) + (totalf - (!periodic)) / ratio; 9862be375d4SJed Brown for (i = 0, startc = 0, startf = 0; i < m; i++) { 9872be375d4SJed Brown PetscInt want = remaining / (m - i) + !!(remaining % (m - i)); 9882be375d4SJed Brown if (i == m - 1) lc[i] = want; 9892be375d4SJed Brown else { 9902be375d4SJed Brown const PetscInt nextf = startf + lf[i]; 9912be375d4SJed Brown /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 9922be375d4SJed Brown * node is within one stencil width. */ 9932be375d4SJed Brown while (nextf / ratio < startc + want - stencil_width) want--; 9942be375d4SJed Brown /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 9952be375d4SJed Brown * fine node is within one stencil width. */ 9962be375d4SJed Brown while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++; 9979371c9d4SSatish Balay if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width)) 9989371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range"); 9992be375d4SJed Brown } 10002be375d4SJed Brown lc[i] = want; 10012be375d4SJed Brown startc += lc[i]; 10022be375d4SJed Brown startf += lf[i]; 10032be375d4SJed Brown remaining -= lc[i]; 10042be375d4SJed Brown } 10053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10062be375d4SJed Brown } 10072be375d4SJed Brown 1008d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref) 1009d71ae5a4SJacob Faibussowitsch { 1010c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 10116858538eSMatthew G. Knepley Vec coordsc, coordsf; 10129a42bb27SBarry Smith DM da2; 101347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data, *dd2; 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith PetscFunctionBegin; 1016a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 10174f572ea9SToby Isaac PetscAssertPointer(daref, 3); 101847c6ae99SBarry Smith 10199566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da, &dim)); 1020bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 102147c6ae99SBarry Smith M = dd->refine_x * dd->M; 102247c6ae99SBarry Smith } else { 102347c6ae99SBarry Smith M = 1 + dd->refine_x * (dd->M - 1); 102447c6ae99SBarry Smith } 1025bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1026c73cfb54SMatthew G. Knepley if (dim > 1) { 102747c6ae99SBarry Smith N = dd->refine_y * dd->N; 102847c6ae99SBarry Smith } else { 10291860e6e9SBarry Smith N = 1; 10301860e6e9SBarry Smith } 10311860e6e9SBarry Smith } else { 103247c6ae99SBarry Smith N = 1 + dd->refine_y * (dd->N - 1); 103347c6ae99SBarry Smith } 1034bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1035c73cfb54SMatthew G. Knepley if (dim > 2) { 103647c6ae99SBarry Smith P = dd->refine_z * dd->P; 103747c6ae99SBarry Smith } else { 10381860e6e9SBarry Smith P = 1; 10391860e6e9SBarry Smith } 10401860e6e9SBarry Smith } else { 104147c6ae99SBarry Smith P = 1 + dd->refine_z * (dd->P - 1); 104247c6ae99SBarry Smith } 10439566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2)); 10449566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix)); 10459566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da2, dim)); 10469566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da2, M, N, P)); 10479566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p)); 10489566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz)); 10499566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da2, dd->w)); 10509566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da2, dd->stencil_type)); 10519566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da2, dd->s)); 1052c73cfb54SMatthew G. Knepley if (dim == 3) { 105347c6ae99SBarry Smith PetscInt *lx, *ly, *lz; 10549566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 10559566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10569566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10579566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz)); 10589566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz)); 10599566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1060c73cfb54SMatthew G. Knepley } else if (dim == 2) { 106147c6ae99SBarry Smith PetscInt *lx, *ly; 10629566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 10639566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10649566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly)); 10659566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL)); 10669566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1067c73cfb54SMatthew G. Knepley } else if (dim == 1) { 106847c6ae99SBarry Smith PetscInt *lx; 10699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 10709566063dSJacob Faibussowitsch PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx)); 10719566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL)); 10729566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 107347c6ae99SBarry Smith } 107447c6ae99SBarry Smith dd2 = (DM_DA *)da2->data; 107547c6ae99SBarry Smith 107647c6ae99SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 107725296bd5SBarry Smith da2->ops->creatematrix = da->ops->creatematrix; 107825296bd5SBarry Smith /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 107947c6ae99SBarry Smith da2->ops->getcoloring = da->ops->getcoloring; 108047c6ae99SBarry Smith dd2->interptype = dd->interptype; 108147c6ae99SBarry Smith 108247c6ae99SBarry Smith /* copy fill information if given */ 108347c6ae99SBarry Smith if (dd->dfill) { 10849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 10859566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 108647c6ae99SBarry Smith } 108747c6ae99SBarry Smith if (dd->ofill) { 10889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 10899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith /* copy the refine information */ 1092397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1093397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1094397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->refine_z; 109547c6ae99SBarry Smith 1096897f7067SBarry Smith if (dd->refine_z_hier) { 1097ad540459SPierre 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]; 1098ad540459SPierre 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]; 1099897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 11009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 11019566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1102897f7067SBarry Smith } 1103897f7067SBarry Smith if (dd->refine_y_hier) { 1104ad540459SPierre 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]; 1105ad540459SPierre 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]; 1106897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 11079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 11089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1109897f7067SBarry Smith } 1110897f7067SBarry Smith if (dd->refine_x_hier) { 1111ad540459SPierre 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]; 1112ad540459SPierre 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]; 1113897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 11149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 11159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1116897f7067SBarry Smith } 1117897f7067SBarry Smith 111847c6ae99SBarry Smith /* copy vector type information */ 11199566063dSJacob Faibussowitsch PetscCall(DMSetVecType(da2, da->vectype)); 1120ddcf8b74SDave May 1121efd51863SBarry Smith dd2->lf = dd->lf; 1122efd51863SBarry Smith dd2->lj = dd->lj; 1123efd51863SBarry Smith 11246e87535bSJed Brown da2->leveldown = da->leveldown; 11256e87535bSJed Brown da2->levelup = da->levelup + 1; 11268865f1eaSKarl Rupp 11279566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2)); 11286e87535bSJed Brown 1129ddcf8b74SDave May /* interpolate coordinates if they are set on the coarse grid */ 11306858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(da, &coordsc)); 11316858538eSMatthew G. Knepley if (coordsc) { 1132ddcf8b74SDave May DM cdaf, cdac; 1133ddcf8b74SDave May Mat II; 1134ddcf8b74SDave May 11359566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cdac)); 11369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da2, &cdaf)); 1137b61d3410SDave May /* force creation of the coordinate vector */ 11389566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 11399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da2, &coordsf)); 11409566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL)); 11419566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 11429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 1143ddcf8b74SDave May } 1144397b6216SJed Brown 1145f3141302SJed Brown for (i = 0; i < da->bs; i++) { 1146f3141302SJed Brown const char *fieldname; 11479566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, i, &fieldname)); 11489566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da2, i, fieldname)); 1149f3141302SJed Brown } 1150397b6216SJed Brown 115147c6ae99SBarry Smith *daref = da2; 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115347c6ae99SBarry Smith } 115447c6ae99SBarry Smith 1155d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc) 1156d71ae5a4SJacob Faibussowitsch { 1157c73cfb54SMatthew G. Knepley PetscInt M, N, P, i, dim; 11586858538eSMatthew G. Knepley Vec coordsc, coordsf; 1159a5bc1bf3SBarry Smith DM dmc2; 1160a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA *)dmf->data, *dd2; 116147c6ae99SBarry Smith 116247c6ae99SBarry Smith PetscFunctionBegin; 1163a5bc1bf3SBarry Smith PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA); 11644f572ea9SToby Isaac PetscAssertPointer(dmc, 3); 116547c6ae99SBarry Smith 11669566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmf, &dim)); 1167bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1168397b6216SJed Brown M = dd->M / dd->coarsen_x; 116947c6ae99SBarry Smith } else { 1170397b6216SJed Brown M = 1 + (dd->M - 1) / dd->coarsen_x; 117147c6ae99SBarry Smith } 1172bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1173c73cfb54SMatthew G. Knepley if (dim > 1) { 1174397b6216SJed Brown N = dd->N / dd->coarsen_y; 117547c6ae99SBarry Smith } else { 11761860e6e9SBarry Smith N = 1; 11771860e6e9SBarry Smith } 11781860e6e9SBarry Smith } else { 1179397b6216SJed Brown N = 1 + (dd->N - 1) / dd->coarsen_y; 118047c6ae99SBarry Smith } 1181bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1182c73cfb54SMatthew G. Knepley if (dim > 2) { 1183397b6216SJed Brown P = dd->P / dd->coarsen_z; 118447c6ae99SBarry Smith } else { 11851860e6e9SBarry Smith P = 1; 11861860e6e9SBarry Smith } 11871860e6e9SBarry Smith } else { 1188397b6216SJed Brown P = 1 + (dd->P - 1) / dd->coarsen_z; 118947c6ae99SBarry Smith } 11909566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2)); 11919566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix)); 11929566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dmc2, dim)); 11939566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(dmc2, M, N, P)); 11949566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p)); 11959566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz)); 11969566063dSJacob Faibussowitsch PetscCall(DMDASetDof(dmc2, dd->w)); 11979566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(dmc2, dd->stencil_type)); 11989566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(dmc2, dd->s)); 1199c73cfb54SMatthew G. Knepley if (dim == 3) { 12002be375d4SJed Brown PetscInt *lx, *ly, *lz; 12019566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz)); 12029566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 12039566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 12049566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz)); 12059566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz)); 12069566063dSJacob Faibussowitsch PetscCall(PetscFree3(lx, ly, lz)); 1207c73cfb54SMatthew G. Knepley } else if (dim == 2) { 12082be375d4SJed Brown PetscInt *lx, *ly; 12099566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly)); 12109566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 12119566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly)); 12129566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL)); 12139566063dSJacob Faibussowitsch PetscCall(PetscFree2(lx, ly)); 1214c73cfb54SMatthew G. Knepley } else if (dim == 1) { 12152be375d4SJed Brown PetscInt *lx; 12169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->m, &lx)); 12179566063dSJacob Faibussowitsch PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx)); 12189566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL)); 12199566063dSJacob Faibussowitsch PetscCall(PetscFree(lx)); 122047c6ae99SBarry Smith } 1221a5bc1bf3SBarry Smith dd2 = (DM_DA *)dmc2->data; 122247c6ae99SBarry Smith 12234dcab191SBarry Smith /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 1224a5bc1bf3SBarry Smith /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 1225a5bc1bf3SBarry Smith dmc2->ops->creatematrix = dmf->ops->creatematrix; 1226a5bc1bf3SBarry Smith dmc2->ops->getcoloring = dmf->ops->getcoloring; 122747c6ae99SBarry Smith dd2->interptype = dd->interptype; 122847c6ae99SBarry Smith 122947c6ae99SBarry Smith /* copy fill information if given */ 123047c6ae99SBarry Smith if (dd->dfill) { 12319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill)); 12329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1)); 123347c6ae99SBarry Smith } 123447c6ae99SBarry Smith if (dd->ofill) { 12359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill)); 12369566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1)); 123747c6ae99SBarry Smith } 123847c6ae99SBarry Smith /* copy the refine information */ 1239397b6216SJed Brown dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1240397b6216SJed Brown dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1241397b6216SJed Brown dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 124247c6ae99SBarry Smith 1243897f7067SBarry Smith if (dd->refine_z_hier) { 1244ad540459SPierre 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]; 1245ad540459SPierre 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]; 1246897f7067SBarry Smith dd2->refine_z_hier_n = dd->refine_z_hier_n; 12479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier)); 12489566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n)); 1249897f7067SBarry Smith } 1250897f7067SBarry Smith if (dd->refine_y_hier) { 1251ad540459SPierre 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]; 1252ad540459SPierre 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]; 1253897f7067SBarry Smith dd2->refine_y_hier_n = dd->refine_y_hier_n; 12549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier)); 12559566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n)); 1256897f7067SBarry Smith } 1257897f7067SBarry Smith if (dd->refine_x_hier) { 1258ad540459SPierre 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]; 1259ad540459SPierre 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]; 1260897f7067SBarry Smith dd2->refine_x_hier_n = dd->refine_x_hier_n; 12619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier)); 12629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n)); 1263897f7067SBarry Smith } 1264897f7067SBarry Smith 126547c6ae99SBarry Smith /* copy vector type information */ 12669566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dmc2, dmf->vectype)); 126747c6ae99SBarry Smith 1268644e2e5bSBarry Smith dd2->lf = dd->lf; 1269644e2e5bSBarry Smith dd2->lj = dd->lj; 1270644e2e5bSBarry Smith 1271a5bc1bf3SBarry Smith dmc2->leveldown = dmf->leveldown + 1; 1272a5bc1bf3SBarry Smith dmc2->levelup = dmf->levelup; 12738865f1eaSKarl Rupp 12749566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmc2)); 12756e87535bSJed Brown 1276ddcf8b74SDave May /* inject coordinates if they are set on the fine grid */ 12776858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dmf, &coordsf)); 12786858538eSMatthew G. Knepley if (coordsf) { 1279ddcf8b74SDave May DM cdaf, cdac; 12806dbf9973SLawrence Mitchell Mat inject; 12816dbf9973SLawrence Mitchell VecScatter vscat; 1282ddcf8b74SDave May 12839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmf, &cdaf)); 12849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmc2, &cdac)); 1285b61d3410SDave May /* force creation of the coordinate vector */ 12869566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 12879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmc2, &coordsc)); 1288ddcf8b74SDave May 12899566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(cdac, cdaf, &inject)); 12909566063dSJacob Faibussowitsch PetscCall(MatScatterGetVecScatter(inject, &vscat)); 12919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12929566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD)); 12939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inject)); 1294ddcf8b74SDave May } 1295f98405f7SJed Brown 1296a5bc1bf3SBarry Smith for (i = 0; i < dmf->bs; i++) { 1297f98405f7SJed Brown const char *fieldname; 12989566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(dmf, i, &fieldname)); 12999566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(dmc2, i, fieldname)); 1300f98405f7SJed Brown } 1301f98405f7SJed Brown 1302a5bc1bf3SBarry Smith *dmc = dmc2; 13033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130447c6ae99SBarry Smith } 130547c6ae99SBarry Smith 1306d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[]) 1307d71ae5a4SJacob Faibussowitsch { 130847c6ae99SBarry Smith PetscInt i, n, *refx, *refy, *refz; 130947c6ae99SBarry Smith 131047c6ae99SBarry Smith PetscFunctionBegin; 131147c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 13127a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 13133ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 13144f572ea9SToby Isaac PetscAssertPointer(daf, 3); 131547c6ae99SBarry Smith 1316aa219208SBarry Smith /* Get refinement factors, defaults taken from the coarse DMDA */ 13179566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz)); 131848a46eb9SPierre Jolivet for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i])); 131947c6ae99SBarry Smith n = nlevels; 13209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL)); 132147c6ae99SBarry Smith n = nlevels; 13229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL)); 132347c6ae99SBarry Smith n = nlevels; 13249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL)); 132547c6ae99SBarry Smith 13269566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0])); 13279566063dSJacob Faibussowitsch PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0])); 132847c6ae99SBarry Smith for (i = 1; i < nlevels; i++) { 13299566063dSJacob Faibussowitsch PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i])); 13309566063dSJacob Faibussowitsch PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i])); 133147c6ae99SBarry Smith } 13329566063dSJacob Faibussowitsch PetscCall(PetscFree3(refx, refy, refz)); 13333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133447c6ae99SBarry Smith } 133547c6ae99SBarry Smith 1336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[]) 1337d71ae5a4SJacob Faibussowitsch { 133847c6ae99SBarry Smith PetscInt i; 133947c6ae99SBarry Smith 134047c6ae99SBarry Smith PetscFunctionBegin; 134147c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 13427a8be351SBarry Smith PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative"); 13433ba16761SJacob Faibussowitsch if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS); 13444f572ea9SToby Isaac PetscAssertPointer(dac, 3); 13459566063dSJacob Faibussowitsch PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0])); 134648a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i])); 13473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134847c6ae99SBarry Smith } 134962197512SBarry Smith 135066976f2fSJacob Faibussowitsch static PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes) 1351d71ae5a4SJacob Faibussowitsch { 13528272889dSSatish Balay PetscInt i, j, xs, xn, q; 135362197512SBarry Smith PetscScalar *xx; 135462197512SBarry Smith PetscReal h; 135562197512SBarry Smith Vec x; 135662197512SBarry Smith DM_DA *da = (DM_DA *)dm->data; 135762197512SBarry Smith 135862197512SBarry Smith PetscFunctionBegin; 135962197512SBarry Smith if (da->bx != DM_BOUNDARY_PERIODIC) { 13609566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 136162197512SBarry Smith q = (q - 1) / (n - 1); /* number of spectral elements */ 136262197512SBarry Smith h = 2.0 / q; 13639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL)); 136462197512SBarry Smith xs = xs / (n - 1); 136562197512SBarry Smith xn = xn / (n - 1); 13669566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.)); 13679566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &x)); 13689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, x, &xx)); 136962197512SBarry Smith 137062197512SBarry Smith /* loop over local spectral elements */ 137162197512SBarry Smith for (j = xs; j < xs + xn; j++) { 137262197512SBarry Smith /* 137362197512SBarry Smith Except for the first process, each process starts on the second GLL point of the first element on that process 137462197512SBarry Smith */ 1375ad540459SPierre 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.; 137662197512SBarry Smith } 13779566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, x, &xx)); 137862197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic"); 13793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138062197512SBarry Smith } 138162197512SBarry Smith 13821fd49c25SBarry Smith /*@ 138362197512SBarry Smith 138462197512SBarry 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 138562197512SBarry Smith 138620f4b53cSBarry Smith Collective 138762197512SBarry Smith 138862197512SBarry Smith Input Parameters: 1389dce8aebaSBarry Smith + da - the `DMDA` object 13902fe279fdSBarry Smith . n - the number of GLL nodes 13918272889dSSatish Balay - nodes - the GLL nodes 139262197512SBarry Smith 1393edc382c3SSatish Balay Level: advanced 1394edc382c3SSatish Balay 1395dce8aebaSBarry Smith Note: 1396dce8aebaSBarry Smith The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 139712b4a537SBarry Smith on each process much be divisible by the number of GLL elements needed per process. This depends on whether the `DMDA` is 1398dce8aebaSBarry Smith periodic or not. 1399dce8aebaSBarry Smith 140012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()` 140162197512SBarry Smith @*/ 1402d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes) 1403d71ae5a4SJacob Faibussowitsch { 140462197512SBarry Smith PetscFunctionBegin; 140562197512SBarry Smith if (da->dim == 1) { 14069566063dSJacob Faibussowitsch PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes)); 140762197512SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d"); 14083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140962197512SBarry Smith } 14107c3cd84eSPatrick Sanan 141166976f2fSJacob Faibussowitsch PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set) 1412d71ae5a4SJacob Faibussowitsch { 14137c3cd84eSPatrick Sanan DM_DA *dd1 = (DM_DA *)da1->data, *dd2; 14147c3cd84eSPatrick Sanan DM da2; 14157c3cd84eSPatrick Sanan DMType dmtype2; 14167c3cd84eSPatrick Sanan PetscBool isda, compatibleLocal; 14177c3cd84eSPatrick Sanan PetscInt i; 14187c3cd84eSPatrick Sanan 14197c3cd84eSPatrick Sanan PetscFunctionBegin; 14207a8be351SBarry Smith PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()"); 14219566063dSJacob Faibussowitsch PetscCall(DMGetType(dm2, &dmtype2)); 14229566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(dmtype2, DMDA, &isda)); 14237c3cd84eSPatrick Sanan if (isda) { 14247c3cd84eSPatrick Sanan da2 = dm2; 14257c3cd84eSPatrick Sanan dd2 = (DM_DA *)da2->data; 14267a8be351SBarry Smith PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()"); 1427110623a0SKarl Rupp compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1428c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1429110623a0SKarl Rupp /* Global size ranks Boundary type */ 1430c790a739SKarl Rupp if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1431c790a739SKarl Rupp if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1432c790a739SKarl Rupp if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 14337c3cd84eSPatrick Sanan if (compatibleLocal) { 14349371c9d4SSatish Balay for (i = 0; i < dd1->m; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ } 14357c3cd84eSPatrick Sanan } 14367c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 1) { 1437ad540459SPierre Jolivet for (i = 0; i < dd1->n; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 14387c3cd84eSPatrick Sanan } 14397c3cd84eSPatrick Sanan if (compatibleLocal && da1->dim > 2) { 1440ad540459SPierre Jolivet for (i = 0; i < dd1->p; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 14417c3cd84eSPatrick Sanan } 14427c3cd84eSPatrick Sanan *compatible = compatibleLocal; 14437c3cd84eSPatrick Sanan *set = PETSC_TRUE; 14447c3cd84eSPatrick Sanan } else { 14457c3cd84eSPatrick Sanan /* Decline to determine compatibility with other DM types */ 14467c3cd84eSPatrick Sanan *set = PETSC_FALSE; 14477c3cd84eSPatrick Sanan } 14483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14497c3cd84eSPatrick Sanan } 1450