xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 37698f3a7f72c8273c9c38edbd354583caf5230e)
1b45d2f2cSJed Brown #include <petsc-private/fortranimpl.h>
2c6db04a5SJed Brown #include <petscts.h>
3e2df7a95SSatish Balay 
4e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
5e2df7a95SSatish Balay #define tssetrhsfunction_                    TSSETRHSFUNCTION
6*37698f3aSJed Brown #define tsgetrhsfunction_                    TSGETRHSFUNCTION
7e2df7a95SSatish Balay #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
8e2df7a95SSatish Balay #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
9*37698f3aSJed Brown #define tssetifunction_                      TSSETIFUNCTION
10*37698f3aSJed Brown #define tsgetifunction_                      TSGETIFUNCTION
11*37698f3aSJed Brown #define tssetijacobian_                      TSSETIJACOBIAN
12*37698f3aSJed Brown #define tsgetijacobian_                      TSGETIJACOBIAN
13e2df7a95SSatish Balay #define tsview_                              TSVIEW
14*37698f3aSJed Brown #define tssetoptionsprefix_                  TSSETOPTIONSPREFIX
15e2df7a95SSatish Balay #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
16*37698f3aSJed Brown #define tsappendoptionsprefix_               TSAPPENDOPTIONSPREFIX
17a6570f20SBarry Smith #define tsmonitorset_                        TSMONITORSET
180e4ef248SJed Brown #define tscomputerhsfunctionlinear_          TSCOMPUTERHSFUNCTIONLINEAR
190e4ef248SJed Brown #define tscomputerhsjacobianconstant_        TSCOMPUTERHSJACOBIANCONSTANT
200fecffdcSJed Brown #define tscomputeifunctionlinear_            TSCOMPUTEIFUNCTIONLINEAR
210fecffdcSJed Brown #define tscomputeijacobianconstant_          TSCOMPUTEIJACOBIANCONSTANT
22a6570f20SBarry Smith #define tsmonitordefault_                    TSMONITORDEFAULT
23dd7ecb2fSBarry Smith #define tssetprestep_                        TSSETPRESTEP
24dd7ecb2fSBarry Smith #define tssetpoststep_                       TSSETPOSTSTEP
25e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
26e2df7a95SSatish Balay #define tssetrhsfunction_                    tssetrhsfunction
27*37698f3aSJed Brown #define tsgetrhsfunction_                    tsgetrhsfunction
28e2df7a95SSatish Balay #define tssetrhsjacobian_                    tssetrhsjacobian
29e2df7a95SSatish Balay #define tsgetrhsjacobian_                    tsgetrhsjacobian
30*37698f3aSJed Brown #define tssetifunction_                      tssetifunction
31*37698f3aSJed Brown #define tsgetifunction_                      tsgetifunction
32*37698f3aSJed Brown #define tssetijacobian_                      tssetijacobian
33*37698f3aSJed Brown #define tsgetijacobian_                      tsgetijacobian
34e2df7a95SSatish Balay #define tsview_                              tsview
35*37698f3aSJed Brown #define tssetoptionsprefix_                  tssetoptionsprefix
36e2df7a95SSatish Balay #define tsgetoptionsprefix_                  tsgetoptionsprefix
37*37698f3aSJed Brown #define tsappendoptionsprefix_               tsappendoptionsprefix
38a6570f20SBarry Smith #define tsmonitorset_                        tsmonitorset
390e4ef248SJed Brown #define tscomputerhsfunctionlinear_          tscomputerhsfunctionlinear
400e4ef248SJed Brown #define tscomputerhsjacobianconstant_        tscomputerhsjacobianconstant
410fecffdcSJed Brown #define tscomputeifunctionlinear_            tscomputeifunctionlinear
420fecffdcSJed Brown #define tscomputeijacobianconstant_          tscomputeijacobianconstant
43a6570f20SBarry Smith #define tsmonitordefault_                    tsmonitordefault
44dd7ecb2fSBarry Smith #define tssetprestep_                        tssetprestep
45dd7ecb2fSBarry Smith #define tssetpoststep_                       tssetpoststep
46e2df7a95SSatish Balay #endif
47e2df7a95SSatish Balay 
480fecffdcSJed Brown enum {OUR_PRESTEP = 0,
490fecffdcSJed Brown       OUR_POSTSTEP,
500fecffdcSJed Brown       OUR_RHSFUNCTION,
510fecffdcSJed Brown       OUR_IFUNCTION,
520fecffdcSJed Brown       OUR_RHSJACOBIAN,
530fecffdcSJed Brown       OUR_IJACOBIAN,
540fecffdcSJed Brown       OUR_MONITOR,
550fecffdcSJed Brown       OUR_MONITORDESTROY,
560fecffdcSJed Brown       OUR_MONITOR_CTX,   /* Casting from function pointer is invalid according to the standard. */
570fecffdcSJed Brown       OUR_COUNT};
580fecffdcSJed Brown 
59dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts)
60dd7ecb2fSBarry Smith {
61dd7ecb2fSBarry Smith   PetscErrorCode ierr = 0;
620fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_PRESTEP]))(&ts,&ierr);
63dd7ecb2fSBarry Smith   return 0;
64dd7ecb2fSBarry Smith }
65dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts)
66dd7ecb2fSBarry Smith {
67dd7ecb2fSBarry Smith   PetscErrorCode ierr = 0;
680fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_POSTSTEP]))(&ts,&ierr);
69dd7ecb2fSBarry Smith   return 0;
70dd7ecb2fSBarry Smith }
710fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
72e2df7a95SSatish Balay {
73e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
740fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSFUNCTION]))(&ts,&d,&x,&f,ctx,&ierr);
75e2df7a95SSatish Balay   return 0;
76e2df7a95SSatish Balay }
770fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,void *ctx)
78e2df7a95SSatish Balay {
79e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
800fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IFUNCTION]))(&ts,&d,&x,&xdot,&f,ctx,&ierr);
810fecffdcSJed Brown   return 0;
820fecffdcSJed Brown }
830fecffdcSJed Brown static PetscErrorCode ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx)
840fecffdcSJed Brown {
850fecffdcSJed Brown   PetscErrorCode ierr = 0;
860fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSJACOBIAN]))(&ts,&d,&x,m,p,type,ctx,&ierr);
870fecffdcSJed Brown   return 0;
880fecffdcSJed Brown }
890fecffdcSJed Brown static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat* m,Mat* p,MatStructure* type,void*ctx)
900fecffdcSJed Brown {
910fecffdcSJed Brown   PetscErrorCode ierr = 0;
920fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IJACOBIAN]))(&ts,&d,&x,&xdot,&shift,m,p,type,ctx,&ierr);
93e2df7a95SSatish Balay   return 0;
94e2df7a95SSatish Balay }
95e2df7a95SSatish Balay 
96c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx)
97e2df7a95SSatish Balay {
98e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
99c2efdce3SBarry Smith   TS          ts = *(TS*)ctx;
1000fecffdcSJed Brown   void        *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
1010fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITORDESTROY]))(mctx,&ierr);
102e2df7a95SSatish Balay   return 0;
103e2df7a95SSatish Balay }
104e2df7a95SSatish Balay 
105e2df7a95SSatish Balay /*
106e2df7a95SSatish Balay    Note ctx is the same as ts so we need to get the Fortran context out of the TS
107e2df7a95SSatish Balay */
1080fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx)
109e2df7a95SSatish Balay {
110e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
1110fecffdcSJed Brown   void           *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
1120fecffdcSJed Brown   (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR]))(&ts,&i,&d,&v,mctx,&ierr);
113e2df7a95SSatish Balay   return 0;
114e2df7a95SSatish Balay }
115e2df7a95SSatish Balay 
116e2df7a95SSatish Balay EXTERN_C_BEGIN
117e2df7a95SSatish Balay 
118dd7ecb2fSBarry Smith void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
119dd7ecb2fSBarry Smith {
1200fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1210fecffdcSJed Brown   ((PetscObject)*ts)->fortran_func_pointers[OUR_PRESTEP] = (PetscVoidFunction)f;
122dd7ecb2fSBarry Smith   *ierr = TSSetPreStep(*ts,ourprestep);
123dd7ecb2fSBarry Smith }
124dd7ecb2fSBarry Smith 
125dd7ecb2fSBarry Smith void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
126dd7ecb2fSBarry Smith {
1270fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1280fecffdcSJed Brown   ((PetscObject)*ts)->fortran_func_pointers[OUR_POSTSTEP] = (PetscVoidFunction)f;
129dd7ecb2fSBarry Smith   *ierr = TSSetPreStep(*ts,ourpoststep);
130dd7ecb2fSBarry Smith }
131dd7ecb2fSBarry Smith 
1320e4ef248SJed Brown void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr)
1330e4ef248SJed Brown {
1340e4ef248SJed Brown   *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx);
1350e4ef248SJed Brown }
136089b2837SJed Brown void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
137e2df7a95SSatish Balay {
1380e4ef248SJed Brown   Vec R;
1390e4ef248SJed Brown   CHKFORTRANNULLOBJECT(r);
1400e4ef248SJed Brown   CHKFORTRANNULLFUNCTION(f);
1410e4ef248SJed Brown   CHKFORTRANNULLOBJECT(fP);
1421db2f0e5SJed Brown   R = r ? *r : (Vec)PETSC_NULL;
1430e4ef248SJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) {
1440e4ef248SJed Brown     *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP);
1450e4ef248SJed Brown   } else {
1460fecffdcSJed Brown     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1470fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSFUNCTION] = (PetscVoidFunction)f;
1480fecffdcSJed Brown     *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,fP);
1490fecffdcSJed Brown   }
1500fecffdcSJed Brown }
151*37698f3aSJed Brown void PETSC_STDCALL tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
152*37698f3aSJed Brown {
153*37698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
154*37698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
155*37698f3aSJed Brown   *ierr = TSGetRHSFunction(*ts,r,PETSC_NULL,ctx);
156*37698f3aSJed Brown }
1570fecffdcSJed Brown 
1580fecffdcSJed Brown void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr)
1590fecffdcSJed Brown {
1600fecffdcSJed Brown   *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx);
1610fecffdcSJed Brown }
1620fecffdcSJed Brown void PETSC_STDCALL tssetifunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
1630fecffdcSJed Brown {
1640fecffdcSJed Brown   Vec R;
1650fecffdcSJed Brown   CHKFORTRANNULLOBJECT(r);
1660fecffdcSJed Brown   CHKFORTRANNULLFUNCTION(f);
1670fecffdcSJed Brown   CHKFORTRANNULLOBJECT(fP);
1681db2f0e5SJed Brown   R = r ? *r : (Vec)PETSC_NULL;
1690fecffdcSJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) {
1700fecffdcSJed Brown     *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP);
1710fecffdcSJed Brown   } else {
1720fecffdcSJed Brown     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1730fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_IFUNCTION] = (PetscVoidFunction)f;
1740fecffdcSJed Brown     *ierr = TSSetIFunction(*ts,R,ourifunction,fP);
1750e4ef248SJed Brown   }
176e2df7a95SSatish Balay }
177*37698f3aSJed Brown void PETSC_STDCALL tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
178*37698f3aSJed Brown {
179*37698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
180*37698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
181*37698f3aSJed Brown   *ierr = TSGetIFunction(*ts,r,PETSC_NULL,ctx);
182*37698f3aSJed Brown }
18326d46c62SHong Zhang 
184e2df7a95SSatish Balay /* ---------------------------------------------------------*/
1850e4ef248SJed Brown void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
1860e4ef248SJed Brown {
1870e4ef248SJed Brown   *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,A,B,flg,ctx);
1880e4ef248SJed Brown }
189e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,
190e2df7a95SSatish Balay                void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
191e2df7a95SSatish Balay {
1920fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
193e2df7a95SSatish Balay   if (FORTRANNULLFUNCTION(f)) {
194e2df7a95SSatish Balay     *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP);
1950e4ef248SJed Brown   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) {
1960e4ef248SJed Brown     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP);
197e2df7a95SSatish Balay   } else {
1980fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSJACOBIAN] = (PetscVoidFunction)f;
1990fecffdcSJed Brown     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,fP);
2000fecffdcSJed Brown   }
2010fecffdcSJed Brown }
2020fecffdcSJed Brown 
2030fecffdcSJed Brown void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
2040fecffdcSJed Brown {
2050fecffdcSJed Brown   *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,A,B,flg,ctx);
2060fecffdcSJed Brown }
2070fecffdcSJed Brown void PETSC_STDCALL tssetijacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,
2080fecffdcSJed Brown                void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
2090fecffdcSJed Brown {
2100fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
2110fecffdcSJed Brown   if (FORTRANNULLFUNCTION(f)) {
2120fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,PETSC_NULL,fP);
2130fecffdcSJed Brown   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) {
2140fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP);
2150fecffdcSJed Brown   } else {
2160fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_IJACOBIAN] = (PetscVoidFunction)f;
2170fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,fP);
218e2df7a95SSatish Balay   }
219e2df7a95SSatish Balay }
220*37698f3aSJed Brown void PETSC_STDCALL tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
221*37698f3aSJed Brown {
222*37698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
223*37698f3aSJed Brown   CHKFORTRANNULLOBJECT(J);
224*37698f3aSJed Brown   CHKFORTRANNULLOBJECT(M);
225*37698f3aSJed Brown   *ierr = TSGetIJacobian(*ts,J,M,0,ctx);
226*37698f3aSJed Brown }
227e2df7a95SSatish Balay 
228e2df7a95SSatish Balay /* ---------------------------------------------------------*/
229e2df7a95SSatish Balay 
230a6570f20SBarry Smith extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
231e2df7a95SSatish Balay 
232a6570f20SBarry Smith void PETSC_STDCALL tsmonitorset_(TS *ts,void (PETSC_STDCALL *func)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*),void (*mctx)(void),void (PETSC_STDCALL *d)(void*,PetscErrorCode*),PetscErrorCode *ierr)
233e2df7a95SSatish Balay {
2340fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
235a6570f20SBarry Smith   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
236a6570f20SBarry Smith     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
237e2df7a95SSatish Balay   } else {
2380fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR]        = (PetscVoidFunction)func;
2390fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITORDESTROY] = (PetscVoidFunction)d;
2400fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR_CTX]    = (PetscVoidFunction)mctx;
241e2df7a95SSatish Balay     if (FORTRANNULLFUNCTION(d)) {
2420fecffdcSJed Brown       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,0);
243e2df7a95SSatish Balay     } else {
2440fecffdcSJed Brown       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy);
245e2df7a95SSatish Balay     }
246e2df7a95SSatish Balay   }
247e2df7a95SSatish Balay }
248e2df7a95SSatish Balay 
249e2df7a95SSatish Balay /* ---------------------------------------------------------*/
250089b2837SJed Brown /*  func is currently ignored from Fortran */
251089b2837SJed Brown void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
252e2df7a95SSatish Balay {
253089b2837SJed Brown   *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx);
254e2df7a95SSatish Balay }
255e2df7a95SSatish Balay 
256e2df7a95SSatish Balay void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
257e2df7a95SSatish Balay {
258e2df7a95SSatish Balay   PetscViewer v;
259e2df7a95SSatish Balay   PetscPatchDefaultViewers_Fortran(viewer,v);
260e2df7a95SSatish Balay   *ierr = TSView(*ts,v);
261e2df7a95SSatish Balay }
262e2df7a95SSatish Balay 
263*37698f3aSJed Brown void PETSC_STDCALL tssetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
264*37698f3aSJed Brown {
265*37698f3aSJed Brown   char *t;
266*37698f3aSJed Brown   FIXCHAR(prefix,len,t);
267*37698f3aSJed Brown   *ierr = TSSetOptionsPrefix(*ts,t);
268*37698f3aSJed Brown   FREECHAR(prefix,t);
269*37698f3aSJed Brown }
270e2df7a95SSatish Balay void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
271e2df7a95SSatish Balay {
272e2df7a95SSatish Balay   const char *tname;
273e2df7a95SSatish Balay 
274e2df7a95SSatish Balay   *ierr = TSGetOptionsPrefix(*ts,&tname);
275e2df7a95SSatish Balay   *ierr = PetscStrncpy(prefix,tname,len);
276e2df7a95SSatish Balay }
277*37698f3aSJed Brown void PETSC_STDCALL tsappendoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
278*37698f3aSJed Brown {
279*37698f3aSJed Brown   char *t;
280*37698f3aSJed Brown   FIXCHAR(prefix,len,t);
281*37698f3aSJed Brown   *ierr = TSAppendOptionsPrefix(*ts,t);
282*37698f3aSJed Brown   FREECHAR(prefix,t);
283*37698f3aSJed Brown }
284e2df7a95SSatish Balay 
285e2df7a95SSatish Balay 
286e2df7a95SSatish Balay EXTERN_C_END
287