1*b45d2f2cSJed Brown #include <petsc-private/fortranimpl.h> 2c6db04a5SJed Brown #include <petscts.h> 3e2df7a95SSatish Balay 4e2df7a95SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS) 5e2df7a95SSatish Balay #define tssetrhsfunction_ TSSETRHSFUNCTION 6e2df7a95SSatish Balay #define tssetrhsjacobian_ TSSETRHSJACOBIAN 7e2df7a95SSatish Balay #define tsgetrhsjacobian_ TSGETRHSJACOBIAN 8e2df7a95SSatish Balay #define tsview_ TSVIEW 9e2df7a95SSatish Balay #define tsgetoptionsprefix_ TSGETOPTIONSPREFIX 10a6570f20SBarry Smith #define tsmonitorset_ TSMONITORSET 110e4ef248SJed Brown #define tscomputerhsfunctionlinear_ TSCOMPUTERHSFUNCTIONLINEAR 120e4ef248SJed Brown #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT 130fecffdcSJed Brown #define tscomputeifunctionlinear_ TSCOMPUTEIFUNCTIONLINEAR 140fecffdcSJed Brown #define tscomputeijacobianconstant_ TSCOMPUTEIJACOBIANCONSTANT 15a6570f20SBarry Smith #define tsmonitordefault_ TSMONITORDEFAULT 16dd7ecb2fSBarry Smith #define tssetprestep_ TSSETPRESTEP 17dd7ecb2fSBarry Smith #define tssetpoststep_ TSSETPOSTSTEP 18e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 19e2df7a95SSatish Balay #define tssetrhsfunction_ tssetrhsfunction 20e2df7a95SSatish Balay #define tssetrhsjacobian_ tssetrhsjacobian 21e2df7a95SSatish Balay #define tsgetrhsjacobian_ tsgetrhsjacobian 22e2df7a95SSatish Balay #define tsview_ tsview 23e2df7a95SSatish Balay #define tsgetoptionsprefix_ tsgetoptionsprefix 24a6570f20SBarry Smith #define tsmonitorset_ tsmonitorset 250e4ef248SJed Brown #define tscomputerhsfunctionlinear_ tscomputerhsfunctionlinear 260e4ef248SJed Brown #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant 270fecffdcSJed Brown #define tscomputeifunctionlinear_ tscomputeifunctionlinear 280fecffdcSJed Brown #define tscomputeijacobianconstant_ tscomputeijacobianconstant 29a6570f20SBarry Smith #define tsmonitordefault_ tsmonitordefault 30dd7ecb2fSBarry Smith #define tssetprestep_ tssetprestep 31dd7ecb2fSBarry Smith #define tssetpoststep_ tssetpoststep 32e2df7a95SSatish Balay #endif 33e2df7a95SSatish Balay 340fecffdcSJed Brown enum {OUR_PRESTEP = 0, 350fecffdcSJed Brown OUR_POSTSTEP, 360fecffdcSJed Brown OUR_RHSFUNCTION, 370fecffdcSJed Brown OUR_IFUNCTION, 380fecffdcSJed Brown OUR_RHSJACOBIAN, 390fecffdcSJed Brown OUR_IJACOBIAN, 400fecffdcSJed Brown OUR_MONITOR, 410fecffdcSJed Brown OUR_MONITORDESTROY, 420fecffdcSJed Brown OUR_MONITOR_CTX, /* Casting from function pointer is invalid according to the standard. */ 430fecffdcSJed Brown OUR_COUNT}; 440fecffdcSJed Brown 45dd7ecb2fSBarry Smith static PetscErrorCode ourprestep(TS ts) 46dd7ecb2fSBarry Smith { 47dd7ecb2fSBarry Smith PetscErrorCode ierr = 0; 480fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_PRESTEP]))(&ts,&ierr); 49dd7ecb2fSBarry Smith return 0; 50dd7ecb2fSBarry Smith } 51dd7ecb2fSBarry Smith static PetscErrorCode ourpoststep(TS ts) 52dd7ecb2fSBarry Smith { 53dd7ecb2fSBarry Smith PetscErrorCode ierr = 0; 540fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_POSTSTEP]))(&ts,&ierr); 55dd7ecb2fSBarry Smith return 0; 56dd7ecb2fSBarry Smith } 570fecffdcSJed Brown static PetscErrorCode ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx) 58e2df7a95SSatish Balay { 59e2df7a95SSatish Balay PetscErrorCode ierr = 0; 600fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSFUNCTION]))(&ts,&d,&x,&f,ctx,&ierr); 61e2df7a95SSatish Balay return 0; 62e2df7a95SSatish Balay } 630fecffdcSJed Brown static PetscErrorCode ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,void *ctx) 64e2df7a95SSatish Balay { 65e2df7a95SSatish Balay PetscErrorCode ierr = 0; 660fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IFUNCTION]))(&ts,&d,&x,&xdot,&f,ctx,&ierr); 670fecffdcSJed Brown return 0; 680fecffdcSJed Brown } 690fecffdcSJed Brown static PetscErrorCode ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx) 700fecffdcSJed Brown { 710fecffdcSJed Brown PetscErrorCode ierr = 0; 720fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_RHSJACOBIAN]))(&ts,&d,&x,m,p,type,ctx,&ierr); 730fecffdcSJed Brown return 0; 740fecffdcSJed Brown } 750fecffdcSJed Brown static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat* m,Mat* p,MatStructure* type,void*ctx) 760fecffdcSJed Brown { 770fecffdcSJed Brown PetscErrorCode ierr = 0; 780fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_IJACOBIAN]))(&ts,&d,&x,&xdot,&shift,m,p,type,ctx,&ierr); 79e2df7a95SSatish Balay return 0; 80e2df7a95SSatish Balay } 81e2df7a95SSatish Balay 82c2efdce3SBarry Smith static PetscErrorCode ourmonitordestroy(void **ctx) 83e2df7a95SSatish Balay { 84e2df7a95SSatish Balay PetscErrorCode ierr = 0; 85c2efdce3SBarry Smith TS ts = *(TS*)ctx; 860fecffdcSJed Brown void *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX]; 870fecffdcSJed Brown (*(void (PETSC_STDCALL *)(void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITORDESTROY]))(mctx,&ierr); 88e2df7a95SSatish Balay return 0; 89e2df7a95SSatish Balay } 90e2df7a95SSatish Balay 91e2df7a95SSatish Balay /* 92e2df7a95SSatish Balay Note ctx is the same as ts so we need to get the Fortran context out of the TS 93e2df7a95SSatish Balay */ 940fecffdcSJed Brown static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx) 95e2df7a95SSatish Balay { 96e2df7a95SSatish Balay PetscErrorCode ierr = 0; 970fecffdcSJed Brown void *mctx = (void*) ((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR_CTX]; 980fecffdcSJed Brown (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[OUR_MONITOR]))(&ts,&i,&d,&v,mctx,&ierr); 99e2df7a95SSatish Balay return 0; 100e2df7a95SSatish Balay } 101e2df7a95SSatish Balay 102e2df7a95SSatish Balay EXTERN_C_BEGIN 103e2df7a95SSatish Balay 104dd7ecb2fSBarry Smith void PETSC_STDCALL tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 105dd7ecb2fSBarry Smith { 1060fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 1070fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_PRESTEP] = (PetscVoidFunction)f; 108dd7ecb2fSBarry Smith *ierr = TSSetPreStep(*ts,ourprestep); 109dd7ecb2fSBarry Smith } 110dd7ecb2fSBarry Smith 111dd7ecb2fSBarry Smith void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 112dd7ecb2fSBarry Smith { 1130fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 1140fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_POSTSTEP] = (PetscVoidFunction)f; 115dd7ecb2fSBarry Smith *ierr = TSSetPreStep(*ts,ourpoststep); 116dd7ecb2fSBarry Smith } 117dd7ecb2fSBarry Smith 1180e4ef248SJed Brown void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr) 1190e4ef248SJed Brown { 1200e4ef248SJed Brown *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx); 1210e4ef248SJed Brown } 122089b2837SJed Brown void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 123e2df7a95SSatish Balay { 1240e4ef248SJed Brown Vec R; 1250e4ef248SJed Brown CHKFORTRANNULLOBJECT(r); 1260e4ef248SJed Brown CHKFORTRANNULLFUNCTION(f); 1270e4ef248SJed Brown CHKFORTRANNULLOBJECT(fP); 1281db2f0e5SJed Brown R = r ? *r : (Vec)PETSC_NULL; 1290e4ef248SJed Brown if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) { 1300e4ef248SJed Brown *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP); 1310e4ef248SJed Brown } else { 1320fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 1330fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSFUNCTION] = (PetscVoidFunction)f; 1340fecffdcSJed Brown *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,fP); 1350fecffdcSJed Brown } 1360fecffdcSJed Brown } 1370fecffdcSJed Brown 1380fecffdcSJed Brown void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr) 1390fecffdcSJed Brown { 1400fecffdcSJed Brown *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx); 1410fecffdcSJed Brown } 1420fecffdcSJed Brown void PETSC_STDCALL tssetifunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 1430fecffdcSJed Brown { 1440fecffdcSJed Brown Vec R; 1450fecffdcSJed Brown CHKFORTRANNULLOBJECT(r); 1460fecffdcSJed Brown CHKFORTRANNULLFUNCTION(f); 1470fecffdcSJed Brown CHKFORTRANNULLOBJECT(fP); 1481db2f0e5SJed Brown R = r ? *r : (Vec)PETSC_NULL; 1490fecffdcSJed Brown if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) { 1500fecffdcSJed Brown *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP); 1510fecffdcSJed Brown } else { 1520fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 1530fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_IFUNCTION] = (PetscVoidFunction)f; 1540fecffdcSJed Brown *ierr = TSSetIFunction(*ts,R,ourifunction,fP); 1550e4ef248SJed Brown } 156e2df7a95SSatish Balay } 15726d46c62SHong Zhang 158e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 1590e4ef248SJed Brown void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr) 1600e4ef248SJed Brown { 1610e4ef248SJed Brown *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,A,B,flg,ctx); 1620e4ef248SJed Brown } 163e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*, 164e2df7a95SSatish Balay void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 165e2df7a95SSatish Balay { 1660fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 167e2df7a95SSatish Balay if (FORTRANNULLFUNCTION(f)) { 168e2df7a95SSatish Balay *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP); 1690e4ef248SJed Brown } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) { 1700e4ef248SJed Brown *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP); 171e2df7a95SSatish Balay } else { 1720fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_RHSJACOBIAN] = (PetscVoidFunction)f; 1730fecffdcSJed Brown *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,fP); 1740fecffdcSJed Brown } 1750fecffdcSJed Brown } 1760fecffdcSJed Brown 1770fecffdcSJed Brown void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,MatStructure *flg,void *ctx,PetscErrorCode *ierr) 1780fecffdcSJed Brown { 1790fecffdcSJed Brown *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,A,B,flg,ctx); 1800fecffdcSJed Brown } 1810fecffdcSJed Brown void PETSC_STDCALL tssetijacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*, 1820fecffdcSJed Brown void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 1830fecffdcSJed Brown { 1840fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 1850fecffdcSJed Brown if (FORTRANNULLFUNCTION(f)) { 1860fecffdcSJed Brown *ierr = TSSetIJacobian(*ts,*A,*B,PETSC_NULL,fP); 1870fecffdcSJed Brown } else if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) { 1880fecffdcSJed Brown *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP); 1890fecffdcSJed Brown } else { 1900fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_IJACOBIAN] = (PetscVoidFunction)f; 1910fecffdcSJed Brown *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,fP); 192e2df7a95SSatish Balay } 193e2df7a95SSatish Balay } 194e2df7a95SSatish Balay 195e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 196e2df7a95SSatish Balay 197a6570f20SBarry Smith extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); 198e2df7a95SSatish Balay 199a6570f20SBarry Smith 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) 200e2df7a95SSatish Balay { 2010fecffdcSJed Brown PetscObjectAllocateFortranPointers(*ts,OUR_COUNT); 202a6570f20SBarry Smith if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) { 203a6570f20SBarry Smith *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0); 204e2df7a95SSatish Balay } else { 2050fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR] = (PetscVoidFunction)func; 2060fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITORDESTROY] = (PetscVoidFunction)d; 2070fecffdcSJed Brown ((PetscObject)*ts)->fortran_func_pointers[OUR_MONITOR_CTX] = (PetscVoidFunction)mctx; 208e2df7a95SSatish Balay if (FORTRANNULLFUNCTION(d)) { 2090fecffdcSJed Brown *ierr = TSMonitorSet(*ts,ourmonitor,*ts,0); 210e2df7a95SSatish Balay } else { 2110fecffdcSJed Brown *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy); 212e2df7a95SSatish Balay } 213e2df7a95SSatish Balay } 214e2df7a95SSatish Balay } 215e2df7a95SSatish Balay 216e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 217089b2837SJed Brown /* func is currently ignored from Fortran */ 218089b2837SJed Brown void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr) 219e2df7a95SSatish Balay { 220089b2837SJed Brown *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx); 221e2df7a95SSatish Balay } 222e2df7a95SSatish Balay 223e2df7a95SSatish Balay void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr) 224e2df7a95SSatish Balay { 225e2df7a95SSatish Balay PetscViewer v; 226e2df7a95SSatish Balay PetscPatchDefaultViewers_Fortran(viewer,v); 227e2df7a95SSatish Balay *ierr = TSView(*ts,v); 228e2df7a95SSatish Balay } 229e2df7a95SSatish Balay 230e2df7a95SSatish Balay void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 231e2df7a95SSatish Balay { 232e2df7a95SSatish Balay const char *tname; 233e2df7a95SSatish Balay 234e2df7a95SSatish Balay *ierr = TSGetOptionsPrefix(*ts,&tname); 235e2df7a95SSatish Balay *ierr = PetscStrncpy(prefix,tname,len); 236e2df7a95SSatish Balay } 237e2df7a95SSatish Balay 238e2df7a95SSatish Balay 239e2df7a95SSatish Balay EXTERN_C_END 240