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