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