xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision fbe7908bcbdc335e0dab871e0aba2fea2d749384)
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,MatStructure *type,void *ctx)
85 {
86   PetscErrorCode ierr = 0;
87   (*(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);
88   return 0;
89 }
90 static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat *m,Mat *p,MatStructure *type,void *ctx)
91 {
92   PetscErrorCode ierr = 0;
93   (*(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);
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*) ((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*) ((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 EXTERN_C_BEGIN
118 
119 void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
120 {
121   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
122   ((PetscObject)*ts)->fortran_func_pointers[OUR_PRESTEP] = (PetscVoidFunction)f;
123 
124   *ierr = TSSetPreStep(*ts,ourprestep);
125 }
126 
127 void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
128 {
129   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
130   ((PetscObject)*ts)->fortran_func_pointers[OUR_POSTSTEP] = (PetscVoidFunction)f;
131 
132   *ierr = TSSetPreStep(*ts,ourpoststep);
133 }
134 
135 void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr)
136 {
137   *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx);
138 }
139 void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
140 {
141   Vec R;
142   CHKFORTRANNULLOBJECT(r);
143   CHKFORTRANNULLFUNCTION(f);
144   CHKFORTRANNULLOBJECT(fP);
145   R = r ? *r : (Vec)NULL;
146   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) {
147     *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP);
148   } else {
149     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
150     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSFUNCTION] = (PetscVoidFunction)f;
151     *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,fP);
152   }
153 }
154 void PETSC_STDCALL tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
155 {
156   CHKFORTRANNULLINTEGER(ctx);
157   CHKFORTRANNULLOBJECT(r);
158   *ierr = TSGetRHSFunction(*ts,r,NULL,ctx);
159 }
160 
161 void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr)
162 {
163   *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx);
164 }
165 void PETSC_STDCALL tssetifunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL*f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr)
166 {
167   Vec R;
168   CHKFORTRANNULLOBJECT(r);
169   CHKFORTRANNULLFUNCTION(f);
170   CHKFORTRANNULLOBJECT(fP);
171   R = r ? *r : (Vec)NULL;
172   if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) {
173     *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP);
174   } else {
175     PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
176     ((PetscObject)*ts)->fortran_func_pointers[OUR_IFUNCTION] = (PetscVoidFunction)f;
177     *ierr = TSSetIFunction(*ts,R,ourifunction,fP);
178   }
179 }
180 void PETSC_STDCALL tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
181 {
182   CHKFORTRANNULLINTEGER(ctx);
183   CHKFORTRANNULLOBJECT(r);
184   *ierr = TSGetIFunction(*ts,r,NULL,ctx);
185 }
186 
187 /* ---------------------------------------------------------*/
188 void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
189 {
190   *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,A,B,flg,ctx);
191 }
192 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)
193 {
194   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
195   if (FORTRANNULLFUNCTION(f)) {
196     *ierr = TSSetRHSJacobian(*ts,*A,*B,NULL,fP);
197   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) {
198     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP);
199   } else {
200     ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSJACOBIAN] = (PetscVoidFunction)f;
201     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,fP);
202   }
203 }
204 
205 void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr)
206 {
207   *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,A,B,flg,ctx);
208 }
209 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)
210 {
211   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
212   if (FORTRANNULLFUNCTION(f)) {
213     *ierr = TSSetIJacobian(*ts,*A,*B,NULL,fP);
214   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) {
215     *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP);
216   } else {
217     ((PetscObject)*ts)->fortran_func_pointers[OUR_IJACOBIAN] = (PetscVoidFunction)f;
218     *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,fP);
219   }
220 }
221 void PETSC_STDCALL tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
222 {
223   CHKFORTRANNULLINTEGER(ctx);
224   CHKFORTRANNULLOBJECT(J);
225   CHKFORTRANNULLOBJECT(M);
226   *ierr = TSGetIJacobian(*ts,J,M,0,ctx);
227 }
228 
229 /* ---------------------------------------------------------*/
230 
231 extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
232 
233 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)
234 {
235   PetscObjectAllocateFortranPointers(*ts,OUR_COUNT);
236   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
237     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
238   } else {
239     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR]        = (PetscVoidFunction)func;
240     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITORDESTROY] = (PetscVoidFunction)d;
241     ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR_CTX]    = (PetscVoidFunction)mctx;
242     if (FORTRANNULLFUNCTION(d)) {
243       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,0);
244     } else {
245       *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy);
246     }
247   }
248 }
249 
250 /* ---------------------------------------------------------*/
251 /*  func is currently ignored from Fortran */
252 void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr)
253 {
254   *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx);
255 }
256 
257 void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
258 {
259   PetscViewer v;
260   PetscPatchDefaultViewers_Fortran(viewer,v);
261   *ierr = TSView(*ts,v);
262 }
263 
264 void PETSC_STDCALL tssetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
265 {
266   char *t;
267   FIXCHAR(prefix,len,t);
268   *ierr = TSSetOptionsPrefix(*ts,t);
269   FREECHAR(prefix,t);
270 }
271 void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
272 {
273   const char *tname;
274 
275   *ierr = TSGetOptionsPrefix(*ts,&tname);
276   *ierr = PetscStrncpy(prefix,tname,len);
277 }
278 void PETSC_STDCALL tsappendoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
279 {
280   char *t;
281   FIXCHAR(prefix,len,t);
282   *ierr = TSAppendOptionsPrefix(*ts,t);
283   FREECHAR(prefix,t);
284 }
285 
286 
287 EXTERN_C_END
288