xref: /petsc/src/tao/linesearch/interface/ftn-custom/ztaolinesearchf.c (revision 5e71baeff2f3138f93cd4f5927dfd596eb8325cc)
1 #include <petsc/private/fortranimpl.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 #define taolinesearchview_                           TAOLINESEARCHVIEW
10 #define taolinesearchsettype_                        TAOLINESEARCHSETTYPE
11 
12 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
13 
14 #define taolinesearchsetobjectiveroutine_            taolinesearchsetobjectiveroutine
15 #define taolinesearchsetgradientroutine_             taolinesearchsetgradientroutine
16 #define taolinesearchsetobjectiveandgradientroutine_ taolinesearchsetobjectiveandgradientroutine
17 #define taolinesearchsetobjectiveandgtsroutine_      taolinesearchsetobjectiveandgtsroutine
18 #define taolinesearchview_                           taolinesearchview
19 #define taolinesearchsettype_                        taolinesearchsettype
20 
21 #endif
22 
23 static int OBJ=0;
24 static int GRAD=1;
25 static int OBJGRAD=2;
26 static int OBJGTS=3;
27 static int NFUNCS=4;
28 
29 static PetscErrorCode ourtaolinesearchobjectiveroutine(TaoLineSearch ls, Vec x, PetscReal *f, void *ctx)
30 {
31     PetscErrorCode ierr = 0;
32     (*(void (PETSC_STDCALL *)(TaoLineSearch*,Vec*,PetscReal*,void*,PetscErrorCode*))
33         (((PetscObject)ls)->fortran_func_pointers[OBJ]))(&ls,&x,f,ctx,&ierr);
34     CHKERRQ(ierr);
35     return 0;
36 }
37 
38 static PetscErrorCode ourtaolinesearchgradientroutine(TaoLineSearch ls, Vec x, Vec g, void *ctx)
39 {
40     PetscErrorCode ierr = 0;
41     (*(void (PETSC_STDCALL *)(TaoLineSearch*,Vec*,Vec*,void*,PetscErrorCode*))
42        (((PetscObject)ls)->fortran_func_pointers[GRAD]))(&ls,&x,&g,ctx,&ierr);
43     CHKERRQ(ierr);
44     return 0;
45 
46 }
47 
48 static PetscErrorCode ourtaolinesearchobjectiveandgradientroutine(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, void* ctx)
49 {
50     PetscErrorCode ierr = 0;
51     (*(void (PETSC_STDCALL *)(TaoLineSearch*,Vec*,PetscReal*,Vec*,void*,PetscErrorCode*))
52      (((PetscObject)ls)->fortran_func_pointers[OBJGRAD]))(&ls,&x,f,&g,ctx,&ierr);
53     CHKERRQ(ierr);
54     return 0;
55 }
56 
57 static PetscErrorCode ourtaolinesearchobjectiveandgtsroutine(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, void* ctx)
58 {
59     PetscErrorCode ierr = 0;
60     (*(void (PETSC_STDCALL *)(TaoLineSearch*,Vec*,Vec*,PetscReal*,PetscReal*,void*,PetscErrorCode*))
61      (((PetscObject)ls)->fortran_func_pointers[OBJGTS]))(&ls,&x,&s,f,gts,ctx,&ierr);
62     CHKERRQ(ierr);
63     return 0;
64 }
65 
66 PETSC_EXTERN void PETSC_STDCALL taolinesearchsetobjectiveroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
67 {
68     CHKFORTRANNULLOBJECT(ctx);
69     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
70     if (!func) {
71         *ierr = TaoLineSearchSetObjectiveRoutine(*ls,0,ctx);
72     } else {
73         ((PetscObject)*ls)->fortran_func_pointers[OBJ] = (PetscVoidFunction)func;
74         *ierr = TaoLineSearchSetObjectiveRoutine(*ls, ourtaolinesearchobjectiveroutine,ctx);
75     }
76 }
77 
78 PETSC_EXTERN void PETSC_STDCALL taolinesearchsetgradientroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
79 {
80     CHKFORTRANNULLOBJECT(ctx);
81     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
82     if (!func) {
83         *ierr = TaoLineSearchSetGradientRoutine(*ls,0,ctx);
84     } else {
85         ((PetscObject)*ls)->fortran_func_pointers[GRAD] = (PetscVoidFunction)func;
86         *ierr = TaoLineSearchSetGradientRoutine(*ls, ourtaolinesearchgradientroutine,ctx);
87     }
88 }
89 
90 PETSC_EXTERN void PETSC_STDCALL taolinesearchsetobjectiveandgradientroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
91 {
92     CHKFORTRANNULLOBJECT(ctx);
93     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
94     if (!func) {
95         *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls,0,ctx);
96     } else {
97         ((PetscObject)*ls)->fortran_func_pointers[OBJGRAD] = (PetscVoidFunction)func;
98         *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, ourtaolinesearchobjectiveandgradientroutine,ctx);
99     }
100 }
101 
102 PETSC_EXTERN void PETSC_STDCALL taolinesearchsetobjectiveandgtsroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, Vec *, PetscReal*, PetscReal*,void*, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
103 {
104     CHKFORTRANNULLOBJECT(ctx);
105     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
106     if (!func) {
107         *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls,0,ctx);
108     } else {
109         ((PetscObject)*ls)->fortran_func_pointers[OBJGTS] = (PetscVoidFunction)func;
110         *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, ourtaolinesearchobjectiveandgtsroutine,ctx);
111     }
112 }
113 
114 PETSC_EXTERN void PETSC_STDCALL taolinesearchsettype_(TaoLineSearch *ls, char* type_name PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
115 
116 {
117     char *t;
118 
119     FIXCHAR(type_name,len,t);
120     *ierr = TaoLineSearchSetType(*ls,t);
121     FREECHAR(type_name,t);
122 
123 }
124 
125 PETSC_EXTERN void PETSC_STDCALL taolinesearchview_(TaoLineSearch *ls, PetscViewer *viewer, PetscErrorCode *ierr)
126 {
127     PetscViewer v;
128     PetscPatchDefaultViewers_Fortran(viewer,v);
129     *ierr = TaoLineSearchView(*ls,v);
130 }
131 
132 PETSC_EXTERN void PETSC_STDCALL taolinesearchgetoptionsprefix_(TaoLineSearch *ls, char* prefix PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
133 {
134   const char *name;
135   *ierr = TaoLineSearchGetOptionsPrefix(*ls,&name);
136   *ierr = PetscStrncpy(prefix,name,len); if (*ierr) return;
137   FIXRETURNCHAR(PETSC_TRUE,prefix,len);
138 
139 }
140 
141 PETSC_EXTERN void PETSC_STDCALL taolinesearchappendoptionsprefix_(TaoLineSearch *ls, char* prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
142 {
143   char *name;
144   FIXCHAR(prefix,len,name);
145   *ierr = TaoLineSearchAppendOptionsPrefix(*ls,name);
146   FREECHAR(prefix,name);
147 }
148 
149 PETSC_EXTERN void PETSC_STDCALL taolinesearchsetoptionsprefix_(TaoLineSearch *ls, char* prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
150 {
151   char *t;
152   FIXCHAR(prefix,len,t);
153   *ierr = TaoLineSearchSetOptionsPrefix(*ls,t);
154   FREECHAR(prefix,t);
155 }
156 
157 PETSC_EXTERN void PETSC_STDCALL taolinesearchgettype_(TaoLineSearch *ls, char* name PETSC_MIXED_LEN(len), PetscErrorCode *ierr  PETSC_END_LEN(len))
158 {
159   const char *tname;
160   *ierr = TaoLineSearchGetType(*ls,&tname);
161   *ierr = PetscStrncpy(name,tname,len); if (*ierr) return;
162   FIXRETURNCHAR(PETSC_TRUE,name,len);
163 
164 }
165