1*6dd63270SBarry 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 160e4ef248SJed Brown #define tscomputerhsfunctionlinear_ TSCOMPUTERHSFUNCTIONLINEAR 170e4ef248SJed Brown #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT 180fecffdcSJed Brown #define tscomputeifunctionlinear_ TSCOMPUTEIFUNCTIONLINEAR 190fecffdcSJed Brown #define tscomputeijacobianconstant_ TSCOMPUTEIJACOBIANCONSTANT 20a6570f20SBarry Smith #define tsmonitordefault_ TSMONITORDEFAULT 21dd7ecb2fSBarry Smith #define tssetprestep_ TSSETPRESTEP 22dd7ecb2fSBarry Smith #define tssetpoststep_ TSSETPOSTSTEP 23e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 24b0bc92c6SBarry Smith #define tsmonitorlgsettransform_ tsmonitorlgsettransform 25e2df7a95SSatish Balay #define tssetrhsfunction_ tssetrhsfunction 2637698f3aSJed Brown #define tsgetrhsfunction_ tsgetrhsfunction 27e2df7a95SSatish Balay #define tssetrhsjacobian_ tssetrhsjacobian 28e2df7a95SSatish Balay #define tsgetrhsjacobian_ tsgetrhsjacobian 2937698f3aSJed Brown #define tssetifunction_ tssetifunction 3037698f3aSJed Brown #define tsgetifunction_ tsgetifunction 3137698f3aSJed Brown #define tssetijacobian_ tssetijacobian 3237698f3aSJed Brown #define tsgetijacobian_ tsgetijacobian 33a6570f20SBarry Smith #define tsmonitorset_ tsmonitorset 340e4ef248SJed Brown #define tscomputerhsfunctionlinear_ tscomputerhsfunctionlinear 350e4ef248SJed Brown #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant 360fecffdcSJed Brown #define tscomputeifunctionlinear_ tscomputeifunctionlinear 370fecffdcSJed Brown #define tscomputeijacobianconstant_ tscomputeijacobianconstant 38a6570f20SBarry Smith #define tsmonitordefault_ tsmonitordefault 39dd7ecb2fSBarry Smith #define tssetprestep_ tssetprestep 40dd7ecb2fSBarry Smith #define tssetpoststep_ tssetpoststep 41e2df7a95SSatish Balay #endif 42e2df7a95SSatish Balay 43109c90ceSBarry Smith static struct { 44109c90ceSBarry Smith PetscFortranCallbackId prestep; 45109c90ceSBarry Smith PetscFortranCallbackId poststep; 46109c90ceSBarry Smith PetscFortranCallbackId rhsfunction; 47109c90ceSBarry Smith PetscFortranCallbackId rhsjacobian; 48109c90ceSBarry Smith PetscFortranCallbackId ifunction; 49109c90ceSBarry Smith PetscFortranCallbackId ijacobian; 50109c90ceSBarry Smith PetscFortranCallbackId monitor; 51109c90ceSBarry Smith PetscFortranCallbackId mondestroy; 52109c90ceSBarry Smith PetscFortranCallbackId transform; 53109c90ceSBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) 54109c90ceSBarry Smith PetscFortranCallbackId function_pgiptr; 55109c90ceSBarry Smith #endif 56109c90ceSBarry Smith } _cb; 570fecffdcSJed Brown 58dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts) 59dd7ecb2fSBarry Smith { 601ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 61109c90ceSBarry Smith void *ptr; 623ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 63109c90ceSBarry Smith #endif 641ab477f3SBarry Smith PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 65dd7ecb2fSBarry Smith } 66dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts) 67dd7ecb2fSBarry Smith { 681ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 69109c90ceSBarry Smith void *ptr; 703ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 71109c90ceSBarry Smith #endif 721ab477f3SBarry Smith PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 73dd7ecb2fSBarry Smith } 740fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, void *ctx) 75e2df7a95SSatish Balay { 761ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 77109c90ceSBarry Smith void *ptr; 783ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 79109c90ceSBarry Smith #endif 801ab477f3SBarry 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) */)); 81e2df7a95SSatish Balay } 820fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, void *ctx) 83e2df7a95SSatish Balay { 841ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 85109c90ceSBarry Smith void *ptr; 863ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 87109c90ceSBarry Smith #endif 881ab477f3SBarry 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) */)); 890fecffdcSJed Brown } 90d1e9a80fSBarry Smith static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, void *ctx) 910fecffdcSJed Brown { 921ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 93109c90ceSBarry Smith void *ptr; 943ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 95109c90ceSBarry Smith #endif 961ab477f3SBarry 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) */)); 970fecffdcSJed Brown } 98d1e9a80fSBarry Smith static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, void *ctx) 990fecffdcSJed Brown { 1001ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 101109c90ceSBarry Smith void *ptr; 1023ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 103109c90ceSBarry Smith #endif 1041ab477f3SBarry 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) */)); 105e2df7a95SSatish Balay } 106e2df7a95SSatish Balay 107c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx) 108e2df7a95SSatish Balay { 109109c90ceSBarry Smith TS ts = (TS)*ctx; 110109c90ceSBarry Smith PetscObjectUseFortranCallback(ts, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr)); 111e2df7a95SSatish Balay } 112e2df7a95SSatish Balay 113e2df7a95SSatish Balay /* 114e2df7a95SSatish Balay Note ctx is the same as ts so we need to get the Fortran context out of the TS 115e2df7a95SSatish Balay */ 1160fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts, PetscInt i, PetscReal d, Vec v, void *ctx) 117e2df7a95SSatish Balay { 118109c90ceSBarry Smith PetscObjectUseFortranCallback(ts, _cb.monitor, (TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), (&ts, &i, &d, &v, _ctx, &ierr)); 119e2df7a95SSatish Balay } 120e2df7a95SSatish Balay 121b0bc92c6SBarry Smith /* 122b0bc92c6SBarry Smith Currently does not handle destroy or context 123b0bc92c6SBarry Smith */ 124b0bc92c6SBarry Smith static PetscErrorCode ourtransform(void *ctx, Vec x, Vec *xout) 125b0bc92c6SBarry Smith { 126109c90ceSBarry Smith PetscObjectUseFortranCallback((TS)ctx, _cb.transform, (void *, Vec *, Vec *, PetscErrorCode *), (_ctx, &x, xout, &ierr)); 127b0bc92c6SBarry Smith } 128b0bc92c6SBarry Smith 12919caf8f3SSatish Balay PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts, void (*f)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode (*d)(void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 130b0bc92c6SBarry Smith { 1315975b3b6SBarry Smith *ierr = TSMonitorLGSetTransform(*ts, ourtransform, NULL, NULL); 1325975b3b6SBarry Smith if (*ierr) return; 1338434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.transform, (PetscVoidFn *)f, ctx); 134b0bc92c6SBarry Smith } 135b0bc92c6SBarry Smith 13619caf8f3SSatish Balay PETSC_EXTERN void tssetprestep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 137dd7ecb2fSBarry Smith { 1385975b3b6SBarry Smith *ierr = TSSetPreStep(*ts, ourprestep); 1395975b3b6SBarry Smith if (*ierr) return; 1408434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.prestep, (PetscVoidFn *)f, NULL); 141dd7ecb2fSBarry Smith } 142dd7ecb2fSBarry Smith 14319caf8f3SSatish Balay PETSC_EXTERN void tssetpoststep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 144dd7ecb2fSBarry Smith { 1455975b3b6SBarry Smith *ierr = TSSetPostStep(*ts, ourpoststep); 1465975b3b6SBarry Smith if (*ierr) return; 1478434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.poststep, (PetscVoidFn *)f, NULL); 148dd7ecb2fSBarry Smith } 149dd7ecb2fSBarry Smith 15069c1e2abSSatish Balay PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *F, void *ctx, PetscErrorCode *ierr) 1510e4ef248SJed Brown { 1520e4ef248SJed Brown *ierr = TSComputeRHSFunctionLinear(*ts, *t, *X, *F, ctx); 1530e4ef248SJed Brown } 15419caf8f3SSatish Balay PETSC_EXTERN void tssetrhsfunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 155e2df7a95SSatish Balay { 1560e4ef248SJed Brown Vec R; 1570e4ef248SJed Brown CHKFORTRANNULLOBJECT(r); 1580e4ef248SJed Brown CHKFORTRANNULLFUNCTION(f); 1590298fd71SBarry Smith R = r ? *r : (Vec)NULL; 1608434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsfunctionlinear_) { 1610e4ef248SJed Brown *ierr = TSSetRHSFunction(*ts, R, TSComputeRHSFunctionLinear, fP); 1620e4ef248SJed Brown } else { 1638434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsfunction, (PetscVoidFn *)f, fP); 164109c90ceSBarry Smith *ierr = TSSetRHSFunction(*ts, R, ourrhsfunction, NULL); 1650fecffdcSJed Brown } 1660fecffdcSJed Brown } 16719caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsfunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 16837698f3aSJed Brown { 16937698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 17037698f3aSJed Brown CHKFORTRANNULLOBJECT(r); 1710298fd71SBarry Smith *ierr = TSGetRHSFunction(*ts, r, NULL, ctx); 17237698f3aSJed Brown } 1730fecffdcSJed Brown 174ce78bad3SBarry Smith PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, Vec *F, void *ctx, PetscErrorCode *ierr); 175ce78bad3SBarry Smith 17619caf8f3SSatish Balay PETSC_EXTERN void tssetifunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 1770fecffdcSJed Brown { 1780fecffdcSJed Brown Vec R; 1790fecffdcSJed Brown CHKFORTRANNULLOBJECT(r); 1800fecffdcSJed Brown CHKFORTRANNULLFUNCTION(f); 1810298fd71SBarry Smith R = r ? *r : (Vec)NULL; 1828434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeifunctionlinear_) { 1830fecffdcSJed Brown *ierr = TSSetIFunction(*ts, R, TSComputeIFunctionLinear, fP); 1840fecffdcSJed Brown } else { 1858434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ifunction, (PetscVoidFn *)f, fP); 186109c90ceSBarry Smith *ierr = TSSetIFunction(*ts, R, ourifunction, NULL); 1870e4ef248SJed Brown } 188e2df7a95SSatish Balay } 18919caf8f3SSatish Balay PETSC_EXTERN void tsgetifunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 19037698f3aSJed Brown { 19137698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 19237698f3aSJed Brown CHKFORTRANNULLOBJECT(r); 1930298fd71SBarry Smith *ierr = TSGetIFunction(*ts, r, NULL, ctx); 19437698f3aSJed Brown } 19526d46c62SHong Zhang 196e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 197d1e9a80fSBarry Smith PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts, PetscReal *t, Vec *X, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr) 1980e4ef248SJed Brown { 199d1e9a80fSBarry Smith *ierr = TSComputeRHSJacobianConstant(*ts, *t, *X, *A, *B, ctx); 2000e4ef248SJed Brown } 20119caf8f3SSatish Balay PETSC_EXTERN void tssetrhsjacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 202e2df7a95SSatish Balay { 203aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(f); 2048434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsjacobianconstant_) { 2050e4ef248SJed Brown *ierr = TSSetRHSJacobian(*ts, *A, *B, TSComputeRHSJacobianConstant, fP); 206e2df7a95SSatish Balay } else { 2078434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobian, (PetscVoidFn *)f, fP); 208109c90ceSBarry Smith *ierr = TSSetRHSJacobian(*ts, *A, *B, ourrhsjacobian, NULL); 2090fecffdcSJed Brown } 2100fecffdcSJed Brown } 2110fecffdcSJed Brown 212ce78bad3SBarry Smith PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, PetscReal *shift, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr); 213ce78bad3SBarry Smith 21419caf8f3SSatish Balay PETSC_EXTERN void tssetijacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 2150fecffdcSJed Brown { 216aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(f); 2178434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeijacobianconstant_) { 2180fecffdcSJed Brown *ierr = TSSetIJacobian(*ts, *A, *B, TSComputeIJacobianConstant, fP); 2190fecffdcSJed Brown } else { 2208434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobian, (PetscVoidFn *)f, fP); 221109c90ceSBarry Smith *ierr = TSSetIJacobian(*ts, *A, *B, ourijacobian, NULL); 222e2df7a95SSatish Balay } 223e2df7a95SSatish Balay } 22419caf8f3SSatish Balay PETSC_EXTERN void tsgetijacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 22537698f3aSJed Brown { 22637698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 22737698f3aSJed Brown CHKFORTRANNULLOBJECT(J); 22837698f3aSJed Brown CHKFORTRANNULLOBJECT(M); 229dfef5ea7SSatish Balay *ierr = TSGetIJacobian(*ts, J, M, NULL, ctx); 23037698f3aSJed Brown } 231e2df7a95SSatish Balay 232ce78bad3SBarry Smith PETSC_EXTERN void tsmonitordefault_(TS *, PetscInt *, PetscReal *, Vec *, PetscViewerAndFormat **, PetscErrorCode *); 23352f0073cSBarry Smith 234e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 235e2df7a95SSatish Balay 23619caf8f3SSatish Balay /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */ 237e2df7a95SSatish Balay 23819caf8f3SSatish Balay PETSC_EXTERN void tsmonitorset_(TS *ts, void (*func)(TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), void *mctx, void (*d)(void *, PetscErrorCode *), PetscErrorCode *ierr) 239e2df7a95SSatish Balay { 240aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(d); 2418434afd1SBarry Smith if ((PetscVoidFn *)func == (PetscVoidFn *)tsmonitordefault_) { 24249abdd8aSBarry Smith *ierr = TSMonitorSet(*ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy); 243e2df7a95SSatish Balay } else { 2448434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscVoidFn *)func, mctx); 2458434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscVoidFn *)d, mctx); 246aecf964fSBarry Smith *ierr = TSMonitorSet(*ts, ourmonitor, *ts, ourmonitordestroy); 247e2df7a95SSatish Balay } 248e2df7a95SSatish Balay } 249e2df7a95SSatish Balay 250e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 251089b2837SJed Brown /* func is currently ignored from Fortran */ 25219caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsjacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 253e2df7a95SSatish Balay { 254dfef5ea7SSatish Balay *ierr = TSGetRHSJacobian(*ts, J, M, NULL, ctx); 255e2df7a95SSatish Balay } 256