xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision 1ebf93c6b7d760d592de6ebe6cdc0debaa3caf75)
1 #include <petsc/private/fortranimpl.h>
2 #include <petscts.h>
3 #include <petscviewer.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6 #define tsmonitorlgsettransform_             TSMONITORLGSETTRANSFORM
7 #define tssetrhsfunction_                    TSSETRHSFUNCTION
8 #define tsgetrhsfunction_                    TSGETRHSFUNCTION
9 #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
10 #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
11 #define tssetifunction_                      TSSETIFUNCTION
12 #define tsgetifunction_                      TSGETIFUNCTION
13 #define tssetijacobian_                      TSSETIJACOBIAN
14 #define tsgetijacobian_                      TSGETIJACOBIAN
15 #define tsview_                              TSVIEW
16 #define tssetoptionsprefix_                  TSSETOPTIONSPREFIX
17 #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
18 #define tsappendoptionsprefix_               TSAPPENDOPTIONSPREFIX
19 #define tsmonitorset_                        TSMONITORSET
20 #define tscomputerhsfunctionlinear_          TSCOMPUTERHSFUNCTIONLINEAR
21 #define tscomputerhsjacobianconstant_        TSCOMPUTERHSJACOBIANCONSTANT
22 #define tscomputeifunctionlinear_            TSCOMPUTEIFUNCTIONLINEAR
23 #define tscomputeijacobianconstant_          TSCOMPUTEIJACOBIANCONSTANT
24 #define tsmonitordefault_                    TSMONITORDEFAULT
25 #define tssetprestep_                        TSSETPRESTEP
26 #define tssetpoststep_                       TSSETPOSTSTEP
27 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
28 #define tsmonitorlgsettransform_             tsmonitorlgsettransform
29 #define tssetrhsfunction_                    tssetrhsfunction
30 #define tsgetrhsfunction_                    tsgetrhsfunction
31 #define tssetrhsjacobian_                    tssetrhsjacobian
32 #define tsgetrhsjacobian_                    tsgetrhsjacobian
33 #define tssetifunction_                      tssetifunction
34 #define tsgetifunction_                      tsgetifunction
35 #define tssetijacobian_                      tssetijacobian
36 #define tsgetijacobian_                      tsgetijacobian
37 #define tsview_                              tsview
38 #define tssetoptionsprefix_                  tssetoptionsprefix
39 #define tsgetoptionsprefix_                  tsgetoptionsprefix
40 #define tsappendoptionsprefix_               tsappendoptionsprefix
41 #define tsmonitorset_                        tsmonitorset
42 #define tscomputerhsfunctionlinear_          tscomputerhsfunctionlinear
43 #define tscomputerhsjacobianconstant_        tscomputerhsjacobianconstant
44 #define tscomputeifunctionlinear_            tscomputeifunctionlinear
45 #define tscomputeijacobianconstant_          tscomputeijacobianconstant
46 #define tsmonitordefault_                    tsmonitordefault
47 #define tssetprestep_                        tssetprestep
48 #define tssetpoststep_                       tssetpoststep
49 #endif
50 
51 enum {OUR_PRESTEP = 0,
52       OUR_POSTSTEP,
53       OUR_RHSFUNCTION,
54       OUR_IFUNCTION,
55       OUR_RHSJACOBIAN,
56       OUR_IJACOBIAN,
57       OUR_MONITOR,
58       OUR_MONITORDESTROY,
59       OUR_MONITOR_CTX,   /* Casting from function pointer is invalid according to the standard. */
60       OUR_COUNT};
61 
62 static PetscErrorCode ourprestep(TS ts)
63 {
64   PetscErrorCode ierr = 0;
65   (*(void (PETSC_STDCALL*)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_PRESTEP]))(&ts,&ierr);
66   return 0;
67 }
68 static PetscErrorCode ourpoststep(TS ts)
69 {
70   PetscErrorCode ierr = 0;
71   (*(void (PETSC_STDCALL*)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_POSTSTEP]))(&ts,&ierr);
72   return 0;
73 }
74 static PetscErrorCode ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
75 {
76   PetscErrorCode ierr = 0;
77   (*(void (PETSC_STDCALL*)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSFUNCTION]))(&ts,&d,&x,&f,ctx,&ierr);
78   return 0;
79 }
80 static PetscErrorCode ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,void *ctx)
81 {
82   PetscErrorCode ierr = 0;
83   (*(void (PETSC_STDCALL*)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IFUNCTION]))(&ts,&d,&x,&xdot,&f,ctx,&ierr);
84   return 0;
85 }
86 static PetscErrorCode ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat m,Mat p,void *ctx)
87 {
88   PetscErrorCode ierr = 0;
89   (*(void (PETSC_STDCALL*)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSJACOBIAN]))(&ts,&d,&x,&m,&p,ctx,&ierr);
90   return 0;
91 }
92 static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,Mat p,void *ctx)
93 {
94   PetscErrorCode ierr = 0;
95   (*(void (PETSC_STDCALL*)(TS*,PetscReal*,Vec*,Vec*,PetscReal*,Mat*,Mat*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IJACOBIAN]))(&ts,&d,&x,&xdot,&shift,&m,&p,ctx,&ierr);
96   return 0;
97 }
98 
99 static PetscErrorCode ourmonitordestroy(void **ctx)
100 {
101   PetscErrorCode ierr  = 0;
102   TS             ts    = *(TS*)ctx;
103   void           *mctx = (void*) (PETSC_UINTPTR_T) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
104   (*(void (PETSC_STDCALL*)(void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITORDESTROY]))(mctx,&ierr);
105   return 0;
106 }
107 
108 /*
109    Note ctx is the same as ts so we need to get the Fortran context out of the TS
110 */
111 static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void *ctx)
112 {
113   PetscErrorCode ierr  = 0;
114   void           *mctx = (void*) (PETSC_UINTPTR_T) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX];
115   (*(void (PETSC_STDCALL*)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR]))(&ts,&i,&d,&v,mctx,&ierr);
116   return 0;
117 }
118 
119 /*
120    Currently does not handle destroy or context
121 */
122 static void (PETSC_STDCALL*yourtransform)(void*ctx,Vec*,Vec*,PetscErrorCode*);
123 static PetscErrorCode ourtransform(void *ctx,Vec x,Vec *xout)
124 {
125   PetscErrorCode ierr = 0;
126   (*yourtransform)(ctx,&x,xout,&ierr);
127   return ierr;
128 }
129 
130 PETSC_EXTERN void PETSC_STDCALL tsmonitorlgsettransform_(TS *ts,void (PETSC_STDCALL*f)(void*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode (PETSC_STDCALL*d)(void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
131 {
132   *ierr = TSMonitorLGSetTransform(*ts,ourtransform,NULL,NULL);
133   yourtransform = f;
134 }
135 
136 PETSC_EXTERN void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
137 {
138   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
139   ((PetscObject)*ts)->fortran_func_pointers[OUR_PRESTEP] = (PetscVoidFunction)f;
140 
141   *ierr = TSSetPreStep(*ts,ourprestep);
142 }
143 
144 PETSC_EXTERN void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
145 {
146   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
147   ((PetscObject)*ts)->fortran_func_pointers[OUR_POSTSTEP] = (PetscVoidFunction)f;
148 
149   *ierr = TSSetPostStep(*ts,ourpoststep);
150 }
151 
152 PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr)
153 {
154   *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx);
155 }
156 PETSC_EXTERN void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
157 {
158   Vec R;
159   CHKFORTRANNULLOBJECT(r);
160   CHKFORTRANNULLFUNCTION(f);
161   CHKFORTRANNULLOBJECT(fP);
162   R = r ? *r : (Vec)NULL;
163   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) {
164     *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP);
165   } else {
166     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
167     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSFUNCTION] = (PetscVoidFunction)f;
168     *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,fP);
169   }
170 }
171 PETSC_EXTERN void PETSC_STDCALL tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
172 {
173   CHKFORTRANNULLINTEGER(ctx);
174   CHKFORTRANNULLOBJECT(r);
175   *ierr = TSGetRHSFunction(*ts,r,NULL,ctx);
176 }
177 
178 PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr)
179 {
180   *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx);
181 }
182 PETSC_EXTERN void PETSC_STDCALL tssetifunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
183 {
184   Vec R;
185   CHKFORTRANNULLOBJECT(r);
186   CHKFORTRANNULLFUNCTION(f);
187   CHKFORTRANNULLOBJECT(fP);
188   R = r ? *r : (Vec)NULL;
189   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) {
190     *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP);
191   } else {
192     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
193     ((PetscObject)*ts)->fortran_func_pointers[OUR_IFUNCTION] = (PetscVoidFunction)f;
194     *ierr = TSSetIFunction(*ts,R,ourifunction,fP);
195   }
196 }
197 PETSC_EXTERN void PETSC_STDCALL tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
198 {
199   CHKFORTRANNULLINTEGER(ctx);
200   CHKFORTRANNULLOBJECT(r);
201   *ierr = TSGetIFunction(*ts,r,NULL,ctx);
202 }
203 
204 /* ---------------------------------------------------------*/
205 PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr)
206 {
207   *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,*A,*B,ctx);
208 }
209 PETSC_EXTERN void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
210 {
211   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
212   if (FORTRANNULLFUNCTION(f)) {
213     *ierr = TSSetRHSJacobian(*ts,*A,*B,NULL,fP);
214   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) {
215     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP);
216   } else {
217     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSJACOBIAN] = (PetscVoidFunction)f;
218     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,fP);
219   }
220 }
221 
222 PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr)
223 {
224   *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,*A,*B,ctx);
225 }
226 PETSC_EXTERN void PETSC_STDCALL tssetijacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
227 {
228   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
229   if (FORTRANNULLFUNCTION(f)) {
230     *ierr = TSSetIJacobian(*ts,*A,*B,NULL,fP);
231   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) {
232     *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP);
233   } else {
234     ((PetscObject)*ts)->fortran_func_pointers[OUR_IJACOBIAN] = (PetscVoidFunction)f;
235     *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,fP);
236   }
237 }
238 PETSC_EXTERN void PETSC_STDCALL tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
239 {
240   CHKFORTRANNULLINTEGER(ctx);
241   CHKFORTRANNULLOBJECT(J);
242   CHKFORTRANNULLOBJECT(M);
243   *ierr = TSGetIJacobian(*ts,J,M,0,ctx);
244 }
245 
246 /* ---------------------------------------------------------*/
247 
248 PETSC_EXTERN void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
249 
250 PETSC_EXTERN 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)
251 {
252   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
253   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
254     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
255   } else {
256     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR]        = (PetscVoidFunction)func;
257     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITORDESTROY] = (PetscVoidFunction)d;
258     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR_CTX]    = (PetscVoidFunction)mctx;
259     if (FORTRANNULLFUNCTION(d)) {
260       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,0);
261     } else {
262       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy);
263     }
264   }
265 }
266 
267 /* ---------------------------------------------------------*/
268 /*  func is currently ignored from Fortran */
269 PETSC_EXTERN void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
270 {
271   *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx);
272 }
273 
274 PETSC_EXTERN void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
275 {
276   PetscViewer v;
277   PetscPatchDefaultViewers_Fortran(viewer,v);
278   *ierr = TSView(*ts,v);
279 }
280 
281 PETSC_EXTERN void PETSC_STDCALL tssetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
282 {
283   char *t;
284   FIXCHAR(prefix,len,t);
285   *ierr = TSSetOptionsPrefix(*ts,t);
286   FREECHAR(prefix,t);
287 }
288 PETSC_EXTERN void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
289 {
290   const char *tname;
291 
292   *ierr = TSGetOptionsPrefix(*ts,&tname);
293   *ierr = PetscStrncpy(prefix,tname,len);
294 }
295 PETSC_EXTERN void PETSC_STDCALL tsappendoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
296 {
297   char *t;
298   FIXCHAR(prefix,len,t);
299   *ierr = TSAppendOptionsPrefix(*ts,t);
300   FREECHAR(prefix,t);
301 }
302 
303