xref: /petsc/src/dm/impls/da/dacreate.c (revision cc85fe4ded5189db5e5e073ce90ef04de0003fdb)
1 
2 #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMSetFromOptions_DA"
6 PetscErrorCode  DMSetFromOptions_DA(DM da)
7 {
8   PetscErrorCode ierr;
9   DM_DA          *dd         = (DM_DA*)da->data;
10   PetscInt       refine      = 0,maxnlevels = 100,refx[100],refy[100],refz[100],n,i;
11   PetscBool      negativeMNP = PETSC_FALSE,bM = PETSC_FALSE,bN = PETSC_FALSE, bP = PETSC_FALSE,flg;
12 
13   PetscFunctionBegin;
14   PetscValidHeaderSpecific(da,DM_CLASSID,1);
15 
16   if (dd->M < 0) {
17     dd->M       = -dd->M;
18     bM          = PETSC_TRUE;
19     negativeMNP = PETSC_TRUE;
20   }
21   if (dd->dim > 1 && dd->N < 0) {
22     dd->N       = -dd->N;
23     bN          = PETSC_TRUE;
24     negativeMNP = PETSC_TRUE;
25   }
26   if (dd->dim > 2 && dd->P < 0) {
27     dd->P       = -dd->P;
28     bP          = PETSC_TRUE;
29     negativeMNP = PETSC_TRUE;
30   }
31 
32   ierr = PetscOptionsHead("DMDA Options");CHKERRQ(ierr);
33   if (bM) {ierr = PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,NULL);CHKERRQ(ierr);}
34   if (bN) {ierr = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL);CHKERRQ(ierr);}
35   if (bP) {ierr = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL);CHKERRQ(ierr);}
36 
37   ierr = PetscOptionsInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg);CHKERRQ(ierr);
38   if (flg) {ierr = DMDASetOverlap(da,dd->xol,dd->xol,dd->xol);CHKERRQ(ierr);}
39   ierr = PetscOptionsInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL);CHKERRQ(ierr);
40   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL);CHKERRQ(ierr);}
41   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL);CHKERRQ(ierr);}
42 
43   ierr = PetscOptionsInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg);CHKERRQ(ierr);
44   if (flg) {ierr = DMDASetNumLocalSubDomains(da,dd->Nsub);CHKERRQ(ierr);}
45 
46   /* Handle DMDA parallel distibution */
47   ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL);CHKERRQ(ierr);
48   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL);CHKERRQ(ierr);}
49   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL);CHKERRQ(ierr);}
50   /* Handle DMDA refinement */
51   ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL);CHKERRQ(ierr);
52   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL);CHKERRQ(ierr);}
53   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL);CHKERRQ(ierr);}
54   dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;
55 
56   /* Get refinement factors, defaults taken from the coarse DMDA */
57   ierr = DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);CHKERRQ(ierr);
58   for (i=1; i<maxnlevels; i++) {
59     refx[i] = refx[0];
60     refy[i] = refy[0];
61     refz[i] = refz[0];
62   }
63   n    = maxnlevels;
64   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
65   if (da->levelup - da->leveldown >= 0) dd->refine_x = refx[da->levelup - da->leveldown];
66   if (da->levelup - da->leveldown >= 1) dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
67   if (dd->dim > 1) {
68     n    = maxnlevels;
69     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
70     if (da->levelup - da->leveldown >= 0) dd->refine_y = refy[da->levelup - da->leveldown];
71     if (da->levelup - da->leveldown >= 1) dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
72   }
73   if (dd->dim > 2) {
74     n    = maxnlevels;
75     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
76     if (da->levelup - da->leveldown >= 0) dd->refine_z = refz[da->levelup - da->leveldown];
77     if (da->levelup - da->leveldown >= 1) dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
78   }
79 
80   if (negativeMNP) {ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,NULL);CHKERRQ(ierr);}
81   ierr = PetscOptionsTail();CHKERRQ(ierr);
82 
83   while (refine--) {
84     if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
85       dd->M = dd->refine_x*dd->M;
86     } else {
87       dd->M = 1 + dd->refine_x*(dd->M - 1);
88     }
89     if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
90       dd->N = dd->refine_y*dd->N;
91     } else {
92       dd->N = 1 + dd->refine_y*(dd->N - 1);
93     }
94     if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
95       dd->P = dd->refine_z*dd->P;
96     } else {
97       dd->P = 1 + dd->refine_z*(dd->P - 1);
98     }
99     da->levelup++;
100     if (da->levelup - da->leveldown >= 0) {
101       dd->refine_x = refx[da->levelup - da->leveldown];
102       dd->refine_y = refy[da->levelup - da->leveldown];
103       dd->refine_z = refz[da->levelup - da->leveldown];
104     }
105     if (da->levelup - da->leveldown >= 1) {
106       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
107       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
108       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
109     }
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
115 extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
116 extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
117 extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
118 extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
119 extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
120 extern PetscErrorCode  DMLocalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
121 extern PetscErrorCode  DMLocalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
122 extern PetscErrorCode  DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
123 extern PetscErrorCode  DMCreateColoring_DA(DM,ISColoringType,ISColoring*);
124 extern PetscErrorCode  DMCreateMatrix_DA(DM,Mat*);
125 extern PetscErrorCode  DMCreateCoordinateDM_DA(DM,DM*);
126 extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
127 extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
128 extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
129 extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
130 extern PetscErrorCode  DMCreateInjection_DA(DM,DM,VecScatter*);
131 extern PetscErrorCode  DMCreateAggregates_DA(DM,DM,Mat*);
132 extern PetscErrorCode  DMView_DA(DM,PetscViewer);
133 extern PetscErrorCode  DMSetUp_DA(DM);
134 extern PetscErrorCode  DMDestroy_DA(DM);
135 extern PetscErrorCode  DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
136 extern PetscErrorCode  DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "DMLoad_DA"
140 PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
141 {
142   PetscErrorCode   ierr;
143   PetscInt         dim,m,n,p,dof,swidth;
144   DMDAStencilType  stencil;
145   DMDABoundaryType bx,by,bz;
146   PetscBool        coors;
147   DM               dac;
148   Vec              c;
149 
150   PetscFunctionBegin;
151   ierr = PetscViewerBinaryRead(viewer,&dim,1,PETSC_INT);CHKERRQ(ierr);
152   ierr = PetscViewerBinaryRead(viewer,&m,1,PETSC_INT);CHKERRQ(ierr);
153   ierr = PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);CHKERRQ(ierr);
154   ierr = PetscViewerBinaryRead(viewer,&p,1,PETSC_INT);CHKERRQ(ierr);
155   ierr = PetscViewerBinaryRead(viewer,&dof,1,PETSC_INT);CHKERRQ(ierr);
156   ierr = PetscViewerBinaryRead(viewer,&swidth,1,PETSC_INT);CHKERRQ(ierr);
157   ierr = PetscViewerBinaryRead(viewer,&bx,1,PETSC_ENUM);CHKERRQ(ierr);
158   ierr = PetscViewerBinaryRead(viewer,&by,1,PETSC_ENUM);CHKERRQ(ierr);
159   ierr = PetscViewerBinaryRead(viewer,&bz,1,PETSC_ENUM);CHKERRQ(ierr);
160   ierr = PetscViewerBinaryRead(viewer,&stencil,1,PETSC_ENUM);CHKERRQ(ierr);
161 
162   ierr = DMDASetDim(da, dim);CHKERRQ(ierr);
163   ierr = DMDASetSizes(da, m,n,p);CHKERRQ(ierr);
164   ierr = DMDASetBoundaryType(da, bx, by, bz);CHKERRQ(ierr);
165   ierr = DMDASetDof(da, dof);CHKERRQ(ierr);
166   ierr = DMDASetStencilType(da, stencil);CHKERRQ(ierr);
167   ierr = DMDASetStencilWidth(da, swidth);CHKERRQ(ierr);
168   ierr = DMSetUp(da);CHKERRQ(ierr);
169   ierr = PetscViewerBinaryRead(viewer,&coors,1,PETSC_ENUM);CHKERRQ(ierr);
170   if (coors) {
171     ierr = DMGetCoordinateDM(da,&dac);CHKERRQ(ierr);
172     ierr = DMCreateGlobalVector(dac,&c);CHKERRQ(ierr);
173     ierr = VecLoad(c,viewer);CHKERRQ(ierr);
174     ierr = DMSetCoordinates(da,c);CHKERRQ(ierr);
175     ierr = VecDestroy(&c);CHKERRQ(ierr);
176   }
177   PetscFunctionReturn(0);
178 }
179 
180 #undef __FUNCT__
181 #define __FUNCT__ "DMCreateSubDM_DA"
182 PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
183 {
184   DM_DA         *da = (DM_DA*) dm->data;
185   PetscSection   section;
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   if (subdm) {
190     ierr = DMClone(dm, subdm);CHKERRQ(ierr);
191     ierr = DMDASetDof(*subdm, numFields);CHKERRQ(ierr);
192   }
193   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
194   if (section) {
195     ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
196   } else {
197     if (is) {
198       PetscInt *indices, cnt = 0, dof = da->w, i, j;
199 
200       ierr = PetscMalloc(da->Nlocal*numFields/dof * sizeof(PetscInt), &indices);CHKERRQ(ierr);
201       for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) {
202         for (j = 0; j < numFields; ++j) {
203           indices[cnt++] = dof*i + fields[j];
204         }
205       }
206       if (cnt != da->Nlocal*numFields/dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %d does not equal expected value %d", cnt, da->Nlocal*numFields/dof);
207       ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
208     }
209   }
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "DMCreateFieldDecomposition_DA"
215 PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
216 {
217   PetscInt       i;
218   PetscErrorCode ierr;
219   DM_DA          *dd = (DM_DA*)dm->data;
220   PetscInt       dof = dd->w;
221 
222   PetscFunctionBegin;
223   if (len) *len = dof;
224   if (islist) {
225     Vec      v;
226     PetscInt rstart,n;
227 
228     ierr = DMGetGlobalVector(dm,&v);CHKERRQ(ierr);
229     ierr = VecGetOwnershipRange(v,&rstart,NULL);CHKERRQ(ierr);
230     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
231     ierr = DMRestoreGlobalVector(dm,&v);CHKERRQ(ierr);
232     ierr = PetscMalloc(dof*sizeof(IS),islist);CHKERRQ(ierr);
233     for (i=0; i<dof; i++) {
234       ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]);CHKERRQ(ierr);
235     }
236   }
237   if (namelist) {
238     ierr = PetscMalloc(dof*sizeof(const char*), namelist);CHKERRQ(ierr);
239     if (dd->fieldname) {
240       for (i=0; i<dof; i++) {
241         ierr = PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);CHKERRQ(ierr);
242       }
243     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
244   }
245   if (dmlist) {
246     DM da;
247 
248     ierr = DMDACreate(PetscObjectComm((PetscObject)dm), &da);CHKERRQ(ierr);
249     ierr = DMDASetDim(da, dd->dim);CHKERRQ(ierr);
250     ierr = DMDASetSizes(da, dd->M, dd->N, dd->P);CHKERRQ(ierr);
251     ierr = DMDASetNumProcs(da, dd->m, dd->n, dd->p);CHKERRQ(ierr);
252     ierr = DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);CHKERRQ(ierr);
253     ierr = DMDASetDof(da, 1);CHKERRQ(ierr);
254     ierr = DMDASetStencilType(da, dd->stencil_type);CHKERRQ(ierr);
255     ierr = DMDASetStencilWidth(da, dd->s);CHKERRQ(ierr);
256     ierr = DMSetUp(da);CHKERRQ(ierr);
257     ierr = PetscMalloc(dof*sizeof(DM),dmlist);CHKERRQ(ierr);
258     for (i=0; i<dof-1; i++) {ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);}
259     for (i=0; i<dof; i++) (*dmlist)[i] = da;
260   }
261   PetscFunctionReturn(0);
262 }
263 
264 #undef __FUNCT__
265 #define __FUNCT__ "DMClone_DA"
266 PetscErrorCode DMClone_DA(DM dm, DM *newdm)
267 {
268   DM_DA         *da = (DM_DA *) dm->data;
269   PetscErrorCode ierr;
270 
271   PetscFunctionBegin;
272   ierr = DMSetType(*newdm, DMDA);CHKERRQ(ierr);
273   ierr = DMDASetDim(*newdm, da->dim);CHKERRQ(ierr);
274   ierr = DMDASetSizes(*newdm, da->M, da->N, da->P);CHKERRQ(ierr);
275   ierr = DMDASetNumProcs(*newdm, da->m, da->n, da->p);CHKERRQ(ierr);
276   ierr = DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz);CHKERRQ(ierr);
277   ierr = DMDASetDof(*newdm, da->w);CHKERRQ(ierr);
278   ierr = DMDASetStencilType(*newdm, da->stencil_type);CHKERRQ(ierr);
279   ierr = DMDASetStencilWidth(*newdm, da->s);CHKERRQ(ierr);
280   ierr = DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz);CHKERRQ(ierr);
281   PetscFunctionReturn(0);
282 }
283 
284 /*MC
285    DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
286          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
287          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width.
288 
289          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
290          vertex centered.
291 
292 
293   Level: intermediate
294 
295 .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType()
296 M*/
297 
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "DMCreate_DA"
301 PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
302 {
303   PetscErrorCode ierr;
304   DM_DA          *dd;
305 
306   PetscFunctionBegin;
307   PetscValidPointer(da,1);
308   ierr     = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr);
309   da->data = dd;
310 
311   dd->dim        = -1;
312   dd->interptype = DMDA_Q1;
313   dd->refine_x   = 2;
314   dd->refine_y   = 2;
315   dd->refine_z   = 2;
316   dd->coarsen_x  = 2;
317   dd->coarsen_y  = 2;
318   dd->coarsen_z  = 2;
319   dd->fieldname  = NULL;
320   dd->nlocal     = -1;
321   dd->Nlocal     = -1;
322   dd->M          = -1;
323   dd->N          = -1;
324   dd->P          = -1;
325   dd->m          = -1;
326   dd->n          = -1;
327   dd->p          = -1;
328   dd->w          = -1;
329   dd->s          = -1;
330 
331   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
332   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
333 
334   dd->Nsub            = 1;
335   dd->xol             = 0;
336   dd->yol             = 0;
337   dd->zol             = 0;
338   dd->xo              = 0;
339   dd->yo              = 0;
340   dd->zo              = 0;
341   dd->Mo              = -1;
342   dd->No              = -1;
343   dd->Po              = -1;
344 
345   dd->gtol         = NULL;
346   dd->ltog         = NULL;
347   dd->ltol         = NULL;
348   dd->ao           = NULL;
349   dd->base         = -1;
350   dd->bx           = DMDA_BOUNDARY_NONE;
351   dd->by           = DMDA_BOUNDARY_NONE;
352   dd->bz           = DMDA_BOUNDARY_NONE;
353   dd->stencil_type = DMDA_STENCIL_BOX;
354   dd->interptype   = DMDA_Q1;
355   dd->idx          = NULL;
356   dd->Nl           = -1;
357   dd->lx           = NULL;
358   dd->ly           = NULL;
359   dd->lz           = NULL;
360 
361   dd->elementtype = DMDA_ELEMENT_Q1;
362 
363   da->ops->globaltolocalbegin          = DMGlobalToLocalBegin_DA;
364   da->ops->globaltolocalend            = DMGlobalToLocalEnd_DA;
365   da->ops->localtoglobalbegin          = DMLocalToGlobalBegin_DA;
366   da->ops->localtoglobalend            = DMLocalToGlobalEnd_DA;
367   da->ops->localtolocalbegin           = DMLocalToLocalBegin_DA;
368   da->ops->localtolocalend             = DMLocalToLocalEnd_DA;
369   da->ops->createglobalvector          = DMCreateGlobalVector_DA;
370   da->ops->createlocalvector           = DMCreateLocalVector_DA;
371   da->ops->createinterpolation         = DMCreateInterpolation_DA;
372   da->ops->getcoloring                 = DMCreateColoring_DA;
373   da->ops->creatematrix                = DMCreateMatrix_DA;
374   da->ops->refine                      = DMRefine_DA;
375   da->ops->coarsen                     = DMCoarsen_DA;
376   da->ops->refinehierarchy             = DMRefineHierarchy_DA;
377   da->ops->coarsenhierarchy            = DMCoarsenHierarchy_DA;
378   da->ops->getinjection                = DMCreateInjection_DA;
379   da->ops->getaggregates               = DMCreateAggregates_DA;
380   da->ops->destroy                     = DMDestroy_DA;
381   da->ops->view                        = 0;
382   da->ops->setfromoptions              = DMSetFromOptions_DA;
383   da->ops->setup                       = DMSetUp_DA;
384   da->ops->clone                       = DMClone_DA;
385   da->ops->load                        = DMLoad_DA;
386   da->ops->createcoordinatedm          = DMCreateCoordinateDM_DA;
387   da->ops->createsubdm                 = DMCreateSubDM_DA;
388   da->ops->createfielddecomposition    = DMCreateFieldDecomposition_DA;
389   da->ops->createdomaindecomposition   = DMCreateDomainDecomposition_DA;
390   da->ops->createddscatters            = DMCreateDomainDecompositionScatters_DA;
391   PetscFunctionReturn(0);
392 }
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "DMDACreate"
396 /*@
397   DMDACreate - Creates a DMDA object.
398 
399   Collective on MPI_Comm
400 
401   Input Parameter:
402 . comm - The communicator for the DMDA object
403 
404   Output Parameter:
405 . da  - The DMDA object
406 
407   Level: advanced
408 
409   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?
410 
411 .keywords: DMDA, create
412 .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
413 @*/
414 PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
415 {
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   PetscValidPointer(da,2);
420   ierr = DMCreate(comm,da);CHKERRQ(ierr);
421   ierr = DMSetType(*da,DMDA);CHKERRQ(ierr);
422   PetscFunctionReturn(0);
423 }
424 
425 
426