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