xref: /petsc/src/snes/interface/ftn-custom/zsnesf.c (revision 4cc2b5b5fe4ffd09e5956b56d7cdc4f43e324103)
1 #include <petsc/private/fortranimpl.h>
2 #include <petscsnes.h>
3 #include <petscviewer.h>
4 #include <petsc/private/f90impl.h>
5 
6 #if defined(PETSC_HAVE_FORTRAN_CAPS)
7   #define snessetpicard_                   SNESSETPICARD
8   #define snessolve_                       SNESSOLVE
9   #define snescomputejacobiandefault_      SNESCOMPUTEJACOBIANDEFAULT
10   #define snescomputejacobiandefaultcolor_ SNESCOMPUTEJACOBIANDEFAULTCOLOR
11   #define snessetjacobian_                 SNESSETJACOBIAN
12   #define snessetjacobian1_                SNESSETJACOBIAN1
13   #define snessetjacobian2_                SNESSETJACOBIAN2
14   #define snessetfunction_                 SNESSETFUNCTION
15   #define snessetobjective_                SNESSETOBJECTIVE
16   #define snessetngs_                      SNESSETNGS
17   #define snessetupdate_                   SNESSETUPDATE
18   #define snesgetfunction_                 SNESGETFUNCTION
19   #define snesgetngs_                      SNESGETNGS
20   #define snessetconvergencetest_          SNESSETCONVERGENCETEST
21   #define snesconvergeddefault_            SNESCONVERGEDDEFAULT
22   #define snesconvergedskip_               SNESCONVERGEDSKIP
23   #define snesgetconvergencehistory_       SNESGETCONVERGENCEHISTORY
24   #define snesgetjacobian_                 SNESGETJACOBIAN
25   #define snesmonitordefault_              SNESMONITORDEFAULT
26   #define snesmonitorsolution_             SNESMONITORSOLUTION
27   #define snesmonitorsolutionupdate_       SNESMONITORSOLUTIONUPDATE
28   #define snesmonitorset_                  SNESMONITORSET
29   #define snesnewtontrsetprecheck_         SNESNEWTONTRSETPRECHECK
30   #define snesnewtontrsetpostcheck_        SNESNEWTONTRSETPOSTCHECK
31   #define snesnewtontrdcsetprecheck_       SNESNEWTONTRDCSETPRECHECK
32   #define snesnewtontrdcsetpostcheck_      SNESNEWTONTRDCSETPOSTCHECK
33   #define matmffdcomputejacobian_          MATMFFDCOMPUTEJACOBIAN
34 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
35   #define snessetpicard_                   snessetpicard
36   #define snessolve_                       snessolve
37   #define snescomputejacobiandefault_      snescomputejacobiandefault
38   #define snescomputejacobiandefaultcolor_ snescomputejacobiandefaultcolor
39   #define snessetjacobian_                 snessetjacobian
40   #define snessetjacobian1_                snessetjacobian1
41   #define snessetjacobian2_                snessetjacobian2
42   #define snessetfunction_                 snessetfunction
43   #define snessetobjective_                snessetobjective
44   #define snessetngs_                      snessetngs
45   #define snessetupdate_                   snessetupdate
46   #define snesgetfunction_                 snesgetfunction
47   #define snesgetngs_                      snesgetngs
48   #define snessetconvergencetest_          snessetconvergencetest
49   #define snesconvergeddefault_            snesconvergeddefault
50   #define snesconvergedskip_               snesconvergedskip
51   #define snesgetjacobian_                 snesgetjacobian
52   #define snesgetconvergencehistory_       snesgetconvergencehistory
53   #define snesmonitordefault_              snesmonitordefault
54   #define snesmonitorsolution_             snesmonitorsolution
55   #define snesmonitorsolutionupdate_       snesmonitorsolutionupdate
56   #define snesmonitorset_                  snesmonitorset
57   #define snesnewtontrsetprecheck_         snesnewtontrsetprecheck
58   #define snesnewtontrsetpostcheck_        snesnewtontrsetpostcheck
59   #define snesnewtontrdcsetprecheck_       snesnewtontrdcsetprecheck
60   #define snesnewtontrdcsetpostcheck_      snesnewtontrdcsetpostcheck
61   #define matmffdcomputejacobian_          matmffdcomputejacobian
62 #endif
63 
64 static struct {
65   PetscFortranCallbackId function;
66   PetscFortranCallbackId objective;
67   PetscFortranCallbackId test;
68   PetscFortranCallbackId destroy;
69   PetscFortranCallbackId jacobian;
70   PetscFortranCallbackId monitor;
71   PetscFortranCallbackId mondestroy;
72   PetscFortranCallbackId ngs;
73   PetscFortranCallbackId update;
74   PetscFortranCallbackId trprecheck;
75   PetscFortranCallbackId trpostcheck;
76 #if defined(PETSC_HAVE_F90_2PTR_ARG)
77   PetscFortranCallbackId function_pgiptr;
78   PetscFortranCallbackId objective_pgiptr;
79   PetscFortranCallbackId trprecheck_pgiptr;
80   PetscFortranCallbackId trpostcheck_pgiptr;
81 #endif
82 } _cb;
83 
84 static PetscErrorCode ourtrprecheckfunction(SNES snes, Vec x, Vec y, PetscBool *changed_y, void *ctx)
85 {
86 #if defined(PETSC_HAVE_F90_2PTR_ARG)
87   void *ptr;
88   PetscCall(PetscObjectGetFortranCallback((PetscObject)snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.trprecheck_pgiptr, NULL, &ptr));
89 #endif
90   PetscObjectUseFortranCallback(snes, _cb.trprecheck, (SNES *, Vec *, Vec *, PetscBool *, void *, PetscErrorCode *PETSC_F90_2PTR_PROTO_NOVAR), (&snes, &x, &y, changed_y, _ctx, &ierr PETSC_F90_2PTR_PARAM(ptr)));
91 }
92 
93 PETSC_EXTERN void snesnewtontrsetprecheck_(SNES *snes, void (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
94 {
95   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trprecheck, (PetscVoidFn *)func, ctx);
96   if (*ierr) return;
97 #if defined(PETSC_HAVE_F90_2PTR_ARG)
98   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trprecheck_pgiptr, NULL, ptr);
99   if (*ierr) return;
100 #endif
101   *ierr = SNESNewtonTRSetPreCheck(*snes, ourtrprecheckfunction, NULL);
102 }
103 
104 PETSC_EXTERN void snesnewtontrdcsetprecheck_(SNES *snes, void (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
105 {
106   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trprecheck, (PetscVoidFn *)func, ctx);
107   if (*ierr) return;
108 #if defined(PETSC_HAVE_F90_2PTR_ARG)
109   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trprecheck_pgiptr, NULL, ptr);
110   if (*ierr) return;
111 #endif
112   *ierr = SNESNewtonTRDCSetPreCheck(*snes, ourtrprecheckfunction, NULL);
113 }
114 
115 static PetscErrorCode ourtrpostcheckfunction(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx)
116 {
117 #if defined(PETSC_HAVE_F90_2PTR_ARG)
118   void *ptr;
119   PetscCall(PetscObjectGetFortranCallback((PetscObject)snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.trpostcheck_pgiptr, NULL, &ptr));
120 #endif
121   PetscObjectUseFortranCallback(snes, _cb.trpostcheck, (SNES *, Vec *, Vec *, Vec *, PetscBool *, PetscBool *, void *, PetscErrorCode *PETSC_F90_2PTR_PROTO_NOVAR), (&snes, &x, &y, &w, changed_y, changed_w, _ctx, &ierr PETSC_F90_2PTR_PARAM(ptr)));
122 }
123 
124 PETSC_EXTERN void snesnewtontrsetpostcheck_(SNES *snes, void (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
125 {
126   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trpostcheck, (PetscVoidFn *)func, ctx);
127   if (*ierr) return;
128 #if defined(PETSC_HAVE_F90_2PTR_ARG)
129   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trpostcheck_pgiptr, NULL, ptr);
130   if (*ierr) return;
131 #endif
132   *ierr = SNESNewtonTRSetPostCheck(*snes, ourtrpostcheckfunction, NULL);
133 }
134 
135 PETSC_EXTERN void snesnewtontrdcsetpostcheck_(SNES *snes, void (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
136 {
137   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trpostcheck, (PetscVoidFn *)func, ctx);
138   if (*ierr) return;
139 #if defined(PETSC_HAVE_F90_2PTR_ARG)
140   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.trpostcheck_pgiptr, NULL, ptr);
141   if (*ierr) return;
142 #endif
143   *ierr = SNESNewtonTRDCSetPostCheck(*snes, ourtrpostcheckfunction, NULL);
144 }
145 
146 static PetscErrorCode oursnesfunction(SNES snes, Vec x, Vec f, void *ctx)
147 {
148 #if defined(PETSC_HAVE_F90_2PTR_ARG)
149   void *ptr;
150   PetscCall(PetscObjectGetFortranCallback((PetscObject)snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
151 #endif
152   PetscObjectUseFortranCallback(snes, _cb.function, (SNES *, Vec *, Vec *, void *, PetscErrorCode *PETSC_F90_2PTR_PROTO_NOVAR), (&snes, &x, &f, _ctx, &ierr PETSC_F90_2PTR_PARAM(ptr)));
153 }
154 
155 static PetscErrorCode oursnesobjective(SNES snes, Vec x, PetscReal *v, void *ctx)
156 {
157 #if defined(PETSC_HAVE_F90_2PTR_ARG)
158   void *ptr;
159   PetscCall(PetscObjectGetFortranCallback((PetscObject)snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.objective_pgiptr, NULL, &ptr));
160 #endif
161   PetscObjectUseFortranCallback(snes, _cb.objective, (SNES *, Vec *, PetscReal *, void *, PetscErrorCode *PETSC_F90_2PTR_PROTO_NOVAR), (&snes, &x, v, _ctx, &ierr PETSC_F90_2PTR_PARAM(ptr)));
162 }
163 
164 static PetscErrorCode oursnestest(SNES snes, PetscInt it, PetscReal a, PetscReal d, PetscReal c, SNESConvergedReason *reason, void *ctx)
165 {
166   PetscObjectUseFortranCallback(snes, _cb.test, (SNES *, PetscInt *, PetscReal *, PetscReal *, PetscReal *, SNESConvergedReason *, void *, PetscErrorCode *), (&snes, &it, &a, &d, &c, reason, _ctx, &ierr));
167 }
168 
169 static PetscErrorCode ourdestroy(void *ctx)
170 {
171   PetscObjectUseFortranCallback(ctx, _cb.destroy, (void *, PetscErrorCode *), (_ctx, &ierr));
172 }
173 
174 static PetscErrorCode oursnesjacobian(SNES snes, Vec x, Mat m, Mat p, void *ctx)
175 {
176   PetscObjectUseFortranCallback(snes, _cb.jacobian, (SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&snes, &x, &m, &p, _ctx, &ierr));
177 }
178 
179 static PetscErrorCode oursnesupdate(SNES snes, PetscInt i)
180 {
181   PetscObjectUseFortranCallback(snes, _cb.update, (SNES *, PetscInt *, PetscErrorCode *), (&snes, &i, &ierr));
182 }
183 static PetscErrorCode oursnesngs(SNES snes, Vec x, Vec b, void *ctx)
184 {
185   PetscObjectUseFortranCallback(snes, _cb.ngs, (SNES *, Vec *, Vec *, void *, PetscErrorCode *), (&snes, &x, &b, _ctx, &ierr));
186 }
187 static PetscErrorCode oursnesmonitor(SNES snes, PetscInt i, PetscReal d, void *ctx)
188 {
189   PetscObjectUseFortranCallback(snes, _cb.monitor, (SNES *, PetscInt *, PetscReal *, void *, PetscErrorCode *), (&snes, &i, &d, _ctx, &ierr));
190 }
191 static PetscErrorCode ourmondestroy(void **ctx)
192 {
193   SNES snes = (SNES)*ctx;
194   PetscObjectUseFortranCallback(snes, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
195 }
196 
197 /* these are generated automatically by bfort */
198 PETSC_EXTERN void snescomputejacobiandefault_(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *);
199 PETSC_EXTERN void snescomputejacobiandefaultcolor_(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *);
200 PETSC_EXTERN void matmffdcomputejacobian_(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *);
201 
202 PETSC_EXTERN void snessetjacobian_(SNES *snes, Mat *A, Mat *B, void (*func)(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
203 {
204   CHKFORTRANNULLFUNCTION(func);
205   if ((PetscVoidFn *)func == (PetscVoidFn *)snescomputejacobiandefault_) {
206     *ierr = SNESSetJacobian(*snes, *A, *B, SNESComputeJacobianDefault, ctx);
207   } else if ((PetscVoidFn *)func == (PetscVoidFn *)snescomputejacobiandefaultcolor_) {
208     if (!ctx) {
209       *ierr = PETSC_ERR_ARG_NULL;
210       return;
211     }
212     *ierr = SNESSetJacobian(*snes, *A, *B, SNESComputeJacobianDefaultColor, *(MatFDColoring *)ctx);
213   } else if ((PetscVoidFn *)func == (PetscVoidFn *)matmffdcomputejacobian_) {
214     *ierr = SNESSetJacobian(*snes, *A, *B, MatMFFDComputeJacobian, ctx);
215   } else {
216     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacobian, (PetscVoidFn *)func, ctx);
217     if (!*ierr) *ierr = SNESSetJacobian(*snes, *A, *B, oursnesjacobian, NULL);
218   }
219 }
220 PETSC_EXTERN void snessetjacobian1_(SNES *snes, Mat *A, Mat *B, void (*func)(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
221 {
222   snessetjacobian_(snes, A, B, func, ctx, ierr);
223 }
224 PETSC_EXTERN void snessetjacobian2_(SNES *snes, Mat *A, Mat *B, void (*func)(SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
225 {
226   snessetjacobian_(snes, A, B, func, ctx, ierr);
227 }
228 
229 static PetscErrorCode oursnespicardfunction(SNES snes, Vec x, Vec f, void *ctx)
230 {
231 #if defined(PETSC_HAVE_F90_2PTR_ARG)
232   void *ptr;
233   PetscCall(PetscObjectGetFortranCallback((PetscObject)snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
234 #endif
235   PetscObjectUseFortranCallback(snes, _cb.function, (SNES *, Vec *, Vec *, void *, PetscErrorCode *PETSC_F90_2PTR_PROTO_NOVAR), (&snes, &x, &f, _ctx, &ierr PETSC_F90_2PTR_PARAM(ptr)));
236 }
237 
238 static PetscErrorCode oursnespicardjacobian(SNES snes, Vec x, Mat m, Mat p, void *ctx)
239 {
240   PetscObjectUseFortranCallback(snes, _cb.jacobian, (SNES *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&snes, &x, &m, &p, _ctx, &ierr));
241 }
242 
243 PETSC_EXTERN void snessetpicard_(SNES *snes, Vec *r, void (*func)(SNES *, Vec *, Vec *, void *, PetscErrorCode *), Mat *A, Mat *B, PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
244 {
245   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.function, (PetscVoidFn *)func, ctx);
246 #if defined(PETSC_HAVE_F90_2PTR_ARG)
247   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.function_pgiptr, NULL, ptr);
248   if (*ierr) return;
249 #endif
250   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacobian, (PetscVoidFn *)J, ctx);
251   if (!*ierr) *ierr = SNESSetPicard(*snes, *r, oursnespicardfunction, *A, *B, oursnespicardjacobian, NULL);
252 }
253 
254 /*
255    These are not usually called from Fortran but allow Fortran users
256    to transparently set these monitors from .F code
257 */
258 
259 PETSC_EXTERN void snessetfunction_(SNES *snes, Vec *r, SNESFunctionFn func, void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
260 {
261   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.function, (PetscVoidFn *)func, ctx);
262   if (*ierr) return;
263 #if defined(PETSC_HAVE_F90_2PTR_ARG)
264   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.function_pgiptr, NULL, ptr);
265   if (*ierr) return;
266 #endif
267   *ierr = SNESSetFunction(*snes, *r, oursnesfunction, NULL);
268 }
269 
270 PETSC_EXTERN void snessetobjective_(SNES *snes, void (*func)(SNES *, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr))
271 {
272   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.objective, (PetscVoidFn *)func, ctx);
273   if (*ierr) return;
274 #if defined(PETSC_HAVE_F90_2PTR_ARG)
275   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.objective_pgiptr, NULL, ptr);
276   if (*ierr) return;
277 #endif
278   *ierr = SNESSetObjective(*snes, oursnesobjective, NULL);
279 }
280 
281 PETSC_EXTERN void snessetngs_(SNES *snes, void (*func)(SNES *, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
282 {
283   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ngs, (PetscVoidFn *)func, ctx);
284   if (*ierr) return;
285   *ierr = SNESSetNGS(*snes, oursnesngs, NULL);
286 }
287 PETSC_EXTERN void snessetupdate_(SNES *snes, void (*func)(SNES *, PetscInt *, PetscErrorCode *), PetscErrorCode *ierr)
288 {
289   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.update, (PetscVoidFn *)func, NULL);
290   if (*ierr) return;
291   *ierr = SNESSetUpdate(*snes, oursnesupdate);
292 }
293 
294 /* the func argument is ignored */
295 PETSC_EXTERN void snesgetfunction_(SNES *snes, Vec *r, SNESFunctionFn func, void **ctx, PetscErrorCode *ierr)
296 {
297   CHKFORTRANNULLOBJECT(r);
298   *ierr = SNESGetFunction(*snes, r, NULL, NULL);
299   if (*ierr) return;
300   if ((PetscVoidFn *)func == (PetscVoidFn *)PETSC_NULL_FUNCTION_Fortran) return;
301   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function, NULL, ctx);
302 }
303 
304 PETSC_EXTERN void snesgetngs_(SNES *snes, void *func, void **ctx, PetscErrorCode *ierr)
305 {
306   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.ngs, NULL, ctx);
307 }
308 
309 PETSC_EXTERN void snesconvergeddefault_(SNES *snes, PetscInt *it, PetscReal *a, PetscReal *b, PetscReal *c, SNESConvergedReason *r, void *ct, PetscErrorCode *ierr)
310 {
311   *ierr = SNESConvergedDefault(*snes, *it, *a, *b, *c, r, ct);
312 }
313 
314 PETSC_EXTERN void snesconvergedskip_(SNES *snes, PetscInt *it, PetscReal *a, PetscReal *b, PetscReal *c, SNESConvergedReason *r, void *ct, PetscErrorCode *ierr)
315 {
316   *ierr = SNESConvergedSkip(*snes, *it, *a, *b, *c, r, ct);
317 }
318 
319 PETSC_EXTERN void snessetconvergencetest_(SNES *snes, void (*func)(SNES *, PetscInt *, PetscReal *, PetscReal *, PetscReal *, SNESConvergedReason *, void *, PetscErrorCode *), void *cctx, void (*destroy)(void *), PetscErrorCode *ierr)
320 {
321   CHKFORTRANNULLFUNCTION(destroy);
322 
323   if ((PetscVoidFn *)func == (PetscVoidFn *)snesconvergeddefault_) {
324     *ierr = SNESSetConvergenceTest(*snes, SNESConvergedDefault, NULL, NULL);
325   } else if ((PetscVoidFn *)func == (PetscVoidFn *)snesconvergedskip_) {
326     *ierr = SNESSetConvergenceTest(*snes, SNESConvergedSkip, NULL, NULL);
327   } else {
328     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.test, (PetscVoidFn *)func, cctx);
329     if (*ierr) return;
330     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.destroy, (PetscVoidFn *)destroy, cctx);
331     if (*ierr) return;
332     *ierr = SNESSetConvergenceTest(*snes, oursnestest, *snes, ourdestroy);
333   }
334 }
335 
336 /*  func is currently ignored from Fortran */
337 PETSC_EXTERN void snesgetjacobian_(SNES *snes, Mat *A, Mat *B, int *func, void **ctx, PetscErrorCode *ierr)
338 {
339   CHKFORTRANNULLINTEGER(ctx);
340   CHKFORTRANNULLOBJECT(A);
341   CHKFORTRANNULLOBJECT(B);
342   *ierr = SNESGetJacobian(*snes, A, B, NULL, NULL);
343   if (*ierr) return;
344   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, _cb.jacobian, NULL, ctx);
345 }
346 
347 PETSC_EXTERN void snesgetconvergencehistory_(SNES *snes, PetscInt *na, PetscErrorCode *ierr)
348 {
349   *ierr = SNESGetConvergenceHistory(*snes, NULL, NULL, na);
350 }
351 
352 PETSC_EXTERN void snesmonitordefault_(SNES *snes, PetscInt *its, PetscReal *fgnorm, PetscViewerAndFormat **dummy, PetscErrorCode *ierr)
353 {
354   *ierr = SNESMonitorDefault(*snes, *its, *fgnorm, *dummy);
355 }
356 
357 PETSC_EXTERN void snesmonitorsolution_(SNES *snes, PetscInt *its, PetscReal *fgnorm, PetscViewerAndFormat **dummy, PetscErrorCode *ierr)
358 {
359   *ierr = SNESMonitorSolution(*snes, *its, *fgnorm, *dummy);
360 }
361 
362 PETSC_EXTERN void snesmonitorsolutionupdate_(SNES *snes, PetscInt *its, PetscReal *fgnorm, PetscViewerAndFormat **dummy, PetscErrorCode *ierr)
363 {
364   *ierr = SNESMonitorSolutionUpdate(*snes, *its, *fgnorm, *dummy);
365 }
366 
367 PETSC_EXTERN void snesmonitorset_(SNES *snes, void (*func)(SNES *, PetscInt *, PetscReal *, void *, PetscErrorCode *), void *mctx, void (*mondestroy)(void *, PetscErrorCode *), PetscErrorCode *ierr)
368 {
369   CHKFORTRANNULLFUNCTION(mondestroy);
370   if ((PetscVoidFn *)func == (PetscVoidFn *)snesmonitordefault_) {
371     *ierr = SNESMonitorSet(*snes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy);
372   } else if ((PetscVoidFn *)func == (PetscVoidFn *)snesmonitorsolution_) {
373     *ierr = SNESMonitorSet(*snes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorSolution, *(PetscViewerAndFormat **)mctx, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy);
374   } else if ((PetscVoidFn *)func == (PetscVoidFn *)snesmonitorsolutionupdate_) {
375     *ierr = SNESMonitorSet(*snes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorSolutionUpdate, *(PetscViewerAndFormat **)mctx, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy);
376   } else {
377     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscVoidFn *)func, mctx);
378     if (*ierr) return;
379     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscVoidFn *)mondestroy, mctx);
380     if (*ierr) return;
381     *ierr = SNESMonitorSet(*snes, oursnesmonitor, *snes, ourmondestroy);
382   }
383 }
384