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