xref: /petsc/src/tao/interface/taosolver_hj.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2 
3 /*@C
4    TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix.
5 
6    Logically collective on Tao
7 
8    Input Parameters:
9 +  tao - the Tao context
10 .  H - Matrix used for the hessian
11 .  Hpre - Matrix that will be used operated on by preconditioner, can be same as H
12 .  hess - Hessian evaluation routine
13 -  ctx - [optional] user-defined context for private data for the
14          Hessian evaluation routine (may be NULL)
15 
16    Calling sequence of hess:
17 $    hess (Tao tao,Vec x,Mat H,Mat Hpre,void *ctx);
18 
19 +  tao - the Tao  context
20 .  x - input vector
21 .  H - Hessian matrix
22 .  Hpre - preconditioner matrix, usually the same as H
23 -  ctx - [optional] user-defined Hessian context
24 
25    Level: beginner
26 
27 @*/
28 PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
29 {
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
34   if (H) {
35     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
36     PetscCheckSameComm(tao,1,H,2);
37   }
38   if (Hpre) {
39     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
40     PetscCheckSameComm(tao,1,Hpre,3);
41   }
42   if (ctx) {
43     tao->user_hessP = ctx;
44   }
45   if (func) {
46     tao->ops->computehessian = func;
47   }
48   if (H) {
49     ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr);
50     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
51     tao->hessian = H;
52   }
53   if (Hpre) {
54     ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr);
55     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
56     tao->hessian_pre = Hpre;
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /*@C
62    TaoComputeHessian - Computes the Hessian matrix that has been
63    set with TaoSetHessianRoutine().
64 
65    Collective on Tao
66 
67    Input Parameters:
68 +  solver - the Tao solver context
69 -  xx - input vector
70 
71    Output Parameters:
72 +  H - Hessian matrix
73 -  Hpre - Preconditioning matrix
74 
75    Notes:
76    Most users should not need to explicitly call this routine, as it
77    is used internally within the minimization solvers.
78 
79    TaoComputeHessian() is typically used within minimization
80    implementations, so most users would not generally call this routine
81    themselves.
82 
83    Level: developer
84 
85 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessian()
86 
87 @*/
88 PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
89 {
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
94   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
95   PetscCheckSameComm(tao,1,X,2);
96 
97   if (!tao->ops->computehessian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessian() first");
98   ++tao->nhess;
99   ierr = PetscLogEventBegin(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
100   PetscStackPush("Tao user Hessian function");
101   ierr = (*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP);CHKERRQ(ierr);
102   PetscStackPop;
103   ierr = PetscLogEventEnd(Tao_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
104   PetscFunctionReturn(0);
105 }
106 
107 /*@C
108    TaoComputeJacobian - Computes the Jacobian matrix that has been
109    set with TaoSetJacobianRoutine().
110 
111    Collective on Tao
112 
113    Input Parameters:
114 +  solver - the Tao solver context
115 -  xx - input vector
116 
117    Output Parameters:
118 +  H - Jacobian matrix
119 -  Hpre - Preconditioning matrix
120 
121    Notes:
122    Most users should not need to explicitly call this routine, as it
123    is used internally within the minimization solvers.
124 
125    TaoComputeJacobian() is typically used within minimization
126    implementations, so most users would not generally call this routine
127    themselves.
128 
129    Level: developer
130 
131 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobian()
132 
133 @*/
134 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
135 {
136   PetscErrorCode ierr;
137 
138   PetscFunctionBegin;
139   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
140   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
141   PetscCheckSameComm(tao,1,X,2);
142 
143   if (!tao->ops->computejacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
144   ++tao->njac;
145   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
146   PetscStackPush("Tao user Jacobian function");
147   ierr = (*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP);CHKERRQ(ierr);
148   PetscStackPop;
149   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
150   PetscFunctionReturn(0);
151 }
152 
153 /*@C
154    TaoComputeJacobianState - Computes the Jacobian matrix that has been
155    set with TaoSetJacobianStateRoutine().
156 
157    Collective on Tao
158 
159    Input Parameters:
160 +  solver - the Tao solver context
161 -  xx - input vector
162 
163    Output Parameters:
164 +  H - Jacobian matrix
165 -  Hpre - Preconditioning matrix
166 
167    Notes:
168    Most users should not need to explicitly call this routine, as it
169    is used internally within the minimization solvers.
170 
171    TaoComputeJacobianState() is typically used within minimization
172    implementations, so most users would not generally call this routine
173    themselves.
174 
175    Level: developer
176 
177 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
178 
179 @*/
180 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
181 {
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
186   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
187   PetscCheckSameComm(tao,1,X,2);
188 
189   if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
190   ++tao->njac_state;
191   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
192   PetscStackPush("Tao user Jacobian(state) function");
193   ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);CHKERRQ(ierr);
194   PetscStackPop;
195   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
196   PetscFunctionReturn(0);
197 }
198 
199 /*@C
200    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
201    set with TaoSetJacobianDesignRoutine().
202 
203    Collective on Tao
204 
205    Input Parameters:
206 +  solver - the Tao solver context
207 -  xx - input vector
208 
209    Output Parameters:
210 .  H - Jacobian matrix
211 
212    Notes:
213    Most users should not need to explicitly call this routine, as it
214    is used internally within the minimization solvers.
215 
216    TaoComputeJacobianDesign() is typically used within minimization
217    implementations, so most users would not generally call this routine
218    themselves.
219 
220    Level: developer
221 
222 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
223 
224 @*/
225 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
226 {
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
231   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
232   PetscCheckSameComm(tao,1,X,2);
233 
234   if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
235   ++tao->njac_design;
236   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
237   PetscStackPush("Tao user Jacobian(design) function");
238   ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);CHKERRQ(ierr);
239   PetscStackPop;
240   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 /*@C
245    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
246 
247    Logically collective on Tao
248 
249    Input Parameters:
250 +  tao - the Tao context
251 .  J - Matrix used for the jacobian
252 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
253 .  jac - Jacobian evaluation routine
254 -  ctx - [optional] user-defined context for private data for the
255          Jacobian evaluation routine (may be NULL)
256 
257    Calling sequence of jac:
258 $    jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
259 
260 +  tao - the Tao  context
261 .  x - input vector
262 .  J - Jacobian matrix
263 .  Jpre - preconditioner matrix, usually the same as J
264 -  ctx - [optional] user-defined Jacobian context
265 
266    Level: intermediate
267 
268 @*/
269 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
270 {
271   PetscErrorCode ierr;
272   PetscFunctionBegin;
273   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
274   if (J) {
275     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
276     PetscCheckSameComm(tao,1,J,2);
277   }
278   if (Jpre) {
279     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
280     PetscCheckSameComm(tao,1,Jpre,3);
281   }
282   if (ctx) {
283     tao->user_jacP = ctx;
284   }
285   if (func) {
286     tao->ops->computejacobian = func;
287   }
288   if (J) {
289     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
290     ierr = MatDestroy(&tao->jacobian);CHKERRQ(ierr);
291     tao->jacobian = J;
292   }
293   if (Jpre) {
294     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
295     ierr = MatDestroy(&tao->jacobian_pre);CHKERRQ(ierr);
296     tao->jacobian_pre=Jpre;
297   }
298   PetscFunctionReturn(0);
299 }
300 
301 /*@C
302    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
303    (and its inverse) of the constraint function with respect to the state variables.
304    Used only for pde-constrained optimization.
305 
306    Logically collective on Tao
307 
308    Input Parameters:
309 +  tao - the Tao context
310 .  J - Matrix used for the jacobian
311 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is NULL
312 .  Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse.
313 .  jac - Jacobian evaluation routine
314 -  ctx - [optional] user-defined context for private data for the
315          Jacobian evaluation routine (may be NULL)
316 
317    Calling sequence of jac:
318 $    jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
319 
320 +  tao - the Tao  context
321 .  x - input vector
322 .  J - Jacobian matrix
323 .  Jpre - preconditioner matrix, usually the same as J
324 .  Jinv - inverse of J
325 -  ctx - [optional] user-defined Jacobian context
326 
327    Level: intermediate
328 .seealse: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
329 @*/
330 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat,void*), void *ctx)
331 {
332   PetscErrorCode ierr;
333   PetscFunctionBegin;
334   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
335   if (J) {
336     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
337     PetscCheckSameComm(tao,1,J,2);
338   }
339   if (Jpre) {
340     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
341     PetscCheckSameComm(tao,1,Jpre,3);
342   }
343   if (Jinv) {
344     PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
345     PetscCheckSameComm(tao,1,Jinv,4);
346   }
347   if (ctx) {
348     tao->user_jac_stateP = ctx;
349   }
350   if (func) {
351     tao->ops->computejacobianstate = func;
352   }
353   if (J) {
354     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
355     ierr = MatDestroy(&tao->jacobian_state);CHKERRQ(ierr);
356     tao->jacobian_state = J;
357   }
358   if (Jpre) {
359     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
360     ierr = MatDestroy(&tao->jacobian_state_pre);CHKERRQ(ierr);
361     tao->jacobian_state_pre=Jpre;
362   }
363   if (Jinv) {
364     ierr = PetscObjectReference((PetscObject)Jinv);CHKERRQ(ierr);
365     ierr = MatDestroy(&tao->jacobian_state_inv);CHKERRQ(ierr);
366     tao->jacobian_state_inv=Jinv;
367   }
368   PetscFunctionReturn(0);
369 }
370 
371 /*@C
372    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
373    the constraint function with respect to the design variables.  Used only for
374    pde-constrained optimization.
375 
376    Logically collective on Tao
377 
378    Input Parameters:
379 +  tao - the Tao context
380 .  J - Matrix used for the jacobian
381 .  jac - Jacobian evaluation routine
382 -  ctx - [optional] user-defined context for private data for the
383          Jacobian evaluation routine (may be NULL)
384 
385    Calling sequence of jac:
386 $    jac (Tao tao,Vec x,Mat *J,void *ctx);
387 
388 +  tao - the Tao  context
389 .  x - input vector
390 .  J - Jacobian matrix
391 -  ctx - [optional] user-defined Jacobian context
392 
393 
394    Notes:
395 
396    The function jac() takes Mat * as the matrix arguments rather than Mat.
397    This allows the Jacobian evaluation routine to replace A and/or B with a
398    completely new new matrix structure (not just different matrix elements)
399    when appropriate, for instance, if the nonzero structure is changing
400    throughout the global iterations.
401 
402    Level: intermediate
403 .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS()
404 @*/
405 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
406 {
407   PetscErrorCode ierr;
408 
409   PetscFunctionBegin;
410   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
411   if (J) {
412     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
413     PetscCheckSameComm(tao,1,J,2);
414   }
415   if (ctx) {
416     tao->user_jac_designP = ctx;
417   }
418   if (func) {
419     tao->ops->computejacobiandesign = func;
420   }
421   if (J) {
422     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
423     ierr = MatDestroy(&tao->jacobian_design);CHKERRQ(ierr);
424     tao->jacobian_design = J;
425   }
426   PetscFunctionReturn(0);
427 }
428 
429 /*@
430    TaoSetStateDesignIS - Indicate to the Tao which variables in the
431    solution vector are state variables and which are design.  Only applies to
432    pde-constrained optimization.
433 
434    Logically Collective on Tao
435 
436    Input Parameters:
437 +  tao - The Tao context
438 .  s_is - the index set corresponding to the state variables
439 -  d_is - the index set corresponding to the design variables
440 
441    Level: intermediate
442 
443 .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine()
444 @*/
445 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
446 {
447   PetscErrorCode ierr;
448 
449   PetscFunctionBegin;
450   ierr = PetscObjectReference((PetscObject)s_is);CHKERRQ(ierr);
451   ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr);
452   tao->state_is = s_is;
453   ierr = PetscObjectReference((PetscObject)(d_is));CHKERRQ(ierr);
454   ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr);
455   tao->design_is = d_is;
456   PetscFunctionReturn(0);
457 }
458 
459 /*@C
460    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
461    set with TaoSetJacobianEqualityRoutine().
462 
463    Collective on Tao
464 
465    Input Parameters:
466 +  solver - the Tao solver context
467 -  xx - input vector
468 
469    Output Parameters:
470 +  H - Jacobian matrix
471 -  Hpre - Preconditioning matrix
472 
473    Notes:
474    Most users should not need to explicitly call this routine, as it
475    is used internally within the minimization solvers.
476 
477    Level: developer
478 
479 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
480 
481 @*/
482 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
483 {
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
488   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
489   PetscCheckSameComm(tao,1,X,2);
490 
491   if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
492   ++tao->njac_equality;
493   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
494   PetscStackPush("Tao user Jacobian(equality) function");
495   ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);CHKERRQ(ierr);
496   PetscStackPop;
497   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 /*@C
502    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
503    set with TaoSetJacobianInequalityRoutine().
504 
505    Collective on Tao
506 
507    Input Parameters:
508 +  solver - the Tao solver context
509 -  xx - input vector
510 
511    Output Parameters:
512 +  H - Jacobian matrix
513 -  Hpre - Preconditioning matrix
514 
515    Notes:
516    Most users should not need to explicitly call this routine, as it
517    is used internally within the minimization solvers.
518 
519    Level: developer
520 
521 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
522 
523 @*/
524 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
525 {
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
530   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
531   PetscCheckSameComm(tao,1,X,2);
532 
533   if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
534   ++tao->njac_inequality;
535   ierr = PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
536   PetscStackPush("Tao user Jacobian(inequality) function");
537   ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);CHKERRQ(ierr);
538   PetscStackPop;
539   ierr = PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 /*@C
544    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
545    (and its inverse) of the constraint function with respect to the equality variables.
546    Used only for pde-constrained optimization.
547 
548    Logically collective on Tao
549 
550    Input Parameters:
551 +  tao - the Tao context
552 .  J - Matrix used for the jacobian
553 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
554 .  jac - Jacobian evaluation routine
555 -  ctx - [optional] user-defined context for private data for the
556          Jacobian evaluation routine (may be NULL)
557 
558    Calling sequence of jac:
559 $    jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
560 
561 +  tao - the Tao  context
562 .  x - input vector
563 .  J - Jacobian matrix
564 .  Jpre - preconditioner matrix, usually the same as J
565 -  ctx - [optional] user-defined Jacobian context
566 
567    Level: intermediate
568 .seealse: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS()
569 @*/
570 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
571 {
572   PetscErrorCode ierr;
573 
574   PetscFunctionBegin;
575   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
576   if (J) {
577     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
578     PetscCheckSameComm(tao,1,J,2);
579   }
580   if (Jpre) {
581     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
582     PetscCheckSameComm(tao,1,Jpre,3);
583   }
584   if (ctx) {
585     tao->user_jac_equalityP = ctx;
586   }
587   if (func) {
588     tao->ops->computejacobianequality = func;
589   }
590   if (J) {
591     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
592     ierr = MatDestroy(&tao->jacobian_equality);CHKERRQ(ierr);
593     tao->jacobian_equality = J;
594   }
595   if (Jpre) {
596     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
597     ierr = MatDestroy(&tao->jacobian_equality_pre);CHKERRQ(ierr);
598     tao->jacobian_equality_pre=Jpre;
599   }
600   PetscFunctionReturn(0);
601 }
602 
603 /*@C
604    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
605    (and its inverse) of the constraint function with respect to the inequality variables.
606    Used only for pde-constrained optimization.
607 
608    Logically collective on Tao
609 
610    Input Parameters:
611 +  tao - the Tao context
612 .  J - Matrix used for the jacobian
613 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
614 .  jac - Jacobian evaluation routine
615 -  ctx - [optional] user-defined context for private data for the
616          Jacobian evaluation routine (may be NULL)
617 
618    Calling sequence of jac:
619 $    jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
620 
621 +  tao - the Tao  context
622 .  x - input vector
623 .  J - Jacobian matrix
624 .  Jpre - preconditioner matrix, usually the same as J
625 -  ctx - [optional] user-defined Jacobian context
626 
627    Level: intermediate
628 .seealse: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS()
629 @*/
630 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
631 {
632   PetscErrorCode ierr;
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
635   if (J) {
636     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
637     PetscCheckSameComm(tao,1,J,2);
638   }
639   if (Jpre) {
640     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
641     PetscCheckSameComm(tao,1,Jpre,3);
642   }
643   if (ctx) {
644     tao->user_jac_inequalityP = ctx;
645   }
646   if (func) {
647     tao->ops->computejacobianinequality = func;
648   }
649   if (J) {
650     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
651     ierr = MatDestroy(&tao->jacobian_inequality);CHKERRQ(ierr);
652     tao->jacobian_inequality = J;
653   }
654   if (Jpre) {
655     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
656     ierr = MatDestroy(&tao->jacobian_inequality_pre);CHKERRQ(ierr);
657     tao->jacobian_inequality_pre=Jpre;
658   }
659   PetscFunctionReturn(0);
660 }
661