xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
1 #include <petsc/private/ftnimpl.h>
2 #include <petscts.h>
3 #include <petscviewer.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6   #define tsmonitorlgsettransform_      TSMONITORLGSETTRANSFORM
7   #define tssetrhsfunction_             TSSETRHSFUNCTION
8   #define tsgetrhsfunction_             TSGETRHSFUNCTION
9   #define tssetrhsjacobian_             TSSETRHSJACOBIAN
10   #define tsgetrhsjacobian_             TSGETRHSJACOBIAN
11   #define tssetifunction_               TSSETIFUNCTION
12   #define tsgetifunction_               TSGETIFUNCTION
13   #define tssetijacobian_               TSSETIJACOBIAN
14   #define tsgetijacobian_               TSGETIJACOBIAN
15   #define tsmonitorset_                 TSMONITORSET
16   #define tssetrhsjacobianp_            TSSETRHSJACOBIANP
17   #define tsgetrhsjacobianp_            TSGETRHSJACOBIANP
18   #define tssetijacobianp_              TSSETIJACOBIANP
19   #define tsgetijacobianp_              TSGETIJACOBIANP
20   #define tscomputerhsfunctionlinear_   TSCOMPUTERHSFUNCTIONLINEAR
21   #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT
22   #define tscomputeifunctionlinear_     TSCOMPUTEIFUNCTIONLINEAR
23   #define tscomputeijacobianconstant_   TSCOMPUTEIJACOBIANCONSTANT
24   #define tsmonitordefault_             TSMONITORDEFAULT
25   #define tssetprestep_                 TSSETPRESTEP
26   #define tssetpoststep_                TSSETPOSTSTEP
27 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
28   #define tsmonitorlgsettransform_      tsmonitorlgsettransform
29   #define tssetrhsfunction_             tssetrhsfunction
30   #define tsgetrhsfunction_             tsgetrhsfunction
31   #define tssetrhsjacobian_             tssetrhsjacobian
32   #define tsgetrhsjacobian_             tsgetrhsjacobian
33   #define tssetifunction_               tssetifunction
34   #define tsgetifunction_               tsgetifunction
35   #define tssetijacobian_               tssetijacobian
36   #define tsgetijacobian_               tsgetijacobian
37   #define tssetijacobianp_              tssetijacobianp
38   #define tsgetijacobianp_              tsgetijacobianp
39   #define tssetrhsjacobianp_            tssetrhsjacobianp
40   #define tsgetrhsjacobianp_            tsgetrhsjacobianp
41   #define tsmonitorset_                 tsmonitorset
42   #define tscomputerhsfunctionlinear_   tscomputerhsfunctionlinear
43   #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant
44   #define tscomputeifunctionlinear_     tscomputeifunctionlinear
45   #define tscomputeijacobianconstant_   tscomputeijacobianconstant
46   #define tsmonitordefault_             tsmonitordefault
47   #define tssetprestep_                 tssetprestep
48   #define tssetpoststep_                tssetpoststep
49 #endif
50 
51 static struct {
52   PetscFortranCallbackId prestep;
53   PetscFortranCallbackId poststep;
54   PetscFortranCallbackId rhsfunction;
55   PetscFortranCallbackId rhsjacobian;
56   PetscFortranCallbackId ifunction;
57   PetscFortranCallbackId ijacobian;
58   PetscFortranCallbackId rhsjacobianp;
59   PetscFortranCallbackId ijacobianp;
60   PetscFortranCallbackId monitor;
61   PetscFortranCallbackId mondestroy;
62   PetscFortranCallbackId transform;
63 #if defined(PETSC_HAVE_F90_2PTR_ARG)
64   PetscFortranCallbackId function_pgiptr;
65 #endif
66 } _cb;
67 
ourprestep(TS ts)68 static PetscErrorCode ourprestep(TS ts)
69 {
70 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
71   void *ptr;
72   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
73 #endif
74   PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
75 }
ourpoststep(TS ts)76 static PetscErrorCode ourpoststep(TS ts)
77 {
78 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
79   void *ptr;
80   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
81 #endif
82   PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
83 }
ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,PetscCtx ctx)84 static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, PetscCtx ctx)
85 {
86 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
87   void *ptr;
88   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
89 #endif
90   PetscObjectUseFortranCallback(ts, _cb.rhsfunction, (TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &f, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
91 }
ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,PetscCtx ctx)92 static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, PetscCtx ctx)
93 {
94 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
95   void *ptr;
96   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
97 #endif
98   PetscObjectUseFortranCallback(ts, _cb.ifunction, (TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &xdot, &f, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
99 }
ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat m,Mat p,PetscCtx ctx)100 static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, PetscCtx ctx)
101 {
102 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
103   void *ptr;
104   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
105 #endif
106   PetscObjectUseFortranCallback(ts, _cb.rhsjacobian, (TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &m, &p, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
107 }
ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,Mat p,PetscCtx ctx)108 static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, PetscCtx ctx)
109 {
110 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
111   void *ptr;
112   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
113 #endif
114   PetscObjectUseFortranCallback(ts, _cb.ijacobian, (TS *, PetscReal *, Vec *, Vec *, PetscReal *, Mat *, Mat *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &xdot, &shift, &m, &p, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
115 }
ourijacobianp(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,PetscCtx ctx)116 static PetscErrorCode ourijacobianp(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, PetscCtx ctx)
117 {
118 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
119   void *ptr;
120   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
121 #endif
122   PetscObjectUseFortranCallback(ts, _cb.ijacobianp, (TS *, PetscReal *, Vec *, Vec *, PetscReal *, Mat *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &xdot, &shift, &m, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
123 }
ourrhsjacobianp(TS ts,PetscReal d,Vec x,Mat m,PetscCtx ctx)124 static PetscErrorCode ourrhsjacobianp(TS ts, PetscReal d, Vec x, Mat m, PetscCtx ctx)
125 {
126 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
127   void *ptr;
128   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
129 #endif
130   PetscObjectUseFortranCallback(ts, _cb.rhsjacobianp, (TS *, PetscReal *, Vec *, Mat *, void *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &d, &x, &m, _ctx, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
131 }
132 
ourmonitordestroy(PetscCtxRt ctx)133 static PetscErrorCode ourmonitordestroy(PetscCtxRt ctx)
134 {
135   TS ts = *(TS *)ctx;
136   PetscObjectUseFortranCallback(ts, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
137 }
138 
139 /*
140    Note ctx is the same as ts so we need to get the Fortran context out of the TS
141 */
ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,PetscCtx ctx)142 static PetscErrorCode ourmonitor(TS ts, PetscInt i, PetscReal d, Vec v, PetscCtx ctx)
143 {
144   PetscObjectUseFortranCallback(ts, _cb.monitor, (TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), (&ts, &i, &d, &v, _ctx, &ierr));
145 }
146 
147 /*
148    Currently does not handle destroy or context
149 */
ourtransform(PetscCtx ctx,Vec x,Vec * xout)150 static PetscErrorCode ourtransform(PetscCtx ctx, Vec x, Vec *xout)
151 {
152   PetscObjectUseFortranCallback((TS)ctx, _cb.transform, (void *, Vec *, Vec *, PetscErrorCode *), (_ctx, &x, xout, &ierr));
153 }
154 
tsmonitorlgsettransform_(TS * ts,void (* f)(void *,Vec *,Vec *,PetscErrorCode *),PetscErrorCode (* d)(void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)155 PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts, void (*f)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode (*d)(void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
156 {
157   *ierr = TSMonitorLGSetTransform(*ts, ourtransform, NULL, NULL);
158   if (*ierr) return;
159   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.transform, (PetscFortranCallbackFn *)f, ctx);
160 }
161 
tssetprestep_(TS * ts,PetscErrorCode (* f)(TS *,PetscErrorCode *),PetscErrorCode * ierr)162 PETSC_EXTERN void tssetprestep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr)
163 {
164   *ierr = TSSetPreStep(*ts, ourprestep);
165   if (*ierr) return;
166   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.prestep, (PetscFortranCallbackFn *)f, NULL);
167 }
168 
tssetpoststep_(TS * ts,PetscErrorCode (* f)(TS *,PetscErrorCode *),PetscErrorCode * ierr)169 PETSC_EXTERN void tssetpoststep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr)
170 {
171   *ierr = TSSetPostStep(*ts, ourpoststep);
172   if (*ierr) return;
173   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.poststep, (PetscFortranCallbackFn *)f, NULL);
174 }
175 
176 PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *);
177 
tssetrhsfunction_(TS * ts,Vec * r,void (* f)(TS *,PetscReal *,Vec *,Vec *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)178 PETSC_EXTERN void tssetrhsfunction_(TS *ts, Vec *r, void (*f)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
179 {
180   Vec R;
181   CHKFORTRANNULLOBJECT(r);
182   CHKFORTRANNULLFUNCTION(f);
183   R = r ? *r : (Vec)NULL;
184   if (f == tscomputerhsfunctionlinear_) {
185     *ierr = TSSetRHSFunction(*ts, R, TSComputeRHSFunctionLinear, fP);
186   } else {
187     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsfunction, (PetscFortranCallbackFn *)f, fP);
188     *ierr = TSSetRHSFunction(*ts, R, ourrhsfunction, NULL);
189   }
190 }
tsgetrhsfunction_(TS * ts,Vec * r,void * func,void ** ctx,PetscErrorCode * ierr)191 PETSC_EXTERN void tsgetrhsfunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr)
192 {
193   CHKFORTRANNULLINTEGER(ctx);
194   CHKFORTRANNULLOBJECT(r);
195   *ierr = TSGetRHSFunction(*ts, r, NULL, ctx);
196 }
197 
198 PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, Vec *F, PetscCtx ctx, PetscErrorCode *ierr);
199 
tssetifunction_(TS * ts,Vec * r,void (* f)(TS *,PetscReal *,Vec *,Vec *,Vec *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)200 PETSC_EXTERN void tssetifunction_(TS *ts, Vec *r, void (*f)(TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
201 {
202   Vec R;
203   CHKFORTRANNULLOBJECT(r);
204   CHKFORTRANNULLFUNCTION(f);
205   R = r ? *r : (Vec)NULL;
206   if (f == tscomputeifunctionlinear_) {
207     *ierr = TSSetIFunction(*ts, R, TSComputeIFunctionLinear, fP);
208   } else {
209     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ifunction, (PetscFortranCallbackFn *)f, fP);
210     *ierr = TSSetIFunction(*ts, R, ourifunction, NULL);
211   }
212 }
tsgetifunction_(TS * ts,Vec * r,void * func,void ** ctx,PetscErrorCode * ierr)213 PETSC_EXTERN void tsgetifunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr)
214 {
215   CHKFORTRANNULLINTEGER(ctx);
216   CHKFORTRANNULLOBJECT(r);
217   *ierr = TSGetIFunction(*ts, r, NULL, ctx);
218 }
219 
220 PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *);
221 
tssetrhsjacobian_(TS * ts,Mat * A,Mat * B,void (* f)(TS *,PetscReal *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)222 PETSC_EXTERN void tssetrhsjacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
223 {
224   CHKFORTRANNULLFUNCTION(f);
225   if (f == tscomputerhsjacobianconstant_) {
226     *ierr = TSSetRHSJacobian(*ts, *A, *B, TSComputeRHSJacobianConstant, fP);
227   } else {
228     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobian, (PetscFortranCallbackFn *)f, fP);
229     *ierr = TSSetRHSJacobian(*ts, *A, *B, ourrhsjacobian, NULL);
230   }
231 }
232 
233 PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, PetscReal *shift, Mat *A, Mat *B, PetscCtx ctx, PetscErrorCode *ierr);
234 
tssetijacobian_(TS * ts,Mat * A,Mat * B,void (* f)(TS *,PetscReal *,Vec *,Vec *,PetscReal *,Mat *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)235 PETSC_EXTERN void tssetijacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Vec *, PetscReal *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
236 {
237   CHKFORTRANNULLFUNCTION(f);
238   if (f == tscomputeijacobianconstant_) {
239     *ierr = TSSetIJacobian(*ts, *A, *B, TSComputeIJacobianConstant, fP);
240   } else {
241     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobian, (PetscFortranCallbackFn *)f, fP);
242     *ierr = TSSetIJacobian(*ts, *A, *B, ourijacobian, NULL);
243   }
244 }
tsgetijacobian_(TS * ts,Mat * J,Mat * M,int * func,void ** ctx,PetscErrorCode * ierr)245 PETSC_EXTERN void tsgetijacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr)
246 {
247   CHKFORTRANNULLINTEGER(ctx);
248   CHKFORTRANNULLOBJECT(J);
249   CHKFORTRANNULLOBJECT(M);
250   *ierr = TSGetIJacobian(*ts, J, M, NULL, ctx);
251 }
tssetijacobianp_(TS * ts,Mat * A,void (* f)(TS *,PetscReal *,Vec *,Vec *,PetscReal,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)252 PETSC_EXTERN void tssetijacobianp_(TS *ts, Mat *A, void (*f)(TS *, PetscReal *, Vec *, Vec *, PetscReal, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
253 {
254   CHKFORTRANNULLFUNCTION(f);
255   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobianp, (PetscFortranCallbackFn *)f, fP);
256   *ierr = TSSetIJacobianP(*ts, *A, ourijacobianp, NULL);
257 }
tsgetijacobianp_(TS * ts,Mat * J,int * func,void ** ctx,PetscErrorCode * ierr)258 PETSC_EXTERN void tsgetijacobianp_(TS *ts, Mat *J, int *func, void **ctx, PetscErrorCode *ierr)
259 {
260   CHKFORTRANNULLINTEGER(ctx);
261   CHKFORTRANNULLOBJECT(J);
262   *ierr = TSGetIJacobianP(*ts, J, NULL, ctx);
263 }
tssetrhsjacobianp_(TS * ts,Mat * A,void (* f)(TS *,PetscReal *,Vec *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)264 PETSC_EXTERN void tssetrhsjacobianp_(TS *ts, Mat *A, void (*f)(TS *, PetscReal *, Vec *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
265 {
266   CHKFORTRANNULLFUNCTION(f);
267   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobianp, (PetscFortranCallbackFn *)f, fP);
268   *ierr = TSSetRHSJacobianP(*ts, *A, ourrhsjacobianp, NULL);
269 }
tsgetrhsjacobianp_(TS * ts,Mat * J,int * func,void ** ctx,PetscErrorCode * ierr)270 PETSC_EXTERN void tsgetrhsjacobianp_(TS *ts, Mat *J, int *func, void **ctx, PetscErrorCode *ierr)
271 {
272   CHKFORTRANNULLINTEGER(ctx);
273   CHKFORTRANNULLOBJECT(J);
274   *ierr = TSGetRHSJacobianP(*ts, J, NULL, ctx);
275 }
276 
277 PETSC_EXTERN void tsmonitordefault_(TS *, PetscInt *, PetscReal *, Vec *, PetscViewerAndFormat **, PetscErrorCode *);
278 
279 /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */
280 
tsmonitorset_(TS * ts,void (* func)(TS *,PetscInt *,PetscReal *,Vec *,void *,PetscErrorCode *),void * mctx,void (* d)(void *,PetscErrorCode *),PetscErrorCode * ierr)281 PETSC_EXTERN void tsmonitorset_(TS *ts, void (*func)(TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), void *mctx, void (*d)(void *, PetscErrorCode *), PetscErrorCode *ierr)
282 {
283   CHKFORTRANNULLFUNCTION(d);
284   if ((PetscFortranCallbackFn *)func == (PetscFortranCallbackFn *)tsmonitordefault_) {
285     *ierr = TSMonitorSet(*ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
286   } else {
287     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscFortranCallbackFn *)func, mctx);
288     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscFortranCallbackFn *)d, mctx);
289     *ierr = TSMonitorSet(*ts, ourmonitor, *ts, ourmonitordestroy);
290   }
291 }
292 
293 /*  func is currently ignored from Fortran */
tsgetrhsjacobian_(TS * ts,Mat * J,Mat * M,int * func,void ** ctx,PetscErrorCode * ierr)294 PETSC_EXTERN void tsgetrhsjacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr)
295 {
296   *ierr = TSGetRHSJacobian(*ts, J, M, NULL, ctx);
297 }
298