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