xref: /petsc/src/dm/interface/dm.c (revision 43842a1e2b93b6377be56a453027dbe636fc393d)
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->creatematrix = dm->ops->creatematrix;
772     (*dmf)->ops->initialguess = dm->ops->initialguess;
773     (*dmf)->ops->function     = dm->ops->function;
774     (*dmf)->ops->functionj    = dm->ops->functionj;
775     if (dm->ops->jacobian != DMComputeJacobianDefault) {
776       (*dmf)->ops->jacobian     = dm->ops->jacobian;
777     }
778     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr);
779     (*dmf)->ctx     = dm->ctx;
780     (*dmf)->levelup = dm->levelup + 1;
781   }
782   PetscFunctionReturn(0);
783 }
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "DMGetRefineLevel"
787 /*@
788     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
789 
790     Not Collective
791 
792     Input Parameter:
793 .   dm - the DM object
794 
795     Output Parameter:
796 .   level - number of refinements
797 
798     Level: developer
799 
800 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
801 
802 @*/
803 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
804 {
805   PetscFunctionBegin;
806   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
807   *level = dm->levelup;
808   PetscFunctionReturn(0);
809 }
810 
811 #undef __FUNCT__
812 #define __FUNCT__ "DMGlobalToLocalBegin"
813 /*@
814     DMGlobalToLocalBegin - Begins updating local vectors from global vector
815 
816     Neighbor-wise Collective on DM
817 
818     Input Parameters:
819 +   dm - the DM object
820 .   g - the global vector
821 .   mode - INSERT_VALUES or ADD_VALUES
822 -   l - the local vector
823 
824 
825     Level: beginner
826 
827 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
828 
829 @*/
830 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
831 {
832   PetscErrorCode ierr;
833 
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
836   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
837   PetscFunctionReturn(0);
838 }
839 
840 #undef __FUNCT__
841 #define __FUNCT__ "DMGlobalToLocalEnd"
842 /*@
843     DMGlobalToLocalEnd - Ends updating local vectors from global vector
844 
845     Neighbor-wise Collective on DM
846 
847     Input Parameters:
848 +   dm - the DM object
849 .   g - the global vector
850 .   mode - INSERT_VALUES or ADD_VALUES
851 -   l - the local vector
852 
853 
854     Level: beginner
855 
856 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
857 
858 @*/
859 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
860 {
861   PetscErrorCode ierr;
862 
863   PetscFunctionBegin;
864   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
865   ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr);
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "DMLocalToGlobalBegin"
871 /*@
872     DMLocalToGlobalBegin - updates global vectors from local vectors
873 
874     Neighbor-wise Collective on DM
875 
876     Input Parameters:
877 +   dm - the DM object
878 .   l - the local vector
879 .   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
880            base point.
881 - - the global vector
882 
883     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
884            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
885            global array to the final global array with VecAXPY().
886 
887     Level: beginner
888 
889 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
890 
891 @*/
892 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
893 {
894   PetscErrorCode ierr;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
898   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
899   PetscFunctionReturn(0);
900 }
901 
902 #undef __FUNCT__
903 #define __FUNCT__ "DMLocalToGlobalEnd"
904 /*@
905     DMLocalToGlobalEnd - updates global vectors from local vectors
906 
907     Neighbor-wise Collective on DM
908 
909     Input Parameters:
910 +   dm - the DM object
911 .   l - the local vector
912 .   mode - INSERT_VALUES or ADD_VALUES
913 -   g - the global vector
914 
915 
916     Level: beginner
917 
918 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
919 
920 @*/
921 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
922 {
923   PetscErrorCode ierr;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
927   ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr);
928   PetscFunctionReturn(0);
929 }
930 
931 #undef __FUNCT__
932 #define __FUNCT__ "DMComputeJacobianDefault"
933 /*@
934     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
935 
936     Collective on DM
937 
938     Input Parameter:
939 +   dm - the DM object
940 .   x - location to compute Jacobian at; may be ignored for linear problems
941 .   A - matrix that defines the operator for the linear solve
942 -   B - the matrix used to construct the preconditioner
943 
944     Level: developer
945 
946 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
947          DMSetFunction()
948 
949 @*/
950 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
951 {
952   PetscErrorCode ierr;
953 
954   PetscFunctionBegin;
955   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
956   *stflag = SAME_NONZERO_PATTERN;
957   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
958   if (A != B) {
959     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
960     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
961   }
962   PetscFunctionReturn(0);
963 }
964 
965 #undef __FUNCT__
966 #define __FUNCT__ "DMCoarsen"
967 /*@
968     DMCoarsen - Coarsens a DM object
969 
970     Collective on DM
971 
972     Input Parameter:
973 +   dm - the DM object
974 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
975 
976     Output Parameter:
977 .   dmc - the coarsened DM
978 
979     Level: developer
980 
981 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
982 
983 @*/
984 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
985 {
986   PetscErrorCode ierr;
987   DMCoarsenHookLink link;
988 
989   PetscFunctionBegin;
990   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
991   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
992   (*dmc)->ops->creatematrix = dm->ops->creatematrix;
993   (*dmc)->ops->initialguess = dm->ops->initialguess;
994   (*dmc)->ops->function     = dm->ops->function;
995   (*dmc)->ops->functionj    = dm->ops->functionj;
996   if (dm->ops->jacobian != DMComputeJacobianDefault) {
997     (*dmc)->ops->jacobian     = dm->ops->jacobian;
998   }
999   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr);
1000   (*dmc)->ctx       = dm->ctx;
1001   (*dmc)->leveldown = dm->leveldown + 1;
1002   for (link=dm->coarsenhook; link; link=link->next) {
1003     if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);}
1004   }
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "DMCoarsenHookAdd"
1010 /*@
1011    DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid
1012 
1013    Logically Collective
1014 
1015    Input Arguments:
1016 +  fine - nonlinear solver context on which to run a hook when restricting to a coarser level
1017 .  coarsenhook - function to run when setting up a coarser level
1018 .  restricthook - function to run to update data on coarser levels (once per SNESSolve())
1019 -  ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL)
1020 
1021    Calling sequence of coarsenhook:
1022 $    coarsenhook(DM fine,DM coarse,void *ctx);
1023 
1024 +  fine - fine level DM
1025 .  coarse - coarse level DM to restrict problem to
1026 -  ctx - optional user-defined function context
1027 
1028    Calling sequence for restricthook:
1029 $    restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx)
1030 
1031 +  fine - fine level DM
1032 .  mrestrict - matrix restricting a fine-level solution to the coarse grid
1033 .  inject - matrix restricting by applying the transpose of injection
1034 .  coarse - coarse level DM to update
1035 -  ctx - optional user-defined function context
1036 
1037    Level: advanced
1038 
1039    Notes:
1040    This function is only needed if auxiliary data needs to be set up on coarse grids.
1041 
1042    If this function is called multiple times, the hooks will be run in the order they are added.
1043 
1044    The
1045 
1046    In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to
1047    extract the finest level information from its context (instead of from the SNES).
1048 
1049 .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate()
1050 @*/
1051 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx)
1052 {
1053   PetscErrorCode ierr;
1054   DMCoarsenHookLink link,*p;
1055 
1056   PetscFunctionBegin;
1057   PetscValidHeaderSpecific(fine,DM_CLASSID,1);
1058   for (p=&fine->coarsenhook; *p; *p=(*p)->next) {} /* Scan to the end of the current list of hooks */
1059   ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr);
1060   link->coarsenhook = coarsenhook;
1061   link->restricthook = restricthook;
1062   link->ctx = ctx;
1063   link->next = PETSC_NULL;
1064   *p = link;
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 #undef __FUNCT__
1069 #define __FUNCT__ "DMRestrict"
1070 /*@
1071    DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd()
1072 
1073    Collective if any hooks are
1074 
1075    Input Arguments:
1076 +  fine - finer DM to use as a base
1077 .  restrct - restriction matrix, apply using MatRestrict()
1078 .  inject - injection matrix, also use MatRestrict()
1079 -  coarse - coarer DM to update
1080 
1081    Level: developer
1082 
1083 .seealso: DMCoarsenHookAdd(), MatRestrict()
1084 @*/
1085 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse)
1086 {
1087   PetscErrorCode ierr;
1088   DMCoarsenHookLink link;
1089 
1090   PetscFunctionBegin;
1091   for (link=fine->coarsenhook; link; link=link->next) {
1092     if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);}
1093   }
1094   PetscFunctionReturn(0);
1095 }
1096 
1097 #undef __FUNCT__
1098 #define __FUNCT__ "DMGetCoarsenLevel"
1099 /*@
1100     DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM.
1101 
1102     Not Collective
1103 
1104     Input Parameter:
1105 .   dm - the DM object
1106 
1107     Output Parameter:
1108 .   level - number of coarsenings
1109 
1110     Level: developer
1111 
1112 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1113 
1114 @*/
1115 PetscErrorCode  DMGetCoarsenLevel(DM dm,PetscInt *level)
1116 {
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1119   *level = dm->leveldown;
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "DMRefineHierarchy"
1127 /*@C
1128     DMRefineHierarchy - Refines a DM object, all levels at once
1129 
1130     Collective on DM
1131 
1132     Input Parameter:
1133 +   dm - the DM object
1134 -   nlevels - the number of levels of refinement
1135 
1136     Output Parameter:
1137 .   dmf - the refined DM hierarchy
1138 
1139     Level: developer
1140 
1141 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1142 
1143 @*/
1144 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
1145 {
1146   PetscErrorCode ierr;
1147 
1148   PetscFunctionBegin;
1149   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1150   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1151   if (nlevels == 0) PetscFunctionReturn(0);
1152   if (dm->ops->refinehierarchy) {
1153     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
1154   } else if (dm->ops->refine) {
1155     PetscInt i;
1156 
1157     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
1158     for (i=1; i<nlevels; i++) {
1159       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
1160     }
1161   } else {
1162     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
1163   }
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 #undef __FUNCT__
1168 #define __FUNCT__ "DMCoarsenHierarchy"
1169 /*@C
1170     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
1171 
1172     Collective on DM
1173 
1174     Input Parameter:
1175 +   dm - the DM object
1176 -   nlevels - the number of levels of coarsening
1177 
1178     Output Parameter:
1179 .   dmc - the coarsened DM hierarchy
1180 
1181     Level: developer
1182 
1183 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation()
1184 
1185 @*/
1186 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
1187 {
1188   PetscErrorCode ierr;
1189 
1190   PetscFunctionBegin;
1191   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1192   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1193   if (nlevels == 0) PetscFunctionReturn(0);
1194   PetscValidPointer(dmc,3);
1195   if (dm->ops->coarsenhierarchy) {
1196     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
1197   } else if (dm->ops->coarsen) {
1198     PetscInt i;
1199 
1200     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
1201     for (i=1; i<nlevels; i++) {
1202       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
1203     }
1204   } else {
1205     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
1206   }
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 #undef __FUNCT__
1211 #define __FUNCT__ "DMCreateAggregates"
1212 /*@
1213    DMCreateAggregates - Gets the aggregates that map between
1214    grids associated with two DMs.
1215 
1216    Collective on DM
1217 
1218    Input Parameters:
1219 +  dmc - the coarse grid DM
1220 -  dmf - the fine grid DM
1221 
1222    Output Parameters:
1223 .  rest - the restriction matrix (transpose of the projection matrix)
1224 
1225    Level: intermediate
1226 
1227 .keywords: interpolation, restriction, multigrid
1228 
1229 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation()
1230 @*/
1231 PetscErrorCode  DMCreateAggregates(DM dmc, DM dmf, Mat *rest)
1232 {
1233   PetscErrorCode ierr;
1234 
1235   PetscFunctionBegin;
1236   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
1237   PetscValidHeaderSpecific(dmf,DM_CLASSID,2);
1238   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNCT__
1243 #define __FUNCT__ "DMSetApplicationContextDestroy"
1244 /*@C
1245     DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed
1246 
1247     Not Collective
1248 
1249     Input Parameters:
1250 +   dm - the DM object
1251 -   destroy - the destroy function
1252 
1253     Level: intermediate
1254 
1255 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1256 
1257 @*/
1258 PetscErrorCode  DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**))
1259 {
1260   PetscFunctionBegin;
1261   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1262   dm->ctxdestroy = destroy;
1263   PetscFunctionReturn(0);
1264 }
1265 
1266 #undef __FUNCT__
1267 #define __FUNCT__ "DMSetApplicationContext"
1268 /*@
1269     DMSetApplicationContext - Set a user context into a DM object
1270 
1271     Not Collective
1272 
1273     Input Parameters:
1274 +   dm - the DM object
1275 -   ctx - the user context
1276 
1277     Level: intermediate
1278 
1279 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1280 
1281 @*/
1282 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
1283 {
1284   PetscFunctionBegin;
1285   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1286   dm->ctx = ctx;
1287   PetscFunctionReturn(0);
1288 }
1289 
1290 #undef __FUNCT__
1291 #define __FUNCT__ "DMGetApplicationContext"
1292 /*@
1293     DMGetApplicationContext - Gets a user context from a DM object
1294 
1295     Not Collective
1296 
1297     Input Parameter:
1298 .   dm - the DM object
1299 
1300     Output Parameter:
1301 .   ctx - the user context
1302 
1303     Level: intermediate
1304 
1305 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext()
1306 
1307 @*/
1308 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1309 {
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1312   *(void**)ctx = dm->ctx;
1313   PetscFunctionReturn(0);
1314 }
1315 
1316 #undef __FUNCT__
1317 #define __FUNCT__ "DMSetInitialGuess"
1318 /*@C
1319     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1320 
1321     Logically Collective on DM
1322 
1323     Input Parameter:
1324 +   dm - the DM object to destroy
1325 -   f - the function to compute the initial guess
1326 
1327     Level: intermediate
1328 
1329 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1330 
1331 @*/
1332 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1333 {
1334   PetscFunctionBegin;
1335   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1336   dm->ops->initialguess = f;
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 #undef __FUNCT__
1341 #define __FUNCT__ "DMSetFunction"
1342 /*@C
1343     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1344 
1345     Logically Collective on DM
1346 
1347     Input Parameter:
1348 +   dm - the DM object
1349 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1350 
1351     Level: intermediate
1352 
1353     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1354            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1355 
1356 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1357          DMSetJacobian()
1358 
1359 @*/
1360 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1361 {
1362   PetscFunctionBegin;
1363   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1364   dm->ops->function = f;
1365   if (f) {
1366     dm->ops->functionj = f;
1367   }
1368   PetscFunctionReturn(0);
1369 }
1370 
1371 #undef __FUNCT__
1372 #define __FUNCT__ "DMSetJacobian"
1373 /*@C
1374     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1375 
1376     Logically Collective on DM
1377 
1378     Input Parameter:
1379 +   dm - the DM object to destroy
1380 -   f - the function to compute the matrix entries
1381 
1382     Level: intermediate
1383 
1384 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1385          DMSetFunction()
1386 
1387 @*/
1388 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1389 {
1390   PetscFunctionBegin;
1391   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1392   dm->ops->jacobian = f;
1393   PetscFunctionReturn(0);
1394 }
1395 
1396 #undef __FUNCT__
1397 #define __FUNCT__ "DMSetVariableBounds"
1398 /*@C
1399     DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI.
1400 
1401     Logically Collective on DM
1402 
1403     Input Parameter:
1404 +   dm - the DM object
1405 -   f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set)
1406 
1407     Level: intermediate
1408 
1409 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1410          DMSetJacobian()
1411 
1412 @*/
1413 PetscErrorCode  DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1414 {
1415   PetscFunctionBegin;
1416   dm->ops->computevariablebounds = f;
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 #undef __FUNCT__
1421 #define __FUNCT__ "DMHasVariableBounds"
1422 /*@
1423     DMHasVariableBounds - does the DM object have a variable bounds function?
1424 
1425     Not Collective
1426 
1427     Input Parameter:
1428 .   dm - the DM object to destroy
1429 
1430     Output Parameter:
1431 .   flg - PETSC_TRUE if the variable bounds function exists
1432 
1433     Level: developer
1434 
1435 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1436 
1437 @*/
1438 PetscErrorCode  DMHasVariableBounds(DM dm,PetscBool  *flg)
1439 {
1440   PetscFunctionBegin;
1441   *flg =  (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE;
1442   PetscFunctionReturn(0);
1443 }
1444 
1445 #undef __FUNCT__
1446 #define __FUNCT__ "DMComputeVariableBounds"
1447 /*@C
1448     DMComputeVariableBounds - compute variable bounds used by SNESVI.
1449 
1450     Logically Collective on DM
1451 
1452     Input Parameters:
1453 +   dm - the DM object to destroy
1454 -   x  - current solution at which the bounds are computed
1455 
1456     Output parameters:
1457 +   xl - lower bound
1458 -   xu - upper bound
1459 
1460     Level: intermediate
1461 
1462 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1463          DMSetFunction(), DMSetVariableBounds()
1464 
1465 @*/
1466 PetscErrorCode  DMComputeVariableBounds(DM dm, Vec xl, Vec xu)
1467 {
1468   PetscErrorCode ierr;
1469   PetscFunctionBegin;
1470   PetscValidHeaderSpecific(xl,VEC_CLASSID,2);
1471   PetscValidHeaderSpecific(xu,VEC_CLASSID,2);
1472   if(dm->ops->computevariablebounds) {
1473     ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr);
1474   }
1475   else {
1476     ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr);
1477     ierr = VecSet(xu,SNES_VI_INF);  CHKERRQ(ierr);
1478   }
1479   PetscFunctionReturn(0);
1480 }
1481 
1482 #undef __FUNCT__
1483 #define __FUNCT__ "DMComputeInitialGuess"
1484 /*@
1485     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1486 
1487     Collective on DM
1488 
1489     Input Parameter:
1490 +   dm - the DM object to destroy
1491 -   x - the vector to hold the initial guess values
1492 
1493     Level: developer
1494 
1495 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1496 
1497 @*/
1498 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1499 {
1500   PetscErrorCode ierr;
1501   PetscFunctionBegin;
1502   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1503   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1504   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1505   PetscFunctionReturn(0);
1506 }
1507 
1508 #undef __FUNCT__
1509 #define __FUNCT__ "DMHasInitialGuess"
1510 /*@
1511     DMHasInitialGuess - does the DM object have an initial guess function
1512 
1513     Not Collective
1514 
1515     Input Parameter:
1516 .   dm - the DM object to destroy
1517 
1518     Output Parameter:
1519 .   flg - PETSC_TRUE if function exists
1520 
1521     Level: developer
1522 
1523 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1524 
1525 @*/
1526 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1527 {
1528   PetscFunctionBegin;
1529   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1530   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1531   PetscFunctionReturn(0);
1532 }
1533 
1534 #undef __FUNCT__
1535 #define __FUNCT__ "DMHasFunction"
1536 /*@
1537     DMHasFunction - does the DM object have a function
1538 
1539     Not Collective
1540 
1541     Input Parameter:
1542 .   dm - the DM object to destroy
1543 
1544     Output Parameter:
1545 .   flg - PETSC_TRUE if function exists
1546 
1547     Level: developer
1548 
1549 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1550 
1551 @*/
1552 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1553 {
1554   PetscFunctionBegin;
1555   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1556   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 #undef __FUNCT__
1561 #define __FUNCT__ "DMHasJacobian"
1562 /*@
1563     DMHasJacobian - does the DM object have a matrix function
1564 
1565     Not Collective
1566 
1567     Input Parameter:
1568 .   dm - the DM object to destroy
1569 
1570     Output Parameter:
1571 .   flg - PETSC_TRUE if function exists
1572 
1573     Level: developer
1574 
1575 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1576 
1577 @*/
1578 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1579 {
1580   PetscFunctionBegin;
1581   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1582   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1583   PetscFunctionReturn(0);
1584 }
1585 
1586 #undef  __FUNCT__
1587 #define __FUNCT__ "DMSetVec"
1588 /*@C
1589     DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear.
1590 
1591     Collective on DM
1592 
1593     Input Parameter:
1594 +   dm - the DM object
1595 -   x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems.
1596 
1597     Level: developer
1598 
1599 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1600          DMSetFunction(), DMSetJacobian(), DMSetVariableBounds()
1601 
1602 @*/
1603 PetscErrorCode  DMSetVec(DM dm,Vec x)
1604 {
1605   PetscErrorCode ierr;
1606   PetscFunctionBegin;
1607   if (x) {
1608     if (!dm->x) {
1609       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1610     }
1611     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1612   }
1613   else if(dm->x) {
1614     ierr = VecDestroy(&dm->x);  CHKERRQ(ierr);
1615   }
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 
1620 #undef __FUNCT__
1621 #define __FUNCT__ "DMComputeFunction"
1622 /*@
1623     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1624 
1625     Collective on DM
1626 
1627     Input Parameter:
1628 +   dm - the DM object to destroy
1629 .   x - the location where the function is evaluationed, may be ignored for linear problems
1630 -   b - the vector to hold the right hand side entries
1631 
1632     Level: developer
1633 
1634 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1635          DMSetJacobian()
1636 
1637 @*/
1638 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1639 {
1640   PetscErrorCode ierr;
1641   PetscFunctionBegin;
1642   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1643   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1644   PetscStackPush("DM user function");
1645   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1646   PetscStackPop;
1647   PetscFunctionReturn(0);
1648 }
1649 
1650 
1651 
1652 #undef __FUNCT__
1653 #define __FUNCT__ "DMComputeJacobian"
1654 /*@
1655     DMComputeJacobian - compute the matrix entries for the solver
1656 
1657     Collective on DM
1658 
1659     Input Parameter:
1660 +   dm - the DM object
1661 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1662 .   A - matrix that defines the operator for the linear solve
1663 -   B - the matrix used to construct the preconditioner
1664 
1665     Level: developer
1666 
1667 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1668          DMSetFunction()
1669 
1670 @*/
1671 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1672 {
1673   PetscErrorCode ierr;
1674 
1675   PetscFunctionBegin;
1676   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1677   if (!dm->ops->jacobian) {
1678     ISColoring     coloring;
1679     MatFDColoring  fd;
1680 
1681     ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
1682     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1683     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1684     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1685     ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1686     ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr);
1687 
1688     dm->fd = fd;
1689     dm->ops->jacobian = DMComputeJacobianDefault;
1690 
1691     /* don't know why this is needed */
1692     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1693   }
1694   if (!x) x = dm->x;
1695   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1696 
1697   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1698   if (x) {
1699     if (!dm->x) {
1700       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1701     }
1702     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1703   }
1704   if (A != B) {
1705     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1706     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1707   }
1708   PetscFunctionReturn(0);
1709 }
1710 
1711 
1712 PetscFList DMList                       = PETSC_NULL;
1713 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1714 
1715 #undef __FUNCT__
1716 #define __FUNCT__ "DMSetType"
1717 /*@C
1718   DMSetType - Builds a DM, for a particular DM implementation.
1719 
1720   Collective on DM
1721 
1722   Input Parameters:
1723 + dm     - The DM object
1724 - method - The name of the DM type
1725 
1726   Options Database Key:
1727 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1728 
1729   Notes:
1730   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1731 
1732   Level: intermediate
1733 
1734 .keywords: DM, set, type
1735 .seealso: DMGetType(), DMCreate()
1736 @*/
1737 PetscErrorCode  DMSetType(DM dm, const DMType method)
1738 {
1739   PetscErrorCode (*r)(DM);
1740   PetscBool      match;
1741   PetscErrorCode ierr;
1742 
1743   PetscFunctionBegin;
1744   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1745   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1746   if (match) PetscFunctionReturn(0);
1747 
1748   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1749   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1750   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1751 
1752   if (dm->ops->destroy) {
1753     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1754     dm->ops->destroy = PETSC_NULL;
1755   }
1756   ierr = (*r)(dm);CHKERRQ(ierr);
1757   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1758   PetscFunctionReturn(0);
1759 }
1760 
1761 #undef __FUNCT__
1762 #define __FUNCT__ "DMGetType"
1763 /*@C
1764   DMGetType - Gets the DM type name (as a string) from the DM.
1765 
1766   Not Collective
1767 
1768   Input Parameter:
1769 . dm  - The DM
1770 
1771   Output Parameter:
1772 . type - The DM type name
1773 
1774   Level: intermediate
1775 
1776 .keywords: DM, get, type, name
1777 .seealso: DMSetType(), DMCreate()
1778 @*/
1779 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1780 {
1781   PetscErrorCode ierr;
1782 
1783   PetscFunctionBegin;
1784   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1785   PetscValidCharPointer(type,2);
1786   if (!DMRegisterAllCalled) {
1787     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1788   }
1789   *type = ((PetscObject)dm)->type_name;
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 #undef __FUNCT__
1794 #define __FUNCT__ "DMConvert"
1795 /*@C
1796   DMConvert - Converts a DM to another DM, either of the same or different type.
1797 
1798   Collective on DM
1799 
1800   Input Parameters:
1801 + dm - the DM
1802 - newtype - new DM type (use "same" for the same type)
1803 
1804   Output Parameter:
1805 . M - pointer to new DM
1806 
1807   Notes:
1808   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1809   the MPI communicator of the generated DM is always the same as the communicator
1810   of the input DM.
1811 
1812   Level: intermediate
1813 
1814 .seealso: DMCreate()
1815 @*/
1816 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1817 {
1818   DM             B;
1819   char           convname[256];
1820   PetscBool      sametype, issame;
1821   PetscErrorCode ierr;
1822 
1823   PetscFunctionBegin;
1824   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1825   PetscValidType(dm,1);
1826   PetscValidPointer(M,3);
1827   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1828   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1829   {
1830     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1831 
1832     /*
1833        Order of precedence:
1834        1) See if a specialized converter is known to the current DM.
1835        2) See if a specialized converter is known to the desired DM class.
1836        3) See if a good general converter is registered for the desired class
1837        4) See if a good general converter is known for the current matrix.
1838        5) Use a really basic converter.
1839     */
1840 
1841     /* 1) See if a specialized converter is known to the current DM and the desired class */
1842     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1843     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1844     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1845     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1846     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1847     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1848     if (conv) goto foundconv;
1849 
1850     /* 2)  See if a specialized converter is known to the desired DM class. */
1851     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1852     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1853     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1854     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1855     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1856     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1857     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1858     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1859     if (conv) {
1860       ierr = DMDestroy(&B);CHKERRQ(ierr);
1861       goto foundconv;
1862     }
1863 
1864 #if 0
1865     /* 3) See if a good general converter is registered for the desired class */
1866     conv = B->ops->convertfrom;
1867     ierr = DMDestroy(&B);CHKERRQ(ierr);
1868     if (conv) goto foundconv;
1869 
1870     /* 4) See if a good general converter is known for the current matrix */
1871     if (dm->ops->convert) {
1872       conv = dm->ops->convert;
1873     }
1874     if (conv) goto foundconv;
1875 #endif
1876 
1877     /* 5) Use a really basic converter. */
1878     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
1879 
1880     foundconv:
1881     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1882     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
1883     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1884   }
1885   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
1886   PetscFunctionReturn(0);
1887 }
1888 
1889 /*--------------------------------------------------------------------------------------------------------------------*/
1890 
1891 #undef __FUNCT__
1892 #define __FUNCT__ "DMRegister"
1893 /*@C
1894   DMRegister - See DMRegisterDynamic()
1895 
1896   Level: advanced
1897 @*/
1898 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1899 {
1900   char fullname[PETSC_MAX_PATH_LEN];
1901   PetscErrorCode ierr;
1902 
1903   PetscFunctionBegin;
1904   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1905   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1906   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1907   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1908   PetscFunctionReturn(0);
1909 }
1910 
1911 
1912 /*--------------------------------------------------------------------------------------------------------------------*/
1913 #undef __FUNCT__
1914 #define __FUNCT__ "DMRegisterDestroy"
1915 /*@C
1916    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1917 
1918    Not Collective
1919 
1920    Level: advanced
1921 
1922 .keywords: DM, register, destroy
1923 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1924 @*/
1925 PetscErrorCode  DMRegisterDestroy(void)
1926 {
1927   PetscErrorCode ierr;
1928 
1929   PetscFunctionBegin;
1930   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1931   DMRegisterAllCalled = PETSC_FALSE;
1932   PetscFunctionReturn(0);
1933 }
1934 
1935 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1936 #include <mex.h>
1937 
1938 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1939 
1940 #undef __FUNCT__
1941 #define __FUNCT__ "DMComputeFunction_Matlab"
1942 /*
1943    DMComputeFunction_Matlab - Calls the function that has been set with
1944                          DMSetFunctionMatlab().
1945 
1946    For linear problems x is null
1947 
1948 .seealso: DMSetFunction(), DMGetFunction()
1949 */
1950 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1951 {
1952   PetscErrorCode    ierr;
1953   DMMatlabContext   *sctx;
1954   int               nlhs = 1,nrhs = 4;
1955   mxArray	    *plhs[1],*prhs[4];
1956   long long int     lx = 0,ly = 0,ls = 0;
1957 
1958   PetscFunctionBegin;
1959   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1960   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1961   PetscCheckSameComm(dm,1,y,3);
1962 
1963   /* call Matlab function in ctx with arguments x and y */
1964   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1965   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1966   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1967   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1968   prhs[0] =  mxCreateDoubleScalar((double)ls);
1969   prhs[1] =  mxCreateDoubleScalar((double)lx);
1970   prhs[2] =  mxCreateDoubleScalar((double)ly);
1971   prhs[3] =  mxCreateString(sctx->funcname);
1972   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1973   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1974   mxDestroyArray(prhs[0]);
1975   mxDestroyArray(prhs[1]);
1976   mxDestroyArray(prhs[2]);
1977   mxDestroyArray(prhs[3]);
1978   mxDestroyArray(plhs[0]);
1979   PetscFunctionReturn(0);
1980 }
1981 
1982 
1983 #undef __FUNCT__
1984 #define __FUNCT__ "DMSetFunctionMatlab"
1985 /*
1986    DMSetFunctionMatlab - Sets the function evaluation routine
1987 
1988 */
1989 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
1990 {
1991   PetscErrorCode    ierr;
1992   DMMatlabContext   *sctx;
1993 
1994   PetscFunctionBegin;
1995   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1996   /* currently sctx is memory bleed */
1997   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1998   if (!sctx) {
1999     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2000   }
2001   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2002   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2003   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
2004   PetscFunctionReturn(0);
2005 }
2006 
2007 #undef __FUNCT__
2008 #define __FUNCT__ "DMComputeJacobian_Matlab"
2009 /*
2010    DMComputeJacobian_Matlab - Calls the function that has been set with
2011                          DMSetJacobianMatlab().
2012 
2013    For linear problems x is null
2014 
2015 .seealso: DMSetFunction(), DMGetFunction()
2016 */
2017 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
2018 {
2019   PetscErrorCode    ierr;
2020   DMMatlabContext   *sctx;
2021   int               nlhs = 2,nrhs = 5;
2022   mxArray	    *plhs[2],*prhs[5];
2023   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
2024 
2025   PetscFunctionBegin;
2026   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2027   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
2028 
2029   /* call MATLAB function in ctx with arguments x, A, and B */
2030   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2031   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
2032   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2033   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
2034   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
2035   prhs[0] =  mxCreateDoubleScalar((double)ls);
2036   prhs[1] =  mxCreateDoubleScalar((double)lx);
2037   prhs[2] =  mxCreateDoubleScalar((double)lA);
2038   prhs[3] =  mxCreateDoubleScalar((double)lB);
2039   prhs[4] =  mxCreateString(sctx->jacname);
2040   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
2041   *str    =  (MatStructure) mxGetScalar(plhs[0]);
2042   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2043   mxDestroyArray(prhs[0]);
2044   mxDestroyArray(prhs[1]);
2045   mxDestroyArray(prhs[2]);
2046   mxDestroyArray(prhs[3]);
2047   mxDestroyArray(prhs[4]);
2048   mxDestroyArray(plhs[0]);
2049   mxDestroyArray(plhs[1]);
2050   PetscFunctionReturn(0);
2051 }
2052 
2053 
2054 #undef __FUNCT__
2055 #define __FUNCT__ "DMSetJacobianMatlab"
2056 /*
2057    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
2058 
2059 */
2060 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
2061 {
2062   PetscErrorCode    ierr;
2063   DMMatlabContext   *sctx;
2064 
2065   PetscFunctionBegin;
2066   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2067   /* currently sctx is memory bleed */
2068   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
2069   if (!sctx) {
2070     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
2071   }
2072   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
2073   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
2074   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
2075   PetscFunctionReturn(0);
2076 }
2077 #endif
2078 
2079 #undef __FUNCT__
2080 #define __FUNCT__ "DMLoad"
2081 /*@C
2082   DMLoad - Loads a DM that has been stored in binary or HDF5 format
2083   with DMView().
2084 
2085   Collective on PetscViewer
2086 
2087   Input Parameters:
2088 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
2089            some related function before a call to DMLoad().
2090 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
2091            HDF5 file viewer, obtained from PetscViewerHDF5Open()
2092 
2093    Level: intermediate
2094 
2095   Notes:
2096   Defaults to the DM DA.
2097 
2098   Notes for advanced users:
2099   Most users should not need to know the details of the binary storage
2100   format, since DMLoad() and DMView() completely hide these details.
2101   But for anyone who's interested, the standard binary matrix storage
2102   format is
2103 .vb
2104      has not yet been determined
2105 .ve
2106 
2107    In addition, PETSc automatically does the byte swapping for
2108 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2109 linux, Windows and the paragon; thus if you write your own binary
2110 read/write routines you have to swap the bytes; see PetscBinaryRead()
2111 and PetscBinaryWrite() to see how this may be done.
2112 
2113   Concepts: vector^loading from file
2114 
2115 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
2116 @*/
2117 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
2118 {
2119   PetscErrorCode ierr;
2120 
2121   PetscFunctionBegin;
2122   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
2123   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2124 
2125   if (!((PetscObject)newdm)->type_name) {
2126     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
2127   }
2128   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
2129   PetscFunctionReturn(0);
2130 }
2131 
2132 /******************************** FEM Support **********************************/
2133 
2134 #undef __FUNCT__
2135 #define __FUNCT__ "DMPrintCellVector"
2136 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) {
2137   PetscInt       f;
2138   PetscErrorCode ierr;
2139 
2140   PetscFunctionBegin;
2141   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2142   for(f = 0; f < len; ++f) {
2143     ierr = PetscPrintf(PETSC_COMM_SELF, "  | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr);
2144   }
2145   PetscFunctionReturn(0);
2146 }
2147 
2148 #undef __FUNCT__
2149 #define __FUNCT__ "DMPrintCellMatrix"
2150 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) {
2151   PetscInt       f, g;
2152   PetscErrorCode ierr;
2153 
2154   PetscFunctionBegin;
2155   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr);
2156   for(f = 0; f < rows; ++f) {
2157     ierr = PetscPrintf(PETSC_COMM_SELF, "  |");CHKERRQ(ierr);
2158     for(g = 0; g < cols; ++g) {
2159       ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr);
2160     }
2161     ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr);
2162   }
2163   PetscFunctionReturn(0);
2164 }
2165