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