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; 7897f7067SBarry Smith PetscBool flg; 847c6ae99SBarry Smith 947c6ae99SBarry Smith PetscFunctionBegin; 1008401ef6SPierre 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"); 1108401ef6SPierre 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"); 1208401ef6SPierre 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"); 13235683edSBarry Smith 14d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMDA Options"); 159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_grid_x", "Number of grid points in x direction", "DMDASetSizes", dd->M, &dd->M, NULL, 1)); 169566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y", "Number of grid points in y direction", "DMDASetSizes", dd->N, &dd->N, NULL, 1)); 179566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z", "Number of grid points in z direction", "DMDASetSizes", dd->P, &dd->P, NULL, 1)); 187ddda789SPeter Brune 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap", "Decomposition overlap in all directions", "DMDASetOverlap", dd->xol, &dd->xol, &flg, 0)); 209566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetOverlap(da, dd->xol, dd->xol, dd->xol)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap_x", "Decomposition overlap in x direction", "DMDASetOverlap", dd->xol, &dd->xol, NULL, 0)); 229566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y", "Decomposition overlap in y direction", "DMDASetOverlap", dd->yol, &dd->yol, NULL, 0)); 239566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z", "Decomposition overlap in z direction", "DMDASetOverlap", dd->zol, &dd->zol, NULL, 0)); 243e7870d2SPeter Brune 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_local_subdomains", "", "DMDASetNumLocalSubdomains", dd->Nsub, &dd->Nsub, &flg, PETSC_DECIDE)); 269566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetNumLocalSubDomains(da, dd->Nsub)); 273e7870d2SPeter Brune 288d07fd27SPatrick Sanan /* Handle DMDA parallel distribution */ 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_processors_x", "Number of processors in x direction", "DMDASetNumProcs", dd->m, &dd->m, NULL, PETSC_DECIDE)); 309566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y", "Number of processors in y direction", "DMDASetNumProcs", dd->n, &dd->n, NULL, PETSC_DECIDE)); 319566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z", "Number of processors in z direction", "DMDASetNumProcs", dd->p, &dd->p, NULL, PETSC_DECIDE)); 32aa219208SBarry Smith /* Handle DMDA refinement */ 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine_x", "Refinement ratio in x direction", "DMDASetRefinementFactor", dd->refine_x, &dd->refine_x, NULL, 1)); 349566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y", "Refinement ratio in y direction", "DMDASetRefinementFactor", dd->refine_y, &dd->refine_y, NULL, 1)); 359566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z", "Refinement ratio in z direction", "DMDASetRefinementFactor", dd->refine_z, &dd->refine_z, NULL, 1)); 369371c9d4SSatish Balay dd->coarsen_x = dd->refine_x; 379371c9d4SSatish Balay dd->coarsen_y = dd->refine_y; 389371c9d4SSatish Balay dd->coarsen_z = dd->refine_z; 3947c6ae99SBarry Smith 40397b6216SJed Brown /* Get refinement factors, defaults taken from the coarse DMDA */ 419566063dSJacob Faibussowitsch PetscCall(DMDAGetRefinementFactor(da, &refx[0], &refy[0], &refz[0])); 42397b6216SJed Brown for (i = 1; i < maxnlevels; i++) { 43397b6216SJed Brown refx[i] = refx[0]; 44397b6216SJed Brown refy[i] = refy[0]; 45397b6216SJed Brown refz[i] = refz[0]; 46397b6216SJed Brown } 47397b6216SJed Brown n = maxnlevels; 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_x", "Refinement factor for each level", "None", refx, &n, &flg)); 49897f7067SBarry Smith if (flg) { 50897f7067SBarry Smith dd->refine_x = refx[0]; 51897f7067SBarry Smith dd->refine_x_hier_n = n; 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_x_hier)); 539566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_x_hier, refx, n)); 54897f7067SBarry Smith } 55c73cfb54SMatthew G. Knepley if (dim > 1) { 56397b6216SJed Brown n = maxnlevels; 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_y", "Refinement factor for each level", "None", refy, &n, &flg)); 58897f7067SBarry Smith if (flg) { 59897f7067SBarry Smith dd->refine_y = refy[0]; 60897f7067SBarry Smith dd->refine_y_hier_n = n; 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_y_hier)); 629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_y_hier, refy, n)); 63897f7067SBarry Smith } 64397b6216SJed Brown } 65c73cfb54SMatthew G. Knepley if (dim > 2) { 66397b6216SJed Brown n = maxnlevels; 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_z", "Refinement factor for each level", "None", refz, &n, &flg)); 68897f7067SBarry Smith if (flg) { 69897f7067SBarry Smith dd->refine_z = refz[0]; 70897f7067SBarry Smith dd->refine_z_hier_n = n; 719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_z_hier)); 729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_z_hier, refz, n)); 73897f7067SBarry Smith } 74397b6216SJed Brown } 75397b6216SJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine", "Uniformly refine DA one or more times", "None", refine, &refine, NULL, 0)); 77d0609cedSBarry Smith PetscOptionsHeadEnd(); 78235683edSBarry Smith 79e0f5d30fSBarry Smith while (refine--) { 80bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 817a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M, &dd->M)); 82e0f5d30fSBarry Smith } else { 837a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M - 1, &dd->M)); 847a3b7fd5SBarry Smith dd->M += 1; 85e0f5d30fSBarry Smith } 86bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 877a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N, &dd->N)); 88e0f5d30fSBarry Smith } else { 897a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N - 1, &dd->N)); 907a3b7fd5SBarry Smith dd->N += 1; 91e0f5d30fSBarry Smith } 92bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 937a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P, &dd->P)); 94e0f5d30fSBarry Smith } else { 957a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P - 1, &dd->P)); 967a3b7fd5SBarry Smith dd->P += 1; 97e0f5d30fSBarry Smith } 98e0f5d30fSBarry Smith da->levelup++; 9985fc4b34SJed Brown if (da->levelup - da->leveldown >= 0) { 10085fc4b34SJed Brown dd->refine_x = refx[da->levelup - da->leveldown]; 10185fc4b34SJed Brown dd->refine_y = refy[da->levelup - da->leveldown]; 10285fc4b34SJed Brown dd->refine_z = refz[da->levelup - da->leveldown]; 10385fc4b34SJed Brown } 10485fc4b34SJed Brown if (da->levelup - da->leveldown >= 1) { 10585fc4b34SJed Brown dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 10685fc4b34SJed Brown dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 10785fc4b34SJed Brown dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 10885fc4b34SJed Brown } 109e0f5d30fSBarry Smith } 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11147c6ae99SBarry Smith } 11247c6ae99SBarry Smith 113ba38deedSJacob Faibussowitsch static PetscErrorCode DMLoad_DA(DM da, PetscViewer viewer) 114d71ae5a4SJacob Faibussowitsch { 115060da220SMatthew G. Knepley PetscInt dim, m, n, p, dof, swidth; 116b859378eSBarry Smith DMDAStencilType stencil; 117bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 118bc2bf880SBarry Smith PetscBool coors; 119bc2bf880SBarry Smith DM dac; 120bc2bf880SBarry Smith Vec c; 121b859378eSBarry Smith 122b859378eSBarry Smith PetscFunctionBegin; 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dim, 1, NULL, PETSC_INT)); 1249566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &m, 1, NULL, PETSC_INT)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &p, 1, NULL, PETSC_INT)); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dof, 1, NULL, PETSC_INT)); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &swidth, 1, NULL, PETSC_INT)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bx, 1, NULL, PETSC_ENUM)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &by, 1, NULL, PETSC_ENUM)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bz, 1, NULL, PETSC_ENUM)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &stencil, 1, NULL, PETSC_ENUM)); 133b859378eSBarry Smith 1349566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1359566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, m, n, p)); 1369566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, bx, by, bz)); 1379566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, dof)); 1389566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, stencil)); 1399566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, swidth)); 1409566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &coors, 1, NULL, PETSC_ENUM)); 142bc2bf880SBarry Smith if (coors) { 1439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &dac)); 1449566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &c)); 1459566063dSJacob Faibussowitsch PetscCall(VecLoad(c, viewer)); 1469566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, c)); 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 148bc2bf880SBarry Smith } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150b859378eSBarry Smith } 151b859378eSBarry Smith 152ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 153d71ae5a4SJacob Faibussowitsch { 154d724dfffSBarry Smith DM_DA *da = (DM_DA *)dm->data; 155d724dfffSBarry Smith 156d724dfffSBarry Smith PetscFunctionBegin; 157d724dfffSBarry Smith if (subdm) { 158c38c1269SMatthew G. Knepley PetscSF sf; 159c38c1269SMatthew G. Knepley Vec coords; 160c38c1269SMatthew G. Knepley void *ctx; 161c38c1269SMatthew G. Knepley /* Cannot use DMClone since the dof stuff is mixed in. Ugh 1629566063dSJacob Faibussowitsch PetscCall(DMClone(dm, subdm)); */ 1639566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm)); 1649566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1659566063dSJacob Faibussowitsch PetscCall(DMSetPointSF(*subdm, sf)); 1669566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 1679566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*subdm, ctx)); 1689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 169c38c1269SMatthew G. Knepley if (coords) { 1709566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*subdm, coords)); 171c38c1269SMatthew G. Knepley } else { 1729566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &coords)); 1739566063dSJacob Faibussowitsch if (coords) PetscCall(DMSetCoordinates(*subdm, coords)); 174c38c1269SMatthew G. Knepley } 175c38c1269SMatthew G. Knepley 1769566063dSJacob Faibussowitsch PetscCall(DMSetType(*subdm, DMDA)); 1779566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*subdm, dm->dim)); 1789566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P)); 1799566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p)); 1809566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz)); 1819566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*subdm, numFields)); 1829566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*subdm, da->stencil_type)); 1839566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*subdm, da->s)); 1849566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz)); 185d724dfffSBarry Smith } 186d724dfffSBarry Smith if (is) { 187d724dfffSBarry Smith PetscInt *indices, cnt = 0, dof = da->w, i, j; 18838221697SMatthew G. Knepley 1899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(da->Nlocal * numFields / dof, &indices)); 19038221697SMatthew G. Knepley for (i = da->base / dof; i < (da->base + da->Nlocal) / dof; ++i) { 191ad540459SPierre Jolivet for (j = 0; j < numFields; ++j) indices[cnt++] = dof * i + fields[j]; 192d724dfffSBarry Smith } 19363a3b9bcSJacob 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); 1949566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), cnt, indices, PETSC_OWN_POINTER, is)); 195d724dfffSBarry Smith } 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197d724dfffSBarry Smith } 198d724dfffSBarry Smith 199ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 200d71ae5a4SJacob Faibussowitsch { 2011c3fb106SBarry Smith PetscInt i; 2021c3fb106SBarry Smith DM_DA *dd = (DM_DA *)dm->data; 2031c3fb106SBarry Smith PetscInt dof = dd->w; 2041c3fb106SBarry Smith 2051c3fb106SBarry Smith PetscFunctionBegin; 206731c8d9eSDmitry Karpeev if (len) *len = dof; 2071c3fb106SBarry Smith if (islist) { 2081c3fb106SBarry Smith Vec v; 2091c3fb106SBarry Smith PetscInt rstart, n; 2101c3fb106SBarry Smith 2119566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(v, &rstart, NULL)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 2149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 2159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, islist)); 21648a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), n / dof, rstart + i, dof, &(*islist)[i])); 2171c3fb106SBarry Smith } 2181c3fb106SBarry Smith if (namelist) { 2199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, namelist)); 2201c3fb106SBarry Smith if (dd->fieldname) { 22148a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(PetscStrallocpy(dd->fieldname[i], &(*namelist)[i])); 2221c3fb106SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently DMDA must have fieldnames"); 2231c3fb106SBarry Smith } 2241c3fb106SBarry Smith if (dmlist) { 2251c3fb106SBarry Smith DM da; 2261c3fb106SBarry Smith 2279566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 2289566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dm->dim)); 2299566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P)); 2309566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p)); 2319566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz)); 2329566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 2339566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, dd->stencil_type)); 2349566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, dd->s)); 2359566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, dmlist)); 2379566063dSJacob Faibussowitsch for (i = 0; i < dof - 1; i++) PetscCall(PetscObjectReference((PetscObject)da)); 2381c3fb106SBarry Smith for (i = 0; i < dof; i++) (*dmlist)[i] = da; 2391c3fb106SBarry Smith } 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2411c3fb106SBarry Smith } 2421c3fb106SBarry Smith 243ba38deedSJacob Faibussowitsch static PetscErrorCode DMClone_DA(DM dm, DM *newdm) 244d71ae5a4SJacob Faibussowitsch { 24538221697SMatthew G. Knepley DM_DA *da = (DM_DA *)dm->data; 24638221697SMatthew G. Knepley 24738221697SMatthew G. Knepley PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(DMSetType(*newdm, DMDA)); 2499566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*newdm, dm->dim)); 2509566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P)); 2519566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p)); 2529566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz)); 2539566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*newdm, da->w)); 2549566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*newdm, da->stencil_type)); 2559566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*newdm, da->s)); 2569566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz)); 2579566063dSJacob Faibussowitsch PetscCall(DMSetUp(*newdm)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25938221697SMatthew G. Knepley } 26038221697SMatthew G. Knepley 261d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg) 262d71ae5a4SJacob Faibussowitsch { 2634a7a4c06SLawrence Mitchell DM_DA *da = (DM_DA *)dm->data; 2644a7a4c06SLawrence Mitchell 2654a7a4c06SLawrence Mitchell PetscFunctionBegin; 2664a7a4c06SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2674f572ea9SToby Isaac PetscAssertPointer(flg, 2); 268131f56ebSLawrence Mitchell *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE; 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2704a7a4c06SLawrence Mitchell } 2714a7a4c06SLawrence Mitchell 272d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 273d71ae5a4SJacob Faibussowitsch { 274793f3fe5SMatthew G. Knepley PetscFunctionBegin; 2759566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd)); 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277793f3fe5SMatthew G. Knepley } 278793f3fe5SMatthew G. Knepley 279d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 280d71ae5a4SJacob Faibussowitsch { 281502a2867SDave May PetscInt dim; 282502a2867SDave May DMDAStencilType st; 283502a2867SDave May 284502a2867SDave May PetscFunctionBegin; 2859566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dm, ranks)); 2869566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &st)); 287502a2867SDave May 288502a2867SDave May switch (dim) { 289502a2867SDave May case 1: 290502a2867SDave May *nranks = 3; 291ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 3; */ 292502a2867SDave May break; 293502a2867SDave May case 2: 294502a2867SDave May *nranks = 9; 295ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 5; */ 296502a2867SDave May break; 297502a2867SDave May case 3: 298502a2867SDave May *nranks = 27; 299ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 7; */ 300502a2867SDave May break; 301d71ae5a4SJacob Faibussowitsch default: 302d71ae5a4SJacob Faibussowitsch break; 303502a2867SDave May } 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 305502a2867SDave May } 306502a2867SDave May 3073efe6655SBarry Smith /*MC 308dce8aebaSBarry Smith DMDA = "da" - A `DM` object that is used to manage data for a structured grid in 1, 2, or 3 dimensions. 3093efe6655SBarry Smith In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points. 31012b4a537SBarry Smith In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width set with `DMDASetStencilWidth()`. 3113efe6655SBarry Smith 3123efe6655SBarry Smith The vectors can be thought of as either cell centered or vertex centered on the mesh. But some variables cannot be cell centered and others 31312b4a537SBarry Smith vertex centered; see the documentation for `DMSTAG`, a similar `DM` implementation which supports more general staggered grids. 3143efe6655SBarry Smith 3153efe6655SBarry Smith Level: intermediate 3163efe6655SBarry Smith 31712b4a537SBarry Smith .seealso: [](sec_struct), `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMDASetStencilWidth()`, `DMDASetStencilType()`, 31812b4a537SBarry Smith `DMDAStencilType` 3193efe6655SBarry Smith M*/ 3203efe6655SBarry Smith 321d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da) 322d71ae5a4SJacob Faibussowitsch { 32347c6ae99SBarry Smith DM_DA *dd; 32447c6ae99SBarry Smith 32547c6ae99SBarry Smith PetscFunctionBegin; 3264f572ea9SToby Isaac PetscAssertPointer(da, 1); 3274dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dd)); 328a4121054SBarry Smith da->data = dd; 32947c6ae99SBarry Smith 330c73cfb54SMatthew G. Knepley da->dim = -1; 331aa219208SBarry Smith dd->interptype = DMDA_Q1; 33247c6ae99SBarry Smith dd->refine_x = 2; 33347c6ae99SBarry Smith dd->refine_y = 2; 33447c6ae99SBarry Smith dd->refine_z = 2; 33581c108dcSJed Brown dd->coarsen_x = 2; 33681c108dcSJed Brown dd->coarsen_y = 2; 33781c108dcSJed Brown dd->coarsen_z = 2; 3380298fd71SBarry Smith dd->fieldname = NULL; 33947c6ae99SBarry Smith dd->nlocal = -1; 34047c6ae99SBarry Smith dd->Nlocal = -1; 34147c6ae99SBarry Smith dd->M = -1; 34247c6ae99SBarry Smith dd->N = -1; 34347c6ae99SBarry Smith dd->P = -1; 34447c6ae99SBarry Smith dd->m = -1; 34547c6ae99SBarry Smith dd->n = -1; 34647c6ae99SBarry Smith dd->p = -1; 34747c6ae99SBarry Smith dd->w = -1; 34847c6ae99SBarry Smith dd->s = -1; 3498865f1eaSKarl Rupp 3509371c9d4SSatish Balay dd->xs = -1; 3519371c9d4SSatish Balay dd->xe = -1; 3529371c9d4SSatish Balay dd->ys = -1; 3539371c9d4SSatish Balay dd->ye = -1; 3549371c9d4SSatish Balay dd->zs = -1; 3559371c9d4SSatish Balay dd->ze = -1; 3569371c9d4SSatish Balay dd->Xs = -1; 3579371c9d4SSatish Balay dd->Xe = -1; 3589371c9d4SSatish Balay dd->Ys = -1; 3599371c9d4SSatish Balay dd->Ye = -1; 3609371c9d4SSatish Balay dd->Zs = -1; 3619371c9d4SSatish Balay dd->Ze = -1; 36247c6ae99SBarry Smith 3633e7870d2SPeter Brune dd->Nsub = 1; 3647ddda789SPeter Brune dd->xol = 0; 3657ddda789SPeter Brune dd->yol = 0; 3667ddda789SPeter Brune dd->zol = 0; 367d886c4f4SPeter Brune dd->xo = 0; 368d886c4f4SPeter Brune dd->yo = 0; 369d886c4f4SPeter Brune dd->zo = 0; 370e30e807fSPeter Brune dd->Mo = -1; 371e30e807fSPeter Brune dd->No = -1; 372e30e807fSPeter Brune dd->Po = -1; 37388661749SPeter Brune 3740298fd71SBarry Smith dd->gtol = NULL; 3750298fd71SBarry Smith dd->ltol = NULL; 3760298fd71SBarry Smith dd->ao = NULL; 3773ba16761SJacob Faibussowitsch PetscCall(PetscStrallocpy(AOBASIC, (char **)&dd->aotype)); 37847c6ae99SBarry Smith dd->base = -1; 379bff4a2f0SMatthew G. Knepley dd->bx = DM_BOUNDARY_NONE; 380bff4a2f0SMatthew G. Knepley dd->by = DM_BOUNDARY_NONE; 381bff4a2f0SMatthew G. Knepley dd->bz = DM_BOUNDARY_NONE; 382aa219208SBarry Smith dd->stencil_type = DMDA_STENCIL_BOX; 383aa219208SBarry Smith dd->interptype = DMDA_Q1; 3840298fd71SBarry Smith dd->lx = NULL; 3850298fd71SBarry Smith dd->ly = NULL; 3860298fd71SBarry Smith dd->lz = NULL; 38747c6ae99SBarry Smith 388454e267fSLisandro Dalcin dd->elementtype = DMDA_ELEMENT_Q1; 389454e267fSLisandro Dalcin 390a4121054SBarry Smith da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA; 391a4121054SBarry Smith da->ops->globaltolocalend = DMGlobalToLocalEnd_DA; 392a4121054SBarry Smith da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA; 393a4121054SBarry Smith da->ops->localtoglobalend = DMLocalToGlobalEnd_DA; 394d78e899eSRichard Tran Mills da->ops->localtolocalbegin = DMLocalToLocalBegin_DA; 395d78e899eSRichard Tran Mills da->ops->localtolocalend = DMLocalToLocalEnd_DA; 396a4121054SBarry Smith da->ops->createglobalvector = DMCreateGlobalVector_DA; 397a4121054SBarry Smith da->ops->createlocalvector = DMCreateLocalVector_DA; 39825296bd5SBarry Smith da->ops->createinterpolation = DMCreateInterpolation_DA; 399e727c939SJed Brown da->ops->getcoloring = DMCreateColoring_DA; 40025296bd5SBarry Smith da->ops->creatematrix = DMCreateMatrix_DA; 401a4121054SBarry Smith da->ops->refine = DMRefine_DA; 402a4121054SBarry Smith da->ops->coarsen = DMCoarsen_DA; 403a4121054SBarry Smith da->ops->refinehierarchy = DMRefineHierarchy_DA; 404a4121054SBarry Smith da->ops->coarsenhierarchy = DMCoarsenHierarchy_DA; 4055a84ad33SLisandro Dalcin da->ops->createinjection = DMCreateInjection_DA; 4064a7a4c06SLawrence Mitchell da->ops->hascreateinjection = DMHasCreateInjection_DA; 407a4121054SBarry Smith da->ops->destroy = DMDestroy_DA; 408ea78f98cSLisandro Dalcin da->ops->view = NULL; 409a4121054SBarry Smith da->ops->setfromoptions = DMSetFromOptions_DA; 410a4121054SBarry Smith da->ops->setup = DMSetUp_DA; 41138221697SMatthew G. Knepley da->ops->clone = DMClone_DA; 412b859378eSBarry Smith da->ops->load = DMLoad_DA; 4136636e97aSMatthew G Knepley da->ops->createcoordinatedm = DMCreateCoordinateDM_DA; 414d724dfffSBarry Smith da->ops->createsubdm = DMCreateSubDM_DA; 41516621825SDmitry Karpeev da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA; 416e30e807fSPeter Brune da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA; 417e30e807fSPeter Brune da->ops->createddscatters = DMCreateDomainDecompositionScatters_DA; 418793f3fe5SMatthew G. Knepley da->ops->getdimpoints = DMGetDimPoints_DA; 419502a2867SDave May da->ops->getneighbors = DMGetNeighbors_DA; 420*6c6a6b79SMatthew G. Knepley da->ops->getlocalboundingbox = DMGetLocalBoundingBox_DA; 42105264a50SDave May da->ops->locatepoints = DMLocatePoints_DA_Regular; 4227c3cd84eSPatrick Sanan da->ops->getcompatibility = DMGetCompatibility_DA; 4239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)da, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_DMDA)); 4243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 425a4121054SBarry Smith } 42647c6ae99SBarry Smith 427a4121054SBarry Smith /*@ 42812b4a537SBarry Smith DMDACreate - Creates a `DMDA` object for managing structured grids. 429a4121054SBarry Smith 4305edebc4bSPatrick Sanan Collective 431a4121054SBarry Smith 432a4121054SBarry Smith Input Parameter: 43312b4a537SBarry Smith . comm - The communicator for the `DMDA` object 434a4121054SBarry Smith 435a4121054SBarry Smith Output Parameter: 436a4121054SBarry Smith . da - The DMDA object 437a4121054SBarry Smith 438e0f5d30fSBarry Smith Level: advanced 439e0f5d30fSBarry Smith 44012b4a537SBarry Smith Developer Note: 4415edebc4bSPatrick Sanan Since there exists DMDACreate1/2/3d() should this routine even exist? 442a4121054SBarry Smith 44312b4a537SBarry Smith .seealso: [](sec_struct), `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 444a4121054SBarry Smith @*/ 445d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate(MPI_Comm comm, DM *da) 446d71ae5a4SJacob Faibussowitsch { 447a4121054SBarry Smith PetscFunctionBegin; 4484f572ea9SToby Isaac PetscAssertPointer(da, 2); 4499566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, da)); 4509566063dSJacob Faibussowitsch PetscCall(DMSetType(*da, DMDA)); 4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45247c6ae99SBarry Smith } 453