xref: /petsc/src/snes/utils/dmsnes.c (revision 0d4d505e9eeaf2031c0432dfe494cb21852b352d)
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__ "SNESDMComputeFunction"
6 /*@C
7   SNESDMComputeFunction - This is a universal function evaluation routine that
8   may be used with SNESSetFunction() as long as the user context has a DM
9   as its first record and the user has called DMSetLocalFunction().
10 
11   Collective on SNES
12 
13   Input Parameters:
14 + snes - the SNES context
15 . X - input vector
16 . F - function vector
17 - ptr - pointer to a structure that must have a DM as its first entry.
18         This ptr must have been passed into SNESDMComputeFunction() as the context.
19 
20   Level: intermediate
21 
22 .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
23 @*/
24 PetscErrorCode SNESDMComputeFunction(SNES snes, Vec X, Vec F, void *ptr)
25 {
26   DM               dm = *(DM*) ptr;
27   PetscErrorCode (*lf)(DM, Vec, Vec, void *);
28   Vec              localX, localF;
29   PetscInt         N, n;
30   PetscErrorCode   ierr;
31 
32   PetscFunctionBegin;
33   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
34   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
35   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
36   if (!dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Looks like you called SNESSetFromFuntion(snes,SNESDMComputeFunction,) without the DM context");
37   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
38 
39   /* determine whether X = localX */
40   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
41   ierr = DMGetLocalVector(dm, &localF);CHKERRQ(ierr);
42   ierr = VecGetSize(X, &N);CHKERRQ(ierr);
43   ierr = VecGetSize(localX, &n);CHKERRQ(ierr);
44 
45   if (n != N){ /* X != localX */
46     /* Scatter ghost points to local vector, using the 2-step process
47         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
48     */
49     ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
50     ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
51   } else {
52     ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
53     localX = X;
54   }
55   ierr = DMGetLocalFunction(dm, &lf);CHKERRQ(ierr);
56   ierr = (*lf)(dm, localX, localF, ptr);CHKERRQ(ierr);
57   if (n != N){
58     ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
59   }
60   ierr = DMLocalToGlobalBegin(dm, localF, ADD_VALUES, F);CHKERRQ(ierr);
61   ierr = DMLocalToGlobalEnd(dm, localF, ADD_VALUES, F);CHKERRQ(ierr);
62   ierr = DMRestoreLocalVector(dm, &localF);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "SNESDMComputeJacobian"
68 /*
69   SNESDMComputeJacobian - This is a universal Jacobian evaluation routine that
70   may be used with SNESSetJacobian() as long as the user context has a DM
71   as its first record and the user has called DMSetLocalJacobian().
72 
73   Collective on SNES
74 
75   Input Parameters:
76 + snes - the SNES context
77 . X - input vector
78 . J - Jacobian
79 . B - Jacobian used in preconditioner (usally same as J)
80 . flag - indicates if the matrix changed its structure
81 - ptr - pointer to a structure that must have a DM as its first entry.
82         This ptr must have been passed into SNESDMComputeFunction() as the context.
83 
84   Level: intermediate
85 
86 .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()
87 */
88 PetscErrorCode SNESDMComputeJacobian(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr)
89 {
90   DM               dm = *(DM*) ptr;
91   PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *);
92   Vec              localX;
93   PetscErrorCode   ierr;
94 
95   PetscFunctionBegin;
96   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
97   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
98   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99   ierr = DMGetLocalJacobian(dm, &lj);CHKERRQ(ierr);
100   ierr = (*lj)(dm, localX, *J, *B, ptr);CHKERRQ(ierr);
101   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
102   /* Assemble true Jacobian; if it is different */
103   if (*J != *B) {
104     ierr  = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105     ierr  = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106   }
107   ierr  = MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr);
108   *flag = SAME_NONZERO_PATTERN;
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "PetscContainerDestroy_SNESDM"
114 static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx)
115 {
116   PetscErrorCode ierr;
117   SNESDM sdm = (SNESDM)ctx;
118 
119   PetscFunctionBegin;
120   if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);}
121   ierr = PetscFree(sdm);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "DMCoarsenHook_SNESDM"
127 /* Attaches the SNESDM to the coarse level.
128  * Under what conditions should we copy versus duplicate?
129  */
130 static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx)
131 {
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "DMRestrictHook_SNESDM"
141 /* This could restrict auxiliary information to the coarse level.
142  */
143 static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
144 {
145 
146   PetscFunctionBegin;
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "DMRefineHook_SNESDM"
152 static PetscErrorCode DMRefineHook_SNESDM(DM dm,DM dmf,void *ctx)
153 {
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   ierr = DMSNESCopyContext(dm,dmf);CHKERRQ(ierr);
158   PetscFunctionReturn(0);
159 }
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "DMInterpolateHook_SNESDM"
163 /* This could restrict auxiliary information to the coarse level.
164  */
165 static PetscErrorCode DMInterpolateHook_SNESDM(DM dm,Mat Interp,DM dmf,void *ctx)
166 {
167 
168   PetscFunctionBegin;
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "DMSNESGetContext"
174 /*@C
175    DMSNESGetContext - get read-only private SNESDM context from a DM
176 
177    Not Collective
178 
179    Input Argument:
180 .  dm - DM to be used with SNES
181 
182    Output Argument:
183 .  snesdm - private SNESDM context
184 
185    Level: developer
186 
187    Notes:
188    Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
189 
190 .seealso: DMSNESGetContextWrite()
191 @*/
192 PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm)
193 {
194   PetscErrorCode ierr;
195   PetscContainer container;
196   SNESDM         sdm;
197 
198   PetscFunctionBegin;
199   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
200   ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
201   if (container) {
202     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
203   } else {
204     ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr);
205     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
206     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
207     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
208     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
209     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
210     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
211     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
212     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "DMSNESGetContextWrite"
219 /*@C
220    DMSNESGetContextWrite - get write access to private SNESDM context from a DM
221 
222    Not Collective
223 
224    Input Argument:
225 .  dm - DM to be used with SNES
226 
227    Output Argument:
228 .  snesdm - private SNESDM context
229 
230    Level: developer
231 
232 .seealso: DMSNESGetContext()
233 @*/
234 PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm)
235 {
236   PetscErrorCode ierr;
237   SNESDM         sdm;
238 
239   PetscFunctionBegin;
240   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
241   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
242   if (!sdm->originaldm) sdm->originaldm = dm;
243   if (sdm->originaldm != dm) {  /* Copy on write */
244     PetscContainer container;
245     SNESDM         oldsdm = sdm;
246     ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr);
247     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
248     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
249     ierr = PetscMemcpy(sdm,oldsdm,sizeof(*sdm));CHKERRQ(ierr);
250     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
251     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
252     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
253     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
254     /* implementation specific copy hooks */
255     ierr = (sdm->duplicate)(oldsdm,dm);CHKERRQ(ierr);
256   }
257   *snesdm = sdm;
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "DMSNESCopyContext"
263 /*@C
264    DMSNESCopyContext - copies a DM context to a new DM
265 
266    Logically Collective
267 
268    Input Arguments:
269 +  dmsrc - DM to obtain context from
270 -  dmdest - DM to add context to
271 
272    Level: developer
273 
274    Note:
275    The context is copied by reference. This function does not ensure that a context exists.
276 
277 .seealso: DMSNESGetContext(), SNESSetDM()
278 @*/
279 PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest)
280 {
281   PetscErrorCode ierr;
282   PetscContainer container;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
286   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
287   ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
288   if (container) {
289     ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
290     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
291     ierr = DMRefineHookAdd(dmdest,DMRefineHook_SNESDM,DMInterpolateHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
292   }
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "DMSNESSetFunction"
298 /*@C
299    DMSNESSetFunction - set SNES residual evaluation function
300 
301    Not Collective
302 
303    Input Arguments:
304 +  dm - DM to be used with SNES
305 .  func - residual evaluation function, see SNESSetFunction() for calling sequence
306 -  ctx - context for residual evaluation
307 
308    Level: advanced
309 
310    Note:
311    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
312    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
313    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
314 
315 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
316 @*/
317 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
318 {
319   PetscErrorCode ierr;
320   SNESDM sdm;
321 
322   PetscFunctionBegin;
323   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
324   if (func || ctx) {
325     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
326   }
327   if (func) sdm->computefunction = func;
328   if (ctx)  sdm->functionctx = ctx;
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "DMSNESGetFunction"
334 /*@C
335    DMSNESGetFunction - get SNES residual evaluation function
336 
337    Not Collective
338 
339    Input Argument:
340 .  dm - DM to be used with SNES
341 
342    Output Arguments:
343 +  func - residual evaluation function, see SNESSetFunction() for calling sequence
344 -  ctx - context for residual evaluation
345 
346    Level: advanced
347 
348    Note:
349    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
350    associated with the DM.
351 
352 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
353 @*/
354 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
355 {
356   PetscErrorCode ierr;
357   SNESDM sdm;
358 
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
361   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
362   if (func) *func = sdm->computefunction;
363   if (ctx)  *ctx = sdm->functionctx;
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "DMSNESSetObjective"
369 /*@C
370    DMSNESSetObjective - set SNES objective evaluation function
371 
372    Not Collective
373 
374    Input Arguments:
375 +  dm - DM to be used with SNES
376 .  func - residual evaluation function, see SNESSetObjective() for calling sequence
377 -  ctx - context for residual evaluation
378 
379    Level: advanced
380 
381 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
382 @*/
383 PetscErrorCode DMSNESSetObjective(DM dm,SNESObjective func,void *ctx)
384 {
385   PetscErrorCode ierr;
386   SNESDM sdm;
387 
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
390   if (func || ctx) {
391     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
392   }
393   if (func) sdm->computeobjective = func;
394   if (ctx)  sdm->objectivectx = ctx;
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "DMSNESGetObjective"
400 /*@C
401    DMSNESGetObjective - get SNES objective evaluation function
402 
403    Not Collective
404 
405    Input Argument:
406 .  dm - DM to be used with SNES
407 
408    Output Arguments:
409 +  func - residual evaluation function, see SNESSetObjective() for calling sequence
410 -  ctx - context for residual evaluation
411 
412    Level: advanced
413 
414    Note:
415    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
416    associated with the DM.
417 
418 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
419 @*/
420 PetscErrorCode DMSNESGetObjective(DM dm,SNESObjective *func,void **ctx)
421 {
422   PetscErrorCode ierr;
423   SNESDM sdm;
424 
425   PetscFunctionBegin;
426   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
427   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
428   if (func) *func = sdm->computeobjective;
429   if (ctx)  *ctx = sdm->objectivectx;
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "DMSNESSetGS"
435 /*@C
436    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
437 
438    Not Collective
439 
440    Input Argument:
441 +  dm - DM to be used with SNES
442 .  func - relaxation function, see SNESSetGS() for calling sequence
443 -  ctx - context for residual evaluation
444 
445    Level: advanced
446 
447    Note:
448    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
449    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
450    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
451 
452 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction()
453 @*/
454 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
455 {
456   PetscErrorCode ierr;
457   SNESDM sdm;
458 
459   PetscFunctionBegin;
460   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
461   if (func || ctx) {
462     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
463   }
464   if (func) sdm->computegs = func;
465   if (ctx)  sdm->gsctx = ctx;
466   PetscFunctionReturn(0);
467 }
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "DMSNESGetGS"
471 /*@C
472    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
473 
474    Not Collective
475 
476    Input Argument:
477 .  dm - DM to be used with SNES
478 
479    Output Arguments:
480 +  func - relaxation function, see SNESSetGS() for calling sequence
481 -  ctx - context for residual evaluation
482 
483    Level: advanced
484 
485    Note:
486    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
487    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
488    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
489 
490 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction()
491 @*/
492 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
493 {
494   PetscErrorCode ierr;
495   SNESDM sdm;
496 
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
499   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
500   if (func) *func = sdm->computegs;
501   if (ctx)  *ctx = sdm->gsctx;
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "DMSNESSetJacobian"
507 /*@C
508    DMSNESSetJacobian - set SNES Jacobian evaluation function
509 
510    Not Collective
511 
512    Input Argument:
513 +  dm - DM to be used with SNES
514 .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
515 -  ctx - context for residual evaluation
516 
517    Level: advanced
518 
519    Note:
520    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
521    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
522    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
523 
524 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
525 @*/
526 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
527 {
528   PetscErrorCode ierr;
529   SNESDM         sdm;
530 
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
533   if (func || ctx) {
534     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
535   }
536   if (func) sdm->computejacobian = func;
537   if (ctx)  sdm->jacobianctx = ctx;
538   PetscFunctionReturn(0);
539 }
540 
541 #undef __FUNCT__
542 #define __FUNCT__ "DMSNESGetJacobian"
543 /*@C
544    DMSNESGetJacobian - get SNES Jacobian evaluation function
545 
546    Not Collective
547 
548    Input Argument:
549 .  dm - DM to be used with SNES
550 
551    Output Arguments:
552 +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
553 -  ctx - context for residual evaluation
554 
555    Level: advanced
556 
557    Note:
558    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
559    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
560    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
561 
562 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
563 @*/
564 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
565 {
566   PetscErrorCode ierr;
567   SNESDM sdm;
568 
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
571   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
572   if (func) *func = sdm->computejacobian;
573   if (ctx)  *ctx = sdm->jacobianctx;
574   PetscFunctionReturn(0);
575 }
576 
577 #undef __FUNCT__
578 #define __FUNCT__ "DMSNESSetPicard"
579 /*@C
580    DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
581 
582    Not Collective
583 
584    Input Argument:
585 +  dm - DM to be used with SNES
586 .  func - RHS evaluation function, see SNESSetFunction() for calling sequence
587 .  pjac - Picard matrix evaluation function, see SNESSetJacobian() for calling sequence
588 -  ctx - context for residual evaluation
589 
590    Level: advanced
591 
592 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
593 @*/
594 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (*pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
595 {
596   PetscErrorCode ierr;
597   SNESDM sdm;
598 
599   PetscFunctionBegin;
600   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
601   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
602   if (pfunc) sdm->computepfunction = pfunc;
603   if (pjac)  sdm->computepjacobian = pjac;
604   if (ctx)   sdm->pctx             = ctx;
605   PetscFunctionReturn(0);
606 }
607 
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "DMSNESGetPicard"
611 /*@C
612    DMSNESGetPicard - get SNES Picard iteration evaluation functions
613 
614    Not Collective
615 
616    Input Argument:
617 .  dm - DM to be used with SNES
618 
619    Output Arguments:
620 +  pfunc - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
621 .  pjac  - RHS evaluation function, see SNESSetFunction() for calling sequence
622 -  ctx - context for residual evaluation
623 
624    Level: advanced
625 
626 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
627 @*/
628 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (**pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
629 {
630   PetscErrorCode ierr;
631   SNESDM sdm;
632 
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
635   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
636   if (pfunc) *pfunc = sdm->computepfunction;
637   if (pjac) *pjac   = sdm->computepjacobian;
638   if (ctx)  *ctx    = sdm->pctx;
639   PetscFunctionReturn(0);
640 }
641 
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy"
645 static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx)
646 {
647   PetscErrorCode ierr;
648   DM             dm;
649 
650   PetscFunctionBegin;
651   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
652   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy"
658 static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
659 {
660   PetscErrorCode ierr;
661   DM             dm;
662 
663   PetscFunctionBegin;
664   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
665   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "DMSNESSetUpLegacy"
671 /* Sets up calling of legacy DM routines */
672 PetscErrorCode DMSNESSetUpLegacy(DM dm)
673 {
674   PetscErrorCode ierr;
675   SNESDM         sdm;
676 
677   PetscFunctionBegin;
678   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
679   if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
680   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
681   if (!sdm->computejacobian) {
682     if (dm->ops->functionj) {
683       ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);
684     } else {
685       ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr);
686     }
687   }
688   PetscFunctionReturn(0);
689 }
690 
691 /* block functions */
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "DMSNESSetBlockFunction"
695 /*@C
696    DMSNESSetBlockFunction - set SNES residual evaluation function
697 
698    Not Collective
699 
700    Input Arguments:
701 +  dm - DM to be used with SNES
702 .  func - residual evaluation function, see SNESSetFunction() for calling sequence
703 -  ctx - context for residual evaluation
704 
705    Level: developer
706 
707    Note:
708    Mostly for use in DM implementations and transferred to a block function rather than being called from here.
709 
710 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
711 @*/
712 PetscErrorCode DMSNESSetBlockFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
713 {
714   PetscErrorCode ierr;
715   SNESDM sdm;
716 
717   PetscFunctionBegin;
718   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
719   if (func || ctx) {
720     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
721   }
722   if (func) sdm->computeblockfunction = func;
723   if (ctx)  sdm->blockfunctionctx = ctx;
724   PetscFunctionReturn(0);
725 }
726 
727 #undef __FUNCT__
728 #define __FUNCT__ "DMSNESGetBlockFunction"
729 /*@C
730    DMSNESGetBlockFunction - get SNES residual evaluation function
731 
732    Not Collective
733 
734    Input Argument:
735 .  dm - DM to be used with SNES
736 
737    Output Arguments:
738 +  func - residual evaluation function, see SNESSetFunction() for calling sequence
739 -  ctx - context for residual evaluation
740 
741    Level: developer
742 
743 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
744 @*/
745 PetscErrorCode DMSNESGetBlockFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
746 {
747   PetscErrorCode ierr;
748   SNESDM sdm;
749 
750   PetscFunctionBegin;
751   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
752   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
753   if (func) *func = sdm->computeblockfunction;
754   if (ctx)  *ctx = sdm->blockfunctionctx;
755   PetscFunctionReturn(0);
756 }
757 
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "DMSNESSetBlockJacobian"
761 /*@C
762    DMSNESSetJacobian - set SNES Jacobian evaluation function
763 
764    Not Collective
765 
766    Input Argument:
767 +  dm - DM to be used with SNES
768 .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
769 -  ctx - context for residual evaluation
770 
771    Level: advanced
772 
773    Note:
774    Mostly for use in DM implementations and transferred to a block function rather than being called from here.
775 
776 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
777 @*/
778 PetscErrorCode DMSNESSetBlockJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
779 {
780   PetscErrorCode ierr;
781   SNESDM sdm;
782 
783   PetscFunctionBegin;
784   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
785   if (func || ctx) {
786     ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
787   }
788   if (func) sdm->computeblockjacobian = func;
789   if (ctx)  sdm->blockjacobianctx = ctx;
790   PetscFunctionReturn(0);
791 }
792 
793 #undef __FUNCT__
794 #define __FUNCT__ "DMSNESGetBlockJacobian"
795 /*@C
796    DMSNESGetBlockJacobian - get SNES Jacobian evaluation function
797 
798    Not Collective
799 
800    Input Argument:
801 .  dm - DM to be used with SNES
802 
803    Output Arguments:
804 +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
805 -  ctx - context for residual evaluation
806 
807    Level: advanced
808 
809 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
810 @*/
811 PetscErrorCode DMSNESGetBlockJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
812 {
813   PetscErrorCode ierr;
814   SNESDM sdm;
815 
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
818   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
819   if (func) *func = sdm->computeblockjacobian;
820   if (ctx)  *ctx = sdm->blockjacobianctx;
821   PetscFunctionReturn(0);
822 }
823