xref: /petsc/src/dm/interface/dm.c (revision 94013140d43fa11e67423739cdc77d55e3c97eb1)
1 #define PETSCDM_DLL
2 
3 #include "private/dmimpl.h"     /*I      "petscda.h"     I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMSetVecType"
7 /*@C
8        DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector()
9 
10    Logically Collective on DA
11 
12    Input Parameter:
13 +  da - initial distributed array
14 .  ctype - the vector type, currently either VECSTANDARD or VECCUDA
15 
16    Options Database:
17 .   -da_vec_type ctype
18 
19    Level: intermediate
20 
21 .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DMDestroy(), DA, DAInterpolationType, VecType
22 @*/
23 PetscErrorCode PETSCDM_DLLEXPORT DMSetVecType(DM da,const VecType ctype)
24 {
25   PetscErrorCode ierr;
26 
27   PetscFunctionBegin;
28   PetscValidHeaderSpecific(da,DM_CLASSID,1);
29   ierr = PetscFree(da->vectype);CHKERRQ(ierr);
30   ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 #undef __FUNCT__
35 #define __FUNCT__ "DMSetOptionsPrefix"
36 /*@C
37    DMSetOptionsPrefix - Sets the prefix used for searching for all
38    DA options in the database.
39 
40    Logically Collective on DA
41 
42    Input Parameter:
43 +  da - the DA context
44 -  prefix - the prefix to prepend to all option names
45 
46    Notes:
47    A hyphen (-) must NOT be given at the beginning of the prefix name.
48    The first character of all runtime options is AUTOMATICALLY the hyphen.
49 
50    Level: advanced
51 
52 .keywords: DA, set, options, prefix, database
53 
54 .seealso: DMSetFromOptions()
55 @*/
56 PetscErrorCode PETSCDM_DLLEXPORT DMSetOptionsPrefix(DM dm,const char prefix[])
57 {
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62   ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "DMDestroy"
68 /*@
69     DMDestroy - Destroys a vector packer or DA.
70 
71     Collective on DM
72 
73     Input Parameter:
74 .   dm - the DM object to destroy
75 
76     Level: developer
77 
78 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
79 
80 @*/
81 PetscErrorCode PETSCDM_DLLEXPORT DMDestroy(DM dm)
82 {
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "DMSetUp"
92 /*@
93     DMSetUp - sets up the data structures inside a DM object
94 
95     Collective on DM
96 
97     Input Parameter:
98 .   dm - the DM object to setup
99 
100     Level: developer
101 
102 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
103 
104 @*/
105 PetscErrorCode PETSCDM_DLLEXPORT DMSetUp(DM dm)
106 {
107   PetscErrorCode ierr;
108 
109   PetscFunctionBegin;
110   if (dm->ops->setup) {
111     ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr);
112   }
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "DMSetFromOptions"
118 /*@
119     DMSetFromOptions - sets parameters in a DM from the options database
120 
121     Collective on DM
122 
123     Input Parameter:
124 .   dm - the DM object to set options for
125 
126     Level: developer
127 
128 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
129 
130 @*/
131 PetscErrorCode PETSCDM_DLLEXPORT DMSetFromOptions(DM dm)
132 {
133   PetscErrorCode ierr;
134 
135   PetscFunctionBegin;
136   if (dm->ops->setfromoptions) {
137     ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr);
138   }
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "DMView"
144 /*@
145     DMView - Views a vector packer or DA.
146 
147     Collective on DM
148 
149     Input Parameter:
150 +   dm - the DM object to view
151 -   v - the viewer
152 
153     Level: developer
154 
155 .seealso DMDestroy(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
156 
157 @*/
158 PetscErrorCode PETSCDM_DLLEXPORT DMView(DM dm,PetscViewer v)
159 {
160   PetscErrorCode ierr;
161 
162   PetscFunctionBegin;
163   if (dm->ops->view) {
164     ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr);
165   }
166   PetscFunctionReturn(0);
167 }
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "DMCreateGlobalVector"
171 /*@
172     DMCreateGlobalVector - Creates a global vector from a DA or DMComposite object
173 
174     Collective on DM
175 
176     Input Parameter:
177 .   dm - the DM object
178 
179     Output Parameter:
180 .   vec - the global vector
181 
182     Level: developer
183 
184 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
185 
186 @*/
187 PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector(DM dm,Vec *vec)
188 {
189   PetscErrorCode ierr;
190 
191   PetscFunctionBegin;
192   ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "DMCreateLocalVector"
198 /*@
199     DMCreateLocalVector - Creates a local vector from a DA or DMComposite object
200 
201     Not Collective
202 
203     Input Parameter:
204 .   dm - the DM object
205 
206     Output Parameter:
207 .   vec - the local vector
208 
209     Level: developer
210 
211 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix()
212 
213 @*/
214 PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector(DM dm,Vec *vec)
215 {
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr);
220   PetscFunctionReturn(0);
221 }
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "DMGetInterpolation"
225 /*@
226     DMGetInterpolation - Gets interpolation matrix between two DA or DMComposite objects
227 
228     Collective on DM
229 
230     Input Parameter:
231 +   dm1 - the DM object
232 -   dm2 - the second, finer DM object
233 
234     Output Parameter:
235 +  mat - the interpolation
236 -  vec - the scaling (optional)
237 
238     Level: developer
239 
240 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix()
241 
242 @*/
243 PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec)
244 {
245   PetscErrorCode ierr;
246 
247   PetscFunctionBegin;
248   ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "DMGetInjection"
254 /*@
255     DMGetInjection - Gets injection matrix between two DA or DMComposite objects
256 
257     Collective on DM
258 
259     Input Parameter:
260 +   dm1 - the DM object
261 -   dm2 - the second, finer DM object
262 
263     Output Parameter:
264 .   ctx - the injection
265 
266     Level: developer
267 
268 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation()
269 
270 @*/
271 PetscErrorCode PETSCDM_DLLEXPORT DMGetInjection(DM dm1,DM dm2,VecScatter *ctx)
272 {
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "DMGetColoring"
282 /*@
283     DMGetColoring - Gets coloring for a DA or DMComposite
284 
285     Collective on DM
286 
287     Input Parameter:
288 +   dm - the DM object
289 .   ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL
290 -   matype - either MATAIJ or MATBAIJ
291 
292     Output Parameter:
293 .   coloring - the coloring
294 
295     Level: developer
296 
297 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
298 
299 @*/
300 PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
301 {
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet");
306   ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "DMGetMatrix"
312 /*@C
313     DMGetMatrix - Gets empty Jacobian for a DA or DMComposite
314 
315     Collective on DM
316 
317     Input Parameter:
318 +   dm - the DM object
319 -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or
320             any type which inherits from one of these (such as MATAIJ)
321 
322     Output Parameter:
323 .   mat - the empty Jacobian
324 
325     Level: developer
326 
327     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
328        do not need to do it yourself.
329 
330        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
331        the nonzero pattern call DASetMatPreallocateOnly()
332 
333        For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
334        internally by PETSc.
335 
336        For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
337        the indices for the global numbering for DAs which is complicated.
338 
339 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
340 
341 @*/
342 PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix(DM dm, const MatType mtype,Mat *mat)
343 {
344   PetscErrorCode ierr;
345 
346   PetscFunctionBegin;
347   ierr = (*dm->ops->getmatrix)(dm,mtype,mat);CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "DMRefine"
353 /*@
354     DMRefine - Refines a DM object
355 
356     Collective on DM
357 
358     Input Parameter:
359 +   dm - the DM object
360 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
361 
362     Output Parameter:
363 .   dmf - the refined DM
364 
365     Level: developer
366 
367 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
368 
369 @*/
370 PetscErrorCode PETSCDM_DLLEXPORT DMRefine(DM dm,MPI_Comm comm,DM *dmf)
371 {
372   PetscErrorCode ierr;
373 
374   PetscFunctionBegin;
375   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
376   PetscFunctionReturn(0);
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "DMGlobalToLocalBegin"
381 /*@
382     DMGlobalToLocalBegin - Begins updating local vectors from global vector
383 
384     Neighbor-wise Collective on DM
385 
386     Input Parameters:
387 +   dm - the DM object
388 .   g - the global vector
389 .   mode - INSERT_VALUES or ADD_VALUES
390 -   l - the local vector
391 
392 
393     Level: beginner
394 
395 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
396 
397 @*/
398 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
399 {
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "DMGlobalToLocalEnd"
409 /*@
410     DMGlobalToLocalEnd - Ends updating local vectors from global vector
411 
412     Neighbor-wise Collective on DM
413 
414     Input Parameters:
415 +   dm - the DM object
416 .   g - the global vector
417 .   mode - INSERT_VALUES or ADD_VALUES
418 -   l - the local vector
419 
420 
421     Level: beginner
422 
423 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
424 
425 @*/
426 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
427 {
428   PetscErrorCode ierr;
429 
430   PetscFunctionBegin;
431   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "DMLocalToGlobalBegin"
437 /*@
438     DMLocalToGlobalBegin - updates global vectors from local vectors
439 
440     Neighbor-wise Collective on DM
441 
442     Input Parameters:
443 +   dm - the DM object
444 .   g - the global vector
445 .   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
446            base point.
447 -   l - the local vector
448 
449     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
450            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
451            global array to the final global array with VecAXPY().
452 
453     Level: beginner
454 
455 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
456 
457 @*/
458 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalBegin(DM dm,Vec g,InsertMode mode,Vec l)
459 {
460   PetscErrorCode ierr;
461 
462   PetscFunctionBegin;
463   ierr = (*dm->ops->localtoglobalbegin)(dm,g,mode,l);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "DMLocalToGlobalEnd"
469 /*@
470     DMLocalToGlobalEnd - updates global vectors from local vectors
471 
472     Neighbor-wise Collective on DM
473 
474     Input Parameters:
475 +   dm - the DM object
476 .   g - the global vector
477 .   mode - INSERT_VALUES or ADD_VALUES
478 -   l - the local vector
479 
480 
481     Level: beginner
482 
483 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
484 
485 @*/
486 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalEnd(DM dm,Vec g,InsertMode mode,Vec l)
487 {
488   PetscErrorCode ierr;
489 
490   PetscFunctionBegin;
491   ierr = (*dm->ops->localtoglobalend)(dm,g,mode,l);CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "DMComputeJacobianDefault"
497 /*@
498     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
499 
500     Collective on DM
501 
502     Input Parameter:
503 +   dm - the DM object
504 .   x - location to compute Jacobian at; may be ignored for linear problems
505 .   A - matrix that defines the operator for the linear solve
506 -   B - the matrix used to construct the preconditioner
507 
508     Level: developer
509 
510 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
511          DMSetFunction()
512 
513 @*/
514 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
515 {
516   PetscErrorCode ierr;
517   PetscFunctionBegin;
518   *stflag = SAME_NONZERO_PATTERN;
519   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
520   if (A != B) {
521     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
522     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523   }
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "DMCoarsen"
529 /*@
530     DMCoarsen - Coarsens a DM object
531 
532     Collective on DM
533 
534     Input Parameter:
535 +   dm - the DM object
536 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
537 
538     Output Parameter:
539 .   dmc - the coarsened DM
540 
541     Level: developer
542 
543 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
544 
545 @*/
546 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
547 {
548   PetscErrorCode ierr;
549 
550   PetscFunctionBegin;
551   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
552   (*dmc)->ops->initialguess = dm->ops->initialguess;
553   (*dmc)->ops->function     = dm->ops->function;
554   (*dmc)->ops->functionj    = dm->ops->functionj;
555   if (dm->ops->jacobian != DMComputeJacobianDefault) {
556     (*dmc)->ops->jacobian     = dm->ops->jacobian;
557   }
558   PetscFunctionReturn(0);
559 }
560 
561 #undef __FUNCT__
562 #define __FUNCT__ "DMRefineHierarchy"
563 /*@C
564     DMRefineHierarchy - Refines a DM object, all levels at once
565 
566     Collective on DM
567 
568     Input Parameter:
569 +   dm - the DM object
570 -   nlevels - the number of levels of refinement
571 
572     Output Parameter:
573 .   dmf - the refined DM hierarchy
574 
575     Level: developer
576 
577 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
578 
579 @*/
580 PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
581 {
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
586   if (nlevels == 0) PetscFunctionReturn(0);
587   if (dm->ops->refinehierarchy) {
588     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
589   } else if (dm->ops->refine) {
590     PetscInt i;
591 
592     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
593     for (i=1; i<nlevels; i++) {
594       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
595     }
596   } else {
597     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
598   }
599   PetscFunctionReturn(0);
600 }
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "DMCoarsenHierarchy"
604 /*@C
605     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
606 
607     Collective on DM
608 
609     Input Parameter:
610 +   dm - the DM object
611 -   nlevels - the number of levels of coarsening
612 
613     Output Parameter:
614 .   dmc - the coarsened DM hierarchy
615 
616     Level: developer
617 
618 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
619 
620 @*/
621 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
622 {
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
627   if (nlevels == 0) PetscFunctionReturn(0);
628   PetscValidPointer(dmc,3);
629   if (dm->ops->coarsenhierarchy) {
630     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
631   } else if (dm->ops->coarsen) {
632     PetscInt i;
633 
634     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
635     for (i=1; i<nlevels; i++) {
636       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
637     }
638   } else {
639     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
640   }
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "DMGetAggregates"
646 /*@
647    DMGetAggregates - Gets the aggregates that map between
648    grids associated with two DMs.
649 
650    Collective on DM
651 
652    Input Parameters:
653 +  dmc - the coarse grid DM
654 -  dmf - the fine grid DM
655 
656    Output Parameters:
657 .  rest - the restriction matrix (transpose of the projection matrix)
658 
659    Level: intermediate
660 
661 .keywords: interpolation, restriction, multigrid
662 
663 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
664 @*/
665 PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates(DM dmc, DM dmf, Mat *rest)
666 {
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "DMSetContext"
676 /*@
677     DMSetContext - Set a user context into a DM object
678 
679     Not Collective
680 
681     Input Parameters:
682 +   dm - the DM object
683 -   ctx - the user context
684 
685     Level: intermediate
686 
687 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
688 
689 @*/
690 PetscErrorCode PETSCDM_DLLEXPORT DMSetContext(DM dm,void *ctx)
691 {
692   PetscFunctionBegin;
693   dm->ctx = ctx;
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "DMGetContext"
699 /*@
700     DMGetContext - Gets a user context from a DM object
701 
702     Not Collective
703 
704     Input Parameter:
705 .   dm - the DM object
706 
707     Output Parameter:
708 .   ctx - the user context
709 
710     Level: intermediate
711 
712 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
713 
714 @*/
715 PetscErrorCode PETSCDM_DLLEXPORT DMGetContext(DM dm,void **ctx)
716 {
717   PetscFunctionBegin;
718   *ctx = dm->ctx;
719   PetscFunctionReturn(0);
720 }
721 
722 #undef __FUNCT__
723 #define __FUNCT__ "DMSetInitialGuess"
724 /*@
725     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
726 
727     Logically Collective on DM
728 
729     Input Parameter:
730 +   dm - the DM object to destroy
731 -   f - the function to compute the initial guess
732 
733     Level: intermediate
734 
735 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
736 
737 @*/
738 PetscErrorCode PETSCDM_DLLEXPORT DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
739 {
740   PetscFunctionBegin;
741   dm->ops->initialguess = f;
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "DMSetFunction"
747 /*@
748     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
749 
750     Logically Collective on DM
751 
752     Input Parameter:
753 +   dm - the DM object
754 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
755 
756     Level: intermediate
757 
758     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
759            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
760 
761 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
762          DMSetJacobian()
763 
764 @*/
765 PetscErrorCode PETSCDM_DLLEXPORT DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
766 {
767   PetscFunctionBegin;
768   dm->ops->function = f;
769   if (f) {
770     dm->ops->functionj = f;
771   }
772   PetscFunctionReturn(0);
773 }
774 
775 #undef __FUNCT__
776 #define __FUNCT__ "DMSetJacobian"
777 /*@
778     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
779 
780     Logically Collective on DM
781 
782     Input Parameter:
783 +   dm - the DM object to destroy
784 -   f - the function to compute the matrix entries
785 
786     Level: intermediate
787 
788 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
789          DMSetFunction()
790 
791 @*/
792 PetscErrorCode PETSCDM_DLLEXPORT DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
793 {
794   PetscFunctionBegin;
795   dm->ops->jacobian = f;
796   PetscFunctionReturn(0);
797 }
798 
799 #undef __FUNCT__
800 #define __FUNCT__ "DMComputeInitialGuess"
801 /*@
802     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
803 
804     Collective on DM
805 
806     Input Parameter:
807 +   dm - the DM object to destroy
808 -   x - the vector to hold the initial guess values
809 
810     Level: developer
811 
812 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetRhs(), DMSetMat()
813 
814 @*/
815 PetscErrorCode PETSCDM_DLLEXPORT DMComputeInitialGuess(DM dm,Vec x)
816 {
817   PetscErrorCode ierr;
818   PetscFunctionBegin;
819   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
820   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "DMHasInitialGuess"
826 /*@
827     DMHasInitialGuess - does the DM object have an initial guess function
828 
829     Not Collective
830 
831     Input Parameter:
832 .   dm - the DM object to destroy
833 
834     Output Parameter:
835 .   flg - PETSC_TRUE if function exists
836 
837     Level: developer
838 
839 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
840 
841 @*/
842 PetscErrorCode PETSCDM_DLLEXPORT DMHasInitialGuess(DM dm,PetscBool  *flg)
843 {
844   PetscFunctionBegin;
845   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
846   PetscFunctionReturn(0);
847 }
848 
849 #undef __FUNCT__
850 #define __FUNCT__ "DMHasFunction"
851 /*@
852     DMHasFunction - does the DM object have a function
853 
854     Not Collective
855 
856     Input Parameter:
857 .   dm - the DM object to destroy
858 
859     Output Parameter:
860 .   flg - PETSC_TRUE if function exists
861 
862     Level: developer
863 
864 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
865 
866 @*/
867 PetscErrorCode PETSCDM_DLLEXPORT DMHasFunction(DM dm,PetscBool  *flg)
868 {
869   PetscFunctionBegin;
870   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
871   PetscFunctionReturn(0);
872 }
873 
874 #undef __FUNCT__
875 #define __FUNCT__ "DMHasJacobian"
876 /*@
877     DMHasJacobian - does the DM object have a matrix function
878 
879     Not Collective
880 
881     Input Parameter:
882 .   dm - the DM object to destroy
883 
884     Output Parameter:
885 .   flg - PETSC_TRUE if function exists
886 
887     Level: developer
888 
889 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
890 
891 @*/
892 PetscErrorCode PETSCDM_DLLEXPORT DMHasJacobian(DM dm,PetscBool  *flg)
893 {
894   PetscFunctionBegin;
895   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
896   PetscFunctionReturn(0);
897 }
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "DMComputeFunction"
901 /*@
902     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
903 
904     Collective on DM
905 
906     Input Parameter:
907 +   dm - the DM object to destroy
908 .   x - the location where the function is evaluationed, may be ignored for linear problems
909 -   b - the vector to hold the right hand side entries
910 
911     Level: developer
912 
913 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
914          DMSetJacobian()
915 
916 @*/
917 PetscErrorCode PETSCDM_DLLEXPORT DMComputeFunction(DM dm,Vec x,Vec b)
918 {
919   PetscErrorCode ierr;
920   PetscFunctionBegin;
921   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
922   if (!x) x = dm->x;
923   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 
928 #undef __FUNCT__
929 #define __FUNCT__ "DMComputeJacobian"
930 /*@
931     DMComputeJacobian - compute the matrix entries for the solver
932 
933     Collective on DM
934 
935     Input Parameter:
936 +   dm - the DM object
937 .   x - location to compute Jacobian at; may be ignored for linear problems
938 .   A - matrix that defines the operator for the linear solve
939 -   B - the matrix used to construct the preconditioner
940 
941     Level: developer
942 
943 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
944          DMSetFunction()
945 
946 @*/
947 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
948 {
949   PetscErrorCode ierr;
950 
951   PetscFunctionBegin;
952   if (!dm->ops->jacobian) {
953     ISColoring    coloring;
954     MatFDColoring fd;
955 
956     ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
957     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
958     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
959     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
960     dm->fd = fd;
961     dm->ops->jacobian = DMComputeJacobianDefault;
962 
963     if (!dm->x) {
964       ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr);
965     }
966   }
967   if (!x) x = dm->x;
968   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
969   PetscFunctionReturn(0);
970 }
971 
972 PetscFList DMList                       = PETSC_NULL;
973 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
974 
975 #undef __FUNCT__
976 #define __FUNCT__ "DMSetType"
977 /*@C
978   DMSetType - Builds a DM, for a particular DM implementation.
979 
980   Collective on DM
981 
982   Input Parameters:
983 + dm     - The DM object
984 - method - The name of the DM type
985 
986   Options Database Key:
987 . -dm_type <type> - Sets the DM type; use -help for a list of available types
988 
989   Notes:
990   See "petsc/include/petscda.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
991 
992   Level: intermediate
993 
994 .keywords: DM, set, type
995 .seealso: DMGetType(), DMCreate()
996 @*/
997 PetscErrorCode PETSCDM_DLLEXPORT DMSetType(DM dm, const DMType method)
998 {
999   PetscErrorCode (*r)(DM);
1000   PetscBool      match;
1001   PetscErrorCode ierr;
1002 
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1005   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
1006   if (match) PetscFunctionReturn(0);
1007 
1008   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1009   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,(void (**)(void)) &r);CHKERRQ(ierr);
1010   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
1011 
1012   if (dm->ops->destroy) {
1013     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1014   }
1015   ierr = (*r)(dm);CHKERRQ(ierr);
1016   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 #undef __FUNCT__
1021 #define __FUNCT__ "DMGetType"
1022 /*@C
1023   DMGetType - Gets the DM type name (as a string) from the DM.
1024 
1025   Not Collective
1026 
1027   Input Parameter:
1028 . dm  - The DM
1029 
1030   Output Parameter:
1031 . type - The DM type name
1032 
1033   Level: intermediate
1034 
1035 .keywords: DM, get, type, name
1036 .seealso: DMSetType(), DMCreate()
1037 @*/
1038 PetscErrorCode PETSCDM_DLLEXPORT DMGetType(DM dm, const DMType *type)
1039 {
1040   PetscErrorCode ierr;
1041 
1042   PetscFunctionBegin;
1043   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1044   PetscValidCharPointer(type,2);
1045   if (!DMRegisterAllCalled) {
1046     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1047   }
1048   *type = ((PetscObject)dm)->type_name;
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 
1053 /*--------------------------------------------------------------------------------------------------------------------*/
1054 
1055 #undef __FUNCT__
1056 #define __FUNCT__ "DMRegister"
1057 /*@C
1058   DMRegister - See DMRegisterDynamic()
1059 
1060   Level: advanced
1061 @*/
1062 PetscErrorCode PETSCDM_DLLEXPORT DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1063 {
1064   char fullname[PETSC_MAX_PATH_LEN];
1065   PetscErrorCode ierr;
1066 
1067   PetscFunctionBegin;
1068   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1069   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1070   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1071   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 
1076 /*--------------------------------------------------------------------------------------------------------------------*/
1077 #undef __FUNCT__
1078 #define __FUNCT__ "DMRegisterDestroy"
1079 /*@C
1080    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1081 
1082    Not Collective
1083 
1084    Level: advanced
1085 
1086 .keywords: DM, register, destroy
1087 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1088 @*/
1089 PetscErrorCode PETSCDM_DLLEXPORT DMRegisterDestroy(void)
1090 {
1091   PetscErrorCode ierr;
1092 
1093   PetscFunctionBegin;
1094   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1095   DMRegisterAllCalled = PETSC_FALSE;
1096   PetscFunctionReturn(0);
1097 }
1098