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