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