xref: /petsc/src/snes/utils/dmsnes.c (revision 047240e14af00aad1ef65e96f6fface8924f7f7e)
1 #include <petsc-private/snesimpl.h>   /*I "petscsnes.h" I*/
2 #include <petsc-private/dmimpl.h>     /*I "petscdm.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMSNESDestroy"
6 static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
7 {
8   PetscErrorCode ierr;
9 
10   PetscFunctionBegin;
11   if (!*kdm) PetscFunctionReturn(0);
12   PetscValidHeaderSpecific((*kdm),DMSNES_CLASSID,1);
13   if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);}
14   if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);}
15   ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr);
16   PetscFunctionReturn(0);
17 }
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "DMSNESLoad"
21 PetscErrorCode DMSNESLoad(DMSNES kdm,PetscViewer viewer)
22 {
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->computefunction,1,PETSC_FUNCTION);CHKERRQ(ierr);
27   ierr = PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,PETSC_FUNCTION);CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "DMSNESView"
33 PetscErrorCode DMSNESView(DMSNES kdm,PetscViewer viewer)
34 {
35   PetscErrorCode ierr;
36   PetscBool      isascii,isbinary;
37 
38   PetscFunctionBegin;
39   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
40   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
41   if (isascii) {
42 #if defined(PETSC_SERIALIZE_FUNCTIONS)
43     const char *fname;
44 
45     ierr = PetscFPTFind(kdm->ops->computefunction,&fname);CHKERRQ(ierr);
46     if (fname) {
47       ierr = PetscViewerASCIIPrintf(viewer,"Function used by SNES: %s\n",fname);CHKERRQ(ierr);
48     }
49     ierr = PetscFPTFind(kdm->ops->computejacobian,&fname);CHKERRQ(ierr);
50     if (fname) {
51       ierr = PetscViewerASCIIPrintf(viewer,"Jacobian function used by SNES: %s\n",fname);CHKERRQ(ierr);
52     }
53 #endif
54   } else if (isbinary) {
55     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->computefunction,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
56     ierr = PetscViewerBinaryWrite(viewer,(void*)kdm->ops->computejacobian,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr);
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "DMSNESCreate"
63 static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm)
64 {
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
69   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
70 #endif
71   ierr = PetscHeaderCreate(*kdm, _p_DMSNES, struct _DMSNESOps, DMSNES_CLASSID, -1, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);CHKERRQ(ierr);
72   ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMSNESOps));CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "DMCoarsenHook_DMSNES"
78 /* Attaches the DMSNES to the coarse level.
79  * Under what conditions should we copy versus duplicate?
80  */
81 static PetscErrorCode DMCoarsenHook_DMSNES(DM dm,DM dmc,void *ctx)
82 {
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr = DMCopyDMSNES(dm,dmc);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "DMRestrictHook_DMSNES"
92 /* This could restrict auxiliary information to the coarse level.
93  */
94 static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
95 {
96 
97   PetscFunctionBegin;
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "DMSubDomainHook_DMSNES"
103 /* Attaches the DMSNES to the subdomain. */
104 static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx)
105 {
106   PetscErrorCode ierr;
107 
108   PetscFunctionBegin;
109   ierr = DMCopyDMSNES(dm,subdm);CHKERRQ(ierr);
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNCT__
114 #define __FUNCT__ "DMSubDomainRestrictHook_DMSNES"
115 /* This could restrict auxiliary information to the coarse level.
116  */
117 static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
118 {
119 
120   PetscFunctionBegin;
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "DMRefineHook_DMSNES"
126 static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx)
127 {
128   PetscErrorCode ierr;
129 
130   PetscFunctionBegin;
131   ierr = DMCopyDMSNES(dm,dmf);CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "DMInterpolateHook_DMSNES"
137 /* This could restrict auxiliary information to the coarse level.
138  */
139 static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx)
140 {
141 
142   PetscFunctionBegin;
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "DMSNESCopy"
148 /*@C
149    DMSNESCopy - copies the information in a DMSNES to another DMSNES
150 
151    Not Collective
152 
153    Input Argument:
154 +  kdm - Original DMSNES
155 -  nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()
156 
157    Level: developer
158 
159 .seealso: DMSNESCreate(), DMSNESDestroy()
160 @*/
161 PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm)
162 {
163   PetscErrorCode ierr;
164 
165   PetscFunctionBegin;
166   PetscValidHeaderSpecific(kdm,DMSNES_CLASSID,1);
167   PetscValidHeaderSpecific(nkdm,DMSNES_CLASSID,2);
168   nkdm->ops->computefunction  = kdm->ops->computefunction;
169   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
170   nkdm->ops->computegs        = kdm->ops->computegs;
171   nkdm->ops->computeobjective = kdm->ops->computeobjective;
172   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
173   nkdm->ops->computepfunction = kdm->ops->computepfunction;
174   nkdm->ops->destroy          = kdm->ops->destroy;
175   nkdm->ops->duplicate        = kdm->ops->duplicate;
176 
177   nkdm->functionctx  = kdm->functionctx;
178   nkdm->gsctx        = kdm->gsctx;
179   nkdm->pctx         = kdm->pctx;
180   nkdm->jacobianctx  = kdm->jacobianctx;
181   nkdm->objectivectx = kdm->objectivectx;
182   nkdm->data         = kdm->data;
183 
184   /*
185   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
186   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
187   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
188   */
189 
190   /* implementation specific copy hooks */
191   if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);}
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "DMGetDMSNES"
197 /*@C
198    DMGetDMSNES - get read-only private DMSNES context from a DM
199 
200    Not Collective
201 
202    Input Argument:
203 .  dm - DM to be used with SNES
204 
205    Output Argument:
206 .  snesdm - private DMSNES context
207 
208    Level: developer
209 
210    Notes:
211    Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
212 
213 .seealso: DMGetDMSNESWrite()
214 @*/
215 PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
216 {
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
221   *snesdm = (DMSNES) dm->dmsnes;
222   if (!*snesdm) {
223     ierr = PetscInfo(dm,"Creating new DMSNES\n");CHKERRQ(ierr);
224     ierr = DMSNESCreate(((PetscObject)dm)->comm,snesdm);CHKERRQ(ierr);
225 
226     dm->dmsnes = (PetscObject) *snesdm;
227 
228     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,PETSC_NULL);CHKERRQ(ierr);
229     ierr = DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,PETSC_NULL);CHKERRQ(ierr);
230     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,PETSC_NULL);CHKERRQ(ierr);
231   }
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "DMGetDMSNESWrite"
237 /*@C
238    DMGetDMSNESWrite - get write access to private DMSNES context from a DM
239 
240    Not Collective
241 
242    Input Argument:
243 .  dm - DM to be used with SNES
244 
245    Output Argument:
246 .  snesdm - private DMSNES context
247 
248    Level: developer
249 
250 .seealso: DMGetDMSNES()
251 @*/
252 PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
253 {
254   PetscErrorCode ierr;
255   DMSNES         sdm;
256 
257   PetscFunctionBegin;
258   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
259   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
260   if (!sdm->originaldm) sdm->originaldm = dm;
261   if (sdm->originaldm != dm) {  /* Copy on write */
262     DMSNES oldsdm = sdm;
263     ierr       = PetscInfo(dm,"Copying DMSNES due to write\n");CHKERRQ(ierr);
264     ierr       = DMSNESCreate(((PetscObject)dm)->comm,&sdm);CHKERRQ(ierr);
265     ierr       = DMSNESCopy(oldsdm,sdm);CHKERRQ(ierr);
266     ierr       = DMSNESDestroy((DMSNES*)&dm->dmsnes);CHKERRQ(ierr);
267     dm->dmsnes = (PetscObject)sdm;
268   }
269   *snesdm = sdm;
270   PetscFunctionReturn(0);
271 }
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "DMCopyDMSNES"
275 /*@C
276    DMCopyDMSNES - copies a DM context to a new DM
277 
278    Logically Collective
279 
280    Input Arguments:
281 +  dmsrc - DM to obtain context from
282 -  dmdest - DM to add context to
283 
284    Level: developer
285 
286    Note:
287    The context is copied by reference. This function does not ensure that a context exists.
288 
289 .seealso: DMGetDMSNES(), SNESSetDM()
290 @*/
291 PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
292 {
293   PetscErrorCode ierr;
294 
295   PetscFunctionBegin;
296   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
297   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
298   ierr = DMSNESDestroy((DMSNES*)&dmdest->dmsnes);CHKERRQ(ierr);
299 
300   dmdest->dmsnes = dmsrc->dmsnes;
301 
302   ierr = PetscObjectReference(dmdest->dmsnes);CHKERRQ(ierr);
303   ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
304   ierr = DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
305   ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,PETSC_NULL);CHKERRQ(ierr);
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "DMSNESSetFunction"
311 /*@C
312    DMSNESSetFunction - set SNES residual evaluation function
313 
314    Not Collective
315 
316    Input Arguments:
317 +  dm - DM to be used with SNES
318 .  SNESFunction - residual evaluation function
319 -  ctx - context for residual evaluation
320 
321    Level: advanced
322 
323    Note:
324    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
325    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
326    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
327 
328 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
329 @*/
330 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*),void *ctx)
331 {
332   PetscErrorCode ierr;
333   DMSNES         sdm;
334 
335   PetscFunctionBegin;
336   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
337   if (SNESFunction || ctx) {
338     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
339   }
340   if (SNESFunction) sdm->ops->computefunction = SNESFunction;
341   if (ctx) sdm->functionctx = ctx;
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "DMSNESGetFunction"
347 /*@C
348    DMSNESGetFunction - get SNES residual evaluation function
349 
350    Not Collective
351 
352    Input Argument:
353 .  dm - DM to be used with SNES
354 
355    Output Arguments:
356 +  SNESFunction - residual evaluation function
357 -  ctx - context for residual evaluation
358 
359    Level: advanced
360 
361    Note:
362    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
363    associated with the DM.
364 
365 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
366 @*/
367 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**SNESFunction)(SNES,Vec,Vec,void*),void **ctx)
368 {
369   PetscErrorCode ierr;
370   DMSNES         sdm;
371 
372   PetscFunctionBegin;
373   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
374   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
375   if (SNESFunction) *SNESFunction = sdm->ops->computefunction;
376   if (ctx) *ctx = sdm->functionctx;
377   PetscFunctionReturn(0);
378 }
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "DMSNESSetObjective"
382 /*@C
383    DMSNESSetObjective - set SNES objective evaluation function
384 
385    Not Collective
386 
387    Input Arguments:
388 +  dm - DM to be used with SNES
389 .  SNESObjectiveFunction - residual evaluation function
390 -  ctx - context for residual evaluation
391 
392    Level: advanced
393 
394 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
395 @*/
396 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*SNESObjectiveFunction)(SNES,Vec,PetscReal*,void*),void *ctx)
397 {
398   PetscErrorCode ierr;
399   DMSNES         sdm;
400 
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
403   if (SNESObjectiveFunction || ctx) {
404     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
405   }
406   if (SNESObjectiveFunction) sdm->ops->computeobjective = SNESObjectiveFunction;
407   if (ctx) sdm->objectivectx = ctx;
408   PetscFunctionReturn(0);
409 }
410 
411 #undef __FUNCT__
412 #define __FUNCT__ "DMSNESGetObjective"
413 /*@C
414    DMSNESGetObjective - get SNES objective evaluation function
415 
416    Not Collective
417 
418    Input Argument:
419 .  dm - DM to be used with SNES
420 
421    Output Arguments:
422 +  SNESObjectiveFunction- residual evaluation function
423 -  ctx - context for residual evaluation
424 
425    Level: advanced
426 
427    Note:
428    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
429    associated with the DM.
430 
431 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
432 @*/
433 PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**SNESObjectiveFunction)(SNES,Vec,PetscReal*,void*),void **ctx)
434 {
435   PetscErrorCode ierr;
436   DMSNES         sdm;
437 
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
440   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
441   if (SNESObjectiveFunction) *SNESObjectiveFunction = sdm->ops->computeobjective;
442   if (ctx) *ctx = sdm->objectivectx;
443   PetscFunctionReturn(0);
444 }
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "DMSNESSetGS"
448 /*@C
449    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
450 
451    Not Collective
452 
453    Input Argument:
454 +  dm - DM to be used with SNES
455 .  SNESGSFunction - relaxation function
456 -  ctx - context for residual evaluation
457 
458    Level: advanced
459 
460    Note:
461    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
462    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
463    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
464 
465 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
466 @*/
467 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*SNESGSFunction)(SNES,Vec,Vec,void*),void *ctx)
468 {
469   PetscErrorCode ierr;
470   DMSNES         sdm;
471 
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
474   if (SNESGSFunction || ctx) {
475     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
476   }
477   if (SNESGSFunction) sdm->ops->computegs = SNESGSFunction;
478   if (ctx) sdm->gsctx = ctx;
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "DMSNESGetGS"
484 /*@C
485    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
486 
487    Not Collective
488 
489    Input Argument:
490 .  dm - DM to be used with SNES
491 
492    Output Arguments:
493 +  SNESGSFunction - relaxation function which performs Gauss-Seidel sweeps
494 -  ctx - context for residual evaluation
495 
496    Level: advanced
497 
498    Note:
499    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
500    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
501    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
502 
503 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESGSFunction
504 @*/
505 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**SNESGSFunction)(SNES,Vec,Vec,void*),void **ctx)
506 {
507   PetscErrorCode ierr;
508   DMSNES         sdm;
509 
510   PetscFunctionBegin;
511   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
512   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
513   if (SNESGSFunction) *SNESGSFunction = sdm->ops->computegs;
514   if (ctx) *ctx = sdm->gsctx;
515   PetscFunctionReturn(0);
516 }
517 
518 #undef __FUNCT__
519 #define __FUNCT__ "DMSNESSetJacobian"
520 /*@C
521    DMSNESSetJacobian - set SNES Jacobian evaluation function
522 
523    Not Collective
524 
525    Input Argument:
526 +  dm - DM to be used with SNES
527 .  SNESJacobianFunction - Jacobian evaluation function
528 -  ctx - context for residual evaluation
529 
530    Level: advanced
531 
532    Note:
533    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
534    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
535    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
536 
537 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
538 @*/
539 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
540 {
541   PetscErrorCode ierr;
542   DMSNES         sdm;
543 
544   PetscFunctionBegin;
545   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
546   if (SNESJacobianFunction || ctx) {
547     ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
548   }
549   if (SNESJacobianFunction) sdm->ops->computejacobian = SNESJacobianFunction;
550   if (ctx) sdm->jacobianctx = ctx;
551   PetscFunctionReturn(0);
552 }
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "DMSNESGetJacobian"
556 /*@C
557    DMSNESGetJacobian - get SNES Jacobian evaluation function
558 
559    Not Collective
560 
561    Input Argument:
562 .  dm - DM to be used with SNES
563 
564    Output Arguments:
565 +  SNESJacobianFunction - Jacobian evaluation function
566 -  ctx - context for residual evaluation
567 
568    Level: advanced
569 
570    Note:
571    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
572    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
573    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
574 
575 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
576 @*/
577 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
578 {
579   PetscErrorCode ierr;
580   DMSNES         sdm;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
584   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
585   if (SNESJacobianFunction) *SNESJacobianFunction = sdm->ops->computejacobian;
586   if (ctx) *ctx = sdm->jacobianctx;
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "DMSNESSetPicard"
592 /*@C
593    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
594 
595    Not Collective
596 
597    Input Argument:
598 +  dm - DM to be used with SNES
599 .  SNESFunction - RHS evaluation function
600 .  SNESJacobianFunction - Picard matrix evaluation function
601 -  ctx - context for residual evaluation
602 
603    Level: advanced
604 
605 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
606 @*/
607 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*),PetscErrorCode (*SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
608 {
609   PetscErrorCode ierr;
610   DMSNES         sdm;
611 
612   PetscFunctionBegin;
613   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
614   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
615   if (SNESFunction) sdm->ops->computepfunction = SNESFunction;
616   if (SNESJacobianFunction) sdm->ops->computepjacobian = SNESJacobianFunction;
617   if (ctx) sdm->pctx = ctx;
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "DMSNESGetPicard"
623 /*@C
624    DMSNESGetPicard - get SNES Picard iteration evaluation functions
625 
626    Not Collective
627 
628    Input Argument:
629 .  dm - DM to be used with SNES
630 
631    Output Arguments:
632 +  SNESFunction - Jacobian evaluation function;
633 .  SNESJacobianFunction  - RHS evaluation function;
634 -  ctx - context for residual evaluation
635 
636    Level: advanced
637 
638 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
639 @*/
640 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**SNESFunction)(SNES,Vec,Vec,void*),PetscErrorCode (**SNESJacobianFunction)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
641 {
642   PetscErrorCode ierr;
643   DMSNES         sdm;
644 
645   PetscFunctionBegin;
646   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
647   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
648   if (SNESFunction) *SNESFunction = sdm->ops->computepfunction;
649   if (SNESJacobianFunction) *SNESJacobianFunction = sdm->ops->computepjacobian;
650   if (ctx) *ctx = sdm->pctx;
651   PetscFunctionReturn(0);
652 }
653