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 tsmonitorset_ TSMONITORSET 17 #define tscomputerhsfunctionlinear_ TSCOMPUTERHSFUNCTIONLINEAR 18 #define tscomputerhsjacobianconstant_ TSCOMPUTERHSJACOBIANCONSTANT 19 #define tscomputeifunctionlinear_ TSCOMPUTEIFUNCTIONLINEAR 20 #define tscomputeijacobianconstant_ TSCOMPUTEIJACOBIANCONSTANT 21 #define tsmonitordefault_ TSMONITORDEFAULT 22 #define tssetprestep_ TSSETPRESTEP 23 #define tssetpoststep_ TSSETPOSTSTEP 24 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 25 #define tsmonitorlgsettransform_ tsmonitorlgsettransform 26 #define tssetrhsfunction_ tssetrhsfunction 27 #define tsgetrhsfunction_ tsgetrhsfunction 28 #define tssetrhsjacobian_ tssetrhsjacobian 29 #define tsgetrhsjacobian_ tsgetrhsjacobian 30 #define tssetifunction_ tssetifunction 31 #define tsgetifunction_ tsgetifunction 32 #define tssetijacobian_ tssetijacobian 33 #define tsgetijacobian_ tsgetijacobian 34 #define tsmonitorset_ tsmonitorset 35 #define tscomputerhsfunctionlinear_ tscomputerhsfunctionlinear 36 #define tscomputerhsjacobianconstant_ tscomputerhsjacobianconstant 37 #define tscomputeifunctionlinear_ tscomputeifunctionlinear 38 #define tscomputeijacobianconstant_ tscomputeijacobianconstant 39 #define tsmonitordefault_ tsmonitordefault 40 #define tssetprestep_ tssetprestep 41 #define tssetpoststep_ tssetpoststep 42 #endif 43 44 static struct { 45 PetscFortranCallbackId prestep; 46 PetscFortranCallbackId poststep; 47 PetscFortranCallbackId rhsfunction; 48 PetscFortranCallbackId rhsjacobian; 49 PetscFortranCallbackId ifunction; 50 PetscFortranCallbackId ijacobian; 51 PetscFortranCallbackId monitor; 52 PetscFortranCallbackId mondestroy; 53 PetscFortranCallbackId transform; 54 #if defined(PETSC_HAVE_F90_2PTR_ARG) 55 PetscFortranCallbackId function_pgiptr; 56 #endif 57 } _cb; 58 59 static PetscErrorCode ourprestep(TS ts) 60 { 61 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 62 void *ptr; 63 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 64 #endif 65 PetscObjectUseFortranCallback(ts, _cb.prestep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 66 } 67 static PetscErrorCode ourpoststep(TS ts) 68 { 69 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 70 void *ptr; 71 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 72 #endif 73 PetscObjectUseFortranCallback(ts, _cb.poststep, (TS *, PetscErrorCode * /* PETSC_F90_2PTR_PROTO_NOVAR */), (&ts, &ierr /* PETSC_F90_2PTR_PARAM(ptr) */)); 74 } 75 static PetscErrorCode ourrhsfunction(TS ts, PetscReal d, Vec x, Vec f, void *ctx) 76 { 77 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 78 void *ptr; 79 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 80 #endif 81 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) */)); 82 } 83 static PetscErrorCode ourifunction(TS ts, PetscReal d, Vec x, Vec xdot, Vec f, void *ctx) 84 { 85 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 86 void *ptr; 87 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 88 #endif 89 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) */)); 90 } 91 static PetscErrorCode ourrhsjacobian(TS ts, PetscReal d, Vec x, Mat m, Mat p, void *ctx) 92 { 93 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 94 void *ptr; 95 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 96 #endif 97 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) */)); 98 } 99 static PetscErrorCode ourijacobian(TS ts, PetscReal d, Vec x, Vec xdot, PetscReal shift, Mat m, Mat p, void *ctx) 100 { 101 #if defined(PETSC_HAVE_F90_2PTR_ARG) && defined(foo) 102 void *ptr; 103 PetscCall(PetscObjectGetFortranCallback((PetscObject)ts, PETSC_FORTRAN_CALLBACK_CLASS, _cb.function_pgiptr, NULL, &ptr)); 104 #endif 105 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) */)); 106 } 107 108 static PetscErrorCode ourmonitordestroy(void **ctx) 109 { 110 TS ts = (TS)*ctx; 111 PetscObjectUseFortranCallback(ts, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr)); 112 } 113 114 /* 115 Note ctx is the same as ts so we need to get the Fortran context out of the TS 116 */ 117 static PetscErrorCode ourmonitor(TS ts, PetscInt i, PetscReal d, Vec v, void *ctx) 118 { 119 PetscObjectUseFortranCallback(ts, _cb.monitor, (TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), (&ts, &i, &d, &v, _ctx, &ierr)); 120 } 121 122 /* 123 Currently does not handle destroy or context 124 */ 125 static PetscErrorCode ourtransform(void *ctx, Vec x, Vec *xout) 126 { 127 PetscObjectUseFortranCallback((TS)ctx, _cb.transform, (void *, Vec *, Vec *, PetscErrorCode *), (_ctx, &x, xout, &ierr)); 128 } 129 130 PETSC_EXTERN void tsmonitorlgsettransform_(TS *ts, void (*f)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode (*d)(void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 131 { 132 *ierr = TSMonitorLGSetTransform(*ts, ourtransform, NULL, NULL); 133 if (*ierr) return; 134 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.transform, (PetscVoidFn *)f, ctx); 135 } 136 137 PETSC_EXTERN void tssetprestep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 138 { 139 *ierr = TSSetPreStep(*ts, ourprestep); 140 if (*ierr) return; 141 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.prestep, (PetscVoidFn *)f, NULL); 142 } 143 144 PETSC_EXTERN void tssetpoststep_(TS *ts, PetscErrorCode (*f)(TS *, PetscErrorCode *), PetscErrorCode *ierr) 145 { 146 *ierr = TSSetPostStep(*ts, ourpoststep); 147 if (*ierr) return; 148 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.poststep, (PetscVoidFn *)f, NULL); 149 } 150 151 PETSC_EXTERN void tscomputerhsfunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *F, void *ctx, PetscErrorCode *ierr) 152 { 153 *ierr = TSComputeRHSFunctionLinear(*ts, *t, *X, *F, ctx); 154 } 155 PETSC_EXTERN void tssetrhsfunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 156 { 157 Vec R; 158 CHKFORTRANNULLOBJECT(r); 159 CHKFORTRANNULLFUNCTION(f); 160 R = r ? *r : (Vec)NULL; 161 if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsfunctionlinear_) { 162 *ierr = TSSetRHSFunction(*ts, R, TSComputeRHSFunctionLinear, fP); 163 } else { 164 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsfunction, (PetscVoidFn *)f, fP); 165 *ierr = TSSetRHSFunction(*ts, R, ourrhsfunction, NULL); 166 } 167 } 168 PETSC_EXTERN void tsgetrhsfunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 169 { 170 CHKFORTRANNULLINTEGER(ctx); 171 CHKFORTRANNULLOBJECT(r); 172 *ierr = TSGetRHSFunction(*ts, r, NULL, ctx); 173 } 174 175 PETSC_EXTERN void tscomputeifunctionlinear_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, Vec *F, void *ctx, PetscErrorCode *ierr); 176 177 PETSC_EXTERN void tssetifunction_(TS *ts, Vec *r, PetscErrorCode (*f)(TS *, PetscReal *, Vec *, Vec *, Vec *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 178 { 179 Vec R; 180 CHKFORTRANNULLOBJECT(r); 181 CHKFORTRANNULLFUNCTION(f); 182 R = r ? *r : (Vec)NULL; 183 if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeifunctionlinear_) { 184 *ierr = TSSetIFunction(*ts, R, TSComputeIFunctionLinear, fP); 185 } else { 186 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ifunction, (PetscVoidFn *)f, fP); 187 *ierr = TSSetIFunction(*ts, R, ourifunction, NULL); 188 } 189 } 190 PETSC_EXTERN void tsgetifunction_(TS *ts, Vec *r, void *func, void **ctx, PetscErrorCode *ierr) 191 { 192 CHKFORTRANNULLINTEGER(ctx); 193 CHKFORTRANNULLOBJECT(r); 194 *ierr = TSGetIFunction(*ts, r, NULL, ctx); 195 } 196 197 /* ---------------------------------------------------------*/ 198 PETSC_EXTERN void tscomputerhsjacobianconstant_(TS *ts, PetscReal *t, Vec *X, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr) 199 { 200 *ierr = TSComputeRHSJacobianConstant(*ts, *t, *X, *A, *B, ctx); 201 } 202 PETSC_EXTERN void tssetrhsjacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 203 { 204 CHKFORTRANNULLFUNCTION(f); 205 if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputerhsjacobianconstant_) { 206 *ierr = TSSetRHSJacobian(*ts, *A, *B, TSComputeRHSJacobianConstant, fP); 207 } else { 208 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.rhsjacobian, (PetscVoidFn *)f, fP); 209 *ierr = TSSetRHSJacobian(*ts, *A, *B, ourrhsjacobian, NULL); 210 } 211 } 212 213 PETSC_EXTERN void tscomputeijacobianconstant_(TS *ts, PetscReal *t, Vec *X, Vec *Xdot, PetscReal *shift, Mat *A, Mat *B, void *ctx, PetscErrorCode *ierr); 214 215 PETSC_EXTERN void tssetijacobian_(TS *ts, Mat *A, Mat *B, void (*f)(TS *, PetscReal *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), void *fP, PetscErrorCode *ierr) 216 { 217 CHKFORTRANNULLFUNCTION(f); 218 if ((PetscVoidFn *)f == (PetscVoidFn *)tscomputeijacobianconstant_) { 219 *ierr = TSSetIJacobian(*ts, *A, *B, TSComputeIJacobianConstant, fP); 220 } else { 221 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.ijacobian, (PetscVoidFn *)f, fP); 222 *ierr = TSSetIJacobian(*ts, *A, *B, ourijacobian, NULL); 223 } 224 } 225 PETSC_EXTERN void tsgetijacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 226 { 227 CHKFORTRANNULLINTEGER(ctx); 228 CHKFORTRANNULLOBJECT(J); 229 CHKFORTRANNULLOBJECT(M); 230 *ierr = TSGetIJacobian(*ts, J, M, NULL, ctx); 231 } 232 233 PETSC_EXTERN void tsmonitordefault_(TS *, PetscInt *, PetscReal *, Vec *, PetscViewerAndFormat **, PetscErrorCode *); 234 235 /* ---------------------------------------------------------*/ 236 237 /* PETSC_EXTERN void tsmonitordefault_(TS*,PetscInt*,PetscReal*,Vec*,void*,PetscErrorCode*); */ 238 239 PETSC_EXTERN void tsmonitorset_(TS *ts, void (*func)(TS *, PetscInt *, PetscReal *, Vec *, void *, PetscErrorCode *), void *mctx, void (*d)(void *, PetscErrorCode *), PetscErrorCode *ierr) 240 { 241 CHKFORTRANNULLFUNCTION(d); 242 if ((PetscVoidFn *)func == (PetscVoidFn *)tsmonitordefault_) { 243 *ierr = TSMonitorSet(*ts, (PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *))TSMonitorDefault, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy); 244 } else { 245 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscVoidFn *)func, mctx); 246 *ierr = PetscObjectSetFortranCallback((PetscObject)*ts, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscVoidFn *)d, mctx); 247 *ierr = TSMonitorSet(*ts, ourmonitor, *ts, ourmonitordestroy); 248 } 249 } 250 251 /* ---------------------------------------------------------*/ 252 /* func is currently ignored from Fortran */ 253 PETSC_EXTERN void tsgetrhsjacobian_(TS *ts, Mat *J, Mat *M, int *func, void **ctx, PetscErrorCode *ierr) 254 { 255 *ierr = TSGetRHSJacobian(*ts, J, M, NULL, ctx); 256 } 257