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