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