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