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 6*76f2fa84SHong Zhang #define tssetmatrices_ TSSETMATRICES 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 12a6570f20SBarry Smith #define tsmonitorset_ TSMONITORSET 13e2df7a95SSatish Balay #define tsdefaultcomputejacobian_ TSDEFAULTCOMPUTEJACOBIAN 14e2df7a95SSatish Balay #define tsdefaultcomputejacobiancolor_ TSDEFAULTCOMPUTEJACOBIANCOLOR 15a6570f20SBarry Smith #define tsmonitordefault_ TSMONITORDEFAULT 16e2df7a95SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 17e2df7a95SSatish Balay #define tssetrhsfunction_ tssetrhsfunction 18*76f2fa84SHong Zhang #define tssetmatrices_ tssetmatrices 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 24a6570f20SBarry Smith #define tsmonitorset_ tsmonitorset 25e2df7a95SSatish Balay #define tsdefaultcomputejacobian_ tsdefaultcomputejacobian 26e2df7a95SSatish Balay #define tsdefaultcomputejacobiancolor_ tsdefaultcomputejacobiancolor 27a6570f20SBarry 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 } 42*76f2fa84SHong Zhang static PetscErrorCode ourtslhsmatrix(TS ts,PetscReal d,Mat* m,Mat* p,MatStructure* type,void*ctx) 43*76f2fa84SHong Zhang { 44*76f2fa84SHong Zhang PetscErrorCode ierr = 0; 45*76f2fa84SHong Zhang (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[7]))(&ts,&d,m,p,type,ctx,&ierr); 46*76f2fa84SHong Zhang return 0; 47*76f2fa84SHong Zhang } 48e2df7a95SSatish Balay static PetscErrorCode ourtsjacobian(TS ts,PetscReal d,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx) 49e2df7a95SSatish Balay { 50e2df7a95SSatish Balay PetscErrorCode ierr = 0; 51e2df7a95SSatish 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); 52e2df7a95SSatish Balay return 0; 53e2df7a95SSatish Balay } 54e2df7a95SSatish Balay 55e2df7a95SSatish Balay static PetscErrorCode ourtsdestroy(void *ctx) 56e2df7a95SSatish Balay { 57e2df7a95SSatish Balay PetscErrorCode ierr = 0; 58e2df7a95SSatish Balay TS ts = (TS)ctx; 59e2df7a95SSatish Balay void (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6]; 60f68b968cSBarry Smith (*(void (PETSC_STDCALL *)(PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[5]))(mctx,&ierr); 61e2df7a95SSatish Balay return 0; 62e2df7a95SSatish Balay } 63e2df7a95SSatish Balay 64e2df7a95SSatish Balay /* 65e2df7a95SSatish Balay Note ctx is the same as ts so we need to get the Fortran context out of the TS 66e2df7a95SSatish Balay */ 67e2df7a95SSatish Balay static PetscErrorCode ourtsmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void*ctx) 68e2df7a95SSatish Balay { 69e2df7a95SSatish Balay PetscErrorCode ierr = 0; 70e2df7a95SSatish Balay void (*mctx)(void) = ((PetscObject)ts)->fortran_func_pointers[6]; 71f68b968cSBarry Smith (*(void (PETSC_STDCALL *)(TS*,PetscInt*,PetscReal*,Vec*,PetscVoidFunction,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[4]))(&ts,&i,&d,&v,mctx,&ierr); 72e2df7a95SSatish Balay return 0; 73e2df7a95SSatish Balay } 74e2df7a95SSatish Balay 75e2df7a95SSatish Balay EXTERN_C_BEGIN 76e2df7a95SSatish Balay 77e2df7a95SSatish Balay 78e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsfunction_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 79e2df7a95SSatish Balay { 80f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[1] = (PetscVoidFunction)f; 81e2df7a95SSatish Balay *ierr = TSSetRHSFunction(*ts,ourtsfunction,fP); 82e2df7a95SSatish Balay } 83*76f2fa84SHong Zhang #ifdef MV 84e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsmatrix_(TS *ts,Mat *A,Mat *B,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Mat*,Mat*,MatStructure*, 85e2df7a95SSatish Balay void*,PetscInt *),void*fP,PetscErrorCode *ierr) 86e2df7a95SSatish Balay { 87e2df7a95SSatish Balay if (FORTRANNULLFUNCTION(f)) { 88e2df7a95SSatish Balay *ierr = TSSetRHSMatrix(*ts,*A,*B,PETSC_NULL,fP); 89e2df7a95SSatish Balay } else { 90f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)f; 91e2df7a95SSatish Balay *ierr = TSSetRHSMatrix(*ts,*A,*B,ourtsmatrix,fP); 92e2df7a95SSatish Balay } 93e2df7a95SSatish Balay } 94*76f2fa84SHong Zhang #endif 95*76f2fa84SHong Zhang void PETSC_STDCALL tssetmatrices_(TS *ts,Mat *Arhs,PetscErrorCode (PETSC_STDCALL *frhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*, 96*76f2fa84SHong Zhang void*,PetscInt *), 97*76f2fa84SHong Zhang Mat *Alhs,PetscErrorCode (PETSC_STDCALL *flhs)(TS*,PetscReal*,Mat*,Mat*,MatStructure*, 98*76f2fa84SHong Zhang void*,PetscInt *), 99*76f2fa84SHong Zhang MatStructure *flag,void*fP,PetscErrorCode *ierr) 100*76f2fa84SHong Zhang { 101*76f2fa84SHong Zhang if (FORTRANNULLFUNCTION(frhs) && FORTRANNULLFUNCTION(flhs)) { 102*76f2fa84SHong Zhang *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,PETSC_NULL,*flag,fP); 103*76f2fa84SHong Zhang } else if (FORTRANNULLFUNCTION(flhs)){ 104*76f2fa84SHong Zhang ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs; 105*76f2fa84SHong Zhang *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,PETSC_NULL,*flag,fP); 106*76f2fa84SHong Zhang } else if (FORTRANNULLFUNCTION(frhs)){ 107*76f2fa84SHong Zhang ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs; 108*76f2fa84SHong Zhang *ierr = TSSetMatrices(*ts,*Arhs,PETSC_NULL,*Alhs,ourtslhsmatrix,*flag,fP); 109*76f2fa84SHong Zhang } else { 110*76f2fa84SHong Zhang ((PetscObject)*ts)->fortran_func_pointers[2] = (PetscVoidFunction)frhs; 111*76f2fa84SHong Zhang ((PetscObject)*ts)->fortran_func_pointers[7] = (PetscVoidFunction)flhs; 112*76f2fa84SHong Zhang *ierr = TSSetMatrices(*ts,*Arhs,ourtsmatrix,*Alhs,ourtslhsmatrix,*flag,fP); 113*76f2fa84SHong Zhang } 114*76f2fa84SHong Zhang } 115e2df7a95SSatish Balay 116e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 117e2df7a95SSatish Balay extern void tsdefaultcomputejacobian_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 118e2df7a95SSatish Balay extern void tsdefaultcomputejacobiancolor_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 119e2df7a95SSatish Balay 120e2df7a95SSatish Balay void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*, 121e2df7a95SSatish Balay void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 122e2df7a95SSatish Balay { 123e2df7a95SSatish Balay if (FORTRANNULLFUNCTION(f)) { 124e2df7a95SSatish Balay *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP); 125f68b968cSBarry Smith } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobian_) { 126e2df7a95SSatish Balay *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobian,fP); 127f68b968cSBarry Smith } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobiancolor_) { 128e2df7a95SSatish Balay *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobianColor,*(MatFDColoring*)fP); 129e2df7a95SSatish Balay } else { 130f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[3] = (PetscVoidFunction)f; 131e2df7a95SSatish Balay *ierr = TSSetRHSJacobian(*ts,*A,*B,ourtsjacobian,fP); 132e2df7a95SSatish Balay } 133e2df7a95SSatish Balay } 134e2df7a95SSatish Balay 135e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 136e2df7a95SSatish Balay 137a6570f20SBarry Smith extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); 138e2df7a95SSatish Balay 139a6570f20SBarry 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) 140e2df7a95SSatish Balay { 141a6570f20SBarry Smith if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) { 142a6570f20SBarry Smith *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0); 143e2df7a95SSatish Balay } else { 144f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[4] = (PetscVoidFunction)func; 145f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[5] = (PetscVoidFunction)d; 146f68b968cSBarry Smith ((PetscObject)*ts)->fortran_func_pointers[6] = (PetscVoidFunction)mctx; 147e2df7a95SSatish Balay if (FORTRANNULLFUNCTION(d)) { 148a6570f20SBarry Smith *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,0); 149e2df7a95SSatish Balay } else { 150a6570f20SBarry Smith *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,ourtsdestroy); 151e2df7a95SSatish Balay } 152e2df7a95SSatish Balay } 153e2df7a95SSatish Balay } 154e2df7a95SSatish Balay 155e2df7a95SSatish Balay /* ---------------------------------------------------------*/ 156e2df7a95SSatish Balay void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr) 157e2df7a95SSatish Balay { 158e2df7a95SSatish Balay *ierr = TSGetRHSJacobian(*ts,J,M,ctx); 159e2df7a95SSatish Balay } 160e2df7a95SSatish Balay 161e2df7a95SSatish Balay void PETSC_STDCALL tsgetrhsmatrix_(TS *ts,Mat *J,Mat *M,void **ctx,PetscErrorCode *ierr) 162e2df7a95SSatish Balay { 163e2df7a95SSatish Balay *ierr = TSGetRHSMatrix(*ts,J,M,ctx); 164e2df7a95SSatish Balay } 165e2df7a95SSatish Balay 166e2df7a95SSatish Balay void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr) 167e2df7a95SSatish Balay { 168e2df7a95SSatish Balay PetscViewer v; 169e2df7a95SSatish Balay PetscPatchDefaultViewers_Fortran(viewer,v); 170e2df7a95SSatish Balay *ierr = TSView(*ts,v); 171e2df7a95SSatish Balay } 172e2df7a95SSatish Balay 173e2df7a95SSatish Balay void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 174e2df7a95SSatish Balay { 175e2df7a95SSatish Balay const char *tname; 176e2df7a95SSatish Balay 177e2df7a95SSatish Balay *ierr = TSGetOptionsPrefix(*ts,&tname); 178e2df7a95SSatish Balay #if defined(PETSC_USES_CPTOFCD) 179e2df7a95SSatish Balay { 180e2df7a95SSatish Balay char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix); 181e2df7a95SSatish Balay *ierr = PetscStrncpy(t,tname,len1); 182e2df7a95SSatish Balay } 183e2df7a95SSatish Balay #else 184e2df7a95SSatish Balay *ierr = PetscStrncpy(prefix,tname,len); 185e2df7a95SSatish Balay #endif 186e2df7a95SSatish Balay } 187e2df7a95SSatish Balay 188e2df7a95SSatish Balay 189e2df7a95SSatish Balay EXTERN_C_END 190