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