Actual source code: dacreate.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <petsc/private/dmdaimpl.h>

  4: PetscErrorCode  DMSetFromOptions_DA(PetscOptionItems *PetscOptionsObject,DM da)
  5: {
  7:   DM_DA          *dd    = (DM_DA*)da->data;
  8:   PetscInt       refine = 0,dim = da->dim,maxnlevels = 100,refx[100],refy[100],refz[100],n,i;
  9:   PetscBool      flg;


 14:   if (dd->M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
 15:   if (dd->N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
 16:   if (dd->P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");

 18:   PetscOptionsHead(PetscOptionsObject,"DMDA Options");
 19:   PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,NULL);
 20:   PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL);
 21:   PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL);

 23:   PetscOptionsInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg);
 24:   if (flg) {DMDASetOverlap(da,dd->xol,dd->xol,dd->xol);}
 25:   PetscOptionsInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL);
 26:   if (dim > 1) {PetscOptionsInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL);}
 27:   if (dim > 2) {PetscOptionsInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL);}

 29:   PetscOptionsInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg);
 30:   if (flg) {DMDASetNumLocalSubDomains(da,dd->Nsub);}

 32:   /* Handle DMDA parallel distibution */
 33:   PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL);
 34:   if (dim > 1) {PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL);}
 35:   if (dim > 2) {PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL);}
 36:   /* Handle DMDA refinement */
 37:   PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL);
 38:   if (dim > 1) {PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL);}
 39:   if (dim > 2) {PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL);}
 40:   dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;

 42:   /* Get refinement factors, defaults taken from the coarse DMDA */
 43:   DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);
 44:   for (i=1; i<maxnlevels; i++) {
 45:     refx[i] = refx[0];
 46:     refy[i] = refy[0];
 47:     refz[i] = refz[0];
 48:   }
 49:   n    = maxnlevels;
 50:   PetscOptionsIntArray("-da_refine_hierarchy_x","Refinement factor for each level","None",refx,&n,&flg);
 51:   if (flg) {
 52:     dd->refine_x = refx[0];
 53:     dd->refine_x_hier_n = n;
 54:     PetscMalloc1(n,&dd->refine_x_hier);
 55:     PetscMemcpy(dd->refine_x_hier,refx,n*sizeof(PetscInt));
 56:   }
 57:   if (dim > 1) {
 58:     n    = maxnlevels;
 59:     PetscOptionsIntArray("-da_refine_hierarchy_y","Refinement factor for each level","None",refy,&n,&flg);
 60:     if (flg) {
 61:       dd->refine_y = refy[0];
 62:       dd->refine_y_hier_n = n;
 63:       PetscMalloc1(n,&dd->refine_y_hier);
 64:       PetscMemcpy(dd->refine_y_hier,refy,n*sizeof(PetscInt));
 65:     }
 66:   }
 67:   if (dim > 2) {
 68:     n    = maxnlevels;
 69:     PetscOptionsIntArray("-da_refine_hierarchy_z","Refinement factor for each level","None",refz,&n,&flg);
 70:     if (flg) {
 71:       dd->refine_z = refz[0];
 72:       dd->refine_z_hier_n = n;
 73:       PetscMalloc1(n,&dd->refine_z_hier);
 74:       PetscMemcpy(dd->refine_z_hier,refz,n*sizeof(PetscInt));
 75:     }
 76:   }

 78:   PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,NULL);
 79:   PetscOptionsTail();

 81:   while (refine--) {
 82:     if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
 83:       dd->M = dd->refine_x*dd->M;
 84:     } else {
 85:       dd->M = 1 + dd->refine_x*(dd->M - 1);
 86:     }
 87:     if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
 88:       dd->N = dd->refine_y*dd->N;
 89:     } else {
 90:       dd->N = 1 + dd->refine_y*(dd->N - 1);
 91:     }
 92:     if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
 93:       dd->P = dd->refine_z*dd->P;
 94:     } else {
 95:       dd->P = 1 + dd->refine_z*(dd->P - 1);
 96:     }
 97:     da->levelup++;
 98:     if (da->levelup - da->leveldown >= 0) {
 99:       dd->refine_x = refx[da->levelup - da->leveldown];
100:       dd->refine_y = refy[da->levelup - da->leveldown];
101:       dd->refine_z = refz[da->levelup - da->leveldown];
102:     }
103:     if (da->levelup - da->leveldown >= 1) {
104:       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
105:       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
106:       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
107:     }
108:   }
109:   return(0);
110: }

