1 #include <petsc/private/fortranimpl.h> 2 #include <petsc/private/f90impl.h> 3 #include <petsc/private/taoimpl.h> 4 5 #if defined(PETSC_HAVE_FORTRAN_CAPS) 6 #define taobrgnsetregularizerobjectiveandgradientroutine_ TAOBRGNSETREGULARIZEROBJECTIVEANDGRADIENTROUTINE 7 #define taobrgnsetregularizerhessianroutine_ TAOBRGNSETREGULARIZERHESSIANROUTINE 8 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 9 #define taobrgnsetregularizerobjectiveandgradientroutine_ taobrgnsetregularizerobjectiveandgradientroutine 10 #define taobrgnsetregularizerhessianroutine_ taobrgnsetregularizerhessianroutine 11 #endif 12 13 static struct { 14 PetscFortranCallbackId objgrad; 15 PetscFortranCallbackId hess; 16 } _cb; 17 18 static PetscErrorCode ourtaobrgnregobjgradroutine(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx) 19 { 20 PetscObjectUseFortranCallback(tao, _cb.objgrad, (Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), (&tao, &x, f, &g, _ctx, &ierr)); 21 } 22 23 static PetscErrorCode ourtaobrgnreghessroutine(Tao tao, Vec x, Mat H, void *ctx) 24 { 25 PetscObjectUseFortranCallback(tao, _cb.hess, (Tao *, Vec *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, _ctx, &ierr)); 26 } 27 28 PETSC_EXTERN void taobrgnsetregularizerobjectiveandgradientroutine_(Tao *tao, void (*func)(Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 29 { 30 CHKFORTRANNULLFUNCTION(func); 31 *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.objgrad, (PetscVoidFn *)func, ctx); 32 if (!*ierr) *ierr = TaoBRGNSetRegularizerObjectiveAndGradientRoutine(*tao, ourtaobrgnregobjgradroutine, ctx); 33 } 34 35 PETSC_EXTERN void taobrgnsetregularizerhessianroutine_(Tao *tao, Mat *H, void (*func)(Tao *, Vec *, Mat *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 36 { 37 CHKFORTRANNULLFUNCTION(func); 38 *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.hess, (PetscVoidFn *)func, ctx); 39 if (!*ierr) *ierr = TaoBRGNSetRegularizerHessianRoutine(*tao, *H, ourtaobrgnreghessroutine, ctx); 40 } 41