Actual source code: taosolver_hj.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/taoimpl.h> /*I "petsctao.h" I*/
5: /*@C
6: TaoSetHessianRoutine - Sets the function to compute the Hessian as well as the location to store the matrix.
8: Logically collective on Tao
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)
18: Calling sequence of hess:
19: $ hess (Tao tao,Vec x,Mat H,Mat Hpre,void *ctx);
21: + tao - the Tao context
22: . x - input vector
23: . H - Hessian matrix
24: . Hpre - preconditioner matrix, usually the same as H
25: - ctx - [optional] user-defined Hessian context
27: Level: beginner
29: @*/
30: PetscErrorCode TaoSetHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
31: {
36: if (H) {
39: }
40: if (Hpre) {
43: }
44: if (ctx) {
45: tao->user_hessP = ctx;
46: }
47: if (func) {
48: tao->ops->computehessian = func;
49: }
50: if (H) {
51: PetscObjectReference((PetscObject)H);
52: MatDestroy(&tao->hessian);
53: tao->hessian = H;
54: }
55: if (Hpre) {
56: PetscObjectReference((PetscObject)Hpre);
57: MatDestroy(&tao->hessian_pre);
58: tao->hessian_pre = Hpre;
59: }
60: return(0);
61: }
65: /*@C
66: TaoComputeHessian - Computes the Hessian matrix that has been
67: set with TaoSetHessianRoutine().
69: Collective on Tao
71: Input Parameters:
72: + solver - the Tao solver context
73: - xx - input vector
75: Output Parameters:
76: + H - Hessian matrix
77: - Hpre - Preconditioning matrix
79: Notes:
80: Most users should not need to explicitly call this routine, as it
81: is used internally within the minimization solvers.
83: TaoComputeHessian() is typically used within minimization
84: implementations, so most users would not generally call this routine
85: themselves.
87: Level: developer
89: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetHessian()
91: @*/
92: PetscErrorCode TaoComputeHessian(Tao tao, Vec X, Mat H, Mat Hpre)
93: {
101: if (!tao->ops->computehessian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetHessian() first");
102: ++tao->nhess;
103: PetscLogEventBegin(Tao_HessianEval,tao,X,H,Hpre);
104: PetscStackPush("Tao user Hessian function");
105: (*tao->ops->computehessian)(tao,X,H,Hpre,tao->user_hessP);
106: PetscStackPop;
107: PetscLogEventEnd(Tao_HessianEval,tao,X,H,Hpre);
108: return(0);
109: }
113: /*@C
114: TaoComputeJacobian - Computes the Jacobian matrix that has been
115: set with TaoSetJacobianRoutine().
117: Collective on Tao
119: Input Parameters:
120: + solver - the Tao solver context
121: - xx - input vector
123: Output Parameters:
124: + H - Jacobian matrix
125: - Hpre - Preconditioning matrix
127: Notes:
128: Most users should not need to explicitly call this routine, as it
129: is used internally within the minimization solvers.
131: TaoComputeJacobian() is typically used within minimization
132: implementations, so most users would not generally call this routine
133: themselves.
135: Level: developer
137: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobian()
139: @*/
140: PetscErrorCode TaoComputeJacobian(Tao tao, Vec X, Mat J, Mat Jpre)
141: {
149: if (!tao->ops->computejacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobian() first");
150: ++tao->njac;
151: PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);
152: PetscStackPush("Tao user Jacobian function");
153: (*tao->ops->computejacobian)(tao,X,J,Jpre,tao->user_jacP);
154: PetscStackPop;
155: PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);
156: return(0);
157: }
161: /*@C
162: TaoComputeJacobianState - Computes the Jacobian matrix that has been
163: set with TaoSetJacobianStateRoutine().
165: Collective on Tao
167: Input Parameters:
168: + solver - the Tao solver context
169: - xx - input vector
171: Output Parameters:
172: + H - Jacobian matrix
173: - Hpre - Preconditioning matrix
175: Notes:
176: Most users should not need to explicitly call this routine, as it
177: is used internally within the minimization solvers.
179: TaoComputeJacobianState() is typically used within minimization
180: implementations, so most users would not generally call this routine
181: themselves.
183: Level: developer
185: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
187: @*/
188: PetscErrorCode TaoComputeJacobianState(Tao tao, Vec X, Mat J, Mat Jpre, Mat Jinv)
189: {
197: if (!tao->ops->computejacobianstate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianState() first");
198: ++tao->njac_state;
199: PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);
200: PetscStackPush("Tao user Jacobian(state) function");
201: (*tao->ops->computejacobianstate)(tao,X,J,Jpre,Jinv,tao->user_jac_stateP);
202: PetscStackPop;
203: PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);
204: return(0);
205: }
209: /*@C
210: TaoComputeJacobianDesign - Computes the Jacobian matrix that has been
211: set with TaoSetJacobianDesignRoutine().
213: Collective on Tao
215: Input Parameters:
216: + solver - the Tao solver context
217: - xx - input vector
219: Output Parameters:
220: . H - Jacobian matrix
222: Notes:
223: Most users should not need to explicitly call this routine, as it
224: is used internally within the minimization solvers.
226: TaoComputeJacobianDesign() is typically used within minimization
227: implementations, so most users would not generally call this routine
228: themselves.
230: Level: developer
232: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianDesignRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
234: @*/
235: PetscErrorCode TaoComputeJacobianDesign(Tao tao, Vec X, Mat J)
236: {
244: if (!tao->ops->computejacobiandesign) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianDesign() first");
245: ++tao->njac_design;
246: PetscLogEventBegin(Tao_JacobianEval,tao,X,J,NULL);
247: PetscStackPush("Tao user Jacobian(design) function");
248: (*tao->ops->computejacobiandesign)(tao,X,J,tao->user_jac_designP);
249: PetscStackPop;
250: PetscLogEventEnd(Tao_JacobianEval,tao,X,J,NULL);
251: return(0);
252: }
256: /*@C
257: TaoSetJacobianRoutine - Sets the function to compute the Jacobian as well as the location to store the matrix.
259: Logically collective on Tao
261: Input Parameters:
262: + tao - the Tao context
263: . J - Matrix used for the jacobian
264: . Jpre - Matrix that will be used operated on by preconditioner, can be same as J
265: . jac - Jacobian evaluation routine
266: - ctx - [optional] user-defined context for private data for the
267: Jacobian evaluation routine (may be NULL)
269: Calling sequence of jac:
270: $ jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
272: + tao - the Tao context
273: . x - input vector
274: . J - Jacobian matrix
275: . Jpre - preconditioner matrix, usually the same as J
276: - ctx - [optional] user-defined Jacobian context
278: Level: intermediate
280: @*/
281: PetscErrorCode TaoSetJacobianRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx)
282: {
286: if (J) {
289: }
290: if (Jpre) {
293: }
294: if (ctx) {
295: tao->user_jacP = ctx;
296: }
297: if (func) {
298: tao->ops->computejacobian = func;
299: }
300: if (J) {
301: PetscObjectReference((PetscObject)J);
302: MatDestroy(&tao->jacobian);
303: tao->jacobian = J;
304: }
305: if (Jpre) {
306: PetscObjectReference((PetscObject)Jpre);
307: MatDestroy(&tao->jacobian_pre);
308: tao->jacobian_pre=Jpre;
309: }
310: return(0);
311: }
315: /*@C
316: TaoSetJacobianStateRoutine - Sets the function to compute the Jacobian
317: (and its inverse) of the constraint function with respect to the state variables.
318: Used only for pde-constrained optimization.
320: Logically collective on Tao
322: Input Parameters:
323: + tao - the Tao context
324: . J - Matrix used for the jacobian
325: . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J. Only used if Jinv is NULL
326: . Jinv - [optional] Matrix used to apply the inverse of the state jacobian. Use NULL to default to PETSc KSP solvers to apply the inverse.
327: . jac - Jacobian evaluation routine
328: - ctx - [optional] user-defined context for private data for the
329: Jacobian evaluation routine (may be NULL)
331: Calling sequence of jac:
332: $ jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
334: + tao - the Tao context
335: . x - input vector
336: . J - Jacobian matrix
337: . Jpre - preconditioner matrix, usually the same as J
338: . Jinv - inverse of J
339: - ctx - [optional] user-defined Jacobian context
341: Level: intermediate
342: .seealse: TaoComputeJacobianState(), TaoSetJacobianDesignRoutine(), TaoSetStateDesignIS()
343: @*/
344: PetscErrorCode TaoSetJacobianStateRoutine(Tao tao, Mat J, Mat Jpre, Mat Jinv, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, Mat,void*), void *ctx)
345: {
349: if (J) {
352: }
353: if (Jpre) {
356: }
357: if (Jinv) {
360: }
361: if (ctx) {
362: tao->user_jac_stateP = ctx;
363: }
364: if (func) {
365: tao->ops->computejacobianstate = func;
366: }
367: if (J) {
368: PetscObjectReference((PetscObject)J);
369: MatDestroy(&tao->jacobian_state);
370: tao->jacobian_state = J;
371: }
372: if (Jpre) {
373: PetscObjectReference((PetscObject)Jpre);
374: MatDestroy(&tao->jacobian_state_pre);
375: tao->jacobian_state_pre=Jpre;
376: }
377: if (Jinv) {
378: PetscObjectReference((PetscObject)Jinv);
379: MatDestroy(&tao->jacobian_state_inv);
380: tao->jacobian_state_inv=Jinv;
381: }
382: return(0);
383: }
387: /*@C
388: TaoSetJacobianDesignRoutine - Sets the function to compute the Jacobian of
389: the constraint function with respect to the design variables. Used only for
390: pde-constrained optimization.
392: Logically collective on Tao
394: Input Parameters:
395: + tao - the Tao context
396: . J - Matrix used for the jacobian
397: . jac - Jacobian evaluation routine
398: - ctx - [optional] user-defined context for private data for the
399: Jacobian evaluation routine (may be NULL)
401: Calling sequence of jac:
402: $ jac (Tao tao,Vec x,Mat *J,void *ctx);
404: + tao - the Tao context
405: . x - input vector
406: . J - Jacobian matrix
407: - ctx - [optional] user-defined Jacobian context
410: Notes:
412: The function jac() takes Mat * as the matrix arguments rather than Mat.
413: This allows the Jacobian evaluation routine to replace A and/or B with a
414: completely new new matrix structure (not just different matrix elements)
415: when appropriate, for instance, if the nonzero structure is changing
416: throughout the global iterations.
418: Level: intermediate
419: .seealso: TaoComputeJacobianDesign(), TaoSetJacobianStateRoutine(), TaoSetStateDesignIS()
420: @*/
421: PetscErrorCode TaoSetJacobianDesignRoutine(Tao tao, Mat J, PetscErrorCode (*func)(Tao, Vec, Mat, void*), void *ctx)
422: {
427: if (J) {
430: }
431: if (ctx) {
432: tao->user_jac_designP = ctx;
433: }
434: if (func) {
435: tao->ops->computejacobiandesign = func;
436: }
437: if (J) {
438: PetscObjectReference((PetscObject)J);
439: MatDestroy(&tao->jacobian_design);
440: tao->jacobian_design = J;
441: }
442: return(0);
443: }
447: /*@
448: TaoSetStateDesignIS - Indicate to the Tao which variables in the
449: solution vector are state variables and which are design. Only applies to
450: pde-constrained optimization.
452: Logically Collective on Tao
454: Input Parameters:
455: + tao - The Tao context
456: . s_is - the index set corresponding to the state variables
457: - d_is - the index set corresponding to the design variables
459: Level: intermediate
461: .seealso: TaoSetJacobianStateRoutine(), TaoSetJacobianDesignRoutine()
462: @*/
463: PetscErrorCode TaoSetStateDesignIS(Tao tao, IS s_is, IS d_is)
464: {
468: PetscObjectReference((PetscObject)s_is);
469: ISDestroy(&tao->state_is);
470: tao->state_is = s_is;
471: PetscObjectReference((PetscObject)(d_is));
472: ISDestroy(&tao->design_is);
473: tao->design_is = d_is;
474: return(0);
475: }
479: /*@C
480: TaoComputeJacobianEquality - Computes the Jacobian matrix that has been
481: set with TaoSetJacobianEqualityRoutine().
483: Collective on Tao
485: Input Parameters:
486: + solver - the Tao solver context
487: - xx - input vector
489: Output Parameters:
490: + H - Jacobian matrix
491: - Hpre - Preconditioning matrix
493: Notes:
494: Most users should not need to explicitly call this routine, as it
495: is used internally within the minimization solvers.
497: Level: developer
499: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
501: @*/
502: PetscErrorCode TaoComputeJacobianEquality(Tao tao, Vec X, Mat J, Mat Jpre)
503: {
511: if (!tao->ops->computejacobianequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianEquality() first");
512: ++tao->njac_equality;
513: PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);
514: PetscStackPush("Tao user Jacobian(equality) function");
515: (*tao->ops->computejacobianequality)(tao,X,J,Jpre,tao->user_jac_equalityP);
516: PetscStackPop;
517: PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);
518: return(0);
519: }
523: /*@C
524: TaoComputeJacobianInequality - Computes the Jacobian matrix that has been
525: set with TaoSetJacobianInequalityRoutine().
527: Collective on Tao
529: Input Parameters:
530: + solver - the Tao solver context
531: - xx - input vector
533: Output Parameters:
534: + H - Jacobian matrix
535: - Hpre - Preconditioning matrix
537: Notes:
538: Most users should not need to explicitly call this routine, as it
539: is used internally within the minimization solvers.
541: Level: developer
543: .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetJacobianStateRoutine(), TaoComputeJacobianDesign(), TaoSetStateDesignIS()
545: @*/
546: PetscErrorCode TaoComputeJacobianInequality(Tao tao, Vec X, Mat J, Mat Jpre)
547: {
555: if (!tao->ops->computejacobianinequality) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetJacobianInequality() first");
556: ++tao->njac_inequality;
557: PetscLogEventBegin(Tao_JacobianEval,tao,X,J,Jpre);
558: PetscStackPush("Tao user Jacobian(inequality) function");
559: (*tao->ops->computejacobianinequality)(tao,X,J,Jpre,tao->user_jac_inequalityP);
560: PetscStackPop;
561: PetscLogEventEnd(Tao_JacobianEval,tao,X,J,Jpre);
562: return(0);
563: }
567: /*@C
568: TaoSetJacobianEqualityRoutine - Sets the function to compute the Jacobian
569: (and its inverse) of the constraint function with respect to the equality variables.
570: Used only for pde-constrained optimization.
572: Logically collective on Tao
574: Input Parameters:
575: + tao - the Tao context
576: . J - Matrix used for the jacobian
577: . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
578: . jac - Jacobian evaluation routine
579: - ctx - [optional] user-defined context for private data for the
580: Jacobian evaluation routine (may be NULL)
582: Calling sequence of jac:
583: $ jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
585: + tao - the Tao context
586: . x - input vector
587: . J - Jacobian matrix
588: . Jpre - preconditioner matrix, usually the same as J
589: - ctx - [optional] user-defined Jacobian context
591: Level: intermediate
592: .seealse: TaoComputeJacobianEquality(), TaoSetJacobianDesignRoutine(), TaoSetEqualityDesignIS()
593: @*/
594: PetscErrorCode TaoSetJacobianEqualityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
595: {
600: if (J) {
603: }
604: if (Jpre) {
607: }
608: if (ctx) {
609: tao->user_jac_equalityP = ctx;
610: }
611: if (func) {
612: tao->ops->computejacobianequality = func;
613: }
614: if (J) {
615: PetscObjectReference((PetscObject)J);
616: MatDestroy(&tao->jacobian_equality);
617: tao->jacobian_equality = J;
618: }
619: if (Jpre) {
620: PetscObjectReference((PetscObject)Jpre);
621: MatDestroy(&tao->jacobian_equality_pre);
622: tao->jacobian_equality_pre=Jpre;
623: }
624: return(0);
625: }
629: /*@C
630: TaoSetJacobianInequalityRoutine - Sets the function to compute the Jacobian
631: (and its inverse) of the constraint function with respect to the inequality variables.
632: Used only for pde-constrained optimization.
634: Logically collective on Tao
636: Input Parameters:
637: + tao - the Tao context
638: . J - Matrix used for the jacobian
639: . Jpre - Matrix that will be used operated on by PETSc preconditioner, can be same as J.
640: . jac - Jacobian evaluation routine
641: - ctx - [optional] user-defined context for private data for the
642: Jacobian evaluation routine (may be NULL)
644: Calling sequence of jac:
645: $ jac (Tao tao,Vec x,Mat *J,Mat *Jpre,void *ctx);
647: + tao - the Tao context
648: . x - input vector
649: . J - Jacobian matrix
650: . Jpre - preconditioner matrix, usually the same as J
651: - ctx - [optional] user-defined Jacobian context
653: Level: intermediate
654: .seealse: TaoComputeJacobianInequality(), TaoSetJacobianDesignRoutine(), TaoSetInequalityDesignIS()
655: @*/
656: PetscErrorCode TaoSetJacobianInequalityRoutine(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat,void*), void *ctx)
657: {
661: if (J) {
664: }
665: if (Jpre) {
668: }
669: if (ctx) {
670: tao->user_jac_inequalityP = ctx;
671: }
672: if (func) {
673: tao->ops->computejacobianinequality = func;
674: }
675: if (J) {
676: PetscObjectReference((PetscObject)J);
677: MatDestroy(&tao->jacobian_inequality);
678: tao->jacobian_inequality = J;
679: }
680: if (Jpre) {
681: PetscObjectReference((PetscObject)Jpre);
682: MatDestroy(&tao->jacobian_inequality_pre);
683: tao->jacobian_inequality_pre=Jpre;
684: }
685: return(0);
686: }