112: extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
113: extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
114: extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
115: extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
116: extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
117: extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
118: extern PetscErrorCode  DMLocalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
119: extern PetscErrorCode  DMLocalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
120: extern PetscErrorCode  DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
121: extern PetscErrorCode  DMCreateColoring_DA(DM,ISColoringType,ISColoring*);
122: extern PetscErrorCode  DMCreateMatrix_DA(DM,Mat*);
123: extern PetscErrorCode  DMCreateCoordinateDM_DA(DM,DM*);
124: extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
125: extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
126: extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
127: extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
128: extern PetscErrorCode  DMCreateInjection_DA(DM,DM,Mat*);
129: extern PetscErrorCode  DMCreateAggregates_DA(DM,DM,Mat*);
130: extern PetscErrorCode  DMView_DA(DM,PetscViewer);
131: extern PetscErrorCode  DMSetUp_DA(DM);
132: extern PetscErrorCode  DMDestroy_DA(DM);
133: extern PetscErrorCode  DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
134: extern PetscErrorCode  DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

136: PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
137: {
138:   PetscErrorCode   ierr;
139:   PetscInt         dim,m,n,p,dof,swidth;
140:   DMDAStencilType  stencil;
141:   DMBoundaryType   bx,by,bz;
142:   PetscBool        coors;
143:   DM               dac;
144:   Vec              c;

147:   PetscViewerBinaryRead(viewer,&dim,1,NULL,PETSC_INT);
148:   PetscViewerBinaryRead(viewer,&m,1,NULL,PETSC_INT);
149:   PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT);
150:   PetscViewerBinaryRead(viewer,&p,1,NULL,PETSC_INT);
151:   PetscViewerBinaryRead(viewer,&dof,1,NULL,PETSC_INT);
152:   PetscViewerBinaryRead(viewer,&swidth,1,NULL,PETSC_INT);
153:   PetscViewerBinaryRead(viewer,&bx,1,NULL,PETSC_ENUM);
154:   PetscViewerBinaryRead(viewer,&by,1,NULL,PETSC_ENUM);
155:   PetscViewerBinaryRead(viewer,&bz,1,NULL,PETSC_ENUM);
156:   PetscViewerBinaryRead(viewer,&stencil,1,NULL,PETSC_ENUM);

158:   DMSetDimension(da, dim);
159:   DMDASetSizes(da, m,n,p);
160:   DMDASetBoundaryType(da, bx, by, bz);
161:   DMDASetDof(da, dof);
162:   DMDASetStencilType(da, stencil);
163:   DMDASetStencilWidth(da, swidth);
164:   DMSetUp(da);
165:   PetscViewerBinaryRead(viewer,&coors,1,NULL,PETSC_ENUM);
166:   if (coors) {
167:     DMGetCoordinateDM(da,&dac);
168:     DMCreateGlobalVector(dac,&c);
169:     VecLoad(c,viewer);
170:     DMSetCoordinates(da,c);
171:     VecDestroy(&c);
172:   }
173:   return(0);
174: }

176: PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
177: {
178:   DM_DA         *da = (DM_DA*) dm->data;
179:   PetscSection   section;

183:   if (subdm) {
184:     PetscSF sf;
185:     Vec     coords;
186:     void   *ctx;
187:     /* Cannot use DMClone since the dof stuff is mixed in. Ugh
188:     DMClone(dm, subdm); */
189:     DMCreate(PetscObjectComm((PetscObject)dm), subdm);
190:     DMGetPointSF(dm, &sf);
191:     DMSetPointSF(*subdm, sf);
192:     DMGetApplicationContext(dm, &ctx);
193:     DMSetApplicationContext(*subdm, ctx);
194:     DMGetCoordinatesLocal(dm, &coords);
195:     if (coords) {
196:       DMSetCoordinatesLocal(*subdm, coords);
197:     } else {
198:       DMGetCoordinates(dm, &coords);
199:       if (coords) {DMSetCoordinates(*subdm, coords);}
200:     }

202:     DMSetType(*subdm, DMDA);
203:     DMSetDimension(*subdm, dm->dim);
204:     DMDASetSizes(*subdm, da->M, da->N, da->P);
205:     DMDASetNumProcs(*subdm, da->m, da->n, da->p);
206:     DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz);
207:     DMDASetDof(*subdm, numFields);
208:     DMDASetStencilType(*subdm, da->stencil_type);
209:     DMDASetStencilWidth(*subdm, da->s);
210:     DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz);
211:   }
212:   DMGetDefaultSection(dm, &section);
213:   if (section) {
214:     DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
215:   } else {
216:     if (is) {
217:       PetscInt *indices, cnt = 0, dof = da->w, i, j;

219:       PetscMalloc1(da->Nlocal*numFields/dof, &indices);
220:       for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) {
221:         for (j = 0; j < numFields; ++j) {
222:           indices[cnt++] = dof*i + fields[j];
223:         }
224:       }
225:       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);
226:       ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is);
227:     }
228:   }
229:   return(0);
230: }

