xref: /petsc/src/ts/interface/ftn-custom/ztsf.c (revision dd7ecb2f948ceb99ce111cd0ccf9220ab4ed17cb)
1 #include <private/fortranimpl.h>
2 #include <petscts.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define tssetrhsfunction_                    TSSETRHSFUNCTION
6 #define tssetmatrices_                       TSSETMATRICES
7 #define tsgetmatrices_                       TSGETMATRICES
8 #define tssetrhsjacobian_                    TSSETRHSJACOBIAN
9 #define tsgetrhsjacobian_                    TSGETRHSJACOBIAN
10 #define tsview_                              TSVIEW
11 #define tsgetoptionsprefix_                  TSGETOPTIONSPREFIX
12 #define tsmonitorset_                        TSMONITORSET
13 #define tsdefaultcomputejacobian_            TSDEFAULTCOMPUTEJACOBIAN
14 #define tsdefaultcomputejacobiancolor_       TSDEFAULTCOMPUTEJACOBIANCOLOR
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 tssetmatrices_                       tssetmatrices
21 #define tsgetmatrices_                       tsgetmatrices
22 #define tssetrhsjacobian_                    tssetrhsjacobian
23 #define tsgetrhsjacobian_                    tsgetrhsjacobian
24 #define tsview_                              tsview
25 #define tsgetoptionsprefix_                  tsgetoptionsprefix
26 #define tsmonitorset_                        tsmonitorset
27 #define tsdefaultcomputejacobian_            tsdefaultcomputejacobian
28 #define tsdefaultcomputejacobiancolor_       tsdefaultcomputejacobiancolor
29 #define tsmonitordefault_                    tsmonitordefault
30 #define tssetprestep_                        tssetprestep
31 #define tssetpoststep_                       tssetpoststep
32 #endif
33 
34 static PetscErrorCode ourprestep(TS ts)
35 {
36   PetscErrorCode ierr = 0;
37   (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[8]))(&ts,&ierr);
38   return 0;
39 }
40 static PetscErrorCode ourpoststep(TS ts)
41 {
42   PetscErrorCode ierr = 0;
43   (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[9]))(&ts,&ierr);
44   return 0;
45 }
46 static PetscErrorCode ourtsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx)
47 {
48   PetscErrorCode ierr = 0;
49   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[1]))(&ts,&d,&x,&f,ctx,&ierr);
50   return 0;
51 }
52 static PetscErrorCode ourtsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
53 {
54   PetscErrorCode ierr = 0;
55   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[2]))(&ts,&d,m,p,type,ctx,&ierr);
56   return 0;
57 }
58 static PetscErrorCode ourtslhsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx)
59 {
60   PetscErrorCode ierr = 0;
61   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[7]))(&ts,&d,m,p,type,ctx,&ierr);
62   return 0;
63 }
64 static PetscErrorCode ourtsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx)
65 {
66   PetscErrorCode ierr = 0;
67   (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[3]))(&ts,&d,&x,m,p,type,ctx,&ierr);
68   return 0;
69 }
70 
71 static PetscErrorCode ourmonitordestroy(void **ctx)
72 {
73   PetscErrorCode ierr = 0;
74   TS          ts = *(TS*)ctx;
75   void        *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[6];
76   (*(void (PETSC_STDCALL *)(void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[5]))(mctx,&ierr);
77   return 0;
78 }
79 
80 /*
81    Note ctx is the same as ts so we need to get the Fortran context out of the TS
82 */
83 static PetscErrorCode ourtsmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx)
84 {
85   PetscErrorCode ierr = 0;
86   void           *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[6];
87   (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[4]))(&ts,&i,&d,&v,mctx,&ierr);
88   return 0;
89 }
90 
91 EXTERN_C_BEGIN
92 
93 void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
94 {
95   PetscObjectAllocateFortranPointers(*ts,10);
96   ((PetscObject)*ts)->fortran_func_pointers[8] = (PetscVoidFunction)f;
97   *ierr = TSSetPreStep(*ts,ourprestep);
98 }
99 
100 void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr)
101 {
102   PetscObjectAllocateFortranPointers(*ts,10);
103   ((PetscObject)*ts)->fortran_func_pointers[9] = (PetscVoidFunction)f;
104   *ierr = TSSetPreStep(*ts,ourpoststep);
105 }
106 
107 void PETSC_STDCALL tssetrhsfunction_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
108 {
109   PetscObjectAllocateFortranPointers(*ts,10);
110   ((PetscObject)*ts)->fortran_func_pointers[1] = (PetscVoidFunction)f;
111   *ierr = TSSetRHSFunction(*ts,ourtsfunction,fP);
112 }
113 
114 void PETSC_STDCALL tssetmatrices_(TS *ts,Mat *Arhs,PetscErrorCode (PETSC_STDCALL *frhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
115                                                    void*,PetscInt *),
116                                          Mat *Alhs,PetscErrorCode (PETSC_STDCALL *flhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,
117                                                    void*,PetscInt *),
118                                          MatStructure *flag,void*fP,PetscErrorCode *ierr)
119 {
120   PetscObjectAllocateFortranPointers(*ts,10);
121   if (FORTRANNULLFUNCTION(frhs) && FORTRANNULLFUNCTION(flhs)) {
122     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,PETSC_NULL,*flag,fP);
123   } else if (FORTRANNULLFUNCTION(flhs)){
124     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
125     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,PETSC_NULL,*flag,fP);
126   } else if (FORTRANNULLFUNCTION(frhs)){
127     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
128     *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,ourtslhsmatrix,*flag,fP);
129   } else {
130     ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs;
131     ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs;
132     *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,ourtslhsmatrix,*flag,fP);
133   }
134 }
135 
136 /* ---------------------------------------------------------*/
137 extern void tsdefaultcomputejacobian_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
138 extern void tsdefaultcomputejacobiancolor_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*);
139 
140 void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,
141                void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr)
142 {
143   PetscObjectAllocateFortranPointers(*ts,10);
144   if (FORTRANNULLFUNCTION(f)) {
145     *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP);
146   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobian_) {
147     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobian,fP);
148   } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobiancolor_) {
149     *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobianColor,*(MatFDColoring*)fP);
150   } else {
151   ((PetscObject)*ts)->fortran_func_pointers[3] = (PetscVoidFunction)f;
152     *ierr = TSSetRHSJacobian(*ts,*A,*B,ourtsjacobian,fP);
153   }
154 }
155 
156 /* ---------------------------------------------------------*/
157 
158 extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*);
159 
160 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)
161 {
162   PetscObjectAllocateFortranPointers(*ts,10);
163   if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) {
164     *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0);
165   } else {
166     ((PetscObject)*ts)->fortran_func_pointers[4] = (PetscVoidFunction)func;
167     ((PetscObject)*ts)->fortran_func_pointers[5] = (PetscVoidFunction)d;
168     ((PetscObject)*ts)->fortran_func_pointers[6] = (PetscVoidFunction)mctx;
169     if (FORTRANNULLFUNCTION(d)) {
170       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,0);
171     } else {
172       *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,ourmonitordestroy);
173     }
174   }
175 }
176 
177 /* ---------------------------------------------------------*/
178 void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr)
179 {
180   *ierr = TSGetRHSJacobian(*ts,J,M,ctx);
181 }
182 
183 void PETSC_STDCALL tsgetmatrices_(TS *ts,Mat *Arhs,Mat *Alhs,void **ctx,PetscErrorCode *ierr)
184 {
185   *ierr = TSGetMatrices(*ts,Arhs,Alhs,ctx);
186 }
187 
188 void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr)
189 {
190   PetscViewer v;
191   PetscPatchDefaultViewers_Fortran(viewer,v);
192   *ierr = TSView(*ts,v);
193 }
194 
195 void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
196 {
197   const char *tname;
198 
199   *ierr = TSGetOptionsPrefix(*ts,&tname);
200   *ierr = PetscStrncpy(prefix,tname,len);
201 }
202 
203 
204 EXTERN_C_END
205