xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2c6db04a5SJed Brown #include <petscts.h>
3665c2dedSJed Brown #include <petscviewer.h>
4e2df7a95SSatish Balay 
5e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
6b0bc92c6SBarry Smith   #define tsmonitorlgsettransform_      TSMONITORLGSETTRANSFORM
7e2df7a95SSatish Balay   #define tssetrhsfunction_             TSSETRHSFUNCTION
837698f3aSJed Brown   #define tsgetrhsfunction_             TSGETRHSFUNCTION
9e2df7a95SSatish Balay   #define tssetrhsjacobian_             TSSETRHSJACOBIAN
10e2df7a95SSatish Balay   #define tsgetrhsjacobian_             TSGETRHSJACOBIAN
1137698f3aSJed Brown   #define tssetifunction_               TSSETIFUNCTION
1237698f3aSJed Brown   #define tsgetifunction_               TSGETIFUNCTION
1337698f3aSJed Brown   #define tssetijacobian_               TSSETIJACOBIAN
1437698f3aSJed Brown   #define tsgetijacobian_               TSGETIJACOBIAN
15a6570f20SBarry Smith   #define tsmonitorset_                 TSMONITORSET
162549da82SBarry Smith   #define tssetrhsjacobianp_            TSSETRHSJACOBIANP
172549da82SBarry Smith   #define tsgetrhsjacobianp_            TSGETRHSJACOBIANP
182549da82SBarry Smith   #define tssetijacobianp_              TSSETIJACOBIANP
192549da82SBarry Smith   #define tsgetijacobianp_              TSGETIJACOBIANP
200e4ef248SJed Brown   #define tscomputerhsfunctionlinear_   TSCOMPUTERHSFUNCTIONLINEAR
210e4ef248SJed Brown   #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT
220fecffdcSJed Brown   #define tscomputeifunctionlinear_     TSCOMPUTEIFUNCTIONLINEAR
230fecffdcSJed Brown   #define tscomputeijacobianconstant_   TSCOMPUTEIJACOBIANCONSTANT
24a6570f20SBarry Smith   #define tsmonitordefault_             TSMONITORDEFAULT
25dd7ecb2fSBarry Smith   #define tssetprestep_                 TSSETPRESTEP
26dd7ecb2fSBarry Smith   #define tssetpoststep_                TSSETPOSTSTEP
27e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
28b0bc92c6SBarry Smith   #define tsmonitorlgsettransform_      tsmonitorlgsettransform
29e2df7a95SSatish Balay   #define tssetrhsfunction_             tssetrhsfunction
3037698f3aSJed Brown   #define tsgetrhsfunction_             tsgetrhsfunction
31e2df7a95SSatish Balay   #define tssetrhsjacobian_             tssetrhsjacobian
32e2df7a95SSatish Balay   #define tsgetrhsjacobian_             tsgetrhsjacobian
3337698f3aSJed Brown   #define tssetifunction_               tssetifunction
3437698f3aSJed Brown   #define tsgetifunction_               tsgetifunction
3537698f3aSJed Brown   #define tssetijacobian_               tssetijacobian
3637698f3aSJed Brown   #define tsgetijacobian_               tsgetijacobian
372549da82SBarry Smith   #define tssetijacobianp_              tssetijacobianp
382549da82SBarry Smith   #define tsgetijacobianp_              tsgetijacobianp
392549da82SBarry Smith   #define tssetrhsjacobianp_            tssetrhsjacobianp
402549da82SBarry Smith   #define tsgetrhsjacobianp_            tsgetrhsjacobianp
41a6570f20SBarry Smith   #define tsmonitorset_                 tsmonitorset
420e4ef248SJed Brown   #define tscomputerhsfunctionlinear_   tscomputerhsfunctionlinear
430e4ef248SJed Brown   #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant
440fecffdcSJed Brown   #define tscomputeifunctionlinear_     tscomputeifunctionlinear
450fecffdcSJed Brown   #define tscomputeijacobianconstant_   tscomputeijacobianconstant
46a6570f20SBarry Smith   #define tsmonitordefault_             tsmonitordefault
47dd7ecb2fSBarry Smith   #define tssetprestep_                 tssetprestep
48dd7ecb2fSBarry Smith   #define tssetpoststep_                tssetpoststep
49e2df7a95SSatish Balay #endif
50e2df7a95SSatish Balay 
51109c90ceSBarry Smith static struct {
52109c90ceSBarry Smith   PetscFortranCallbackId prestep;
53109c90ceSBarry Smith   PetscFortranCallbackId poststep;
54109c90ceSBarry Smith   PetscFortranCallbackId rhsfunction;
55109c90ceSBarry Smith   PetscFortranCallbackId rhsjacobian;
56109c90ceSBarry Smith   PetscFortranCallbackId ifunction;
57109c90ceSBarry Smith   PetscFortranCallbackId ijacobian;
582549da82SBarry Smith   PetscFortranCallbackId rhsjacobianp;
592549da82SBarry Smith   PetscFortranCallbackId ijacobianp;
60109c90ceSBarry Smith   PetscFortranCallbackId monitor;
61109c90ceSBarry Smith   PetscFortranCallbackId mondestroy;
62109c90ceSBarry Smith   PetscFortranCallbackId transform;
63109c90ceSBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG)
64109c90ceSBarry Smith   PetscFortranCallbackId function_pgiptr;
65109c90ceSBarry Smith #endif
66109c90ceSBarry Smith } _cb;
670fecffdcSJed Brown 
ourprestep(TS ts)68dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts)
69dd7ecb2fSBarry Smith {
701ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
71109c90ceSBarry Smith   void *ptr;
723ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
73109c90ceSBarry Smith #endif
741ab477f3SBarry Smith   PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
75dd7ecb2fSBarry Smith }
ourpoststep(TS ts)76dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts)
77dd7ecb2fSBarry Smith {
781ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
79109c90ceSBarry Smith   void *ptr;
803ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
81109c90ceSBarry Smith #endif
821ab477f3SBarry Smith   PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
83dd7ecb2fSBarry Smith }
ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,PetscCtx ctx)84*2a8381b2SBarry Smith static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, PetscCtx ctx)
85e2df7a95SSatish Balay {
861ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
87109c90ceSBarry Smith   void *ptr;
883ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
89109c90ceSBarry Smith #endif
901ab477f3SBarry Smith   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) */));
91e2df7a95SSatish Balay }
ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,PetscCtx ctx)92*2a8381b2SBarry Smith static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, PetscCtx ctx)
93e2df7a95SSatish Balay {
941ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
95109c90ceSBarry Smith   void *ptr;
963ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
97109c90ceSBarry Smith #endif
981ab477f3SBarry Smith   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) */));
990fecffdcSJed Brown }
ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat m,Mat p,PetscCtx ctx)100*2a8381b2SBarry Smith static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, PetscCtx ctx)
1010fecffdcSJed Brown {
1021ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
103109c90ceSBarry Smith   void *ptr;
1043ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
105109c90ceSBarry Smith #endif
1061ab477f3SBarry Smith   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) */));
1070fecffdcSJed Brown }
ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,Mat p,PetscCtx ctx)108*2a8381b2SBarry Smith static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, PetscCtx ctx)
1090fecffdcSJed Brown {
1101ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
111109c90ceSBarry Smith   void *ptr;
1123ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
113109c90ceSBarry Smith #endif
1141ab477f3SBarry Smith   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) */));
115e2df7a95SSatish Balay }
ourijacobianp(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,PetscCtx ctx)116*2a8381b2SBarry Smith static PetscErrorCode ourijacobianp(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, PetscCtx ctx)
1172549da82SBarry Smith {
1182549da82SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
1192549da82SBarry Smith   void *ptr;
1202549da82SBarry Smith   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
1212549da82SBarry Smith #endif
1222549da82SBarry Smith   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) */));
1232549da82SBarry Smith }
ourrhsjacobianp(TS ts,PetscReal d,Vec x,Mat m,PetscCtx ctx)124*2a8381b2SBarry Smith static PetscErrorCode ourrhsjacobianp(TS ts, PetscReal d, Vec x, Mat m, PetscCtx ctx)
1252549da82SBarry Smith {
1262549da82SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
1272549da82SBarry Smith   void *ptr;
1282549da82SBarry Smith   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
1292549da82SBarry Smith #endif
1302549da82SBarry Smith   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) */));
1312549da82SBarry Smith }
132e2df7a95SSatish Balay 
ourmonitordestroy(PetscCtxRt ctx)133*2a8381b2SBarry Smith static PetscErrorCode ourmonitordestroy(PetscCtxRt ctx)
134e2df7a95SSatish Balay {
135*2a8381b2SBarry Smith   TS ts = *(TS *)ctx;
136109c90ceSBarry Smith   PetscObjectUseFortranCallback(ts, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
137e2df7a95SSatish Balay }
138e2df7a95SSatish Balay 
139e2df7a95SSatish Balay /*
140e2df7a95SSatish Balay    Note ctx is the same as ts so we need to get the Fortran context out of the TS
141e2df7a95SSatish Balay */
ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,PetscCtx ctx)142*2a8381b2SBarry Smith static PetscErrorCode ourmonitor(TS ts, PetscInt i, PetscReal d, Vec v, PetscCtx ctx)
143e2df7a95SSatish Balay {
144109c90ceSBarry Smith   PetscObjectUseFortranCallback(ts, _cb.monitor, (TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), (&ts, &i, &d, &v, _ctx, &ierr));
145e2df7a95SSatish Balay }
146e2df7a95SSatish Balay 
147b0bc92c6SBarry Smith /*
148b0bc92c6SBarry Smith    Currently does not handle destroy or context
149b0bc92c6SBarry Smith */
ourtransform(PetscCtx ctx,Vec x,Vec * xout)150*2a8381b2SBarry Smith static PetscErrorCode ourtransform(PetscCtx ctx, Vec x, Vec *xout)
151b0bc92c6SBarry Smith {
152109c90ceSBarry Smith   PetscObjectUseFortranCallback((TS)ctx, _cb.transform, (void *, Vec *, Vec *, PetscErrorCode *), (_ctx, &x, xout, &ierr));
153b0bc92c6SBarry Smith }
154b0bc92c6SBarry Smith 
tsmonitorlgsettransform_(TS * ts,void (* f)(void *,Vec *,Vec *,PetscErrorCode *),PetscErrorCode (* d)(void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)155*2a8381b2SBarry Smith PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts, void (*f)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode (*d)(void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
156b0bc92c6SBarry Smith {
1575975b3b6SBarry Smith   *ierr = TSMonitorLGSetTransform(*ts, ourtransform, NULL, NULL);
1585975b3b6SBarry Smith   if (*ierr) return;
1595ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.transform, (PetscFortranCallbackFn *)f, ctx);
160b0bc92c6SBarry Smith }
161b0bc92c6SBarry Smith 
tssetprestep_(TS * ts,PetscErrorCode (* f)(TS *,PetscErrorCode *),PetscErrorCode * ierr)16219caf8f3SSatish Balay PETSC_EXTERN void tssetprestep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr)
163dd7ecb2fSBarry Smith {
1645975b3b6SBarry Smith   *ierr = TSSetPreStep(*ts, ourprestep);
1655975b3b6SBarry Smith   if (*ierr) return;
1665ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.prestep, (PetscFortranCallbackFn *)f, NULL);
167dd7ecb2fSBarry Smith }
168dd7ecb2fSBarry Smith 
tssetpoststep_(TS * ts,PetscErrorCode (* f)(TS *,PetscErrorCode *),PetscErrorCode * ierr)16919caf8f3SSatish Balay PETSC_EXTERN void tssetpoststep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr)
170dd7ecb2fSBarry Smith {
1715975b3b6SBarry Smith   *ierr = TSSetPostStep(*ts, ourpoststep);
1725975b3b6SBarry Smith   if (*ierr) return;
1735ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.poststep, (PetscFortranCallbackFn *)f, NULL);
174dd7ecb2fSBarry Smith }
175dd7ecb2fSBarry Smith 
1765b669ad3SBarry Smith PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *);
1775b669ad3SBarry Smith 
tssetrhsfunction_(TS * ts,Vec * r,void (* f)(TS *,PetscReal *,Vec *,Vec *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)1785ebfa9e9SBarry Smith PETSC_EXTERN void tssetrhsfunction_(TS *ts, Vec *r, void (*f)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
179e2df7a95SSatish Balay {
1800e4ef248SJed Brown   Vec R;
1810e4ef248SJed Brown   CHKFORTRANNULLOBJECT(r);
1820e4ef248SJed Brown   CHKFORTRANNULLFUNCTION(f);
1830298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
1845ebfa9e9SBarry Smith   if (f == tscomputerhsfunctionlinear_) {
1850e4ef248SJed Brown     *ierr = TSSetRHSFunction(*ts, R, TSComputeRHSFunctionLinear, fP);
1860e4ef248SJed Brown   } else {
1875ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsfunction, (PetscFortranCallbackFn *)f, fP);
188109c90ceSBarry Smith     *ierr = TSSetRHSFunction(*ts, R, ourrhsfunction, NULL);
1890fecffdcSJed Brown   }
1900fecffdcSJed Brown }
tsgetrhsfunction_(TS * ts,Vec * r,void * func,void ** ctx,PetscErrorCode * ierr)19119caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsfunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr)
19237698f3aSJed Brown {
19337698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
19437698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
1950298fd71SBarry Smith   *ierr = TSGetRHSFunction(*ts, r, NULL, ctx);
19637698f3aSJed Brown }
1970fecffdcSJed Brown 
198*2a8381b2SBarry Smith PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, Vec *F, PetscCtx ctx, PetscErrorCode *ierr);
199ce78bad3SBarry Smith 
tssetifunction_(TS * ts,Vec * r,void (* f)(TS *,PetscReal *,Vec *,Vec *,Vec *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)2005ebfa9e9SBarry Smith PETSC_EXTERN void tssetifunction_(TS *ts, Vec *r, void (*f)(TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
2010fecffdcSJed Brown {
2020fecffdcSJed Brown   Vec R;
2030fecffdcSJed Brown   CHKFORTRANNULLOBJECT(r);
2040fecffdcSJed Brown   CHKFORTRANNULLFUNCTION(f);
2050298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
2065ebfa9e9SBarry Smith   if (f == tscomputeifunctionlinear_) {
2070fecffdcSJed Brown     *ierr = TSSetIFunction(*ts, R, TSComputeIFunctionLinear, fP);
2080fecffdcSJed Brown   } else {
2095ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ifunction, (PetscFortranCallbackFn *)f, fP);
210109c90ceSBarry Smith     *ierr = TSSetIFunction(*ts, R, ourifunction, NULL);
2110e4ef248SJed Brown   }
212e2df7a95SSatish Balay }
tsgetifunction_(TS * ts,Vec * r,void * func,void ** ctx,PetscErrorCode * ierr)21319caf8f3SSatish Balay PETSC_EXTERN void tsgetifunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr)
21437698f3aSJed Brown {
21537698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
21637698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
2170298fd71SBarry Smith   *ierr = TSGetIFunction(*ts, r, NULL, ctx);
21837698f3aSJed Brown }
21926d46c62SHong Zhang 
2205b669ad3SBarry Smith PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *);
2215b669ad3SBarry Smith 
tssetrhsjacobian_(TS * ts,Mat * A,Mat * B,void (* f)(TS *,PetscReal *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)22219caf8f3SSatish Balay PETSC_EXTERN void tssetrhsjacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
223e2df7a95SSatish Balay {
224aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(f);
2255ebfa9e9SBarry Smith   if (f == tscomputerhsjacobianconstant_) {
2260e4ef248SJed Brown     *ierr = TSSetRHSJacobian(*ts, *A, *B, TSComputeRHSJacobianConstant, fP);
227e2df7a95SSatish Balay   } else {
2285ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobian, (PetscFortranCallbackFn *)f, fP);
229109c90ceSBarry Smith     *ierr = TSSetRHSJacobian(*ts, *A, *B, ourrhsjacobian, NULL);
2300fecffdcSJed Brown   }
2310fecffdcSJed Brown }
2320fecffdcSJed Brown 
233*2a8381b2SBarry Smith PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, PetscReal *shift, Mat *A, Mat *B, PetscCtx ctx, PetscErrorCode *ierr);
234ce78bad3SBarry Smith 
tssetijacobian_(TS * ts,Mat * A,Mat * B,void (* f)(TS *,PetscReal *,Vec *,Vec *,PetscReal *,Mat *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)2355ebfa9e9SBarry Smith 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)
2360fecffdcSJed Brown {
237aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(f);
2385ebfa9e9SBarry Smith   if (f == tscomputeijacobianconstant_) {
2390fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts, *A, *B, TSComputeIJacobianConstant, fP);
2400fecffdcSJed Brown   } else {
2415ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobian, (PetscFortranCallbackFn *)f, fP);
242109c90ceSBarry Smith     *ierr = TSSetIJacobian(*ts, *A, *B, ourijacobian, NULL);
243e2df7a95SSatish Balay   }
244e2df7a95SSatish Balay }
tsgetijacobian_(TS * ts,Mat * J,Mat * M,int * func,void ** ctx,PetscErrorCode * ierr)24519caf8f3SSatish Balay PETSC_EXTERN void tsgetijacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr)
24637698f3aSJed Brown {
24737698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
24837698f3aSJed Brown   CHKFORTRANNULLOBJECT(J);
24937698f3aSJed Brown   CHKFORTRANNULLOBJECT(M);
250dfef5ea7SSatish Balay   *ierr = TSGetIJacobian(*ts, J, M, NULL, ctx);
25137698f3aSJed Brown }
tssetijacobianp_(TS * ts,Mat * A,void (* f)(TS *,PetscReal *,Vec *,Vec *,PetscReal,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)2522549da82SBarry Smith PETSC_EXTERN void tssetijacobianp_(TS *ts, Mat *A, void (*f)(TS *, PetscReal *, Vec *, Vec *, PetscReal, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
2532549da82SBarry Smith {
2542549da82SBarry Smith   CHKFORTRANNULLFUNCTION(f);
2555ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobianp, (PetscFortranCallbackFn *)f, fP);
2562549da82SBarry Smith   *ierr = TSSetIJacobianP(*ts, *A, ourijacobianp, NULL);
2572549da82SBarry Smith }
tsgetijacobianp_(TS * ts,Mat * J,int * func,void ** ctx,PetscErrorCode * ierr)2582549da82SBarry Smith PETSC_EXTERN void tsgetijacobianp_(TS *ts, Mat *J, int *func, void **ctx, PetscErrorCode *ierr)
2592549da82SBarry Smith {
2602549da82SBarry Smith   CHKFORTRANNULLINTEGER(ctx);
2612549da82SBarry Smith   CHKFORTRANNULLOBJECT(J);
2622549da82SBarry Smith   *ierr = TSGetIJacobianP(*ts, J, NULL, ctx);
2632549da82SBarry Smith }
tssetrhsjacobianp_(TS * ts,Mat * A,void (* f)(TS *,PetscReal *,Vec *,Mat *,void *,PetscErrorCode *),void * fP,PetscErrorCode * ierr)2642549da82SBarry Smith PETSC_EXTERN void tssetrhsjacobianp_(TS *ts, Mat *A, void (*f)(TS *, PetscReal *, Vec *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr)
2652549da82SBarry Smith {
2662549da82SBarry Smith   CHKFORTRANNULLFUNCTION(f);
2675ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobianp, (PetscFortranCallbackFn *)f, fP);
2682549da82SBarry Smith   *ierr = TSSetRHSJacobianP(*ts, *A, ourrhsjacobianp, NULL);
2692549da82SBarry Smith }
tsgetrhsjacobianp_(TS * ts,Mat * J,int * func,void ** ctx,PetscErrorCode * ierr)2702549da82SBarry Smith PETSC_EXTERN void tsgetrhsjacobianp_(TS *ts, Mat *J, int *func, void **ctx, PetscErrorCode *ierr)
2712549da82SBarry Smith {
2722549da82SBarry Smith   CHKFORTRANNULLINTEGER(ctx);
2732549da82SBarry Smith   CHKFORTRANNULLOBJECT(J);
2742549da82SBarry Smith   *ierr = TSGetRHSJacobianP(*ts, J, NULL, ctx);
2752549da82SBarry Smith }
276e2df7a95SSatish Balay 
277ce78bad3SBarry Smith PETSC_EXTERN void tsmonitordefault_(TS *, PetscInt *, PetscReal *, Vec *, PetscViewerAndFormat **, PetscErrorCode *);
27852f0073cSBarry Smith 
27919caf8f3SSatish Balay /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */
280e2df7a95SSatish Balay 
tsmonitorset_(TS * ts,void (* func)(TS *,PetscInt *,PetscReal *,Vec *,void *,PetscErrorCode *),void * mctx,void (* d)(void *,PetscErrorCode *),PetscErrorCode * ierr)28119caf8f3SSatish Balay PETSC_EXTERN void tsmonitorset_(TS *ts, void (*func)(TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), void *mctx, void (*d)(void *, PetscErrorCode *), PetscErrorCode *ierr)
282e2df7a95SSatish Balay {
283aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(d);
2845ebfa9e9SBarry Smith   if ((PetscFortranCallbackFn *)func == (PetscFortranCallbackFn *)tsmonitordefault_) {
28549abdd8aSBarry Smith     *ierr = TSMonitorSet(*ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
286e2df7a95SSatish Balay   } else {
2875ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscFortranCallbackFn *)func, mctx);
2885ebfa9e9SBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscFortranCallbackFn *)d, mctx);
289aecf964fSBarry Smith     *ierr = TSMonitorSet(*ts, ourmonitor, *ts, ourmonitordestroy);
290e2df7a95SSatish Balay   }
291e2df7a95SSatish Balay }
292e2df7a95SSatish Balay 
293089b2837SJed Brown /*  func is currently ignored from Fortran */
tsgetrhsjacobian_(TS * ts,Mat * J,Mat * M,int * func,void ** ctx,PetscErrorCode * ierr)29419caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsjacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr)
295e2df7a95SSatish Balay {
296dfef5ea7SSatish Balay   *ierr = TSGetRHSJacobian(*ts, J, M, NULL, ctx);
297e2df7a95SSatish Balay }
298