xref: /petsc/src/tao/interface/ftn-custom/ztaosolverf.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1 #include <petsc/private/ftnimpl.h>
2 #include <petsc/private/taoimpl.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define taosetobjective_                    TAOSETOBJECTIVE
6   #define taosetgradient_                     TAOSETGRADIENT
7   #define taosetobjectiveandgradient_         TAOSETOBJECTIVEANDGRADIENT
8   #define taosethessian_                      TAOSETHESSIAN
9   #define taosetresidualroutine_              TAOSETRESIDUALROUTINE
10   #define taosetjacobianresidualroutine_      TAOSETJACOBIANRESIDUALROUTINE
11   #define taosetjacobianroutine_              TAOSETJACOBIANROUTINE
12   #define taosetjacobianstateroutine_         TAOSETJACOBIANSTATEROUTINE
13   #define taosetjacobiandesignroutine_        TAOSETJACOBIANDESIGNROUTINE
14   #define taosetjacobianinequalityroutine_    TAOSETJACOBIANINEQUALITYROUTINE
15   #define taosetjacobianequalityroutine_      TAOSETJACOBIANEQUALITYROUTINE
16   #define taosetinequalityconstraintsroutine_ TAOSETINEQUALITYCONSTRAINTSROUTINE
17   #define taosetequalityconstraintsroutine_   TAOSETEQUALITYCONSTRAINTSROUTINE
18   #define taosetvariableboundsroutine_        TAOSETVARIABLEBOUNDSROUTINE
19   #define taosetconstraintsroutine_           TAOSETCONSTRAINTSROUTINE
20   #define taomonitorset_                      TAOMONITORSET
21   #define taogetconvergencehistory_           TAOGETCONVERGENCEHISTORY
22   #define taosetconvergencetest_              TAOSETCONVERGENCETEST
23   #define taosetupdate_                       TAOSETUPDATE
24 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
25   #define taosetobjective_                    taosetobjective
26   #define taosetgradient_                     taosetgradient
27   #define taosetobjectiveandgradient_         taosetobjectiveandgradient
28   #define taosethessian_                      taosethessian
29   #define taosetresidualroutine_              taosetresidualroutine
30   #define taosetjacobianresidualroutine_      taosetjacobianresidualroutine
31   #define taosetjacobianroutine_              taosetjacobianroutine
32   #define taosetjacobianstateroutine_         taosetjacobianstateroutine
33   #define taosetjacobiandesignroutine_        taosetjacobiandesignroutine
34   #define taosetjacobianinequalityroutine_    taosetjacobianinequalityroutine
35   #define taosetjacobianequalityroutine_      taosetjacobianequalityroutine
36   #define taosetinequalityconstraintsroutine_ taosetinequalityconstraintsroutine
37   #define taosetequalityconstraintsroutine_   taosetequalityconstraintsroutine
38   #define taosetvariableboundsroutine_        taosetvariableboundsroutine
39   #define taosetconstraintsroutine_           taosetconstraintsroutine
40   #define taomonitorset_                      taomonitorset
41   #define taogetconvergencehistory_           taogetconvergencehistory
42   #define taosetconvergencetest_              taosetconvergencetest
43   #define taosetupdate_                       taosetupdate
44 #endif
45 
46 static struct {
47   PetscFortranCallbackId obj;
48   PetscFortranCallbackId grad;
49   PetscFortranCallbackId objgrad;
50   PetscFortranCallbackId hess;
51   PetscFortranCallbackId lsres;
52   PetscFortranCallbackId lsjac;
53   PetscFortranCallbackId jac;
54   PetscFortranCallbackId jacstate;
55   PetscFortranCallbackId jacdesign;
56   PetscFortranCallbackId bounds;
57   PetscFortranCallbackId mon;
58   PetscFortranCallbackId mondestroy;
59   PetscFortranCallbackId convtest;
60   PetscFortranCallbackId constraints;
61   PetscFortranCallbackId jacineq;
62   PetscFortranCallbackId jaceq;
63   PetscFortranCallbackId conineq;
64   PetscFortranCallbackId coneq;
65   PetscFortranCallbackId nfuncs;
66   PetscFortranCallbackId update;
67 #if defined(PETSC_HAVE_F90_2PTR_ARG)
68   PetscFortranCallbackId function_pgiptr;
69 #endif
70 } _cb;
71 
ourtaoobjectiveroutine(Tao tao,Vec x,PetscReal * f,PetscCtx ctx)72 static PetscErrorCode ourtaoobjectiveroutine(Tao tao, Vec x, PetscReal *f, PetscCtx ctx)
73 {
74   PetscObjectUseFortranCallback(tao, _cb.obj, (Tao *, Vec *, PetscReal *, void *, PetscErrorCode *), (&tao, &x, f, _ctx, &ierr));
75 }
76 
ourtaogradientroutine(Tao tao,Vec x,Vec g,PetscCtx ctx)77 static PetscErrorCode ourtaogradientroutine(Tao tao, Vec x, Vec g, PetscCtx ctx)
78 {
79   PetscObjectUseFortranCallback(tao, _cb.grad, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &g, _ctx, &ierr));
80 }
81 
ourtaoobjectiveandgradientroutine(Tao tao,Vec x,PetscReal * f,Vec g,PetscCtx ctx)82 static PetscErrorCode ourtaoobjectiveandgradientroutine(Tao tao, Vec x, PetscReal *f, Vec g, PetscCtx ctx)
83 {
84   PetscObjectUseFortranCallback(tao, _cb.objgrad, (Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), (&tao, &x, f, &g, _ctx, &ierr));
85 }
86 
ourtaohessianroutine(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx)87 static PetscErrorCode ourtaohessianroutine(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx)
88 {
89   PetscObjectUseFortranCallback(tao, _cb.hess, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, _ctx, &ierr));
90 }
91 
ourtaojacobianroutine(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx)92 static PetscErrorCode ourtaojacobianroutine(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx)
93 {
94   PetscObjectUseFortranCallback(tao, _cb.jac, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, _ctx, &ierr));
95 }
96 
ourtaojacobianstateroutine(Tao tao,Vec x,Mat H,Mat Hpre,Mat Hinv,PetscCtx ctx)97 static PetscErrorCode ourtaojacobianstateroutine(Tao tao, Vec x, Mat H, Mat Hpre, Mat Hinv, PetscCtx ctx)
98 {
99   PetscObjectUseFortranCallback(tao, _cb.jacstate, (Tao *, Vec *, Mat *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, &Hinv, _ctx, &ierr));
100 }
101 
ourtaojacobiandesignroutine(Tao tao,Vec x,Mat H,PetscCtx ctx)102 static PetscErrorCode ourtaojacobiandesignroutine(Tao tao, Vec x, Mat H, PetscCtx ctx)
103 {
104   PetscObjectUseFortranCallback(tao, _cb.jacdesign, (Tao *, Vec *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, _ctx, &ierr));
105 }
106 
ourtaoboundsroutine(Tao tao,Vec xl,Vec xu,PetscCtx ctx)107 static PetscErrorCode ourtaoboundsroutine(Tao tao, Vec xl, Vec xu, PetscCtx ctx)
108 {
109   PetscObjectUseFortranCallback(tao, _cb.bounds, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &xl, &xu, _ctx, &ierr));
110 }
ourtaoresidualroutine(Tao tao,Vec x,Vec f,PetscCtx ctx)111 static PetscErrorCode ourtaoresidualroutine(Tao tao, Vec x, Vec f, PetscCtx ctx)
112 {
113   PetscObjectUseFortranCallback(tao, _cb.lsres, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &f, _ctx, &ierr));
114 }
115 
ourtaojacobianresidualroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)116 static PetscErrorCode ourtaojacobianresidualroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
117 {
118   PetscObjectUseFortranCallback(tao, _cb.lsjac, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
119 }
120 
ourtaomonitor(Tao tao,PetscCtx ctx)121 static PetscErrorCode ourtaomonitor(Tao tao, PetscCtx ctx)
122 {
123   PetscObjectUseFortranCallback(tao, _cb.mon, (Tao *, void *, PetscErrorCode *), (&tao, _ctx, &ierr));
124 }
125 
ourtaomondestroy(PetscCtxRt ctx)126 static PetscErrorCode ourtaomondestroy(PetscCtxRt ctx)
127 {
128   Tao tao = *(Tao *)ctx;
129   PetscObjectUseFortranCallback(tao, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
130 }
ourtaoconvergencetest(Tao tao,PetscCtx ctx)131 static PetscErrorCode ourtaoconvergencetest(Tao tao, PetscCtx ctx)
132 {
133   PetscObjectUseFortranCallback(tao, _cb.convtest, (Tao *, void *, PetscErrorCode *), (&tao, _ctx, &ierr));
134 }
135 
ourtaoconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)136 static PetscErrorCode ourtaoconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
137 {
138   PetscObjectUseFortranCallback(tao, _cb.constraints, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
139 }
140 
ourtaojacobianinequalityroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)141 static PetscErrorCode ourtaojacobianinequalityroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
142 {
143   PetscObjectUseFortranCallback(tao, _cb.jacineq, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
144 }
145 
ourtaojacobianequalityroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)146 static PetscErrorCode ourtaojacobianequalityroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
147 {
148   PetscObjectUseFortranCallback(tao, _cb.jaceq, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
149 }
150 
ourtaoinequalityconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)151 static PetscErrorCode ourtaoinequalityconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
152 {
153   PetscObjectUseFortranCallback(tao, _cb.conineq, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
154 }
155 
ourtaoequalityconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)156 static PetscErrorCode ourtaoequalityconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
157 {
158   PetscObjectUseFortranCallback(tao, _cb.coneq, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
159 }
160 
ourtaoupdateroutine(Tao tao,PetscInt iter,PetscCtx ctx)161 static PetscErrorCode ourtaoupdateroutine(Tao tao, PetscInt iter, PetscCtx ctx)
162 {
163   PetscObjectUseFortranCallback(tao, _cb.update, (Tao *, PetscInt *, void *), (&tao, &iter, _ctx));
164 }
165 
taosetobjective_(Tao * tao,void (* func)(Tao *,Vec *,PetscReal *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)166 PETSC_EXTERN void taosetobjective_(Tao *tao, void (*func)(Tao *, Vec *, PetscReal *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
167 {
168   CHKFORTRANNULLFUNCTION(func);
169   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.obj, (PetscFortranCallbackFn *)func, ctx);
170   if (!*ierr) *ierr = TaoSetObjective(*tao, ourtaoobjectiveroutine, ctx);
171 }
172 
taosetgradient_(Tao * tao,Vec * g,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)173 PETSC_EXTERN void taosetgradient_(Tao *tao, Vec *g, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
174 {
175   CHKFORTRANNULLFUNCTION(func);
176   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.grad, (PetscFortranCallbackFn *)func, ctx);
177   if (!*ierr) *ierr = TaoSetGradient(*tao, *g, ourtaogradientroutine, ctx);
178 }
179 
taosetobjectiveandgradient_(Tao * tao,Vec * g,void (* func)(Tao *,Vec *,PetscReal *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)180 PETSC_EXTERN void taosetobjectiveandgradient_(Tao *tao, Vec *g, void (*func)(Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
181 {
182   CHKFORTRANNULLFUNCTION(func);
183   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.objgrad, (PetscFortranCallbackFn *)func, ctx);
184   if (!*ierr) *ierr = TaoSetObjectiveAndGradient(*tao, *g, ourtaoobjectiveandgradientroutine, ctx);
185 }
186 
taosethessian_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)187 PETSC_EXTERN void taosethessian_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
188 {
189   CHKFORTRANNULLFUNCTION(func);
190   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.hess, (PetscFortranCallbackFn *)func, ctx);
191   if (!*ierr) *ierr = TaoSetHessian(*tao, *J, *Jp, ourtaohessianroutine, ctx);
192 }
193 
taosetresidualroutine_(Tao * tao,Vec * F,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)194 PETSC_EXTERN void taosetresidualroutine_(Tao *tao, Vec *F, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
195 {
196   CHKFORTRANNULLFUNCTION(func);
197   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.lsres, (PetscFortranCallbackFn *)func, ctx);
198   if (!*ierr) *ierr = TaoSetResidualRoutine(*tao, *F, ourtaoresidualroutine, ctx);
199 }
200 
taosetjacobianresidualroutine_(Tao * tao,Mat * J,Mat * Jpre,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)201 PETSC_EXTERN void taosetjacobianresidualroutine_(Tao *tao, Mat *J, Mat *Jpre, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
202 {
203   CHKFORTRANNULLFUNCTION(func);
204   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.lsjac, (PetscFortranCallbackFn *)func, ctx);
205   if (!*ierr) *ierr = TaoSetJacobianResidualRoutine(*tao, *J, *Jpre, ourtaojacobianresidualroutine, ctx);
206 }
207 
taosetjacobianroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)208 PETSC_EXTERN void taosetjacobianroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
209 {
210   CHKFORTRANNULLFUNCTION(func);
211   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jac, (PetscFortranCallbackFn *)func, ctx);
212   if (!*ierr) *ierr = TaoSetJacobianRoutine(*tao, *J, *Jp, ourtaojacobianroutine, ctx);
213 }
214 
taosetjacobianstateroutine_(Tao * tao,Mat * J,Mat * Jp,Mat * Jinv,void (* func)(Tao *,Vec *,Mat *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)215 PETSC_EXTERN void taosetjacobianstateroutine_(Tao *tao, Mat *J, Mat *Jp, Mat *Jinv, void (*func)(Tao *, Vec *, Mat *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
216 {
217   CHKFORTRANNULLFUNCTION(func);
218   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacstate, (PetscFortranCallbackFn *)func, ctx);
219   if (!*ierr) *ierr = TaoSetJacobianStateRoutine(*tao, *J, *Jp, *Jinv, ourtaojacobianstateroutine, ctx);
220 }
221 
taosetjacobiandesignroutine_(Tao * tao,Mat * J,void (* func)(Tao *,Vec *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)222 PETSC_EXTERN void taosetjacobiandesignroutine_(Tao *tao, Mat *J, void (*func)(Tao *, Vec *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
223 {
224   CHKFORTRANNULLFUNCTION(func);
225   *ierr = PetscObjectSetFortranCallback((PetscObject)tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacdesign, (PetscFortranCallbackFn *)func, ctx);
226   if (!*ierr) *ierr = TaoSetJacobianDesignRoutine(*tao, *J, ourtaojacobiandesignroutine, ctx);
227 }
228 
taosetvariableboundsroutine_(Tao * tao,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)229 PETSC_EXTERN void taosetvariableboundsroutine_(Tao *tao, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
230 {
231   CHKFORTRANNULLFUNCTION(func);
232   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.bounds, (PetscFortranCallbackFn *)func, ctx);
233   if (!*ierr) *ierr = TaoSetVariableBoundsRoutine(*tao, ourtaoboundsroutine, ctx);
234 }
235 
taomonitorset_(Tao * tao,void (* func)(Tao *,void *,PetscErrorCode *),PetscCtx ctx,void (* mondestroy)(void *,PetscErrorCode *),PetscErrorCode * ierr)236 PETSC_EXTERN void taomonitorset_(Tao *tao, void (*func)(Tao *, void *, PetscErrorCode *), PetscCtx ctx, void (*mondestroy)(void *, PetscErrorCode *), PetscErrorCode *ierr)
237 {
238   CHKFORTRANNULLFUNCTION(mondestroy);
239   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mon, (PetscFortranCallbackFn *)func, ctx);
240   if (*ierr) return;
241   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscFortranCallbackFn *)mondestroy, ctx);
242   if (*ierr) return;
243   *ierr = TaoMonitorSet(*tao, ourtaomonitor, *tao, ourtaomondestroy);
244 }
245 
taosetconvergencetest_(Tao * tao,void (* func)(Tao *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)246 PETSC_EXTERN void taosetconvergencetest_(Tao *tao, void (*func)(Tao *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
247 {
248   CHKFORTRANNULLFUNCTION(func);
249   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.convtest, (PetscFortranCallbackFn *)func, ctx);
250   if (!*ierr) *ierr = TaoSetConvergenceTest(*tao, ourtaoconvergencetest, ctx);
251 }
252 
taosetconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)253 PETSC_EXTERN void taosetconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
254 {
255   CHKFORTRANNULLFUNCTION(func);
256   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.constraints, (PetscFortranCallbackFn *)func, ctx);
257   if (!*ierr) *ierr = TaoSetConstraintsRoutine(*tao, *C, ourtaoconstraintsroutine, ctx);
258 }
259 
taogetconvergencehistory_(Tao * tao,PetscInt * nhist,PetscErrorCode * ierr)260 PETSC_EXTERN void taogetconvergencehistory_(Tao *tao, PetscInt *nhist, PetscErrorCode *ierr)
261 {
262   *ierr = TaoGetConvergenceHistory(*tao, NULL, NULL, NULL, NULL, nhist);
263 }
264 
taosetjacobianinequalityroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)265 PETSC_EXTERN void taosetjacobianinequalityroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
266 {
267   CHKFORTRANNULLFUNCTION(func);
268   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacineq, (PetscFortranCallbackFn *)func, ctx);
269   if (!*ierr) *ierr = TaoSetJacobianInequalityRoutine(*tao, *J, *Jp, ourtaojacobianinequalityroutine, ctx);
270 }
271 
taosetjacobianequalityroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)272 PETSC_EXTERN void taosetjacobianequalityroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
273 {
274   CHKFORTRANNULLFUNCTION(func);
275   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jaceq, (PetscFortranCallbackFn *)func, ctx);
276   if (!*ierr) *ierr = TaoSetJacobianEqualityRoutine(*tao, *J, *Jp, ourtaojacobianequalityroutine, ctx);
277 }
278 
taosetinequalityconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)279 PETSC_EXTERN void taosetinequalityconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
280 {
281   CHKFORTRANNULLFUNCTION(func);
282   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.conineq, (PetscFortranCallbackFn *)func, ctx);
283   if (!*ierr) *ierr = TaoSetInequalityConstraintsRoutine(*tao, *C, ourtaoinequalityconstraintsroutine, ctx);
284 }
285 
taosetequalityconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)286 PETSC_EXTERN void taosetequalityconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
287 {
288   CHKFORTRANNULLFUNCTION(func);
289   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.coneq, (PetscFortranCallbackFn *)func, ctx);
290   if (!*ierr) *ierr = TaoSetEqualityConstraintsRoutine(*tao, *C, ourtaoequalityconstraintsroutine, ctx);
291 }
292 
taosetupdate_(Tao * tao,void (* func)(Tao *,PetscInt *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)293 PETSC_EXTERN void taosetupdate_(Tao *tao, void (*func)(Tao *, PetscInt *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
294 {
295   CHKFORTRANNULLFUNCTION(func);
296   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.update, (PetscFortranCallbackFn *)func, ctx);
297   if (!*ierr) *ierr = TaoSetUpdate(*tao, ourtaoupdateroutine, ctx);
298 }
299