xref: /petsc/src/tao/interface/taosolver_hj.c (revision b0a7d7e7f246badab30cb8ce3f95dd1540bfb513)
1 #include "tao-private/taosolver_impl.h" /*I "taosolver.h" I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "TaoSetHessianRoutine"
5 /*@C
6    TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix.
7 
8    Logically collective on TaoSolver
9 
10    Input Parameters:
11 +  tao - the TaoSolver context
12 .  H - Matrix used for the hessian
13 .  Hpre - Matrix that will be used operated on by preconditioner, can be same as H
14 .  hess - Hessian evaluation routine
15 -  ctx - [optional] user-defined context for private data for the
16          Hessian evaluation routine (may be PETSC_NULL)
17 
18    Calling sequence of hess:
19 $    hess (TaoSolver tao,Vec x,Mat *H,Mat *Hpre,MatStructure *flag,void *ctx);
20 
21 +  tao - the TaoSolver  context
22 .  x - input vector
23 .  H - Hessian matrix
24 .  Hpre - preconditioner matrix, usually the same as H
25 .  flag - flag indicating information about the preconditioner matrix
26    structure (see below)
27 -  ctx - [optional] user-defined Hessian context
28 
29 
30    Notes:
31 
32    The function hess() takes Mat * as the matrix arguments rather than Mat.
33    This allows the Hessian evaluation routine to replace A and/or B with a
34    completely new new matrix structure (not just different matrix elements)
35    when appropriate, for instance, if the nonzero structure is changing
36    throughout the global iterations.
37 
38    The flag can be used to eliminate unnecessary work in the preconditioner
39    during the repeated solution of linear systems of the same size.  The
40    available options are
41 $    SAME_PRECONDITIONER -
42 $      Hpre is identical during successive linear solves.
43 $      This option is intended for folks who are using
44 $      different Amat and Pmat matrices and want to reuse the
45 $      same preconditioner matrix.  For example, this option
46 $      saves work by not recomputing incomplete factorization
47 $      for ILU/ICC preconditioners.
48 $    SAME_NONZERO_PATTERN -
49 $      Hpre has the same nonzero structure during
50 $      successive linear solves.
51 $    DIFFERENT_NONZERO_PATTERN -
52 $      Hpre does not have the same nonzero structure.
53 
54    Caution:
55    If you specify SAME_NONZERO_PATTERN, the software believes your assertion
56    and does not check the structure of the matrix.  If you erroneously
57    claim that the structure is the same when it actually is not, the new
58    preconditioner will not function correctly.  Thus, use this optimization
59    feature carefully!
60 
61    If in doubt about whether your preconditioner matrix has changed
62    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
63 
64    Level: beginner
65 
66 @*/
67 PetscErrorCode TaoSetHessianRoutine(TaoSolver tao, Mat H, Mat Hpre, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, Mat *, MatStructure *, void*), void *ctx)
68 {
69     PetscErrorCode ierr;
70     PetscFunctionBegin;
71     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
72     if (H) {
73 	PetscValidHeaderSpecific(H,MAT_CLASSID,2);
74 	PetscCheckSameComm(tao,1,H,2);
75     }
76     if (Hpre) {
77 	PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3);
78 	PetscCheckSameComm(tao,1,Hpre,3);
79     }
80     if (ctx) {
81 	tao->user_hessP = ctx;
82     }
83     if (func) {
84 	tao->ops->computehessian = func;
85     }
86 
87 
88     if (H) {
89 	ierr = PetscObjectReference((PetscObject)H); CHKERRQ(ierr);
90 	if (tao->hessian) {   ierr = MatDestroy(&tao->hessian); CHKERRQ(ierr);}
91 	tao->hessian = H;
92     }
93     if (Hpre) {
94 	ierr = PetscObjectReference((PetscObject)Hpre); CHKERRQ(ierr);
95 	if (tao->hessian_pre) { ierr = MatDestroy(&tao->hessian_pre); CHKERRQ(ierr);}
96 	tao->hessian_pre=Hpre;
97     }
98     PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "TaoComputeHessian"
103 /*@C
104    TaoComputeHessian - Computes the Hessian matrix that has been
105    set with TaoSetHessianRoutine().
106 
107    Collective on TaoSolver
108 
109    Input Parameters:
110 +  solver - the TaoSolver solver context
111 -  xx - input vector
112 
113    Output Parameters:
114 +  H - Hessian matrix
115 .  Hpre - Preconditioning matrix
116 -  flag - flag indicating matrix structure (SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, or SAME_PRECONDITIONER)
117 
118    Notes:
119    Most users should not need to explicitly call this routine, as it
120    is used internally within the minimization solvers.
121 
122    TaoComputeHessian() is typically used within minimization
123    implementations, so most users would not generally call this routine
124    themselves.
125 
126    Level: developer
127 
128 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessian()
129 
130 @*/
131 PetscErrorCode TaoComputeHessian(TaoSolver tao, Vec X, Mat *H, Mat *Hpre, MatStructure *flg)
132 {
133     PetscErrorCode ierr;
134     PetscFunctionBegin;
135     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
136     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
137     PetscValidPointer(flg,5);
138     PetscCheckSameComm(tao,1,X,2);
139 
140     if (!tao->ops->computehessian) {
141 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessian() first");
142     }
143     *flg = DIFFERENT_NONZERO_PATTERN;
144     ++tao->nhess;
145     ierr = PetscLogEventBegin(TaoSolver_HessianEval,tao,X,*H,*Hpre); CHKERRQ(ierr);
146     PetscStackPush("TaoSolver user Hessian function");
147     CHKMEMQ;
148     ierr = (*tao->ops->computehessian)(tao,X,H,Hpre,flg,tao->user_hessP); CHKERRQ(ierr);
149     CHKMEMQ;
150     PetscStackPop;
151     ierr = PetscLogEventEnd(TaoSolver_HessianEval,tao,X,*H,*Hpre); CHKERRQ(ierr);
152 
153     PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "TaoComputeJacobian"
158 /*@C
159    TaoComputeJacobian - Computes the Jacobian matrix that has been
160    set with TaoSetJacobianRoutine().
161 
162    Collective on TaoSolver
163 
164    Input Parameters:
165 +  solver - the TaoSolver solver context
166 -  xx - input vector
167 
168    Output Parameters:
169 +  H - Jacobian matrix
170 .  Hpre - Preconditioning matrix
171 -  flag - flag indicating matrix structure (SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, or SAME_PRECONDITIONER)
172 
173    Notes:
174    Most users should not need to explicitly call this routine, as it
175    is used internally within the minimization solvers.
176 
177    TaoComputeJacobian() is typically used within minimization
178    implementations, so most users would not generally call this routine
179    themselves.
180 
181    Level: developer
182 
183 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobian()
184 
185 @*/
186 PetscErrorCode TaoComputeJacobian(TaoSolver tao, Vec X, Mat *J, Mat *Jpre, MatStructure *flg)
187 {
188     PetscErrorCode ierr;
189     PetscFunctionBegin;
190     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
191     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
192     PetscValidPointer(flg,5);
193     PetscCheckSameComm(tao,1,X,2);
194 
195     if (!tao->ops->computejacobian) {
196 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
197     }
198     *flg = DIFFERENT_NONZERO_PATTERN;
199     ++tao->njac;
200     ierr = PetscLogEventBegin(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
201     PetscStackPush("TaoSolver user Jacobian function");
202     CHKMEMQ;
203     ierr = (*tao->ops->computejacobian)(tao,X,J,Jpre,flg,tao->user_jacP); CHKERRQ(ierr);
204     CHKMEMQ;
205     PetscStackPop;
206     ierr = PetscLogEventEnd(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
207 
208     PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "TaoComputeJacobianState"
213 /*@C
214    TaoComputeJacobianState - Computes the Jacobian matrix that has been
215    set with TaoSetJacobianStateRoutine().
216 
217    Collective on TaoSolver
218 
219    Input Parameters:
220 +  solver - the TaoSolver solver context
221 -  xx - input vector
222 
223    Output Parameters:
224 +  H - Jacobian matrix
225 .  Hpre - Preconditioning matrix
226 -  flag - flag indicating matrix structure (SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, or SAME_PRECONDITIONER)
227 
228    Notes:
229    Most users should not need to explicitly call this routine, as it
230    is used internally within the minimization solvers.
231 
232    TaoComputeJacobianState() is typically used within minimization
233    implementations, so most users would not generally call this routine
234    themselves.
235 
236    Level: developer
237 
238 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
239 
240 @*/
241 PetscErrorCode TaoComputeJacobianState(TaoSolver tao, Vec X, Mat *J, Mat *Jpre, Mat *Jinv, MatStructure *flg)
242 {
243     PetscErrorCode ierr;
244     PetscFunctionBegin;
245     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
246     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
247     PetscValidPointer(flg,5);
248     PetscCheckSameComm(tao,1,X,2);
249 
250     if (!tao->ops->computejacobianstate) {
251 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
252     }
253     *flg = DIFFERENT_NONZERO_PATTERN;
254     ++tao->njac_state;
255     ierr = PetscLogEventBegin(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
256     PetscStackPush("TaoSolver user Jacobian(state) function");
257     CHKMEMQ;
258     ierr = (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,flg,tao->user_jac_stateP); CHKERRQ(ierr);
259     CHKMEMQ;
260     PetscStackPop;
261     ierr = PetscLogEventEnd(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
262 
263     PetscFunctionReturn(0);
264 }
265 
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "TaoComputeJacobianDesign"
269 /*@C
270    TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
271    set with TaoSetJacobianDesignRoutine().
272 
273    Collective on TaoSolver
274 
275    Input Parameters:
276 +  solver - the TaoSolver solver context
277 -  xx - input vector
278 
279    Output Parameters:
280 .  H - Jacobian matrix
281 
282    Notes:
283    Most users should not need to explicitly call this routine, as it
284    is used internally within the minimization solvers.
285 
286    TaoComputeJacobianDesign() is typically used within minimization
287    implementations, so most users would not generally call this routine
288    themselves.
289 
290    Level: developer
291 
292 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
293 
294 @*/
295 PetscErrorCode TaoComputeJacobianDesign(TaoSolver tao, Vec X, Mat *J)
296 {
297     PetscErrorCode ierr;
298     PetscFunctionBegin;
299     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
300     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
301     PetscCheckSameComm(tao,1,X,2);
302 
303     if (!tao->ops->computejacobiandesign) {
304 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
305     }
306     ++tao->njac_design;
307     ierr = PetscLogEventBegin(TaoSolver_JacobianEval,tao,X,*J,PETSC_NULL); CHKERRQ(ierr);
308     PetscStackPush("TaoSolver user Jacobian(design) function");
309     CHKMEMQ;
310     ierr = (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP); CHKERRQ(ierr);
311     CHKMEMQ;
312     PetscStackPop;
313     ierr = PetscLogEventEnd(TaoSolver_JacobianEval,tao,X,*J,PETSC_NULL); CHKERRQ(ierr);
314 
315     PetscFunctionReturn(0);
316 }
317 
318 
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "TaoSetJacobianRoutine"
322 /*@C
323    TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
324 
325    Logically collective on TaoSolver
326 
327    Input Parameters:
328 +  tao - the TaoSolver context
329 .  J - Matrix used for the jacobian
330 .  Jpre - Matrix that will be used operated on by preconditioner, can be same as J
331 .  jac - Jacobian evaluation routine
332 -  ctx - [optional] user-defined context for private data for the
333          Jacobian evaluation routine (may be PETSC_NULL)
334 
335    Calling sequence of jac:
336 $    jac (TaoSolver tao,Vec x,Mat *J,Mat *Jpre,MatStructure *flag,void *ctx);
337 
338 +  tao - the TaoSolver  context
339 .  x - input vector
340 .  J - Jacobian matrix
341 .  Jpre - preconditioner matrix, usually the same as J
342 .  flag - flag indicating information about the preconditioner matrix
343    structure (see below)
344 -  ctx - [optional] user-defined Jacobian context
345 
346 
347    Notes:
348 
349    The function jac() takes Mat * as the matrix arguments rather than Mat.
350    This allows the Jacobian evaluation routine to replace A and/or B with a
351    completely new new matrix structure (not just different matrix elements)
352    when appropriate, for instance, if the nonzero structure is changing
353    throughout the global iterations.
354 
355    The flag can be used to eliminate unnecessary work in the preconditioner
356    during the repeated solution of linear systems of the same size.  The
357    available options are
358 $    SAME_PRECONDITIONER -
359 $      Jpre is identical during successive linear solves.
360 $      This option is intended for folks who are using
361 $      different Amat and Pmat matrices and want to reuse the
362 $      same preconditioner matrix.  For example, this option
363 $      saves work by not recomputing incomplete factorization
364 $      for ILU/ICC preconditioners.
365 $    SAME_NONZERO_PATTERN -
366 $      Jpre has the same nonzero structure during
367 $      successive linear solves.
368 $    DIFFERENT_NONZERO_PATTERN -
369 $      Jpre does not have the same nonzero structure.
370 
371    Caution:
372    If you specify SAME_NONZERO_PATTERN, the software believes your assertion
373    and does not check the structure of the matrix.  If you erroneously
374    claim that the structure is the same when it actually is not, the new
375    preconditioner will not function correctly.  Thus, use this optimization
376    feature carefully!
377 
378    If in doubt about whether your preconditioner matrix has changed
379    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
380 
381    Level: intermediate
382 
383 @*/
384 PetscErrorCode TaoSetJacobianRoutine(TaoSolver tao, Mat J, Mat Jpre, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, Mat *, MatStructure *, void*), void *ctx)
385 {
386     PetscErrorCode ierr;
387     PetscFunctionBegin;
388     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
389     if (J) {
390 	PetscValidHeaderSpecific(J,MAT_CLASSID,2);
391 	PetscCheckSameComm(tao,1,J,2);
392     }
393     if (Jpre) {
394 	PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
395 	PetscCheckSameComm(tao,1,Jpre,3);
396     }
397     if (ctx) {
398 	tao->user_jacP = ctx;
399     }
400     if (func) {
401 	tao->ops->computejacobian = func;
402     }
403 
404 
405     if (J) {
406 	ierr = PetscObjectReference((PetscObject)J); CHKERRQ(ierr);
407 	if (tao->jacobian) {   ierr = MatDestroy(&tao->jacobian); CHKERRQ(ierr);}
408 	tao->jacobian = J;
409     }
410     if (Jpre) {
411 	ierr = PetscObjectReference((PetscObject)Jpre); CHKERRQ(ierr);
412 	if (tao->jacobian_pre) { ierr = MatDestroy(&tao->jacobian_pre); CHKERRQ(ierr);}
413 	tao->jacobian_pre=Jpre;
414     }
415     PetscFunctionReturn(0);
416 }
417 
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "TaoSetJacobianStateRoutine"
421 /*@C
422    TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
423    (and its inverse) of the constraint function with respect to the state variables.
424    Used only for pde-constrained optimization.
425 
426    Logically collective on TaoSolver
427 
428    Input Parameters:
429 +  tao - the TaoSolver context
430 .  J - Matrix used for the jacobian
431 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.  Only used if Jinv is PETSC_NULL
432 .  Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use PETSC_NULL to default to PETSc KSP solvers to apply the inverse.
433 .  jac - Jacobian evaluation routine
434 -  ctx - [optional] user-defined context for private data for the
435          Jacobian evaluation routine (may be PETSC_NULL)
436 
437    Calling sequence of jac:
438 $    jac (TaoSolver tao,Vec x,Mat *J,Mat *Jpre,MatStructure *flag,void *ctx);
439 
440 +  tao - the TaoSolver  context
441 .  x - input vector
442 .  J - Jacobian matrix
443 .  Jpre - preconditioner matrix, usually the same as J
444 .  Jinv - inverse of J
445 .  flag - flag indicating information about the preconditioner matrix
446    structure (see below)
447 -  ctx - [optional] user-defined Jacobian context
448 
449 
450    Notes:
451    Because of the structure of the jacobian matrix,
452    It may be more efficient for a pde-constrained application to provide
453    its own Jinv matrix.
454 
455    The function jac() takes Mat * as the matrix arguments rather than Mat.
456    This allows the Jacobian evaluation routine to replace A and/or B with a
457    completely new new maitrix structure (not just different matrix elements)
458    when appropriate, for instance, if the nonzero structure is changing
459    throughout the global iterations.
460 
461    The flag can be used to eliminate unnecessary work in the preconditioner
462    during the repeated solution of linear systems of the same size.  The
463    available options are
464 $    SAME_PRECONDITIONER -
465 $      Jpre is identical during successive linear solves.
466 $      This option is intended for folks who are using
467 $      different Amat and Pmat matrices and want to reuse the
468 $      same preconditioner matrix.  For example, this option
469 $      saves work by not recomputing incomplete factorization
470 $      for ILU/ICC preconditioners.
471 $    SAME_NONZERO_PATTERN -
472 $      Jpre has the same nonzero structure during
473 $      successive linear solves.
474 $    DIFFERENT_NONZERO_PATTERN -
475 $      Jpre does not have the same nonzero structure.
476 
477    Caution:
478    If you specify SAME_NONZERO_PATTERN, the software believes your assertion
479    and does not check the structure of the matrix.  If you erroneously
480    claim that the structure is the same when it actually is not, the new
481    preconditioner will not function correctly.  Thus, use this optimization
482    feature carefully!
483 
484    If in doubt about whether your preconditioner matrix has changed
485    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
486 
487    Level: intermediate
488 .seealse: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
489 @*/
490 PetscErrorCode TaoSetJacobianStateRoutine(TaoSolver tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, Mat *, Mat *, MatStructure *, void*), void *ctx)
491 {
492     PetscErrorCode ierr;
493     PetscFunctionBegin;
494     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
495     if (J) {
496 	PetscValidHeaderSpecific(J,MAT_CLASSID,2);
497 	PetscCheckSameComm(tao,1,J,2);
498     }
499     if (Jpre) {
500 	PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
501 	PetscCheckSameComm(tao,1,Jpre,3);
502     }
503     if (Jinv) {
504       PetscValidHeaderSpecific(Jinv,MAT_CLASSID,4);
505       PetscCheckSameComm(tao,1,Jinv,4);
506     }
507     if (ctx) {
508 	tao->user_jac_stateP = ctx;
509     }
510     if (func) {
511 	tao->ops->computejacobianstate = func;
512     }
513 
514 
515     if (J) {
516 	ierr = PetscObjectReference((PetscObject)J); CHKERRQ(ierr);
517 	if (tao->jacobian_state) {   ierr = MatDestroy(&tao->jacobian_state); CHKERRQ(ierr);}
518 	tao->jacobian_state = J;
519     }
520     if (Jpre) {
521 	ierr = PetscObjectReference((PetscObject)Jpre); CHKERRQ(ierr);
522 	if (tao->jacobian_state_pre) { ierr = MatDestroy(&tao->jacobian_state_pre); CHKERRQ(ierr);}
523 	tao->jacobian_state_pre=Jpre;
524     }
525     if (Jinv) {
526       ierr = PetscObjectReference((PetscObject)Jinv); CHKERRQ(ierr);
527       if (tao->jacobian_state_inv) {ierr = MatDestroy(&tao->jacobian_state_inv); CHKERRQ(ierr);}
528       tao->jacobian_state_inv=Jinv;
529     }
530     PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "TaoSetJacobianDesignRoutine"
535 /*@C
536    TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
537    the constraint function with respect to the design variables.  Used only for
538    pde-constrained optimization.
539 
540    Logically collective on TaoSolver
541 
542    Input Parameters:
543 +  tao - the TaoSolver context
544 .  J - Matrix used for the jacobian
545 .  jac - Jacobian evaluation routine
546 -  ctx - [optional] user-defined context for private data for the
547          Jacobian evaluation routine (may be PETSC_NULL)
548 
549    Calling sequence of jac:
550 $    jac (TaoSolver tao,Vec x,Mat *J,void *ctx);
551 
552 +  tao - the TaoSolver  context
553 .  x - input vector
554 .  J - Jacobian matrix
555 -  ctx - [optional] user-defined Jacobian context
556 
557 
558    Notes:
559 
560    The function jac() takes Mat * as the matrix arguments rather than Mat.
561    This allows the Jacobian evaluation routine to replace A and/or B with a
562    completely new new matrix structure (not just different matrix elements)
563    when appropriate, for instance, if the nonzero structure is changing
564    throughout the global iterations.
565 
566    Level: intermediate
567 .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS()
568 @*/
569 PetscErrorCode TaoSetJacobianDesignRoutine(TaoSolver tao, Mat J, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, void*), void *ctx)
570 {
571     PetscErrorCode ierr;
572     PetscFunctionBegin;
573     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
574     if (J) {
575 	PetscValidHeaderSpecific(J,MAT_CLASSID,2);
576 	PetscCheckSameComm(tao,1,J,2);
577     }
578     if (ctx) {
579 	tao->user_jac_designP = ctx;
580     }
581     if (func) {
582 	tao->ops->computejacobiandesign = func;
583     }
584 
585 
586     if (J) {
587 	ierr = PetscObjectReference((PetscObject)J); CHKERRQ(ierr);
588 	if (tao->jacobian_design) {   ierr = MatDestroy(&tao->jacobian_design); CHKERRQ(ierr);}
589 	tao->jacobian_design = J;
590     }
591     PetscFunctionReturn(0);
592 }
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "TaoSetStateDesignIS"
596 /*@
597    TaoSetStateDesignIS - Indicate to the TaoSolver which variables in the
598    solution vector are state variables and which are design.  Only applies to
599    pde-constrained optimization.
600 
601    Logically Collective on TaoSolver
602 
603    Input Parameters:
604 +  tao - The TaoSolver context
605 .  s_is - the index set corresponding to the state variables
606 -  d_is - the index set corresponding to the design variables
607 
608    Level: intermediate
609 
610 .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine()
611 @*/
612 PetscErrorCode TaoSetStateDesignIS(TaoSolver tao, IS s_is, IS d_is)
613 {
614   PetscErrorCode ierr;
615   if (tao->state_is) {
616     ierr = PetscObjectDereference((PetscObject)(tao->state_is)); CHKERRQ(ierr);
617   }
618   if (tao->design_is) {
619     ierr = PetscObjectDereference((PetscObject)(tao->design_is)); CHKERRQ(ierr);
620   }
621   tao->state_is = s_is;
622   tao->design_is = d_is;
623   if (s_is) {
624     ierr = PetscObjectReference((PetscObject)(tao->state_is)); CHKERRQ(ierr);
625   }
626   if (d_is) {
627     ierr = PetscObjectReference((PetscObject)(tao->design_is)); CHKERRQ(ierr);
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "TaoComputeJacobianEquality"
635 /*@C
636    TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
637    set with TaoSetJacobianEqualityRoutine().
638 
639    Collective on TaoSolver
640 
641    Input Parameters:
642 +  solver - the TaoSolver solver context
643 -  xx - input vector
644 
645    Output Parameters:
646 +  H - Jacobian matrix
647 .  Hpre - Preconditioning matrix
648 -  flag - flag indicating matrix structure (SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, or SAME_PRECONDITIONER)
649 
650    Notes:
651    Most users should not need to explicitly call this routine, as it
652    is used internally within the minimization solvers.
653 
654    Level: developer
655 
656 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
657 
658 @*/
659 PetscErrorCode TaoComputeJacobianEquality(TaoSolver tao, Vec X, Mat *J, Mat *Jpre, MatStructure *flg)
660 {
661     PetscErrorCode ierr;
662     PetscFunctionBegin;
663     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
664     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
665     PetscValidPointer(flg,5);
666     PetscCheckSameComm(tao,1,X,2);
667 
668     if (!tao->ops->computejacobianequality) {
669 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
670     }
671     *flg = DIFFERENT_NONZERO_PATTERN;
672     ++tao->njac_equality;
673     ierr = PetscLogEventBegin(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
674     PetscStackPush("TaoSolver user Jacobian(equality) function");
675     CHKMEMQ;
676     ierr = (*tao->ops->computejacobianequality)(tao,X,J,Jpre,flg,tao->user_jac_equalityP); CHKERRQ(ierr);
677     CHKMEMQ;
678     PetscStackPop;
679     ierr = PetscLogEventEnd(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
680 
681     PetscFunctionReturn(0);
682 }
683 
684 #undef __FUNCT__
685 #define __FUNCT__ "TaoComputeJacobianInequality"
686 /*@C
687    TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
688    set with TaoSetJacobianInequalityRoutine().
689 
690    Collective on TaoSolver
691 
692    Input Parameters:
693 +  solver - the TaoSolver solver context
694 -  xx - input vector
695 
696    Output Parameters:
697 +  H - Jacobian matrix
698 .  Hpre - Preconditioning matrix
699 -  flag - flag indicating matrix structure (SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN, or SAME_PRECONDITIONER)
700 
701    Notes:
702    Most users should not need to explicitly call this routine, as it
703    is used internally within the minimization solvers.
704 
705    Level: developer
706 
707 .seealso:  TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
708 
709 @*/
710 PetscErrorCode TaoComputeJacobianInequality(TaoSolver tao, Vec X, Mat *J, Mat *Jpre, MatStructure *flg)
711 {
712     PetscErrorCode ierr;
713     PetscFunctionBegin;
714     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
715     PetscValidHeaderSpecific(X, VEC_CLASSID,2);
716     PetscValidPointer(flg,5);
717     PetscCheckSameComm(tao,1,X,2);
718 
719     if (!tao->ops->computejacobianinequality) {
720 	SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
721     }
722     *flg = DIFFERENT_NONZERO_PATTERN;
723     ++tao->njac_inequality;
724     ierr = PetscLogEventBegin(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
725     PetscStackPush("TaoSolver user Jacobian(inequality) function");
726     CHKMEMQ;
727     ierr = (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,flg,tao->user_jac_inequalityP); CHKERRQ(ierr);
728     CHKMEMQ;
729     PetscStackPop;
730     ierr = PetscLogEventEnd(TaoSolver_JacobianEval,tao,X,*J,*Jpre); CHKERRQ(ierr);
731 
732     PetscFunctionReturn(0);
733 }
734 
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "TaoSetJacobianEqualityRoutine"
738 /*@C
739    TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
740    (and its inverse) of the constraint function with respect to the equality variables.
741    Used only for pde-constrained optimization.
742 
743    Logically collective on TaoSolver
744 
745    Input Parameters:
746 +  tao - the TaoSolver context
747 .  J - Matrix used for the jacobian
748 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
749 .  jac - Jacobian evaluation routine
750 -  ctx - [optional] user-defined context for private data for the
751          Jacobian evaluation routine (may be PETSC_NULL)
752 
753    Calling sequence of jac:
754 $    jac (TaoSolver tao,Vec x,Mat *J,Mat *Jpre,MatStructure *flag,void *ctx);
755 
756 +  tao - the TaoSolver  context
757 .  x - input vector
758 .  J - Jacobian matrix
759 .  Jpre - preconditioner matrix, usually the same as J
760 .  flag - flag indicating information about the preconditioner matrix
761    structure (see below)
762 -  ctx - [optional] user-defined Jacobian context
763 
764 
765    Notes:
766    Because of the structure of the jacobian matrix,
767    It may be more efficient for a pde-constrained application to provide
768    its own Jinv matrix.
769 
770    The function jac() takes Mat * as the matrix arguments rather than Mat.
771    This allows the Jacobian evaluation routine to replace A and/or B with a
772    completely new new maitrix structure (not just different matrix elements)
773    when appropriate, for instance, if the nonzero structure is changing
774    throughout the global iterations.
775 
776    The flag can be used to eliminate unnecessary work in the preconditioner
777    during the repeated solution of linear systems of the same size.  The
778    available options are
779 $    SAME_PRECONDITIONER -
780 $      Jpre is identical during successive linear solves.
781 $      This option is intended for folks who are using
782 $      different Amat and Pmat matrices and want to reuse the
783 $      same preconditioner matrix.  For example, this option
784 $      saves work by not recomputing incomplete factorization
785 $      for ILU/ICC preconditioners.
786 $    SAME_NONZERO_PATTERN -
787 $      Jpre has the same nonzero structure during
788 $      successive linear solves.
789 $    DIFFERENT_NONZERO_PATTERN -
790 $      Jpre does not have the same nonzero structure.
791 
792    Caution:
793    If you specify SAME_NONZERO_PATTERN, the software believes your assertion
794    and does not check the structure of the matrix.  If you erroneously
795    claim that the structure is the same when it actually is not, the new
796    preconditioner will not function correctly.  Thus, use this optimization
797    feature carefully!
798 
799    If in doubt about whether your preconditioner matrix has changed
800    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
801 
802    Level: intermediate
803 .seealse: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS()
804 @*/
805 PetscErrorCode TaoSetJacobianEqualityRoutine(TaoSolver tao, Mat J, Mat Jpre, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, Mat *, MatStructure *, void*), void *ctx)
806 {
807     PetscErrorCode ierr;
808     PetscFunctionBegin;
809     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
810     if (J) {
811 	PetscValidHeaderSpecific(J,MAT_CLASSID,2);
812 	PetscCheckSameComm(tao,1,J,2);
813     }
814     if (Jpre) {
815 	PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
816 	PetscCheckSameComm(tao,1,Jpre,3);
817     }
818     if (ctx) {
819 	tao->user_jac_equalityP = ctx;
820     }
821     if (func) {
822 	tao->ops->computejacobianequality = func;
823     }
824 
825 
826     if (J) {
827 	ierr = PetscObjectReference((PetscObject)J); CHKERRQ(ierr);
828 	if (tao->jacobian_equality) {   ierr = MatDestroy(&tao->jacobian_equality); CHKERRQ(ierr);}
829 	tao->jacobian_equality = J;
830     }
831     if (Jpre) {
832 	ierr = PetscObjectReference((PetscObject)Jpre); CHKERRQ(ierr);
833 	if (tao->jacobian_equality_pre) { ierr = MatDestroy(&tao->jacobian_equality_pre); CHKERRQ(ierr);}
834 	tao->jacobian_equality_pre=Jpre;
835     }
836     PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "TaoSetJacobianInequalityRoutine"
841 /*@C
842    TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
843    (and its inverse) of the constraint function with respect to the inequality variables.
844    Used only for pde-constrained optimization.
845 
846    Logically collective on TaoSolver
847 
848    Input Parameters:
849 +  tao - the TaoSolver context
850 .  J - Matrix used for the jacobian
851 .  Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
852 .  jac - Jacobian evaluation routine
853 -  ctx - [optional] user-defined context for private data for the
854          Jacobian evaluation routine (may be PETSC_NULL)
855 
856    Calling sequence of jac:
857 $    jac (TaoSolver tao,Vec x,Mat *J,Mat *Jpre,MatStructure *flag,void *ctx);
858 
859 +  tao - the TaoSolver  context
860 .  x - input vector
861 .  J - Jacobian matrix
862 .  Jpre - preconditioner matrix, usually the same as J
863 .  flag - flag indicating information about the preconditioner matrix
864    structure (see below)
865 -  ctx - [optional] user-defined Jacobian context
866 
867 
868    Notes:
869    Because of the structure of the jacobian matrix,
870    It may be more efficient for a pde-constrained application to provide
871    its own Jinv matrix.
872 
873    The function jac() takes Mat * as the matrix arguments rather than Mat.
874    This allows the Jacobian evaluation routine to replace A and/or B with a
875    completely new new maitrix structure (not just different matrix elements)
876    when appropriate, for instance, if the nonzero structure is changing
877    throughout the global iterations.
878 
879    The flag can be used to eliminate unnecessary work in the preconditioner
880    during the repeated solution of linear systems of the same size.  The
881    available options are
882 $    SAME_PRECONDITIONER -
883 $      Jpre is identical during successive linear solves.
884 $      This option is intended for folks who are using
885 $      different Amat and Pmat matrices and want to reuse the
886 $      same preconditioner matrix.  For example, this option
887 $      saves work by not recomputing incomplete factorization
888 $      for ILU/ICC preconditioners.
889 $    SAME_NONZERO_PATTERN -
890 $      Jpre has the same nonzero structure during
891 $      successive linear solves.
892 $    DIFFERENT_NONZERO_PATTERN -
893 $      Jpre does not have the same nonzero structure.
894 
895    Caution:
896    If you specify SAME_NONZERO_PATTERN, the software believes your assertion
897    and does not check the structure of the matrix.  If you erroneously
898    claim that the structure is the same when it actually is not, the new
899    preconditioner will not function correctly.  Thus, use this optimization
900    feature carefully!
901 
902    If in doubt about whether your preconditioner matrix has changed
903    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
904 
905    Level: intermediate
906 .seealse: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS()
907 @*/
908 PetscErrorCode TaoSetJacobianInequalityRoutine(TaoSolver tao, Mat J, Mat Jpre, PetscErrorCode (*func)(TaoSolver, Vec, Mat*, Mat *, MatStructure *, void*), void *ctx)
909 {
910     PetscErrorCode ierr;
911     PetscFunctionBegin;
912     PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
913     if (J) {
914 	PetscValidHeaderSpecific(J,MAT_CLASSID,2);
915 	PetscCheckSameComm(tao,1,J,2);
916     }
917     if (Jpre) {
918 	PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3);
919 	PetscCheckSameComm(tao,1,Jpre,3);
920     }
921     if (ctx) {
922 	tao->user_jac_inequalityP = ctx;
923     }
924     if (func) {
925 	tao->ops->computejacobianinequality = func;
926     }
927 
928 
929     if (J) {
930 	ierr = PetscObjectReference((PetscObject)J); CHKERRQ(ierr);
931 	if (tao->jacobian_inequality) {   ierr = MatDestroy(&tao->jacobian_inequality); CHKERRQ(ierr);}
932 	tao->jacobian_inequality = J;
933     }
934     if (Jpre) {
935 	ierr = PetscObjectReference((PetscObject)Jpre); CHKERRQ(ierr);
936 	if (tao->jacobian_inequality_pre) { ierr = MatDestroy(&tao->jacobian_inequality_pre); CHKERRQ(ierr);}
937 	tao->jacobian_inequality_pre=Jpre;
938     }
939     PetscFunctionReturn(0);
940 }
941