xref: /petsc/src/tao/interface/taosolver_hj.c (revision e0ffd71f47fcc4029f6759c19bf7e02dd09d54e9)
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 .  func - 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 func:
17 $    func(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 PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
28 {
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
33   if (H) {
34     PetscValidHeaderSpecific(H,MAT_CLASSID,2);
35     PetscCheckSameComm(tao,1,H,2);
36   }
37   if (Hpre) {
38     PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
39     PetscCheckSameComm(tao,1,Hpre,3);
40   }
41   if (ctx) {
42     tao->user_hessP = ctx;
43   }
44   if (func) {
45     tao->ops->computehessian = func;
46   }
47   if (H) {
48     ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr);
49     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
50     tao->hessian = H;
51   }
52   if (Hpre) {
53     ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr);
54     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
55     tao->hessian_pre = Hpre;
56   }
57   PetscFunctionReturn(0);
58 }
59 
60 PetscErrorCode TaoTestHessian(Tao tao)
61 {
62   Mat               A,B,C,D,hessian;
63   Vec               x = tao->solution;
64   PetscErrorCode    ierr;
65   PetscReal         nrm,gnorm;
66   PetscReal         threshold = 1.e-5;
67   PetscInt          m,n,M,N;
68   PetscBool         complete_print = PETSC_FALSE,test = PETSC_FALSE,flg;
69   PetscViewer       viewer,mviewer;
70   MPI_Comm          comm;
71   PetscInt          tabs;
72   static PetscBool  directionsprinted = PETSC_FALSE;
73   PetscViewerFormat format;
74 
75   PetscFunctionBegin;
76   ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr);
77   ierr = PetscOptionsName("-tao_test_hessian","Compare hand-coded and finite difference Hessians","None",&test);CHKERRQ(ierr);
78   ierr = PetscOptionsReal("-tao_test_hessian", "Threshold for element difference between hand-coded and finite difference being meaningful","None",threshold,&threshold,NULL);CHKERRQ(ierr);
79   ierr = PetscOptionsViewer("-tao_test_hessian_view","View difference between hand-coded and finite difference Hessians element entries","None",&mviewer,&format,&complete_print);CHKERRQ(ierr);
80   ierr = PetscOptionsEnd();CHKERRQ(ierr);
81   if (!test) PetscFunctionReturn(0);
82 
83   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
84   ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
85   ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr);
86   ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel);CHKERRQ(ierr);
87   ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Hessian -------------\n");CHKERRQ(ierr);
88   if (!complete_print && !directionsprinted) {
89     ierr = PetscViewerASCIIPrintf(viewer,"  Run with -tao_test_hessian_view and optionally -tao_test_hessian <threshold> to show difference\n");CHKERRQ(ierr);
90     ierr = PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference Hessian entries greater than <threshold>.\n");CHKERRQ(ierr);
91   }
92   if (!directionsprinted) {
93     ierr = PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Hessian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");CHKERRQ(ierr);
94     ierr = PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Hessian is probably correct.\n");CHKERRQ(ierr);
95     directionsprinted = PETSC_TRUE;
96   }
97   if (complete_print) {
98     ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr);
99   }
100 
101   ierr = PetscObjectTypeCompare((PetscObject)tao->hessian,MATMFFD,&flg);CHKERRQ(ierr);
102   if (!flg) hessian = tao->hessian;
103   else hessian = tao->hessian_pre;
104 
105   while (hessian) {
106     ierr = PetscObjectTypeCompareAny((PetscObject)hessian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr);
107     if (flg) {
108       A    = hessian;
109       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
110     } else {
111       ierr = MatComputeExplicitOperator(hessian,&A);CHKERRQ(ierr);
112     }
113 
114     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
115     ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
116     ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
117     ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr);
118     ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
119     ierr = MatSetUp(B);CHKERRQ(ierr);
120     ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
121 
122     ierr = TaoDefaultComputeHessian(tao,x,B,B,NULL);CHKERRQ(ierr);
123 
124     ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr);
125     ierr = MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
126     ierr = MatNorm(D,NORM_FROBENIUS,&nrm);CHKERRQ(ierr);
127     ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr);
128     ierr = MatDestroy(&D);CHKERRQ(ierr);
129     if (!gnorm) gnorm = 1; /* just in case */
130     ierr = PetscViewerASCIIPrintf(viewer,"  ||H - Hfd||_F/||H||_F = %g, ||H - Hfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);CHKERRQ(ierr);
131 
132     if (complete_print) {
133       ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded Hessian ----------\n");CHKERRQ(ierr);
134       ierr = MatView(hessian,mviewer);CHKERRQ(ierr);
135       ierr = PetscViewerASCIIPrintf(viewer,"  Finite difference Hessian ----------\n");CHKERRQ(ierr);
136       ierr = MatView(B,mviewer);CHKERRQ(ierr);
137     }
138 
139     if (complete_print) {
140       PetscInt          Istart, Iend, *ccols, bncols, cncols, j, row;
141       PetscScalar       *cvals;
142       const PetscInt    *bcols;
143       const PetscScalar *bvals;
144 
145       ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
146       ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
147       ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
148       ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
149       ierr = MatSetUp(C);CHKERRQ(ierr);
150       ierr = MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
151       ierr = MatGetOwnershipRange(B,&Istart,&Iend);CHKERRQ(ierr);
152 
153       for (row = Istart; row < Iend; row++) {
154         ierr = MatGetRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr);
155         ierr = PetscMalloc2(bncols,&ccols,bncols,&cvals);CHKERRQ(ierr);
156         for (j = 0, cncols = 0; j < bncols; j++) {
157           if (PetscAbsScalar(bvals[j]) > threshold) {
158             ccols[cncols] = bcols[j];
159             cvals[cncols] = bvals[j];
160             cncols += 1;
161           }
162         }
163         if (cncols) {
164           ierr = MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);CHKERRQ(ierr);
165         }
166         ierr = MatRestoreRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr);
167         ierr = PetscFree2(ccols,cvals);CHKERRQ(ierr);
168       }
169       ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170       ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171       ierr = PetscViewerASCIIPrintf(viewer,"  Finite-difference minus hand-coded Hessian with tolerance %g ----------\n",(double)threshold);CHKERRQ(ierr);
172       ierr = MatView(C,mviewer);CHKERRQ(ierr);
173       ierr = MatDestroy(&C);CHKERRQ(ierr);
174     }
175     ierr = MatDestroy(&A);CHKERRQ(ierr);
176     ierr = MatDestroy(&B);CHKERRQ(ierr);
177 
178     if (hessian != tao->hessian_pre) {
179       hessian = tao->hessian_pre;
180       ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Hessian for preconditioner -------------\n");CHKERRQ(ierr);
181     } else hessian = NULL;
182   }
183   if (complete_print) {
184     ierr = PetscViewerPopFormat(mviewer);CHKERRQ(ierr);
185   }
186   ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
187   PetscFunctionReturn(0);
188 }
189 
190 /*@C
191    TaoComputeHessian - Computes the Hessian matrix that has been
192    set with TaoSetHessianRoutine().
193 
194    Collective on Tao
195 
196    Input Parameters:
197 +  tao - the Tao solver context
198 -  X   - input vector
199 
200    Output Parameters:
201 +  H    - Hessian matrix
202 -  Hpre - Preconditioning matrix
203 
204    Options Database Keys:
205 +     -tao_test_hessian - compare the user provided Hessian with one compute via finite differences to check for errors
206 .     -tao_test_hessian <numerical value>  - display entries in the difference between the user provided Hessian and finite difference Hessian that are greater than a certain value to help users detect errors
207 -     -tao_test_hessian_view - display the user provided Hessian, the finite difference Hessian and the difference between them to help users detect the location of errors in the user provided Hessian
208 
209    Notes:
210    Most users should not need to explicitly call this routine, as it
211    is used internally within the minimization solvers.
212 
213    TaoComputeHessian() is typically used within minimization
214    implementations, so most users would not generally call this routine
215    themselves.
216 
217    Developer Notes:
218    The Hessian test mechanism follows SNESTestJacobian().
219 
220    Level: developer
221 
222 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessianRoutine()
223 @*/
224 PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
225 {
226   PetscErrorCode ierr;
227 
228   PetscFunctionBegin;
229   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
230   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
231   PetscCheckSameComm(tao,1,X,2);
232   if (!tao->ops->computehessian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessianRoutine() first");
233   ++tao->nhess;
234   ierr = VecLockPush(X);CHKERRQ(ierr);
235   ierr = PetscLogEventBegin(TAO_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
236   PetscStackPush("Tao user Hessian function");
237   ierr = (*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP);CHKERRQ(ierr);
238   PetscStackPop;
239   ierr = PetscLogEventEnd(TAO_HessianEval,tao,X,H,Hpre);CHKERRQ(ierr);
240   ierr = VecLockPop(X);CHKERRQ(ierr);
241 
242   ierr = TaoTestHessian(tao);CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 /*@C
247    TaoComputeJacobian - Computes the Jacobian matrix that has been
248    set with TaoSetJacobianRoutine().
249 
250    Collective on Tao
251 
252    Input Parameters:
253 +  tao - the Tao solver context
254 -  X   - input vector
255 
256    Output Parameters:
257 +  J    - Jacobian matrix
258 -  Jpre - Preconditioning matrix
259 
260    Notes:
261    Most users should not need to explicitly call this routine, as it
262    is used internally within the minimization solvers.
263 
264    TaoComputeJacobian() is typically used within minimization
265    implementations, so most users would not generally call this routine
266    themselves.
267 
268    Level: developer
269 
270 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianRoutine()
271 @*/
272 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
273 {
274   PetscErrorCode ierr;
275 
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
278   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
279   PetscCheckSameComm(tao,1,X,2);
280   if (!tao->ops->computejacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
281   ++tao->njac;
282   ierr = VecLockPush(X);CHKERRQ(ierr);
283   ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
284   PetscStackPush("Tao user Jacobian function");
285   ierr = (*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP);CHKERRQ(ierr);
286   PetscStackPop;
287   ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
288   ierr = VecLockPop(X);CHKERRQ(ierr);
289   PetscFunctionReturn(0);
290 }
291 
292 /*@C
293    TaoComputeJacobianState - Computes the Jacobian matrix that has been
294    set with TaoSetJacobianStateRoutine().
295 
296    Collective on Tao
297 
298    Input Parameters:
299 +  tao - the Tao solver context
300 -  X   - input vector
301 
302    Output Parameters:
303 +  Jpre - Jacobian matrix
304 -  Jinv - Preconditioning matrix
305 
306    Notes:
307    Most users should not need to explicitly call this routine, as it
308    is used internally within the minimization solvers.
309 
310    TaoComputeJacobianState() is typically used within minimization
311    implementations, so most users would not generally call this routine
312    themselves.
313 
314    Level: developer
315 
316 .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
317 @*/
318 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
319 {
320   PetscErrorCode ierr;
321 
322   PetscFunctionBegin;
323   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
324   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
325   PetscCheckSameComm(tao,1,X,2);
326   if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
327   ++tao->njac_state;
328   ierr = VecLockPush(X);CHKERRQ(ierr);
329   ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
330   PetscStackPush("Tao user Jacobian(state) function");
331   ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);CHKERRQ(ierr);
332   PetscStackPop;
333   ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
334   ierr = VecLockPop(X);CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 /*@C
339    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
340    set with TaoSetJacobianDesignRoutine().
341 
342    Collective on Tao
343 
344    Input Parameters:
345 +  tao - the Tao solver context
346 -  X   - input vector
347 
348    Output Parameters:
349 .  J - Jacobian matrix
350 
351    Notes:
352    Most users should not need to explicitly call this routine, as it
353    is used internally within the minimization solvers.
354 
355    TaoComputeJacobianDesign() is typically used within minimization
356    implementations, so most users would not generally call this routine
357    themselves.
358 
359    Level: developer
360 
361 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
362 @*/
363 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
364 {
365   PetscErrorCode ierr;
366 
367   PetscFunctionBegin;
368   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
369   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
370   PetscCheckSameComm(tao,1,X,2);
371   if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
372   ++tao->njac_design;
373   ierr = VecLockPush(X);CHKERRQ(ierr);
374   ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
375   PetscStackPush("Tao user Jacobian(design) function");
376   ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);CHKERRQ(ierr);
377   PetscStackPop;
378   ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL);CHKERRQ(ierr);
379   ierr = VecLockPop(X);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 /*@C
384    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
385 
386    Logically collective on Tao
387 
388    Input Parameters:
389 +  tao  - the Tao context
390 .  J    - Matrix used for the jacobian
391 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
392 .  func - Jacobian evaluation routine
393 -  ctx  - [optional] user-defined context for private data for the
394           Jacobian evaluation routine (may be NULL)
395 
396    Calling sequence of func:
397 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
398 
399 +  tao  - the Tao  context
400 .  x    - input vector
401 .  J    - Jacobian matrix
402 .  Jpre - preconditioning matrix, usually the same as J
403 -  ctx  - [optional] user-defined Jacobian context
404 
405    Level: intermediate
406 @*/
407 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
408 {
409   PetscErrorCode ierr;
410 
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
413   if (J) {
414     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
415     PetscCheckSameComm(tao,1,J,2);
416   }
417   if (Jpre) {
418     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
419     PetscCheckSameComm(tao,1,Jpre,3);
420   }
421   if (ctx) {
422     tao->user_jacP = ctx;
423   }
424   if (func) {
425     tao->ops->computejacobian = func;
426   }
427   if (J) {
428     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
429     ierr = MatDestroy(&tao->jacobian);CHKERRQ(ierr);
430     tao->jacobian = J;
431   }
432   if (Jpre) {
433     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
434     ierr = MatDestroy(&tao->jacobian_pre);CHKERRQ(ierr);
435     tao->jacobian_pre=Jpre;
436   }
437   PetscFunctionReturn(0);
438 }
439 
440 /*@C
441    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
442    (and its inverse) of the constraint function with respect to the state variables.
443    Used only for pde-constrained optimization.
444 
445    Logically collective on Tao
446 
447    Input Parameters:
448 +  tao  - the Tao context
449 .  J    - Matrix used for the jacobian
450 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is NULL
451 .  Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse.
452 .  func - Jacobian evaluation routine
453 -  ctx  - [optional] user-defined context for private data for the
454           Jacobian evaluation routine (may be NULL)
455 
456    Calling sequence of func:
457 $    func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx);
458 
459 +  tao  - the Tao  context
460 .  x    - input vector
461 .  J    - Jacobian matrix
462 .  Jpre - preconditioner matrix, usually the same as J
463 .  Jinv - inverse of J
464 -  ctx  - [optional] user-defined Jacobian context
465 
466    Level: intermediate
467 .seealso: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
468 @*/
469 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx)
470 {
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
475   if (J) {
476     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
477     PetscCheckSameComm(tao,1,J,2);
478   }
479   if (Jpre) {
480     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
481     PetscCheckSameComm(tao,1,Jpre,3);
482   }
483   if (Jinv) {
484     PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
485     PetscCheckSameComm(tao,1,Jinv,4);
486   }
487   if (ctx) {
488     tao->user_jac_stateP = ctx;
489   }
490   if (func) {
491     tao->ops->computejacobianstate = func;
492   }
493   if (J) {
494     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
495     ierr = MatDestroy(&tao->jacobian_state);CHKERRQ(ierr);
496     tao->jacobian_state = J;
497   }
498   if (Jpre) {
499     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
500     ierr = MatDestroy(&tao->jacobian_state_pre);CHKERRQ(ierr);
501     tao->jacobian_state_pre=Jpre;
502   }
503   if (Jinv) {
504     ierr = PetscObjectReference((PetscObject)Jinv);CHKERRQ(ierr);
505     ierr = MatDestroy(&tao->jacobian_state_inv);CHKERRQ(ierr);
506     tao->jacobian_state_inv=Jinv;
507   }
508   PetscFunctionReturn(0);
509 }
510 
511 /*@C
512    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
513    the constraint function with respect to the design variables.  Used only for
514    pde-constrained optimization.
515 
516    Logically collective on Tao
517 
518    Input Parameters:
519 +  tao  - the Tao context
520 .  J    - Matrix used for the jacobian
521 .  func - Jacobian evaluation routine
522 -  ctx  - [optional] user-defined context for private data for the
523           Jacobian evaluation routine (may be NULL)
524 
525    Calling sequence of func:
526 $    func(Tao tao,Vec x,Mat J,void *ctx);
527 
528 +  tao - the Tao  context
529 .  x   - input vector
530 .  J   - Jacobian matrix
531 -  ctx - [optional] user-defined Jacobian context
532 
533    Level: intermediate
534 
535 .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS()
536 @*/
537 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
538 {
539   PetscErrorCode ierr;
540 
541   PetscFunctionBegin;
542   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
543   if (J) {
544     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
545     PetscCheckSameComm(tao,1,J,2);
546   }
547   if (ctx) {
548     tao->user_jac_designP = ctx;
549   }
550   if (func) {
551     tao->ops->computejacobiandesign = func;
552   }
553   if (J) {
554     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
555     ierr = MatDestroy(&tao->jacobian_design);CHKERRQ(ierr);
556     tao->jacobian_design = J;
557   }
558   PetscFunctionReturn(0);
559 }
560 
561 /*@
562    TaoSetStateDesignIS - Indicate to the Tao which variables in the
563    solution vector are state variables and which are design.  Only applies to
564    pde-constrained optimization.
565 
566    Logically Collective on Tao
567 
568    Input Parameters:
569 +  tao  - The Tao context
570 .  s_is - the index set corresponding to the state variables
571 -  d_is - the index set corresponding to the design variables
572 
573    Level: intermediate
574 
575 .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine()
576 @*/
577 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
578 {
579   PetscErrorCode ierr;
580 
581   PetscFunctionBegin;
582   ierr = PetscObjectReference((PetscObject)s_is);CHKERRQ(ierr);
583   ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr);
584   tao->state_is = s_is;
585   ierr = PetscObjectReference((PetscObject)(d_is));CHKERRQ(ierr);
586   ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr);
587   tao->design_is = d_is;
588   PetscFunctionReturn(0);
589 }
590 
591 /*@C
592    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
593    set with TaoSetJacobianEqualityRoutine().
594 
595    Collective on Tao
596 
597    Input Parameters:
598 +  tao - the Tao solver context
599 -  X   - input vector
600 
601    Output Parameters:
602 +  J    - Jacobian matrix
603 -  Jpre - Preconditioning matrix
604 
605    Notes:
606    Most users should not need to explicitly call this routine, as it
607    is used internally within the minimization solvers.
608 
609    Level: developer
610 
611 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
612 @*/
613 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
614 {
615   PetscErrorCode ierr;
616 
617   PetscFunctionBegin;
618   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
619   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
620   PetscCheckSameComm(tao,1,X,2);
621   if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
622   ++tao->njac_equality;
623   ierr = VecLockPush(X);CHKERRQ(ierr);
624   ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
625   PetscStackPush("Tao user Jacobian(equality) function");
626   ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);CHKERRQ(ierr);
627   PetscStackPop;
628   ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
629   ierr = VecLockPop(X);CHKERRQ(ierr);
630   PetscFunctionReturn(0);
631 }
632 
633 /*@C
634    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
635    set with TaoSetJacobianInequalityRoutine().
636 
637    Collective on Tao
638 
639    Input Parameters:
640 +  tao - the Tao solver context
641 -  X   - input vector
642 
643    Output Parameters:
644 +  J    - Jacobian matrix
645 -  Jpre - Preconditioning matrix
646 
647    Notes:
648    Most users should not need to explicitly call this routine, as it
649    is used internally within the minimization solvers.
650 
651    Level: developer
652 
653 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
654 @*/
655 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
656 {
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
661   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
662   PetscCheckSameComm(tao,1,X,2);
663   if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
664   ++tao->njac_inequality;
665   ierr = VecLockPush(X);CHKERRQ(ierr);
666   ierr = PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
667   PetscStackPush("Tao user Jacobian(inequality) function");
668   ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);CHKERRQ(ierr);
669   PetscStackPop;
670   ierr = PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre);CHKERRQ(ierr);
671   ierr = VecLockPop(X);CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 
675 /*@C
676    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
677    (and its inverse) of the constraint function with respect to the equality variables.
678    Used only for pde-constrained optimization.
679 
680    Logically collective on Tao
681 
682    Input Parameters:
683 +  tao  - the Tao context
684 .  J    - Matrix used for the jacobian
685 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
686 .  func - Jacobian evaluation routine
687 -  ctx  - [optional] user-defined context for private data for the
688           Jacobian evaluation routine (may be NULL)
689 
690    Calling sequence of func:
691 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
692 
693 +  tao  - the Tao  context
694 .  x    - input vector
695 .  J    - Jacobian matrix
696 .  Jpre - preconditioner matrix, usually the same as J
697 -  ctx  - [optional] user-defined Jacobian context
698 
699    Level: intermediate
700 
701 .seealso: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS()
702 @*/
703 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
704 {
705   PetscErrorCode ierr;
706 
707   PetscFunctionBegin;
708   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
709   if (J) {
710     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
711     PetscCheckSameComm(tao,1,J,2);
712   }
713   if (Jpre) {
714     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
715     PetscCheckSameComm(tao,1,Jpre,3);
716   }
717   if (ctx) {
718     tao->user_jac_equalityP = ctx;
719   }
720   if (func) {
721     tao->ops->computejacobianequality = func;
722   }
723   if (J) {
724     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
725     ierr = MatDestroy(&tao->jacobian_equality);CHKERRQ(ierr);
726     tao->jacobian_equality = J;
727   }
728   if (Jpre) {
729     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
730     ierr = MatDestroy(&tao->jacobian_equality_pre);CHKERRQ(ierr);
731     tao->jacobian_equality_pre=Jpre;
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 /*@C
737    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
738    (and its inverse) of the constraint function with respect to the inequality variables.
739    Used only for pde-constrained optimization.
740 
741    Logically collective on Tao
742 
743    Input Parameters:
744 +  tao  - the Tao context
745 .  J    - Matrix used for the jacobian
746 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
747 .  func - Jacobian evaluation routine
748 -  ctx  - [optional] user-defined context for private data for the
749           Jacobian evaluation routine (may be NULL)
750 
751    Calling sequence of func:
752 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
753 
754 +  tao  - the Tao  context
755 .  x    - input vector
756 .  J    - Jacobian matrix
757 .  Jpre - preconditioner matrix, usually the same as J
758 -  ctx  - [optional] user-defined Jacobian context
759 
760    Level: intermediate
761 
762 .seealso: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS()
763 @*/
764 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
765 {
766   PetscErrorCode ierr;
767 
768   PetscFunctionBegin;
769   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
770   if (J) {
771     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
772     PetscCheckSameComm(tao,1,J,2);
773   }
774   if (Jpre) {
775     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
776     PetscCheckSameComm(tao,1,Jpre,3);
777   }
778   if (ctx) {
779     tao->user_jac_inequalityP = ctx;
780   }
781   if (func) {
782     tao->ops->computejacobianinequality = func;
783   }
784   if (J) {
785     ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
786     ierr = MatDestroy(&tao->jacobian_inequality);CHKERRQ(ierr);
787     tao->jacobian_inequality = J;
788   }
789   if (Jpre) {
790     ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr);
791     ierr = MatDestroy(&tao->jacobian_inequality_pre);CHKERRQ(ierr);
792     tao->jacobian_inequality_pre=Jpre;
793   }
794   PetscFunctionReturn(0);
795 }
796