1 #include <petsc/private/ftnimpl.h> 2 #include <petsc/private/taolinesearchimpl.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define taolinesearchsetobjectiveroutine_ TAOLINESEARCHSETOBJECTIVEROUTINE 6 #define taolinesearchsetgradientroutine_ TAOLINESEARCHSETGRADIENTROUTINE 7 #define taolinesearchsetobjectiveandgradientroutine_ TAOLINESEARCHSETOBJECTIVEANDGRADIENTROUTINE 8 #define taolinesearchsetobjectiveandgtsroutine_ TAOLINESEARCHSETOBJECTIVEANDGTSROUTINE 9 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 10 #define taolinesearchsetobjectiveroutine_ taolinesearchsetobjectiveroutine 11 #define taolinesearchsetgradientroutine_ taolinesearchsetgradientroutine 12 #define taolinesearchsetobjectiveandgradientroutine_ taolinesearchsetobjectiveandgradientroutine 13 #define taolinesearchsetobjectiveandgtsroutine_ taolinesearchsetobjectiveandgtsroutine 14 #endif 15 16 static int OBJ = 0; 17 static int GRAD = 1; 18 static int OBJGRAD = 2; 19 static int OBJGTS = 3; 20 static size_t NFUNCS = 4; 21 22 static PetscErrorCode ourtaolinesearchobjectiveroutine(TaoLineSearch ls, Vec x, PetscReal *f, void *ctx) 23 { 24 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch *, Vec *, PetscReal *, void *, PetscErrorCode *))(((PetscObject)ls)->fortran_func_pointers[OBJ]))(&ls, &x, f, ctx, &ierr)); 25 return PETSC_SUCCESS; 26 } 27 28 static PetscErrorCode ourtaolinesearchgradientroutine(TaoLineSearch ls, Vec x, Vec g, void *ctx) 29 { 30 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch *, Vec *, Vec *, void *, PetscErrorCode *))(((PetscObject)ls)->fortran_func_pointers[GRAD]))(&ls, &x, &g, ctx, &ierr)); 31 return PETSC_SUCCESS; 32 } 33 34 static PetscErrorCode ourtaolinesearchobjectiveandgradientroutine(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, void *ctx) 35 { 36 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *))(((PetscObject)ls)->fortran_func_pointers[OBJGRAD]))(&ls, &x, f, &g, ctx, &ierr)); 37 return PETSC_SUCCESS; 38 } 39 40 static PetscErrorCode ourtaolinesearchobjectiveandgtsroutine(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, void *ctx) 41 { 42 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch *, Vec *, Vec *, PetscReal *, PetscReal *, void *, PetscErrorCode *))(((PetscObject)ls)->fortran_func_pointers[OBJGTS]))(&ls, &x, &s, f, gts, ctx, &ierr)); 43 return PETSC_SUCCESS; 44 } 45 46 PETSC_EXTERN void taolinesearchsetobjectiveroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch *, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 47 { 48 PetscObjectAllocateFortranPointers(*ls, NFUNCS); 49 if (!func) { 50 *ierr = TaoLineSearchSetObjectiveRoutine(*ls, NULL, ctx); 51 } else { 52 ((PetscObject)*ls)->fortran_func_pointers[OBJ] = (PetscFortranCallbackFn *)func; 53 *ierr = TaoLineSearchSetObjectiveRoutine(*ls, ourtaolinesearchobjectiveroutine, ctx); 54 } 55 } 56 57 PETSC_EXTERN void taolinesearchsetgradientroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch *, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 58 { 59 PetscObjectAllocateFortranPointers(*ls, NFUNCS); 60 if (!func) { 61 *ierr = TaoLineSearchSetGradientRoutine(*ls, NULL, ctx); 62 } else { 63 ((PetscObject)*ls)->fortran_func_pointers[GRAD] = (PetscFortranCallbackFn *)func; 64 *ierr = TaoLineSearchSetGradientRoutine(*ls, ourtaolinesearchgradientroutine, ctx); 65 } 66 } 67 68 PETSC_EXTERN void taolinesearchsetobjectiveandgradientroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 69 { 70 PetscObjectAllocateFortranPointers(*ls, NFUNCS); 71 if (!func) { 72 *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, NULL, ctx); 73 } else { 74 ((PetscObject)*ls)->fortran_func_pointers[OBJGRAD] = (PetscFortranCallbackFn *)func; 75 *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, ourtaolinesearchobjectiveandgradientroutine, ctx); 76 } 77 } 78 79 PETSC_EXTERN void taolinesearchsetobjectiveandgtsroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch *, Vec *, Vec *, PetscReal *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 80 { 81 PetscObjectAllocateFortranPointers(*ls, NFUNCS); 82 if (!func) { 83 *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, NULL, ctx); 84 } else { 85 ((PetscObject)*ls)->fortran_func_pointers[OBJGTS] = (PetscFortranCallbackFn *)func; 86 *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, ourtaolinesearchobjectiveandgtsroutine, ctx); 87 } 88 } 89