xref: /petsc/src/dm/impls/da/dacreate.c (revision 3e7870d263ac1815e57c04a2ae8ebb11bfa82906)
17d0a6c19SBarry Smith 
24035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
59a42bb27SBarry Smith #define __FUNCT__ "DMSetFromOptions_DA"
67087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions_DA(DM da)
747c6ae99SBarry Smith {
847c6ae99SBarry Smith   PetscErrorCode ierr;
947c6ae99SBarry Smith   DM_DA          *dd         = (DM_DA*)da->data;
10397b6216SJed Brown   PetscInt       refine      = 0,maxnlevels = 100,*refx,*refy,*refz,n,i;
117ddda789SPeter Brune   PetscBool      negativeMNP = PETSC_FALSE,bM = PETSC_FALSE,bN = PETSC_FALSE, bP = PETSC_FALSE,flg;
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith   PetscFunctionBegin;
1447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   if (dd->M < 0) {
17235683edSBarry Smith     dd->M       = -dd->M;
18235683edSBarry Smith     bM          = PETSC_TRUE;
19e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
2047c6ae99SBarry Smith   }
2147c6ae99SBarry Smith   if (dd->dim > 1 && dd->N < 0) {
22235683edSBarry Smith     dd->N       = -dd->N;
23235683edSBarry Smith     bN          = PETSC_TRUE;
24e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
2547c6ae99SBarry Smith   }
2647c6ae99SBarry Smith   if (dd->dim > 2 && dd->P < 0) {
27235683edSBarry Smith     dd->P       = -dd->P;
28235683edSBarry Smith     bP          = PETSC_TRUE;
29e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
3047c6ae99SBarry Smith   }
31235683edSBarry Smith 
32f9ba7244SBarry Smith   ierr = PetscOptionsHead("DMDA Options");CHKERRQ(ierr);
330298fd71SBarry Smith   if (bM) {ierr = PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,NULL);CHKERRQ(ierr);}
340298fd71SBarry Smith   if (bN) {ierr = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL);CHKERRQ(ierr);}
350298fd71SBarry Smith   if (bP) {ierr = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL);CHKERRQ(ierr);}
367ddda789SPeter Brune 
377ddda789SPeter Brune   ierr = PetscOptionsInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg);CHKERRQ(ierr);
387ddda789SPeter Brune   if (flg) {ierr = DMDASetOverlap(da,dd->xol,dd->xol,dd->xol);CHKERRQ(ierr);}
397ddda789SPeter Brune   ierr = PetscOptionsInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL);CHKERRQ(ierr);
407ddda789SPeter Brune   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL);CHKERRQ(ierr);}
417ddda789SPeter Brune   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL);CHKERRQ(ierr);}
42*3e7870d2SPeter Brune 
43*3e7870d2SPeter Brune   ierr = PetscOptionsInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg);CHKERRQ(ierr);
44*3e7870d2SPeter Brune   if (flg) {ierr = DMDASetNumLocalSubDomains(da,dd->Nsub);CHKERRQ(ierr);}
45*3e7870d2SPeter Brune 
46aa219208SBarry Smith   /* Handle DMDA parallel distibution */
470298fd71SBarry Smith   ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL);CHKERRQ(ierr);
480298fd71SBarry Smith   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL);CHKERRQ(ierr);}
490298fd71SBarry Smith   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL);CHKERRQ(ierr);}
50aa219208SBarry Smith   /* Handle DMDA refinement */
510298fd71SBarry Smith   ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL);CHKERRQ(ierr);
520298fd71SBarry Smith   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL);CHKERRQ(ierr);}
530298fd71SBarry Smith   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL);CHKERRQ(ierr);}
5438fb6c91SJed Brown   dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;
5547c6ae99SBarry Smith 
56397b6216SJed Brown   /* Get refinement factors, defaults taken from the coarse DMDA */
57397b6216SJed Brown   ierr = PetscMalloc3(maxnlevels,PetscInt,&refx,maxnlevels,PetscInt,&refy,maxnlevels,PetscInt,&refz);CHKERRQ(ierr);
58397b6216SJed Brown   ierr = DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);CHKERRQ(ierr);
59397b6216SJed Brown   for (i=1; i<maxnlevels; i++) {
60397b6216SJed Brown     refx[i] = refx[0];
61397b6216SJed Brown     refy[i] = refy[0];
62397b6216SJed Brown     refz[i] = refz[0];
63397b6216SJed Brown   }
64397b6216SJed Brown   n    = maxnlevels;
650298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
66397b6216SJed Brown   if (da->levelup - da->leveldown >= 0) dd->refine_x = refx[da->levelup - da->leveldown];
67397b6216SJed Brown   if (da->levelup - da->leveldown >= 1) dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
68397b6216SJed Brown   if (dd->dim > 1) {
69397b6216SJed Brown     n    = maxnlevels;
700298fd71SBarry Smith     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
71397b6216SJed Brown     if (da->levelup - da->leveldown >= 0) dd->refine_y = refy[da->levelup - da->leveldown];
72397b6216SJed Brown     if (da->levelup - da->leveldown >= 1) dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
73397b6216SJed Brown   }
74397b6216SJed Brown   if (dd->dim > 2) {
75397b6216SJed Brown     n    = maxnlevels;
760298fd71SBarry Smith     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
77397b6216SJed Brown     if (da->levelup - da->leveldown >= 0) dd->refine_z = refz[da->levelup - da->leveldown];
78397b6216SJed Brown     if (da->levelup - da->leveldown >= 1) dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
79397b6216SJed Brown   }
80397b6216SJed Brown 
810298fd71SBarry Smith   if (negativeMNP) {ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,NULL);CHKERRQ(ierr);}
82f9ba7244SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
83235683edSBarry Smith 
84e0f5d30fSBarry Smith   while (refine--) {
85e0f5d30fSBarry Smith     if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
86e0f5d30fSBarry Smith       dd->M = dd->refine_x*dd->M;
87e0f5d30fSBarry Smith     } else {
88e0f5d30fSBarry Smith       dd->M = 1 + dd->refine_x*(dd->M - 1);
89e0f5d30fSBarry Smith     }
90e0f5d30fSBarry Smith     if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
91e0f5d30fSBarry Smith       dd->N = dd->refine_y*dd->N;
92e0f5d30fSBarry Smith     } else {
93e0f5d30fSBarry Smith       dd->N = 1 + dd->refine_y*(dd->N - 1);
94e0f5d30fSBarry Smith     }
95e0f5d30fSBarry Smith     if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
96e0f5d30fSBarry Smith       dd->P = dd->refine_z*dd->P;
97e0f5d30fSBarry Smith     } else {
98e0f5d30fSBarry Smith       dd->P = 1 + dd->refine_z*(dd->P - 1);
99e0f5d30fSBarry Smith     }
100e0f5d30fSBarry Smith     da->levelup++;
10185fc4b34SJed Brown     if (da->levelup - da->leveldown >= 0) {
10285fc4b34SJed Brown       dd->refine_x = refx[da->levelup - da->leveldown];
10385fc4b34SJed Brown       dd->refine_y = refy[da->levelup - da->leveldown];
10485fc4b34SJed Brown       dd->refine_z = refz[da->levelup - da->leveldown];
10585fc4b34SJed Brown     }
10685fc4b34SJed Brown     if (da->levelup - da->leveldown >= 1) {
10785fc4b34SJed Brown       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
10885fc4b34SJed Brown       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
10985fc4b34SJed Brown       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
11085fc4b34SJed Brown     }
111e0f5d30fSBarry Smith   }
112397b6216SJed Brown   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
11347c6ae99SBarry Smith   PetscFunctionReturn(0);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
1167087cfbeSBarry Smith extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
1177087cfbeSBarry Smith extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
1187087cfbeSBarry Smith extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
1197087cfbeSBarry Smith extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
1207087cfbeSBarry Smith extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
1217087cfbeSBarry Smith extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
122e727c939SJed Brown extern PetscErrorCode  DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
12319fd82e9SBarry Smith extern PetscErrorCode  DMCreateColoring_DA(DM,ISColoringType,MatType,ISColoring*);
12419fd82e9SBarry Smith extern PetscErrorCode  DMCreateMatrix_DA(DM,MatType,Mat*);
1256636e97aSMatthew G Knepley extern PetscErrorCode  DMCreateCoordinateDM_DA(DM,DM*);
1267087cfbeSBarry Smith extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
1277087cfbeSBarry Smith extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
1287087cfbeSBarry Smith extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
1297087cfbeSBarry Smith extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
130e727c939SJed Brown extern PetscErrorCode  DMCreateInjection_DA(DM,DM,VecScatter*);
131e727c939SJed Brown extern PetscErrorCode  DMCreateAggregates_DA(DM,DM,Mat*);
1327087cfbeSBarry Smith extern PetscErrorCode  DMView_DA(DM,PetscViewer);
1337087cfbeSBarry Smith extern PetscErrorCode  DMSetUp_DA(DM);
1347087cfbeSBarry Smith extern PetscErrorCode  DMDestroy_DA(DM);
135e30e807fSPeter Brune extern PetscErrorCode  DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
136e30e807fSPeter Brune extern PetscErrorCode  DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
1379a42bb27SBarry Smith 
138b859378eSBarry Smith #undef __FUNCT__
139b859378eSBarry Smith #define __FUNCT__ "DMLoad_DA"
140b859378eSBarry Smith PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
141b859378eSBarry Smith {
142b859378eSBarry Smith   PetscErrorCode   ierr;
143bc2bf880SBarry Smith   PetscInt         dim,m,n,p,dof,swidth;
144b859378eSBarry Smith   DMDAStencilType  stencil;
145b859378eSBarry Smith   DMDABoundaryType bx,by,bz;
146bc2bf880SBarry Smith   PetscBool        coors;
147bc2bf880SBarry Smith   DM               dac;
148bc2bf880SBarry Smith   Vec              c;
149b859378eSBarry Smith 
150b859378eSBarry Smith   PetscFunctionBegin;
151b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&dim,1,PETSC_INT);CHKERRQ(ierr);
152b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&m,1,PETSC_INT);CHKERRQ(ierr);
153b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);CHKERRQ(ierr);
154b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&p,1,PETSC_INT);CHKERRQ(ierr);
155b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&dof,1,PETSC_INT);CHKERRQ(ierr);
156b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&swidth,1,PETSC_INT);CHKERRQ(ierr);
157b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&bx,1,PETSC_ENUM);CHKERRQ(ierr);
158b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&by,1,PETSC_ENUM);CHKERRQ(ierr);
159b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&bz,1,PETSC_ENUM);CHKERRQ(ierr);
160b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&stencil,1,PETSC_ENUM);CHKERRQ(ierr);
161b859378eSBarry Smith 
162bc2bf880SBarry Smith   ierr = DMDASetDim(da, dim);CHKERRQ(ierr);
163bc2bf880SBarry Smith   ierr = DMDASetSizes(da, m,n,p);CHKERRQ(ierr);
164b859378eSBarry Smith   ierr = DMDASetBoundaryType(da, bx, by, bz);CHKERRQ(ierr);
165b859378eSBarry Smith   ierr = DMDASetDof(da, dof);CHKERRQ(ierr);
166b859378eSBarry Smith   ierr = DMDASetStencilType(da, stencil);CHKERRQ(ierr);
167b859378eSBarry Smith   ierr = DMDASetStencilWidth(da, swidth);CHKERRQ(ierr);
168b859378eSBarry Smith   ierr = DMSetUp(da);CHKERRQ(ierr);
169bc2bf880SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&coors,1,PETSC_ENUM);CHKERRQ(ierr);
170bc2bf880SBarry Smith   if (coors) {
1716636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&dac);CHKERRQ(ierr);
172bc2bf880SBarry Smith     ierr = DMCreateGlobalVector(dac,&c);CHKERRQ(ierr);
173bc2bf880SBarry Smith     ierr = VecLoad(c,viewer);CHKERRQ(ierr);
1746636e97aSMatthew G Knepley     ierr = DMSetCoordinates(da,c);CHKERRQ(ierr);
175bc2bf880SBarry Smith     ierr = VecDestroy(&c);CHKERRQ(ierr);
176bc2bf880SBarry Smith   }
177b859378eSBarry Smith   PetscFunctionReturn(0);
178b859378eSBarry Smith }
179b859378eSBarry Smith 
1801c3fb106SBarry Smith #undef __FUNCT__
181d724dfffSBarry Smith #define __FUNCT__ "DMCreateSubDM_DA"
182d724dfffSBarry Smith PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
183d724dfffSBarry Smith {
184d724dfffSBarry Smith   PetscErrorCode ierr;
185d724dfffSBarry Smith   DM_DA          *da = (DM_DA*)dm->data;
186d724dfffSBarry Smith 
187d724dfffSBarry Smith   PetscFunctionBegin;
188ce94432eSBarry Smith   if (da->dim != 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support only implemented for 2d");
189d724dfffSBarry Smith   if (subdm) {
190ce94432eSBarry Smith     ierr = DMDACreate2d(PetscObjectComm((PetscObject)dm),da->bx,da->by,da->stencil_type,da->M,da->N,da->m,da->n,numFields,da->s,da->lx,da->ly,subdm);CHKERRQ(ierr);
191d724dfffSBarry Smith   }
192d724dfffSBarry Smith   if (is) {
193d724dfffSBarry Smith     PetscInt *indices,cnt = 0, dof = da->w,i,j;
194d724dfffSBarry Smith     ierr = PetscMalloc(sizeof(PetscInt)*da->Nlocal*numFields/dof,&indices);CHKERRQ(ierr);
195d724dfffSBarry Smith     for (i=da->base/dof; i<(da->base+da->Nlocal)/dof; i++) {
196d724dfffSBarry Smith       for (j=0; j<numFields; j++) {
197d724dfffSBarry Smith         indices[cnt++] = dof*i + fields[j];
198d724dfffSBarry Smith       }
199d724dfffSBarry Smith     }
200d724dfffSBarry Smith     if (cnt != da->Nlocal*numFields/dof) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not equal expected value");
201ce94432eSBarry Smith     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),da->Nlocal*numFields/dof,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr);
202d724dfffSBarry Smith   }
203d724dfffSBarry Smith   PetscFunctionReturn(0);
204d724dfffSBarry Smith }
205d724dfffSBarry Smith 
206d724dfffSBarry Smith #undef __FUNCT__
20716621825SDmitry Karpeev #define __FUNCT__ "DMCreateFieldDecomposition_DA"
20816621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
2091c3fb106SBarry Smith {
2101c3fb106SBarry Smith   PetscInt       i;
2111c3fb106SBarry Smith   PetscErrorCode ierr;
2121c3fb106SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
2131c3fb106SBarry Smith   PetscInt       dof = dd->w;
2141c3fb106SBarry Smith 
2151c3fb106SBarry Smith   PetscFunctionBegin;
216731c8d9eSDmitry Karpeev   if (len) *len = dof;
2171c3fb106SBarry Smith   if (islist) {
2181c3fb106SBarry Smith     Vec      v;
2191c3fb106SBarry Smith     PetscInt rstart,n;
2201c3fb106SBarry Smith 
2211c3fb106SBarry Smith     ierr = DMGetGlobalVector(dm,&v);CHKERRQ(ierr);
2220298fd71SBarry Smith     ierr = VecGetOwnershipRange(v,&rstart,NULL);CHKERRQ(ierr);
2231c3fb106SBarry Smith     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2241c3fb106SBarry Smith     ierr = DMRestoreGlobalVector(dm,&v);CHKERRQ(ierr);
2251c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(IS),islist);CHKERRQ(ierr);
2261c3fb106SBarry Smith     for (i=0; i<dof; i++) {
227ce94432eSBarry Smith       ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]);CHKERRQ(ierr);
2281c3fb106SBarry Smith     }
2291c3fb106SBarry Smith   }
2301c3fb106SBarry Smith   if (namelist) {
2311c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(const char*), namelist);CHKERRQ(ierr);
2321c3fb106SBarry Smith     if (dd->fieldname) {
2331c3fb106SBarry Smith       for (i=0; i<dof; i++) {
2341c3fb106SBarry Smith         ierr = PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);CHKERRQ(ierr);
2351c3fb106SBarry Smith       }
2361c3fb106SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
2371c3fb106SBarry Smith   }
2381c3fb106SBarry Smith   if (dmlist) {
2391c3fb106SBarry Smith     DM da;
2401c3fb106SBarry Smith 
241ce94432eSBarry Smith     ierr = DMDACreate(PetscObjectComm((PetscObject)dm), &da);CHKERRQ(ierr);
2421c3fb106SBarry Smith     ierr = DMDASetDim(da, dd->dim);CHKERRQ(ierr);
2431c3fb106SBarry Smith     ierr = DMDASetSizes(da, dd->M, dd->N, dd->P);CHKERRQ(ierr);
2441c3fb106SBarry Smith     ierr = DMDASetNumProcs(da, dd->m, dd->n, dd->p);CHKERRQ(ierr);
2451c3fb106SBarry Smith     ierr = DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);CHKERRQ(ierr);
2461c3fb106SBarry Smith     ierr = DMDASetDof(da, 1);CHKERRQ(ierr);
2471c3fb106SBarry Smith     ierr = DMDASetStencilType(da, dd->stencil_type);CHKERRQ(ierr);
2481c3fb106SBarry Smith     ierr = DMDASetStencilWidth(da, dd->s);CHKERRQ(ierr);
2491c3fb106SBarry Smith     ierr = DMSetUp(da);CHKERRQ(ierr);
2501c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(DM),dmlist);CHKERRQ(ierr);
251c893992eSBarry Smith     for (i=0; i<dof-1; i++) {ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);}
2521c3fb106SBarry Smith     for (i=0; i<dof; i++) (*dmlist)[i] = da;
2531c3fb106SBarry Smith   }
2541c3fb106SBarry Smith   PetscFunctionReturn(0);
2551c3fb106SBarry Smith }
2561c3fb106SBarry Smith 
2573efe6655SBarry Smith /*MC
2583efe6655SBarry Smith    DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
2593efe6655SBarry Smith          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
2603efe6655SBarry Smith          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width.
2613efe6655SBarry Smith 
2623efe6655SBarry 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
2633efe6655SBarry Smith          vertex centered.
2643efe6655SBarry Smith 
2653efe6655SBarry Smith 
2663efe6655SBarry Smith   Level: intermediate
2673efe6655SBarry Smith 
2683efe6655SBarry Smith .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType()
2693efe6655SBarry Smith M*/
2703efe6655SBarry Smith 
2711c3fb106SBarry Smith 
27247c6ae99SBarry Smith #undef __FUNCT__
273a4121054SBarry Smith #define __FUNCT__ "DMCreate_DA"
2748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
27547c6ae99SBarry Smith {
27647c6ae99SBarry Smith   PetscErrorCode ierr;
27747c6ae99SBarry Smith   DM_DA          *dd;
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   PetscFunctionBegin;
280a4121054SBarry Smith   PetscValidPointer(da,1);
281a4121054SBarry Smith   ierr     = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr);
282a4121054SBarry Smith   da->data = dd;
28347c6ae99SBarry Smith 
28447c6ae99SBarry Smith   dd->dim        = -1;
285aa219208SBarry Smith   dd->interptype = DMDA_Q1;
28647c6ae99SBarry Smith   dd->refine_x   = 2;
28747c6ae99SBarry Smith   dd->refine_y   = 2;
28847c6ae99SBarry Smith   dd->refine_z   = 2;
28981c108dcSJed Brown   dd->coarsen_x  = 2;
29081c108dcSJed Brown   dd->coarsen_y  = 2;
29181c108dcSJed Brown   dd->coarsen_z  = 2;
2920298fd71SBarry Smith   dd->fieldname  = NULL;
29347c6ae99SBarry Smith   dd->nlocal     = -1;
29447c6ae99SBarry Smith   dd->Nlocal     = -1;
29547c6ae99SBarry Smith   dd->M          = -1;
29647c6ae99SBarry Smith   dd->N          = -1;
29747c6ae99SBarry Smith   dd->P          = -1;
29847c6ae99SBarry Smith   dd->m          = -1;
29947c6ae99SBarry Smith   dd->n          = -1;
30047c6ae99SBarry Smith   dd->p          = -1;
30147c6ae99SBarry Smith   dd->w          = -1;
30247c6ae99SBarry Smith   dd->s          = -1;
3038865f1eaSKarl Rupp 
30447c6ae99SBarry Smith   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
30547c6ae99SBarry Smith   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
30647c6ae99SBarry Smith 
307*3e7870d2SPeter Brune   dd->Nsub            = 1;
3087ddda789SPeter Brune   dd->xol             = 0;
3097ddda789SPeter Brune   dd->yol             = 0;
3107ddda789SPeter Brune   dd->zol             = 0;
311d886c4f4SPeter Brune   dd->xo              = 0;
312d886c4f4SPeter Brune   dd->yo              = 0;
313d886c4f4SPeter Brune   dd->zo              = 0;
314e30e807fSPeter Brune   dd->Mo              = -1;
315e30e807fSPeter Brune   dd->No              = -1;
316e30e807fSPeter Brune   dd->Po              = -1;
31788661749SPeter Brune 
3180298fd71SBarry Smith   dd->gtol         = NULL;
3190298fd71SBarry Smith   dd->ltog         = NULL;
3200298fd71SBarry Smith   dd->ltol         = NULL;
3210298fd71SBarry Smith   dd->ao           = NULL;
32247c6ae99SBarry Smith   dd->base         = -1;
3231321219cSEthan Coon   dd->bx           = DMDA_BOUNDARY_NONE;
3241321219cSEthan Coon   dd->by           = DMDA_BOUNDARY_NONE;
3251321219cSEthan Coon   dd->bz           = DMDA_BOUNDARY_NONE;
326aa219208SBarry Smith   dd->stencil_type = DMDA_STENCIL_BOX;
327aa219208SBarry Smith   dd->interptype   = DMDA_Q1;
3280298fd71SBarry Smith   dd->idx          = NULL;
32947c6ae99SBarry Smith   dd->Nl           = -1;
3300298fd71SBarry Smith   dd->lx           = NULL;
3310298fd71SBarry Smith   dd->ly           = NULL;
3320298fd71SBarry Smith   dd->lz           = NULL;
33347c6ae99SBarry Smith 
334454e267fSLisandro Dalcin   dd->elementtype = DMDA_ELEMENT_Q1;
335454e267fSLisandro Dalcin 
33619fd82e9SBarry Smith   ierr = PetscStrallocpy(VECSTANDARD,(char**)&da->vectype);CHKERRQ(ierr);
3378865f1eaSKarl Rupp 
338a4121054SBarry Smith   da->ops->globaltolocalbegin          = DMGlobalToLocalBegin_DA;
339a4121054SBarry Smith   da->ops->globaltolocalend            = DMGlobalToLocalEnd_DA;
340a4121054SBarry Smith   da->ops->localtoglobalbegin          = DMLocalToGlobalBegin_DA;
341a4121054SBarry Smith   da->ops->localtoglobalend            = DMLocalToGlobalEnd_DA;
342a4121054SBarry Smith   da->ops->createglobalvector          = DMCreateGlobalVector_DA;
343a4121054SBarry Smith   da->ops->createlocalvector           = DMCreateLocalVector_DA;
34425296bd5SBarry Smith   da->ops->createinterpolation         = DMCreateInterpolation_DA;
345e727c939SJed Brown   da->ops->getcoloring                 = DMCreateColoring_DA;
34625296bd5SBarry Smith   da->ops->creatematrix                = DMCreateMatrix_DA;
347a4121054SBarry Smith   da->ops->refine                      = DMRefine_DA;
348a4121054SBarry Smith   da->ops->coarsen                     = DMCoarsen_DA;
349a4121054SBarry Smith   da->ops->refinehierarchy             = DMRefineHierarchy_DA;
350a4121054SBarry Smith   da->ops->coarsenhierarchy            = DMCoarsenHierarchy_DA;
351e727c939SJed Brown   da->ops->getinjection                = DMCreateInjection_DA;
352e727c939SJed Brown   da->ops->getaggregates               = DMCreateAggregates_DA;
353a4121054SBarry Smith   da->ops->destroy                     = DMDestroy_DA;
354a4121054SBarry Smith   da->ops->view                        = 0;
355a4121054SBarry Smith   da->ops->setfromoptions              = DMSetFromOptions_DA;
356a4121054SBarry Smith   da->ops->setup                       = DMSetUp_DA;
357b859378eSBarry Smith   da->ops->load                        = DMLoad_DA;
3586636e97aSMatthew G Knepley   da->ops->createcoordinatedm          = DMCreateCoordinateDM_DA;
359d724dfffSBarry Smith   da->ops->createsubdm                 = DMCreateSubDM_DA;
36016621825SDmitry Karpeev   da->ops->createfielddecomposition    = DMCreateFieldDecomposition_DA;
361e30e807fSPeter Brune   da->ops->createdomaindecomposition   = DMCreateDomainDecomposition_DA;
362e30e807fSPeter Brune   da->ops->createddscatters            = DMCreateDomainDecompositionScatters_DA;
363a4121054SBarry Smith   PetscFunctionReturn(0);
364a4121054SBarry Smith }
36547c6ae99SBarry Smith 
366a4121054SBarry Smith #undef __FUNCT__
367a4121054SBarry Smith #define __FUNCT__ "DMDACreate"
368a4121054SBarry Smith /*@
369a4121054SBarry Smith   DMDACreate - Creates a DMDA object.
370a4121054SBarry Smith 
371a4121054SBarry Smith   Collective on MPI_Comm
372a4121054SBarry Smith 
373a4121054SBarry Smith   Input Parameter:
374a4121054SBarry Smith . comm - The communicator for the DMDA object
375a4121054SBarry Smith 
376a4121054SBarry Smith   Output Parameter:
377a4121054SBarry Smith . da  - The DMDA object
378a4121054SBarry Smith 
379e0f5d30fSBarry Smith   Level: advanced
380e0f5d30fSBarry Smith 
381e0f5d30fSBarry Smith   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?
382a4121054SBarry Smith 
383a4121054SBarry Smith .keywords: DMDA, create
384e0f5d30fSBarry Smith .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
385a4121054SBarry Smith @*/
3867087cfbeSBarry Smith PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
387a4121054SBarry Smith {
388a4121054SBarry Smith   PetscErrorCode ierr;
389a4121054SBarry Smith 
390a4121054SBarry Smith   PetscFunctionBegin;
391a4121054SBarry Smith   PetscValidPointer(da,2);
392a4121054SBarry Smith   ierr = DMCreate(comm,da);CHKERRQ(ierr);
393a4121054SBarry Smith   ierr = DMSetType(*da,DMDA);CHKERRQ(ierr);
39447c6ae99SBarry Smith   PetscFunctionReturn(0);
39547c6ae99SBarry Smith }
396b859378eSBarry Smith 
397b859378eSBarry Smith 
398