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