xref: /petsc/src/tao/interface/ftn-custom/ztaosolverf.c (revision b0a7d7e7f246badab30cb8ce3f95dd1540bfb513)
1 #include "petsc-private/fortranimpl.h"
2 #include "tao-private/taosolver_impl.h"
3 
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6 #define taosetobjectiveroutine_             TAOSETOBJECTIVEROUTINE
7 #define taosetgradientroutine_              TAOSETGRADIENTROUTINE
8 #define taosetobjectiveandgradientroutine_  TAOSETOBJECTIVEANDGRADIENTROUTINE
9 #define taosethessianroutine_               TAOSETHESSIANROUTINE
10 #define taosetseparableobjectiveroutine_    TAOSETSEPARABLEOBJECTIVEROUTINE
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 taosetmonitor_                      TAOSETMONITOR
21 #define taosettype_                         TAOSETTYPE
22 #define taoview_                            TAOVIEW
23 #define taogethistory_                      TAOGETHISTORY
24 #define taosetconvergencetest_              TAOSETCONVERGENCETEST
25 #define taogetoptionsprefix_                TAOGETOPTIONSPREFIX
26 #define taosetoptionsprefix_                TAOSETOPTIONSPREFIX
27 #define taoappendoptionsprefix_             TAOAPPENDOPTIONSPREFIX
28 #define taogettype_                         TAOGETTYPE
29 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
30 
31 #define taosetobjectiveroutine_             taosetobjectiveroutine
32 #define taosetgradientroutine_              taosetgradientroutine
33 #define taosetobjectiveandgradientroutine_  taosetobjectiveandgradientroutine
34 #define taosethessianroutine_               taosethessianroutine
35 #define taosetseparableobjectiveroutine_    taosetseparableobjectiveroutine
36 #define taosetjacobianroutine_              taosetjacobianroutine
37 #define taosetjacobianstateroutine_         taosetjacobianstateroutine
38 #define taosetjacobiandesignroutine_        taosetjacobiandesignroutine
39 #define taosetjacobianinequalityroutine_    taosetjacobianinequalityroutine
40 #define taosetjacobianequalityroutine_      taosetjacobianequalityroutine
41 #define taosetinequalityconstraintsroutine_ taosetinequalityconstraintsroutine
42 #define taosetequalityconstraintsroutine_   taosetequalityconstraintsroutine
43 #define taosetvariableboundsroutine_        taosetvariableboundsroutine
44 #define taosetconstraintsroutine_           taosetconstraintsroutine
45 #define taosetmonitor_                      taosetmonitor
46 #define taosettype_                         taosettype
47 #define taoview_                            taoview
48 #define taogethistory_                      taogethistory
49 #define taosetconvergencetest_              taosetconvergencetest
50 #define taogetoptionsprefix_                taogetoptionsprefix
51 #define taosetoptionsprefix_                taosetoptionsprefix
52 #define taoappendoptionsprefix_             taoappendoptionsprefix
53 #define taogettype_                         taogettype
54 #endif
55 
56 static int OBJ=0;       /*  objective routine index */
57 static int GRAD=1;      /*  gradient routine index */
58 static int OBJGRAD=2;   /*  objective and gradient routine */
59 static int HESS=3;      /*  hessian routine index */
60 static int SEPOBJ=4;    /*  separable objective routine index */
61 static int JAC=5;       /*  jacobian routine index */
62 static int JACSTATE=6;  /*  jacobian state routine index */
63 static int JACDESIGN=7; /*  jacobian design routine index */
64 static int BOUNDS=8;
65 static int MON=9;       /*  monitor routine index */
66 static int MONCTX=10;       /*  monitor routine index */
67 static int MONDESTROY=11; /*  monitor destroy index */
68 static int CONVTEST=12;  /*  */
69 static int CONSTRAINTS=13;
70 static int JACINEQ=14;
71 static int JACEQ=15;
72 static int CONINEQ=16;
73 static int CONEQ=17;
74 static int NFUNCS=18;
75 
76 static PetscErrorCode ourtaoobjectiveroutine(TaoSolver tao, Vec x, PetscReal *f, void *ctx)
77 {
78     PetscErrorCode ierr = 0;
79     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,PetscReal*,void*,PetscErrorCode*))
80 	(((PetscObject)tao)->fortran_func_pointers[OBJ]))(&tao,&x,f,ctx,&ierr);
81     CHKERRQ(ierr);
82     return 0;
83 }
84 
85 static PetscErrorCode ourtaogradientroutine(TaoSolver tao, Vec x, Vec g, void *ctx)
86 {
87     PetscErrorCode ierr = 0;
88     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
89        (((PetscObject)tao)->fortran_func_pointers[GRAD]))(&tao,&x,&g,ctx,&ierr);
90     CHKERRQ(ierr);
91     return 0;
92 
93 }
94 
95 static PetscErrorCode ourtaoobjectiveandgradientroutine(TaoSolver tao, Vec x, PetscReal *f, Vec g, void* ctx)
96 {
97     PetscErrorCode ierr = 0;
98     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,PetscReal*,Vec*,void*,PetscErrorCode*))
99      (((PetscObject)tao)->fortran_func_pointers[OBJGRAD]))(&tao,&x,f,&g,ctx,&ierr);
100     CHKERRQ(ierr);
101     return 0;
102 }
103 
104 static PetscErrorCode ourtaohessianroutine(TaoSolver tao, Vec x, Mat *H, Mat *Hpre, MatStructure *type, void *ctx)
105 {
106     PetscErrorCode ierr = 0;
107     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))
108      (((PetscObject)tao)->fortran_func_pointers[HESS]))(&tao,&x,H,Hpre,type,ctx,&ierr); CHKERRQ(ierr);
109     return 0;
110 }
111 
112 static PetscErrorCode ourtaojacobianroutine(TaoSolver tao, Vec x, Mat *H, Mat *Hpre, MatStructure *type, void *ctx)
113 {
114     PetscErrorCode ierr = 0;
115     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))
116      (((PetscObject)tao)->fortran_func_pointers[JAC]))(&tao,&x,H,Hpre,type,ctx,&ierr); CHKERRQ(ierr);
117     return 0;
118 }
119 
120 static PetscErrorCode ourtaojacobianstateroutine(TaoSolver tao, Vec x, Mat *H, Mat *Hpre, Mat *Hinv, MatStructure *type, void *ctx)
121 {
122     PetscErrorCode ierr = 0;
123     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))
124      (((PetscObject)tao)->fortran_func_pointers[JACSTATE]))(&tao,&x,H,Hpre,Hinv,type,ctx,&ierr); CHKERRQ(ierr);
125     return 0;
126 }
127 
128 static PetscErrorCode ourtaojacobiandesignroutine(TaoSolver tao, Vec x, Mat *H, void *ctx)
129 {
130     PetscErrorCode ierr = 0;
131     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,void*,PetscErrorCode*))
132      (((PetscObject)tao)->fortran_func_pointers[JACDESIGN]))(&tao,&x,H,ctx,&ierr); CHKERRQ(ierr);
133     return 0;
134 }
135 
136 static PetscErrorCode ourtaoboundsroutine(TaoSolver tao, Vec xl, Vec xu, void *ctx)
137 {
138     PetscErrorCode ierr = 0;
139     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
140      (((PetscObject)tao)->fortran_func_pointers[BOUNDS]))(&tao,&xl,&xu,ctx,&ierr); CHKERRQ(ierr);
141     return 0;
142 }
143 static PetscErrorCode ourtaoseparableobjectiveroutine(TaoSolver tao, Vec x, Vec f, void *ctx)
144 {
145     PetscErrorCode ierr = 0;
146     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
147      (((PetscObject)tao)->fortran_func_pointers[SEPOBJ]))(&tao,&x,&f,ctx,&ierr);
148     return 0;
149 }
150 
151 static PetscErrorCode ourtaomonitor(TaoSolver tao, void *ctx)
152 {
153     PetscErrorCode ierr = 0;
154     (*(void (PETSC_STDCALL *)(TaoSolver *, void*, PetscErrorCode*))
155      (((PetscObject)tao)->fortran_func_pointers[MON]))(&tao,ctx,&ierr);
156     CHKERRQ(ierr);
157     return 0;
158 }
159 
160 static PetscErrorCode ourtaomondestroy(void **ctx)
161 {
162     PetscErrorCode ierr = 0;
163     TaoSolver tao = *(TaoSolver*)ctx;
164     void *mctx = (void*)((PetscObject)tao)->fortran_func_pointers[MONCTX];
165     (*(void (PETSC_STDCALL *)(void*,PetscErrorCode*))(((PetscObject)tao)->fortran_func_pointers[MONDESTROY]))(mctx,&ierr); CHKERRQ(ierr);
166     return 0;
167 }
168 static PetscErrorCode ourtaoconvergencetest(TaoSolver tao, void *ctx)
169 {
170     PetscErrorCode ierr = 0;
171     (*(void (PETSC_STDCALL *)(TaoSolver *, void*, PetscErrorCode*))
172      (((PetscObject)tao)->fortran_func_pointers[CONVTEST]))(&tao,ctx,&ierr);
173     CHKERRQ(ierr);
174     return 0;
175 }
176 
177 
178 static PetscErrorCode ourtaoconstraintsroutine(TaoSolver tao, Vec x, Vec c, void *ctx)
179 {
180     PetscErrorCode ierr = 0;
181     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
182        (((PetscObject)tao)->fortran_func_pointers[CONSTRAINTS]))(&tao,&x,&c,ctx,&ierr);
183     CHKERRQ(ierr);
184     return 0;
185 
186 }
187 
188 static PetscErrorCode ourtaojacobianinequalityroutine(TaoSolver tao, Vec x, Mat *J, Mat *Jpre, MatStructure *type, void *ctx)
189 {
190     PetscErrorCode ierr = 0;
191     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))
192      (((PetscObject)tao)->fortran_func_pointers[JACINEQ]))(&tao,&x,J,Jpre,type,ctx,&ierr); CHKERRQ(ierr);
193     return 0;
194 }
195 
196 static PetscErrorCode ourtaojacobianequalityroutine(TaoSolver tao, Vec x, Mat *J, Mat *Jpre, MatStructure *type, void *ctx)
197 {
198     PetscErrorCode ierr = 0;
199     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*))
200      (((PetscObject)tao)->fortran_func_pointers[JACEQ]))(&tao,&x,J,Jpre,type,ctx,&ierr); CHKERRQ(ierr);
201     return 0;
202 }
203 
204 static PetscErrorCode ourtaoinequalityconstraintsroutine(TaoSolver tao, Vec x, Vec c, void *ctx)
205 {
206     PetscErrorCode ierr = 0;
207     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
208        (((PetscObject)tao)->fortran_func_pointers[CONINEQ]))(&tao,&x,&c,ctx,&ierr);
209     CHKERRQ(ierr);
210     return 0;
211 
212 }
213 
214 static PetscErrorCode ourtaoequalityconstraintsroutine(TaoSolver tao, Vec x, Vec c, void *ctx)
215 {
216     PetscErrorCode ierr = 0;
217     (*(void (PETSC_STDCALL *)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*))
218        (((PetscObject)tao)->fortran_func_pointers[CONEQ]))(&tao,&x,&c,ctx,&ierr);
219     CHKERRQ(ierr);
220     return 0;
221 
222 }
223 
224 
225 EXTERN_C_BEGIN
226 
227 
228 void PETSC_STDCALL taosetobjectiveroutine_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
229 {
230     CHKFORTRANNULLOBJECT(ctx);
231     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
232     if (!func) {
233 	*ierr = TaoSetObjectiveRoutine(*tao,0,ctx);
234     } else {
235 	((PetscObject)*tao)->fortran_func_pointers[OBJ] = (PetscVoidFunction)func;
236 	*ierr = TaoSetObjectiveRoutine(*tao, ourtaoobjectiveroutine,ctx);
237     }
238 }
239 
240 void PETSC_STDCALL taosetgradientroutine_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
241 {
242     CHKFORTRANNULLOBJECT(ctx);
243     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
244     if (!func) {
245 	*ierr = TaoSetGradientRoutine(*tao,0,ctx);
246     } else {
247 	((PetscObject)*tao)->fortran_func_pointers[GRAD] = (PetscVoidFunction)func;
248 	*ierr = TaoSetGradientRoutine(*tao, ourtaogradientroutine,ctx);
249     }
250 }
251 
252 void PETSC_STDCALL taosetobjectiveandgradientroutine_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
253 {
254     CHKFORTRANNULLOBJECT(ctx);
255     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
256     if (!func) {
257 	*ierr = TaoSetObjectiveAndGradientRoutine(*tao,0,ctx);
258     } else {
259 	((PetscObject)*tao)->fortran_func_pointers[OBJGRAD] = (PetscVoidFunction)func;
260 	*ierr = TaoSetObjectiveAndGradientRoutine(*tao, ourtaoobjectiveandgradientroutine,ctx);
261     }
262 }
263 
264 
265 
266 
267 void PETSC_STDCALL taosetseparableobjectiveroutine_(TaoSolver *tao, Vec *F, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
268 {
269     CHKFORTRANNULLOBJECT(ctx);
270     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
271     if (!func) {
272 	*ierr = TaoSetSeparableObjectiveRoutine(*tao,*F,0,ctx);
273     } else {
274 	((PetscObject)*tao)->fortran_func_pointers[SEPOBJ] = (PetscVoidFunction)func;
275 	*ierr = TaoSetSeparableObjectiveRoutine(*tao,*F, ourtaoseparableobjectiveroutine,ctx);
276     }
277 }
278 
279 
280 
281 void PETSC_STDCALL taosetjacobianroutine_(TaoSolver *tao, Mat *J, Mat *Jp, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, Mat *, MatStructure *,void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
282 {
283     CHKFORTRANNULLOBJECT(ctx);
284     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
285     if (!func) {
286 	*ierr = TaoSetJacobianRoutine(*tao,*J,*Jp,0,ctx);
287     } else {
288 	((PetscObject)*tao)->fortran_func_pointers[JAC] = (PetscVoidFunction)func;
289 	*ierr = TaoSetJacobianRoutine(*tao,*J, *Jp, ourtaojacobianroutine,ctx);
290     }
291 }
292 
293 void PETSC_STDCALL taosetjacobianstateroutine_(TaoSolver *tao, Mat *J, Mat *Jp, Mat*Jinv, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, Mat *, Mat*, MatStructure *,void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
294 {
295     CHKFORTRANNULLOBJECT(ctx);
296     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
297     if (!func) {
298       *ierr = TaoSetJacobianStateRoutine(*tao,*J,*Jp,*Jinv,0,ctx);
299     } else {
300       ((PetscObject)*tao)->fortran_func_pointers[JACSTATE] = (PetscVoidFunction)func;
301       *ierr = TaoSetJacobianStateRoutine(*tao,*J, *Jp, *Jinv, ourtaojacobianstateroutine,ctx);
302     }
303 }
304 
305 void PETSC_STDCALL taosetjacobiandesignroutine_(TaoSolver *tao, Mat *J, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
306 {
307     CHKFORTRANNULLOBJECT(ctx);
308     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
309     if (!func) {
310 	*ierr = TaoSetJacobianDesignRoutine(*tao,*J,0,ctx);
311     } else {
312 	((PetscObject)*tao)->fortran_func_pointers[JACDESIGN] = (PetscVoidFunction)func;
313 	*ierr = TaoSetJacobianDesignRoutine(*tao,*J, ourtaojacobiandesignroutine,ctx);
314     }
315 }
316 
317 
318 void PETSC_STDCALL taosethessianroutine_(TaoSolver *tao, Mat *J, Mat *Jp, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, Mat *, MatStructure *,void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
319 {
320     CHKFORTRANNULLOBJECT(ctx);
321     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
322     if (!func) {
323 	*ierr = TaoSetHessianRoutine(*tao,*J,*Jp,0,ctx);
324     } else {
325 	((PetscObject)*tao)->fortran_func_pointers[HESS] = (PetscVoidFunction)func;
326 	*ierr = TaoSetHessianRoutine(*tao,*J, *Jp, ourtaohessianroutine,ctx);
327     }
328 }
329 
330 void PETSC_STDCALL taosetvariableboundsroutine_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx, PetscErrorCode *ierr)
331 {
332     CHKFORTRANNULLOBJECT(ctx);
333     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
334     if (func) {
335 	((PetscObject)*tao)->fortran_func_pointers[BOUNDS] = (PetscVoidFunction)func;
336 	*ierr = TaoSetVariableBoundsRoutine(*tao,ourtaoboundsroutine,ctx);
337     } else {
338 	*ierr = TaoSetVariableBoundsRoutine(*tao,0,ctx);
339     }
340 
341 }
342 void PETSC_STDCALL taosetmonitor_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*,void*,PetscErrorCode*),void *ctx, void (PETSC_STDCALL *mondestroy)(void*,PetscErrorCode*),PetscErrorCode *ierr)
343 {
344     CHKFORTRANNULLOBJECT(ctx);
345     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
346     if (func) {
347 	((PetscObject)*tao)->fortran_func_pointers[MON] = (PetscVoidFunction)func;
348 	if (FORTRANNULLFUNCTION(mondestroy)){
349 	  *ierr = TaoSetMonitor(*tao,ourtaomonitor,*tao,PETSC_NULL);
350 	} else {
351 	  *ierr = TaoSetMonitor(*tao,ourtaomonitor,*tao,ourtaomondestroy);
352 	}
353     }
354 }
355 
356 void PETSC_STDCALL taosetconvergencetest_(TaoSolver *tao, void (PETSC_STDCALL *func)(TaoSolver*,void*,PetscErrorCode*),void *ctx, PetscErrorCode *ierr)
357 {
358     CHKFORTRANNULLOBJECT(ctx);
359     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
360     if (!func) {
361 	*ierr = TaoSetConvergenceTest(*tao,0,ctx);
362     } else {
363 	((PetscObject)*tao)->fortran_func_pointers[CONVTEST] = (PetscVoidFunction)func;
364 	*ierr = TaoSetConvergenceTest(*tao,ourtaoconvergencetest,ctx);
365     }
366 }
367 
368 
369 void PETSC_STDCALL taosetconstraintsroutine_(TaoSolver *tao, Vec *C, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
370 {
371     CHKFORTRANNULLOBJECT(ctx);
372     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
373     if (!func) {
374       *ierr = TaoSetConstraintsRoutine(*tao,*C,0,ctx);
375     } else {
376 	((PetscObject)*tao)->fortran_func_pointers[CONSTRAINTS] = (PetscVoidFunction)func;
377 	*ierr = TaoSetConstraintsRoutine(*tao, *C, ourtaoconstraintsroutine,ctx);
378     }
379 }
380 
381 
382 void PETSC_STDCALL taosettype_(TaoSolver *tao, CHAR type_name PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
383 
384 {
385     char *t;
386 
387     FIXCHAR(type_name,len,t);
388     *ierr = TaoSetType(*tao,t);
389     FREECHAR(type_name,t);
390 
391 }
392 
393 void PETSC_STDCALL taoview_(TaoSolver *tao, PetscViewer *viewer, PetscErrorCode *ierr)
394 {
395     PetscViewer v;
396     PetscPatchDefaultViewers_Fortran(viewer,v);
397     *ierr = TaoView(*tao,v);
398 }
399 
400 void PETSC_STDCALL taogethistory_(TaoSolver *tao, PetscInt *nhist, PetscErrorCode *ierr)
401 {
402   *nhist  = (*tao)->hist_len;
403   *ierr = 0;
404 }
405 
406 void PETSC_STDCALL taogetoptionsprefix_(TaoSolver *tao, CHAR prefix PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
407 {
408   const char *name;
409   *ierr = TaoGetOptionsPrefix(*tao,&name);
410   *ierr = PetscStrncpy(prefix,name,len); if (*ierr) return;
411   FIXRETURNCHAR(PETSC_TRUE,prefix,len);
412 
413 }
414 
415 void PETSC_STDCALL taoappendoptionsprefix_(TaoSolver *tao, CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
416 {
417   char *name;
418   FIXCHAR(prefix,len,name);
419   *ierr = TaoAppendOptionsPrefix(*tao,name);
420   FREECHAR(prefix,name);
421 }
422 
423 void PETSC_STDCALL taosetoptionsprefix_(TaoSolver *tao, CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
424 {
425   char *t;
426   FIXCHAR(prefix,len,t);
427   *ierr = TaoSetOptionsPrefix(*tao,t);
428   FREECHAR(prefix,t);
429 }
430 
431 void PETSC_STDCALL taogettype_(TaoSolver *tao, CHAR name PETSC_MIXED_LEN(len), PetscErrorCode *ierr  PETSC_END_LEN(len))
432 {
433   const char *tname;
434   *ierr = TaoGetType(*tao,&tname);
435   *ierr = PetscStrncpy(name,tname,len); if (*ierr) return;
436   FIXRETURNCHAR(PETSC_TRUE,name,len);
437 
438 }
439 
440 
441 void PETSC_STDCALL taosetjacobianinequalityroutine_(TaoSolver *tao, Mat *J, Mat *Jp, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, Mat *, MatStructure *,void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
442 {
443     CHKFORTRANNULLOBJECT(ctx);
444     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
445     if (!func) {
446 	*ierr = TaoSetJacobianInequalityRoutine(*tao,*J,*Jp,0,ctx);
447     } else {
448 	((PetscObject)*tao)->fortran_func_pointers[JACINEQ] = (PetscVoidFunction)func;
449 	*ierr = TaoSetJacobianInequalityRoutine(*tao,*J, *Jp, ourtaojacobianinequalityroutine,ctx);
450     }
451 }
452 
453 void PETSC_STDCALL taosetjacobianequalityroutine_(TaoSolver *tao, Mat *J, Mat *Jp, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Mat *, Mat *, MatStructure *,void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
454 {
455     CHKFORTRANNULLOBJECT(ctx);
456     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
457     if (!func) {
458 	*ierr = TaoSetJacobianEqualityRoutine(*tao,*J,*Jp,0,ctx);
459     } else {
460 	((PetscObject)*tao)->fortran_func_pointers[JACEQ] = (PetscVoidFunction)func;
461 	*ierr = TaoSetJacobianEqualityRoutine(*tao,*J, *Jp, ourtaojacobianequalityroutine,ctx);
462     }
463 }
464 
465 
466 void PETSC_STDCALL taosetinequalityconstraintsroutine_(TaoSolver *tao, Vec *C, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
467 {
468     CHKFORTRANNULLOBJECT(ctx);
469     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
470     if (!func) {
471       *ierr = TaoSetInequalityConstraintsRoutine(*tao,*C,0,ctx);
472     } else {
473 	((PetscObject)*tao)->fortran_func_pointers[CONINEQ] = (PetscVoidFunction)func;
474 	*ierr = TaoSetInequalityConstraintsRoutine(*tao, *C, ourtaoinequalityconstraintsroutine,ctx);
475     }
476 }
477 
478 void PETSC_STDCALL taosetequalityconstraintsroutine_(TaoSolver *tao, Vec *C, void (PETSC_STDCALL *func)(TaoSolver*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr)
479 {
480     CHKFORTRANNULLOBJECT(ctx);
481     PetscObjectAllocateFortranPointers(*tao,NFUNCS);
482     if (!func) {
483       *ierr = TaoSetEqualityConstraintsRoutine(*tao,*C,0,ctx);
484     } else {
485 	((PetscObject)*tao)->fortran_func_pointers[CONEQ] = (PetscVoidFunction)func;
486 	*ierr = TaoSetEqualityConstraintsRoutine(*tao, *C, ourtaoequalityconstraintsroutine,ctx);
487     }
488 }
489 
490 
491 EXTERN_C_END
492 
493 
494