17d0a6c19SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 347c6ae99SBarry Smith 44416b707SBarry Smith PetscErrorCode DMSetFromOptions_DA(PetscOptionItems *PetscOptionsObject,DM da) 547c6ae99SBarry Smith { 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; 11064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(da,DM_CLASSID,2); 1247c6ae99SBarry Smith 13*08401ef6SPierre 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"); 14*08401ef6SPierre 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"); 15*08401ef6SPierre 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"); 16235683edSBarry Smith 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"DMDA Options")); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,NULL,1)); 199566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL,1)); 209566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL,1)); 217ddda789SPeter Brune 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg,0)); 239566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetOverlap(da,dd->xol,dd->xol,dd->xol)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL,0)); 259566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL,0)); 269566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL,0)); 273e7870d2SPeter Brune 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg,PETSC_DECIDE)); 299566063dSJacob Faibussowitsch if (flg) PetscCall(DMDASetNumLocalSubDomains(da,dd->Nsub)); 303e7870d2SPeter Brune 318d07fd27SPatrick Sanan /* Handle DMDA parallel distribution */ 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL,PETSC_DECIDE)); 339566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL,PETSC_DECIDE)); 349566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL,PETSC_DECIDE)); 35aa219208SBarry Smith /* Handle DMDA refinement */ 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL,1)); 379566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL,1)); 389566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL,1)); 3938fb6c91SJed Brown dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; 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)); 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 79235683edSBarry Smith 80e0f5d30fSBarry Smith while (refine--) { 81bff4a2f0SMatthew G. Knepley if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 82e0f5d30fSBarry Smith dd->M = dd->refine_x*dd->M; 83e0f5d30fSBarry Smith } else { 84e0f5d30fSBarry Smith dd->M = 1 + dd->refine_x*(dd->M - 1); 85e0f5d30fSBarry Smith } 86bff4a2f0SMatthew G. Knepley if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 87e0f5d30fSBarry Smith dd->N = dd->refine_y*dd->N; 88e0f5d30fSBarry Smith } else { 89e0f5d30fSBarry Smith dd->N = 1 + dd->refine_y*(dd->N - 1); 90e0f5d30fSBarry Smith } 91bff4a2f0SMatthew G. Knepley if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 92e0f5d30fSBarry Smith dd->P = dd->refine_z*dd->P; 93e0f5d30fSBarry Smith } else { 94e0f5d30fSBarry Smith dd->P = 1 + dd->refine_z*(dd->P - 1); 95e0f5d30fSBarry Smith } 96e0f5d30fSBarry Smith da->levelup++; 9785fc4b34SJed Brown if (da->levelup - da->leveldown >= 0) { 9885fc4b34SJed Brown dd->refine_x = refx[da->levelup - da->leveldown]; 9985fc4b34SJed Brown dd->refine_y = refy[da->levelup - da->leveldown]; 10085fc4b34SJed Brown dd->refine_z = refz[da->levelup - da->leveldown]; 10185fc4b34SJed Brown } 10285fc4b34SJed Brown if (da->levelup - da->leveldown >= 1) { 10385fc4b34SJed Brown dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 10485fc4b34SJed Brown dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 10585fc4b34SJed Brown dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 10685fc4b34SJed Brown } 107e0f5d30fSBarry Smith } 10847c6ae99SBarry Smith PetscFunctionReturn(0); 10947c6ae99SBarry Smith } 11047c6ae99SBarry Smith 1117087cfbeSBarry Smith extern PetscErrorCode DMCreateGlobalVector_DA(DM,Vec*); 1127087cfbeSBarry Smith extern PetscErrorCode DMCreateLocalVector_DA(DM,Vec*); 1137087cfbeSBarry Smith extern PetscErrorCode DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec); 1147087cfbeSBarry Smith extern PetscErrorCode DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec); 1157087cfbeSBarry Smith extern PetscErrorCode DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec); 1167087cfbeSBarry Smith extern PetscErrorCode DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec); 117d78e899eSRichard Tran Mills extern PetscErrorCode DMLocalToLocalBegin_DA(DM,Vec,InsertMode,Vec); 118d78e899eSRichard Tran Mills extern PetscErrorCode DMLocalToLocalEnd_DA(DM,Vec,InsertMode,Vec); 119e727c939SJed Brown extern PetscErrorCode DMCreateInterpolation_DA(DM,DM,Mat*,Vec*); 120b412c318SBarry Smith extern PetscErrorCode DMCreateColoring_DA(DM,ISColoringType,ISColoring*); 121b412c318SBarry Smith extern PetscErrorCode DMCreateMatrix_DA(DM,Mat*); 1226636e97aSMatthew G Knepley extern PetscErrorCode DMCreateCoordinateDM_DA(DM,DM*); 1237087cfbeSBarry Smith extern PetscErrorCode DMRefine_DA(DM,MPI_Comm,DM*); 1247087cfbeSBarry Smith extern PetscErrorCode DMCoarsen_DA(DM,MPI_Comm,DM*); 1257087cfbeSBarry Smith extern PetscErrorCode DMRefineHierarchy_DA(DM,PetscInt,DM[]); 1267087cfbeSBarry Smith extern PetscErrorCode DMCoarsenHierarchy_DA(DM,PetscInt,DM[]); 1276dbf9973SLawrence Mitchell extern PetscErrorCode DMCreateInjection_DA(DM,DM,Mat*); 1287087cfbeSBarry Smith extern PetscErrorCode DMView_DA(DM,PetscViewer); 1297087cfbeSBarry Smith extern PetscErrorCode DMSetUp_DA(DM); 1307087cfbeSBarry Smith extern PetscErrorCode DMDestroy_DA(DM); 131e30e807fSPeter Brune extern PetscErrorCode DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**); 132e30e807fSPeter Brune extern PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**); 1337c3cd84eSPatrick Sanan PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM,DM,PetscBool*,PetscBool*); 1349a42bb27SBarry Smith 135b859378eSBarry Smith PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer) 136b859378eSBarry Smith { 137060da220SMatthew G. Knepley PetscInt dim,m,n,p,dof,swidth; 138b859378eSBarry Smith DMDAStencilType stencil; 139bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 140bc2bf880SBarry Smith PetscBool coors; 141bc2bf880SBarry Smith DM dac; 142bc2bf880SBarry Smith Vec c; 143b859378eSBarry Smith 144b859378eSBarry Smith PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&dim,1,NULL,PETSC_INT)); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&m,1,NULL,PETSC_INT)); 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&p,1,NULL,PETSC_INT)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&dof,1,NULL,PETSC_INT)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&swidth,1,NULL,PETSC_INT)); 1519566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&bx,1,NULL,PETSC_ENUM)); 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&by,1,NULL,PETSC_ENUM)); 1539566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&bz,1,NULL,PETSC_ENUM)); 1549566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&stencil,1,NULL,PETSC_ENUM)); 155b859378eSBarry Smith 1569566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1579566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, m,n,p)); 1589566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, bx, by, bz)); 1599566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, dof)); 1609566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, stencil)); 1619566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, swidth)); 1629566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1639566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer,&coors,1,NULL,PETSC_ENUM)); 164bc2bf880SBarry Smith if (coors) { 1659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&dac)); 1669566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac,&c)); 1679566063dSJacob Faibussowitsch PetscCall(VecLoad(c,viewer)); 1689566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da,c)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 170bc2bf880SBarry Smith } 171b859378eSBarry Smith PetscFunctionReturn(0); 172b859378eSBarry Smith } 173b859378eSBarry Smith 174276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 175d724dfffSBarry Smith { 176d724dfffSBarry Smith DM_DA *da = (DM_DA*) dm->data; 177d724dfffSBarry Smith 178d724dfffSBarry Smith PetscFunctionBegin; 179d724dfffSBarry Smith if (subdm) { 180c38c1269SMatthew G. Knepley PetscSF sf; 181c38c1269SMatthew G. Knepley Vec coords; 182c38c1269SMatthew G. Knepley void *ctx; 183c38c1269SMatthew G. Knepley /* Cannot use DMClone since the dof stuff is mixed in. Ugh 1849566063dSJacob Faibussowitsch PetscCall(DMClone(dm, subdm)); */ 1859566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm)); 1869566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1879566063dSJacob Faibussowitsch PetscCall(DMSetPointSF(*subdm, sf)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 1899566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*subdm, ctx)); 1909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 191c38c1269SMatthew G. Knepley if (coords) { 1929566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*subdm, coords)); 193c38c1269SMatthew G. Knepley } else { 1949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &coords)); 1959566063dSJacob Faibussowitsch if (coords) PetscCall(DMSetCoordinates(*subdm, coords)); 196c38c1269SMatthew G. Knepley } 197c38c1269SMatthew G. Knepley 1989566063dSJacob Faibussowitsch PetscCall(DMSetType(*subdm, DMDA)); 1999566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*subdm, dm->dim)); 2009566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P)); 2019566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p)); 2029566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz)); 2039566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*subdm, numFields)); 2049566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*subdm, da->stencil_type)); 2059566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*subdm, da->s)); 2069566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz)); 207d724dfffSBarry Smith } 208d724dfffSBarry Smith if (is) { 209d724dfffSBarry Smith PetscInt *indices, cnt = 0, dof = da->w, i, j; 21038221697SMatthew G. Knepley 2119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(da->Nlocal*numFields/dof, &indices)); 21238221697SMatthew G. Knepley for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) { 21338221697SMatthew G. Knepley for (j = 0; j < numFields; ++j) { 214d724dfffSBarry Smith indices[cnt++] = dof*i + fields[j]; 215d724dfffSBarry Smith } 216d724dfffSBarry Smith } 217*08401ef6SPierre Jolivet PetscCheck(cnt == da->Nlocal*numFields/dof,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %D does not equal expected value %D", cnt, da->Nlocal*numFields/dof); 2189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is)); 219d724dfffSBarry Smith } 220d724dfffSBarry Smith PetscFunctionReturn(0); 221d724dfffSBarry Smith } 222d724dfffSBarry Smith 22316621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist) 2241c3fb106SBarry Smith { 2251c3fb106SBarry Smith PetscInt i; 2261c3fb106SBarry Smith DM_DA *dd = (DM_DA*)dm->data; 2271c3fb106SBarry Smith PetscInt dof = dd->w; 2281c3fb106SBarry Smith 2291c3fb106SBarry Smith PetscFunctionBegin; 230731c8d9eSDmitry Karpeev if (len) *len = dof; 2311c3fb106SBarry Smith if (islist) { 2321c3fb106SBarry Smith Vec v; 2331c3fb106SBarry Smith PetscInt rstart,n; 2341c3fb106SBarry Smith 2359566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm,&v)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(v,&rstart,NULL)); 2379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v,&n)); 2389566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm,&v)); 2399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,islist)); 2401c3fb106SBarry Smith for (i=0; i<dof; i++) { 2419566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i])); 2421c3fb106SBarry Smith } 2431c3fb106SBarry Smith } 2441c3fb106SBarry Smith if (namelist) { 2459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, namelist)); 2461c3fb106SBarry Smith if (dd->fieldname) { 2471c3fb106SBarry Smith for (i=0; i<dof; i++) { 2489566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(dd->fieldname[i],&(*namelist)[i])); 2491c3fb106SBarry Smith } 2501c3fb106SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames"); 2511c3fb106SBarry Smith } 2521c3fb106SBarry Smith if (dmlist) { 2531c3fb106SBarry Smith DM da; 2541c3fb106SBarry Smith 2559566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 2569566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dm->dim)); 2579566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P)); 2589566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p)); 2599566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz)); 2609566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 2619566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da, dd->stencil_type)); 2629566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, dd->s)); 2639566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 2649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,dmlist)); 2659566063dSJacob Faibussowitsch for (i=0; i<dof-1; i++) PetscCall(PetscObjectReference((PetscObject)da)); 2661c3fb106SBarry Smith for (i=0; i<dof; i++) (*dmlist)[i] = da; 2671c3fb106SBarry Smith } 2681c3fb106SBarry Smith PetscFunctionReturn(0); 2691c3fb106SBarry Smith } 2701c3fb106SBarry Smith 27138221697SMatthew G. Knepley PetscErrorCode DMClone_DA(DM dm, DM *newdm) 27238221697SMatthew G. Knepley { 27338221697SMatthew G. Knepley DM_DA *da = (DM_DA *) dm->data; 27438221697SMatthew G. Knepley 27538221697SMatthew G. Knepley PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(DMSetType(*newdm, DMDA)); 2779566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*newdm, dm->dim)); 2789566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P)); 2799566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p)); 2809566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz)); 2819566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*newdm, da->w)); 2829566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*newdm, da->stencil_type)); 2839566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*newdm, da->s)); 2849566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz)); 2859566063dSJacob Faibussowitsch PetscCall(DMSetUp(*newdm)); 28638221697SMatthew G. Knepley PetscFunctionReturn(0); 28738221697SMatthew G. Knepley } 28838221697SMatthew G. Knepley 2894a7a4c06SLawrence Mitchell static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg) 2904a7a4c06SLawrence Mitchell { 2914a7a4c06SLawrence Mitchell DM_DA *da = (DM_DA *)dm->data; 2924a7a4c06SLawrence Mitchell 2934a7a4c06SLawrence Mitchell PetscFunctionBegin; 2944a7a4c06SLawrence Mitchell PetscValidHeaderSpecific(dm,DM_CLASSID,1); 295534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 296131f56ebSLawrence Mitchell *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE; 2974a7a4c06SLawrence Mitchell PetscFunctionReturn(0); 2984a7a4c06SLawrence Mitchell } 2994a7a4c06SLawrence Mitchell 300793f3fe5SMatthew G. Knepley static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 301793f3fe5SMatthew G. Knepley { 302793f3fe5SMatthew G. Knepley PetscFunctionBegin; 3039566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd)); 304793f3fe5SMatthew G. Knepley PetscFunctionReturn(0); 305793f3fe5SMatthew G. Knepley } 306793f3fe5SMatthew G. Knepley 307502a2867SDave May static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[]) 308502a2867SDave May { 309502a2867SDave May PetscInt dim; 310502a2867SDave May DMDAStencilType st; 311502a2867SDave May 312502a2867SDave May PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(DMDAGetNeighbors(dm,ranks)); 3149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&st)); 315502a2867SDave May 316502a2867SDave May switch (dim) { 317502a2867SDave May case 1: 318502a2867SDave May *nranks = 3; 319d2ebda5dSSatish Balay /* if (st == DMDA_STENCIL_STAR) { *nranks = 3; } */ 320502a2867SDave May break; 321502a2867SDave May case 2: 322502a2867SDave May *nranks = 9; 323d2ebda5dSSatish Balay /* if (st == DMDA_STENCIL_STAR) { *nranks = 5; } */ 324502a2867SDave May break; 325502a2867SDave May case 3: 326502a2867SDave May *nranks = 27; 327d2ebda5dSSatish Balay /* if (st == DMDA_STENCIL_STAR) { *nranks = 7; } */ 328502a2867SDave May break; 329502a2867SDave May default: 330502a2867SDave May break; 331502a2867SDave May } 332502a2867SDave May PetscFunctionReturn(0); 333502a2867SDave May } 334502a2867SDave May 3353efe6655SBarry Smith /*MC 3363efe6655SBarry Smith DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions. 3373efe6655SBarry Smith In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points. 3383efe6655SBarry Smith In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width. 3393efe6655SBarry Smith 3403efe6655SBarry 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 3413bd220d7SPatrick Sanan vertex centered; see the documentation for DMSTAG, a similar DM implementation which supports these staggered grids. 3423efe6655SBarry Smith 3433efe6655SBarry Smith Level: intermediate 3443efe6655SBarry Smith 3453bd220d7SPatrick Sanan .seealso: DMType, DMCOMPOSITE, DMSTAG, DMDACreate(), DMCreate(), DMSetType() 3463efe6655SBarry Smith M*/ 3473efe6655SBarry Smith 34805264a50SDave May extern PetscErrorCode DMLocatePoints_DA_Regular(DM,Vec,DMPointLocationType,PetscSF); 3498135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject,PetscViewer); 3501c3fb106SBarry Smith 3518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da) 35247c6ae99SBarry Smith { 35347c6ae99SBarry Smith DM_DA *dd; 35447c6ae99SBarry Smith 35547c6ae99SBarry Smith PetscFunctionBegin; 356a4121054SBarry Smith PetscValidPointer(da,1); 3579566063dSJacob Faibussowitsch PetscCall(PetscNewLog(da,&dd)); 358a4121054SBarry Smith da->data = dd; 35947c6ae99SBarry Smith 360c73cfb54SMatthew G. Knepley da->dim = -1; 361aa219208SBarry Smith dd->interptype = DMDA_Q1; 36247c6ae99SBarry Smith dd->refine_x = 2; 36347c6ae99SBarry Smith dd->refine_y = 2; 36447c6ae99SBarry Smith dd->refine_z = 2; 36581c108dcSJed Brown dd->coarsen_x = 2; 36681c108dcSJed Brown dd->coarsen_y = 2; 36781c108dcSJed Brown dd->coarsen_z = 2; 3680298fd71SBarry Smith dd->fieldname = NULL; 36947c6ae99SBarry Smith dd->nlocal = -1; 37047c6ae99SBarry Smith dd->Nlocal = -1; 37147c6ae99SBarry Smith dd->M = -1; 37247c6ae99SBarry Smith dd->N = -1; 37347c6ae99SBarry Smith dd->P = -1; 37447c6ae99SBarry Smith dd->m = -1; 37547c6ae99SBarry Smith dd->n = -1; 37647c6ae99SBarry Smith dd->p = -1; 37747c6ae99SBarry Smith dd->w = -1; 37847c6ae99SBarry Smith dd->s = -1; 3798865f1eaSKarl Rupp 38047c6ae99SBarry Smith dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1; 38147c6ae99SBarry Smith dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1; 38247c6ae99SBarry Smith 3833e7870d2SPeter Brune dd->Nsub = 1; 3847ddda789SPeter Brune dd->xol = 0; 3857ddda789SPeter Brune dd->yol = 0; 3867ddda789SPeter Brune dd->zol = 0; 387d886c4f4SPeter Brune dd->xo = 0; 388d886c4f4SPeter Brune dd->yo = 0; 389d886c4f4SPeter Brune dd->zo = 0; 390e30e807fSPeter Brune dd->Mo = -1; 391e30e807fSPeter Brune dd->No = -1; 392e30e807fSPeter Brune dd->Po = -1; 39388661749SPeter Brune 3940298fd71SBarry Smith dd->gtol = NULL; 3950298fd71SBarry Smith dd->ltol = NULL; 3960298fd71SBarry Smith dd->ao = NULL; 3979db3d8bcSStefano Zampini PetscStrallocpy(AOBASIC,(char**)&dd->aotype); 39847c6ae99SBarry Smith dd->base = -1; 399bff4a2f0SMatthew G. Knepley dd->bx = DM_BOUNDARY_NONE; 400bff4a2f0SMatthew G. Knepley dd->by = DM_BOUNDARY_NONE; 401bff4a2f0SMatthew G. Knepley dd->bz = DM_BOUNDARY_NONE; 402aa219208SBarry Smith dd->stencil_type = DMDA_STENCIL_BOX; 403aa219208SBarry Smith dd->interptype = DMDA_Q1; 4040298fd71SBarry Smith dd->lx = NULL; 4050298fd71SBarry Smith dd->ly = NULL; 4060298fd71SBarry Smith dd->lz = NULL; 40747c6ae99SBarry Smith 408454e267fSLisandro Dalcin dd->elementtype = DMDA_ELEMENT_Q1; 409454e267fSLisandro Dalcin 410a4121054SBarry Smith da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA; 411a4121054SBarry Smith da->ops->globaltolocalend = DMGlobalToLocalEnd_DA; 412a4121054SBarry Smith da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA; 413a4121054SBarry Smith da->ops->localtoglobalend = DMLocalToGlobalEnd_DA; 414d78e899eSRichard Tran Mills da->ops->localtolocalbegin = DMLocalToLocalBegin_DA; 415d78e899eSRichard Tran Mills da->ops->localtolocalend = DMLocalToLocalEnd_DA; 416a4121054SBarry Smith da->ops->createglobalvector = DMCreateGlobalVector_DA; 417a4121054SBarry Smith da->ops->createlocalvector = DMCreateLocalVector_DA; 41825296bd5SBarry Smith da->ops->createinterpolation = DMCreateInterpolation_DA; 419e727c939SJed Brown da->ops->getcoloring = DMCreateColoring_DA; 42025296bd5SBarry Smith da->ops->creatematrix = DMCreateMatrix_DA; 421a4121054SBarry Smith da->ops->refine = DMRefine_DA; 422a4121054SBarry Smith da->ops->coarsen = DMCoarsen_DA; 423a4121054SBarry Smith da->ops->refinehierarchy = DMRefineHierarchy_DA; 424a4121054SBarry Smith da->ops->coarsenhierarchy = DMCoarsenHierarchy_DA; 4255a84ad33SLisandro Dalcin da->ops->createinjection = DMCreateInjection_DA; 4264a7a4c06SLawrence Mitchell da->ops->hascreateinjection = DMHasCreateInjection_DA; 427a4121054SBarry Smith da->ops->destroy = DMDestroy_DA; 428ea78f98cSLisandro Dalcin da->ops->view = NULL; 429a4121054SBarry Smith da->ops->setfromoptions = DMSetFromOptions_DA; 430a4121054SBarry Smith da->ops->setup = DMSetUp_DA; 43138221697SMatthew G. Knepley da->ops->clone = DMClone_DA; 432b859378eSBarry Smith da->ops->load = DMLoad_DA; 4336636e97aSMatthew G Knepley da->ops->createcoordinatedm = DMCreateCoordinateDM_DA; 434d724dfffSBarry Smith da->ops->createsubdm = DMCreateSubDM_DA; 43516621825SDmitry Karpeev da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA; 436e30e807fSPeter Brune da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA; 437e30e807fSPeter Brune da->ops->createddscatters = DMCreateDomainDecompositionScatters_DA; 438793f3fe5SMatthew G. Knepley da->ops->getdimpoints = DMGetDimPoints_DA; 439502a2867SDave May da->ops->getneighbors = DMGetNeighbors_DA; 44005264a50SDave May da->ops->locatepoints = DMLocatePoints_DA_Regular; 4417c3cd84eSPatrick Sanan da->ops->getcompatibility = DMGetCompatibility_DA; 4429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)da,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_DMDA)); 443a4121054SBarry Smith PetscFunctionReturn(0); 444a4121054SBarry Smith } 44547c6ae99SBarry Smith 446a4121054SBarry Smith /*@ 447a4121054SBarry Smith DMDACreate - Creates a DMDA object. 448a4121054SBarry Smith 4495edebc4bSPatrick Sanan Collective 450a4121054SBarry Smith 451a4121054SBarry Smith Input Parameter: 452a4121054SBarry Smith . comm - The communicator for the DMDA object 453a4121054SBarry Smith 454a4121054SBarry Smith Output Parameter: 455a4121054SBarry Smith . da - The DMDA object 456a4121054SBarry Smith 457e0f5d30fSBarry Smith Level: advanced 458e0f5d30fSBarry Smith 4595edebc4bSPatrick Sanan Developers Note: 4605edebc4bSPatrick Sanan Since there exists DMDACreate1/2/3d() should this routine even exist? 461a4121054SBarry Smith 4625edebc4bSPatrick Sanan .seealso: DMDASetSizes(), DMClone(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 463a4121054SBarry Smith @*/ 4647087cfbeSBarry Smith PetscErrorCode DMDACreate(MPI_Comm comm, DM *da) 465a4121054SBarry Smith { 466a4121054SBarry Smith PetscFunctionBegin; 467a4121054SBarry Smith PetscValidPointer(da,2); 4689566063dSJacob Faibussowitsch PetscCall(DMCreate(comm,da)); 4699566063dSJacob Faibussowitsch PetscCall(DMSetType(*da,DMDA)); 47047c6ae99SBarry Smith PetscFunctionReturn(0); 47147c6ae99SBarry Smith } 472