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