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