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