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