1 #include <private/fortranimpl.h> 2 #include <petscts.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define tssetrhsfunction_ TSSETRHSFUNCTION 6 #define tssetrhsjacobian_ TSSETRHSJACOBIAN 7 #define tsgetrhsjacobian_ TSGETRHSJACOBIAN 8 #define tsview_ TSVIEW 9 #define tsgetoptionsprefix_ TSGETOPTIONSPREFIX 10 #define tsmonitorset_ TSMONITORSET 11 #define tsdefaultcomputejacobian_ TSDEFAULTCOMPUTEJACOBIAN 12 #define tsdefaultcomputejacobiancolor_ TSDEFAULTCOMPUTEJACOBIANCOLOR 13 #define tsmonitordefault_ TSMONITORDEFAULT 14 #define tssetprestep_ TSSETPRESTEP 15 #define tssetpoststep_ TSSETPOSTSTEP 16 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 17 #define tssetrhsfunction_ tssetrhsfunction 18 #define tssetrhsjacobian_ tssetrhsjacobian 19 #define tsgetrhsjacobian_ tsgetrhsjacobian 20 #define tsview_ tsview 21 #define tsgetoptionsprefix_ tsgetoptionsprefix 22 #define tsmonitorset_ tsmonitorset 23 #define tsdefaultcomputejacobian_ tsdefaultcomputejacobian 24 #define tsdefaultcomputejacobiancolor_ tsdefaultcomputejacobiancolor 25 #define tsmonitordefault_ tsmonitordefault 26 #define tssetprestep_ tssetprestep 27 #define tssetpoststep_ tssetpoststep 28 #endif 29 30 static PetscErrorCode ourprestep(TS ts) 31 { 32 PetscErrorCode ierr = 0; 33 (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[8]))(&ts,&ierr); 34 return 0; 35 } 36 static PetscErrorCode ourpoststep(TS ts) 37 { 38 PetscErrorCode ierr = 0; 39 (*(void (PETSC_STDCALL *)(TS*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[9]))(&ts,&ierr); 40 return 0; 41 } 42 static PetscErrorCode ourtsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx) 43 { 44 PetscErrorCode ierr = 0; 45 (*(void (PETSC_STDCALL *)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)ts)->fortran_func_pointers[1]))(&ts,&d,&x,&f,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 ourmonitordestroy(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 *)(void*,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*,void*,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 tssetprestep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 78 { 79 PetscObjectAllocateFortranPointers(*ts,10); 80 ((PetscObject)*ts)->fortran_func_pointers[8] = (PetscVoidFunction)f; 81 *ierr = TSSetPreStep(*ts,ourprestep); 82 } 83 84 void PETSC_STDCALL tssetpoststep_(TS *ts,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 85 { 86 PetscObjectAllocateFortranPointers(*ts,10); 87 ((PetscObject)*ts)->fortran_func_pointers[9] = (PetscVoidFunction)f; 88 *ierr = TSSetPreStep(*ts,ourpoststep); 89 } 90 91 void PETSC_STDCALL tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 92 { 93 PetscObjectAllocateFortranPointers(*ts,10); 94 ((PetscObject)*ts)->fortran_func_pointers[1] = (PetscVoidFunction)f; 95 *ierr = TSSetRHSFunction(*ts,*r,ourtsfunction,fP); 96 } 97 98 /* ---------------------------------------------------------*/ 99 extern void tsdefaultcomputejacobian_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 100 extern void tsdefaultcomputejacobiancolor_(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 101 102 void PETSC_STDCALL tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (PETSC_STDCALL *f)(TS*,PetscReal*,Vec*,Mat*,Mat*,MatStructure*, 103 void*,PetscErrorCode*),void*fP,PetscErrorCode *ierr) 104 { 105 PetscObjectAllocateFortranPointers(*ts,10); 106 if (FORTRANNULLFUNCTION(f)) { 107 *ierr = TSSetRHSJacobian(*ts,*A,*B,PETSC_NULL,fP); 108 } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobian_) { 109 *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobian,fP); 110 } else if ((PetscVoidFunction)f == (PetscVoidFunction)tsdefaultcomputejacobiancolor_) { 111 *ierr = TSSetRHSJacobian(*ts,*A,*B,TSDefaultComputeJacobianColor,*(MatFDColoring*)fP); 112 } else { 113 ((PetscObject)*ts)->fortran_func_pointers[3] = (PetscVoidFunction)f; 114 *ierr = TSSetRHSJacobian(*ts,*A,*B,ourtsjacobian,fP); 115 } 116 } 117 118 /* ---------------------------------------------------------*/ 119 120 extern void PETSC_STDCALL tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); 121 122 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) 123 { 124 PetscObjectAllocateFortranPointers(*ts,10); 125 if ((PetscVoidFunction)func == (PetscVoidFunction)tsmonitordefault_) { 126 *ierr = TSMonitorSet(*ts,TSMonitorDefault,0,0); 127 } else { 128 ((PetscObject)*ts)->fortran_func_pointers[4] = (PetscVoidFunction)func; 129 ((PetscObject)*ts)->fortran_func_pointers[5] = (PetscVoidFunction)d; 130 ((PetscObject)*ts)->fortran_func_pointers[6] = (PetscVoidFunction)mctx; 131 if (FORTRANNULLFUNCTION(d)) { 132 *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,0); 133 } else { 134 *ierr = TSMonitorSet(*ts,ourtsmonitor,*ts,ourmonitordestroy); 135 } 136 } 137 } 138 139 /* ---------------------------------------------------------*/ 140 /* func is currently ignored from Fortran */ 141 void PETSC_STDCALL tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr) 142 { 143 *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx); 144 } 145 146 void PETSC_STDCALL tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr) 147 { 148 PetscViewer v; 149 PetscPatchDefaultViewers_Fortran(viewer,v); 150 *ierr = TSView(*ts,v); 151 } 152 153 void PETSC_STDCALL tsgetoptionsprefix_(TS *ts,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 154 { 155 const char *tname; 156 157 *ierr = TSGetOptionsPrefix(*ts,&tname); 158 *ierr = PetscStrncpy(prefix,tname,len); 159 } 160 161 162 EXTERN_C_END 163