232: PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
233: {
234:   PetscInt       i;
236:   DM_DA          *dd = (DM_DA*)dm->data;
237:   PetscInt       dof = dd->w;

240:   if (len) *len = dof;
241:   if (islist) {
242:     Vec      v;
243:     PetscInt rstart,n;

245:     DMGetGlobalVector(dm,&v);
246:     VecGetOwnershipRange(v,&rstart,NULL);
247:     VecGetLocalSize(v,&n);
248:     DMRestoreGlobalVector(dm,&v);
249:     PetscMalloc1(dof,islist);
250:     for (i=0; i<dof; i++) {
251:       ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]);
252:     }
253:   }
254:   if (namelist) {
255:     PetscMalloc1(dof, namelist);
256:     if (dd->fieldname) {
257:       for (i=0; i<dof; i++) {
258:         PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);
259:       }
260:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
261:   }
262:   if (dmlist) {
263:     DM da;

265:     DMDACreate(PetscObjectComm((PetscObject)dm), &da);
266:     DMSetDimension(da, dm->dim);
267:     DMDASetSizes(da, dd->M, dd->N, dd->P);
268:     DMDASetNumProcs(da, dd->m, dd->n, dd->p);
269:     DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);
270:     DMDASetDof(da, 1);
271:     DMDASetStencilType(da, dd->stencil_type);
272:     DMDASetStencilWidth(da, dd->s);
273:     DMSetUp(da);
274:     PetscMalloc1(dof,dmlist);
275:     for (i=0; i<dof-1; i++) {PetscObjectReference((PetscObject)da);}
276:     for (i=0; i<dof; i++) (*dmlist)[i] = da;
277:   }
278:   return(0);
279: }

281: PetscErrorCode DMClone_DA(DM dm, DM *newdm)
282: {
283:   DM_DA         *da = (DM_DA *) dm->data;

287:   DMSetType(*newdm, DMDA);
288:   DMSetDimension(*newdm, dm->dim);
289:   DMDASetSizes(*newdm, da->M, da->N, da->P);
290:   DMDASetNumProcs(*newdm, da->m, da->n, da->p);
291:   DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz);
292:   DMDASetDof(*newdm, da->w);
293:   DMDASetStencilType(*newdm, da->stencil_type);
294:   DMDASetStencilWidth(*newdm, da->s);
295:   DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz);
296:   DMSetUp(*newdm);
297:   return(0);
298: }

300: static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
301: {

305:   DMDAGetDepthStratum(dm, dim, pStart, pEnd);
306:   return(0);
307: }

309: static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
310: {
312:   PetscInt dim;
313:   DMDAStencilType st;
314: 
316:   DMDAGetNeighbors(dm,ranks);
317:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&st);

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:   return(0);
336: }

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.

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.

346:   Level: intermediate

348: .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType()
349: M*/

351: extern PetscErrorCode DMProjectFunctionLocal_DA(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
352: extern PetscErrorCode DMComputeL2Diff_DA(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
353: extern PetscErrorCode DMComputeL2GradientDiff_DA(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [],PetscInt, PetscScalar *, void *), void **, Vec,const PetscReal [], PetscReal *);
354: extern PetscErrorCode DMLocatePoints_DA_Regular(DM,Vec,DMPointLocationType,PetscSF);
355: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject,PetscViewer);

