xref: /petsc/src/tao/interface/taosolver_hj.c (revision 2e16c0ce58b3a4ec287cbc0a0807bfb0a0fa5ac9)
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) {
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 optimization algorithms,
246    so most users would not generally call this routine
247    themselves.
248 
249    Developer Note:
250    The Hessian test mechanism follows `SNESTestJacobian()`.
251 
252    Level: developer
253 
254 .seealso: `Tao`, `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: `Tao`, `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 - unknown
378 
379    Notes:
380    Most users should not need to explicitly call this routine, as it
381    is used internally within the optimization algorithms.
382 
383    Level: developer
384 
385 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
386 @*/
387 PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
388 {
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
391   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
392   PetscCheckSameComm(tao,1,X,2);
393   PetscCheck(tao->ops->computejacobianstate,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
394   ++tao->njac_state;
395   PetscCall(VecLockReadPush(X));
396   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
397   PetscStackPush("Tao user Jacobian(state) function");
398   PetscCall((*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP));
399   PetscStackPop;
400   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
401   PetscCall(VecLockReadPop(X));
402   PetscFunctionReturn(0);
403 }
404 
405 /*@C
406    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
407    set with `TaoSetJacobianDesignRoutine()`.
408 
409    Collective on tao
410 
411    Input Parameters:
412 +  tao - the Tao solver context
413 -  X   - input vector
414 
415    Output Parameters:
416 .  J - Jacobian matrix
417 
418    Notes:
419    Most users should not need to explicitly call this routine, as it
420    is used internally within the optimization algorithms.
421 
422    Level: developer
423 
424 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianDesignRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
425 @*/
426 PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
427 {
428   PetscFunctionBegin;
429   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
430   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
431   PetscCheckSameComm(tao,1,X,2);
432   PetscCheck(tao->ops->computejacobiandesign,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
433   ++tao->njac_design;
434   PetscCall(VecLockReadPush(X));
435   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,NULL));
436   PetscStackPush("Tao user Jacobian(design) function");
437   PetscCall((*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP));
438   PetscStackPop;
439   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,NULL));
440   PetscCall(VecLockReadPop(X));
441   PetscFunctionReturn(0);
442 }
443 
444 /*@C
445    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
446 
447    Logically collective on tao
448 
449    Input Parameters:
450 +  tao  - the Tao context
451 .  J    - Matrix used for the jacobian
452 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
453 .  func - Jacobian evaluation routine
454 -  ctx  - [optional] user-defined context for private data for the
455           Jacobian evaluation routine (may be NULL)
456 
457    Calling sequence of func:
458 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
459 
460 +  tao  - the Tao  context
461 .  x    - input vector
462 .  J    - Jacobian matrix
463 .  Jpre - preconditioning matrix, usually the same as J
464 -  ctx  - [optional] user-defined Jacobian context
465 
466    Level: intermediate
467 
468 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()`
469 @*/
470 PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
471 {
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
474   if (J) {
475     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
476     PetscCheckSameComm(tao,1,J,2);
477   }
478   if (Jpre) {
479     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
480     PetscCheckSameComm(tao,1,Jpre,3);
481   }
482   if (ctx) {
483     tao->user_jacP = ctx;
484   }
485   if (func) {
486     tao->ops->computejacobian = func;
487   }
488   if (J) {
489     PetscCall(PetscObjectReference((PetscObject)J));
490     PetscCall(MatDestroy(&tao->jacobian));
491     tao->jacobian = J;
492   }
493   if (Jpre) {
494     PetscCall(PetscObjectReference((PetscObject)Jpre));
495     PetscCall(MatDestroy(&tao->jacobian_pre));
496     tao->jacobian_pre=Jpre;
497   }
498   PetscFunctionReturn(0);
499 }
500 
501 /*@C
502    TaoSetJacobianResidualRoutine - Sets the function to compute the least-squares residual Jacobian as well as the
503    location to store the matrix.
504 
505    Logically collective on tao
506 
507    Input Parameters:
508 +  tao  - the Tao context
509 .  J    - Matrix used for the jacobian
510 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
511 .  func - Jacobian evaluation routine
512 -  ctx  - [optional] user-defined context for private data for the
513           Jacobian evaluation routine (may be NULL)
514 
515    Calling sequence of func:
516 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
517 
518 +  tao  - the Tao  context
519 .  x    - input vector
520 .  J    - Jacobian matrix
521 .  Jpre - preconditioning matrix, usually the same as J
522 -  ctx  - [optional] user-defined Jacobian context
523 
524    Level: intermediate
525 
526 .seealso: `Tao`, `TaoSetGradient()`, `TaoSetObjective()`
527 @*/
528 PetscErrorCode TaoSetJacobianResidualRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
529 {
530   PetscFunctionBegin;
531   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
532   if (J) {
533     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
534     PetscCheckSameComm(tao,1,J,2);
535   }
536   if (Jpre) {
537     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
538     PetscCheckSameComm(tao,1,Jpre,3);
539   }
540   if (ctx) {
541     tao->user_lsjacP = ctx;
542   }
543   if (func) {
544     tao->ops->computeresidualjacobian = func;
545   }
546   if (J) {
547     PetscCall(PetscObjectReference((PetscObject)J));
548     PetscCall(MatDestroy(&tao->ls_jac));
549     tao->ls_jac = J;
550   }
551   if (Jpre) {
552     PetscCall(PetscObjectReference((PetscObject)Jpre));
553     PetscCall(MatDestroy(&tao->ls_jac_pre));
554     tao->ls_jac_pre=Jpre;
555   }
556   PetscFunctionReturn(0);
557 }
558 
559 /*@C
560    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
561    (and its inverse) of the constraint function with respect to the state variables.
562    Used only for PDE-constrained optimization.
563 
564    Logically collective on tao
565 
566    Input Parameters:
567 +  tao  - the Tao context
568 .  J    - Matrix used for the jacobian
569 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is NULL
570 .  Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse.
571 .  func - Jacobian evaluation routine
572 -  ctx  - [optional] user-defined context for private data for the
573           Jacobian evaluation routine (may be NULL)
574 
575    Calling sequence of func:
576 $    func(Tao tao,Vec x,Mat J,Mat Jpre,Mat Jinv,void *ctx);
577 
578 +  tao  - the Tao  context
579 .  x    - input vector
580 .  J    - Jacobian matrix
581 .  Jpre - preconditioner matrix, usually the same as J
582 .  Jinv - inverse of J
583 -  ctx  - [optional] user-defined Jacobian context
584 
585    Level: intermediate
586 
587 .seealso: `Tao`, `TaoComputeJacobianState()`, `TaoSetJacobianDesignRoutine()`, `TaoSetStateDesignIS()`
588 @*/
589 PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat, void*), void *ctx)
590 {
591   PetscFunctionBegin;
592   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
593   if (J) {
594     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
595     PetscCheckSameComm(tao,1,J,2);
596   }
597   if (Jpre) {
598     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
599     PetscCheckSameComm(tao,1,Jpre,3);
600   }
601   if (Jinv) {
602     PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
603     PetscCheckSameComm(tao,1,Jinv,4);
604   }
605   if (ctx) {
606     tao->user_jac_stateP = ctx;
607   }
608   if (func) {
609     tao->ops->computejacobianstate = func;
610   }
611   if (J) {
612     PetscCall(PetscObjectReference((PetscObject)J));
613     PetscCall(MatDestroy(&tao->jacobian_state));
614     tao->jacobian_state = J;
615   }
616   if (Jpre) {
617     PetscCall(PetscObjectReference((PetscObject)Jpre));
618     PetscCall(MatDestroy(&tao->jacobian_state_pre));
619     tao->jacobian_state_pre=Jpre;
620   }
621   if (Jinv) {
622     PetscCall(PetscObjectReference((PetscObject)Jinv));
623     PetscCall(MatDestroy(&tao->jacobian_state_inv));
624     tao->jacobian_state_inv=Jinv;
625   }
626   PetscFunctionReturn(0);
627 }
628 
629 /*@C
630    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
631    the constraint function with respect to the design variables.  Used only for
632    PDE-constrained optimization.
633 
634    Logically collective on tao
635 
636    Input Parameters:
637 +  tao  - the Tao context
638 .  J    - Matrix used for the jacobian
639 .  func - Jacobian evaluation routine
640 -  ctx  - [optional] user-defined context for private data for the
641           Jacobian evaluation routine (may be NULL)
642 
643    Calling sequence of func:
644 $    func(Tao tao,Vec x,Mat J,void *ctx);
645 
646 +  tao - the Tao  context
647 .  x   - input vector
648 .  J   - Jacobian matrix
649 -  ctx - [optional] user-defined Jacobian context
650 
651    Level: intermediate
652 
653 .seealso: `Tao`, `TaoComputeJacobianDesign()`, `TaoSetJacobianStateRoutine()`, `TaoSetStateDesignIS()`
654 @*/
655 PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
656 {
657   PetscFunctionBegin;
658   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
659   if (J) {
660     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
661     PetscCheckSameComm(tao,1,J,2);
662   }
663   if (ctx) {
664     tao->user_jac_designP = ctx;
665   }
666   if (func) {
667     tao->ops->computejacobiandesign = func;
668   }
669   if (J) {
670     PetscCall(PetscObjectReference((PetscObject)J));
671     PetscCall(MatDestroy(&tao->jacobian_design));
672     tao->jacobian_design = J;
673   }
674   PetscFunctionReturn(0);
675 }
676 
677 /*@
678    TaoSetStateDesignIS - Indicate to the Tao which variables in the
679    solution vector are state variables and which are design.  Only applies to
680    PDE-constrained optimization.
681 
682    Logically Collective on Tao
683 
684    Input Parameters:
685 +  tao  - The Tao context
686 .  s_is - the index set corresponding to the state variables
687 -  d_is - the index set corresponding to the design variables
688 
689    Level: intermediate
690 
691 .seealso: `Tao`, `TaoSetJacobianStateRoutine()`, `TaoSetJacobianDesignRoutine()`
692 @*/
693 PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
694 {
695   PetscFunctionBegin;
696   PetscCall(PetscObjectReference((PetscObject)s_is));
697   PetscCall(ISDestroy(&tao->state_is));
698   tao->state_is = s_is;
699   PetscCall(PetscObjectReference((PetscObject)(d_is)));
700   PetscCall(ISDestroy(&tao->design_is));
701   tao->design_is = d_is;
702   PetscFunctionReturn(0);
703 }
704 
705 /*@C
706    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
707    set with `TaoSetJacobianEqualityRoutine()`.
708 
709    Collective on tao
710 
711    Input Parameters:
712 +  tao - the Tao solver context
713 -  X   - input vector
714 
715    Output Parameters:
716 +  J    - Jacobian matrix
717 -  Jpre - Preconditioning matrix
718 
719    Notes:
720    Most users should not need to explicitly call this routine, as it
721    is used internally within the optimization algorithms.
722 
723    Level: developer
724 
725 .seealso: `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
726 @*/
727 PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
728 {
729   PetscFunctionBegin;
730   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
731   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
732   PetscCheckSameComm(tao,1,X,2);
733   PetscCheck(tao->ops->computejacobianequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
734   ++tao->njac_equality;
735   PetscCall(VecLockReadPush(X));
736   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
737   PetscStackPush("Tao user Jacobian(equality) function");
738   PetscCall((*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP));
739   PetscStackPop;
740   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
741   PetscCall(VecLockReadPop(X));
742   PetscFunctionReturn(0);
743 }
744 
745 /*@C
746    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
747    set with `TaoSetJacobianInequalityRoutine()`.
748 
749    Collective on tao
750 
751    Input Parameters:
752 +  tao - the Tao solver context
753 -  X   - input vector
754 
755    Output Parameters:
756 +  J    - Jacobian matrix
757 -  Jpre - Preconditioning matrix
758 
759    Notes:
760    Most users should not need to explicitly call this routine, as it
761    is used internally within the minimization solvers.
762 
763    Level: developer
764 
765 .seealso: `Tao`, `TaoComputeObjective()`, `TaoComputeObjectiveAndGradient()`, `TaoSetJacobianStateRoutine()`, `TaoComputeJacobianDesign()`, `TaoSetStateDesignIS()`
766 @*/
767 PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
768 {
769   PetscFunctionBegin;
770   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
771   PetscValidHeaderSpecific(X, VEC_CLASSID,2);
772   PetscCheckSameComm(tao,1,X,2);
773   PetscCheck(tao->ops->computejacobianinequality,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
774   ++tao->njac_inequality;
775   PetscCall(VecLockReadPush(X));
776   PetscCall(PetscLogEventBegin(TAO_JacobianEval,tao,X,J,Jpre));
777   PetscStackPush("Tao user Jacobian(inequality) function");
778   PetscCall((*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP));
779   PetscStackPop;
780   PetscCall(PetscLogEventEnd(TAO_JacobianEval,tao,X,J,Jpre));
781   PetscCall(VecLockReadPop(X));
782   PetscFunctionReturn(0);
783 }
784 
785 /*@C
786    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
787    (and its inverse) of the constraint function with respect to the equality 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`, `TaoComputeJacobianEquality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetEqualityDesignIS()`
812 @*/
813 PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
814 {
815   PetscFunctionBegin;
816   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
817   if (J) {
818     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
819     PetscCheckSameComm(tao,1,J,2);
820   }
821   if (Jpre) {
822     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
823     PetscCheckSameComm(tao,1,Jpre,3);
824   }
825   if (ctx) {
826     tao->user_jac_equalityP = ctx;
827   }
828   if (func) {
829     tao->ops->computejacobianequality = func;
830   }
831   if (J) {
832     PetscCall(PetscObjectReference((PetscObject)J));
833     PetscCall(MatDestroy(&tao->jacobian_equality));
834     tao->jacobian_equality = J;
835   }
836   if (Jpre) {
837     PetscCall(PetscObjectReference((PetscObject)Jpre));
838     PetscCall(MatDestroy(&tao->jacobian_equality_pre));
839     tao->jacobian_equality_pre=Jpre;
840   }
841   PetscFunctionReturn(0);
842 }
843 
844 /*@C
845    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
846    (and its inverse) of the constraint function with respect to the inequality variables.
847    Used only for PDE-constrained optimization.
848 
849    Logically collective on tao
850 
851    Input Parameters:
852 +  tao  - the Tao context
853 .  J    - Matrix used for the jacobian
854 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
855 .  func - Jacobian evaluation routine
856 -  ctx  - [optional] user-defined context for private data for the
857           Jacobian evaluation routine (may be NULL)
858 
859    Calling sequence of func:
860 $    func(Tao tao,Vec x,Mat J,Mat Jpre,void *ctx);
861 
862 +  tao  - the Tao  context
863 .  x    - input vector
864 .  J    - Jacobian matrix
865 .  Jpre - preconditioner matrix, usually the same as J
866 -  ctx  - [optional] user-defined Jacobian context
867 
868    Level: intermediate
869 
870 .seealso: `Tao`, `TaoComputeJacobianInequality()`, `TaoSetJacobianDesignRoutine()`, `TaoSetInequalityDesignIS()`
871 @*/
872 PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
873 {
874   PetscFunctionBegin;
875   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
876   if (J) {
877     PetscValidHeaderSpecific(J,MAT_CLASSID,2);
878     PetscCheckSameComm(tao,1,J,2);
879   }
880   if (Jpre) {
881     PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
882     PetscCheckSameComm(tao,1,Jpre,3);
883   }
884   if (ctx) {
885     tao->user_jac_inequalityP = ctx;
886   }
887   if (func) {
888     tao->ops->computejacobianinequality = func;
889   }
890   if (J) {
891     PetscCall(PetscObjectReference((PetscObject)J));
892     PetscCall(MatDestroy(&tao->jacobian_inequality));
893     tao->jacobian_inequality = J;
894   }
895   if (Jpre) {
896     PetscCall(PetscObjectReference((PetscObject)Jpre));
897     PetscCall(MatDestroy(&tao->jacobian_inequality_pre));
898     tao->jacobian_inequality_pre=Jpre;
899   }
900   PetscFunctionReturn(0);
901 }
902