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