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