xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 665c2ded495bb9782a7454dcfef3abf1536c3670)
1b45d2f2cSJed Brown #include <petsc-private/fortranimpl.h>
2c6db04a5SJed Brown #include <petscts.h>
3*665c2dedSJed Brown #include <petscviewer.h>
4e2df7a95SSatish Balay 
5e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
6e2df7a95SSatish Balay #define tssetrhsfunction_                    TSSETRHSFUNCTION
737698f3aSJed Brown #define tsgetrhsfunction_                    TSGETRHSFUNCTION
8e2df7a95SSatish Balay #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
9e2df7a95SSatish Balay #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
1037698f3aSJed Brown #define tssetifunction_                      TSSETIFUNCTION
1137698f3aSJed Brown #define tsgetifunction_                      TSGETIFUNCTION
1237698f3aSJed Brown #define tssetijacobian_                      TSSETIJACOBIAN
1337698f3aSJed Brown #define tsgetijacobian_                      TSGETIJACOBIAN
14e2df7a95SSatish Balay #define tsview_                              TSVIEW
1537698f3aSJed Brown #define tssetoptionsprefix_                  TSSETOPTIONSPREFIX
16e2df7a95SSatish Balay #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
1737698f3aSJed Brown #define tsappendoptionsprefix_               TSAPPENDOPTIONSPREFIX
18a6570f20SBarry Smith #define tsmonitorset_                        TSMONITORSET
190e4ef248SJed Brown #define tscomputerhsfunctionlinear_          TSCOMPUTERHSFUNCTIONLINEAR
200e4ef248SJed Brown #define tscomputerhsjacobianconstant_        TSCOMPUTERHSJACOBIANCONSTANT
210fecffdcSJed Brown #define tscomputeifunctionlinear_            TSCOMPUTEIFUNCTIONLINEAR
220fecffdcSJed Brown #define tscomputeijacobianconstant_          TSCOMPUTEIJACOBIANCONSTANT
23a6570f20SBarry Smith #define tsmonitordefault_                    TSMONITORDEFAULT
24dd7ecb2fSBarry Smith #define tssetprestep_                        TSSETPRESTEP
25dd7ecb2fSBarry Smith #define tssetpoststep_                       TSSETPOSTSTEP
26e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
27e2df7a95SSatish Balay #define tssetrhsfunction_                    tssetrhsfunction
2837698f3aSJed Brown #define tsgetrhsfunction_                    tsgetrhsfunction
29e2df7a95SSatish Balay #define tssetrhsjacobian_                    tssetrhsjacobian
30e2df7a95SSatish Balay #define tsgetrhsjacobian_                    tsgetrhsjacobian
3137698f3aSJed Brown #define tssetifunction_                      tssetifunction
3237698f3aSJed Brown #define tsgetifunction_                      tsgetifunction
3337698f3aSJed Brown #define tssetijacobian_                      tssetijacobian
3437698f3aSJed Brown #define tsgetijacobian_                      tsgetijacobian
35e2df7a95SSatish Balay #define tsview_                              tsview
3637698f3aSJed Brown #define tssetoptionsprefix_                  tssetoptionsprefix
37e2df7a95SSatish Balay #define tsgetoptionsprefix_                  tsgetoptionsprefix
3837698f3aSJed Brown #define tsappendoptionsprefix_               tsappendoptionsprefix
39a6570f20SBarry Smith #define tsmonitorset_                        tsmonitorset
400e4ef248SJed Brown #define tscomputerhsfunctionlinear_          tscomputerhsfunctionlinear
410e4ef248SJed Brown #define tscomputerhsjacobianconstant_        tscomputerhsjacobianconstant
420fecffdcSJed Brown #define tscomputeifunctionlinear_            tscomputeifunctionlinear
430fecffdcSJed Brown #define tscomputeijacobianconstant_          tscomputeijacobianconstant
44a6570f20SBarry Smith #define tsmonitordefault_                    tsmonitordefault
45dd7ecb2fSBarry Smith #define tssetprestep_                        tssetprestep
46dd7ecb2fSBarry Smith #define tssetpoststep_                       tssetpoststep
47e2df7a95SSatish Balay #endif
48e2df7a95SSatish Balay 
490fecffdcSJed Brown enum {OUR_PRESTEP = 0,
500fecffdcSJed Brown       OUR_POSTSTEP,
510fecffdcSJed Brown       OUR_RHSFUNCTION,
520fecffdcSJed Brown       OUR_IFUNCTION,
530fecffdcSJed Brown       OUR_RHSJACOBIAN,
540fecffdcSJed Brown       OUR_IJACOBIAN,
550fecffdcSJed Brown       OUR_MONITOR,
560fecffdcSJed Brown       OUR_MONITORDESTROY,
570fecffdcSJed Brown       OUR_MONITOR_CTX,   /* Casting from function pointer is invalid according to the standard. */
580fecffdcSJed Brown       OUR_COUNT};
590fecffdcSJed Brown 
60dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts)
61dd7ecb2fSBarry Smith {
62dd7ecb2fSBarry Smith   PetscErrorCode ierr = 0;
630fecffdcSJed Brown   (*(void (PETSC_STDCALL*)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_PRESTEP]))(&ts,&ierr);
64dd7ecb2fSBarry Smith   return 0;
65dd7ecb2fSBarry Smith }
66dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts)
67dd7ecb2fSBarry Smith {
68dd7ecb2fSBarry Smith   PetscErrorCode ierr = 0;
690fecffdcSJed Brown   (*(void (PETSC_STDCALL*)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_POSTSTEP]))(&ts,&ierr);
70dd7ecb2fSBarry Smith   return 0;
71dd7ecb2fSBarry Smith }
720fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
73e2df7a95SSatish Balay {
74e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
750fecffdcSJed Brown   (*(void (PETSC_STDCALL*)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSFUNCTION]))(&ts,&d,&x,&f,ctx,&ierr);
76e2df7a95SSatish Balay   return 0;
77e2df7a95SSatish Balay }
780fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,void *ctx)
79e2df7a95SSatish Balay {
80e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
810fecffdcSJed 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);
820fecffdcSJed Brown   return 0;
830fecffdcSJed Brown }
840fecffdcSJed Brown static PetscErrorCode ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat *m,Mat *p,MatStructure *type,void *ctx)
850fecffdcSJed Brown {
860fecffdcSJed Brown   PetscErrorCode ierr = 0;
870fecffdcSJed 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);
880fecffdcSJed Brown   return 0;
890fecffdcSJed Brown }
900fecffdcSJed Brown static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat *m,Mat *p,MatStructure *type,void *ctx)
910fecffdcSJed Brown {
920fecffdcSJed Brown   PetscErrorCode ierr = 0;
930fecffdcSJed 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);
94e2df7a95SSatish Balay   return 0;
95e2df7a95SSatish Balay }
96e2df7a95SSatish Balay 
97c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx)
98e2df7a95SSatish Balay {
99e2df7a95SSatish Balay   PetscErrorCode ierr  = 0;
100c2efdce3SBarry Smith   TS             ts    = *(TS*)ctx;
1010fecffdcSJed Brown   void           *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
1020fecffdcSJed Brown   (*(void (PETSC_STDCALL*)(void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITORDESTROY]))(mctx,&ierr);
103e2df7a95SSatish Balay   return 0;
104e2df7a95SSatish Balay }
105e2df7a95SSatish Balay 
106e2df7a95SSatish Balay /*
107e2df7a95SSatish Balay    Note ctx is the same as ts so we need to get the Fortran context out of the TS
108e2df7a95SSatish Balay */
1090fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void *ctx)
110e2df7a95SSatish Balay {
111e2df7a95SSatish Balay   PetscErrorCode ierr  = 0;
1120fecffdcSJed Brown   void           *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
1130fecffdcSJed Brown   (*(void (PETSC_STDCALL*)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR]))(&ts,&i,&d,&v,mctx,&ierr);
114e2df7a95SSatish Balay   return 0;
115e2df7a95SSatish Balay }
116e2df7a95SSatish Balay 
117e2df7a95SSatish Balay EXTERN_C_BEGIN
118e2df7a95SSatish Balay 
119dd7ecb2fSBarry Smith void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
120dd7ecb2fSBarry Smith {
1210fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1220fecffdcSJed Brown   ((PetscObject)*ts)->fortran_func_pointers[OUR_PRESTEP] = (PetscVoidFunction)f;
123bbd56ea5SKarl Rupp 
124dd7ecb2fSBarry Smith   *ierr = TSSetPreStep(*ts,ourprestep);
125dd7ecb2fSBarry Smith }
126dd7ecb2fSBarry Smith 
127dd7ecb2fSBarry Smith void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
128dd7ecb2fSBarry Smith {
1290fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1300fecffdcSJed Brown   ((PetscObject)*ts)->fortran_func_pointers[OUR_POSTSTEP] = (PetscVoidFunction)f;
131bbd56ea5SKarl Rupp 
132dd7ecb2fSBarry Smith   *ierr = TSSetPreStep(*ts,ourpoststep);
133dd7ecb2fSBarry Smith }
134dd7ecb2fSBarry Smith 
1350e4ef248SJed Brown void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr)
1360e4ef248SJed Brown {
1370e4ef248SJed Brown   *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx);
1380e4ef248SJed Brown }
139089b2837SJed Brown void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
140e2df7a95SSatish Balay {
1410e4ef248SJed Brown   Vec R;
1420e4ef248SJed Brown   CHKFORTRANNULLOBJECT(r);
1430e4ef248SJed Brown   CHKFORTRANNULLFUNCTION(f);
1440e4ef248SJed Brown   CHKFORTRANNULLOBJECT(fP);
1450298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
1460e4ef248SJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) {
1470e4ef248SJed Brown     *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP);
1480e4ef248SJed Brown   } else {
1490fecffdcSJed Brown     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1500fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSFUNCTION] = (PetscVoidFunction)f;
1510fecffdcSJed Brown     *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,fP);
1520fecffdcSJed Brown   }
1530fecffdcSJed Brown }
15437698f3aSJed Brown void PETSC_STDCALL tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
15537698f3aSJed Brown {
15637698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
15737698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
1580298fd71SBarry Smith   *ierr = TSGetRHSFunction(*ts,r,NULL,ctx);
15937698f3aSJed Brown }
1600fecffdcSJed Brown 
1610fecffdcSJed Brown void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr)
1620fecffdcSJed Brown {
1630fecffdcSJed Brown   *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx);
1640fecffdcSJed Brown }
1650fecffdcSJed Brown void PETSC_STDCALL tssetifunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
1660fecffdcSJed Brown {
1670fecffdcSJed Brown   Vec R;
1680fecffdcSJed Brown   CHKFORTRANNULLOBJECT(r);
1690fecffdcSJed Brown   CHKFORTRANNULLFUNCTION(f);
1700fecffdcSJed Brown   CHKFORTRANNULLOBJECT(fP);
1710298fd71SBarry Smith   R = r ? *r : (Vec)NULL;
1720fecffdcSJed Brown   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) {
1730fecffdcSJed Brown     *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP);
1740fecffdcSJed Brown   } else {
1750fecffdcSJed Brown     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
1760fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_IFUNCTION] = (PetscVoidFunction)f;
1770fecffdcSJed Brown     *ierr = TSSetIFunction(*ts,R,ourifunction,fP);
1780e4ef248SJed Brown   }
179e2df7a95SSatish Balay }
18037698f3aSJed Brown void PETSC_STDCALL tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
18137698f3aSJed Brown {
18237698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
18337698f3aSJed Brown   CHKFORTRANNULLOBJECT(r);
1840298fd71SBarry Smith   *ierr = TSGetIFunction(*ts,r,NULL,ctx);
18537698f3aSJed Brown }
18626d46c62SHong Zhang 
187e2df7a95SSatish Balay /* ---------------------------------------------------------*/
1880e4ef248SJed Brown void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
1890e4ef248SJed Brown {
1900e4ef248SJed Brown   *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,A,B,flg,ctx);
1910e4ef248SJed Brown }
192bbd56ea5SKarl Rupp void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
193e2df7a95SSatish Balay {
1940fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
195e2df7a95SSatish Balay   if (FORTRANNULLFUNCTION(f)) {
1960298fd71SBarry Smith     *ierr = TSSetRHSJacobian(*ts,*A,*B,NULL,fP);
1970e4ef248SJed Brown   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) {
1980e4ef248SJed Brown     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP);
199e2df7a95SSatish Balay   } else {
2000fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSJACOBIAN] = (PetscVoidFunction)f;
2010fecffdcSJed Brown     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,fP);
2020fecffdcSJed Brown   }
2030fecffdcSJed Brown }
2040fecffdcSJed Brown 
2050fecffdcSJed Brown void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
2060fecffdcSJed Brown {
2070fecffdcSJed Brown   *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,A,B,flg,ctx);
2080fecffdcSJed Brown }
209bbd56ea5SKarl Rupp void PETSC_STDCALL tssetijacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
2100fecffdcSJed Brown {
2110fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
2120fecffdcSJed Brown   if (FORTRANNULLFUNCTION(f)) {
2130298fd71SBarry Smith     *ierr = TSSetIJacobian(*ts,*A,*B,NULL,fP);
2140fecffdcSJed Brown   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) {
2150fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP);
2160fecffdcSJed Brown   } else {
2170fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_IJACOBIAN] = (PetscVoidFunction)f;
2180fecffdcSJed Brown     *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,fP);
219e2df7a95SSatish Balay   }
220e2df7a95SSatish Balay }
22137698f3aSJed Brown void PETSC_STDCALL tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
22237698f3aSJed Brown {
22337698f3aSJed Brown   CHKFORTRANNULLINTEGER(ctx);
22437698f3aSJed Brown   CHKFORTRANNULLOBJECT(J);
22537698f3aSJed Brown   CHKFORTRANNULLOBJECT(M);
22637698f3aSJed Brown   *ierr = TSGetIJacobian(*ts,J,M,0,ctx);
22737698f3aSJed Brown }
228e2df7a95SSatish Balay 
229e2df7a95SSatish Balay /* ---------------------------------------------------------*/
230e2df7a95SSatish Balay 
231a6570f20SBarry Smith extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
232e2df7a95SSatish Balay 
233a6570f20SBarry 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)
234e2df7a95SSatish Balay {
2350fecffdcSJed Brown   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
236a6570f20SBarry Smith   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
237a6570f20SBarry Smith     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
238e2df7a95SSatish Balay   } else {
2390fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR]        = (PetscVoidFunction)func;
2400fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITORDESTROY] = (PetscVoidFunction)d;
2410fecffdcSJed Brown     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR_CTX]    = (PetscVoidFunction)mctx;
242e2df7a95SSatish Balay     if (FORTRANNULLFUNCTION(d)) {
2430fecffdcSJed Brown       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,0);
244e2df7a95SSatish Balay     } else {
2450fecffdcSJed Brown       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy);
246e2df7a95SSatish Balay     }
247e2df7a95SSatish Balay   }
248e2df7a95SSatish Balay }
249e2df7a95SSatish Balay 
250e2df7a95SSatish Balay /* ---------------------------------------------------------*/
251089b2837SJed Brown /*  func is currently ignored from Fortran */
252089b2837SJed Brown void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
253e2df7a95SSatish Balay {
254089b2837SJed Brown   *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx);
255e2df7a95SSatish Balay }
256e2df7a95SSatish Balay 
257e2df7a95SSatish Balay void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
258e2df7a95SSatish Balay {
259e2df7a95SSatish Balay   PetscViewer v;
260e2df7a95SSatish Balay   PetscPatchDefaultViewers_Fortran(viewer,v);
261e2df7a95SSatish Balay   *ierr = TSView(*ts,v);
262e2df7a95SSatish Balay }
263e2df7a95SSatish Balay 
26437698f3aSJed Brown void PETSC_STDCALL tssetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
26537698f3aSJed Brown {
26637698f3aSJed Brown   char *t;
26737698f3aSJed Brown   FIXCHAR(prefix,len,t);
26837698f3aSJed Brown   *ierr = TSSetOptionsPrefix(*ts,t);
26937698f3aSJed Brown   FREECHAR(prefix,t);
27037698f3aSJed Brown }
271e2df7a95SSatish Balay void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
272e2df7a95SSatish Balay {
273e2df7a95SSatish Balay   const char *tname;
274e2df7a95SSatish Balay 
275e2df7a95SSatish Balay   *ierr = TSGetOptionsPrefix(*ts,&tname);
276e2df7a95SSatish Balay   *ierr = PetscStrncpy(prefix,tname,len);
277e2df7a95SSatish Balay }
27837698f3aSJed Brown void PETSC_STDCALL tsappendoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
27937698f3aSJed Brown {
28037698f3aSJed Brown   char *t;
28137698f3aSJed Brown   FIXCHAR(prefix,len,t);
28237698f3aSJed Brown   *ierr = TSAppendOptionsPrefix(*ts,t);
28337698f3aSJed Brown   FREECHAR(prefix,t);
28437698f3aSJed Brown }
285e2df7a95SSatish Balay 
286e2df7a95SSatish Balay 
287e2df7a95SSatish Balay EXTERN_C_END
288