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