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