xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 26d46c62b9c7bef6a87dcb65c60e507b04e8055b)
1e2df7a95SSatish Balay #include "zpetsc.h"
2e2df7a95SSatish Balay #include "petscts.h"
3e2df7a95SSatish Balay 
4e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
5e2df7a95SSatish Balay #define tssetrhsfunction_                    TSSETRHSFUNCTION
676f2fa84SHong Zhang #define tssetmatrices_                       TSSETMATRICES
7*26d46c62SHong Zhang #define tsgetmatrices_                       TSGETMATRICES
8e2df7a95SSatish Balay #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
9e2df7a95SSatish Balay #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
10e2df7a95SSatish Balay #define tsview_                              TSVIEW
11e2df7a95SSatish Balay #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
12a6570f20SBarry Smith #define tsmonitorset_                        TSMONITORSET
13e2df7a95SSatish Balay #define tsdefaultcomputejacobian_            TSDEFAULTCOMPUTEJACOBIAN
14e2df7a95SSatish Balay #define tsdefaultcomputejacobiancolor_       TSDEFAULTCOMPUTEJACOBIANCOLOR
15a6570f20SBarry Smith #define tsmonitordefault_                    TSMONITORDEFAULT
16e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
17e2df7a95SSatish Balay #define tssetrhsfunction_                    tssetrhsfunction
1876f2fa84SHong Zhang #define tssetmatrices_                       tssetmatrices
19*26d46c62SHong Zhang #define tsgetmatrices_                       tsgetmatrices
20e2df7a95SSatish Balay #define tssetrhsjacobian_                    tssetrhsjacobian
21e2df7a95SSatish Balay #define tsgetrhsjacobian_                    tsgetrhsjacobian
22e2df7a95SSatish Balay #define tsview_                              tsview
23e2df7a95SSatish Balay #define tsgetoptionsprefix_                  tsgetoptionsprefix
24a6570f20SBarry Smith #define tsmonitorset_                        tsmonitorset
25e2df7a95SSatish Balay #define tsdefaultcomputejacobian_            tsdefaultcomputejacobian
26e2df7a95SSatish Balay #define tsdefaultcomputejacobiancolor_       tsdefaultcomputejacobiancolor
27a6570f20SBarry Smith #define tsmonitordefault_                    tsmonitordefault
28e2df7a95SSatish Balay #endif
29e2df7a95SSatish Balay 
30e2df7a95SSatish Balay static PetscErrorCode ourtsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
31e2df7a95SSatish Balay {
32e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
33e2df7a95SSatish Balay   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[1]))(&ts,&d,&x,&f,ctx,&ierr);
34e2df7a95SSatish Balay   return 0;
35e2df7a95SSatish Balay }
36e2df7a95SSatish Balay static PetscErrorCode ourtsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
37e2df7a95SSatish Balay {
38e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
39e2df7a95SSatish Balay   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[2]))(&ts,&d,m,p,type,ctx,&ierr);
40e2df7a95SSatish Balay   return 0;
41e2df7a95SSatish Balay }
4276f2fa84SHong Zhang static PetscErrorCode ourtslhsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
4376f2fa84SHong Zhang {
4476f2fa84SHong Zhang   PetscErrorCode ierr = 0;
4576f2fa84SHong Zhang   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[7]))(&ts,&d,m,p,type,ctx,&ierr);
4676f2fa84SHong Zhang   return 0;
4776f2fa84SHong Zhang }
48e2df7a95SSatish Balay static PetscErrorCode ourtsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx)
49e2df7a95SSatish Balay {
50e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
51e2df7a95SSatish Balay   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[3]))(&ts,&d,&x,m,p,type,ctx,&ierr);
52e2df7a95SSatish Balay   return 0;
53e2df7a95SSatish Balay }
54e2df7a95SSatish Balay 
55e2df7a95SSatish Balay static PetscErrorCode ourtsdestroy(void *ctx)
56e2df7a95SSatish Balay {
57e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
58e2df7a95SSatish Balay   TS          ts = (TS)ctx;
59e2df7a95SSatish Balay   void        (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6];
60f68b968cSBarry Smith   (*(void (PETSC_STDCALL *)(PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[5]))(mctx,&ierr);
61e2df7a95SSatish Balay   return 0;
62e2df7a95SSatish Balay }
63e2df7a95SSatish Balay 
64e2df7a95SSatish Balay /*
65e2df7a95SSatish Balay    Note ctx is the same as ts so we need to get the Fortran context out of the TS
66e2df7a95SSatish Balay */
67e2df7a95SSatish Balay static PetscErrorCode ourtsmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx)
68e2df7a95SSatish Balay {
69e2df7a95SSatish Balay   PetscErrorCode ierr = 0;
70e2df7a95SSatish Balay   void       (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6];
71f68b968cSBarry Smith   (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[4]))(&ts,&i,&d,&v,mctx,&ierr);
72e2df7a95SSatish Balay   return 0;
73e2df7a95SSatish Balay }
74e2df7a95SSatish Balay 
75e2df7a95SSatish Balay EXTERN_C_BEGIN
76e2df7a95SSatish Balay 
77e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsfunction_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
78e2df7a95SSatish Balay {
79f68b968cSBarry Smith   ((PetscObject)*ts)->fortran_func_pointers[1] = (PetscVoidFunction)f;
80e2df7a95SSatish Balay   *ierr = TSSetRHSFunction(*ts,ourtsfunction,fP);
81e2df7a95SSatish Balay }
82*26d46c62SHong Zhang 
8376f2fa84SHong Zhang void PETSC_STDCALL tssetmatrices_(TS *ts,Mat *Arhs,PetscErrorCode (PETSC_STDCALL *frhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
8476f2fa84SHong Zhang                                                    void*,PetscInt *),
8576f2fa84SHong Zhang                                          Mat *Alhs,PetscErrorCode (PETSC_STDCALL *flhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
8676f2fa84SHong Zhang                                                    void*,PetscInt *),
8776f2fa84SHong Zhang                                          MatStructure *flag,void*fP,PetscErrorCode *ierr)
8876f2fa84SHong Zhang {
8976f2fa84SHong Zhang   if (FORTRANNULLFUNCTION(frhs) && FORTRANNULLFUNCTION(flhs)) {
9076f2fa84SHong Zhang     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,PETSC_NULL,*flag,fP);
9176f2fa84SHong Zhang   } else if (FORTRANNULLFUNCTION(flhs)){
9276f2fa84SHong Zhang     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
9376f2fa84SHong Zhang     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,PETSC_NULL,*flag,fP);
9476f2fa84SHong Zhang   } else if (FORTRANNULLFUNCTION(frhs)){
9576f2fa84SHong Zhang     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
9676f2fa84SHong Zhang     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,ourtslhsmatrix,*flag,fP);
9776f2fa84SHong Zhang   } else {
9876f2fa84SHong Zhang     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
9976f2fa84SHong Zhang     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
10076f2fa84SHong Zhang     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,ourtslhsmatrix,*flag,fP);
10176f2fa84SHong Zhang   }
10276f2fa84SHong Zhang }
103e2df7a95SSatish Balay 
104e2df7a95SSatish Balay /* ---------------------------------------------------------*/
105e2df7a95SSatish Balay extern void tsdefaultcomputejacobian_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
106e2df7a95SSatish Balay extern void tsdefaultcomputejacobiancolor_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
107e2df7a95SSatish Balay 
108e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,
109e2df7a95SSatish Balay                void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
110e2df7a95SSatish Balay {
111e2df7a95SSatish Balay   if (FORTRANNULLFUNCTION(f)) {
112e2df7a95SSatish Balay     *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP);
113f68b968cSBarry Smith   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobian_) {
114e2df7a95SSatish Balay     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobian,fP);
115f68b968cSBarry Smith   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobiancolor_) {
116e2df7a95SSatish Balay     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobianColor,*(MatFDColoring*)fP);
117e2df7a95SSatish Balay   } else {
118f68b968cSBarry Smith   ((PetscObject)*ts)->fortran_func_pointers[3] = (PetscVoidFunction)f;
119e2df7a95SSatish Balay     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourtsjacobian,fP);
120e2df7a95SSatish Balay   }
121e2df7a95SSatish Balay }
122e2df7a95SSatish Balay 
123e2df7a95SSatish Balay /* ---------------------------------------------------------*/
124e2df7a95SSatish Balay 
125a6570f20SBarry Smith extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
126e2df7a95SSatish Balay 
127a6570f20SBarry 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)
128e2df7a95SSatish Balay {
129a6570f20SBarry Smith   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
130a6570f20SBarry Smith     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
131e2df7a95SSatish Balay   } else {
132f68b968cSBarry Smith     ((PetscObject)*ts)->fortran_func_pointers[4] = (PetscVoidFunction)func;
133f68b968cSBarry Smith     ((PetscObject)*ts)->fortran_func_pointers[5] = (PetscVoidFunction)d;
134f68b968cSBarry Smith     ((PetscObject)*ts)->fortran_func_pointers[6] = (PetscVoidFunction)mctx;
135e2df7a95SSatish Balay     if (FORTRANNULLFUNCTION(d)) {
136a6570f20SBarry Smith       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,0);
137e2df7a95SSatish Balay     } else {
138a6570f20SBarry Smith       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,ourtsdestroy);
139e2df7a95SSatish Balay     }
140e2df7a95SSatish Balay   }
141e2df7a95SSatish Balay }
142e2df7a95SSatish Balay 
143e2df7a95SSatish Balay /* ---------------------------------------------------------*/
144e2df7a95SSatish Balay void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr)
145e2df7a95SSatish Balay {
146e2df7a95SSatish Balay   *ierr = TSGetRHSJacobian(*ts,J,M,ctx);
147e2df7a95SSatish Balay }
148e2df7a95SSatish Balay 
149*26d46c62SHong Zhang void PETSC_STDCALL tsgetmatrices_(TS *ts,Mat *Arhs,Mat *Alhs,void **ctx,PetscErrorCode *ierr)
150e2df7a95SSatish Balay {
151*26d46c62SHong Zhang   *ierr = TSGetMatrices(*ts,Arhs,Alhs,ctx);
152e2df7a95SSatish Balay }
153e2df7a95SSatish Balay 
154e2df7a95SSatish Balay void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
155e2df7a95SSatish Balay {
156e2df7a95SSatish Balay   PetscViewer v;
157e2df7a95SSatish Balay   PetscPatchDefaultViewers_Fortran(viewer,v);
158e2df7a95SSatish Balay   *ierr = TSView(*ts,v);
159e2df7a95SSatish Balay }
160e2df7a95SSatish Balay 
161e2df7a95SSatish Balay void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
162e2df7a95SSatish Balay {
163e2df7a95SSatish Balay   const char *tname;
164e2df7a95SSatish Balay 
165e2df7a95SSatish Balay   *ierr = TSGetOptionsPrefix(*ts,&tname);
166e2df7a95SSatish Balay #if defined(PETSC_USES_CPTOFCD)
167e2df7a95SSatish Balay   {
168e2df7a95SSatish Balay     char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix);
169e2df7a95SSatish Balay     *ierr = PetscStrncpy(t,tname,len1);
170e2df7a95SSatish Balay   }
171e2df7a95SSatish Balay #else
172e2df7a95SSatish Balay   *ierr = PetscStrncpy(prefix,tname,len);
173e2df7a95SSatish Balay #endif
174e2df7a95SSatish Balay }
175e2df7a95SSatish Balay 
176e2df7a95SSatish Balay 
177e2df7a95SSatish Balay EXTERN_C_END
178