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