17d0a6c19SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 347c6ae99SBarry Smith 4*ba38deedSJacob Faibussowitsch static PetscErrorCode DMSetFromOptions_DA(DM da, PetscOptionItems *PetscOptionsObject) 5d71ae5a4SJacob Faibussowitsch { 647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 7c73cfb54SMatthew G. Knepley PetscInt refine = 0, dim = da->dim, maxnlevels = 100, refx[100], refy[100], refz[100], n, i; 8897f7067SBarry Smith PetscBool flg; 947c6ae99SBarry Smith 1047c6ae99SBarry Smith PetscFunctionBegin; 1108401ef6SPierre 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"); 1208401ef6SPierre 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"); 1308401ef6SPierre 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"); 14235683edSBarry Smith 15d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMDA Options"); 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_grid_x", "Number of grid points in x direction", "DMDASetSizes", dd->M, &dd->M, NULL, 1)); 179566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y", "Number of grid points in y direction", "DMDASetSizes", dd->N, &dd->N, NULL, 1)); 189566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z", "Number of grid points in z direction", "DMDASetSizes", dd->P, &dd->P, NULL, 1)); 197ddda789SPeter Brune 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap", "Decomposition overlap in all directions", "DMDASetOverlap", dd->xol, &dd->xol, &flg, 0)); 219566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetOverlap(da, dd->xol, dd->xol, dd->xol)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap_x", "Decomposition overlap in x direction", "DMDASetOverlap", dd->xol, &dd->xol, NULL, 0)); 239566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y", "Decomposition overlap in y direction", "DMDASetOverlap", dd->yol, &dd->yol, NULL, 0)); 249566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z", "Decomposition overlap in z direction", "DMDASetOverlap", dd->zol, &dd->zol, NULL, 0)); 253e7870d2SPeter Brune 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_local_subdomains", "", "DMDASetNumLocalSubdomains", dd->Nsub, &dd->Nsub, &flg, PETSC_DECIDE)); 279566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetNumLocalSubDomains(da, dd->Nsub)); 283e7870d2SPeter Brune 298d07fd27SPatrick Sanan /* Handle DMDA parallel distribution */ 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_processors_x", "Number of processors in x direction", "DMDASetNumProcs", dd->m, &dd->m, NULL, PETSC_DECIDE)); 319566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y", "Number of processors in y direction", "DMDASetNumProcs", dd->n, &dd->n, NULL, PETSC_DECIDE)); 329566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z", "Number of processors in z direction", "DMDASetNumProcs", dd->p, &dd->p, NULL, PETSC_DECIDE)); 33aa219208SBarry Smith /* Handle DMDA refinement */ 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine_x", "Refinement ratio in x direction", "DMDASetRefinementFactor", dd->refine_x, &dd->refine_x, NULL, 1)); 359566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y", "Refinement ratio in y direction", "DMDASetRefinementFactor", dd->refine_y, &dd->refine_y, NULL, 1)); 369566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z", "Refinement ratio in z direction", "DMDASetRefinementFactor", dd->refine_z, &dd->refine_z, NULL, 1)); 379371c9d4SSatish Balay dd->coarsen_x = dd->refine_x; 389371c9d4SSatish Balay dd->coarsen_y = dd->refine_y; 399371c9d4SSatish Balay dd->coarsen_z = dd->refine_z; 4047c6ae99SBarry Smith 41397b6216SJed Brown /* Get refinement factors, defaults taken from the coarse DMDA */ 429566063dSJacob Faibussowitsch PetscCall(DMDAGetRefinementFactor(da, &refx[0], &refy[0], &refz[0])); 43397b6216SJed Brown for (i = 1; i < maxnlevels; i++) { 44397b6216SJed Brown refx[i] = refx[0]; 45397b6216SJed Brown refy[i] = refy[0]; 46397b6216SJed Brown refz[i] = refz[0]; 47397b6216SJed Brown } 48397b6216SJed Brown n = maxnlevels; 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_x", "Refinement factor for each level", "None", refx, &n, &flg)); 50897f7067SBarry Smith if (flg) { 51897f7067SBarry Smith dd->refine_x = refx[0]; 52897f7067SBarry Smith dd->refine_x_hier_n = n; 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_x_hier)); 549566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_x_hier, refx, n)); 55897f7067SBarry Smith } 56c73cfb54SMatthew G. Knepley if (dim > 1) { 57397b6216SJed Brown n = maxnlevels; 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_y", "Refinement factor for each level", "None", refy, &n, &flg)); 59897f7067SBarry Smith if (flg) { 60897f7067SBarry Smith dd->refine_y = refy[0]; 61897f7067SBarry Smith dd->refine_y_hier_n = n; 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_y_hier)); 639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_y_hier, refy, n)); 64897f7067SBarry Smith } 65397b6216SJed Brown } 66c73cfb54SMatthew G. Knepley if (dim > 2) { 67397b6216SJed Brown n = maxnlevels; 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_z", "Refinement factor for each level", "None", refz, &n, &flg)); 69897f7067SBarry Smith if (flg) { 70897f7067SBarry Smith dd->refine_z = refz[0]; 71897f7067SBarry Smith dd->refine_z_hier_n = n; 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->refine_z_hier)); 739566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dd->refine_z_hier, refz, n)); 74897f7067SBarry Smith } 75397b6216SJed Brown } 76397b6216SJed Brown 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine", "Uniformly refine DA one or more times", "None", refine, &refine, NULL, 0)); 78d0609cedSBarry Smith PetscOptionsHeadEnd(); 79235683edSBarry Smith 80e0f5d30fSBarry Smith while (refine--) { 81bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 827a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M, &dd->M)); 83e0f5d30fSBarry Smith } else { 847a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_x, dd->M - 1, &dd->M)); 857a3b7fd5SBarry Smith dd->M += 1; 86e0f5d30fSBarry Smith } 87bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 887a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N, &dd->N)); 89e0f5d30fSBarry Smith } else { 907a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_y, dd->N - 1, &dd->N)); 917a3b7fd5SBarry Smith dd->N += 1; 92e0f5d30fSBarry Smith } 93bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 947a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P, &dd->P)); 95e0f5d30fSBarry Smith } else { 967a3b7fd5SBarry Smith PetscCall(PetscIntMultError(dd->refine_z, dd->P - 1, &dd->P)); 977a3b7fd5SBarry Smith dd->P += 1; 98e0f5d30fSBarry Smith } 99e0f5d30fSBarry Smith da->levelup++; 10085fc4b34SJed Brown if (da->levelup - da->leveldown >= 0) { 10185fc4b34SJed Brown dd->refine_x = refx[da->levelup - da->leveldown]; 10285fc4b34SJed Brown dd->refine_y = refy[da->levelup - da->leveldown]; 10385fc4b34SJed Brown dd->refine_z = refz[da->levelup - da->leveldown]; 10485fc4b34SJed Brown } 10585fc4b34SJed Brown if (da->levelup - da->leveldown >= 1) { 10685fc4b34SJed Brown dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 10785fc4b34SJed Brown dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 10885fc4b34SJed Brown dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 10985fc4b34SJed Brown } 110e0f5d30fSBarry Smith } 1113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11247c6ae99SBarry Smith } 11347c6ae99SBarry Smith 114*ba38deedSJacob Faibussowitsch static PetscErrorCode DMLoad_DA(DM da, PetscViewer viewer) 115d71ae5a4SJacob Faibussowitsch { 116060da220SMatthew G. Knepley PetscInt dim, m, n, p, dof, swidth; 117b859378eSBarry Smith DMDAStencilType stencil; 118bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz; 119bc2bf880SBarry Smith PetscBool coors; 120bc2bf880SBarry Smith DM dac; 121bc2bf880SBarry Smith Vec c; 122b859378eSBarry Smith 123b859378eSBarry Smith PetscFunctionBegin; 1249566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dim, 1, NULL, PETSC_INT)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &m, 1, NULL, PETSC_INT)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT)); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &p, 1, NULL, PETSC_INT)); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &dof, 1, NULL, PETSC_INT)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &swidth, 1, NULL, PETSC_INT)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bx, 1, NULL, PETSC_ENUM)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &by, 1, NULL, PETSC_ENUM)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &bz, 1, NULL, PETSC_ENUM)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &stencil, 1, NULL, PETSC_ENUM)); 134b859378eSBarry Smith 1359566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1369566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, m, n, p)); 1379566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, bx, by, bz)); 1389566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, dof)); 1399566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, stencil)); 1409566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, swidth)); 1419566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &coors, 1, NULL, PETSC_ENUM)); 143bc2bf880SBarry Smith if (coors) { 1449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &dac)); 1459566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &c)); 1469566063dSJacob Faibussowitsch PetscCall(VecLoad(c, viewer)); 1479566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, c)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 149bc2bf880SBarry Smith } 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151b859378eSBarry Smith } 152b859378eSBarry Smith 153*ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 154d71ae5a4SJacob Faibussowitsch { 155d724dfffSBarry Smith DM_DA *da = (DM_DA *)dm->data; 156d724dfffSBarry Smith 157d724dfffSBarry Smith PetscFunctionBegin; 158d724dfffSBarry Smith if (subdm) { 159c38c1269SMatthew G. Knepley PetscSF sf; 160c38c1269SMatthew G. Knepley Vec coords; 161c38c1269SMatthew G. Knepley void *ctx; 162c38c1269SMatthew G. Knepley /* Cannot use DMClone since the dof stuff is mixed in. Ugh 1639566063dSJacob Faibussowitsch PetscCall(DMClone(dm, subdm)); */ 1649566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm)); 1659566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1669566063dSJacob Faibussowitsch PetscCall(DMSetPointSF(*subdm, sf)); 1679566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 1689566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*subdm, ctx)); 1699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 170c38c1269SMatthew G. Knepley if (coords) { 1719566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*subdm, coords)); 172c38c1269SMatthew G. Knepley } else { 1739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &coords)); 1749566063dSJacob Faibussowitsch if (coords) PetscCall(DMSetCoordinates(*subdm, coords)); 175c38c1269SMatthew G. Knepley } 176c38c1269SMatthew G. Knepley 1779566063dSJacob Faibussowitsch PetscCall(DMSetType(*subdm, DMDA)); 1789566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*subdm, dm->dim)); 1799566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P)); 1809566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p)); 1819566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz)); 1829566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*subdm, numFields)); 1839566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*subdm, da->stencil_type)); 1849566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*subdm, da->s)); 1859566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz)); 186d724dfffSBarry Smith } 187d724dfffSBarry Smith if (is) { 188d724dfffSBarry Smith PetscInt *indices, cnt = 0, dof = da->w, i, j; 18938221697SMatthew G. Knepley 1909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(da->Nlocal * numFields / dof, &indices)); 19138221697SMatthew G. Knepley for (i = da->base / dof; i < (da->base + da->Nlocal) / dof; ++i) { 192ad540459SPierre Jolivet for (j = 0; j < numFields; ++j) indices[cnt++] = dof * i + fields[j]; 193d724dfffSBarry Smith } 19463a3b9bcSJacob 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); 1959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), cnt, indices, PETSC_OWN_POINTER, is)); 196d724dfffSBarry Smith } 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 198d724dfffSBarry Smith } 199d724dfffSBarry Smith 200*ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 201d71ae5a4SJacob Faibussowitsch { 2021c3fb106SBarry Smith PetscInt i; 2031c3fb106SBarry Smith DM_DA *dd = (DM_DA *)dm->data; 2041c3fb106SBarry Smith PetscInt dof = dd->w; 2051c3fb106SBarry Smith 2061c3fb106SBarry Smith PetscFunctionBegin; 207731c8d9eSDmitry Karpeev if (len) *len = dof; 2081c3fb106SBarry Smith if (islist) { 2091c3fb106SBarry Smith Vec v; 2101c3fb106SBarry Smith PetscInt rstart, n; 2111c3fb106SBarry Smith 2129566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(v, &rstart, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 2159566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 2169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, islist)); 21748a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), n / dof, rstart + i, dof, &(*islist)[i])); 2181c3fb106SBarry Smith } 2191c3fb106SBarry Smith if (namelist) { 2209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, namelist)); 2211c3fb106SBarry Smith if (dd->fieldname) { 22248a46eb9SPierre Jolivet for (i = 0; i < dof; i++) PetscCall(PetscStrallocpy(dd->fieldname[i], &(*namelist)[i])); 2231c3fb106SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently DMDA must have fieldnames"); 2241c3fb106SBarry Smith } 2251c3fb106SBarry Smith if (dmlist) { 2261c3fb106SBarry Smith DM da; 2271c3fb106SBarry Smith 2289566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 2299566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dm->dim)); 2309566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P)); 2319566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p)); 2329566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz)); 2339566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 2349566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, dd->stencil_type)); 2359566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, dd->s)); 2369566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 2379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, dmlist)); 2389566063dSJacob Faibussowitsch for (i = 0; i < dof - 1; i++) PetscCall(PetscObjectReference((PetscObject)da)); 2391c3fb106SBarry Smith for (i = 0; i < dof; i++) (*dmlist)[i] = da; 2401c3fb106SBarry Smith } 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2421c3fb106SBarry Smith } 2431c3fb106SBarry Smith 244*ba38deedSJacob Faibussowitsch static PetscErrorCode DMClone_DA(DM dm, DM *newdm) 245d71ae5a4SJacob Faibussowitsch { 24638221697SMatthew G. Knepley DM_DA *da = (DM_DA *)dm->data; 24738221697SMatthew G. Knepley 24838221697SMatthew G. Knepley PetscFunctionBegin; 2499566063dSJacob Faibussowitsch PetscCall(DMSetType(*newdm, DMDA)); 2509566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*newdm, dm->dim)); 2519566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P)); 2529566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p)); 2539566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz)); 2549566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*newdm, da->w)); 2559566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*newdm, da->stencil_type)); 2569566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*newdm, da->s)); 2579566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz)); 2589566063dSJacob Faibussowitsch PetscCall(DMSetUp(*newdm)); 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26038221697SMatthew G. Knepley } 26138221697SMatthew G. Knepley 262d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg) 263d71ae5a4SJacob Faibussowitsch { 2644a7a4c06SLawrence Mitchell DM_DA *da = (DM_DA *)dm->data; 2654a7a4c06SLawrence Mitchell 2664a7a4c06SLawrence Mitchell PetscFunctionBegin; 2674a7a4c06SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2684f572ea9SToby Isaac PetscAssertPointer(flg, 2); 269131f56ebSLawrence Mitchell *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE; 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2714a7a4c06SLawrence Mitchell } 2724a7a4c06SLawrence Mitchell 273d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 274d71ae5a4SJacob Faibussowitsch { 275793f3fe5SMatthew G. Knepley PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd)); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278793f3fe5SMatthew G. Knepley } 279793f3fe5SMatthew G. Knepley 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 281d71ae5a4SJacob Faibussowitsch { 282502a2867SDave May PetscInt dim; 283502a2867SDave May DMDAStencilType st; 284502a2867SDave May 285502a2867SDave May PetscFunctionBegin; 2869566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dm, ranks)); 2879566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &st)); 288502a2867SDave May 289502a2867SDave May switch (dim) { 290502a2867SDave May case 1: 291502a2867SDave May *nranks = 3; 292ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 3; */ 293502a2867SDave May break; 294502a2867SDave May case 2: 295502a2867SDave May *nranks = 9; 296ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 5; */ 297502a2867SDave May break; 298502a2867SDave May case 3: 299502a2867SDave May *nranks = 27; 300ad540459SPierre Jolivet /* if (st == DMDA_STENCIL_STAR) *nranks = 7; */ 301502a2867SDave May break; 302d71ae5a4SJacob Faibussowitsch default: 303d71ae5a4SJacob Faibussowitsch break; 304502a2867SDave May } 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306502a2867SDave May } 307502a2867SDave May 3083efe6655SBarry Smith /*MC 309dce8aebaSBarry Smith DMDA = "da" - A `DM` object that is used to manage data for a structured grid in 1, 2, or 3 dimensions. 3103efe6655SBarry Smith In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points. 3113efe6655SBarry Smith In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width. 3123efe6655SBarry Smith 3133efe6655SBarry 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 314dce8aebaSBarry Smith vertex centered; see the documentation for `DMSTAG`, a similar DM implementation which supports these staggered grids. 3153efe6655SBarry Smith 3163efe6655SBarry Smith Level: intermediate 3173efe6655SBarry Smith 318db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()` 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; 42005264a50SDave May da->ops->locatepoints = DMLocatePoints_DA_Regular; 4217c3cd84eSPatrick Sanan da->ops->getcompatibility = DMGetCompatibility_DA; 4229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)da, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_DMDA)); 4233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 424a4121054SBarry Smith } 42547c6ae99SBarry Smith 426a4121054SBarry Smith /*@ 427a4121054SBarry Smith DMDACreate - Creates a DMDA object. 428a4121054SBarry Smith 4295edebc4bSPatrick Sanan Collective 430a4121054SBarry Smith 431a4121054SBarry Smith Input Parameter: 432a4121054SBarry Smith . comm - The communicator for the DMDA object 433a4121054SBarry Smith 434a4121054SBarry Smith Output Parameter: 435a4121054SBarry Smith . da - The DMDA object 436a4121054SBarry Smith 437e0f5d30fSBarry Smith Level: advanced 438e0f5d30fSBarry Smith 43960225df5SJacob Faibussowitsch Developer Notes: 4405edebc4bSPatrick Sanan Since there exists DMDACreate1/2/3d() should this routine even exist? 441a4121054SBarry Smith 442db781477SPatrick Sanan .seealso: `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 443a4121054SBarry Smith @*/ 444d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate(MPI_Comm comm, DM *da) 445d71ae5a4SJacob Faibussowitsch { 446a4121054SBarry Smith PetscFunctionBegin; 4474f572ea9SToby Isaac PetscAssertPointer(da, 2); 4489566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, da)); 4499566063dSJacob Faibussowitsch PetscCall(DMSetType(*da, DMDA)); 4503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45147c6ae99SBarry Smith } 452