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