xref: /petsc/src/dm/interface/dm.c (revision f3e3d7df08ad49fe665df0df50a5f980ce6310e9)
1 
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 vector object. The type can then be set with DMetType().
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 = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
37 #endif
38 
39   ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", comm, DMDestroy, DMView);CHKERRQ(ierr);
40   ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr);
41 
42   v->ltogmap      = PETSC_NULL;
43   v->ltogmapb     = PETSC_NULL;
44   v->bs           = 1;
45 
46   *dm = v;
47   PetscFunctionReturn(0);
48 }
49 
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "DMSetVecType"
53 /*@C
54        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
55 
56    Logically Collective on DMDA
57 
58    Input Parameter:
59 +  da - initial distributed array
60 .  ctype - the vector type, currently either VECSTANDARD or VECCUSP
61 
62    Options Database:
63 .   -da_vec_type ctype
64 
65    Level: intermediate
66 
67 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType
68 @*/
69 PetscErrorCode  DMSetVecType(DM da,const VecType ctype)
70 {
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   PetscValidHeaderSpecific(da,DM_CLASSID,1);
75   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
76   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "DMSetOptionsPrefix"
82 /*@C
83    DMSetOptionsPrefix - Sets the prefix used for searching for all
84    DMDA options in the database.
85 
86    Logically Collective on DMDA
87 
88    Input Parameter:
89 +  da - the DMDA context
90 -  prefix - the prefix to prepend to all option names
91 
92    Notes:
93    A hyphen (-) must NOT be given at the beginning of the prefix name.
94    The first character of all runtime options is AUTOMATICALLY the hyphen.
95 
96    Level: advanced
97 
98 .keywords: DMDA, set, options, prefix, database
99 
100 .seealso: DMSetFromOptions()
101 @*/
102 PetscErrorCode  DMSetOptionsPrefix(DM dm,const char prefix[])
103 {
104   PetscErrorCode ierr;
105 
106   PetscFunctionBegin;
107   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
108   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "DMDestroy"
114 /*@
115     DMDestroy - Destroys a vector packer or DMDA.
116 
117     Collective on DM
118 
119     Input Parameter:
120 .   dm - the DM object to destroy
121 
122     Level: developer
123 
124 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
125 
126 @*/
127 PetscErrorCode  DMDestroy(DM *dm)
128 {
129   PetscInt       i, cnt = 0;
130   PetscErrorCode ierr;
131 
132   PetscFunctionBegin;
133   if (!*dm) PetscFunctionReturn(0);
134   PetscValidHeaderSpecific((*dm),DM_CLASSID,1);
135 
136   /* count all the circular references of DM and its contained Vecs */
137   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
138     if ((*dm)->localin[i])  {cnt++;}
139     if ((*dm)->globalin[i]) {cnt++;}
140   }
141   if ((*dm)->x) {
142     PetscObject obj;
143     ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr);
144     if (obj == (PetscObject)*dm) cnt++;
145   }
146 
147   if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);}
148   /*
149      Need this test because the dm references the vectors that
150      reference the dm, so destroying the dm calls destroy on the
151      vectors that cause another destroy on the dm
152   */
153   if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0);
154   ((PetscObject) (*dm))->refct = 0;
155   for (i=0; i<DM_MAX_WORK_VECTORS; i++) {
156     if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()");
157     ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr);
158   }
159   ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr);
160   ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr);
161   ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr);
162   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr);
163   ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr);
164   ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr);
165 
166   /* if memory was published with AMS then destroy it */
167   ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr);
168 
169   ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr);
170   ierr = PetscFree((*dm)->data);CHKERRQ(ierr);
171   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "DMSetUp"
177 /*@
178     DMSetUp - sets up the data structures inside a DM object
179 
180     Collective on DM
181 
182     Input Parameter:
183 .   dm - the DM object to setup
184 
185     Level: developer
186 
187 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
188 
189 @*/
190 PetscErrorCode  DMSetUp(DM dm)
191 {
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   if (dm->setupcalled) PetscFunctionReturn(0);
196   if (dm->ops->setup) {
197     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
198   }
199   dm->setupcalled = PETSC_TRUE;
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "DMSetFromOptions"
205 /*@
206     DMSetFromOptions - sets parameters in a DM from the options database
207 
208     Collective on DM
209 
210     Input Parameter:
211 .   dm - the DM object to set options for
212 
213     Options Database:
214 .   -dm_preallocate_only: Only preallocate the matrix for DMGetMatrix(), but do not fill it with zeros
215 
216     Level: developer
217 
218 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
219 
220 @*/
221 PetscErrorCode  DMSetFromOptions(DM dm)
222 {
223   PetscBool      flg1 = PETSC_FALSE;
224   PetscErrorCode ierr;
225 
226   PetscFunctionBegin;
227   ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_preallocate_only", &dm->prealloc_only, PETSC_NULL);CHKERRQ(ierr);
228   if (dm->ops->setfromoptions) {
229     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
230   }
231   ierr = PetscOptionsBegin(((PetscObject)dm)->comm,((PetscObject)dm)->prefix,"DM Options","DM");CHKERRQ(ierr);
232     ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr);
233   ierr = PetscOptionsEnd();CHKERRQ(ierr);
234   if (flg1) {
235     ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
236   }
237   PetscFunctionReturn(0);
238 }
239 
240 #undef __FUNCT__
241 #define __FUNCT__ "DMView"
242 /*@C
243     DMView - Views a vector packer or DMDA.
244 
245     Collective on DM
246 
247     Input Parameter:
248 +   dm - the DM object to view
249 -   v - the viewer
250 
251     Level: developer
252 
253 .seealso DMDestroy(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
254 
255 @*/
256 PetscErrorCode  DMView(DM dm,PetscViewer v)
257 {
258   PetscErrorCode ierr;
259 
260   PetscFunctionBegin;
261  if (!v) {
262     ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr);
263   }
264   if (dm->ops->view) {
265     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
266   }
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "DMCreateGlobalVector"
272 /*@
273     DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object
274 
275     Collective on DM
276 
277     Input Parameter:
278 .   dm - the DM object
279 
280     Output Parameter:
281 .   vec - the global vector
282 
283     Level: developer
284 
285 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
286 
287 @*/
288 PetscErrorCode  DMCreateGlobalVector(DM dm,Vec *vec)
289 {
290   PetscErrorCode ierr;
291 
292   PetscFunctionBegin;
293   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "DMCreateLocalVector"
299 /*@
300     DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object
301 
302     Not Collective
303 
304     Input Parameter:
305 .   dm - the DM object
306 
307     Output Parameter:
308 .   vec - the local vector
309 
310     Level: developer
311 
312 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
313 
314 @*/
315 PetscErrorCode  DMCreateLocalVector(DM dm,Vec *vec)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "DMGetLocalToGlobalMapping"
326 /*@
327    DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM.
328 
329    Collective on DM
330 
331    Input Parameter:
332 .  dm - the DM that provides the mapping
333 
334    Output Parameter:
335 .  ltog - the mapping
336 
337    Level: intermediate
338 
339    Notes:
340    This mapping can then be used by VecSetLocalToGlobalMapping() or
341    MatSetLocalToGlobalMapping().
342 
343 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock()
344 @*/
345 PetscErrorCode  DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog)
346 {
347   PetscErrorCode ierr;
348 
349   PetscFunctionBegin;
350   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
351   PetscValidPointer(ltog,2);
352   if (!dm->ltogmap) {
353     if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping");
354     ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr);
355   }
356   *ltog = dm->ltogmap;
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock"
362 /*@
363    DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM.
364 
365    Collective on DM
366 
367    Input Parameter:
368 .  da - the distributed array that provides the mapping
369 
370    Output Parameter:
371 .  ltog - the block mapping
372 
373    Level: intermediate
374 
375    Notes:
376    This mapping can then be used by VecSetLocalToGlobalMappingBlock() or
377    MatSetLocalToGlobalMappingBlock().
378 
379 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize()
380 @*/
381 PetscErrorCode  DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
387   PetscValidPointer(ltog,2);
388   if (!dm->ltogmapb) {
389     PetscInt bs;
390     ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr);
391     if (bs > 1) {
392       if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock");
393       ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr);
394     } else {
395       ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr);
396       ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr);
397     }
398   }
399   *ltog = dm->ltogmapb;
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "DMGetBlockSize"
405 /*@
406    DMGetBlockSize - Gets the inherent block size associated with a DM
407 
408    Not Collective
409 
410    Input Parameter:
411 .  dm - the DM with block structure
412 
413    Output Parameter:
414 .  bs - the block size, 1 implies no exploitable block structure
415 
416    Level: intermediate
417 
418 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock()
419 @*/
420 PetscErrorCode  DMGetBlockSize(DM dm,PetscInt *bs)
421 {
422   PetscFunctionBegin;
423   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
424   PetscValidPointer(bs,2);
425   if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet");
426   *bs = dm->bs;
427   PetscFunctionReturn(0);
428 }
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "DMGetInterpolation"
432 /*@
433     DMGetInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects
434 
435     Collective on DM
436 
437     Input Parameter:
438 +   dm1 - the DM object
439 -   dm2 - the second, finer DM object
440 
441     Output Parameter:
442 +  mat - the interpolation
443 -  vec - the scaling (optional)
444 
445     Level: developer
446 
447     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
448         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation.
449 
450         For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors
451         EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic.
452 
453 
454 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMRefine(), DMCoarsen()
455 
456 @*/
457 PetscErrorCode  DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
458 {
459   PetscErrorCode ierr;
460 
461   PetscFunctionBegin;
462   ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 
466 #undef __FUNCT__
467 #define __FUNCT__ "DMGetInjection"
468 /*@
469     DMGetInjection - Gets injection matrix between two DMDA or DMComposite objects
470 
471     Collective on DM
472 
473     Input Parameter:
474 +   dm1 - the DM object
475 -   dm2 - the second, finer DM object
476 
477     Output Parameter:
478 .   ctx - the injection
479 
480     Level: developer
481 
482    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
483         DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection.
484 
485 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation()
486 
487 @*/
488 PetscErrorCode  DMGetInjection(DM dm1,DM dm2,VecScatter *ctx)
489 {
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
494   PetscFunctionReturn(0);
495 }
496 
497 #undef __FUNCT__
498 #define __FUNCT__ "DMGetColoring"
499 /*@C
500     DMGetColoring - Gets coloring for a DMDA or DMComposite
501 
502     Collective on DM
503 
504     Input Parameter:
505 +   dm - the DM object
506 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
507 -   matype - either MATAIJ or MATBAIJ
508 
509     Output Parameter:
510 .   coloring - the coloring
511 
512     Level: developer
513 
514 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
515 
516 @*/
517 PetscErrorCode  DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
518 {
519   PetscErrorCode ierr;
520 
521   PetscFunctionBegin;
522   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
523   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "DMGetMatrix"
529 /*@C
530     DMGetMatrix - Gets empty Jacobian for a DMDA or DMComposite
531 
532     Collective on DM
533 
534     Input Parameter:
535 +   dm - the DM object
536 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
537             any type which inherits from one of these (such as MATAIJ)
538 
539     Output Parameter:
540 .   mat - the empty Jacobian
541 
542     Level: developer
543 
544     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
545        do not need to do it yourself.
546 
547        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
548        the nonzero pattern call DMDASetMatPreallocateOnly()
549 
550        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
551        internally by PETSc.
552 
553        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
554        the indices for the global numbering for DMDAs which is complicated.
555 
556 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
557 
558 @*/
559 PetscErrorCode  DMGetMatrix(DM dm,const MatType mtype,Mat *mat)
560 {
561   PetscErrorCode ierr;
562   char           ttype[256];
563   PetscBool      flg;
564 
565   PetscFunctionBegin;
566 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
567   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
568 #endif
569   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
570   PetscValidPointer(mat,3);
571   ierr = PetscStrncpy(ttype,mtype,sizeof(ttype));CHKERRQ(ierr);
572   ttype[sizeof(ttype)-1] = 0;
573   ierr = PetscOptionsBegin(((PetscObject)dm)->comm,((PetscObject)dm)->prefix,"DM options","Mat");CHKERRQ(ierr);
574   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,ttype,ttype,sizeof(ttype),&flg);CHKERRQ(ierr);
575   ierr = PetscOptionsEnd();
576   if (flg || mtype) {
577     ierr = (*dm->ops->getmatrix)(dm,ttype,mat);CHKERRQ(ierr);
578   } else {                      /* Let the implementation decide */
579     ierr = (*dm->ops->getmatrix)(dm,PETSC_NULL,mat);CHKERRQ(ierr);
580   }
581   PetscFunctionReturn(0);
582 }
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "DMSetMatrixPreallocateOnly"
586 /*@
587   DMSetMatrixPreallocateOnly - When DMGetMatrix() is called the matrix will be properly
588     preallocated but the nonzero structure and zero values will not be set.
589 
590   Logically Collective on DMDA
591 
592   Input Parameter:
593 + dm - the DM
594 - only - PETSC_TRUE if only want preallocation
595 
596   Level: developer
597 .seealso DMGetMatrix()
598 @*/
599 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only)
600 {
601   PetscFunctionBegin;
602   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
603   dm->prealloc_only = only;
604   PetscFunctionReturn(0);
605 }
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "DMRefine"
609 /*@
610     DMRefine - Refines a DM object
611 
612     Collective on DM
613 
614     Input Parameter:
615 +   dm - the DM object
616 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
617 
618     Output Parameter:
619 .   dmf - the refined DM
620 
621     Level: developer
622 
623 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
624 
625 @*/
626 PetscErrorCode  DMRefine(DM dm,MPI_Comm comm,DM *dmf)
627 {
628   PetscErrorCode ierr;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
632   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
633   (*dmf)->ops->initialguess = dm->ops->initialguess;
634   (*dmf)->ops->function     = dm->ops->function;
635   (*dmf)->ops->functionj    = dm->ops->functionj;
636   if (dm->ops->jacobian != DMComputeJacobianDefault) {
637     (*dmf)->ops->jacobian     = dm->ops->jacobian;
638   }
639   (*dmf)->ctx     = dm->ctx;
640   (*dmf)->levelup = dm->levelup + 1;
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "DMGetRefineLevel"
646 /*@
647     DMGetRefineLevel - Get's the number of refinements that have generated this DM.
648 
649     Not Collective
650 
651     Input Parameter:
652 .   dm - the DM object
653 
654     Output Parameter:
655 .   level - number of refinements
656 
657     Level: developer
658 
659 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
660 
661 @*/
662 PetscErrorCode  DMGetRefineLevel(DM dm,PetscInt *level)
663 {
664   PetscFunctionBegin;
665   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
666   *level = dm->levelup;
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "DMGlobalToLocalBegin"
672 /*@
673     DMGlobalToLocalBegin - Begins updating local vectors from global vector
674 
675     Neighbor-wise Collective on DM
676 
677     Input Parameters:
678 +   dm - the DM object
679 .   g - the global vector
680 .   mode - INSERT_VALUES or ADD_VALUES
681 -   l - the local vector
682 
683 
684     Level: beginner
685 
686 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
687 
688 @*/
689 PetscErrorCode  DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
690 {
691   PetscErrorCode ierr;
692 
693   PetscFunctionBegin;
694   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
695   PetscFunctionReturn(0);
696 }
697 
698 #undef __FUNCT__
699 #define __FUNCT__ "DMGlobalToLocalEnd"
700 /*@
701     DMGlobalToLocalEnd - Ends updating local vectors from global vector
702 
703     Neighbor-wise Collective on DM
704 
705     Input Parameters:
706 +   dm - the DM object
707 .   g - the global vector
708 .   mode - INSERT_VALUES or ADD_VALUES
709 -   l - the local vector
710 
711 
712     Level: beginner
713 
714 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
715 
716 @*/
717 PetscErrorCode  DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
718 {
719   PetscErrorCode ierr;
720 
721   PetscFunctionBegin;
722   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNCT__
727 #define __FUNCT__ "DMLocalToGlobalBegin"
728 /*@
729     DMLocalToGlobalBegin - updates global vectors from local vectors
730 
731     Neighbor-wise Collective on DM
732 
733     Input Parameters:
734 +   dm - the DM object
735 .   l - the local vector
736 .   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
737            base point.
738 - - the global vector
739 
740     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
741            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
742            global array to the final global array with VecAXPY().
743 
744     Level: beginner
745 
746 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
747 
748 @*/
749 PetscErrorCode  DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g)
750 {
751   PetscErrorCode ierr;
752 
753   PetscFunctionBegin;
754   ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr);
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "DMLocalToGlobalEnd"
760 /*@
761     DMLocalToGlobalEnd - updates global vectors from local vectors
762 
763     Neighbor-wise Collective on DM
764 
765     Input Parameters:
766 +   dm - the DM object
767 .   l - the local vector
768 .   mode - INSERT_VALUES or ADD_VALUES
769 -   g - the global vector
770 
771 
772     Level: beginner
773 
774 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
775 
776 @*/
777 PetscErrorCode  DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g)
778 {
779   PetscErrorCode ierr;
780 
781   PetscFunctionBegin;
782   ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr);
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "DMComputeJacobianDefault"
788 /*@
789     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
790 
791     Collective on DM
792 
793     Input Parameter:
794 +   dm - the DM object
795 .   x - location to compute Jacobian at; may be ignored for linear problems
796 .   A - matrix that defines the operator for the linear solve
797 -   B - the matrix used to construct the preconditioner
798 
799     Level: developer
800 
801 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
802          DMSetFunction()
803 
804 @*/
805 PetscErrorCode  DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
806 {
807   PetscErrorCode ierr;
808   PetscFunctionBegin;
809   *stflag = SAME_NONZERO_PATTERN;
810   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
811   if (A != B) {
812     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
813     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
814   }
815   PetscFunctionReturn(0);
816 }
817 
818 #undef __FUNCT__
819 #define __FUNCT__ "DMCoarsen"
820 /*@
821     DMCoarsen - Coarsens a DM object
822 
823     Collective on DM
824 
825     Input Parameter:
826 +   dm - the DM object
827 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
828 
829     Output Parameter:
830 .   dmc - the coarsened DM
831 
832     Level: developer
833 
834 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
835 
836 @*/
837 PetscErrorCode  DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
838 {
839   PetscErrorCode ierr;
840 
841   PetscFunctionBegin;
842   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
843   (*dmc)->ops->initialguess = dm->ops->initialguess;
844   (*dmc)->ops->function     = dm->ops->function;
845   (*dmc)->ops->functionj    = dm->ops->functionj;
846   if (dm->ops->jacobian != DMComputeJacobianDefault) {
847     (*dmc)->ops->jacobian     = dm->ops->jacobian;
848   }
849   (*dmc)->ctx       = dm->ctx;
850   (*dmc)->leveldown = dm->leveldown + 1;
851   PetscFunctionReturn(0);
852 }
853 
854 #undef __FUNCT__
855 #define __FUNCT__ "DMRefineHierarchy"
856 /*@C
857     DMRefineHierarchy - Refines a DM object, all levels at once
858 
859     Collective on DM
860 
861     Input Parameter:
862 +   dm - the DM object
863 -   nlevels - the number of levels of refinement
864 
865     Output Parameter:
866 .   dmf - the refined DM hierarchy
867 
868     Level: developer
869 
870 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
871 
872 @*/
873 PetscErrorCode  DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
874 {
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
879   if (nlevels == 0) PetscFunctionReturn(0);
880   if (dm->ops->refinehierarchy) {
881     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
882   } else if (dm->ops->refine) {
883     PetscInt i;
884 
885     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
886     for (i=1; i<nlevels; i++) {
887       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
888     }
889   } else {
890     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
891   }
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "DMCoarsenHierarchy"
897 /*@C
898     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
899 
900     Collective on DM
901 
902     Input Parameter:
903 +   dm - the DM object
904 -   nlevels - the number of levels of coarsening
905 
906     Output Parameter:
907 .   dmc - the coarsened DM hierarchy
908 
909     Level: developer
910 
911 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
912 
913 @*/
914 PetscErrorCode  DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
915 {
916   PetscErrorCode ierr;
917 
918   PetscFunctionBegin;
919   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
920   if (nlevels == 0) PetscFunctionReturn(0);
921   PetscValidPointer(dmc,3);
922   if (dm->ops->coarsenhierarchy) {
923     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
924   } else if (dm->ops->coarsen) {
925     PetscInt i;
926 
927     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
928     for (i=1; i<nlevels; i++) {
929       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
930     }
931   } else {
932     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
933   }
934   PetscFunctionReturn(0);
935 }
936 
937 #undef __FUNCT__
938 #define __FUNCT__ "DMGetAggregates"
939 /*@
940    DMGetAggregates - Gets the aggregates that map between
941    grids associated with two DMs.
942 
943    Collective on DM
944 
945    Input Parameters:
946 +  dmc - the coarse grid DM
947 -  dmf - the fine grid DM
948 
949    Output Parameters:
950 .  rest - the restriction matrix (transpose of the projection matrix)
951 
952    Level: intermediate
953 
954 .keywords: interpolation, restriction, multigrid
955 
956 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
957 @*/
958 PetscErrorCode  DMGetAggregates(DM dmc, DM dmf, Mat *rest)
959 {
960   PetscErrorCode ierr;
961 
962   PetscFunctionBegin;
963   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
964   PetscFunctionReturn(0);
965 }
966 
967 #undef __FUNCT__
968 #define __FUNCT__ "DMSetApplicationContext"
969 /*@
970     DMSetApplicationContext - Set a user context into a DM object
971 
972     Not Collective
973 
974     Input Parameters:
975 +   dm - the DM object
976 -   ctx - the user context
977 
978     Level: intermediate
979 
980 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext()
981 
982 @*/
983 PetscErrorCode  DMSetApplicationContext(DM dm,void *ctx)
984 {
985   PetscFunctionBegin;
986   dm->ctx = ctx;
987   PetscFunctionReturn(0);
988 }
989 
990 #undef __FUNCT__
991 #define __FUNCT__ "DMGetApplicationContext"
992 /*@
993     DMGetApplicationContext - Gets a user context from a DM object
994 
995     Not Collective
996 
997     Input Parameter:
998 .   dm - the DM object
999 
1000     Output Parameter:
1001 .   ctx - the user context
1002 
1003     Level: intermediate
1004 
1005 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext()
1006 
1007 @*/
1008 PetscErrorCode  DMGetApplicationContext(DM dm,void *ctx)
1009 {
1010   PetscFunctionBegin;
1011   *(void**)ctx = dm->ctx;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "DMSetInitialGuess"
1017 /*@
1018     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
1019 
1020     Logically Collective on DM
1021 
1022     Input Parameter:
1023 +   dm - the DM object to destroy
1024 -   f - the function to compute the initial guess
1025 
1026     Level: intermediate
1027 
1028 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1029 
1030 @*/
1031 PetscErrorCode  DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
1032 {
1033   PetscFunctionBegin;
1034   dm->ops->initialguess = f;
1035   PetscFunctionReturn(0);
1036 }
1037 
1038 #undef __FUNCT__
1039 #define __FUNCT__ "DMSetFunction"
1040 /*@
1041     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
1042 
1043     Logically Collective on DM
1044 
1045     Input Parameter:
1046 +   dm - the DM object
1047 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
1048 
1049     Level: intermediate
1050 
1051     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
1052            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
1053 
1054 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1055          DMSetJacobian()
1056 
1057 @*/
1058 PetscErrorCode  DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
1059 {
1060   PetscFunctionBegin;
1061   dm->ops->function = f;
1062   if (f) {
1063     dm->ops->functionj = f;
1064   }
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 #undef __FUNCT__
1069 #define __FUNCT__ "DMSetJacobian"
1070 /*@
1071     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
1072 
1073     Logically Collective on DM
1074 
1075     Input Parameter:
1076 +   dm - the DM object to destroy
1077 -   f - the function to compute the matrix entries
1078 
1079     Level: intermediate
1080 
1081 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1082          DMSetFunction()
1083 
1084 @*/
1085 PetscErrorCode  DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
1086 {
1087   PetscFunctionBegin;
1088   dm->ops->jacobian = f;
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "DMComputeInitialGuess"
1094 /*@
1095     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
1096 
1097     Collective on DM
1098 
1099     Input Parameter:
1100 +   dm - the DM object to destroy
1101 -   x - the vector to hold the initial guess values
1102 
1103     Level: developer
1104 
1105 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat()
1106 
1107 @*/
1108 PetscErrorCode  DMComputeInitialGuess(DM dm,Vec x)
1109 {
1110   PetscErrorCode ierr;
1111   PetscFunctionBegin;
1112   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
1113   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 #undef __FUNCT__
1118 #define __FUNCT__ "DMHasInitialGuess"
1119 /*@
1120     DMHasInitialGuess - does the DM object have an initial guess function
1121 
1122     Not Collective
1123 
1124     Input Parameter:
1125 .   dm - the DM object to destroy
1126 
1127     Output Parameter:
1128 .   flg - PETSC_TRUE if function exists
1129 
1130     Level: developer
1131 
1132 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1133 
1134 @*/
1135 PetscErrorCode  DMHasInitialGuess(DM dm,PetscBool  *flg)
1136 {
1137   PetscFunctionBegin;
1138   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMHasFunction"
1144 /*@
1145     DMHasFunction - does the DM object have a function
1146 
1147     Not Collective
1148 
1149     Input Parameter:
1150 .   dm - the DM object to destroy
1151 
1152     Output Parameter:
1153 .   flg - PETSC_TRUE if function exists
1154 
1155     Level: developer
1156 
1157 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1158 
1159 @*/
1160 PetscErrorCode  DMHasFunction(DM dm,PetscBool  *flg)
1161 {
1162   PetscFunctionBegin;
1163   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 #undef __FUNCT__
1168 #define __FUNCT__ "DMHasJacobian"
1169 /*@
1170     DMHasJacobian - does the DM object have a matrix function
1171 
1172     Not Collective
1173 
1174     Input Parameter:
1175 .   dm - the DM object to destroy
1176 
1177     Output Parameter:
1178 .   flg - PETSC_TRUE if function exists
1179 
1180     Level: developer
1181 
1182 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian()
1183 
1184 @*/
1185 PetscErrorCode  DMHasJacobian(DM dm,PetscBool  *flg)
1186 {
1187   PetscFunctionBegin;
1188   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 #undef __FUNCT__
1193 #define __FUNCT__ "DMComputeFunction"
1194 /*@
1195     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
1196 
1197     Collective on DM
1198 
1199     Input Parameter:
1200 +   dm - the DM object to destroy
1201 .   x - the location where the function is evaluationed, may be ignored for linear problems
1202 -   b - the vector to hold the right hand side entries
1203 
1204     Level: developer
1205 
1206 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1207          DMSetJacobian()
1208 
1209 @*/
1210 PetscErrorCode  DMComputeFunction(DM dm,Vec x,Vec b)
1211 {
1212   PetscErrorCode ierr;
1213   PetscFunctionBegin;
1214   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
1215   PetscStackPush("DM user function");
1216   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
1217   PetscStackPop;
1218   PetscFunctionReturn(0);
1219 }
1220 
1221 
1222 #undef __FUNCT__
1223 #define __FUNCT__ "DMComputeJacobian"
1224 /*@
1225     DMComputeJacobian - compute the matrix entries for the solver
1226 
1227     Collective on DM
1228 
1229     Input Parameter:
1230 +   dm - the DM object
1231 .   x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM
1232 .   A - matrix that defines the operator for the linear solve
1233 -   B - the matrix used to construct the preconditioner
1234 
1235     Level: developer
1236 
1237 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(),
1238          DMSetFunction()
1239 
1240 @*/
1241 PetscErrorCode  DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
1242 {
1243   PetscErrorCode ierr;
1244 
1245   PetscFunctionBegin;
1246   if (!dm->ops->jacobian) {
1247     ISColoring     coloring;
1248     MatFDColoring  fd;
1249 
1250     ierr = DMGetColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&coloring);CHKERRQ(ierr);
1251     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
1252     ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1253     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
1254 
1255     dm->fd = fd;
1256     dm->ops->jacobian = DMComputeJacobianDefault;
1257 
1258     /* don't know why this is needed */
1259     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1260   }
1261   if (!x) x = dm->x;
1262   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
1263 
1264   /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */
1265   if (x) {
1266     if (!dm->x) {
1267       ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr);
1268     }
1269     ierr = VecCopy(x,dm->x);CHKERRQ(ierr);
1270   }
1271   if (A != B) {
1272     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1273     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1274   }
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 
1279 PetscFList DMList                       = PETSC_NULL;
1280 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
1281 
1282 #undef __FUNCT__
1283 #define __FUNCT__ "DMSetType"
1284 /*@C
1285   DMSetType - Builds a DM, for a particular DM implementation.
1286 
1287   Collective on DM
1288 
1289   Input Parameters:
1290 + dm     - The DM object
1291 - method - The name of the DM type
1292 
1293   Options Database Key:
1294 . -dm_type <type> - Sets the DM type; use -help for a list of available types
1295 
1296   Notes:
1297   See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
1298 
1299   Level: intermediate
1300 
1301 .keywords: DM, set, type
1302 .seealso: DMGetType(), DMCreate()
1303 @*/
1304 PetscErrorCode  DMSetType(DM dm, const DMType method)
1305 {
1306   PetscErrorCode (*r)(DM);
1307   PetscBool      match;
1308   PetscErrorCode ierr;
1309 
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1312   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1313   if (match) PetscFunctionReturn(0);
1314 
1315   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1316   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
1317   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1318 
1319   if (dm->ops->destroy) {
1320     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1321   }
1322   ierr = (*r)(dm);CHKERRQ(ierr);
1323   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 #undef __FUNCT__
1328 #define __FUNCT__ "DMGetType"
1329 /*@C
1330   DMGetType - Gets the DM type name (as a string) from the DM.
1331 
1332   Not Collective
1333 
1334   Input Parameter:
1335 . dm  - The DM
1336 
1337   Output Parameter:
1338 . type - The DM type name
1339 
1340   Level: intermediate
1341 
1342 .keywords: DM, get, type, name
1343 .seealso: DMSetType(), DMCreate()
1344 @*/
1345 PetscErrorCode  DMGetType(DM dm, const DMType *type)
1346 {
1347   PetscErrorCode ierr;
1348 
1349   PetscFunctionBegin;
1350   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1351   PetscValidCharPointer(type,2);
1352   if (!DMRegisterAllCalled) {
1353     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1354   }
1355   *type = ((PetscObject)dm)->type_name;
1356   PetscFunctionReturn(0);
1357 }
1358 
1359 #undef __FUNCT__
1360 #define __FUNCT__ "DMConvert"
1361 /*@C
1362   DMConvert - Converts a DM to another DM, either of the same or different type.
1363 
1364   Collective on DM
1365 
1366   Input Parameters:
1367 + dm - the DM
1368 - newtype - new DM type (use "same" for the same type)
1369 
1370   Output Parameter:
1371 . M - pointer to new DM
1372 
1373   Notes:
1374   Cannot be used to convert a sequential DM to parallel or parallel to sequential,
1375   the MPI communicator of the generated DM is always the same as the communicator
1376   of the input DM.
1377 
1378   Level: intermediate
1379 
1380 .seealso: DMCreate()
1381 @*/
1382 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M)
1383 {
1384   DM             B;
1385   char           convname[256];
1386   PetscBool      sametype, issame;
1387   PetscErrorCode ierr;
1388 
1389   PetscFunctionBegin;
1390   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1391   PetscValidType(dm,1);
1392   PetscValidPointer(M,3);
1393   ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr);
1394   ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr);
1395   {
1396     PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL;
1397 
1398     /*
1399        Order of precedence:
1400        1) See if a specialized converter is known to the current DM.
1401        2) See if a specialized converter is known to the desired DM class.
1402        3) See if a good general converter is registered for the desired class
1403        4) See if a good general converter is known for the current matrix.
1404        5) Use a really basic converter.
1405     */
1406 
1407     /* 1) See if a specialized converter is known to the current DM and the desired class */
1408     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1409     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1410     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1411     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1412     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1413     ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1414     if (conv) goto foundconv;
1415 
1416     /* 2)  See if a specialized converter is known to the desired DM class. */
1417     ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr);
1418     ierr = DMSetType(B, newtype);CHKERRQ(ierr);
1419     ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr);
1420     ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr);
1421     ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
1422     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
1423     ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
1424     ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr);
1425     if (conv) {
1426       ierr = DMDestroy(&B);CHKERRQ(ierr);
1427       goto foundconv;
1428     }
1429 
1430 #if 0
1431     /* 3) See if a good general converter is registered for the desired class */
1432     conv = B->ops->convertfrom;
1433     ierr = DMDestroy(&B);CHKERRQ(ierr);
1434     if (conv) goto foundconv;
1435 
1436     /* 4) See if a good general converter is known for the current matrix */
1437     if (dm->ops->convert) {
1438       conv = dm->ops->convert;
1439     }
1440     if (conv) goto foundconv;
1441 #endif
1442 
1443     /* 5) Use a really basic converter. */
1444     SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype);
1445 
1446     foundconv:
1447     ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1448     ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr);
1449     ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr);
1450   }
1451   ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr);
1452   PetscFunctionReturn(0);
1453 }
1454 
1455 /*--------------------------------------------------------------------------------------------------------------------*/
1456 
1457 #undef __FUNCT__
1458 #define __FUNCT__ "DMRegister"
1459 /*@C
1460   DMRegister - See DMRegisterDynamic()
1461 
1462   Level: advanced
1463 @*/
1464 PetscErrorCode  DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1465 {
1466   char fullname[PETSC_MAX_PATH_LEN];
1467   PetscErrorCode ierr;
1468 
1469   PetscFunctionBegin;
1470   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1471   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1472   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1473   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1474   PetscFunctionReturn(0);
1475 }
1476 
1477 
1478 /*--------------------------------------------------------------------------------------------------------------------*/
1479 #undef __FUNCT__
1480 #define __FUNCT__ "DMRegisterDestroy"
1481 /*@C
1482    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1483 
1484    Not Collective
1485 
1486    Level: advanced
1487 
1488 .keywords: DM, register, destroy
1489 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1490 @*/
1491 PetscErrorCode  DMRegisterDestroy(void)
1492 {
1493   PetscErrorCode ierr;
1494 
1495   PetscFunctionBegin;
1496   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1497   DMRegisterAllCalled = PETSC_FALSE;
1498   PetscFunctionReturn(0);
1499 }
1500 
1501 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1502 #include <mex.h>
1503 
1504 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext;
1505 
1506 #undef __FUNCT__
1507 #define __FUNCT__ "DMComputeFunction_Matlab"
1508 /*
1509    DMComputeFunction_Matlab - Calls the function that has been set with
1510                          DMSetFunctionMatlab().
1511 
1512    For linear problems x is null
1513 
1514 .seealso: DMSetFunction(), DMGetFunction()
1515 */
1516 PetscErrorCode  DMComputeFunction_Matlab(DM dm,Vec x,Vec y)
1517 {
1518   PetscErrorCode    ierr;
1519   DMMatlabContext   *sctx;
1520   int               nlhs = 1,nrhs = 4;
1521   mxArray	    *plhs[1],*prhs[4];
1522   long long int     lx = 0,ly = 0,ls = 0;
1523 
1524   PetscFunctionBegin;
1525   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1526   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1527   PetscCheckSameComm(dm,1,y,3);
1528 
1529   /* call Matlab function in ctx with arguments x and y */
1530   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1531   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1532   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1533   ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr);
1534   prhs[0] =  mxCreateDoubleScalar((double)ls);
1535   prhs[1] =  mxCreateDoubleScalar((double)lx);
1536   prhs[2] =  mxCreateDoubleScalar((double)ly);
1537   prhs[3] =  mxCreateString(sctx->funcname);
1538   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr);
1539   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
1540   mxDestroyArray(prhs[0]);
1541   mxDestroyArray(prhs[1]);
1542   mxDestroyArray(prhs[2]);
1543   mxDestroyArray(prhs[3]);
1544   mxDestroyArray(plhs[0]);
1545   PetscFunctionReturn(0);
1546 }
1547 
1548 
1549 #undef __FUNCT__
1550 #define __FUNCT__ "DMSetFunctionMatlab"
1551 /*
1552    DMSetFunctionMatlab - Sets the function evaluation routine
1553 
1554 */
1555 PetscErrorCode  DMSetFunctionMatlab(DM dm,const char *func)
1556 {
1557   PetscErrorCode    ierr;
1558   DMMatlabContext   *sctx;
1559 
1560   PetscFunctionBegin;
1561   /* currently sctx is memory bleed */
1562   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1563   if (!sctx) {
1564     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1565   }
1566   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
1567   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1568   ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr);
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 #undef __FUNCT__
1573 #define __FUNCT__ "DMComputeJacobian_Matlab"
1574 /*
1575    DMComputeJacobian_Matlab - Calls the function that has been set with
1576                          DMSetJacobianMatlab().
1577 
1578    For linear problems x is null
1579 
1580 .seealso: DMSetFunction(), DMGetFunction()
1581 */
1582 PetscErrorCode  DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str)
1583 {
1584   PetscErrorCode    ierr;
1585   DMMatlabContext   *sctx;
1586   int               nlhs = 2,nrhs = 5;
1587   mxArray	    *plhs[2],*prhs[5];
1588   long long int     lx = 0,lA = 0,lB = 0,ls = 0;
1589 
1590   PetscFunctionBegin;
1591   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1592   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
1593 
1594   /* call MATLAB function in ctx with arguments x, A, and B */
1595   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1596   ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr);
1597   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
1598   ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr);
1599   ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr);
1600   prhs[0] =  mxCreateDoubleScalar((double)ls);
1601   prhs[1] =  mxCreateDoubleScalar((double)lx);
1602   prhs[2] =  mxCreateDoubleScalar((double)lA);
1603   prhs[3] =  mxCreateDoubleScalar((double)lB);
1604   prhs[4] =  mxCreateString(sctx->jacname);
1605   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr);
1606   *str    =  (MatStructure) mxGetScalar(plhs[0]);
1607   ierr    =  (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr);
1608   mxDestroyArray(prhs[0]);
1609   mxDestroyArray(prhs[1]);
1610   mxDestroyArray(prhs[2]);
1611   mxDestroyArray(prhs[3]);
1612   mxDestroyArray(prhs[4]);
1613   mxDestroyArray(plhs[0]);
1614   mxDestroyArray(plhs[1]);
1615   PetscFunctionReturn(0);
1616 }
1617 
1618 
1619 #undef __FUNCT__
1620 #define __FUNCT__ "DMSetJacobianMatlab"
1621 /*
1622    DMSetJacobianMatlab - Sets the Jacobian function evaluation routine
1623 
1624 */
1625 PetscErrorCode  DMSetJacobianMatlab(DM dm,const char *func)
1626 {
1627   PetscErrorCode    ierr;
1628   DMMatlabContext   *sctx;
1629 
1630   PetscFunctionBegin;
1631   /* currently sctx is memory bleed */
1632   ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr);
1633   if (!sctx) {
1634     ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr);
1635   }
1636   ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr);
1637   ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr);
1638   ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr);
1639   PetscFunctionReturn(0);
1640 }
1641 #endif
1642 
1643 #undef __FUNCT__
1644 #define __FUNCT__ "DMLoad"
1645 /*@C
1646   DMLoad - Loads a DM that has been stored in binary or HDF5 format
1647   with DMView().
1648 
1649   Collective on PetscViewer
1650 
1651   Input Parameters:
1652 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or
1653            some related function before a call to DMLoad().
1654 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
1655            HDF5 file viewer, obtained from PetscViewerHDF5Open()
1656 
1657    Level: intermediate
1658 
1659   Notes:
1660   Defaults to the DM DA.
1661 
1662   Notes for advanced users:
1663   Most users should not need to know the details of the binary storage
1664   format, since DMLoad() and DMView() completely hide these details.
1665   But for anyone who's interested, the standard binary matrix storage
1666   format is
1667 .vb
1668      has not yet been determined
1669 .ve
1670 
1671    In addition, PETSc automatically does the byte swapping for
1672 machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1673 linux, Windows and the paragon; thus if you write your own binary
1674 read/write routines you have to swap the bytes; see PetscBinaryRead()
1675 and PetscBinaryWrite() to see how this may be done.
1676 
1677   Concepts: vector^loading from file
1678 
1679 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad()
1680 @*/
1681 PetscErrorCode  DMLoad(DM newdm, PetscViewer viewer)
1682 {
1683   PetscErrorCode ierr;
1684 
1685   PetscFunctionBegin;
1686   PetscValidHeaderSpecific(newdm,DM_CLASSID,1);
1687   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1688 
1689   if (!((PetscObject)newdm)->type_name) {
1690     ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr);
1691   }
1692   ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);
1693   PetscFunctionReturn(0);
1694 }
1695 
1696