xref: /petsc/src/dm/interface/dm.c (revision 9a42bb27a39f0cdf3306a1e22d33cd9809484eaa)
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 DACreateLocalVector() and DACreateGlobalVector()
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, MATLUSOL, etc.).
321 
322     Output Parameter:
323 .   mat - the empty Jacobian
324 
325     Level: developer
326 
327 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix()
328 
329 @*/
330 PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix(DM dm, const MatType mtype,Mat *mat)
331 {
332   PetscErrorCode ierr;
333 
334   PetscFunctionBegin;
335   ierr = (*dm->ops->getmatrix)(dm,mtype,mat);CHKERRQ(ierr);
336   PetscFunctionReturn(0);
337 }
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "DMRefine"
341 /*@
342     DMRefine - Refines a DM object
343 
344     Collective on DM
345 
346     Input Parameter:
347 +   dm - the DM object
348 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
349 
350     Output Parameter:
351 .   dmf - the refined DM
352 
353     Level: developer
354 
355 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
356 
357 @*/
358 PetscErrorCode PETSCDM_DLLEXPORT DMRefine(DM dm,MPI_Comm comm,DM *dmf)
359 {
360   PetscErrorCode ierr;
361 
362   PetscFunctionBegin;
363   ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "DMGlobalToLocalBegin"
369 /*@
370     DMGlobalToLocalBegin - Begins updating local vectors from global vector
371 
372     Neighbor-wise Collective on DM
373 
374     Input Parameters:
375 +   dm - the DM object
376 .   g - the global vector
377 .   mode - INSERT_VALUES or ADD_VALUES
378 -   l - the local vector
379 
380 
381     Level: beginner
382 
383 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
384 
385 @*/
386 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l)
387 {
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr);
392   PetscFunctionReturn(0);
393 }
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "DMGlobalToLocalEnd"
397 /*@
398     DMGlobalToLocalEnd - Ends updating local vectors from global vector
399 
400     Neighbor-wise Collective on DM
401 
402     Input Parameters:
403 +   dm - the DM object
404 .   g - the global vector
405 .   mode - INSERT_VALUES or ADD_VALUES
406 -   l - the local vector
407 
408 
409     Level: beginner
410 
411 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin()
412 
413 @*/
414 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l)
415 {
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr);
420   PetscFunctionReturn(0);
421 }
422 
423 #undef __FUNCT__
424 #define __FUNCT__ "DMLocalToGlobalBegin"
425 /*@
426     DMLocalToGlobalBegin - updates global vectors from local vectors
427 
428     Neighbor-wise Collective on DM
429 
430     Input Parameters:
431 +   dm - the DM object
432 .   g - the global vector
433 .   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
434            base point.
435 -   l - the local vector
436 
437     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
438            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
439            global array to the final global array with VecAXPY().
440 
441     Level: beginner
442 
443 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin()
444 
445 @*/
446 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalBegin(DM dm,Vec g,InsertMode mode,Vec l)
447 {
448   PetscErrorCode ierr;
449 
450   PetscFunctionBegin;
451   ierr = (*dm->ops->localtoglobalbegin)(dm,g,mode,l);CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "DMLocalToGlobalEnd"
457 /*@
458     DMLocalToGlobalEnd - updates global vectors from local vectors
459 
460     Neighbor-wise Collective on DM
461 
462     Input Parameters:
463 +   dm - the DM object
464 .   g - the global vector
465 .   mode - INSERT_VALUES or ADD_VALUES
466 -   l - the local vector
467 
468 
469     Level: beginner
470 
471 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd()
472 
473 @*/
474 PetscErrorCode PETSCDM_DLLEXPORT DMLocalToGlobalEnd(DM dm,Vec g,InsertMode mode,Vec l)
475 {
476   PetscErrorCode ierr;
477 
478   PetscFunctionBegin;
479   ierr = (*dm->ops->localtoglobalend)(dm,g,mode,l);CHKERRQ(ierr);
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "DMComputeJacobianDefault"
485 /*@
486     DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided
487 
488     Collective on DM
489 
490     Input Parameter:
491 +   dm - the DM object
492 .   x - location to compute Jacobian at; may be ignored for linear problems
493 .   A - matrix that defines the operator for the linear solve
494 -   B - the matrix used to construct the preconditioner
495 
496     Level: developer
497 
498 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
499          DMSetFunction()
500 
501 @*/
502 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
503 {
504   PetscErrorCode ierr;
505   PetscFunctionBegin;
506   *stflag = SAME_NONZERO_PATTERN;
507   ierr  = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr);
508   if (A != B) {
509     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
510     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
511   }
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "DMCoarsen"
517 /*@
518     DMCoarsen - Coarsens a DM object
519 
520     Collective on DM
521 
522     Input Parameter:
523 +   dm - the DM object
524 -   comm - the communicator to contain the new DM object (or PETSC_NULL)
525 
526     Output Parameter:
527 .   dmc - the coarsened DM
528 
529     Level: developer
530 
531 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
532 
533 @*/
534 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen(DM dm, MPI_Comm comm, DM *dmc)
535 {
536   PetscErrorCode ierr;
537 
538   PetscFunctionBegin;
539   ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr);
540   (*dmc)->ops->initialguess = dm->ops->initialguess;
541   (*dmc)->ops->function     = dm->ops->function;
542   (*dmc)->ops->functionj    = dm->ops->functionj;
543   if (dm->ops->jacobian != DMComputeJacobianDefault) {
544     (*dmc)->ops->jacobian     = dm->ops->jacobian;
545   }
546   PetscFunctionReturn(0);
547 }
548 
549 #undef __FUNCT__
550 #define __FUNCT__ "DMRefineHierarchy"
551 /*@C
552     DMRefineHierarchy - Refines a DM object, all levels at once
553 
554     Collective on DM
555 
556     Input Parameter:
557 +   dm - the DM object
558 -   nlevels - the number of levels of refinement
559 
560     Output Parameter:
561 .   dmf - the refined DM hierarchy
562 
563     Level: developer
564 
565 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
566 
567 @*/
568 PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[])
569 {
570   PetscErrorCode ierr;
571 
572   PetscFunctionBegin;
573   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
574   if (nlevels == 0) PetscFunctionReturn(0);
575   if (dm->ops->refinehierarchy) {
576     ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr);
577   } else if (dm->ops->refine) {
578     PetscInt i;
579 
580     ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr);
581     for (i=1; i<nlevels; i++) {
582       ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr);
583     }
584   } else {
585     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet");
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "DMCoarsenHierarchy"
592 /*@C
593     DMCoarsenHierarchy - Coarsens a DM object, all levels at once
594 
595     Collective on DM
596 
597     Input Parameter:
598 +   dm - the DM object
599 -   nlevels - the number of levels of coarsening
600 
601     Output Parameter:
602 .   dmc - the coarsened DM hierarchy
603 
604     Level: developer
605 
606 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation()
607 
608 @*/
609 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[])
610 {
611   PetscErrorCode ierr;
612 
613   PetscFunctionBegin;
614   if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
615   if (nlevels == 0) PetscFunctionReturn(0);
616   PetscValidPointer(dmc,3);
617   if (dm->ops->coarsenhierarchy) {
618     ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr);
619   } else if (dm->ops->coarsen) {
620     PetscInt i;
621 
622     ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr);
623     for (i=1; i<nlevels; i++) {
624       ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr);
625     }
626   } else {
627     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet");
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "DMGetAggregates"
634 /*@
635    DMGetAggregates - Gets the aggregates that map between
636    grids associated with two DMs.
637 
638    Collective on DM
639 
640    Input Parameters:
641 +  dmc - the coarse grid DM
642 -  dmf - the fine grid DM
643 
644    Output Parameters:
645 .  rest - the restriction matrix (transpose of the projection matrix)
646 
647    Level: intermediate
648 
649 .keywords: interpolation, restriction, multigrid
650 
651 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation()
652 @*/
653 PetscErrorCode PETSCDM_DLLEXPORT DMGetAggregates(DM dmc, DM dmf, Mat *rest)
654 {
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "DMSetContext"
664 /*@
665     DMSetContext - Set a user context into a DM object
666 
667     Not Collective
668 
669     Input Parameters:
670 +   dm - the DM object
671 -   ctx - the user context
672 
673     Level: intermediate
674 
675 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
676 
677 @*/
678 PetscErrorCode PETSCDM_DLLEXPORT DMSetContext(DM dm,void *ctx)
679 {
680   PetscFunctionBegin;
681   dm->ctx = ctx;
682   PetscFunctionReturn(0);
683 }
684 
685 #undef __FUNCT__
686 #define __FUNCT__ "DMGetContext"
687 /*@
688     DMGetContext - Gets a user context from a DM object
689 
690     Not Collective
691 
692     Input Parameter:
693 .   dm - the DM object
694 
695     Output Parameter:
696 .   ctx - the user context
697 
698     Level: intermediate
699 
700 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext()
701 
702 @*/
703 PetscErrorCode PETSCDM_DLLEXPORT DMGetContext(DM dm,void **ctx)
704 {
705   PetscFunctionBegin;
706   *ctx = dm->ctx;
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "DMSetInitialGuess"
712 /*@
713     DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers
714 
715     Logically Collective on DM
716 
717     Input Parameter:
718 +   dm - the DM object to destroy
719 -   f - the function to compute the initial guess
720 
721     Level: intermediate
722 
723 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
724 
725 @*/
726 PetscErrorCode PETSCDM_DLLEXPORT DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec))
727 {
728   PetscFunctionBegin;
729   dm->ops->initialguess = f;
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "DMSetFunction"
735 /*@
736     DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES
737 
738     Logically Collective on DM
739 
740     Input Parameter:
741 +   dm - the DM object
742 -   f - the function to compute (use PETSC_NULL to cancel a previous function that was set)
743 
744     Level: intermediate
745 
746     Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian
747            computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian.
748 
749 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
750          DMSetJacobian()
751 
752 @*/
753 PetscErrorCode PETSCDM_DLLEXPORT DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec))
754 {
755   PetscFunctionBegin;
756   dm->ops->function = f;
757   if (f) {
758     dm->ops->functionj = f;
759   }
760   PetscFunctionReturn(0);
761 }
762 
763 #undef __FUNCT__
764 #define __FUNCT__ "DMSetJacobian"
765 /*@
766     DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES
767 
768     Logically Collective on DM
769 
770     Input Parameter:
771 +   dm - the DM object to destroy
772 -   f - the function to compute the matrix entries
773 
774     Level: intermediate
775 
776 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
777          DMSetFunction()
778 
779 @*/
780 PetscErrorCode PETSCDM_DLLEXPORT DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*))
781 {
782   PetscFunctionBegin;
783   dm->ops->jacobian = f;
784   PetscFunctionReturn(0);
785 }
786 
787 #undef __FUNCT__
788 #define __FUNCT__ "DMComputeInitialGuess"
789 /*@
790     DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers
791 
792     Collective on DM
793 
794     Input Parameter:
795 +   dm - the DM object to destroy
796 -   x - the vector to hold the initial guess values
797 
798     Level: developer
799 
800 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetRhs(), DMSetMat()
801 
802 @*/
803 PetscErrorCode PETSCDM_DLLEXPORT DMComputeInitialGuess(DM dm,Vec x)
804 {
805   PetscErrorCode ierr;
806   PetscFunctionBegin;
807   if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()");
808   ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr);
809   PetscFunctionReturn(0);
810 }
811 
812 #undef __FUNCT__
813 #define __FUNCT__ "DMHasInitialGuess"
814 /*@
815     DMHasInitialGuess - does the DM object have an initial guess function
816 
817     Not Collective
818 
819     Input Parameter:
820 .   dm - the DM object to destroy
821 
822     Output Parameter:
823 .   flg - PETSC_TRUE if function exists
824 
825     Level: developer
826 
827 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
828 
829 @*/
830 PetscErrorCode PETSCDM_DLLEXPORT DMHasInitialGuess(DM dm,PetscBool  *flg)
831 {
832   PetscFunctionBegin;
833   *flg =  (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE;
834   PetscFunctionReturn(0);
835 }
836 
837 #undef __FUNCT__
838 #define __FUNCT__ "DMHasFunction"
839 /*@
840     DMHasFunction - does the DM object have a function
841 
842     Not Collective
843 
844     Input Parameter:
845 .   dm - the DM object to destroy
846 
847     Output Parameter:
848 .   flg - PETSC_TRUE if function exists
849 
850     Level: developer
851 
852 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
853 
854 @*/
855 PetscErrorCode PETSCDM_DLLEXPORT DMHasFunction(DM dm,PetscBool  *flg)
856 {
857   PetscFunctionBegin;
858   *flg =  (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE;
859   PetscFunctionReturn(0);
860 }
861 
862 #undef __FUNCT__
863 #define __FUNCT__ "DMHasJacobian"
864 /*@
865     DMHasJacobian - does the DM object have a matrix function
866 
867     Not Collective
868 
869     Input Parameter:
870 .   dm - the DM object to destroy
871 
872     Output Parameter:
873 .   flg - PETSC_TRUE if function exists
874 
875     Level: developer
876 
877 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian()
878 
879 @*/
880 PetscErrorCode PETSCDM_DLLEXPORT DMHasJacobian(DM dm,PetscBool  *flg)
881 {
882   PetscFunctionBegin;
883   *flg =  (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE;
884   PetscFunctionReturn(0);
885 }
886 
887 #undef __FUNCT__
888 #define __FUNCT__ "DMComputeFunction"
889 /*@
890     DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES
891 
892     Collective on DM
893 
894     Input Parameter:
895 +   dm - the DM object to destroy
896 .   x - the location where the function is evaluationed, may be ignored for linear problems
897 -   b - the vector to hold the right hand side entries
898 
899     Level: developer
900 
901 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
902          DMSetJacobian()
903 
904 @*/
905 PetscErrorCode PETSCDM_DLLEXPORT DMComputeFunction(DM dm,Vec x,Vec b)
906 {
907   PetscErrorCode ierr;
908   PetscFunctionBegin;
909   if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()");
910   if (!x) x = dm->x;
911   ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr);
912   PetscFunctionReturn(0);
913 }
914 
915 
916 #undef __FUNCT__
917 #define __FUNCT__ "DMComputeJacobian"
918 /*@
919     DMComputeJacobian - compute the matrix entries for the solver
920 
921     Collective on DM
922 
923     Input Parameter:
924 +   dm - the DM object
925 .   x - location to compute Jacobian at; may be ignored for linear problems
926 .   A - matrix that defines the operator for the linear solve
927 -   B - the matrix used to construct the preconditioner
928 
929     Level: developer
930 
931 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(),
932          DMSetFunction()
933 
934 @*/
935 PetscErrorCode PETSCDM_DLLEXPORT DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag)
936 {
937   PetscErrorCode ierr;
938 
939   PetscFunctionBegin;
940   if (!dm->ops->jacobian) {
941     ISColoring    coloring;
942     MatFDColoring fd;
943 
944     ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
945     ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr);
946     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
947     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr);
948     dm->fd = fd;
949     dm->ops->jacobian = DMComputeJacobianDefault;
950 
951     if (!dm->x) {
952       ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr);
953     }
954   }
955   if (!x) x = dm->x;
956   ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr);
957   PetscFunctionReturn(0);
958 }
959 
960 PetscFList DMList                       = PETSC_NULL;
961 PetscBool  DMRegisterAllCalled          = PETSC_FALSE;
962 
963 #undef __FUNCT__
964 #define __FUNCT__ "DMSetType"
965 /*@C
966   DMSetType - Builds a DM, for a particular DM implementation.
967 
968   Collective on DM
969 
970   Input Parameters:
971 + dm     - The DM object
972 - method - The name of the DM type
973 
974   Options Database Key:
975 . -dm_type <type> - Sets the DM type; use -help for a list of available types
976 
977   Notes:
978   See "petsc/include/petscda.h" for available DM types (for instance, DM1D, DM2D, or DM3D).
979 
980   Level: intermediate
981 
982 .keywords: DM, set, type
983 .seealso: DMGetType(), DMCreate()
984 @*/
985 PetscErrorCode PETSCDM_DLLEXPORT DMSetType(DM dm, const DMType method)
986 {
987   PetscErrorCode (*r)(DM);
988   PetscBool      match;
989   PetscErrorCode ierr;
990 
991   PetscFunctionBegin;
992   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
993   ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr);
994   if (match) PetscFunctionReturn(0);
995 
996   if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
997   ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,(void (**)(void)) &r);CHKERRQ(ierr);
998   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method);
999 
1000   if (dm->ops->destroy) {
1001     ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr);
1002   }
1003   ierr = (*r)(dm);CHKERRQ(ierr);
1004   ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "DMGetType"
1010 /*@C
1011   DMGetType - Gets the DM type name (as a string) from the DM.
1012 
1013   Not Collective
1014 
1015   Input Parameter:
1016 . dm  - The DM
1017 
1018   Output Parameter:
1019 . type - The DM type name
1020 
1021   Level: intermediate
1022 
1023 .keywords: DM, get, type, name
1024 .seealso: DMSetType(), DMCreate()
1025 @*/
1026 PetscErrorCode PETSCDM_DLLEXPORT DMGetType(DM dm, const DMType *type)
1027 {
1028   PetscErrorCode ierr;
1029 
1030   PetscFunctionBegin;
1031   PetscValidHeaderSpecific(dm, DM_CLASSID,1);
1032   PetscValidCharPointer(type,2);
1033   if (!DMRegisterAllCalled) {
1034     ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1035   }
1036   *type = ((PetscObject)dm)->type_name;
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 
1041 /*--------------------------------------------------------------------------------------------------------------------*/
1042 
1043 #undef __FUNCT__
1044 #define __FUNCT__ "DMRegister"
1045 /*@C
1046   DMRegister - See DMRegisterDynamic()
1047 
1048   Level: advanced
1049 @*/
1050 PetscErrorCode PETSCDM_DLLEXPORT DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM))
1051 {
1052   char fullname[PETSC_MAX_PATH_LEN];
1053   PetscErrorCode ierr;
1054 
1055   PetscFunctionBegin;
1056   ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr);
1057   ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr);
1058   ierr = PetscStrcat(fullname, name);CHKERRQ(ierr);
1059   ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 
1064 /*--------------------------------------------------------------------------------------------------------------------*/
1065 #undef __FUNCT__
1066 #define __FUNCT__ "DMRegisterDestroy"
1067 /*@C
1068    DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic().
1069 
1070    Not Collective
1071 
1072    Level: advanced
1073 
1074 .keywords: DM, register, destroy
1075 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic()
1076 @*/
1077 PetscErrorCode PETSCDM_DLLEXPORT DMRegisterDestroy(void)
1078 {
1079   PetscErrorCode ierr;
1080 
1081   PetscFunctionBegin;
1082   ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr);
1083   DMRegisterAllCalled = PETSC_FALSE;
1084   PetscFunctionReturn(0);
1085 }
1086