357: PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
358: {
360:   DM_DA          *dd;

364:   PetscNewLog(da,&dd);
365:   da->data = dd;

367:   da->dim        = -1;
368:   dd->interptype = DMDA_Q1;
369:   dd->refine_x   = 2;
370:   dd->refine_y   = 2;
371:   dd->refine_z   = 2;
372:   dd->coarsen_x  = 2;
373:   dd->coarsen_y  = 2;
374:   dd->coarsen_z  = 2;
375:   dd->fieldname  = NULL;
376:   dd->nlocal     = -1;
377:   dd->Nlocal     = -1;
378:   dd->M          = -1;
379:   dd->N          = -1;
380:   dd->P          = -1;
381:   dd->m          = -1;
382:   dd->n          = -1;
383:   dd->p          = -1;
384:   dd->w          = -1;
385:   dd->s          = -1;

387:   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
388:   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;

390:   dd->Nsub            = 1;
391:   dd->xol             = 0;
392:   dd->yol             = 0;
393:   dd->zol             = 0;
394:   dd->xo              = 0;
395:   dd->yo              = 0;
396:   dd->zo              = 0;
397:   dd->Mo              = -1;
398:   dd->No              = -1;
399:   dd->Po              = -1;

401:   dd->gtol         = NULL;
402:   dd->ltol         = NULL;
403:   dd->ao           = NULL;
404:   PetscStrallocpy(AOBASIC,(char**)&dd->aotype);
405:   dd->base         = -1;
406:   dd->bx           = DM_BOUNDARY_NONE;
407:   dd->by           = DM_BOUNDARY_NONE;
408:   dd->bz           = DM_BOUNDARY_NONE;
409:   dd->stencil_type = DMDA_STENCIL_BOX;
410:   dd->interptype   = DMDA_Q1;
411:   dd->lx           = NULL;
412:   dd->ly           = NULL;
413:   dd->lz           = NULL;

415:   dd->elementtype = DMDA_ELEMENT_Q1;

417:   da->ops->globaltolocalbegin          = DMGlobalToLocalBegin_DA;
418:   da->ops->globaltolocalend            = DMGlobalToLocalEnd_DA;
419:   da->ops->localtoglobalbegin          = DMLocalToGlobalBegin_DA;
420:   da->ops->localtoglobalend            = DMLocalToGlobalEnd_DA;
421:   da->ops->localtolocalbegin           = DMLocalToLocalBegin_DA;
422:   da->ops->localtolocalend             = DMLocalToLocalEnd_DA;
423:   da->ops->createglobalvector          = DMCreateGlobalVector_DA;
424:   da->ops->createlocalvector           = DMCreateLocalVector_DA;
425:   da->ops->createinterpolation         = DMCreateInterpolation_DA;
426:   da->ops->getcoloring                 = DMCreateColoring_DA;
427:   da->ops->creatematrix                = DMCreateMatrix_DA;
428:   da->ops->refine                      = DMRefine_DA;
429:   da->ops->coarsen                     = DMCoarsen_DA;
430:   da->ops->refinehierarchy             = DMRefineHierarchy_DA;
431:   da->ops->coarsenhierarchy            = DMCoarsenHierarchy_DA;
432:   da->ops->getinjection                = DMCreateInjection_DA;
433:   da->ops->getaggregates               = DMCreateAggregates_DA;
434:   da->ops->destroy                     = DMDestroy_DA;
435:   da->ops->view                        = 0;
436:   da->ops->setfromoptions              = DMSetFromOptions_DA;
437:   da->ops->setup                       = DMSetUp_DA;
438:   da->ops->clone                       = DMClone_DA;
439:   da->ops->load                        = DMLoad_DA;
440:   da->ops->createcoordinatedm          = DMCreateCoordinateDM_DA;
441:   da->ops->createsubdm                 = DMCreateSubDM_DA;
442:   da->ops->createfielddecomposition    = DMCreateFieldDecomposition_DA;
443:   da->ops->createdomaindecomposition   = DMCreateDomainDecomposition_DA;
444:   da->ops->createddscatters            = DMCreateDomainDecompositionScatters_DA;
445:   da->ops->getdimpoints                = DMGetDimPoints_DA;
446:   da->ops->projectfunctionlocal        = DMProjectFunctionLocal_DA;
447:   da->ops->computel2diff               = DMComputeL2Diff_DA;
448:   da->ops->computel2gradientdiff       = DMComputeL2GradientDiff_DA;
449:   da->ops->getneighbors                = DMGetNeighbors_DA;
450:   da->ops->locatepoints                = DMLocatePoints_DA_Regular;
451:   PetscObjectComposeFunction((PetscObject)da,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_DMDA);
452:   return(0);
453: }

455: /*@
456:   DMDACreate - Creates a DMDA object.

458:   Collective on MPI_Comm

460:   Input Parameter:
461: . comm - The communicator for the DMDA object

463:   Output Parameter:
464: . da  - The DMDA object

466:   Level: advanced

468:   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?

470: .keywords: DMDA, create
471: .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
472: @*/
473: PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
474: {

479:   DMCreate(comm,da);
480:   DMSetType(*da,DMDA);
481:   return(0);
482: }