xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision b52f573bf7dd30aec2f45ea71ed0f119e1fe824d)
1 #include "zpetsc.h"
2 #include "petscts.h"
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define tssetrhsfunction_                    TSSETRHSFUNCTION
6 #define tssetmatrices_                       TSSETMATRICES
7 #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
8 #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
9 #define tsgetrhsmatrix_                      TSGETRHSMATRIX
10 #define tsview_                              TSVIEW
11 #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
12 #define tsmonitorset_                        TSMONITORSET
13 #define tsdefaultcomputejacobian_            TSDEFAULTCOMPUTEJACOBIAN
14 #define tsdefaultcomputejacobiancolor_       TSDEFAULTCOMPUTEJACOBIANCOLOR
15 #define tsmonitordefault_                    TSMONITORDEFAULT
16 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
17 #define tssetrhsfunction_                    tssetrhsfunction
18 #define tssetmatrices_                       tssetmatrices
19 #define tssetrhsjacobian_                    tssetrhsjacobian
20 #define tsgetrhsjacobian_                    tsgetrhsjacobian
21 #define tsgetrhsmatrix_                      tsgetrhsmatrix
22 #define tsview_                              tsview
23 #define tsgetoptionsprefix_                  tsgetoptionsprefix
24 #define tsmonitorset_                        tsmonitorset
25 #define tsdefaultcomputejacobian_            tsdefaultcomputejacobian
26 #define tsdefaultcomputejacobiancolor_       tsdefaultcomputejacobiancolor
27 #define tsmonitordefault_                    tsmonitordefault
28 #endif
29 
30 static PetscErrorCode ourtsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
31 {
32   PetscErrorCode ierr = 0;
33   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[1]))(&ts,&d,&x,&f,ctx,&ierr);
34   return 0;
35 }
36 static PetscErrorCode ourtsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
37 {
38   PetscErrorCode ierr = 0;
39   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[2]))(&ts,&d,m,p,type,ctx,&ierr);
40   return 0;
41 }
42 static PetscErrorCode ourtslhsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
43 {
44   PetscErrorCode ierr = 0;
45   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[7]))(&ts,&d,m,p,type,ctx,&ierr);
46   return 0;
47 }
48 static PetscErrorCode ourtsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx)
49 {
50   PetscErrorCode ierr = 0;
51   (*(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);
52   return 0;
53 }
54 
55 static PetscErrorCode ourtsdestroy(void *ctx)
56 {
57   PetscErrorCode ierr = 0;
58   TS          ts = (TS)ctx;
59   void        (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6];
60   (*(void (PETSC_STDCALL *)(PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[5]))(mctx,&ierr);
61   return 0;
62 }
63 
64 /*
65    Note ctx is the same as ts so we need to get the Fortran context out of the TS
66 */
67 static PetscErrorCode ourtsmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx)
68 {
69   PetscErrorCode ierr = 0;
70   void       (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6];
71   (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[4]))(&ts,&i,&d,&v,mctx,&ierr);
72   return 0;
73 }
74 
75 EXTERN_C_BEGIN
76 
77 
78 void PETSC_STDCALL tssetrhsfunction_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
79 {
80   ((PetscObject)*ts)->fortran_func_pointers[1] = (PetscVoidFunction)f;
81   *ierr = TSSetRHSFunction(*ts,ourtsfunction,fP);
82 }
83 #ifdef MV
84 void PETSC_STDCALL tssetrhsmatrix_(TS *ts,Mat *A,Mat *B,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
85                                                    void*,PetscInt *),void*fP,PetscErrorCode *ierr)
86 {
87   if (FORTRANNULLFUNCTION(f)) {
88     *ierr = TSSetRHSMatrix(*ts,*A,*B,PETSC_NULL,fP);
89   } else {
90     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)f;
91     *ierr = TSSetRHSMatrix(*ts,*A,*B,ourtsmatrix,fP);
92   }
93 }
94 #endif
95 void PETSC_STDCALL tssetmatrices_(TS *ts,Mat *Arhs,PetscErrorCode (PETSC_STDCALL *frhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
96                                                    void*,PetscInt *),
97                                          Mat *Alhs,PetscErrorCode (PETSC_STDCALL *flhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
98                                                    void*,PetscInt *),
99                                          MatStructure *flag,void*fP,PetscErrorCode *ierr)
100 {
101   if (FORTRANNULLFUNCTION(frhs) && FORTRANNULLFUNCTION(flhs)) {
102     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,PETSC_NULL,*flag,fP);
103   } else if (FORTRANNULLFUNCTION(flhs)){
104     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
105     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,PETSC_NULL,*flag,fP);
106   } else if (FORTRANNULLFUNCTION(frhs)){
107     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
108     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,ourtslhsmatrix,*flag,fP);
109   } else {
110     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
111     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
112     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,ourtslhsmatrix,*flag,fP);
113   }
114 }
115 
116 /* ---------------------------------------------------------*/
117 extern void tsdefaultcomputejacobian_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
118 extern void tsdefaultcomputejacobiancolor_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
119 
120 void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,
121                void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
122 {
123   if (FORTRANNULLFUNCTION(f)) {
124     *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP);
125   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobian_) {
126     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobian,fP);
127   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobiancolor_) {
128     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobianColor,*(MatFDColoring*)fP);
129   } else {
130   ((PetscObject)*ts)->fortran_func_pointers[3] = (PetscVoidFunction)f;
131     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourtsjacobian,fP);
132   }
133 }
134 
135 /* ---------------------------------------------------------*/
136 
137 extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
138 
139 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)
140 {
141   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
142     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
143   } else {
144     ((PetscObject)*ts)->fortran_func_pointers[4] = (PetscVoidFunction)func;
145     ((PetscObject)*ts)->fortran_func_pointers[5] = (PetscVoidFunction)d;
146     ((PetscObject)*ts)->fortran_func_pointers[6] = (PetscVoidFunction)mctx;
147     if (FORTRANNULLFUNCTION(d)) {
148       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,0);
149     } else {
150       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,ourtsdestroy);
151     }
152   }
153 }
154 
155 /* ---------------------------------------------------------*/
156 void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr)
157 {
158   *ierr = TSGetRHSJacobian(*ts,J,M,ctx);
159 }
160 
161 void PETSC_STDCALL tsgetrhsmatrix_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr)
162 {
163   *ierr = TSGetRHSMatrix(*ts,J,M,ctx);
164 }
165 
166 void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
167 {
168   PetscViewer v;
169   PetscPatchDefaultViewers_Fortran(viewer,v);
170   *ierr = TSView(*ts,v);
171 }
172 
173 void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
174 {
175   const char *tname;
176 
177   *ierr = TSGetOptionsPrefix(*ts,&tname);
178 #if defined(PETSC_USES_CPTOFCD)
179   {
180     char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix);
181     *ierr = PetscStrncpy(t,tname,len1);
182   }
183 #else
184   *ierr = PetscStrncpy(prefix,tname,len);
185 #endif
186 }
187 
188 
189 EXTERN_C_END
190