xref: /petsc/src/tao/interface/taosolver_hj.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
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: `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: `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) {
195           PetscCall(MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES));
196         }
197         PetscCall(MatRestoreRow(B,row,&bncols,&bcols,&bvals));
198         PetscCall(PetscFree2(ccols,cvals));
199       }
200       PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
201       PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
202       PetscCall(PetscViewerASCIIPrintf(viewer,"  Finite-difference minus hand-coded Hessian with tolerance %g ----------\n",(double)threshold));
203       PetscCall(MatView(C,mviewer));
204       PetscCall(MatDestroy(&C));
205     }
206     PetscCall(MatDestroy(&A));
207     PetscCall(MatDestroy(&B));
208 
209     if (hessian != tao->hessian_pre) {
210       hessian = tao->hessian_pre;
211       PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------- Testing Hessian for preconditioner -------------\n"));
212     } else hessian = NULL;
213   }
214   if (complete_print) {
215     PetscCall(PetscViewerPopFormat(mviewer));
216     PetscCall(PetscViewerDestroy(&mviewer));
217   }
218   PetscCall(PetscViewerASCIISetTab(viewer,tabs));
219   PetscFunctionReturn(0);
220 }
221 
222 /*@C
223    TaoComputeHessian - Computes the Hessian matrix that has been
224    set with TaoSetHessian().
225 
226    Collective on Tao
227 
228    Input Parameters:
229 +  tao - the Tao solver context
230 -  X   - input vector
231 
232    Output Parameters:
233 +  H    - Hessian matrix
234 -  Hpre - Preconditioning matrix
235 
236    Options Database Keys:
237 +     -tao_test_hessian - compare the user provided Hessian with one compute via finite differences to check for errors
238 .     -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
239 -     -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
240 
241    Notes:
242    Most users should not need to explicitly call this routine, as it
243    is used internally within the minimization solvers.
244 
245    TaoComputeHessian() is typically used within minimization
246    implementations, so most users would not generally call this routine
247    themselves.
248 
249    Developer Notes:
250    The Hessian test mechanism follows SNESTestJacobian().
251 
252    Level: developer
253 
254 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetHessian()`
255 @*/
256 PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
257 {
258   PetscFunctionBegin;
259   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
260   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
261   PetscCheckSameComm(tao,1,X,2);
262   PetscCheck(tao->ops->computehessian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessian() first");
263   ++tao->nhess;
264   PetscCall(VecLockReadPush(X));
265   PetscCall(PetscLogEventBegin(TAO_HessianEval,tao,X,H,Hpre));
266   PetscStackPush("Tao user Hessian function");
267   PetscCall((*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP));
268   PetscStackPop;
269   PetscCall(PetscLogEventEnd(TAO_HessianEval,tao,X,H,Hpre));
270   PetscCall(VecLockReadPop(X));
271 
272   PetscCall(TaoTestHessian(tao));
273   PetscFunctionReturn(0);
274 }
275 
276 /*@C
277    TaoComputeJacobian - Computes the Jacobian matrix that has been
278    set with TaoSetJacobianRoutine().
279 
280    Collective on Tao
281 
282    Input Parameters:
283 +  tao - the Tao solver context
284 -  X   - input vector
285 
286    Output Parameters:
287 +  J    - Jacobian matrix
288 -  Jpre - Preconditioning matrix
289 
290    Notes:
291    Most users should not need to explicitly call this routine, as it
292    is used internally within the minimization solvers.
293 
294    TaoComputeJacobian() is typically used within minimization
295    implementations, so most users would not generally call this routine
296    themselves.
297 
298    Level: developer
299 
300 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianRoutine()`
301 @*/
302 PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
303 {
304   PetscFunctionBegin;
305   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
306   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
307   PetscCheckSameComm(tao,1,X,2);
308   PetscCheck(tao->ops->computejacobian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
309   ++tao->njac;
310   PetscCall(VecLockReadPush(X));
311   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
312   PetscStackPush("Tao user Jacobian function");
313   PetscCall((*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP));
314   PetscStackPop;
315   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
316   PetscCall(VecLockReadPop(X));
317   PetscFunctionReturn(0);
318 }
319 
320 /*@C
321    TaoComputeResidualJacobian - Computes the least-squares residual Jacobian matrix that has been
322    set with TaoSetJacobianResidual().
323 
324    Collective on Tao
325 
326    Input Parameters:
327 +  tao - the Tao solver context
328 -  X   - input vector
329 
330    Output Parameters:
331 +  J    - Jacobian matrix
332 -  Jpre - Preconditioning matrix
333 
334    Notes:
335    Most users should not need to explicitly call this routine, as it
336    is used internally within the minimization solvers.
337 
338    TaoComputeResidualJacobian() is typically used within least-squares
339    implementations, so most users would not generally call this routine
340    themselves.
341 
342    Level: developer
343 
344 .seealso: `TaoComputeResidual()`, `TaoSetJacobianResidual()`
345 @*/
346 PetscErrorCode TaoComputeResidualJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
347 {
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
350   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
351   PetscCheckSameComm(tao,1,X,2);
352   PetscCheck(tao->ops->computeresidualjacobian,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetResidualJacobian() first");
353   ++tao->njac;
354   PetscCall(VecLockReadPush(X));
355   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
356   PetscStackPush("Tao user least-squares residual Jacobian function");
357   PetscCall((*tao->ops->computeresidualjacobian)(tao,X,J,Jpre,tao->user_lsjacP));
358   PetscStackPop;
359   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
360   PetscCall(VecLockReadPop(X));
361   PetscFunctionReturn(0);
362 }
363 
364 /*@C
365    TaoComputeJacobianState - Computes the Jacobian matrix that has been
366    set with TaoSetJacobianStateRoutine().
367 
368    Collective on Tao
369 
370    Input Parameters:
371 +  tao - the Tao solver context
372 -  X   - input vector
373 
374    Output Parameters:
375 +  J    - Jacobian matrix
376 .  Jpre - Preconditioning matrix
377 -  Jinv -
378 
379    Notes:
380    Most users should not need to explicitly call this routine, as it
381    is used internally within the minimization solvers.
382 
383    TaoComputeJacobianState() is typically used within minimization
384    implementations, so most users would not generally call this routine
385    themselves.
386 
387    Level: developer
388 
389 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
390 @*/
391 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
392 {
393   PetscFunctionBegin;
394   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
395   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
396   PetscCheckSameComm(tao,1,X,2);
397   PetscCheck(tao->ops->computejacobianstate,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
398   ++tao->njac_state;
399   PetscCall(VecLockReadPush(X));
400   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
401   PetscStackPush("Tao user Jacobian(state) function");
402   PetscCall((*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP));
403   PetscStackPop;
404   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
405   PetscCall(VecLockReadPop(X));
406   PetscFunctionReturn(0);
407 }
408 
409 /*@C
410    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
411    set with TaoSetJacobianDesignRoutine().
412 
413    Collective on Tao
414 
415    Input Parameters:
416 +  tao - the Tao solver context
417 -  X   - input vector
418 
419    Output Parameters:
420 .  J - Jacobian matrix
421 
422    Notes:
423    Most users should not need to explicitly call this routine, as it
424    is used internally within the minimization solvers.
425 
426    TaoComputeJacobianDesign() is typically used within minimization
427    implementations, so most users would not generally call this routine
428    themselves.
429 
430    Level: developer
431 
432 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianDesignRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
433 @*/
434 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
435 {
436   PetscFunctionBegin;
437   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
438   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
439   PetscCheckSameComm(tao,1,X,2);
440   PetscCheck(tao->ops->computejacobiandesign,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
441   ++tao->njac_design;
442   PetscCall(VecLockReadPush(X));
443   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL));
444   PetscStackPush("Tao user Jacobian(design) function");
445   PetscCall((*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP));
446   PetscStackPop;
447   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL));
448   PetscCall(VecLockReadPop(X));
449   PetscFunctionReturn(0);
450 }
451 
452 /*@C
453    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
454 
455    Logically collective on Tao
456 
457    Input Parameters:
458 +  tao  - the Tao context
459 .  J    - Matrix used for the jacobian
460 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
461 .  func - Jacobian evaluation routine
462 -  ctx  - [optional] user-defined context for private data for the
463           Jacobian evaluation routine (may be NULL)
464 
465    Calling sequence of func:
466 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
467 
468 +  tao  - the Tao  context
469 .  x    - input vector
470 .  J    - Jacobian matrix
471 .  Jpre - preconditioning matrix, usually the same as J
472 -  ctx  - [optional] user-defined Jacobian context
473 
474    Level: intermediate
475 @*/
476 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
477 {
478   PetscFunctionBegin;
479   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
480   if (J) {
481     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
482     PetscCheckSameComm(tao,1,J,2);
483   }
484   if (Jpre) {
485     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
486     PetscCheckSameComm(tao,1,Jpre,3);
487   }
488   if (ctx) {
489     tao->user_jacP = ctx;
490   }
491   if (func) {
492     tao->ops->computejacobian = func;
493   }
494   if (J) {
495     PetscCall(PetscObjectReference((PetscObject)J));
496     PetscCall(MatDestroy(&tao->jacobian));
497     tao->jacobian = J;
498   }
499   if (Jpre) {
500     PetscCall(PetscObjectReference((PetscObject)Jpre));
501     PetscCall(MatDestroy(&tao->jacobian_pre));
502     tao->jacobian_pre=Jpre;
503   }
504   PetscFunctionReturn(0);
505 }
506 
507 /*@C
508    TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the
509    location to store the matrix.
510 
511    Logically collective on Tao
512 
513    Input Parameters:
514 +  tao  - the Tao context
515 .  J    - Matrix used for the jacobian
516 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
517 .  func - Jacobian evaluation routine
518 -  ctx  - [optional] user-defined context for private data for the
519           Jacobian evaluation routine (may be NULL)
520 
521    Calling sequence of func:
522 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
523 
524 +  tao  - the Tao  context
525 .  x    - input vector
526 .  J    - Jacobian matrix
527 .  Jpre - preconditioning matrix, usually the same as J
528 -  ctx  - [optional] user-defined Jacobian context
529 
530    Level: intermediate
531 @*/
532 PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
533 {
534   PetscFunctionBegin;
535   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
536   if (J) {
537     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
538     PetscCheckSameComm(tao,1,J,2);
539   }
540   if (Jpre) {
541     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
542     PetscCheckSameComm(tao,1,Jpre,3);
543   }
544   if (ctx) {
545     tao->user_lsjacP = ctx;
546   }
547   if (func) {
548     tao->ops->computeresidualjacobian = func;
549   }
550   if (J) {
551     PetscCall(PetscObjectReference((PetscObject)J));
552     PetscCall(MatDestroy(&tao->ls_jac));
553     tao->ls_jac = J;
554   }
555   if (Jpre) {
556     PetscCall(PetscObjectReference((PetscObject)Jpre));
557     PetscCall(MatDestroy(&tao->ls_jac_pre));
558     tao->ls_jac_pre=Jpre;
559   }
560   PetscFunctionReturn(0);
561 }
562 
563 /*@C
564    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
565    (and its inverse) of the constraint function with respect to the state variables.
566    Used only for pde-constrained optimization.
567 
568    Logically collective on Tao
569 
570    Input Parameters:
571 +  tao  - the Tao context
572 .  J    - Matrix used for the jacobian
573 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is NULL
574 .  Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse.
575 .  func - Jacobian evaluation routine
576 -  ctx  - [optional] user-defined context for private data for the
577           Jacobian evaluation routine (may be NULL)
578 
579    Calling sequence of func:
580 $    func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx);
581 
582 +  tao  - the Tao  context
583 .  x    - input vector
584 .  J    - Jacobian matrix
585 .  Jpre - preconditioner matrix, usually the same as J
586 .  Jinv - inverse of J
587 -  ctx  - [optional] user-defined Jacobian context
588 
589    Level: intermediate
590 .seealso: `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()`
591 @*/
592 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx)
593 {
594   PetscFunctionBegin;
595   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
596   if (J) {
597     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
598     PetscCheckSameComm(tao,1,J,2);
599   }
600   if (Jpre) {
601     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
602     PetscCheckSameComm(tao,1,Jpre,3);
603   }
604   if (Jinv) {
605     PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
606     PetscCheckSameComm(tao,1,Jinv,4);
607   }
608   if (ctx) {
609     tao->user_jac_stateP = ctx;
610   }
611   if (func) {
612     tao->ops->computejacobianstate = func;
613   }
614   if (J) {
615     PetscCall(PetscObjectReference((PetscObject)J));
616     PetscCall(MatDestroy(&tao->jacobian_state));
617     tao->jacobian_state = J;
618   }
619   if (Jpre) {
620     PetscCall(PetscObjectReference((PetscObject)Jpre));
621     PetscCall(MatDestroy(&tao->jacobian_state_pre));
622     tao->jacobian_state_pre=Jpre;
623   }
624   if (Jinv) {
625     PetscCall(PetscObjectReference((PetscObject)Jinv));
626     PetscCall(MatDestroy(&tao->jacobian_state_inv));
627     tao->jacobian_state_inv=Jinv;
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 /*@C
633    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
634    the constraint function with respect to the design variables.  Used only for
635    pde-constrained optimization.
636 
637    Logically collective on Tao
638 
639    Input Parameters:
640 +  tao  - the Tao context
641 .  J    - Matrix used for the jacobian
642 .  func - Jacobian evaluation routine
643 -  ctx  - [optional] user-defined context for private data for the
644           Jacobian evaluation routine (may be NULL)
645 
646    Calling sequence of func:
647 $    func(Tao tao,Vec x,Mat J,void *ctx);
648 
649 +  tao - the Tao  context
650 .  x   - input vector
651 .  J   - Jacobian matrix
652 -  ctx - [optional] user-defined Jacobian context
653 
654    Level: intermediate
655 
656 .seealso: `TaoComputeJacobianDesign()`, `TaoSetJacobianStateRoutine()`, `TaoSetStateDesignIS()`
657 @*/
658 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
659 {
660   PetscFunctionBegin;
661   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
662   if (J) {
663     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
664     PetscCheckSameComm(tao,1,J,2);
665   }
666   if (ctx) {
667     tao->user_jac_designP = ctx;
668   }
669   if (func) {
670     tao->ops->computejacobiandesign = func;
671   }
672   if (J) {
673     PetscCall(PetscObjectReference((PetscObject)J));
674     PetscCall(MatDestroy(&tao->jacobian_design));
675     tao->jacobian_design = J;
676   }
677   PetscFunctionReturn(0);
678 }
679 
680 /*@
681    TaoSetStateDesignIS - Indicate to the Tao which variables in the
682    solution vector are state variables and which are design.  Only applies to
683    pde-constrained optimization.
684 
685    Logically Collective on Tao
686 
687    Input Parameters:
688 +  tao  - The Tao context
689 .  s_is - the index set corresponding to the state variables
690 -  d_is - the index set corresponding to the design variables
691 
692    Level: intermediate
693 
694 .seealso: `TaoSetJacobianStateRoutine()`, `TaoSetJacobianDesignRoutine()`
695 @*/
696 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
697 {
698   PetscFunctionBegin;
699   PetscCall(PetscObjectReference((PetscObject)s_is));
700   PetscCall(ISDestroy(&tao->state_is));
701   tao->state_is = s_is;
702   PetscCall(PetscObjectReference((PetscObject)(d_is)));
703   PetscCall(ISDestroy(&tao->design_is));
704   tao->design_is = d_is;
705   PetscFunctionReturn(0);
706 }
707 
708 /*@C
709    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
710    set with TaoSetJacobianEqualityRoutine().
711 
712    Collective on Tao
713 
714    Input Parameters:
715 +  tao - the Tao solver context
716 -  X   - input vector
717 
718    Output Parameters:
719 +  J    - Jacobian matrix
720 -  Jpre - Preconditioning matrix
721 
722    Notes:
723    Most users should not need to explicitly call this routine, as it
724    is used internally within the minimization solvers.
725 
726    Level: developer
727 
728 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
729 @*/
730 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
731 {
732   PetscFunctionBegin;
733   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
734   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
735   PetscCheckSameComm(tao,1,X,2);
736   PetscCheck(tao->ops->computejacobianequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
737   ++tao->njac_equality;
738   PetscCall(VecLockReadPush(X));
739   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
740   PetscStackPush("Tao user Jacobian(equality) function");
741   PetscCall((*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP));
742   PetscStackPop;
743   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
744   PetscCall(VecLockReadPop(X));
745   PetscFunctionReturn(0);
746 }
747 
748 /*@C
749    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
750    set with TaoSetJacobianInequalityRoutine().
751 
752    Collective on Tao
753 
754    Input Parameters:
755 +  tao - the Tao solver context
756 -  X   - input vector
757 
758    Output Parameters:
759 +  J    - Jacobian matrix
760 -  Jpre - Preconditioning matrix
761 
762    Notes:
763    Most users should not need to explicitly call this routine, as it
764    is used internally within the minimization solvers.
765 
766    Level: developer
767 
768 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
769 @*/
770 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
771 {
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
774   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
775   PetscCheckSameComm(tao,1,X,2);
776   PetscCheck(tao->ops->computejacobianinequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
777   ++tao->njac_inequality;
778   PetscCall(VecLockReadPush(X));
779   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
780   PetscStackPush("Tao user Jacobian(inequality) function");
781   PetscCall((*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP));
782   PetscStackPop;
783   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
784   PetscCall(VecLockReadPop(X));
785   PetscFunctionReturn(0);
786 }
787 
788 /*@C
789    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
790    (and its inverse) of the constraint function with respect to the equality variables.
791    Used only for pde-constrained optimization.
792 
793    Logically collective on Tao
794 
795    Input Parameters:
796 +  tao  - the Tao context
797 .  J    - Matrix used for the jacobian
798 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
799 .  func - Jacobian evaluation routine
800 -  ctx  - [optional] user-defined context for private data for the
801           Jacobian evaluation routine (may be NULL)
802 
803    Calling sequence of func:
804 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
805 
806 +  tao  - the Tao  context
807 .  x    - input vector
808 .  J    - Jacobian matrix
809 .  Jpre - preconditioner matrix, usually the same as J
810 -  ctx  - [optional] user-defined Jacobian context
811 
812    Level: intermediate
813 
814 .seealso: `TaoComputeJacobianEquality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetEqualityDesignIS()`
815 @*/
816 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
817 {
818   PetscFunctionBegin;
819   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
820   if (J) {
821     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
822     PetscCheckSameComm(tao,1,J,2);
823   }
824   if (Jpre) {
825     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
826     PetscCheckSameComm(tao,1,Jpre,3);
827   }
828   if (ctx) {
829     tao->user_jac_equalityP = ctx;
830   }
831   if (func) {
832     tao->ops->computejacobianequality = func;
833   }
834   if (J) {
835     PetscCall(PetscObjectReference((PetscObject)J));
836     PetscCall(MatDestroy(&tao->jacobian_equality));
837     tao->jacobian_equality = J;
838   }
839   if (Jpre) {
840     PetscCall(PetscObjectReference((PetscObject)Jpre));
841     PetscCall(MatDestroy(&tao->jacobian_equality_pre));
842     tao->jacobian_equality_pre=Jpre;
843   }
844   PetscFunctionReturn(0);
845 }
846 
847 /*@C
848    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
849    (and its inverse) of the constraint function with respect to the inequality variables.
850    Used only for pde-constrained optimization.
851 
852    Logically collective on Tao
853 
854    Input Parameters:
855 +  tao  - the Tao context
856 .  J    - Matrix used for the jacobian
857 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
858 .  func - Jacobian evaluation routine
859 -  ctx  - [optional] user-defined context for private data for the
860           Jacobian evaluation routine (may be NULL)
861 
862    Calling sequence of func:
863 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
864 
865 +  tao  - the Tao  context
866 .  x    - input vector
867 .  J    - Jacobian matrix
868 .  Jpre - preconditioner matrix, usually the same as J
869 -  ctx  - [optional] user-defined Jacobian context
870 
871    Level: intermediate
872 
873 .seealso: `TaoComputeJacobianInequality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetInequalityDesignIS()`
874 @*/
875 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
876 {
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
879   if (J) {
880     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
881     PetscCheckSameComm(tao,1,J,2);
882   }
883   if (Jpre) {
884     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
885     PetscCheckSameComm(tao,1,Jpre,3);
886   }
887   if (ctx) {
888     tao->user_jac_inequalityP = ctx;
889   }
890   if (func) {
891     tao->ops->computejacobianinequality = func;
892   }
893   if (J) {
894     PetscCall(PetscObjectReference((PetscObject)J));
895     PetscCall(MatDestroy(&tao->jacobian_inequality));
896     tao->jacobian_inequality = J;
897   }
898   if (Jpre) {
899     PetscCall(PetscObjectReference((PetscObject)Jpre));
900     PetscCall(MatDestroy(&tao->jacobian_inequality_pre));
901     tao->jacobian_inequality_pre=Jpre;
902   }
903   PetscFunctionReturn(0);
904 }
905