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