xref: /petsc/src/tao/interface/ftn-custom/ztaosolverf.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2af0996ceSBarry Smith #include <petsc/private/taoimpl.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
5a82e8c82SStefano Zampini   #define taosetobjective_                    TAOSETOBJECTIVE
6a82e8c82SStefano Zampini   #define taosetgradient_                     TAOSETGRADIENT
7a82e8c82SStefano Zampini   #define taosetobjectiveandgradient_         TAOSETOBJECTIVEANDGRADIENT
8a82e8c82SStefano Zampini   #define taosethessian_                      TAOSETHESSIAN
9737f463aSAlp Dener   #define taosetresidualroutine_              TAOSETRESIDUALROUTINE
104ffbe8acSAlp Dener   #define taosetjacobianresidualroutine_      TAOSETJACOBIANRESIDUALROUTINE
11a7e14dcfSSatish Balay   #define taosetjacobianroutine_              TAOSETJACOBIANROUTINE
12a7e14dcfSSatish Balay   #define taosetjacobianstateroutine_         TAOSETJACOBIANSTATEROUTINE
13a7e14dcfSSatish Balay   #define taosetjacobiandesignroutine_        TAOSETJACOBIANDESIGNROUTINE
14a7e14dcfSSatish Balay   #define taosetjacobianinequalityroutine_    TAOSETJACOBIANINEQUALITYROUTINE
15a7e14dcfSSatish Balay   #define taosetjacobianequalityroutine_      TAOSETJACOBIANEQUALITYROUTINE
16a7e14dcfSSatish Balay   #define taosetinequalityconstraintsroutine_ TAOSETINEQUALITYCONSTRAINTSROUTINE
17a7e14dcfSSatish Balay   #define taosetequalityconstraintsroutine_   TAOSETEQUALITYCONSTRAINTSROUTINE
18a7e14dcfSSatish Balay   #define taosetvariableboundsroutine_        TAOSETVARIABLEBOUNDSROUTINE
19a7e14dcfSSatish Balay   #define taosetconstraintsroutine_           TAOSETCONSTRAINTSROUTINE
2010978b7dSBarry Smith   #define taomonitorset_                      TAOMONITORSET
21ae93cb3cSJason Sarich   #define taogetconvergencehistory_           TAOGETCONVERGENCEHISTORY
22a7e14dcfSSatish Balay   #define taosetconvergencetest_              TAOSETCONVERGENCETEST
23e1e80dc8SAlp Dener   #define taosetupdate_                       TAOSETUPDATE
24a7e14dcfSSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
25a82e8c82SStefano Zampini   #define taosetobjective_                    taosetobjective
26a82e8c82SStefano Zampini   #define taosetgradient_                     taosetgradient
27a82e8c82SStefano Zampini   #define taosetobjectiveandgradient_         taosetobjectiveandgradient
28a82e8c82SStefano Zampini   #define taosethessian_                      taosethessian
29737f463aSAlp Dener   #define taosetresidualroutine_              taosetresidualroutine
304ffbe8acSAlp Dener   #define taosetjacobianresidualroutine_      taosetjacobianresidualroutine
31a7e14dcfSSatish Balay   #define taosetjacobianroutine_              taosetjacobianroutine
32a7e14dcfSSatish Balay   #define taosetjacobianstateroutine_         taosetjacobianstateroutine
33a7e14dcfSSatish Balay   #define taosetjacobiandesignroutine_        taosetjacobiandesignroutine
34a7e14dcfSSatish Balay   #define taosetjacobianinequalityroutine_    taosetjacobianinequalityroutine
35a7e14dcfSSatish Balay   #define taosetjacobianequalityroutine_      taosetjacobianequalityroutine
36a7e14dcfSSatish Balay   #define taosetinequalityconstraintsroutine_ taosetinequalityconstraintsroutine
37a7e14dcfSSatish Balay   #define taosetequalityconstraintsroutine_   taosetequalityconstraintsroutine
38a7e14dcfSSatish Balay   #define taosetvariableboundsroutine_        taosetvariableboundsroutine
39a7e14dcfSSatish Balay   #define taosetconstraintsroutine_           taosetconstraintsroutine
4010978b7dSBarry Smith   #define taomonitorset_                      taomonitorset
41ae93cb3cSJason Sarich   #define taogetconvergencehistory_           taogetconvergencehistory
42a7e14dcfSSatish Balay   #define taosetconvergencetest_              taosetconvergencetest
43e1e80dc8SAlp Dener   #define taosetupdate_                       taosetupdate
44a7e14dcfSSatish Balay #endif
45a7e14dcfSSatish Balay 
468732526dSAlp Dener static struct {
478732526dSAlp Dener   PetscFortranCallbackId obj;
488732526dSAlp Dener   PetscFortranCallbackId grad;
498732526dSAlp Dener   PetscFortranCallbackId objgrad;
508732526dSAlp Dener   PetscFortranCallbackId hess;
51737f463aSAlp Dener   PetscFortranCallbackId lsres;
52737f463aSAlp Dener   PetscFortranCallbackId lsjac;
538732526dSAlp Dener   PetscFortranCallbackId jac;
548732526dSAlp Dener   PetscFortranCallbackId jacstate;
558732526dSAlp Dener   PetscFortranCallbackId jacdesign;
568732526dSAlp Dener   PetscFortranCallbackId bounds;
578732526dSAlp Dener   PetscFortranCallbackId mon;
588732526dSAlp Dener   PetscFortranCallbackId mondestroy;
598732526dSAlp Dener   PetscFortranCallbackId convtest;
608732526dSAlp Dener   PetscFortranCallbackId constraints;
618732526dSAlp Dener   PetscFortranCallbackId jacineq;
628732526dSAlp Dener   PetscFortranCallbackId jaceq;
638732526dSAlp Dener   PetscFortranCallbackId conineq;
648732526dSAlp Dener   PetscFortranCallbackId coneq;
658732526dSAlp Dener   PetscFortranCallbackId nfuncs;
66e1e80dc8SAlp Dener   PetscFortranCallbackId update;
678732526dSAlp Dener #if defined(PETSC_HAVE_F90_2PTR_ARG)
688732526dSAlp Dener   PetscFortranCallbackId function_pgiptr;
698732526dSAlp Dener #endif
708732526dSAlp Dener } _cb;
71a7e14dcfSSatish Balay 
ourtaoobjectiveroutine(Tao tao,Vec x,PetscReal * f,PetscCtx ctx)72*2a8381b2SBarry Smith static PetscErrorCode ourtaoobjectiveroutine(Tao tao, Vec x, PetscReal *f, PetscCtx ctx)
73a7e14dcfSSatish Balay {
748732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.obj, (Tao *, Vec *, PetscReal *, void *, PetscErrorCode *), (&tao, &x, f, _ctx, &ierr));
75a7e14dcfSSatish Balay }
76a7e14dcfSSatish Balay 
ourtaogradientroutine(Tao tao,Vec x,Vec g,PetscCtx ctx)77*2a8381b2SBarry Smith static PetscErrorCode ourtaogradientroutine(Tao tao, Vec x, Vec g, PetscCtx ctx)
78a7e14dcfSSatish Balay {
798732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.grad, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &g, _ctx, &ierr));
80a7e14dcfSSatish Balay }
81a7e14dcfSSatish Balay 
ourtaoobjectiveandgradientroutine(Tao tao,Vec x,PetscReal * f,Vec g,PetscCtx ctx)82*2a8381b2SBarry Smith static PetscErrorCode ourtaoobjectiveandgradientroutine(Tao tao, Vec x, PetscReal *f, Vec g, PetscCtx ctx)
83a7e14dcfSSatish Balay {
848732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.objgrad, (Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), (&tao, &x, f, &g, _ctx, &ierr));
85a7e14dcfSSatish Balay }
86a7e14dcfSSatish Balay 
ourtaohessianroutine(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx)87*2a8381b2SBarry Smith static PetscErrorCode ourtaohessianroutine(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx)
88a7e14dcfSSatish Balay {
898732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.hess, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, _ctx, &ierr));
90a7e14dcfSSatish Balay }
91a7e14dcfSSatish Balay 
ourtaojacobianroutine(Tao tao,Vec x,Mat H,Mat Hpre,PetscCtx ctx)92*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobianroutine(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx)
93a7e14dcfSSatish Balay {
948732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.jac, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, _ctx, &ierr));
95a7e14dcfSSatish Balay }
96a7e14dcfSSatish Balay 
ourtaojacobianstateroutine(Tao tao,Vec x,Mat H,Mat Hpre,Mat Hinv,PetscCtx ctx)97*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobianstateroutine(Tao tao, Vec x, Mat H, Mat Hpre, Mat Hinv, PetscCtx ctx)
98a7e14dcfSSatish Balay {
998732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.jacstate, (Tao *, Vec *, Mat *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, &Hpre, &Hinv, _ctx, &ierr));
100a7e14dcfSSatish Balay }
101a7e14dcfSSatish Balay 
ourtaojacobiandesignroutine(Tao tao,Vec x,Mat H,PetscCtx ctx)102*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobiandesignroutine(Tao tao, Vec x, Mat H, PetscCtx ctx)
103a7e14dcfSSatish Balay {
1048732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.jacdesign, (Tao *, Vec *, Mat *, void *, PetscErrorCode *), (&tao, &x, &H, _ctx, &ierr));
105a7e14dcfSSatish Balay }
106a7e14dcfSSatish Balay 
ourtaoboundsroutine(Tao tao,Vec xl,Vec xu,PetscCtx ctx)107*2a8381b2SBarry Smith static PetscErrorCode ourtaoboundsroutine(Tao tao, Vec xl, Vec xu, PetscCtx ctx)
108a7e14dcfSSatish Balay {
1098732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.bounds, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &xl, &xu, _ctx, &ierr));
110a7e14dcfSSatish Balay }
ourtaoresidualroutine(Tao tao,Vec x,Vec f,PetscCtx ctx)111*2a8381b2SBarry Smith static PetscErrorCode ourtaoresidualroutine(Tao tao, Vec x, Vec f, PetscCtx ctx)
112a7e14dcfSSatish Balay {
113737f463aSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.lsres, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &f, _ctx, &ierr));
114737f463aSAlp Dener }
115737f463aSAlp Dener 
ourtaojacobianresidualroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)116*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobianresidualroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
117737f463aSAlp Dener {
118737f463aSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.lsjac, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
119a7e14dcfSSatish Balay }
120a7e14dcfSSatish Balay 
ourtaomonitor(Tao tao,PetscCtx ctx)121*2a8381b2SBarry Smith static PetscErrorCode ourtaomonitor(Tao tao, PetscCtx ctx)
122a7e14dcfSSatish Balay {
1238732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.mon, (Tao *, void *, PetscErrorCode *), (&tao, _ctx, &ierr));
124a7e14dcfSSatish Balay }
125a7e14dcfSSatish Balay 
ourtaomondestroy(PetscCtxRt ctx)126*2a8381b2SBarry Smith static PetscErrorCode ourtaomondestroy(PetscCtxRt ctx)
127a7e14dcfSSatish Balay {
128*2a8381b2SBarry Smith   Tao tao = *(Tao *)ctx;
1298732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.mondestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
130a7e14dcfSSatish Balay }
ourtaoconvergencetest(Tao tao,PetscCtx ctx)131*2a8381b2SBarry Smith static PetscErrorCode ourtaoconvergencetest(Tao tao, PetscCtx ctx)
132a7e14dcfSSatish Balay {
1338732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.convtest, (Tao *, void *, PetscErrorCode *), (&tao, _ctx, &ierr));
134a7e14dcfSSatish Balay }
135a7e14dcfSSatish Balay 
ourtaoconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)136*2a8381b2SBarry Smith static PetscErrorCode ourtaoconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
137a7e14dcfSSatish Balay {
1388732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.constraints, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
139a7e14dcfSSatish Balay }
140a7e14dcfSSatish Balay 
ourtaojacobianinequalityroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)141*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobianinequalityroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
142a7e14dcfSSatish Balay {
1438732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.jacineq, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
144a7e14dcfSSatish Balay }
145a7e14dcfSSatish Balay 
ourtaojacobianequalityroutine(Tao tao,Vec x,Mat J,Mat Jpre,PetscCtx ctx)146*2a8381b2SBarry Smith static PetscErrorCode ourtaojacobianequalityroutine(Tao tao, Vec x, Mat J, Mat Jpre, PetscCtx ctx)
147a7e14dcfSSatish Balay {
1488732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.jaceq, (Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), (&tao, &x, &J, &Jpre, _ctx, &ierr));
149a7e14dcfSSatish Balay }
150a7e14dcfSSatish Balay 
ourtaoinequalityconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)151*2a8381b2SBarry Smith static PetscErrorCode ourtaoinequalityconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
152a7e14dcfSSatish Balay {
1538732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.conineq, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
154a7e14dcfSSatish Balay }
155a7e14dcfSSatish Balay 
ourtaoequalityconstraintsroutine(Tao tao,Vec x,Vec c,PetscCtx ctx)156*2a8381b2SBarry Smith static PetscErrorCode ourtaoequalityconstraintsroutine(Tao tao, Vec x, Vec c, PetscCtx ctx)
157a7e14dcfSSatish Balay {
1588732526dSAlp Dener   PetscObjectUseFortranCallback(tao, _cb.coneq, (Tao *, Vec *, Vec *, void *, PetscErrorCode *), (&tao, &x, &c, _ctx, &ierr));
159a7e14dcfSSatish Balay }
160a7e14dcfSSatish Balay 
ourtaoupdateroutine(Tao tao,PetscInt iter,PetscCtx ctx)161*2a8381b2SBarry Smith static PetscErrorCode ourtaoupdateroutine(Tao tao, PetscInt iter, PetscCtx ctx)
162e1e80dc8SAlp Dener {
1638fcddce6SStefano Zampini   PetscObjectUseFortranCallback(tao, _cb.update, (Tao *, PetscInt *, void *), (&tao, &iter, _ctx));
164e1e80dc8SAlp Dener }
165e1e80dc8SAlp Dener 
taosetobjective_(Tao * tao,void (* func)(Tao *,Vec *,PetscReal *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)166*2a8381b2SBarry Smith PETSC_EXTERN void taosetobjective_(Tao *tao, void (*func)(Tao *, Vec *, PetscReal *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
167a7e14dcfSSatish Balay {
1688732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
1695ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.obj, (PetscFortranCallbackFn *)func, ctx);
170a82e8c82SStefano Zampini   if (!*ierr) *ierr = TaoSetObjective(*tao, ourtaoobjectiveroutine, ctx);
171a7e14dcfSSatish Balay }
172a7e14dcfSSatish Balay 
taosetgradient_(Tao * tao,Vec * g,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)173*2a8381b2SBarry Smith PETSC_EXTERN void taosetgradient_(Tao *tao, Vec *g, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
174a7e14dcfSSatish Balay {
1758732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
1765ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.grad, (PetscFortranCallbackFn *)func, ctx);
177a82e8c82SStefano Zampini   if (!*ierr) *ierr = TaoSetGradient(*tao, *g, ourtaogradientroutine, ctx);
178a7e14dcfSSatish Balay }
179a7e14dcfSSatish Balay 
taosetobjectiveandgradient_(Tao * tao,Vec * g,void (* func)(Tao *,Vec *,PetscReal *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)180*2a8381b2SBarry Smith PETSC_EXTERN void taosetobjectiveandgradient_(Tao *tao, Vec *g, void (*func)(Tao *, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
181a7e14dcfSSatish Balay {
1828732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
1835ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.objgrad, (PetscFortranCallbackFn *)func, ctx);
184a82e8c82SStefano Zampini   if (!*ierr) *ierr = TaoSetObjectiveAndGradient(*tao, *g, ourtaoobjectiveandgradientroutine, ctx);
185a82e8c82SStefano Zampini }
186a82e8c82SStefano Zampini 
taosethessian_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)187*2a8381b2SBarry Smith PETSC_EXTERN void taosethessian_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
188a82e8c82SStefano Zampini {
189a82e8c82SStefano Zampini   CHKFORTRANNULLFUNCTION(func);
1905ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.hess, (PetscFortranCallbackFn *)func, ctx);
191a82e8c82SStefano Zampini   if (!*ierr) *ierr = TaoSetHessian(*tao, *J, *Jp, ourtaohessianroutine, ctx);
192a7e14dcfSSatish Balay }
193a7e14dcfSSatish Balay 
taosetresidualroutine_(Tao * tao,Vec * F,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)194*2a8381b2SBarry Smith PETSC_EXTERN void taosetresidualroutine_(Tao *tao, Vec *F, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
195a7e14dcfSSatish Balay {
1968732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
1975ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.lsres, (PetscFortranCallbackFn *)func, ctx);
198737f463aSAlp Dener   if (!*ierr) *ierr = TaoSetResidualRoutine(*tao, *F, ourtaoresidualroutine, ctx);
199737f463aSAlp Dener }
200737f463aSAlp Dener 
taosetjacobianresidualroutine_(Tao * tao,Mat * J,Mat * Jpre,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)201*2a8381b2SBarry Smith PETSC_EXTERN void taosetjacobianresidualroutine_(Tao *tao, Mat *J, Mat *Jpre, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
202737f463aSAlp Dener {
203737f463aSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2045ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.lsjac, (PetscFortranCallbackFn *)func, ctx);
2054ffbe8acSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianResidualRoutine(*tao, *J, *Jpre, ourtaojacobianresidualroutine, ctx);
206a7e14dcfSSatish Balay }
207a7e14dcfSSatish Balay 
taosetjacobianroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)208*2a8381b2SBarry Smith PETSC_EXTERN void taosetjacobianroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
209a7e14dcfSSatish Balay {
2108732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2115ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jac, (PetscFortranCallbackFn *)func, ctx);
2128732526dSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianRoutine(*tao, *J, *Jp, ourtaojacobianroutine, ctx);
213a7e14dcfSSatish Balay }
214a7e14dcfSSatish Balay 
taosetjacobianstateroutine_(Tao * tao,Mat * J,Mat * Jp,Mat * Jinv,void (* func)(Tao *,Vec *,Mat *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)215*2a8381b2SBarry Smith 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)
216a7e14dcfSSatish Balay {
2178732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2185ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacstate, (PetscFortranCallbackFn *)func, ctx);
2198732526dSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianStateRoutine(*tao, *J, *Jp, *Jinv, ourtaojacobianstateroutine, ctx);
220a7e14dcfSSatish Balay }
221a7e14dcfSSatish Balay 
taosetjacobiandesignroutine_(Tao * tao,Mat * J,void (* func)(Tao *,Vec *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)222*2a8381b2SBarry Smith PETSC_EXTERN void taosetjacobiandesignroutine_(Tao *tao, Mat *J, void (*func)(Tao *, Vec *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
223a7e14dcfSSatish Balay {
2248732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2255ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacdesign, (PetscFortranCallbackFn *)func, ctx);
2268732526dSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianDesignRoutine(*tao, *J, ourtaojacobiandesignroutine, ctx);
227a7e14dcfSSatish Balay }
228a7e14dcfSSatish Balay 
taosetvariableboundsroutine_(Tao * tao,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)229*2a8381b2SBarry Smith PETSC_EXTERN void taosetvariableboundsroutine_(Tao *tao, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
230a7e14dcfSSatish Balay {
2318732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2325ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.bounds, (PetscFortranCallbackFn *)func, ctx);
2338732526dSAlp Dener   if (!*ierr) *ierr = TaoSetVariableBoundsRoutine(*tao, ourtaoboundsroutine, ctx);
234a7e14dcfSSatish Balay }
235a7e14dcfSSatish Balay 
taomonitorset_(Tao * tao,void (* func)(Tao *,void *,PetscErrorCode *),PetscCtx ctx,void (* mondestroy)(void *,PetscErrorCode *),PetscErrorCode * ierr)236*2a8381b2SBarry Smith PETSC_EXTERN void taomonitorset_(Tao *tao, void (*func)(Tao *, void *, PetscErrorCode *), PetscCtx ctx, void (*mondestroy)(void *, PetscErrorCode *), PetscErrorCode *ierr)
237a7e14dcfSSatish Balay {
238aecf964fSBarry Smith   CHKFORTRANNULLFUNCTION(mondestroy);
2395ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mon, (PetscFortranCallbackFn *)func, ctx);
2405975b3b6SBarry Smith   if (*ierr) return;
2415ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.mondestroy, (PetscFortranCallbackFn *)mondestroy, ctx);
2425975b3b6SBarry Smith   if (*ierr) return;
24310978b7dSBarry Smith   *ierr = TaoMonitorSet(*tao, ourtaomonitor, *tao, ourtaomondestroy);
244a7e14dcfSSatish Balay }
245a7e14dcfSSatish Balay 
taosetconvergencetest_(Tao * tao,void (* func)(Tao *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)246*2a8381b2SBarry Smith PETSC_EXTERN void taosetconvergencetest_(Tao *tao, void (*func)(Tao *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
247a7e14dcfSSatish Balay {
2488732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2495ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.convtest, (PetscFortranCallbackFn *)func, ctx);
2508732526dSAlp Dener   if (!*ierr) *ierr = TaoSetConvergenceTest(*tao, ourtaoconvergencetest, ctx);
251a7e14dcfSSatish Balay }
252a7e14dcfSSatish Balay 
taosetconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)253*2a8381b2SBarry Smith PETSC_EXTERN void taosetconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
254a7e14dcfSSatish Balay {
2558732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2565ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.constraints, (PetscFortranCallbackFn *)func, ctx);
2578732526dSAlp Dener   if (!*ierr) *ierr = TaoSetConstraintsRoutine(*tao, *C, ourtaoconstraintsroutine, ctx);
258a7e14dcfSSatish Balay }
259a7e14dcfSSatish Balay 
taogetconvergencehistory_(Tao * tao,PetscInt * nhist,PetscErrorCode * ierr)26019caf8f3SSatish Balay PETSC_EXTERN void taogetconvergencehistory_(Tao *tao, PetscInt *nhist, PetscErrorCode *ierr)
261a7e14dcfSSatish Balay {
262ae93cb3cSJason Sarich   *ierr = TaoGetConvergenceHistory(*tao, NULL, NULL, NULL, NULL, nhist);
263a7e14dcfSSatish Balay }
264a7e14dcfSSatish Balay 
taosetjacobianinequalityroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)265*2a8381b2SBarry Smith PETSC_EXTERN void taosetjacobianinequalityroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
266a7e14dcfSSatish Balay {
2678732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2685ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jacineq, (PetscFortranCallbackFn *)func, ctx);
2698732526dSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianInequalityRoutine(*tao, *J, *Jp, ourtaojacobianinequalityroutine, ctx);
270a7e14dcfSSatish Balay }
271a7e14dcfSSatish Balay 
taosetjacobianequalityroutine_(Tao * tao,Mat * J,Mat * Jp,void (* func)(Tao *,Vec *,Mat *,Mat *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)272*2a8381b2SBarry Smith PETSC_EXTERN void taosetjacobianequalityroutine_(Tao *tao, Mat *J, Mat *Jp, void (*func)(Tao *, Vec *, Mat *, Mat *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
273a7e14dcfSSatish Balay {
2748732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2755ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.jaceq, (PetscFortranCallbackFn *)func, ctx);
2768732526dSAlp Dener   if (!*ierr) *ierr = TaoSetJacobianEqualityRoutine(*tao, *J, *Jp, ourtaojacobianequalityroutine, ctx);
277a7e14dcfSSatish Balay }
278a7e14dcfSSatish Balay 
taosetinequalityconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)279*2a8381b2SBarry Smith PETSC_EXTERN void taosetinequalityconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
280a7e14dcfSSatish Balay {
2818732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2825ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.conineq, (PetscFortranCallbackFn *)func, ctx);
2838732526dSAlp Dener   if (!*ierr) *ierr = TaoSetInequalityConstraintsRoutine(*tao, *C, ourtaoinequalityconstraintsroutine, ctx);
284a7e14dcfSSatish Balay }
285a7e14dcfSSatish Balay 
taosetequalityconstraintsroutine_(Tao * tao,Vec * C,void (* func)(Tao *,Vec *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)286*2a8381b2SBarry Smith PETSC_EXTERN void taosetequalityconstraintsroutine_(Tao *tao, Vec *C, void (*func)(Tao *, Vec *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
287a7e14dcfSSatish Balay {
2888732526dSAlp Dener   CHKFORTRANNULLFUNCTION(func);
2895ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.coneq, (PetscFortranCallbackFn *)func, ctx);
2908732526dSAlp Dener   if (!*ierr) *ierr = TaoSetEqualityConstraintsRoutine(*tao, *C, ourtaoequalityconstraintsroutine, ctx);
291a7e14dcfSSatish Balay }
292a7e14dcfSSatish Balay 
taosetupdate_(Tao * tao,void (* func)(Tao *,PetscInt *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)293*2a8381b2SBarry Smith PETSC_EXTERN void taosetupdate_(Tao *tao, void (*func)(Tao *, PetscInt *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
294e1e80dc8SAlp Dener {
295e1e80dc8SAlp Dener   CHKFORTRANNULLFUNCTION(func);
2965ebfa9e9SBarry Smith   *ierr = PetscObjectSetFortranCallback((PetscObject)*tao, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.update, (PetscFortranCallbackFn *)func, ctx);
297e1e80dc8SAlp Dener   if (!*ierr) *ierr = TaoSetUpdate(*tao, ourtaoupdateroutine, ctx);
298e1e80dc8SAlp Dener }
299