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