xref: /petsc/src/dm/impls/da/dacreate.c (revision 2e16c0ce58b3a4ec287cbc0a0807bfb0a0fa5ac9)
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   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   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   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");
16 
17   PetscOptionsHeadBegin(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   PetscOptionsHeadEnd();
79 
80   while (refine--) {
81     if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
82       PetscCall(PetscIntMultError(dd->refine_x,dd->M,&dd->M));
83     } else {
84       PetscCall(PetscIntMultError(dd->refine_x,dd->M-1,&dd->M));
85       dd->M += 1;
86     }
87     if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
88       PetscCall(PetscIntMultError(dd->refine_y,dd->N,&dd->N));
89     } else {
90       PetscCall(PetscIntMultError(dd->refine_y,dd->N-1,&dd->N));
91       dd->N += 1;
92     }
93     if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
94       PetscCall(PetscIntMultError(dd->refine_z,dd->P,&dd->P));
95     } else {
96       PetscCall(PetscIntMultError(dd->refine_z,dd->P-1,&dd->P));
97       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,Mat*);
131 extern PetscErrorCode  DMView_DA(DM,PetscViewer);
132 extern PetscErrorCode  DMSetUp_DA(DM);
133 extern PetscErrorCode  DMDestroy_DA(DM);
134 extern PetscErrorCode  DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
135 extern PetscErrorCode  DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
136 PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM,DM,PetscBool*,PetscBool*);
137 
138 PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
139 {
140   PetscInt         dim,m,n,p,dof,swidth;
141   DMDAStencilType  stencil;
142   DMBoundaryType   bx,by,bz;
143   PetscBool        coors;
144   DM               dac;
145   Vec              c;
146 
147   PetscFunctionBegin;
148   PetscCall(PetscViewerBinaryRead(viewer,&dim,1,NULL,PETSC_INT));
149   PetscCall(PetscViewerBinaryRead(viewer,&m,1,NULL,PETSC_INT));
150   PetscCall(PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT));
151   PetscCall(PetscViewerBinaryRead(viewer,&p,1,NULL,PETSC_INT));
152   PetscCall(PetscViewerBinaryRead(viewer,&dof,1,NULL,PETSC_INT));
153   PetscCall(PetscViewerBinaryRead(viewer,&swidth,1,NULL,PETSC_INT));
154   PetscCall(PetscViewerBinaryRead(viewer,&bx,1,NULL,PETSC_ENUM));
155   PetscCall(PetscViewerBinaryRead(viewer,&by,1,NULL,PETSC_ENUM));
156   PetscCall(PetscViewerBinaryRead(viewer,&bz,1,NULL,PETSC_ENUM));
157   PetscCall(PetscViewerBinaryRead(viewer,&stencil,1,NULL,PETSC_ENUM));
158 
159   PetscCall(DMSetDimension(da, dim));
160   PetscCall(DMDASetSizes(da, m,n,p));
161   PetscCall(DMDASetBoundaryType(da, bx, by, bz));
162   PetscCall(DMDASetDof(da, dof));
163   PetscCall(DMDASetStencilType(da, stencil));
164   PetscCall(DMDASetStencilWidth(da, swidth));
165   PetscCall(DMSetUp(da));
166   PetscCall(PetscViewerBinaryRead(viewer,&coors,1,NULL,PETSC_ENUM));
167   if (coors) {
168     PetscCall(DMGetCoordinateDM(da,&dac));
169     PetscCall(DMCreateGlobalVector(dac,&c));
170     PetscCall(VecLoad(c,viewer));
171     PetscCall(DMSetCoordinates(da,c));
172     PetscCall(VecDestroy(&c));
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
178 {
179   DM_DA         *da = (DM_DA*) dm->data;
180 
181   PetscFunctionBegin;
182   if (subdm) {
183     PetscSF sf;
184     Vec     coords;
185     void   *ctx;
186     /* Cannot use DMClone since the dof stuff is mixed in. Ugh
187     PetscCall(DMClone(dm, subdm)); */
188     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm));
189     PetscCall(DMGetPointSF(dm, &sf));
190     PetscCall(DMSetPointSF(*subdm, sf));
191     PetscCall(DMGetApplicationContext(dm, &ctx));
192     PetscCall(DMSetApplicationContext(*subdm, ctx));
193     PetscCall(DMGetCoordinatesLocal(dm, &coords));
194     if (coords) {
195       PetscCall(DMSetCoordinatesLocal(*subdm, coords));
196     } else {
197       PetscCall(DMGetCoordinates(dm, &coords));
198       if (coords) PetscCall(DMSetCoordinates(*subdm, coords));
199     }
200 
201     PetscCall(DMSetType(*subdm, DMDA));
202     PetscCall(DMSetDimension(*subdm, dm->dim));
203     PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P));
204     PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p));
205     PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz));
206     PetscCall(DMDASetDof(*subdm, numFields));
207     PetscCall(DMDASetStencilType(*subdm, da->stencil_type));
208     PetscCall(DMDASetStencilWidth(*subdm, da->s));
209     PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz));
210   }
211   if (is) {
212     PetscInt *indices, cnt = 0, dof = da->w, i, j;
213 
214     PetscCall(PetscMalloc1(da->Nlocal*numFields/dof, &indices));
215     for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) {
216       for (j = 0; j < numFields; ++j) {
217         indices[cnt++] = dof*i + fields[j];
218       }
219     }
220     PetscCheck(cnt == da->Nlocal*numFields/dof,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %" PetscInt_FMT " does not equal expected value %" PetscInt_FMT, cnt, da->Nlocal*numFields/dof);
221     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is));
222   }
223   PetscFunctionReturn(0);
224 }
225 
226 PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
227 {
228   PetscInt       i;
229   DM_DA          *dd = (DM_DA*)dm->data;
230   PetscInt       dof = dd->w;
231 
232   PetscFunctionBegin;
233   if (len) *len = dof;
234   if (islist) {
235     Vec      v;
236     PetscInt rstart,n;
237 
238     PetscCall(DMGetGlobalVector(dm,&v));
239     PetscCall(VecGetOwnershipRange(v,&rstart,NULL));
240     PetscCall(VecGetLocalSize(v,&n));
241     PetscCall(DMRestoreGlobalVector(dm,&v));
242     PetscCall(PetscMalloc1(dof,islist));
243     for (i=0; i<dof; i++) {
244       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]));
245     }
246   }
247   if (namelist) {
248     PetscCall(PetscMalloc1(dof, namelist));
249     if (dd->fieldname) {
250       for (i=0; i<dof; i++) {
251         PetscCall(PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]));
252       }
253     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
254   }
255   if (dmlist) {
256     DM da;
257 
258     PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da));
259     PetscCall(DMSetDimension(da, dm->dim));
260     PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P));
261     PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p));
262     PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz));
263     PetscCall(DMDASetDof(da, 1));
264     PetscCall(DMDASetStencilType(da, dd->stencil_type));
265     PetscCall(DMDASetStencilWidth(da, dd->s));
266     PetscCall(DMSetUp(da));
267     PetscCall(PetscMalloc1(dof,dmlist));
268     for (i=0; i<dof-1; i++) PetscCall(PetscObjectReference((PetscObject)da));
269     for (i=0; i<dof; i++) (*dmlist)[i] = da;
270   }
271   PetscFunctionReturn(0);
272 }
273 
274 PetscErrorCode DMClone_DA(DM dm, DM *newdm)
275 {
276   DM_DA         *da = (DM_DA *) dm->data;
277 
278   PetscFunctionBegin;
279   PetscCall(DMSetType(*newdm, DMDA));
280   PetscCall(DMSetDimension(*newdm, dm->dim));
281   PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P));
282   PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p));
283   PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz));
284   PetscCall(DMDASetDof(*newdm, da->w));
285   PetscCall(DMDASetStencilType(*newdm, da->stencil_type));
286   PetscCall(DMDASetStencilWidth(*newdm, da->s));
287   PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz));
288   PetscCall(DMSetUp(*newdm));
289   PetscFunctionReturn(0);
290 }
291 
292 static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg)
293 {
294   DM_DA          *da = (DM_DA *)dm->data;
295 
296   PetscFunctionBegin;
297   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
298   PetscValidBoolPointer(flg,2);
299   *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE;
300   PetscFunctionReturn(0);
301 }
302 
303 static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
304 {
305   PetscFunctionBegin;
306   PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd));
307   PetscFunctionReturn(0);
308 }
309 
310 static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
311 {
312   PetscInt dim;
313   DMDAStencilType st;
314 
315   PetscFunctionBegin;
316   PetscCall(DMDAGetNeighbors(dm,ranks));
317   PetscCall(DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&st));
318 
319   switch (dim) {
320     case 1:
321       *nranks = 3;
322       /* if (st == DMDA_STENCIL_STAR) { *nranks = 3; } */
323       break;
324     case 2:
325       *nranks = 9;
326       /* if (st == DMDA_STENCIL_STAR) { *nranks = 5; } */
327       break;
328     case 3:
329       *nranks = 27;
330       /* if (st == DMDA_STENCIL_STAR) { *nranks = 7; } */
331       break;
332     default:
333       break;
334   }
335   PetscFunctionReturn(0);
336 }
337 
338 /*MC
339    DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
340          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
341          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width.
342 
343          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
344          vertex centered; see the documentation for DMSTAG, a similar DM implementation which supports these staggered grids.
345 
346   Level: intermediate
347 
348 .seealso: `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()`
349 M*/
350 
351 extern PetscErrorCode DMLocatePoints_DA_Regular(DM,Vec,DMPointLocationType,PetscSF);
352 PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject,PetscViewer);
353 
354 PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
355 {
356   DM_DA          *dd;
357 
358   PetscFunctionBegin;
359   PetscValidPointer(da,1);
360   PetscCall(PetscNewLog(da,&dd));
361   da->data = dd;
362 
363   da->dim        = -1;
364   dd->interptype = DMDA_Q1;
365   dd->refine_x   = 2;
366   dd->refine_y   = 2;
367   dd->refine_z   = 2;
368   dd->coarsen_x  = 2;
369   dd->coarsen_y  = 2;
370   dd->coarsen_z  = 2;
371   dd->fieldname  = NULL;
372   dd->nlocal     = -1;
373   dd->Nlocal     = -1;
374   dd->M          = -1;
375   dd->N          = -1;
376   dd->P          = -1;
377   dd->m          = -1;
378   dd->n          = -1;
379   dd->p          = -1;
380   dd->w          = -1;
381   dd->s          = -1;
382 
383   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
384   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
385 
386   dd->Nsub            = 1;
387   dd->xol             = 0;
388   dd->yol             = 0;
389   dd->zol             = 0;
390   dd->xo              = 0;
391   dd->yo              = 0;
392   dd->zo              = 0;
393   dd->Mo              = -1;
394   dd->No              = -1;
395   dd->Po              = -1;
396 
397   dd->gtol         = NULL;
398   dd->ltol         = NULL;
399   dd->ao           = NULL;
400   PetscStrallocpy(AOBASIC,(char**)&dd->aotype);
401   dd->base         = -1;
402   dd->bx           = DM_BOUNDARY_NONE;
403   dd->by           = DM_BOUNDARY_NONE;
404   dd->bz           = DM_BOUNDARY_NONE;
405   dd->stencil_type = DMDA_STENCIL_BOX;
406   dd->interptype   = DMDA_Q1;
407   dd->lx           = NULL;
408   dd->ly           = NULL;
409   dd->lz           = NULL;
410 
411   dd->elementtype = DMDA_ELEMENT_Q1;
412 
413   da->ops->globaltolocalbegin          = DMGlobalToLocalBegin_DA;
414   da->ops->globaltolocalend            = DMGlobalToLocalEnd_DA;
415   da->ops->localtoglobalbegin          = DMLocalToGlobalBegin_DA;
416   da->ops->localtoglobalend            = DMLocalToGlobalEnd_DA;
417   da->ops->localtolocalbegin           = DMLocalToLocalBegin_DA;
418   da->ops->localtolocalend             = DMLocalToLocalEnd_DA;
419   da->ops->createglobalvector          = DMCreateGlobalVector_DA;
420   da->ops->createlocalvector           = DMCreateLocalVector_DA;
421   da->ops->createinterpolation         = DMCreateInterpolation_DA;
422   da->ops->getcoloring                 = DMCreateColoring_DA;
423   da->ops->creatematrix                = DMCreateMatrix_DA;
424   da->ops->refine                      = DMRefine_DA;
425   da->ops->coarsen                     = DMCoarsen_DA;
426   da->ops->refinehierarchy             = DMRefineHierarchy_DA;
427   da->ops->coarsenhierarchy            = DMCoarsenHierarchy_DA;
428   da->ops->createinjection             = DMCreateInjection_DA;
429   da->ops->hascreateinjection          = DMHasCreateInjection_DA;
430   da->ops->destroy                     = DMDestroy_DA;
431   da->ops->view                        = NULL;
432   da->ops->setfromoptions              = DMSetFromOptions_DA;
433   da->ops->setup                       = DMSetUp_DA;
434   da->ops->clone                       = DMClone_DA;
435   da->ops->load                        = DMLoad_DA;
436   da->ops->createcoordinatedm          = DMCreateCoordinateDM_DA;
437   da->ops->createsubdm                 = DMCreateSubDM_DA;
438   da->ops->createfielddecomposition    = DMCreateFieldDecomposition_DA;
439   da->ops->createdomaindecomposition   = DMCreateDomainDecomposition_DA;
440   da->ops->createddscatters            = DMCreateDomainDecompositionScatters_DA;
441   da->ops->getdimpoints                = DMGetDimPoints_DA;
442   da->ops->getneighbors                = DMGetNeighbors_DA;
443   da->ops->locatepoints                = DMLocatePoints_DA_Regular;
444   da->ops->getcompatibility            = DMGetCompatibility_DA;
445   PetscCall(PetscObjectComposeFunction((PetscObject)da,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_DMDA));
446   PetscFunctionReturn(0);
447 }
448 
449 /*@
450   DMDACreate - Creates a DMDA object.
451 
452   Collective
453 
454   Input Parameter:
455 . comm - The communicator for the DMDA object
456 
457   Output Parameter:
458 . da  - The DMDA object
459 
460   Level: advanced
461 
462   Developers Note:
463   Since there exists DMDACreate1/2/3d() should this routine even exist?
464 
465 .seealso: `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
466 @*/
467 PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
468 {
469   PetscFunctionBegin;
470   PetscValidPointer(da,2);
471   PetscCall(DMCreate(comm,da));
472   PetscCall(DMSetType(*da,DMDA));
473   PetscFunctionReturn(0);
474 }
475