xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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
16e2df7a95SSatish Balay #define tsview_                              TSVIEW
1737698f3aSJed Brown #define tssetoptionsprefix_                  TSSETOPTIONSPREFIX
18e2df7a95SSatish Balay #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
1937698f3aSJed Brown #define tsappendoptionsprefix_               TSAPPENDOPTIONSPREFIX
20a6570f20SBarry Smith #define tsmonitorset_                        TSMONITORSET
210e4ef248SJed Brown #define tscomputerhsfunctionlinear_          TSCOMPUTERHSFUNCTIONLINEAR
220e4ef248SJed Brown #define tscomputerhsjacobianconstant_        TSCOMPUTERHSJACOBIANCONSTANT
230fecffdcSJed Brown #define tscomputeifunctionlinear_            TSCOMPUTEIFUNCTIONLINEAR
240fecffdcSJed Brown #define tscomputeijacobianconstant_          TSCOMPUTEIJACOBIANCONSTANT
25a6570f20SBarry Smith #define tsmonitordefault_                    TSMONITORDEFAULT
26dd7ecb2fSBarry Smith #define tssetprestep_                        TSSETPRESTEP
27dd7ecb2fSBarry Smith #define tssetpoststep_                       TSSETPOSTSTEP
28fe2efc57SMark #define tsviewfromoptions_                   TSVIEWFROMOPTIONS
29e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
30b0bc92c6SBarry Smith #define tsmonitorlgsettransform_             tsmonitorlgsettransform
31e2df7a95SSatish Balay #define tssetrhsfunction_                    tssetrhsfunction
3237698f3aSJed Brown #define tsgetrhsfunction_                    tsgetrhsfunction
33e2df7a95SSatish Balay #define tssetrhsjacobian_                    tssetrhsjacobian
34e2df7a95SSatish Balay #define tsgetrhsjacobian_                    tsgetrhsjacobian
3537698f3aSJed Brown #define tssetifunction_                      tssetifunction
3637698f3aSJed Brown #define tsgetifunction_                      tsgetifunction
3737698f3aSJed Brown #define tssetijacobian_                      tssetijacobian
3837698f3aSJed Brown #define tsgetijacobian_                      tsgetijacobian
39e2df7a95SSatish Balay #define tsview_                              tsview
4037698f3aSJed Brown #define tssetoptionsprefix_                  tssetoptionsprefix
41e2df7a95SSatish Balay #define tsgetoptionsprefix_                  tsgetoptionsprefix
4237698f3aSJed Brown #define tsappendoptionsprefix_               tsappendoptionsprefix
43a6570f20SBarry Smith #define tsmonitorset_                        tsmonitorset
440e4ef248SJed Brown #define tscomputerhsfunctionlinear_          tscomputerhsfunctionlinear
450e4ef248SJed Brown #define tscomputerhsjacobianconstant_        tscomputerhsjacobianconstant
460fecffdcSJed Brown #define tscomputeifunctionlinear_            tscomputeifunctionlinear
470fecffdcSJed Brown #define tscomputeijacobianconstant_          tscomputeijacobianconstant
48a6570f20SBarry Smith #define tsmonitordefault_                    tsmonitordefault
49dd7ecb2fSBarry Smith #define tssetprestep_                        tssetprestep
50dd7ecb2fSBarry Smith #define tssetpoststep_                       tssetpoststep
51fe2efc57SMark #define tsviewfromoptions_                   tsviewfromoptions
52e2df7a95SSatish Balay #endif
53e2df7a95SSatish Balay 
54109c90ceSBarry Smith static struct {
55109c90ceSBarry Smith   PetscFortranCallbackId prestep;
56109c90ceSBarry Smith   PetscFortranCallbackId poststep;
57109c90ceSBarry Smith   PetscFortranCallbackId rhsfunction;
58109c90ceSBarry Smith   PetscFortranCallbackId rhsjacobian;
59109c90ceSBarry Smith   PetscFortranCallbackId ifunction;
60109c90ceSBarry Smith   PetscFortranCallbackId ijacobian;
61109c90ceSBarry Smith   PetscFortranCallbackId monitor;
62109c90ceSBarry Smith   PetscFortranCallbackId mondestroy;
63109c90ceSBarry Smith   PetscFortranCallbackId transform;
64109c90ceSBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG)
65109c90ceSBarry Smith   PetscFortranCallbackId function_pgiptr;
66109c90ceSBarry Smith #endif
67109c90ceSBarry Smith } _cb;
680fecffdcSJed Brown 
69dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts)
70dd7ecb2fSBarry Smith {
711ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
72109c90ceSBarry Smith   void *ptr;
73*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
74109c90ceSBarry Smith #endif
751ab477f3SBarry Smith   PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
76dd7ecb2fSBarry Smith }
77dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts)
78dd7ecb2fSBarry Smith {
791ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
80109c90ceSBarry Smith   void *ptr;
81*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
82109c90ceSBarry Smith #endif
831ab477f3SBarry Smith   PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */));
84dd7ecb2fSBarry Smith }
850fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, void *ctx)
86e2df7a95SSatish Balay {
871ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
88109c90ceSBarry Smith   void *ptr;
89*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
90109c90ceSBarry Smith #endif
911ab477f3SBarry 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) */));
92e2df7a95SSatish Balay }
930fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, void *ctx)
94e2df7a95SSatish Balay {
951ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
96109c90ceSBarry Smith   void *ptr;
97*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
98109c90ceSBarry Smith #endif
991ab477f3SBarry 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) */));
1000fecffdcSJed Brown }
101d1e9a80fSBarry Smith static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, void *ctx)
1020fecffdcSJed Brown {
1031ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
104109c90ceSBarry Smith   void *ptr;
105*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
106109c90ceSBarry Smith #endif
1071ab477f3SBarry 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) */));
1080fecffdcSJed Brown }
109d1e9a80fSBarry Smith static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, void *ctx)
1100fecffdcSJed Brown {
1111ab477f3SBarry Smith #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo)
112109c90ceSBarry Smith   void *ptr;
113*3ba16761SJacob Faibussowitsch   PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr));
114109c90ceSBarry Smith #endif
1151ab477f3SBarry 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) */));
116e2df7a95SSatish Balay }
117e2df7a95SSatish Balay 
118c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx)
119e2df7a95SSatish Balay {
120109c90ceSBarry Smith   TS ts = (TS)*ctx;
121109c90ceSBarry Smith   PetscObjectUseFortranCallback(ts,_cb.mondestroy,(void*,PetscErrorCode*),(_ctx,&ierr));
122e2df7a95SSatish Balay }
123e2df7a95SSatish Balay 
124e2df7a95SSatish Balay /*
125e2df7a95SSatish Balay    Note ctx is the same as ts so we need to get the Fortran context out of the TS
126e2df7a95SSatish Balay */
1270fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void *ctx)
128e2df7a95SSatish Balay {
129109c90ceSBarry Smith   PetscObjectUseFortranCallback(ts,_cb.monitor,(TS*,PetscInt*,PetscReal*,Vec *,void*,PetscErrorCode*),(&ts,&i,&d,&v,_ctx,&ierr));
130e2df7a95SSatish Balay }
131e2df7a95SSatish Balay 
132b0bc92c6SBarry Smith /*
133b0bc92c6SBarry Smith    Currently does not handle destroy or context
134b0bc92c6SBarry Smith */
135b0bc92c6SBarry Smith static PetscErrorCode ourtransform(void *ctx,Vec x,Vec *xout)
136b0bc92c6SBarry Smith {
137109c90ceSBarry Smith   PetscObjectUseFortranCallback((TS)ctx,_cb.transform,(void*,Vec *,Vec *,PetscErrorCode*),(_ctx,&x,xout,&ierr));
138b0bc92c6SBarry Smith }
139b0bc92c6SBarry Smith 
14019caf8f3SSatish Balay PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts,void (*f)(void*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode (*d)(void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
141b0bc92c6SBarry Smith {
142109c90ceSBarry Smith   *ierr = TSMonitorLGSetTransform(*ts,ourtransform,NULL,NULL); if (*ierr) return;
143109c90ceSBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.transform,(PetscVoidFunction)f,ctx);
144b0bc92c6SBarry Smith }
145b0bc92c6SBarry Smith 
14619caf8f3SSatish Balay PETSC_EXTERN void tssetprestep_(TS *ts,PetscErrorCode (*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
147dd7ecb2fSBarry Smith {
148109c90ceSBarry Smith   *ierr = TSSetPreStep(*ts,ourprestep);if (*ierr) return;
149109c90ceSBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.prestep,(PetscVoidFunction)f,NULL);
150dd7ecb2fSBarry Smith }
151dd7ecb2fSBarry Smith 
15219caf8f3SSatish Balay PETSC_EXTERN void tssetpoststep_(TS *ts,PetscErrorCode (*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
153dd7ecb2fSBarry Smith {
154109c90ceSBarry Smith   *ierr = TSSetPostStep(*ts,ourpoststep);if (*ierr) return;
155109c90ceSBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.poststep,(PetscVoidFunction)f,NULL);
156dd7ecb2fSBarry Smith }
157dd7ecb2fSBarry Smith 
15869c1e2abSSatish Balay PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr)
1590e4ef248SJed Brown {
1600e4ef248SJed Brown   *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx);
1610e4ef248SJed Brown }
16219caf8f3SSatish Balay PETSC_EXTERN void tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (*f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
163e2df7a95SSatish Balay {
1640e4ef248SJed Brown   Vec R;
1650e4ef248SJed Brown   CHKFORTRANNULLOBJECT(r);
1660e4ef248SJed Brown   CHKFORTRANNULLFUNCTION(f);
1670298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
1680e4ef248SJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) {
1690e4ef248SJed Brown     *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP);
1700e4ef248SJed Brown   } else {
171109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.rhsfunction,(PetscVoidFunction)f,fP);
172109c90ceSBarry Smith     *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,NULL);
1730fecffdcSJed Brown   }
1740fecffdcSJed Brown }
17519caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
17637698f3aSJed Brown {
17737698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
17837698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
1790298fd71SBarry Smith   *ierr = TSGetRHSFunction(*ts,r,NULL,ctx);
18037698f3aSJed Brown }
1810fecffdcSJed Brown 
18269c1e2abSSatish Balay PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr)
1830fecffdcSJed Brown {
1840fecffdcSJed Brown   *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx);
1850fecffdcSJed Brown }
18619caf8f3SSatish Balay PETSC_EXTERN void tssetifunction_(TS *ts,Vec *r,PetscErrorCode (*f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
1870fecffdcSJed Brown {
1880fecffdcSJed Brown   Vec R;
1890fecffdcSJed Brown   CHKFORTRANNULLOBJECT(r);
1900fecffdcSJed Brown   CHKFORTRANNULLFUNCTION(f);
1910298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
1920fecffdcSJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) {
1930fecffdcSJed Brown     *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP);
1940fecffdcSJed Brown   } else {
195109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ifunction,(PetscVoidFunction)f,fP);
196109c90ceSBarry Smith     *ierr = TSSetIFunction(*ts,R,ourifunction,NULL);
1970e4ef248SJed Brown   }
198e2df7a95SSatish Balay }
19919caf8f3SSatish Balay PETSC_EXTERN void tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
20037698f3aSJed Brown {
20137698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
20237698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
2030298fd71SBarry Smith   *ierr = TSGetIFunction(*ts,r,NULL,ctx);
20437698f3aSJed Brown }
20526d46c62SHong Zhang 
206e2df7a95SSatish Balay /* ---------------------------------------------------------*/
207d1e9a80fSBarry Smith PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr)
2080e4ef248SJed Brown {
209d1e9a80fSBarry Smith   *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,*A,*B,ctx);
2100e4ef248SJed Brown }
21119caf8f3SSatish Balay PETSC_EXTERN void tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
212e2df7a95SSatish Balay {
213aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(f);
214aecf964fSBarry Smith   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) {
2150e4ef248SJed Brown     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP);
216e2df7a95SSatish Balay   } else {
217109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.rhsjacobian,(PetscVoidFunction)f,fP);
218109c90ceSBarry Smith     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,NULL);
2190fecffdcSJed Brown   }
2200fecffdcSJed Brown }
2210fecffdcSJed Brown 
222d1e9a80fSBarry Smith PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr)
2230fecffdcSJed Brown {
224d1e9a80fSBarry Smith   *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,*A,*B,ctx);
2250fecffdcSJed Brown }
22619caf8f3SSatish Balay PETSC_EXTERN void tssetijacobian_(TS *ts,Mat *A,Mat *B,void (*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
2270fecffdcSJed Brown {
228aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(f);
229aecf964fSBarry Smith   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) {
2300fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP);
2310fecffdcSJed Brown   } else {
232109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ijacobian,(PetscVoidFunction)f,fP);
233109c90ceSBarry Smith     *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,NULL);
234e2df7a95SSatish Balay   }
235e2df7a95SSatish Balay }
23619caf8f3SSatish Balay PETSC_EXTERN void tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
23737698f3aSJed Brown {
23837698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
23937698f3aSJed Brown   CHKFORTRANNULLOBJECT(J);
24037698f3aSJed Brown   CHKFORTRANNULLOBJECT(M);
24137698f3aSJed Brown   *ierr = TSGetIJacobian(*ts,J,M,0,ctx);
24237698f3aSJed Brown }
243e2df7a95SSatish Balay 
24452f0073cSBarry Smith PETSC_EXTERN void tsmonitordefault_(TS *ts,PetscInt *its,PetscReal *fgnorm,Vec *u,PetscViewerAndFormat **dummy,PetscErrorCode *ierr)
24552f0073cSBarry Smith {
246410efd14SBarry Smith   *ierr = TSMonitorDefault(*ts,*its,*fgnorm,*u,*dummy);
24752f0073cSBarry Smith }
24852f0073cSBarry Smith 
249e2df7a95SSatish Balay /* ---------------------------------------------------------*/
250e2df7a95SSatish Balay 
25119caf8f3SSatish Balay /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */
252e2df7a95SSatish Balay 
25319caf8f3SSatish Balay PETSC_EXTERN void tsmonitorset_(TS *ts,void (*func)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*),void *mctx,void (*d)(void*,PetscErrorCode*),PetscErrorCode *ierr)
254e2df7a95SSatish Balay {
255aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(d);
25652f0073cSBarry Smith   if ((PetscVoidFunction)func == (PetscVoidFunction) tsmonitordefault_) {
2571cb03803SBarry Smith     *ierr = TSMonitorSet(*ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,void*))TSMonitorDefault,*(PetscViewerAndFormat**)mctx,(PetscErrorCode (*)(void **))PetscViewerAndFormatDestroy);
258e2df7a95SSatish Balay   } else {
259109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.monitor,(PetscVoidFunction)func,mctx);
260109c90ceSBarry Smith     *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.mondestroy,(PetscVoidFunction)d,mctx);
261aecf964fSBarry Smith     *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy);
262e2df7a95SSatish Balay   }
263e2df7a95SSatish Balay }
264e2df7a95SSatish Balay 
265e2df7a95SSatish Balay /* ---------------------------------------------------------*/
266089b2837SJed Brown /*  func is currently ignored from Fortran */
26719caf8f3SSatish Balay PETSC_EXTERN void tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
268e2df7a95SSatish Balay {
269089b2837SJed Brown   *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx);
270e2df7a95SSatish Balay }
271e2df7a95SSatish Balay 
27219caf8f3SSatish Balay PETSC_EXTERN void tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
273e2df7a95SSatish Balay {
274e2df7a95SSatish Balay   PetscViewer v;
275e2df7a95SSatish Balay   PetscPatchDefaultViewers_Fortran(viewer,v);
276e2df7a95SSatish Balay   *ierr = TSView(*ts,v);
277e2df7a95SSatish Balay }
278e2df7a95SSatish Balay 
27919caf8f3SSatish Balay PETSC_EXTERN void tssetoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
28037698f3aSJed Brown {
28137698f3aSJed Brown   char *t;
28237698f3aSJed Brown   FIXCHAR(prefix,len,t);
283d49bb8f9SBarry Smith   *ierr = TSSetOptionsPrefix(*ts,t);if (*ierr) return;
28437698f3aSJed Brown   FREECHAR(prefix,t);
28537698f3aSJed Brown }
28619caf8f3SSatish Balay PETSC_EXTERN void tsgetoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
287e2df7a95SSatish Balay {
288e2df7a95SSatish Balay   const char *tname;
289e2df7a95SSatish Balay 
290e2df7a95SSatish Balay   *ierr = TSGetOptionsPrefix(*ts,&tname);
291e2df7a95SSatish Balay   *ierr = PetscStrncpy(prefix,tname,len);
292d6a8cea5SBarry Smith   FIXRETURNCHAR(PETSC_TRUE,prefix,len);
293e2df7a95SSatish Balay }
29419caf8f3SSatish Balay PETSC_EXTERN void tsappendoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
29537698f3aSJed Brown {
29637698f3aSJed Brown   char *t;
29737698f3aSJed Brown   FIXCHAR(prefix,len,t);
298d49bb8f9SBarry Smith   *ierr = TSAppendOptionsPrefix(*ts,t);if (*ierr) return;
29937698f3aSJed Brown   FREECHAR(prefix,t);
30037698f3aSJed Brown }
301e2df7a95SSatish Balay 
30219caf8f3SSatish Balay PETSC_EXTERN void tsviewfromoptions_(TS *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
303fe2efc57SMark {
304fe2efc57SMark   char *t;
305fe2efc57SMark 
306fe2efc57SMark   FIXCHAR(type,len,t);
307b14c0cbaSBlaise Bourdin   CHKFORTRANNULLOBJECT(obj);
308fe2efc57SMark   *ierr = TSViewFromOptions(*ao,obj,t);if (*ierr) return;
309fe2efc57SMark   FREECHAR(type,t);
310fe2efc57SMark }
311