xref: /petsc/src/dm/interface/dm.c (revision a5057860b72a2a574eaf7af7a4da4e3775d71143)
1 #include <petscsnes.h>
2 #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/
3 
4 PetscClassId  DM_CLASSID;
5 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal;
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "DMCreate"
9 /*@
10   DMCreate - Creates an empty DM object. The type can then be set with DMSetType().
11 
12    If you never  call DMSetType()  it will generate an
13    error when you try to use the vector.
14 
15   Collective on MPI_Comm
16 
17   Input Parameter:
18 . comm - The communicator for the DM object
19 
20   Output Parameter:
21 . dm - The DM object
22 
23   Level: beginner
24 
25 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE
26 @*/
27 PetscErrorCode  DMCreate(MPI_Comm comm,DM *dm)
28 {
29   DM             v;
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   PetscValidPointer(dm,2);
34   *dm = PETSC_NULL;
35 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
36   ierr = VecInitializePackage(PETSC_NULL);CHKERRQ(ierr);
37   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
38   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
39 #endif
40 
41   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
42   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
43 
44 
45   v->workSize     = 0;
46   v->workArray    = PETSC_NULL;
47   v->ltogmap      = PETSC_NULL;
48   v->ltogmapb     = PETSC_NULL;
49   v->bs           = 1;
50   v->coloringtype = IS_COLORING_GLOBAL;
51   v->lf           = PETSC_NULL;
52   v->lj           = PETSC_NULL;
53 
54   *dm = v;
55   PetscFunctionReturn(0);
56 }
57 
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "DMSetVecType"
61 /*@C
62        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
63 
64    Logically Collective on DMDA
65 
66    Input Parameter:
67 +  da - initial distributed array
68 .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
69 
70    Options Database:
71 .   -dm_vec_type ctype
72 
73    Level: intermediate
74 
75 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
76 @*/
77 PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
78 {
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(da,DM_CLASSID,1);
83   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
84   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "DMSetMatType"
90 /*@C
91        DMSetMatType - Sets the type of matrix created with DMCreateMatrix()
92 
93    Logically Collective on DM
94 
95    Input Parameter:
96 +  dm - the DM context
97 .  ctype - the matrix type
98 
99    Options Database:
100 .   -dm_mat_type ctype
101 
102    Level: intermediate
103 
104 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType
105 @*/
106 PetscErrorCode  DMSetMatType(DM dm,const MatType ctype)
107 {
108   PetscErrorCode ierr;
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
111   ierr = PetscFree(dm->mattype);CHKERRQ(ierr);
112   ierr = PetscStrallocpy(ctype,&dm->mattype);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "DMSetOptionsPrefix"
118 /*@C
119    DMSetOptionsPrefix - Sets the prefix used for searching for all
120    DMDA options in the database.
121 
122    Logically Collective on DMDA
123 
124    Input Parameter:
125 +  da - the DMDA context
126 -  prefix - the prefix to prepend to all option names
127 
128    Notes:
129    A hyphen (-) must NOT be given at the beginning of the prefix name.
130    The first character of all runtime options is AUTOMATICALLY the hyphen.
131 
132    Level: advanced
133 
134 .keywords: DMDA, set, options, prefix, database
135 
136 .seealso: DMSetFromOptions()
137 @*/
138 PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
139 {
140   PetscErrorCode ierr;
141 
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
144   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "DMDestroy"
150 /*@
151     DMDestroy - Destroys a vector packer or DMDA.
152 
153     Collective on DM
154 
155     Input Parameter:
156 .   dm - the DM object to destroy
157 
158     Level: developer
159 
160 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
161 
162 @*/
163 PetscErrorCode  DMDestroy(DM *dm)
164 {
165   PetscInt       i, cnt = 0;
166   DMCoarsenHookLink link,next;
167   DMNamedVecLink nlink,nnext;
168   PetscErrorCode ierr;
169 
170   PetscFunctionBegin;
171   if (!*dm) PetscFunctionReturn(0);
172   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
173 
174   /* count all the circular references of DM and its contained Vecs */
175   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
176     if ((*dm)->localin[i])  {cnt++;}
177     if ((*dm)->globalin[i]) {cnt++;}
178   }
179   for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++;
180   if ((*dm)->x) {
181     PetscObject obj;
182     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
183     if (obj == (PetscObject)*dm) cnt++;
184   }
185 
186   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
187   /*
188      Need this test because the dm references the vectors that
189      reference the dm, so destroying the dm calls destroy on the
190      vectors that cause another destroy on the dm
191   */
192   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
193   ((PetscObject) (*dm))->refct = 0;
194   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
195     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
196     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
197   }
198   for (nlink=(*dm)->namedglobal; nlink; nlink=nnext) { /* Destroy the named vectors */
199     nnext = nlink->next;
200     if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name);
201     ierr = PetscFree(nlink->name);CHKERRQ(ierr);
202     ierr = VecDestroy(&nlink->X);CHKERRQ(ierr);
203     ierr = PetscFree(nlink);CHKERRQ(ierr);
204   }
205   (*dm)->namedglobal = PETSC_NULL;
206 
207   /* Destroy the list of hooks */
208   for (link=(*dm)->coarsenhook; link; link=next) {
209     next = link->next;
210     ierr = PetscFree(link);CHKERRQ(ierr);
211   }
212   (*dm)->coarsenhook = PETSC_NULL;
213 
214   if ((*dm)->ctx && (*dm)->ctxdestroy) {
215     ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr);
216   }
217   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
218   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
219   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
220   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
221   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
222   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
223   ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr);
224   ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr);
225   /* if memory was published with AMS then destroy it */
226   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
227 
228   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
229   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
230   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "DMSetUp"
236 /*@
237     DMSetUp - sets up the data structures inside a DM object
238 
239     Collective on DM
240 
241     Input Parameter:
242 .   dm - the DM object to setup
243 
244     Level: developer
245 
246 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
247 
248 @*/
249 PetscErrorCode  DMSetUp(DM dm)
250 {
251   PetscErrorCode ierr;
252 
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
255   if (dm->setupcalled) PetscFunctionReturn(0);
256   if (dm->ops->setup) {
257     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
258   }
259   dm->setupcalled = PETSC_TRUE;
260   PetscFunctionReturn(0);
261 }
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "DMSetFromOptions"
265 /*@
266     DMSetFromOptions - sets parameters in a DM from the options database
267 
268     Collective on DM
269 
270     Input Parameter:
271 .   dm - the DM object to set options for
272 
273     Options Database:
274 +   -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros
275 .   -dm_vec_type <type>  type of vector to create inside DM
276 .   -dm_mat_type <type>  type of matrix to create inside DM
277 -   -dm_coloring_type <global or ghosted>
278 
279     Level: developer
280 
281 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
282 
283 @*/
284 PetscErrorCode  DMSetFromOptions(DM dm)
285 {
286   PetscBool      flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg;
287   PetscErrorCode ierr;
288   char           typeName[256] = MATAIJ;
289 
290   PetscFunctionBegin;
291   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
292   ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr);
293     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
294     ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr);
295     ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr);
296     ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr);
297     ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,PETSC_NULL);CHKERRQ(ierr);
298     ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr);
299     if (flg) {
300       ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr);
301     }
302     ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype?dm->mattype:typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr);
303     if (flg) {
304       ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr);
305     }
306     ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr);
307     if (dm->ops->setfromoptions) {
308       ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
309     }
310     /* process any options handlers added with PetscObjectAddOptionsHandler() */
311     ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr);
312   ierr = PetscOptionsEnd();CHKERRQ(ierr);
313   if (flg1) {
314     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
315   }
316   if (flg2) {
317     PetscViewer viewer;
318 
319     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
320     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
321     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr);
322     ierr = DMView(dm, viewer);CHKERRQ(ierr);
323     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
324   }
325   if (flg3) {
326     PetscViewer viewer;
327 
328     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
329     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
330     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
331     ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr);
332     ierr = DMView(dm, viewer);CHKERRQ(ierr);
333     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
334   }
335   if (flg4) {
336     PetscViewer viewer;
337 
338     ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr);
339     ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
340     ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr);
341     ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr);
342     ierr = DMView(dm, viewer);CHKERRQ(ierr);
343     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
344   }
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "DMView"
350 /*@C
351     DMView - Views a vector packer or DMDA.
352 
353     Collective on DM
354 
355     Input Parameter:
356 +   dm - the DM object to view
357 -   v - the viewer
358 
359     Level: developer
360 
361 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
362 
363 @*/
364 PetscErrorCode  DMView(DM dm,PetscViewer v)
365 {
366   PetscErrorCode ierr;
367 
368   PetscFunctionBegin;
369   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
370  if (!v) {
371     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
372   }
373   if (dm->ops->view) {
374     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "DMCreateGlobalVector"
381 /*@
382     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
383 
384     Collective on DM
385 
386     Input Parameter:
387 .   dm - the DM object
388 
389     Output Parameter:
390 .   vec - the global vector
391 
392     Level: beginner
393 
394 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
395 
396 @*/
397 PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
398 {
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
403   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "DMCreateLocalVector"
409 /*@
410     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
411 
412     Not Collective
413 
414     Input Parameter:
415 .   dm - the DM object
416 
417     Output Parameter:
418 .   vec - the local vector
419 
420     Level: beginner
421 
422 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
423 
424 @*/
425 PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
426 {
427   PetscErrorCode ierr;
428 
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
431   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "DMGetLocalToGlobalMapping"
437 /*@
438    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
439 
440    Collective on DM
441 
442    Input Parameter:
443 .  dm - the DM that provides the mapping
444 
445    Output Parameter:
446 .  ltog - the mapping
447 
448    Level: intermediate
449 
450    Notes:
451    This mapping can then be used by VecSetLocalToGlobalMapping() or
452    MatSetLocalToGlobalMapping().
453 
454 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
455 @*/
456 PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
457 {
458   PetscErrorCode ierr;
459 
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
462   PetscValidPointer(ltog,2);
463   if (!dm->ltogmap) {
464     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
465     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
466   }
467   *ltog = dm->ltogmap;
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
473 /*@
474    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
475 
476    Collective on DM
477 
478    Input Parameter:
479 .  da - the distributed array that provides the mapping
480 
481    Output Parameter:
482 .  ltog - the block mapping
483 
484    Level: intermediate
485 
486    Notes:
487    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
488    MatSetLocalToGlobalMappingBlock().
489 
490 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
491 @*/
492 PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
493 {
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
498   PetscValidPointer(ltog,2);
499   if (!dm->ltogmapb) {
500     PetscInt bs;
501     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
502     if (bs > 1) {
503       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
504       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
505     } else {
506       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
507       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
508     }
509   }
510   *ltog = dm->ltogmapb;
511   PetscFunctionReturn(0);
512 }
513 
514 #undef __FUNCT__
515 #define __FUNCT__ "DMGetBlockSize"
516 /*@
517    DMGetBlockSize - Gets the inherent block size associated with a DM
518 
519    Not Collective
520 
521    Input Parameter:
522 .  dm - the DM with block structure
523 
524    Output Parameter:
525 .  bs - the block size, 1 implies no exploitable block structure
526 
527    Level: intermediate
528 
529 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
530 @*/
531 PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
532 {
533   PetscFunctionBegin;
534   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
535   PetscValidPointer(bs,2);
536   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
537   *bs = dm->bs;
538   PetscFunctionReturn(0);
539 }
540 
541 #undef __FUNCT__
542 #define __FUNCT__ "DMCreateInterpolation"
543 /*@
544     DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
545 
546     Collective on DM
547 
548     Input Parameter:
549 +   dm1 - the DM object
550 -   dm2 - the second, finer DM object
551 
552     Output Parameter:
553 +  mat - the interpolation
554 -  vec - the scaling (optional)
555 
556     Level: developer
557 
558     Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
559         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
560 
561         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
562         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
563 
564 
565 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen()
566 
567 @*/
568 PetscErrorCode  DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
569 {
570   PetscErrorCode ierr;
571 
572   PetscFunctionBegin;
573   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
574   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
575   ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
576   PetscFunctionReturn(0);
577 }
578 
579 #undef __FUNCT__
580 #define __FUNCT__ "DMCreateInjection"
581 /*@
582     DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects
583 
584     Collective on DM
585 
586     Input Parameter:
587 +   dm1 - the DM object
588 -   dm2 - the second, finer DM object
589 
590     Output Parameter:
591 .   ctx - the injection
592 
593     Level: developer
594 
595    Notes:  For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by
596         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
597 
598 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation()
599 
600 @*/
601 PetscErrorCode  DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx)
602 {
603   PetscErrorCode ierr;
604 
605   PetscFunctionBegin;
606   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
607   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
608   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
609   PetscFunctionReturn(0);
610 }
611 
612 #undef __FUNCT__
613 #define __FUNCT__ "DMCreateColoring"
614 /*@C
615     DMCreateColoring - Gets coloring for a DMDA or DMComposite
616 
617     Collective on DM
618 
619     Input Parameter:
620 +   dm - the DM object
621 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
622 -   matype - either MATAIJ or MATBAIJ
623 
624     Output Parameter:
625 .   coloring - the coloring
626 
627     Level: developer
628 
629 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix()
630 
631 @*/
632 PetscErrorCode  DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
633 {
634   PetscErrorCode ierr;
635 
636   PetscFunctionBegin;
637   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
638   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
639   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "DMCreateMatrix"
645 /*@C
646     DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite
647 
648     Collective on DM
649 
650     Input Parameter:
651 +   dm - the DM object
652 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
653             any type which inherits from one of these (such as MATAIJ)
654 
655     Output Parameter:
656 .   mat - the empty Jacobian
657 
658     Level: beginner
659 
660     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
661        do not need to do it yourself.
662 
663        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
664        the nonzero pattern call DMDASetMatPreallocateOnly()
665 
666        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
667        internally by PETSc.
668 
669        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
670        the indices for the global numbering for DMDAs which is complicated.
671 
672 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
673 
674 @*/
675 PetscErrorCode  DMCreateMatrix(DM dm,const MatType mtype,Mat *mat)
676 {
677   PetscErrorCode ierr;
678 
679   PetscFunctionBegin;
680   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
681 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
682   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
683 #endif
684   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
685   PetscValidPointer(mat,3);
686   if (dm->mattype) {
687     ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr);
688   } else {
689     ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr);
690   }
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
696 /*@
697   DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly
698     preallocated but the nonzero structure and zero values will not be set.
699 
700   Logically Collective on DMDA
701 
702   Input Parameter:
703 + dm - the DM
704 - only - PETSC_TRUE if only want preallocation
705 
706   Level: developer
707 .seealso DMCreateMatrix()
708 @*/
709 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
710 {
711   PetscFunctionBegin;
712   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
713   dm->prealloc_only = only;
714   PetscFunctionReturn(0);
715 }
716 
717 #undef __FUNCT__
718 #define __FUNCT__ "DMGetWorkArray"
719 /*@C
720   DMGetWorkArray - Gets a work array guaranteed to be at least the input size
721 
722   Not Collective
723 
724   Input Parameters:
725 + dm - the DM object
726 - size - The minium size
727 
728   Output Parameter:
729 . array - the work array
730 
731   Level: developer
732 
733 .seealso DMDestroy(), DMCreate()
734 @*/
735 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array)
736 {
737   PetscErrorCode ierr;
738 
739   PetscFunctionBegin;
740   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
741   PetscValidPointer(array,3);
742   if (size > dm->workSize) {
743     dm->workSize = size;
744     ierr = PetscFree(dm->workArray);CHKERRQ(ierr);
745     ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr);
746   }
747   *array = dm->workArray;
748   PetscFunctionReturn(0);
749 }
750 
751 #undef __FUNCT__
752 #define __FUNCT__ "DMCreateDecompositionDM"
753 /*@C
754   DMCreateDecompositionDM - creates a DM that encapsulates a decomposition of the original DM.
755 
756   Not Collective
757 
758   Input Parameters:
759 + dm   - the DM object
760 - name - the name of the decomposition
761 
762   Output Parameter:
763 . ddm  - the decomposition DM (PETSC_NULL, if no such decomposition is known)
764 
765   Level: advanced
766 
767 .seealso DMDestroy(), DMCreate(), DMCreateDecomposition()
768 @*/
769 PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm)
770 {
771   PetscErrorCode ierr;
772 
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
775   PetscValidCharPointer(name,2);
776   PetscValidPointer(ddm,3);
777   if(!dm->ops->createdecompositiondm) {
778     *ddm = PETSC_NULL;
779   }
780   else {
781     ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr);
782   }
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "DMCreateFieldIS"
788 /*@C
789   DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field
790 
791   Not collective
792 
793   Input Parameter:
794 . dm - the DM object
795 
796   Output Parameters:
797 + numFields - The number of fields (or PETSC_NULL if not requested)
798 . names     - The name for each field (or PETSC_NULL if not requested)
799 - fields    - The global indices for each field (or PETSC_NULL if not requested)
800 
801   Level: intermediate
802 
803   Notes:
804   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
805   PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with
806   PetscFree().
807 
808 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix()
809 @*/
810 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields)
811 {
812   PetscErrorCode ierr;
813 
814   PetscFunctionBegin;
815   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
816   if (numFields) {
817     PetscValidPointer(numFields,2);
818     *numFields = 0;
819   }
820   if (names) {
821     PetscValidPointer(names,3);
822     *names = PETSC_NULL;
823   }
824   if (fields) {
825     PetscValidPointer(fields,4);
826     *fields = PETSC_NULL;
827   }
828   if(dm->ops->createfieldis) {
829     ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr);
830   }
831   PetscFunctionReturn(0);
832 }
833 
834 
835 #undef __FUNCT__
836 #define __FUNCT__ "DMCreateDecomposition"
837 /*@C
838   DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems:
839                           each IS contains the global indices of the dofs of the corresponding subproblem.
840                           The optional list of DMs define the DM for each subproblem.
841                           Generalizes DMCreateFieldIS().
842 
843   Not collective
844 
845   Input Parameter:
846 . dm - the DM object
847 
848   Output Parameters:
849 + len       - The number of subproblems in the decomposition (or PETSC_NULL if not requested)
850 . namelist  - The name for each subproblem (or PETSC_NULL if not requested)
851 . islist    - The global indices for each subproblem (or PETSC_NULL if not requested)
852 - dmlist    - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined)
853 
854   Level: intermediate
855 
856   Notes:
857   The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with
858   PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(),
859   and all of the arrays should be freed with PetscFree().
860 
861 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS()
862 @*/
863 PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
864 {
865   PetscErrorCode ierr;
866 
867   PetscFunctionBegin;
868   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
869   if (len) {PetscValidPointer(len,2);}
870   if (namelist) {PetscValidPointer(namelist,3);}
871   if (islist) {PetscValidPointer(islist,4);}
872   if (dmlist) {PetscValidPointer(dmlist,5);}
873   if(!dm->ops->createdecomposition) {
874     ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr);
875     /* By default there are no DMs associated with subproblems. */
876     if(dmlist) *dmlist = PETSC_NULL;
877   }
878   else {
879     ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr);
880   }
881   PetscFunctionReturn(0);
882 }
883 
884 
885 #undef __FUNCT__
886 #define __FUNCT__ "DMRefine"
887 /*@
888   DMRefine - Refines a DM object
889 
890   Collective on DM
891 
892   Input Parameter:
893 + dm   - the DM object
894 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
895 
896   Output Parameter:
897 . dmf - the refined DM, or PETSC_NULL
898 
899   Note: If no refinement was done, the return value is PETSC_NULL
900 
901   Level: developer
902 
903 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
904 @*/
905 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
906 {
907   PetscErrorCode ierr;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
911   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
912   if (*dmf) {
913     (*dmf)->ops->creatematrix = dm->ops->creatematrix;
914     (*dmf)->ops->initialguess = dm->ops->initialguess;
915     (*dmf)->ops->function     = dm->ops->function;
916     (*dmf)->ops->functionj    = dm->ops->functionj;
917     if (dm->ops->jacobian != DMComputeJacobianDefault) {
918       (*dmf)->ops->jacobian     = dm->ops->jacobian;
919     }
920     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
921     (*dmf)->ctx     = dm->ctx;
922     (*dmf)->levelup = dm->levelup + 1;
923   }
924   PetscFunctionReturn(0);
925 }
926 
927 #undef __FUNCT__
928 #define __FUNCT__ "DMGetRefineLevel"
929 /*@
930     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
931 
932     Not Collective
933 
934     Input Parameter:
935 .   dm - the DM object
936 
937     Output Parameter:
938 .   level - number of refinements
939 
940     Level: developer
941 
942 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
943 
944 @*/
945 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
946 {
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
949   *level = dm->levelup;
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "DMGlobalToLocalBegin"
955 /*@
956     DMGlobalToLocalBegin - Begins updating local vectors from global vector
957 
958     Neighbor-wise Collective on DM
959 
960     Input Parameters:
961 +   dm - the DM object
962 .   g - the global vector
963 .   mode - INSERT_VALUES or ADD_VALUES
964 -   l - the local vector
965 
966 
967     Level: beginner
968 
969 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
970 
971 @*/
972 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
973 {
974   PetscErrorCode ierr;
975 
976   PetscFunctionBegin;
977   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
978   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
979   PetscFunctionReturn(0);
980 }
981 
982 #undef __FUNCT__
983 #define __FUNCT__ "DMGlobalToLocalEnd"
984 /*@
985     DMGlobalToLocalEnd - Ends updating local vectors from global vector
986 
987     Neighbor-wise Collective on DM
988 
989     Input Parameters:
990 +   dm - the DM object
991 .   g - the global vector
992 .   mode - INSERT_VALUES or ADD_VALUES
993 -   l - the local vector
994 
995 
996     Level: beginner
997 
998 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
999 
1000 @*/
1001 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
1002 {
1003   PetscErrorCode ierr;
1004 
1005   PetscFunctionBegin;
1006   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1007   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "DMLocalToGlobalBegin"
1013 /*@
1014     DMLocalToGlobalBegin - updates global vectors from local vectors
1015 
1016     Neighbor-wise Collective on DM
1017 
1018     Input Parameters:
1019 +   dm - the DM object
1020 .   l - the local vector
1021 .   mode - if INSERT_VALUES then no parallel communication is used, if ADD_VALUES then all ghost points from the same base point accumulate into that
1022            base point.
1023 - - the global vector
1024 
1025     Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. If you would like to simply add the non-ghosted values in the local
1026            array into the global array you need to either (1) zero the ghosted locations and use ADD_VALUES or (2) use INSERT_VALUES into a work global array and then add the work
1027            global array to the final global array with VecAXPY().
1028 
1029     Level: beginner
1030 
1031 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
1032 
1033 @*/
1034 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
1035 {
1036   PetscErrorCode ierr;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1040   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 #undef __FUNCT__
1045 #define __FUNCT__ "DMLocalToGlobalEnd"
1046 /*@
1047     DMLocalToGlobalEnd - updates global vectors from local vectors
1048 
1049     Neighbor-wise Collective on DM
1050 
1051     Input Parameters:
1052 +   dm - the DM object
1053 .   l - the local vector
1054 .   mode - INSERT_VALUES or ADD_VALUES
1055 -   g - the global vector
1056 
1057 
1058     Level: beginner
1059 
1060 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
1061 
1062 @*/
1063 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
1064 {
1065   PetscErrorCode ierr;
1066 
1067   PetscFunctionBegin;
1068   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1069   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
1070   PetscFunctionReturn(0);
1071 }
1072 
1073 #undef __FUNCT__
1074 #define __FUNCT__ "DMComputeJacobianDefault"
1075 /*@
1076     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
1077 
1078     Collective on DM
1079 
1080     Input Parameter:
1081 +   dm - the DM object
1082 .   x - location to compute Jacobian at; may be ignored for linear problems
1083 .   A - matrix that defines the operator for the linear solve
1084 -   B - the matrix used to construct the preconditioner
1085 
1086     Level: developer
1087 
1088 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1089          DMSetFunction()
1090 
1091 @*/
1092 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1093 {
1094   PetscErrorCode ierr;
1095 
1096   PetscFunctionBegin;
1097   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1098   *stflag = SAME_NONZERO_PATTERN;
1099   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
1100   if (A != B) {
1101     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1102     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1103   }
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 #undef __FUNCT__
1108 #define __FUNCT__ "DMCoarsen"
1109 /*@
1110     DMCoarsen - Coarsens a DM object
1111 
1112     Collective on DM
1113 
1114     Input Parameter:
1115 +   dm - the DM object
1116 -   comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
1117 
1118     Output Parameter:
1119 .   dmc - the coarsened DM
1120 
1121     Level: developer
1122 
1123 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1124 
1125 @*/
1126 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
1127 {
1128   PetscErrorCode ierr;
1129   DMCoarsenHookLink link;
1130 
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1133   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
1134   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
1135   (*dmc)->ops->initialguess = dm->ops->initialguess;
1136   (*dmc)->ops->function     = dm->ops->function;
1137   (*dmc)->ops->functionj    = dm->ops->functionj;
1138   if (dm->ops->jacobian != DMComputeJacobianDefault) {
1139     (*dmc)->ops->jacobian     = dm->ops->jacobian;
1140   }
1141   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1142   (*dmc)->ctx       = dm->ctx;
1143   (*dmc)->leveldown = dm->leveldown + 1;
1144   for (link=dm->coarsenhook; link; link=link->next) {
1145     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1146   }
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 #undef __FUNCT__
1151 #define __FUNCT__ "DMCoarsenHookAdd"
1152 /*@
1153    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1154 
1155    Logically Collective
1156 
1157    Input Arguments:
1158 +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1159 .  coarsenhook - function to run when setting up a coarser level
1160 .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1161 -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
1162 
1163    Calling sequence of coarsenhook:
1164 $    coarsenhook(DM fine,DM coarse,void *ctx);
1165 
1166 +  fine - fine level DM
1167 .  coarse - coarse level DM to restrict problem to
1168 -  ctx - optional user-defined function context
1169 
1170    Calling sequence for restricthook:
1171 $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
1172 
1173 +  fine - fine level DM
1174 .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1175 .  inject - matrix restricting by applying the transpose of injection
1176 .  coarse - coarse level DM to update
1177 -  ctx - optional user-defined function context
1178 
1179    Level: advanced
1180 
1181    Notes:
1182    This function is only needed if auxiliary data needs to be set up on coarse grids.
1183 
1184    If this function is called multiple times, the hooks will be run in the order they are added.
1185 
1186    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1187    extract the finest level information from its context (instead of from the SNES).
1188 
1189 .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1190 @*/
1191 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1192 {
1193   PetscErrorCode ierr;
1194   DMCoarsenHookLink link,*p;
1195 
1196   PetscFunctionBegin;
1197   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1198   for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */
1199   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1200   link->coarsenhook = coarsenhook;
1201   link->restricthook = restricthook;
1202   link->ctx = ctx;
1203   link->next = PETSC_NULL;
1204   *p = link;
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "DMRestrict"
1210 /*@
1211    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1212 
1213    Collective if any hooks are
1214 
1215    Input Arguments:
1216 +  fine - finer DM to use as a base
1217 .  restrct - restriction matrix, apply using MatRestrict()
1218 .  inject - injection matrix, also use MatRestrict()
1219 -  coarse - coarer DM to update
1220 
1221    Level: developer
1222 
1223 .seealso: DMCoarsenHookAdd(), MatRestrict()
1224 @*/
1225 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1226 {
1227   PetscErrorCode ierr;
1228   DMCoarsenHookLink link;
1229 
1230   PetscFunctionBegin;
1231   for (link=fine->coarsenhook; link; link=link->next) {
1232     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1233   }
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 #undef __FUNCT__
1238 #define __FUNCT__ "DMGetCoarsenLevel"
1239 /*@
1240     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
1241 
1242     Not Collective
1243 
1244     Input Parameter:
1245 .   dm - the DM object
1246 
1247     Output Parameter:
1248 .   level - number of coarsenings
1249 
1250     Level: developer
1251 
1252 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1253 
1254 @*/
1255 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
1256 {
1257   PetscFunctionBegin;
1258   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1259   *level = dm->leveldown;
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 
1264 
1265 #undef __FUNCT__
1266 #define __FUNCT__ "DMRefineHierarchy"
1267 /*@C
1268     DMRefineHierarchy - Refines a DM object, all levels at once
1269 
1270     Collective on DM
1271 
1272     Input Parameter:
1273 +   dm - the DM object
1274 -   nlevels - the number of levels of refinement
1275 
1276     Output Parameter:
1277 .   dmf - the refined DM hierarchy
1278 
1279     Level: developer
1280 
1281 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1282 
1283 @*/
1284 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
1285 {
1286   PetscErrorCode ierr;
1287 
1288   PetscFunctionBegin;
1289   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1290   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1291   if (nlevels == 0) PetscFunctionReturn(0);
1292   if (dm->ops->refinehierarchy) {
1293     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
1294   } else if (dm->ops->refine) {
1295     PetscInt i;
1296 
1297     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
1298     for (i=1; i<nlevels; i++) {
1299       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
1300     }
1301   } else {
1302     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
1303   }
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 #undef __FUNCT__
1308 #define __FUNCT__ "DMCoarsenHierarchy"
1309 /*@C
1310     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
1311 
1312     Collective on DM
1313 
1314     Input Parameter:
1315 +   dm - the DM object
1316 -   nlevels - the number of levels of coarsening
1317 
1318     Output Parameter:
1319 .   dmc - the coarsened DM hierarchy
1320 
1321     Level: developer
1322 
1323 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1324 
1325 @*/
1326 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
1327 {
1328   PetscErrorCode ierr;
1329 
1330   PetscFunctionBegin;
1331   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1332   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1333   if (nlevels == 0) PetscFunctionReturn(0);
1334   PetscValidPointer(dmc,3);
1335   if (dm->ops->coarsenhierarchy) {
1336     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
1337   } else if (dm->ops->coarsen) {
1338     PetscInt i;
1339 
1340     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
1341     for (i=1; i<nlevels; i++) {
1342       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
1343     }
1344   } else {
1345     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
1346   }
1347   PetscFunctionReturn(0);
1348 }
1349 
1350 #undef __FUNCT__
1351 #define __FUNCT__ "DMCreateAggregates"
1352 /*@
1353    DMCreateAggregates - Gets the aggregates that map between
1354    grids associated with two DMs.
1355 
1356    Collective on DM
1357 
1358    Input Parameters:
1359 +  dmc - the coarse grid DM
1360 -  dmf - the fine grid DM
1361 
1362    Output Parameters:
1363 .  rest - the restriction matrix (transpose of the projection matrix)
1364 
1365    Level: intermediate
1366 
1367 .keywords: interpolation, restriction, multigrid
1368 
1369 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
1370 @*/
1371 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
1372 {
1373   PetscErrorCode ierr;
1374 
1375   PetscFunctionBegin;
1376   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1377   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
1378   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
1379   PetscFunctionReturn(0);
1380 }
1381 
1382 #undef __FUNCT__
1383 #define __FUNCT__ "DMSetApplicationContextDestroy"
1384 /*@C
1385     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
1386 
1387     Not Collective
1388 
1389     Input Parameters:
1390 +   dm - the DM object
1391 -   destroy - the destroy function
1392 
1393     Level: intermediate
1394 
1395 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1396 
1397 @*/
1398 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
1399 {
1400   PetscFunctionBegin;
1401   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1402   dm->ctxdestroy = destroy;
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "DMSetApplicationContext"
1408 /*@
1409     DMSetApplicationContext - Set a user context into a DM object
1410 
1411     Not Collective
1412 
1413     Input Parameters:
1414 +   dm - the DM object
1415 -   ctx - the user context
1416 
1417     Level: intermediate
1418 
1419 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1420 
1421 @*/
1422 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
1423 {
1424   PetscFunctionBegin;
1425   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1426   dm->ctx = ctx;
1427   PetscFunctionReturn(0);
1428 }
1429 
1430 #undef __FUNCT__
1431 #define __FUNCT__ "DMGetApplicationContext"
1432 /*@
1433     DMGetApplicationContext - Gets a user context from a DM object
1434 
1435     Not Collective
1436 
1437     Input Parameter:
1438 .   dm - the DM object
1439 
1440     Output Parameter:
1441 .   ctx - the user context
1442 
1443     Level: intermediate
1444 
1445 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1446 
1447 @*/
1448 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1449 {
1450   PetscFunctionBegin;
1451   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1452   *(void**)ctx = dm->ctx;
1453   PetscFunctionReturn(0);
1454 }
1455 
1456 #undef __FUNCT__
1457 #define __FUNCT__ "DMSetInitialGuess"
1458 /*@C
1459     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1460 
1461     Logically Collective on DM
1462 
1463     Input Parameter:
1464 +   dm - the DM object to destroy
1465 -   f - the function to compute the initial guess
1466 
1467     Level: intermediate
1468 
1469 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1470 
1471 @*/
1472 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1473 {
1474   PetscFunctionBegin;
1475   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1476   dm->ops->initialguess = f;
1477   PetscFunctionReturn(0);
1478 }
1479 
1480 #undef __FUNCT__
1481 #define __FUNCT__ "DMSetFunction"
1482 /*@C
1483     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1484 
1485     Logically Collective on DM
1486 
1487     Input Parameter:
1488 +   dm - the DM object
1489 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1490 
1491     Level: intermediate
1492 
1493     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1494            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1495 
1496 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1497          DMSetJacobian()
1498 
1499 @*/
1500 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1501 {
1502   PetscFunctionBegin;
1503   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1504   dm->ops->function = f;
1505   if (f) {
1506     dm->ops->functionj = f;
1507   }
1508   PetscFunctionReturn(0);
1509 }
1510 
1511 #undef __FUNCT__
1512 #define __FUNCT__ "DMSetJacobian"
1513 /*@C
1514     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1515 
1516     Logically Collective on DM
1517 
1518     Input Parameter:
1519 +   dm - the DM object to destroy
1520 -   f - the function to compute the matrix entries
1521 
1522     Level: intermediate
1523 
1524 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1525          DMSetFunction()
1526 
1527 @*/
1528 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1529 {
1530   PetscFunctionBegin;
1531   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1532   dm->ops->jacobian = f;
1533   PetscFunctionReturn(0);
1534 }
1535 
1536 #undef __FUNCT__
1537 #define __FUNCT__ "DMSetVariableBounds"
1538 /*@C
1539     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
1540 
1541     Logically Collective on DM
1542 
1543     Input Parameter:
1544 +   dm - the DM object
1545 -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
1546 
1547     Level: intermediate
1548 
1549 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1550          DMSetJacobian()
1551 
1552 @*/
1553 PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1554 {
1555   PetscFunctionBegin;
1556   dm->ops->computevariablebounds = f;
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 #undef __FUNCT__
1561 #define __FUNCT__ "DMHasVariableBounds"
1562 /*@
1563     DMHasVariableBounds - does the DM object have a variable bounds function?
1564 
1565     Not Collective
1566 
1567     Input Parameter:
1568 .   dm - the DM object to destroy
1569 
1570     Output Parameter:
1571 .   flg - PETSC_TRUE if the variable bounds function exists
1572 
1573     Level: developer
1574 
1575 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1576 
1577 @*/
1578 PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
1579 {
1580   PetscFunctionBegin;
1581   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
1582   PetscFunctionReturn(0);
1583 }
1584 
1585 #undef __FUNCT__
1586 #define __FUNCT__ "DMComputeVariableBounds"
1587 /*@C
1588     DMComputeVariableBounds - compute variable bounds used by SNESVI.
1589 
1590     Logically Collective on DM
1591 
1592     Input Parameters:
1593 +   dm - the DM object to destroy
1594 -   x  - current solution at which the bounds are computed
1595 
1596     Output parameters:
1597 +   xl - lower bound
1598 -   xu - upper bound
1599 
1600     Level: intermediate
1601 
1602 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1603          DMSetFunction(), DMSetVariableBounds()
1604 
1605 @*/
1606 PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
1607 {
1608   PetscErrorCode ierr;
1609   PetscFunctionBegin;
1610   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1611   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
1612   if(dm->ops->computevariablebounds) {
1613     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
1614   }
1615   else {
1616     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
1617     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
1618   }
1619   PetscFunctionReturn(0);
1620 }
1621 
1622 #undef __FUNCT__
1623 #define __FUNCT__ "DMComputeInitialGuess"
1624 /*@
1625     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1626 
1627     Collective on DM
1628 
1629     Input Parameter:
1630 +   dm - the DM object to destroy
1631 -   x - the vector to hold the initial guess values
1632 
1633     Level: developer
1634 
1635 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1636 
1637 @*/
1638 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1639 {
1640   PetscErrorCode ierr;
1641   PetscFunctionBegin;
1642   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1643   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1644   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1645   PetscFunctionReturn(0);
1646 }
1647 
1648 #undef __FUNCT__
1649 #define __FUNCT__ "DMHasInitialGuess"
1650 /*@
1651     DMHasInitialGuess - does the DM object have an initial guess function
1652 
1653     Not Collective
1654 
1655     Input Parameter:
1656 .   dm - the DM object to destroy
1657 
1658     Output Parameter:
1659 .   flg - PETSC_TRUE if function exists
1660 
1661     Level: developer
1662 
1663 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1664 
1665 @*/
1666 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1667 {
1668   PetscFunctionBegin;
1669   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1670   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1671   PetscFunctionReturn(0);
1672 }
1673 
1674 #undef __FUNCT__
1675 #define __FUNCT__ "DMHasFunction"
1676 /*@
1677     DMHasFunction - does the DM object have a function
1678 
1679     Not Collective
1680 
1681     Input Parameter:
1682 .   dm - the DM object to destroy
1683 
1684     Output Parameter:
1685 .   flg - PETSC_TRUE if function exists
1686 
1687     Level: developer
1688 
1689 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1690 
1691 @*/
1692 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1693 {
1694   PetscFunctionBegin;
1695   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1696   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1697   PetscFunctionReturn(0);
1698 }
1699 
1700 #undef __FUNCT__
1701 #define __FUNCT__ "DMHasJacobian"
1702 /*@
1703     DMHasJacobian - does the DM object have a matrix function
1704 
1705     Not Collective
1706 
1707     Input Parameter:
1708 .   dm - the DM object to destroy
1709 
1710     Output Parameter:
1711 .   flg - PETSC_TRUE if function exists
1712 
1713     Level: developer
1714 
1715 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1716 
1717 @*/
1718 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1719 {
1720   PetscFunctionBegin;
1721   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1722   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1723   PetscFunctionReturn(0);
1724 }
1725 
1726 #undef  __FUNCT__
1727 #define __FUNCT__ "DMSetVec"
1728 /*@C
1729     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
1730 
1731     Collective on DM
1732 
1733     Input Parameter:
1734 +   dm - the DM object
1735 -   x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems.
1736 
1737     Level: developer
1738 
1739 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1740          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
1741 
1742 @*/
1743 PetscErrorCode  DMSetVec(DM dm,Vec x)
1744 {
1745   PetscErrorCode ierr;
1746   PetscFunctionBegin;
1747   if (x) {
1748     if (!dm->x) {
1749       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1750     }
1751     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1752   }
1753   else if(dm->x) {
1754     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
1755   }
1756   PetscFunctionReturn(0);
1757 }
1758 
1759 
1760 #undef __FUNCT__
1761 #define __FUNCT__ "DMComputeFunction"
1762 /*@
1763     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1764 
1765     Collective on DM
1766 
1767     Input Parameter:
1768 +   dm - the DM object to destroy
1769 .   x - the location where the function is evaluationed, may be ignored for linear problems
1770 -   b - the vector to hold the right hand side entries
1771 
1772     Level: developer
1773 
1774 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1775          DMSetJacobian()
1776 
1777 @*/
1778 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1779 {
1780   PetscErrorCode ierr;
1781   PetscFunctionBegin;
1782   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1783   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1784   PetscStackPush("DM user function");
1785   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1786   PetscStackPop;
1787   PetscFunctionReturn(0);
1788 }
1789 
1790 
1791 
1792 #undef __FUNCT__
1793 #define __FUNCT__ "DMComputeJacobian"
1794 /*@
1795     DMComputeJacobian - compute the matrix entries for the solver
1796 
1797     Collective on DM
1798 
1799     Input Parameter:
1800 +   dm - the DM object
1801 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1802 .   A - matrix that defines the operator for the linear solve
1803 -   B - the matrix used to construct the preconditioner
1804 
1805     Level: developer
1806 
1807 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1808          DMSetFunction()
1809 
1810 @*/
1811 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1812 {
1813   PetscErrorCode ierr;
1814 
1815   PetscFunctionBegin;
1816   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1817   if (!dm->ops->jacobian) {
1818     ISColoring     coloring;
1819     MatFDColoring  fd;
1820     const MatType  mtype;
1821 
1822     ierr = PetscObjectGetType((PetscObject)B,&mtype);CHKERRQ(ierr);
1823     ierr = DMCreateColoring(dm,dm->coloringtype,mtype,&coloring);CHKERRQ(ierr);
1824     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1825     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1826     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1827     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1828     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
1829 
1830     dm->fd = fd;
1831     dm->ops->jacobian = DMComputeJacobianDefault;
1832 
1833     /* don't know why this is needed */
1834     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1835   }
1836   if (!x) x = dm->x;
1837   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1838 
1839   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1840   if (x) {
1841     if (!dm->x) {
1842       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1843     }
1844     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1845   }
1846   if (A != B) {
1847     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1848     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1849   }
1850   PetscFunctionReturn(0);
1851 }
1852 
1853 
1854 PetscFList DMList                       = PETSC_NULL;
1855 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1856 
1857 #undef __FUNCT__
1858 #define __FUNCT__ "DMSetType"
1859 /*@C
1860   DMSetType - Builds a DM, for a particular DM implementation.
1861 
1862   Collective on DM
1863 
1864   Input Parameters:
1865 + dm     - The DM object
1866 - method - The name of the DM type
1867 
1868   Options Database Key:
1869 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1870 
1871   Notes:
1872   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1873 
1874   Level: intermediate
1875 
1876 .keywords: DM, set, type
1877 .seealso: DMGetType(), DMCreate()
1878 @*/
1879 PetscErrorCode  DMSetType(DM dm, const DMType method)
1880 {
1881   PetscErrorCode (*r)(DM);
1882   PetscBool      match;
1883   PetscErrorCode ierr;
1884 
1885   PetscFunctionBegin;
1886   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1887   ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1888   if (match) PetscFunctionReturn(0);
1889 
1890   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1891   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1892   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1893 
1894   if (dm->ops->destroy) {
1895     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1896     dm->ops->destroy = PETSC_NULL;
1897   }
1898   ierr = (*r)(dm);CHKERRQ(ierr);
1899   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1900   PetscFunctionReturn(0);
1901 }
1902 
1903 #undef __FUNCT__
1904 #define __FUNCT__ "DMGetType"
1905 /*@C
1906   DMGetType - Gets the DM type name (as a string) from the DM.
1907 
1908   Not Collective
1909 
1910   Input Parameter:
1911 . dm  - The DM
1912 
1913   Output Parameter:
1914 . type - The DM type name
1915 
1916   Level: intermediate
1917 
1918 .keywords: DM, get, type, name
1919 .seealso: DMSetType(), DMCreate()
1920 @*/
1921 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1922 {
1923   PetscErrorCode ierr;
1924 
1925   PetscFunctionBegin;
1926   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1927   PetscValidCharPointer(type,2);
1928   if (!DMRegisterAllCalled) {
1929     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1930   }
1931   *type = ((PetscObject)dm)->type_name;
1932   PetscFunctionReturn(0);
1933 }
1934 
1935 #undef __FUNCT__
1936 #define __FUNCT__ "DMConvert"
1937 /*@C
1938   DMConvert - Converts a DM to another DM, either of the same or different type.
1939 
1940   Collective on DM
1941 
1942   Input Parameters:
1943 + dm - the DM
1944 - newtype - new DM type (use "same" for the same type)
1945 
1946   Output Parameter:
1947 . M - pointer to new DM
1948 
1949   Notes:
1950   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1951   the MPI communicator of the generated DM is always the same as the communicator
1952   of the input DM.
1953 
1954   Level: intermediate
1955 
1956 .seealso: DMCreate()
1957 @*/
1958 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1959 {
1960   DM             B;
1961   char           convname[256];
1962   PetscBool      sametype, issame;
1963   PetscErrorCode ierr;
1964 
1965   PetscFunctionBegin;
1966   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1967   PetscValidType(dm,1);
1968   PetscValidPointer(M,3);
1969   ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1970   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1971   {
1972     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1973 
1974     /*
1975        Order of precedence:
1976        1) See if a specialized converter is known to the current DM.
1977        2) See if a specialized converter is known to the desired DM class.
1978        3) See if a good general converter is registered for the desired class
1979        4) See if a good general converter is known for the current matrix.
1980        5) Use a really basic converter.
1981     */
1982 
1983     /* 1) See if a specialized converter is known to the current DM and the desired class */
1984     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1985     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1986     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1987     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1988     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1989     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1990     if (conv) goto foundconv;
1991 
1992     /* 2)  See if a specialized converter is known to the desired DM class. */
1993     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1994     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1995     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1996     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1997     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1998     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1999     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2000     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
2001     if (conv) {
2002       ierr = DMDestroy(&B);CHKERRQ(ierr);
2003       goto foundconv;
2004     }
2005 
2006 #if 0
2007     /* 3) See if a good general converter is registered for the desired class */
2008     conv = B->ops->convertfrom;
2009     ierr = DMDestroy(&B);CHKERRQ(ierr);
2010     if (conv) goto foundconv;
2011 
2012     /* 4) See if a good general converter is known for the current matrix */
2013     if (dm->ops->convert) {
2014       conv = dm->ops->convert;
2015     }
2016     if (conv) goto foundconv;
2017 #endif
2018 
2019     /* 5) Use a really basic converter. */
2020     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
2021 
2022     foundconv:
2023     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2024     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
2025     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
2026   }
2027   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
2028   PetscFunctionReturn(0);
2029 }
2030 
2031 /*--------------------------------------------------------------------------------------------------------------------*/
2032 
2033 #undef __FUNCT__
2034 #define __FUNCT__ "DMRegister"
2035 /*@C
2036   DMRegister - See DMRegisterDynamic()
2037 
2038   Level: advanced
2039 @*/
2040 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
2041 {
2042   char fullname[PETSC_MAX_PATH_LEN];
2043   PetscErrorCode ierr;
2044 
2045   PetscFunctionBegin;
2046   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
2047   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
2048   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
2049   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
2050   PetscFunctionReturn(0);
2051 }
2052 
2053 
2054 /*--------------------------------------------------------------------------------------------------------------------*/
2055 #undef __FUNCT__
2056 #define __FUNCT__ "DMRegisterDestroy"
2057 /*@C
2058    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
2059 
2060    Not Collective
2061 
2062    Level: advanced
2063 
2064 .keywords: DM, register, destroy
2065 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
2066 @*/
2067 PetscErrorCode  DMRegisterDestroy(void)
2068 {
2069   PetscErrorCode ierr;
2070 
2071   PetscFunctionBegin;
2072   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
2073   DMRegisterAllCalled = PETSC_FALSE;
2074   PetscFunctionReturn(0);
2075 }
2076 
2077 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2078 #include <mex.h>
2079 
2080 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
2081 
2082 #undef __FUNCT__
2083 #define __FUNCT__ "DMComputeFunction_Matlab"
2084 /*
2085    DMComputeFunction_Matlab - Calls the function that has been set with
2086                          DMSetFunctionMatlab().
2087 
2088    For linear problems x is null
2089 
2090 .seealso: DMSetFunction(), DMGetFunction()
2091 */
2092 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
2093 {
2094   PetscErrorCode    ierr;
2095   DMMatlabContext   *sctx;
2096   int               nlhs = 1,nrhs = 4;
2097   mxArray	    *plhs[1],*prhs[4];
2098   long long int     lx = 0,ly = 0,ls = 0;
2099 
2100   PetscFunctionBegin;
2101   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2102   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
2103   PetscCheckSameComm(dm,1,y,3);
2104 
2105   /* call Matlab function in ctx with arguments x and y */
2106   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2107   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
2108   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2109   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
2110   prhs[0] =  mxCreateDoubleScalar((double)ls);
2111   prhs[1] =  mxCreateDoubleScalar((double)lx);
2112   prhs[2] =  mxCreateDoubleScalar((double)ly);
2113   prhs[3] =  mxCreateString(sctx->funcname);
2114   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
2115   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2116   mxDestroyArray(prhs[0]);
2117   mxDestroyArray(prhs[1]);
2118   mxDestroyArray(prhs[2]);
2119   mxDestroyArray(prhs[3]);
2120   mxDestroyArray(plhs[0]);
2121   PetscFunctionReturn(0);
2122 }
2123 
2124 
2125 #undef __FUNCT__
2126 #define __FUNCT__ "DMSetFunctionMatlab"
2127 /*
2128    DMSetFunctionMatlab - Sets the function evaluation routine
2129 
2130 */
2131 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
2132 {
2133   PetscErrorCode    ierr;
2134   DMMatlabContext   *sctx;
2135 
2136   PetscFunctionBegin;
2137   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2138   /* currently sctx is memory bleed */
2139   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2140   if (!sctx) {
2141     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2142   }
2143   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2144   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2145   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
2146   PetscFunctionReturn(0);
2147 }
2148 
2149 #undef __FUNCT__
2150 #define __FUNCT__ "DMComputeJacobian_Matlab"
2151 /*
2152    DMComputeJacobian_Matlab - Calls the function that has been set with
2153                          DMSetJacobianMatlab().
2154 
2155    For linear problems x is null
2156 
2157 .seealso: DMSetFunction(), DMGetFunction()
2158 */
2159 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
2160 {
2161   PetscErrorCode    ierr;
2162   DMMatlabContext   *sctx;
2163   int               nlhs = 2,nrhs = 5;
2164   mxArray	    *plhs[2],*prhs[5];
2165   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
2166 
2167   PetscFunctionBegin;
2168   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2169   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
2170 
2171   /* call MATLAB function in ctx with arguments x, A, and B */
2172   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2173   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
2174   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2175   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
2176   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
2177   prhs[0] =  mxCreateDoubleScalar((double)ls);
2178   prhs[1] =  mxCreateDoubleScalar((double)lx);
2179   prhs[2] =  mxCreateDoubleScalar((double)lA);
2180   prhs[3] =  mxCreateDoubleScalar((double)lB);
2181   prhs[4] =  mxCreateString(sctx->jacname);
2182   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2183   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2184   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2185   mxDestroyArray(prhs[0]);
2186   mxDestroyArray(prhs[1]);
2187   mxDestroyArray(prhs[2]);
2188   mxDestroyArray(prhs[3]);
2189   mxDestroyArray(prhs[4]);
2190   mxDestroyArray(plhs[0]);
2191   mxDestroyArray(plhs[1]);
2192   PetscFunctionReturn(0);
2193 }
2194 
2195 
2196 #undef __FUNCT__
2197 #define __FUNCT__ "DMSetJacobianMatlab"
2198 /*
2199    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
2200 
2201 */
2202 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
2203 {
2204   PetscErrorCode    ierr;
2205   DMMatlabContext   *sctx;
2206 
2207   PetscFunctionBegin;
2208   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2209   /* currently sctx is memory bleed */
2210   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2211   if (!sctx) {
2212     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2213   }
2214   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
2215   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2216   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
2217   PetscFunctionReturn(0);
2218 }
2219 #endif
2220 
2221 #undef __FUNCT__
2222 #define __FUNCT__ "DMLoad"
2223 /*@C
2224   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2225   with DMView().
2226 
2227   Collective on PetscViewer
2228 
2229   Input Parameters:
2230 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2231            some related function before a call to DMLoad().
2232 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2233            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2234 
2235    Level: intermediate
2236 
2237   Notes:
2238   Defaults to the DM DA.
2239 
2240   Notes for advanced users:
2241   Most users should not need to know the details of the binary storage
2242   format, since DMLoad() and DMView() completely hide these details.
2243   But for anyone who's interested, the standard binary matrix storage
2244   format is
2245 .vb
2246      has not yet been determined
2247 .ve
2248 
2249    In addition, PETSc automatically does the byte swapping for
2250 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2251 linux, Windows and the paragon; thus if you write your own binary
2252 read/write routines you have to swap the bytes; see PetscBinaryRead()
2253 and PetscBinaryWrite() to see how this may be done.
2254 
2255   Concepts: vector^loading from file
2256 
2257 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2258 @*/
2259 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2260 {
2261   PetscErrorCode ierr;
2262 
2263   PetscFunctionBegin;
2264   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2265   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2266 
2267   if (!((PetscObject)newdm)->type_name) {
2268     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2269   }
2270   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2271   PetscFunctionReturn(0);
2272 }
2273 
2274 /******************************** FEM Support **********************************/
2275 
2276 #undef __FUNCT__
2277 #define __FUNCT__ "DMPrintCellVector"
2278 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
2279   PetscInt       f;
2280   PetscErrorCode ierr;
2281 
2282   PetscFunctionBegin;
2283   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2284   for(f = 0; f < len; ++f) {
2285     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
2286   }
2287   PetscFunctionReturn(0);
2288 }
2289 
2290 #undef __FUNCT__
2291 #define __FUNCT__ "DMPrintCellMatrix"
2292 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
2293   PetscInt       f, g;
2294   PetscErrorCode ierr;
2295 
2296   PetscFunctionBegin;
2297   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2298   for(f = 0; f < rows; ++f) {
2299     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
2300     for(g = 0; g < cols; ++g) {
2301       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
2302     }
2303     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
2304   }
2305   PetscFunctionReturn(0);
2306 }
2307 
2308 #undef __FUNCT__
2309 #define __FUNCT__ "DMGetLocalFunction"
2310 /*@C
2311   DMGetLocalFunction - Get the local residual function from this DM
2312 
2313   Not collective
2314 
2315   Input Parameter:
2316 . dm - The DM
2317 
2318   Output Parameter:
2319 . lf - The local residual function
2320 
2321    Calling sequence of lf:
2322 $    lf (SNES snes, Vec x, Vec f, void *ctx);
2323 
2324 +  snes - the SNES context
2325 .  x - local vector with the state at which to evaluate residual
2326 .  f - local vector to put residual in
2327 -  ctx - optional user-defined function context
2328 
2329   Level: intermediate
2330 
2331 .seealso DMSetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian()
2332 @*/
2333 PetscErrorCode DMGetLocalFunction(DM dm, PetscErrorCode (**lf)(DM, Vec, Vec, void *))
2334 {
2335   PetscFunctionBegin;
2336   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2337   if (lf) *lf = dm->lf;
2338   PetscFunctionReturn(0);
2339 }
2340 
2341 #undef __FUNCT__
2342 #define __FUNCT__ "DMSetLocalFunction"
2343 /*@C
2344   DMSetLocalFunction - Set the local residual function from this DM
2345 
2346   Not collective
2347 
2348   Input Parameters:
2349 + dm - The DM
2350 - lf - The local residual function
2351 
2352    Calling sequence of lf:
2353 $    lf (SNES snes, Vec x, Vec f, void *ctx);
2354 
2355 +  snes - the SNES context
2356 .  x - local vector with the state at which to evaluate residual
2357 .  f - local vector to put residual in
2358 -  ctx - optional user-defined function context
2359 
2360   Level: intermediate
2361 
2362 .seealso DMGetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian()
2363 @*/
2364 PetscErrorCode DMSetLocalFunction(DM dm, PetscErrorCode (*lf)(DM, Vec, Vec, void *))
2365 {
2366   PetscFunctionBegin;
2367   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2368   dm->lf = lf;
2369   PetscFunctionReturn(0);
2370 }
2371 
2372 #undef __FUNCT__
2373 #define __FUNCT__ "DMGetLocalJacobian"
2374 /*@C
2375   DMGetLocalJacobian - Get the local Jacobian function from this DM
2376 
2377   Not collective
2378 
2379   Input Parameter:
2380 . dm - The DM
2381 
2382   Output Parameter:
2383 . lj - The local Jacobian function
2384 
2385    Calling sequence of lj:
2386 $    lj (SNES snes, Vec x, Mat J, Mat M, void *ctx);
2387 
2388 +  snes - the SNES context
2389 .  x - local vector with the state at which to evaluate residual
2390 .  J - matrix to put Jacobian in
2391 .  M - matrix to use for defining Jacobian preconditioner
2392 -  ctx - optional user-defined function context
2393 
2394   Level: intermediate
2395 
2396 .seealso DMSetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction()
2397 @*/
2398 PetscErrorCode DMGetLocalJacobian(DM dm, PetscErrorCode (**lj)(DM, Vec, Mat, Mat, void *))
2399 {
2400   PetscFunctionBegin;
2401   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2402   if (lj) *lj = dm->lj;
2403   PetscFunctionReturn(0);
2404 }
2405 
2406 #undef __FUNCT__
2407 #define __FUNCT__ "DMSetLocalJacobian"
2408 /*@C
2409   DMSetLocalJacobian - Set the local Jacobian function from this DM
2410 
2411   Not collective
2412 
2413   Input Parameters:
2414 + dm - The DM
2415 - lj - The local Jacobian function
2416 
2417    Calling sequence of lj:
2418 $    lj (SNES snes, Vec x, Mat J, Mat M, void *ctx);
2419 
2420 +  snes - the SNES context
2421 .  x - local vector with the state at which to evaluate residual
2422 .  J - matrix to put Jacobian in
2423 .  M - matrix to use for defining Jacobian preconditioner
2424 -  ctx - optional user-defined function context
2425 
2426   Level: intermediate
2427 
2428 .seealso DMGetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction()
2429 @*/
2430 PetscErrorCode DMSetLocalJacobian(DM dm, PetscErrorCode (*lj)(DM, Vec, Mat,  Mat, void *))
2431 {
2432   PetscFunctionBegin;
2433   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2434   dm->lj = lj;
2435   PetscFunctionReturn(0);
2436 }
2437