xref: /petsc/src/tao/linesearch/interface/ftn-custom/ztaolinesearchf.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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 EXTERN_C_BEGIN
67 
68 
69 
70 void PETSC_STDCALL taolinesearchsetobjectiveroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
71 {
72     CHKFORTRANNULLOBJECT(ctx);
73     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
74     if (!func) {
75         *ierr = TaoLineSearchSetObjectiveRoutine(*ls,0,ctx);
76     } else {
77         ((PetscObject)*ls)->fortran_func_pointers[OBJ] = (PetscVoidFunction)func;
78         *ierr = TaoLineSearchSetObjectiveRoutine(*ls, ourtaolinesearchobjectiveroutine,ctx);
79     }
80 }
81 
82 void PETSC_STDCALL taolinesearchsetgradientroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
83 {
84     CHKFORTRANNULLOBJECT(ctx);
85     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
86     if (!func) {
87         *ierr = TaoLineSearchSetGradientRoutine(*ls,0,ctx);
88     } else {
89         ((PetscObject)*ls)->fortran_func_pointers[GRAD] = (PetscVoidFunction)func;
90         *ierr = TaoLineSearchSetGradientRoutine(*ls, ourtaolinesearchgradientroutine,ctx);
91     }
92 }
93 
94 void PETSC_STDCALL taolinesearchsetobjectiveandgradientroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
95 {
96     CHKFORTRANNULLOBJECT(ctx);
97     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
98     if (!func) {
99         *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls,0,ctx);
100     } else {
101         ((PetscObject)*ls)->fortran_func_pointers[OBJGRAD] = (PetscVoidFunction)func;
102         *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, ourtaolinesearchobjectiveandgradientroutine,ctx);
103     }
104 }
105 
106 
107 
108 
109 void PETSC_STDCALL taolinesearchsetobjectiveandgtsroutine_(TaoLineSearch *ls, void (PETSC_STDCALL *func)(TaoLineSearch*, Vec *, Vec *, PetscReal*, PetscReal*,void*, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
110 {
111     CHKFORTRANNULLOBJECT(ctx);
112     PetscObjectAllocateFortranPointers(*ls,NFUNCS);
113     if (!func) {
114         *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls,0,ctx);
115     } else {
116         ((PetscObject)*ls)->fortran_func_pointers[OBJGTS] = (PetscVoidFunction)func;
117         *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, ourtaolinesearchobjectiveandgtsroutine,ctx);
118     }
119 }
120 
121 
122 
123 void PETSC_STDCALL taolinesearchsettype_(TaoLineSearch *ls, CHAR type_name PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
124 
125 {
126     char *t;
127 
128     FIXCHAR(type_name,len,t);
129     *ierr = TaoLineSearchSetType(*ls,t);
130     FREECHAR(type_name,t);
131 
132 }
133 
134 void PETSC_STDCALL taolinesearchview_(TaoLineSearch *ls, PetscViewer *viewer, PetscErrorCode *ierr)
135 {
136     PetscViewer v;
137     PetscPatchDefaultViewers_Fortran(viewer,v);
138     *ierr = TaoLineSearchView(*ls,v);
139 }
140 
141 void PETSC_STDCALL taolinesearchgetoptionsprefix_(TaoLineSearch *ls, CHAR prefix PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
142 {
143   const char *name;
144   *ierr = TaoLineSearchGetOptionsPrefix(*ls,&name);
145   *ierr = PetscStrncpy(prefix,name,len); if (*ierr) return;
146   FIXRETURNCHAR(PETSC_TRUE,prefix,len);
147 
148 }
149 
150 void PETSC_STDCALL taolinesearchappendoptionsprefix_(TaoLineSearch *ls, CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
151 {
152   char *name;
153   FIXCHAR(prefix,len,name);
154   *ierr = TaoLineSearchAppendOptionsPrefix(*ls,name);
155   FREECHAR(prefix,name);
156 }
157 
158 void PETSC_STDCALL taolinesearchsetoptionsprefix_(TaoLineSearch *ls, CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
159 {
160   char *t;
161   FIXCHAR(prefix,len,t);
162   *ierr = TaoLineSearchSetOptionsPrefix(*ls,t);
163   FREECHAR(prefix,t);
164 }
165 
166 void PETSC_STDCALL taolinesearchgettype_(TaoLineSearch *ls, CHAR name PETSC_MIXED_LEN(len), PetscErrorCode *ierr  PETSC_END_LEN(len))
167 {
168   const char *tname;
169   *ierr = TaoLineSearchGetType(*ls,&tname);
170   *ierr = PetscStrncpy(name,tname,len); if (*ierr) return;
171   FIXRETURNCHAR(PETSC_TRUE,name,len);
172 
173 }
174 EXTERN_C_END
175