1 #include <petsc/private/fortranimpl.h> 2 #include <petscts.h> 3 #include <petscviewer.h> 4 #include <petsc/private/f90impl.h> 5 6 #if defined(PETSC_HAVE_FORTRAN_CAPS) 7 #define tsmonitorlgsettransform_ TSMONITORLGSETTRANSFORM 8 #define tssetrhsfunction_ TSSETRHSFUNCTION 9 #define tsgetrhsfunction_ TSGETRHSFUNCTION 10 #define tssetrhsjacobian_ TSSETRHSJACOBIAN 11 #define tsgetrhsjacobian_ TSGETRHSJACOBIAN 12 #define tssetifunction_ TSSETIFUNCTION 13 #define tsgetifunction_ TSGETIFUNCTION 14 #define tssetijacobian_ TSSETIJACOBIAN 15 #define tsgetijacobian_ TSGETIJACOBIAN 16 #define tsview_ TSVIEW 17 #define tssetoptionsprefix_ TSSETOPTIONSPREFIX 18 #define tsgetoptionsprefix_ TSGETOPTIONSPREFIX 19 #define tsappendoptionsprefix_ TSAPPENDOPTIONSPREFIX 20 #define tsmonitorset_ TSMONITORSET 21 #define tscomputerhsfunctionlinear_ TSCOMPUTERHSFUNCTIONLINEAR 22 #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT 23 #define tscomputeifunctionlinear_ TSCOMPUTEIFUNCTIONLINEAR 24 #define tscomputeijacobianconstant_ TSCOMPUTEIJACOBIANCONSTANT 25 #define tsmonitordefault_ TSMONITORDEFAULT 26 #define tssetprestep_ TSSETPRESTEP 27 #define tssetpoststep_ TSSETPOSTSTEP 28 #define tsviewfromoptions_ TSVIEWFROMOPTIONS 29 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 30 #define tsmonitorlgsettransform_ tsmonitorlgsettransform 31 #define tssetrhsfunction_ tssetrhsfunction 32 #define tsgetrhsfunction_ tsgetrhsfunction 33 #define tssetrhsjacobian_ tssetrhsjacobian 34 #define tsgetrhsjacobian_ tsgetrhsjacobian 35 #define tssetifunction_ tssetifunction 36 #define tsgetifunction_ tsgetifunction 37 #define tssetijacobian_ tssetijacobian 38 #define tsgetijacobian_ tsgetijacobian 39 #define tsview_ tsview 40 #define tssetoptionsprefix_ tssetoptionsprefix 41 #define tsgetoptionsprefix_ tsgetoptionsprefix 42 #define tsappendoptionsprefix_ tsappendoptionsprefix 43 #define tsmonitorset_ tsmonitorset 44 #define tscomputerhsfunctionlinear_ tscomputerhsfunctionlinear 45 #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant 46 #define tscomputeifunctionlinear_ tscomputeifunctionlinear 47 #define tscomputeijacobianconstant_ tscomputeijacobianconstant 48 #define tsmonitordefault_ tsmonitordefault 49 #define tssetprestep_ tssetprestep 50 #define tssetpoststep_ tssetpoststep 51 #define tsviewfromoptions_ tsviewfromoptions 52 #endif 53 54 static struct { 55 PetscFortranCallbackId prestep; 56 PetscFortranCallbackId poststep; 57 PetscFortranCallbackId rhsfunction; 58 PetscFortranCallbackId rhsjacobian; 59 PetscFortranCallbackId ifunction; 60 PetscFortranCallbackId ijacobian; 61 PetscFortranCallbackId monitor; 62 PetscFortranCallbackId mondestroy; 63 PetscFortranCallbackId transform; 64 #if defined(PETSC_HAVE_F90_2PTR_ARG) 65 PetscFortranCallbackId function_pgiptr; 66 #endif 67 } _cb; 68 69 static PetscErrorCode ourprestep(TS ts) 70 { 71 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 72 void* ptr; 73 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 74 #endif 75 PetscObjectUseFortranCallback(ts,_cb.prestep,(TS*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 76 return 0; 77 } 78 static PetscErrorCode ourpoststep(TS ts) 79 { 80 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 81 void* ptr; 82 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 83 #endif 84 PetscObjectUseFortranCallback(ts,_cb.poststep,(TS*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 85 return 0; 86 } 87 static PetscErrorCode ourrhsfunction(TS ts,PetscReal d,Vec x,Vec f,void *ctx) 88 { 89 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 90 void* ptr; 91 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 92 #endif 93 PetscObjectUseFortranCallback(ts,_cb.rhsfunction,(TS*,PetscReal*,Vec*, Vec*, void*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&d,&x,&f,_ctx,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 94 return 0; 95 } 96 static PetscErrorCode ourifunction(TS ts,PetscReal d,Vec x,Vec xdot,Vec f,void *ctx) 97 { 98 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 99 void* ptr; 100 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 101 #endif 102 PetscObjectUseFortranCallback(ts,_cb.ifunction,(TS*,PetscReal*,Vec*, Vec*, Vec*, void*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&d,&x,&xdot,&f,_ctx,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 103 return 0; 104 } 105 static PetscErrorCode ourrhsjacobian(TS ts,PetscReal d,Vec x,Mat m,Mat p,void *ctx) 106 { 107 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 108 void* ptr; 109 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 110 #endif 111 PetscObjectUseFortranCallback(ts,_cb.rhsjacobian,(TS*,PetscReal*, Vec*, Mat*, Mat*, void*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&d,&x,&m,&p,_ctx,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 112 return 0; 113 } 114 static PetscErrorCode ourijacobian(TS ts,PetscReal d,Vec x,Vec xdot,PetscReal shift,Mat m,Mat p,void *ctx) 115 { 116 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 117 void* ptr; 118 PetscObjectGetFortranCallback((PetscObject)ts,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 119 #endif 120 PetscObjectUseFortranCallback(ts,_cb.ijacobian,(TS*,PetscReal*,Vec*, Vec*, PetscReal *,Mat*, Mat*, void*,PetscErrorCode* /* PETSC_F90_2PTR_PROTO_NOVAR */),(&ts,&d,&x,&xdot,&shift,&m,&p,_ctx,&ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 121 return 0; 122 } 123 124 static PetscErrorCode ourmonitordestroy(void **ctx) 125 { 126 TS ts = (TS)*ctx; 127 PetscObjectUseFortranCallback(ts,_cb.mondestroy,(void*,PetscErrorCode*),(_ctx,&ierr)); 128 } 129 130 /* 131 Note ctx is the same as ts so we need to get the Fortran context out of the TS 132 */ 133 static PetscErrorCode ourmonitor(TS ts,PetscInt i,PetscReal d,Vec v,void *ctx) 134 { 135 PetscObjectUseFortranCallback(ts,_cb.monitor,(TS*,PetscInt*,PetscReal*,Vec *,void*,PetscErrorCode*),(&ts,&i,&d,&v,_ctx,&ierr)); 136 return 0; 137 } 138 139 /* 140 Currently does not handle destroy or context 141 */ 142 static PetscErrorCode ourtransform(void *ctx,Vec x,Vec *xout) 143 { 144 PetscObjectUseFortranCallback((TS)ctx,_cb.transform,(void*,Vec *,Vec *,PetscErrorCode*),(_ctx,&x,xout,&ierr)); 145 } 146 147 PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts,void (*f)(void*,Vec*,Vec*,PetscErrorCode*),PetscErrorCode (*d)(void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr) 148 { 149 *ierr = TSMonitorLGSetTransform(*ts,ourtransform,NULL,NULL); if (*ierr) return; 150 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.transform,(PetscVoidFunction)f,ctx); 151 } 152 153 PETSC_EXTERN void tssetprestep_(TS *ts,PetscErrorCode (*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 154 { 155 *ierr = TSSetPreStep(*ts,ourprestep);if (*ierr) return; 156 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.prestep,(PetscVoidFunction)f,NULL); 157 } 158 159 PETSC_EXTERN void tssetpoststep_(TS *ts,PetscErrorCode (*f)(TS*,PetscErrorCode*),PetscErrorCode *ierr) 160 { 161 *ierr = TSSetPostStep(*ts,ourpoststep);if (*ierr) return; 162 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.poststep,(PetscVoidFunction)f,NULL); 163 } 164 165 PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *F,void *ctx,PetscErrorCode *ierr) 166 { 167 *ierr = TSComputeRHSFunctionLinear(*ts,*t,*X,*F,ctx); 168 } 169 PETSC_EXTERN void tssetrhsfunction_(TS *ts,Vec *r,PetscErrorCode (*f)(TS*,PetscReal*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr) 170 { 171 Vec R; 172 CHKFORTRANNULLOBJECT(r); 173 CHKFORTRANNULLFUNCTION(f); 174 R = r ? *r : (Vec)NULL; 175 if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsfunctionlinear_) { 176 *ierr = TSSetRHSFunction(*ts,R,TSComputeRHSFunctionLinear,fP); 177 } else { 178 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.rhsfunction,(PetscVoidFunction)f,fP); 179 *ierr = TSSetRHSFunction(*ts,R,ourrhsfunction,NULL); 180 } 181 } 182 PETSC_EXTERN void tsgetrhsfunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 183 { 184 CHKFORTRANNULLINTEGER(ctx); 185 CHKFORTRANNULLOBJECT(r); 186 *ierr = TSGetRHSFunction(*ts,r,NULL,ctx); 187 } 188 189 PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,Vec *F,void *ctx,PetscErrorCode *ierr) 190 { 191 *ierr = TSComputeIFunctionLinear(*ts,*t,*X,*Xdot,*F,ctx); 192 } 193 PETSC_EXTERN void tssetifunction_(TS *ts,Vec *r,PetscErrorCode (*f)(TS*,PetscReal*,Vec*,Vec*,Vec*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr) 194 { 195 Vec R; 196 CHKFORTRANNULLOBJECT(r); 197 CHKFORTRANNULLFUNCTION(f); 198 R = r ? *r : (Vec)NULL; 199 if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeifunctionlinear_) { 200 *ierr = TSSetIFunction(*ts,R,TSComputeIFunctionLinear,fP); 201 } else { 202 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ifunction,(PetscVoidFunction)f,fP); 203 *ierr = TSSetIFunction(*ts,R,ourifunction,NULL); 204 } 205 } 206 PETSC_EXTERN void tsgetifunction_(TS *ts,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 207 { 208 CHKFORTRANNULLINTEGER(ctx); 209 CHKFORTRANNULLOBJECT(r); 210 *ierr = TSGetIFunction(*ts,r,NULL,ctx); 211 } 212 213 /* ---------------------------------------------------------*/ 214 PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts,PetscReal *t,Vec *X,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr) 215 { 216 *ierr = TSComputeRHSJacobianConstant(*ts,*t,*X,*A,*B,ctx); 217 } 218 PETSC_EXTERN void tssetrhsjacobian_(TS *ts,Mat *A,Mat *B,void (*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr) 219 { 220 CHKFORTRANNULLFUNCTION(f); 221 if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputerhsjacobianconstant_) { 222 *ierr = TSSetRHSJacobian(*ts,*A,*B,TSComputeRHSJacobianConstant,fP); 223 } else { 224 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.rhsjacobian,(PetscVoidFunction)f,fP); 225 *ierr = TSSetRHSJacobian(*ts,*A,*B,ourrhsjacobian,NULL); 226 } 227 } 228 229 PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts,PetscReal *t,Vec *X,Vec *Xdot,PetscReal *shift,Mat *A,Mat *B,void *ctx,PetscErrorCode *ierr) 230 { 231 *ierr = TSComputeIJacobianConstant(*ts,*t,*X,*Xdot,*shift,*A,*B,ctx); 232 } 233 PETSC_EXTERN void tssetijacobian_(TS *ts,Mat *A,Mat *B,void (*f)(TS*,PetscReal*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),void *fP,PetscErrorCode *ierr) 234 { 235 CHKFORTRANNULLFUNCTION(f); 236 if ((PetscVoidFunction)f == (PetscVoidFunction)tscomputeijacobianconstant_) { 237 *ierr = TSSetIJacobian(*ts,*A,*B,TSComputeIJacobianConstant,fP); 238 } else { 239 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ijacobian,(PetscVoidFunction)f,fP); 240 *ierr = TSSetIJacobian(*ts,*A,*B,ourijacobian,NULL); 241 } 242 } 243 PETSC_EXTERN void tsgetijacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr) 244 { 245 CHKFORTRANNULLINTEGER(ctx); 246 CHKFORTRANNULLOBJECT(J); 247 CHKFORTRANNULLOBJECT(M); 248 *ierr = TSGetIJacobian(*ts,J,M,0,ctx); 249 } 250 251 PETSC_EXTERN void tsmonitordefault_(TS *ts,PetscInt *its,PetscReal *fgnorm,Vec *u,PetscViewerAndFormat **dummy,PetscErrorCode *ierr) 252 { 253 *ierr = TSMonitorDefault(*ts,*its,*fgnorm,*u,*dummy); 254 } 255 256 /* ---------------------------------------------------------*/ 257 258 /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */ 259 260 PETSC_EXTERN void tsmonitorset_(TS *ts,void (*func)(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*),void *mctx,void (*d)(void*,PetscErrorCode*),PetscErrorCode *ierr) 261 { 262 CHKFORTRANNULLFUNCTION(d); 263 if ((PetscVoidFunction)func == (PetscVoidFunction) tsmonitordefault_) { 264 *ierr = TSMonitorSet(*ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,void*))TSMonitorDefault,*(PetscViewerAndFormat**)mctx,(PetscErrorCode (*)(void **))PetscViewerAndFormatDestroy); 265 } else { 266 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.monitor,(PetscVoidFunction)func,mctx); 267 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.mondestroy,(PetscVoidFunction)d,mctx); 268 *ierr = TSMonitorSet(*ts,ourmonitor,*ts,ourmonitordestroy); 269 } 270 } 271 272 /* ---------------------------------------------------------*/ 273 /* func is currently ignored from Fortran */ 274 PETSC_EXTERN void tsgetrhsjacobian_(TS *ts,Mat *J,Mat *M,int *func,void **ctx,PetscErrorCode *ierr) 275 { 276 *ierr = TSGetRHSJacobian(*ts,J,M,0,ctx); 277 } 278 279 PETSC_EXTERN void tsview_(TS *ts,PetscViewer *viewer, PetscErrorCode *ierr) 280 { 281 PetscViewer v; 282 PetscPatchDefaultViewers_Fortran(viewer,v); 283 *ierr = TSView(*ts,v); 284 } 285 286 PETSC_EXTERN void tssetoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 287 { 288 char *t; 289 FIXCHAR(prefix,len,t); 290 *ierr = TSSetOptionsPrefix(*ts,t);if (*ierr) return; 291 FREECHAR(prefix,t); 292 } 293 PETSC_EXTERN void tsgetoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 294 { 295 const char *tname; 296 297 *ierr = TSGetOptionsPrefix(*ts,&tname); 298 *ierr = PetscStrncpy(prefix,tname,len); 299 FIXRETURNCHAR(PETSC_TRUE,prefix,len); 300 } 301 PETSC_EXTERN void tsappendoptionsprefix_(TS *ts,char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 302 { 303 char *t; 304 FIXCHAR(prefix,len,t); 305 *ierr = TSAppendOptionsPrefix(*ts,t);if (*ierr) return; 306 FREECHAR(prefix,t); 307 } 308 309 PETSC_EXTERN void tsviewfromoptions_(TS *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 310 { 311 char *t; 312 313 FIXCHAR(type,len,t); 314 *ierr = TSViewFromOptions(*ao,obj,t);if (*ierr) return; 315 FREECHAR(type,t); 316 } 317