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