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