1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h> 2c6db04a5SJed Brown #include <petscts.h> 3665c2dedSJed Brown #include <petscviewer.h> 48ad9143bSBarry Smith #include <petsc/private/f90impl.h> 5e2df7a95SSatish Balay 6e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS) 7b0bc92c6SBarry Smith #define tsmonitorlgsettransform_ TSMONITORLGSETTRANSFORM 8e2df7a95SSatish Balay #define tssetrhsfunction_ TSSETRHSFUNCTION 937698f3aSJed Brown #define tsgetrhsfunction_ TSGETRHSFUNCTION 10e2df7a95SSatish Balay #define tssetrhsjacobian_ TSSETRHSJACOBIAN 11e2df7a95SSatish Balay #define tsgetrhsjacobian_ TSGETRHSJACOBIAN 1237698f3aSJed Brown #define tssetifunction_ TSSETIFUNCTION 1337698f3aSJed Brown #define tsgetifunction_ TSGETIFUNCTION 1437698f3aSJed Brown #define tssetijacobian_ TSSETIJACOBIAN 1537698f3aSJed Brown #define tsgetijacobian_ TSGETIJACOBIAN 16a6570f20SBarry Smith #define tsmonitorset_ TSMONITORSET 170e4ef248SJed Brown #define tscomputerhsfunctionlinear_ TSCOMPUTERHSFUNCTIONLINEAR 180e4ef248SJed Brown #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT 190fecffdcSJed Brown #define tscomputeifunctionlinear_ TSCOMPUTEIFUNCTIONLINEAR 200fecffdcSJed Brown #define tscomputeijacobianconstant_ TSCOMPUTEIJACOBIANCONSTANT 21a6570f20SBarry Smith #define tsmonitordefault_ TSMONITORDEFAULT 22dd7ecb2fSBarry Smith #define tssetprestep_ TSSETPRESTEP 23dd7ecb2fSBarry Smith #define tssetpoststep_ TSSETPOSTSTEP 24e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 25b0bc92c6SBarry Smith #define tsmonitorlgsettransform_ tsmonitorlgsettransform 26e2df7a95SSatish Balay #define tssetrhsfunction_ tssetrhsfunction 2737698f3aSJed Brown #define tsgetrhsfunction_ tsgetrhsfunction 28e2df7a95SSatish Balay #define tssetrhsjacobian_ tssetrhsjacobian 29e2df7a95SSatish Balay #define tsgetrhsjacobian_ tsgetrhsjacobian 3037698f3aSJed Brown #define tssetifunction_ tssetifunction 3137698f3aSJed Brown #define tsgetifunction_ tsgetifunction 3237698f3aSJed Brown #define tssetijacobian_ tssetijacobian 3337698f3aSJed Brown #define tsgetijacobian_ tsgetijacobian 34a6570f20SBarry Smith #define tsmonitorset_ tsmonitorset 350e4ef248SJed Brown #define tscomputerhsfunctionlinear_ tscomputerhsfunctionlinear 360e4ef248SJed Brown #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant 370fecffdcSJed Brown #define tscomputeifunctionlinear_ tscomputeifunctionlinear 380fecffdcSJed Brown #define tscomputeijacobianconstant_ tscomputeijacobianconstant 39a6570f20SBarry Smith #define tsmonitordefault_ tsmonitordefault 40dd7ecb2fSBarry Smith #define tssetprestep_ tssetprestep 41dd7ecb2fSBarry Smith #define tssetpoststep_ tssetpoststep 42e2df7a95SSatish Balay #endif 43e2df7a95SSatish Balay 44109c90ceSBarry Smith static struct { 45109c90ceSBarry Smith PetscFortranCallbackId prestep; 46109c90ceSBarry Smith PetscFortranCallbackId poststep; 47109c90ceSBarry Smith PetscFortranCallbackId rhsfunction; 48109c90ceSBarry Smith PetscFortranCallbackId rhsjacobian; 49109c90ceSBarry Smith PetscFortranCallbackId ifunction; 50109c90ceSBarry Smith PetscFortranCallbackId ijacobian; 51109c90ceSBarry Smith PetscFortranCallbackId monitor; 52109c90ceSBarry Smith PetscFortranCallbackId mondestroy; 53109c90ceSBarry Smith PetscFortranCallbackId transform; 54109c90ceSBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) 55109c90ceSBarry Smith PetscFortranCallbackId function_pgiptr; 56109c90ceSBarry Smith #endif 57109c90ceSBarry Smith } _cb; 580fecffdcSJed Brown 59dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts) 60dd7ecb2fSBarry Smith { 611ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 62109c90ceSBarry Smith void *ptr; 633ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 64109c90ceSBarry Smith #endif 651ab477f3SBarry Smith PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 66dd7ecb2fSBarry Smith } 67dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts) 68dd7ecb2fSBarry Smith { 691ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 70109c90ceSBarry Smith void *ptr; 713ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 72109c90ceSBarry Smith #endif 731ab477f3SBarry Smith PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 74dd7ecb2fSBarry Smith } 750fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, void *ctx) 76e2df7a95SSatish Balay { 771ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 78109c90ceSBarry Smith void *ptr; 793ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 80109c90ceSBarry Smith #endif 811ab477f3SBarry 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) */)); 82e2df7a95SSatish Balay } 830fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, void *ctx) 84e2df7a95SSatish Balay { 851ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 86109c90ceSBarry Smith void *ptr; 873ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 88109c90ceSBarry Smith #endif 891ab477f3SBarry 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) */)); 900fecffdcSJed Brown } 91d1e9a80fSBarry Smith static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, void *ctx) 920fecffdcSJed Brown { 931ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 94109c90ceSBarry Smith void *ptr; 953ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 96109c90ceSBarry Smith #endif 971ab477f3SBarry 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) */)); 980fecffdcSJed Brown } 99d1e9a80fSBarry Smith static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, void *ctx) 1000fecffdcSJed Brown { 1011ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 102109c90ceSBarry Smith void *ptr; 1033ba16761SJacob Faibussowitsch PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 104109c90ceSBarry Smith #endif 1051ab477f3SBarry 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) */)); 106e2df7a95SSatish Balay } 107e2df7a95SSatish Balay 108c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx) 109e2df7a95SSatish Balay { 110109c90ceSBarry Smith TS ts = (TS)*ctx; 111109c90ceSBarry Smith PetscObjectUseFortranCallback(ts, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr)); 112e2df7a95SSatish Balay } 113e2df7a95SSatish Balay 114e2df7a95SSatish Balay /* 115e2df7a95SSatish Balay Note ctx is the same as ts so we need to get the Fortran context out of the TS 116e2df7a95SSatish Balay */ 1170fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts, PetscInt i, PetscReal d, Vec v, void *ctx) 118e2df7a95SSatish Balay { 119109c90ceSBarry Smith PetscObjectUseFortranCallback(ts, _cb.monitor, (TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), (&ts, &i, &d, &v, _ctx, &ierr)); 120e2df7a95SSatish Balay } 121e2df7a95SSatish Balay 122b0bc92c6SBarry Smith /* 123b0bc92c6SBarry Smith Currently does not handle destroy or context 124b0bc92c6SBarry Smith */ 125b0bc92c6SBarry Smith static PetscErrorCode ourtransform(void *ctx, Vec x, Vec *xout) 126b0bc92c6SBarry Smith { 127109c90ceSBarry Smith PetscObjectUseFortranCallback((TS)ctx, _cb.transform, (void *, Vec *, Vec *, PetscErrorCode *), (_ctx, &x, xout, &ierr)); 128b0bc92c6SBarry Smith } 129b0bc92c6SBarry Smith 13019caf8f3SSatish Balay PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts, void (*f)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode (*d)(void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 131b0bc92c6SBarry Smith { 1325975b3b6SBarry Smith *ierr = TSMonitorLGSetTransform(*ts, ourtransform, NULL, NULL); 1335975b3b6SBarry Smith if (*ierr) return; 1348434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.transform, (PetscVoidFn *)f, ctx); 135b0bc92c6SBarry Smith } 136b0bc92c6SBarry Smith 13719caf8f3SSatish Balay PETSC_EXTERN void tssetprestep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 138dd7ecb2fSBarry Smith { 1395975b3b6SBarry Smith *ierr = TSSetPreStep(*ts, ourprestep); 1405975b3b6SBarry Smith if (*ierr) return; 1418434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.prestep, (PetscVoidFn *)f, NULL); 142dd7ecb2fSBarry Smith } 143dd7ecb2fSBarry Smith 14419caf8f3SSatish Balay PETSC_EXTERN void tssetpoststep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 145dd7ecb2fSBarry Smith { 1465975b3b6SBarry Smith *ierr = TSSetPostStep(*ts, ourpoststep); 1475975b3b6SBarry Smith if (*ierr) return; 1488434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.poststep, (PetscVoidFn *)f, NULL); 149dd7ecb2fSBarry Smith } 150dd7ecb2fSBarry Smith 15169c1e2abSSatish Balay PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *F, void *ctx, PetscErrorCode *ierr) 1520e4ef248SJed Brown { 1530e4ef248SJed Brown *ierr = TSComputeRHSFunctionLinear(*ts, *t, *X, *F, ctx); 1540e4ef248SJed Brown } 15519caf8f3SSatish Balay PETSC_EXTERN void tssetrhsfunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 156e2df7a95SSatish Balay { 1570e4ef248SJed Brown Vec R; 1580e4ef248SJed Brown CHKFORTRANNULLOBJECT(r); 1590e4ef248SJed Brown CHKFORTRANNULLFUNCTION(f); 1600298fd71SBarry Smith R = r ? *r : (Vec)NULL; 1618434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsfunctionlinear_) { 1620e4ef248SJed Brown *ierr = TSSetRHSFunction(*ts, R, TSComputeRHSFunctionLinear, fP); 1630e4ef248SJed Brown } else { 1648434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsfunction, (PetscVoidFn *)f, fP); 165109c90ceSBarry Smith *ierr = TSSetRHSFunction(*ts, R, ourrhsfunction, NULL); 1660fecffdcSJed Brown } 1670fecffdcSJed Brown } 16819caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsfunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 16937698f3aSJed Brown { 17037698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 17137698f3aSJed Brown CHKFORTRANNULLOBJECT(r); 1720298fd71SBarry Smith *ierr = TSGetRHSFunction(*ts, r, NULL, ctx); 17337698f3aSJed Brown } 1740fecffdcSJed Brown 175*ce78bad3SBarry Smith PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, Vec *F, void *ctx, PetscErrorCode *ierr); 176*ce78bad3SBarry Smith 17719caf8f3SSatish Balay PETSC_EXTERN void tssetifunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 1780fecffdcSJed Brown { 1790fecffdcSJed Brown Vec R; 1800fecffdcSJed Brown CHKFORTRANNULLOBJECT(r); 1810fecffdcSJed Brown CHKFORTRANNULLFUNCTION(f); 1820298fd71SBarry Smith R = r ? *r : (Vec)NULL; 1838434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeifunctionlinear_) { 1840fecffdcSJed Brown *ierr = TSSetIFunction(*ts, R, TSComputeIFunctionLinear, fP); 1850fecffdcSJed Brown } else { 1868434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ifunction, (PetscVoidFn *)f, fP); 187109c90ceSBarry Smith *ierr = TSSetIFunction(*ts, R, ourifunction, NULL); 1880e4ef248SJed Brown } 189e2df7a95SSatish Balay } 19019caf8f3SSatish Balay PETSC_EXTERN void tsgetifunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 19137698f3aSJed Brown { 19237698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 19337698f3aSJed Brown CHKFORTRANNULLOBJECT(r); 1940298fd71SBarry Smith *ierr = TSGetIFunction(*ts, r, NULL, ctx); 19537698f3aSJed Brown } 19626d46c62SHong Zhang 197e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 198d1e9a80fSBarry Smith PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts, PetscReal *t, Vec *X, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr) 1990e4ef248SJed Brown { 200d1e9a80fSBarry Smith *ierr = TSComputeRHSJacobianConstant(*ts, *t, *X, *A, *B, ctx); 2010e4ef248SJed Brown } 20219caf8f3SSatish Balay PETSC_EXTERN void tssetrhsjacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 203e2df7a95SSatish Balay { 204aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(f); 2058434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsjacobianconstant_) { 2060e4ef248SJed Brown *ierr = TSSetRHSJacobian(*ts, *A, *B, TSComputeRHSJacobianConstant, fP); 207e2df7a95SSatish Balay } else { 2088434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobian, (PetscVoidFn *)f, fP); 209109c90ceSBarry Smith *ierr = TSSetRHSJacobian(*ts, *A, *B, ourrhsjacobian, NULL); 2100fecffdcSJed Brown } 2110fecffdcSJed Brown } 2120fecffdcSJed Brown 213*ce78bad3SBarry Smith PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, PetscReal *shift, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr); 214*ce78bad3SBarry Smith 21519caf8f3SSatish Balay PETSC_EXTERN void tssetijacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 2160fecffdcSJed Brown { 217aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(f); 2188434afd1SBarry Smith if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeijacobianconstant_) { 2190fecffdcSJed Brown *ierr = TSSetIJacobian(*ts, *A, *B, TSComputeIJacobianConstant, fP); 2200fecffdcSJed Brown } else { 2218434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobian, (PetscVoidFn *)f, fP); 222109c90ceSBarry Smith *ierr = TSSetIJacobian(*ts, *A, *B, ourijacobian, NULL); 223e2df7a95SSatish Balay } 224e2df7a95SSatish Balay } 22519caf8f3SSatish Balay PETSC_EXTERN void tsgetijacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 22637698f3aSJed Brown { 22737698f3aSJed Brown CHKFORTRANNULLINTEGER(ctx); 22837698f3aSJed Brown CHKFORTRANNULLOBJECT(J); 22937698f3aSJed Brown CHKFORTRANNULLOBJECT(M); 230dfef5ea7SSatish Balay *ierr = TSGetIJacobian(*ts, J, M, NULL, ctx); 23137698f3aSJed Brown } 232e2df7a95SSatish Balay 233*ce78bad3SBarry Smith PETSC_EXTERN void tsmonitordefault_(TS *, PetscInt *, PetscReal *, Vec *, PetscViewerAndFormat **, PetscErrorCode *); 23452f0073cSBarry Smith 235e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 236e2df7a95SSatish Balay 23719caf8f3SSatish Balay /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */ 238e2df7a95SSatish Balay 23919caf8f3SSatish Balay PETSC_EXTERN void tsmonitorset_(TS *ts, void (*func)(TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), void *mctx, void (*d)(void *, PetscErrorCode *), PetscErrorCode *ierr) 240e2df7a95SSatish Balay { 241aecf964fSBarry Smith CHKFORTRANNULLFUNCTION(d); 2428434afd1SBarry Smith if ((PetscVoidFn *)func == (PetscVoidFn *)tsmonitordefault_) { 24349abdd8aSBarry Smith *ierr = TSMonitorSet(*ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy); 244e2df7a95SSatish Balay } else { 2458434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscVoidFn *)func, mctx); 2468434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscVoidFn *)d, mctx); 247aecf964fSBarry Smith *ierr = TSMonitorSet(*ts, ourmonitor, *ts, ourmonitordestroy); 248e2df7a95SSatish Balay } 249e2df7a95SSatish Balay } 250e2df7a95SSatish Balay 251e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 252089b2837SJed Brown /* func is currently ignored from Fortran */ 25319caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsjacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 254e2df7a95SSatish Balay { 255dfef5ea7SSatish Balay *ierr = TSGetRHSJacobian(*ts, J, M, NULL, ctx); 256e2df7a95SSatish Balay } 257