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