1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 3ba38deedSJacob Faibussowitsch static PetscErrorCode DMSetFromOptions_DA(DM da, PetscOptionItems *PetscOptionsObject) 4d71ae5a4SJacob Faibussowitsch { 547c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 6c73cfb54SMatthew G. Knepley PetscInt refine = 0, dim = da->dim, maxnlevels = 100, refx[100], refy[100], refz[100], n, i; 7fe58071bSMatthew G. Knepley DMBoundaryType bt = DM_BOUNDARY_NONE; 8897f7067SBarry Smith PetscBool flg; 947c6ae99SBarry Smith 1047c6ae99SBarry Smith PetscFunctionBegin; 118e704042SBarry Smith PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot call DMSetFromOptions() after DMSetUp() for DMDA"); 1208401ef6SPierre Jolivet PetscCheck(dd->M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime"); 1308401ef6SPierre Jolivet PetscCheck(dd->N >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime"); 1408401ef6SPierre Jolivet PetscCheck(dd->P >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime"); 15235683edSBarry Smith 16d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMDA Options"); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_grid_x", "Number of grid points in x direction", "DMDASetSizes", dd->M, &dd->M, NULL, 1)); 189566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y", "Number of grid points in y direction", "DMDASetSizes", dd->N, &dd->N, NULL, 1)); 199566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z", "Number of grid points in z direction", "DMDASetSizes", dd->P, &dd->P, NULL, 1)); 207ddda789SPeter Brune 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap", "Decomposition overlap in all directions", "DMDASetOverlap", dd->xol, &dd->xol, &flg, 0)); 229566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetOverlap(da, dd->xol, dd->xol, dd->xol)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap_x", "Decomposition overlap in x direction", "DMDASetOverlap", dd->xol, &dd->xol, NULL, 0)); 249566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y", "Decomposition overlap in y direction", "DMDASetOverlap", dd->yol, &dd->yol, NULL, 0)); 259566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z", "Decomposition overlap in z direction", "DMDASetOverlap", dd->zol, &dd->zol, NULL, 0)); 263e7870d2SPeter Brune 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_local_subdomains", "", "DMDASetNumLocalSubdomains", dd->Nsub, &dd->Nsub, &flg, PETSC_DECIDE)); 289566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetNumLocalSubDomains(da, dd->Nsub)); 293e7870d2SPeter Brune 308d07fd27SPatrick Sanan /* Handle DMDA parallel distribution */ 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_processors_x", "Number of processors in x direction", "DMDASetNumProcs", dd->m, &dd->m, NULL, PETSC_DECIDE)); 329566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y", "Number of processors in y direction", "DMDASetNumProcs", dd->n, &dd->n, NULL, PETSC_DECIDE)); 339566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z", "Number of processors in z direction", "DMDASetNumProcs", dd->p, &dd->p, NULL, PETSC_DECIDE)); 34fe58071bSMatthew G. Knepley // Handle boundaries 35fe58071bSMatthew G. Knepley PetscCall(PetscOptionsEnum("-da_bd_x", "Boundary type for x direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bx, (PetscEnum *)&dd->bx, NULL)); 36fe58071bSMatthew G. Knepley if (dim > 1) PetscCall(PetscOptionsEnum("-da_bd_y", "Boundary type for y direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->by, (PetscEnum *)&dd->by, NULL)); 37fe58071bSMatthew G. Knepley if (dim > 2) PetscCall(PetscOptionsEnum("-da_bd_z", "Boundary type for z direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bz, (PetscEnum *)&dd->bz, NULL)); 38fe58071bSMatthew G. Knepley PetscCall(PetscOptionsEnum("-da_bd_all", "Boundary type for every direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)bt, (PetscEnum *)&bt, &flg)); 39fe58071bSMatthew G. Knepley if (flg) PetscCall(DMDASetBoundaryType(da, bt, bt, bt)); 40aa219208SBarry Smith /* Handle DMDA refinement */ 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine_x", "Refinement ratio in x direction", "DMDASetRefinementFactor", dd->refine_x, &dd->refine_x, NULL, 1)); 429566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y", "Refinement ratio in y direction", "DMDASetRefinementFactor", dd->refine_y, &dd->refine_y, NULL, 1)); 439566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z", "Refinement ratio in z direction", "DMDASetRefinementFactor", dd->refine_z, &dd->refine_z, NULL, 1)); 449371c9d4SSatish Balay dd->coarsen_x = dd->refine_x; 459371c9d4SSatish Balay dd->coarsen_y = dd->refine_y; 469371c9d4SSatish Balay dd->coarsen_z = dd->refine_z; 4747c6ae99SBarry Smith 48397b6216SJed Brown /* Get refinement factors, defaults taken from the coarse DMDA */ 499566063dSJacob Faibussowitsch PetscCall(DMDAGetRefinementFactor(da, &refx[0], &refy[0], &refz[0])); 50397b6216SJed Brown for (i = 1; i < maxnlevels; i++) { 51397b6216SJed Brown refx[i] = refx[0]; 52397b6216SJed Brown refy[i] = refy[0]; 53397b6216SJed Brown refz[i] = refz[0]; 54397b6216SJed Brown } 55397b6216SJed Brown n = maxnlevels; 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_x", "Refinement factor for each level", "None", refx, &n, &flg)); 57897f7067SBarry Smith if (flg) { 58897f7067SBarry Smith dd->refine_x = refx[0]; 59897f7067SBarry Smith dd->refine_x_hier_n = n; 609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_x_hier)); 619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_x_hier, refx, n)); 62897f7067SBarry Smith } 63c73cfb54SMatthew G. Knepley if (dim > 1) { 64397b6216SJed Brown n = maxnlevels; 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_y", "Refinement factor for each level", "None", refy, &n, &flg)); 66897f7067SBarry Smith if (flg) { 67897f7067SBarry Smith dd->refine_y = refy[0]; 68897f7067SBarry Smith dd->refine_y_hier_n = n; 699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_y_hier)); 709566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_y_hier, refy, n)); 71897f7067SBarry Smith } 72397b6216SJed Brown } 73c73cfb54SMatthew G. Knepley if (dim > 2) { 74397b6216SJed Brown n = maxnlevels; 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_z", "Refinement factor for each level", "None", refz, &n, &flg)); 76897f7067SBarry Smith if (flg) { 77897f7067SBarry Smith dd->refine_z = refz[0]; 78897f7067SBarry Smith dd->refine_z_hier_n = n; 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_z_hier)); 809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_z_hier, refz, n)); 81897f7067SBarry Smith } 82397b6216SJed Brown } 83397b6216SJed Brown 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine", "Uniformly refine DA one or more times", "None", refine, &refine, NULL, 0)); 85d0609cedSBarry Smith PetscOptionsHeadEnd(); 86235683edSBarry Smith 87e0f5d30fSBarry Smith while (refine--) { 88bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 897a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M, &dd->M)); 90e0f5d30fSBarry Smith } else { 917a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M - 1, &dd->M)); 927a3b7fd5SBarry Smith dd->M += 1; 93e0f5d30fSBarry Smith } 9442e65cf2SBarry Smith if (dim > 1 && (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0)) { 957a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N, &dd->N)); 9642e65cf2SBarry Smith } else if (dim > 1) { 977a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N - 1, &dd->N)); 987a3b7fd5SBarry Smith dd->N += 1; 99e0f5d30fSBarry Smith } 10042e65cf2SBarry Smith if (dim > 2 && (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0)) { 1017a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P, &dd->P)); 10242e65cf2SBarry Smith } else if (dim > 2) { 1037a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P - 1, &dd->P)); 1047a3b7fd5SBarry Smith dd->P += 1; 105e0f5d30fSBarry Smith } 106e0f5d30fSBarry Smith da->levelup++; 10785fc4b34SJed Brown if (da->levelup - da->leveldown >= 0) { 10885fc4b34SJed Brown dd->refine_x = refx[da->levelup - da->leveldown]; 10985fc4b34SJed Brown dd->refine_y = refy[da->levelup - da->leveldown]; 11085fc4b34SJed Brown dd->refine_z = refz[da->levelup - da->leveldown]; 11185fc4b34SJed Brown } 11285fc4b34SJed Brown if (da->levelup - da->leveldown >= 1) { 11385fc4b34SJed Brown dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 11485fc4b34SJed Brown dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 11585fc4b34SJed Brown dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 11685fc4b34SJed Brown } 117e0f5d30fSBarry Smith } 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11947c6ae99SBarry Smith } 12047c6ae99SBarry Smith 121ba38deedSJacob Faibussowitsch static PetscErrorCode DMLoad_DA(DM da, PetscViewer viewer) 122d71ae5a4SJacob Faibussowitsch { 123060da220SMatthew G. Knepley PetscInt dim, m, n, p, dof, swidth; 124b859378eSBarry Smith DMDAStencilType stencil; 125bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 126bc2bf880SBarry Smith PetscBool coors; 127bc2bf880SBarry Smith DM dac; 128bc2bf880SBarry Smith Vec c; 129b859378eSBarry Smith 130b859378eSBarry Smith PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dim, 1, NULL, PETSC_INT)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &m, 1, NULL, PETSC_INT)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &p, 1, NULL, PETSC_INT)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dof, 1, NULL, PETSC_INT)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &swidth, 1, NULL, PETSC_INT)); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bx, 1, NULL, PETSC_ENUM)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &by, 1, NULL, PETSC_ENUM)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bz, 1, NULL, PETSC_ENUM)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &stencil, 1, NULL, PETSC_ENUM)); 141b859378eSBarry Smith 1429566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1439566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, m, n, p)); 1449566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, bx, by, bz)); 1459566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, dof)); 1469566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, stencil)); 1479566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, swidth)); 1489566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &coors, 1, NULL, PETSC_ENUM)); 150bc2bf880SBarry Smith if (coors) { 1519566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &dac)); 1529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &c)); 1539566063dSJacob Faibussowitsch PetscCall(VecLoad(c, viewer)); 1549566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, c)); 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 156bc2bf880SBarry Smith } 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158b859378eSBarry Smith } 159b859378eSBarry Smith 160ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 161d71ae5a4SJacob Faibussowitsch { 162d724dfffSBarry Smith DM_DA *da = (DM_DA *)dm->data; 163d724dfffSBarry Smith 164d724dfffSBarry Smith PetscFunctionBegin; 165d724dfffSBarry Smith if (subdm) { 166c38c1269SMatthew G. Knepley PetscSF sf; 167c38c1269SMatthew G. Knepley Vec coords; 168c38c1269SMatthew G. Knepley void *ctx; 169c38c1269SMatthew G. Knepley /* Cannot use DMClone since the dof stuff is mixed in. Ugh 1709566063dSJacob Faibussowitsch PetscCall(DMClone(dm, subdm)); */ 1719566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm)); 1729566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1739566063dSJacob Faibussowitsch PetscCall(DMSetPointSF(*subdm, sf)); 1749566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 1759566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*subdm, ctx)); 1769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 177c38c1269SMatthew G. Knepley if (coords) { 1789566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*subdm, coords)); 179c38c1269SMatthew G. Knepley } else { 1809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &coords)); 1819566063dSJacob Faibussowitsch if (coords) PetscCall(DMSetCoordinates(*subdm, coords)); 182c38c1269SMatthew G. Knepley } 183c38c1269SMatthew G. Knepley 1849566063dSJacob Faibussowitsch PetscCall(DMSetType(*subdm, DMDA)); 1859566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*subdm, dm->dim)); 1869566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P)); 1879566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p)); 1889566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz)); 1899566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*subdm, numFields)); 1909566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*subdm, da->stencil_type)); 1919566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*subdm, da->s)); 1929566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz)); 193d724dfffSBarry Smith } 194d724dfffSBarry Smith if (is) { 195d724dfffSBarry Smith PetscInt *indices, cnt = 0, dof = da->w, i, j; 19638221697SMatthew G. Knepley 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(da->Nlocal * numFields / dof, &indices)); 19838221697SMatthew G. Knepley for (i = da->base / dof; i < (da->base + da->Nlocal) / dof; ++i) { 199ad540459SPierre Jolivet for (j = 0; j < numFields; ++j) indices[cnt++] = dof * i + fields[j]; 200d724dfffSBarry Smith } 20163a3b9bcSJacob Faibussowitsch PetscCheck(cnt == da->Nlocal * numFields / dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %" PetscInt_FMT " does not equal expected value %" PetscInt_FMT, cnt, da->Nlocal * numFields / dof); 2029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), cnt, indices, PETSC_OWN_POINTER, is)); 203d724dfffSBarry Smith } 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205d724dfffSBarry Smith } 206d724dfffSBarry Smith 207ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 208d71ae5a4SJacob Faibussowitsch { 2091c3fb106SBarry Smith PetscInt i; 2101c3fb106SBarry Smith DM_DA *dd = (DM_DA *)dm->data; 2111c3fb106SBarry Smith PetscInt dof = dd->w; 2121c3fb106SBarry Smith 2131c3fb106SBarry Smith PetscFunctionBegin; 214731c8d9eSDmitry Karpeev if (len) *len = dof; 2151c3fb106SBarry Smith if (islist) { 2161c3fb106SBarry Smith Vec v; 2171c3fb106SBarry Smith PetscInt rstart, n; 2181c3fb106SBarry Smith 2199566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 2209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(v, &rstart, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 2229566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 2239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, islist)); 22448a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), n / dof, rstart + i, dof, &(*islist)[i])); 2251c3fb106SBarry Smith } 2261c3fb106SBarry Smith if (namelist) { 2279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, namelist)); 2281c3fb106SBarry Smith if (dd->fieldname) { 22948a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(PetscStrallocpy(dd->fieldname[i], &(*namelist)[i])); 2301c3fb106SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently DMDA must have fieldnames"); 2311c3fb106SBarry Smith } 2321c3fb106SBarry Smith if (dmlist) { 2331c3fb106SBarry Smith DM da; 2341c3fb106SBarry Smith 2359566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 2369566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dm->dim)); 2379566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P)); 2389566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p)); 2399566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz)); 2409566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 2419566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, dd->stencil_type)); 2429566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, dd->s)); 2439566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, dmlist)); 2459566063dSJacob Faibussowitsch for (i = 0; i < dof - 1; i++) PetscCall(PetscObjectReference((PetscObject)da)); 2461c3fb106SBarry Smith for (i = 0; i < dof; i++) (*dmlist)[i] = da; 2471c3fb106SBarry Smith } 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2491c3fb106SBarry Smith } 2501c3fb106SBarry Smith 251ba38deedSJacob Faibussowitsch static PetscErrorCode DMClone_DA(DM dm, DM *newdm) 252d71ae5a4SJacob Faibussowitsch { 25338221697SMatthew G. Knepley DM_DA *da = (DM_DA *)dm->data; 25438221697SMatthew G. Knepley 25538221697SMatthew G. Knepley PetscFunctionBegin; 2569566063dSJacob Faibussowitsch PetscCall(DMSetType(*newdm, DMDA)); 2579566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*newdm, dm->dim)); 2589566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P)); 2599566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p)); 2609566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz)); 2619566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*newdm, da->w)); 2629566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*newdm, da->stencil_type)); 2639566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*newdm, da->s)); 2649566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz)); 2659566063dSJacob Faibussowitsch PetscCall(DMSetUp(*newdm)); 2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26738221697SMatthew G. Knepley } 26838221697SMatthew G. Knepley 269d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg) 270d71ae5a4SJacob Faibussowitsch { 2714a7a4c06SLawrence Mitchell DM_DA *da = (DM_DA *)dm->data; 2724a7a4c06SLawrence Mitchell 2734a7a4c06SLawrence Mitchell PetscFunctionBegin; 2744a7a4c06SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2754f572ea9SToby Isaac PetscAssertPointer(flg, 2); 276131f56ebSLawrence Mitchell *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE; 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2784a7a4c06SLawrence Mitchell } 2794a7a4c06SLawrence Mitchell 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 281d71ae5a4SJacob Faibussowitsch { 282793f3fe5SMatthew G. Knepley PetscFunctionBegin; 2839566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd)); 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285793f3fe5SMatthew G. Knepley } 286793f3fe5SMatthew G. Knepley 287d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 288d71ae5a4SJacob Faibussowitsch { 289502a2867SDave May PetscInt dim; 290502a2867SDave May DMDAStencilType st; 291502a2867SDave May 292502a2867SDave May PetscFunctionBegin; 2939566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dm, ranks)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &st)); 295502a2867SDave May 296502a2867SDave May switch (dim) { 297502a2867SDave May case 1: 298502a2867SDave May *nranks = 3; 299ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 3; */ 300502a2867SDave May break; 301502a2867SDave May case 2: 302502a2867SDave May *nranks = 9; 303ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 5; */ 304502a2867SDave May break; 305502a2867SDave May case 3: 306502a2867SDave May *nranks = 27; 307ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 7; */ 308502a2867SDave May break; 309d71ae5a4SJacob Faibussowitsch default: 310d71ae5a4SJacob Faibussowitsch break; 311502a2867SDave May } 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313502a2867SDave May } 314502a2867SDave May 3153efe6655SBarry Smith /*MC 316*0b4b7b1cSBarry Smith DMDA = "da" - A `DM` object that is used to help solve PDEs on a structured grid (or mesh) in 1, 2, or 3 dimensions. 3173efe6655SBarry Smith 3183efe6655SBarry Smith Level: intermediate 3193efe6655SBarry Smith 320*0b4b7b1cSBarry Smith Notes: 321*0b4b7b1cSBarry Smith In the global representation of the vectors each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points. 322*0b4b7b1cSBarry Smith In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width set with `DMDASetStencilWidth()`. 323*0b4b7b1cSBarry Smith 324*0b4b7b1cSBarry Smith The vectors can be thought of as either cell centered or vertex centered on the grid (or mesh). But some variables cannot be cell centered and others 325*0b4b7b1cSBarry Smith vertex centered; see the documentation for `DMSTAG`, a similar `DM` implementation which supports more general staggered grids. 326*0b4b7b1cSBarry Smith 327*0b4b7b1cSBarry Smith Periodic boundary conditions can be handled by using a `DMBoundaryType` of `DM_BOUNDARY_PERIODIC` provided with `DMDASetBoundaryType()`. 328*0b4b7b1cSBarry Smith Other `DMBoundaryType`values allow for different handling of terms along the boundary of the grid (or mesh). 329*0b4b7b1cSBarry Smith 33012b4a537SBarry Smith .seealso: [](sec_struct), `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMDASetStencilWidth()`, `DMDASetStencilType()`, 33112b4a537SBarry Smith `DMDAStencilType` 3323efe6655SBarry Smith M*/ 3333efe6655SBarry Smith 334d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da) 335d71ae5a4SJacob Faibussowitsch { 33647c6ae99SBarry Smith DM_DA *dd; 33747c6ae99SBarry Smith 33847c6ae99SBarry Smith PetscFunctionBegin; 3394f572ea9SToby Isaac PetscAssertPointer(da, 1); 3404dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dd)); 341a4121054SBarry Smith da->data = dd; 34247c6ae99SBarry Smith 343c73cfb54SMatthew G. Knepley da->dim = -1; 344aa219208SBarry Smith dd->interptype = DMDA_Q1; 34547c6ae99SBarry Smith dd->refine_x = 2; 34647c6ae99SBarry Smith dd->refine_y = 2; 34747c6ae99SBarry Smith dd->refine_z = 2; 34881c108dcSJed Brown dd->coarsen_x = 2; 34981c108dcSJed Brown dd->coarsen_y = 2; 35081c108dcSJed Brown dd->coarsen_z = 2; 3510298fd71SBarry Smith dd->fieldname = NULL; 35247c6ae99SBarry Smith dd->nlocal = -1; 35347c6ae99SBarry Smith dd->Nlocal = -1; 35447c6ae99SBarry Smith dd->M = -1; 35547c6ae99SBarry Smith dd->N = -1; 35647c6ae99SBarry Smith dd->P = -1; 35747c6ae99SBarry Smith dd->m = -1; 35847c6ae99SBarry Smith dd->n = -1; 35947c6ae99SBarry Smith dd->p = -1; 36047c6ae99SBarry Smith dd->w = -1; 36147c6ae99SBarry Smith dd->s = -1; 3628865f1eaSKarl Rupp 3639371c9d4SSatish Balay dd->xs = -1; 3649371c9d4SSatish Balay dd->xe = -1; 3659371c9d4SSatish Balay dd->ys = -1; 3669371c9d4SSatish Balay dd->ye = -1; 3679371c9d4SSatish Balay dd->zs = -1; 3689371c9d4SSatish Balay dd->ze = -1; 3699371c9d4SSatish Balay dd->Xs = -1; 3709371c9d4SSatish Balay dd->Xe = -1; 3719371c9d4SSatish Balay dd->Ys = -1; 3729371c9d4SSatish Balay dd->Ye = -1; 3739371c9d4SSatish Balay dd->Zs = -1; 3749371c9d4SSatish Balay dd->Ze = -1; 37547c6ae99SBarry Smith 3763e7870d2SPeter Brune dd->Nsub = 1; 3777ddda789SPeter Brune dd->xol = 0; 3787ddda789SPeter Brune dd->yol = 0; 3797ddda789SPeter Brune dd->zol = 0; 380d886c4f4SPeter Brune dd->xo = 0; 381d886c4f4SPeter Brune dd->yo = 0; 382d886c4f4SPeter Brune dd->zo = 0; 383e30e807fSPeter Brune dd->Mo = -1; 384e30e807fSPeter Brune dd->No = -1; 385e30e807fSPeter Brune dd->Po = -1; 38688661749SPeter Brune 3870298fd71SBarry Smith dd->gtol = NULL; 3880298fd71SBarry Smith dd->ltol = NULL; 3890298fd71SBarry Smith dd->ao = NULL; 3903ba16761SJacob Faibussowitsch PetscCall(PetscStrallocpy(AOBASIC, (char **)&dd->aotype)); 39147c6ae99SBarry Smith dd->base = -1; 392bff4a2f0SMatthew G. Knepley dd->bx = DM_BOUNDARY_NONE; 393bff4a2f0SMatthew G. Knepley dd->by = DM_BOUNDARY_NONE; 394bff4a2f0SMatthew G. Knepley dd->bz = DM_BOUNDARY_NONE; 395aa219208SBarry Smith dd->stencil_type = DMDA_STENCIL_BOX; 396aa219208SBarry Smith dd->interptype = DMDA_Q1; 3970298fd71SBarry Smith dd->lx = NULL; 3980298fd71SBarry Smith dd->ly = NULL; 3990298fd71SBarry Smith dd->lz = NULL; 40047c6ae99SBarry Smith 401454e267fSLisandro Dalcin dd->elementtype = DMDA_ELEMENT_Q1; 402454e267fSLisandro Dalcin 403a4121054SBarry Smith da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA; 404a4121054SBarry Smith da->ops->globaltolocalend = DMGlobalToLocalEnd_DA; 405a4121054SBarry Smith da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA; 406a4121054SBarry Smith da->ops->localtoglobalend = DMLocalToGlobalEnd_DA; 407d78e899eSRichard Tran Mills da->ops->localtolocalbegin = DMLocalToLocalBegin_DA; 408d78e899eSRichard Tran Mills da->ops->localtolocalend = DMLocalToLocalEnd_DA; 409a4121054SBarry Smith da->ops->createglobalvector = DMCreateGlobalVector_DA; 410a4121054SBarry Smith da->ops->createlocalvector = DMCreateLocalVector_DA; 41125296bd5SBarry Smith da->ops->createinterpolation = DMCreateInterpolation_DA; 412e727c939SJed Brown da->ops->getcoloring = DMCreateColoring_DA; 41325296bd5SBarry Smith da->ops->creatematrix = DMCreateMatrix_DA; 414a4121054SBarry Smith da->ops->refine = DMRefine_DA; 415a4121054SBarry Smith da->ops->coarsen = DMCoarsen_DA; 416a4121054SBarry Smith da->ops->refinehierarchy = DMRefineHierarchy_DA; 417a4121054SBarry Smith da->ops->coarsenhierarchy = DMCoarsenHierarchy_DA; 4185a84ad33SLisandro Dalcin da->ops->createinjection = DMCreateInjection_DA; 4194a7a4c06SLawrence Mitchell da->ops->hascreateinjection = DMHasCreateInjection_DA; 420a4121054SBarry Smith da->ops->destroy = DMDestroy_DA; 421ea78f98cSLisandro Dalcin da->ops->view = NULL; 422a4121054SBarry Smith da->ops->setfromoptions = DMSetFromOptions_DA; 423a4121054SBarry Smith da->ops->setup = DMSetUp_DA; 42438221697SMatthew G. Knepley da->ops->clone = DMClone_DA; 425b859378eSBarry Smith da->ops->load = DMLoad_DA; 4266636e97aSMatthew G Knepley da->ops->createcoordinatedm = DMCreateCoordinateDM_DA; 427d724dfffSBarry Smith da->ops->createsubdm = DMCreateSubDM_DA; 42816621825SDmitry Karpeev da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA; 429e30e807fSPeter Brune da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA; 430e30e807fSPeter Brune da->ops->createddscatters = DMCreateDomainDecompositionScatters_DA; 431793f3fe5SMatthew G. Knepley da->ops->getdimpoints = DMGetDimPoints_DA; 432502a2867SDave May da->ops->getneighbors = DMGetNeighbors_DA; 4336c6a6b79SMatthew G. Knepley da->ops->getlocalboundingbox = DMGetLocalBoundingBox_DA; 43405264a50SDave May da->ops->locatepoints = DMLocatePoints_DA_Regular; 4357c3cd84eSPatrick Sanan da->ops->getcompatibility = DMGetCompatibility_DA; 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)da, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_DMDA)); 4373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 438a4121054SBarry Smith } 43947c6ae99SBarry Smith 440a4121054SBarry Smith /*@ 44112b4a537SBarry Smith DMDACreate - Creates a `DMDA` object for managing structured grids. 442a4121054SBarry Smith 4435edebc4bSPatrick Sanan Collective 444a4121054SBarry Smith 445a4121054SBarry Smith Input Parameter: 44612b4a537SBarry Smith . comm - The communicator for the `DMDA` object 447a4121054SBarry Smith 448a4121054SBarry Smith Output Parameter: 44972fd0fbdSBarry Smith . da - the `DMDA` object 450a4121054SBarry Smith 451e0f5d30fSBarry Smith Level: advanced 452e0f5d30fSBarry Smith 45372fd0fbdSBarry Smith Notes: 45472fd0fbdSBarry Smith See [](sec_struct) for details on the construction of a `DMDA` 455a4121054SBarry Smith 45672fd0fbdSBarry Smith `DMDACreate1d()`, `DMDACreate2d()`, and `DMDACreate3d()` are convenience routines to quickly completely create a `DMDA` 45772fd0fbdSBarry Smith 45872fd0fbdSBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMSetUp()`, `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 459a4121054SBarry Smith @*/ 460d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate(MPI_Comm comm, DM *da) 461d71ae5a4SJacob Faibussowitsch { 462a4121054SBarry Smith PetscFunctionBegin; 4634f572ea9SToby Isaac PetscAssertPointer(da, 2); 4649566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, da)); 4659566063dSJacob Faibussowitsch PetscCall(DMSetType(*da, DMDA)); 4663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46747c6ae99SBarry Smith } 468