Actual source code: taoimpl.h

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: #ifndef __TAO_IMPL_H

  4: #include <petsctaolinesearch.h>
  5: #include <petsc/private/petscimpl.h>
  6: #include <petscksp.h>

  8: PETSC_EXTERN PetscErrorCode TaoRegisterAll(void);

 10: typedef struct _TaoOps *TaoOps;

 12: struct _TaoOps {
 13:   /* Methods set by application */
 14:     PetscErrorCode (*computeobjective)(Tao, Vec, PetscReal*, void*);
 15:     PetscErrorCode (*computeobjectiveandgradient)(Tao, Vec, PetscReal*, Vec, void*);
 16:     PetscErrorCode (*computegradient)(Tao, Vec, Vec, void*);
 17:     PetscErrorCode (*computehessian)(Tao, Vec, Mat, Mat,  void*);
 18:     PetscErrorCode (*computeseparableobjective)(Tao, Vec, Vec, void*);
 19:     PetscErrorCode (*computeconstraints)(Tao, Vec, Vec, void*);
 20:     PetscErrorCode (*computeinequalityconstraints)(Tao, Vec, Vec, void*);
 21:     PetscErrorCode (*computeequalityconstraints)(Tao, Vec, Vec, void*);
 22:     PetscErrorCode (*computejacobian)(Tao, Vec, Mat, Mat,  void*);
 23:     PetscErrorCode (*computejacobianstate)(Tao, Vec, Mat, Mat, Mat,  void*);
 24:     PetscErrorCode (*computejacobiandesign)(Tao, Vec, Mat, void*);
 25:     PetscErrorCode (*computejacobianinequality)(Tao, Vec, Mat, Mat,  void*);
 26:     PetscErrorCode (*computejacobianequality)(Tao, Vec, Mat, Mat,  void*);
 27:     PetscErrorCode (*computebounds)(Tao, Vec, Vec, void*);

 29:     PetscErrorCode (*convergencetest)(Tao,void*);
 30:     PetscErrorCode (*convergencedestroy)(void*);

 32:   /* Methods set by solver */
 33:     PetscErrorCode (*computedual)(Tao, Vec, Vec);
 34:     PetscErrorCode (*setup)(Tao);
 35:     PetscErrorCode (*solve)(Tao);
 36:     PetscErrorCode (*view)(Tao, PetscViewer);
 37:     PetscErrorCode (*setfromoptions)(PetscOptions*,Tao);
 38:     PetscErrorCode (*destroy)(Tao);
 39: };

 41: #define MAXTAOMONITORS 10

 43: struct _p_Tao {
 44:     PETSCHEADER(struct _TaoOps);
 45:     void *user;
 46:     void *user_objP;
 47:     void *user_objgradP;
 48:     void *user_gradP;
 49:     void *user_hessP;
 50:     void *user_sepobjP;
 51:     void *user_conP;
 52:     void *user_con_equalityP;
 53:     void *user_con_inequalityP;
 54:     void *user_jacP;
 55:     void *user_jac_equalityP;
 56:     void *user_jac_inequalityP;
 57:     void *user_jac_stateP;
 58:     void *user_jac_designP;
 59:     void *user_boundsP;

 61:     PetscErrorCode (*monitor[MAXTAOMONITORS])(Tao,void*);
 62:     PetscErrorCode (*monitordestroy[MAXTAOMONITORS])(void**);
 63:     void *monitorcontext[MAXTAOMONITORS];
 64:     PetscInt numbermonitors;
 65:     void *cnvP;
 66:     TaoConvergedReason reason;

 68:     PetscBool setupcalled;
 69:     void *data;

 71:     Vec solution;
 72:     Vec gradient;
 73:     Vec stepdirection;
 74:     Vec XL;
 75:     Vec XU;
 76:     Vec IL;
 77:     Vec IU;
 78:     Vec DI;
 79:     Vec DE;
 80:     Mat hessian;
 81:     Mat hessian_pre;
 82:     Vec sep_objective;
 83:     Vec constraints;
 84:     Vec constraints_equality;
 85:     Vec constraints_inequality;
 86:     Mat jacobian;
 87:     Mat jacobian_pre;
 88:     Mat jacobian_inequality;
 89:     Mat jacobian_inequality_pre;
 90:     Mat jacobian_equality;
 91:     Mat jacobian_equality_pre;
 92:     Mat jacobian_state;
 93:     Mat jacobian_state_inv;
 94:     Mat jacobian_design;
 95:     Mat jacobian_state_pre;
 96:     Mat jacobian_design_pre;
 97:     IS state_is;
 98:     IS design_is;
 99:     PetscReal step;
100:     PetscReal residual;
101:     PetscReal gnorm0;
102:     PetscReal cnorm;
103:     PetscReal cnorm0;
104:     PetscReal fc;


107:     PetscInt  max_it;
108:     PetscInt  max_funcs;
109:     PetscInt  max_constraints;
110:     PetscInt  nfuncs;
111:     PetscInt  ngrads;
112:     PetscInt  nfuncgrads;
113:     PetscInt  nhess;
114:     PetscInt  niter;
115:     PetscInt  ntotalits;
116:     PetscInt  nconstraints;
117:     PetscInt  niconstraints;
118:     PetscInt  neconstraints;
119:     PetscInt  njac;
120:     PetscInt  njac_equality;
121:     PetscInt  njac_inequality;
122:     PetscInt  njac_state;
123:     PetscInt  njac_design;

125:     PetscInt  ksp_its; /* KSP iterations for this solver iteration */
126:     PetscInt  ksp_tot_its; /* Total (cumulative) KSP iterations */


129:     TaoLineSearch linesearch;
130:     PetscBool lsflag; /* goes up when line search fails */
131:     KSP ksp;
132:     PetscReal trust0; /* initial trust region radius */
133:     PetscReal trust;  /* Current trust region */

135:     PetscReal fatol;
136:     PetscReal frtol;
137:     PetscReal gatol;
138:     PetscReal grtol;
139:     PetscReal gttol;
140:     PetscReal catol;
141:     PetscReal crtol;
142:     PetscReal steptol;
143:     PetscReal fmin;
144:     PetscBool max_funcs_changed;
145:     PetscBool max_it_changed;
146:     PetscBool fatol_changed;
147:     PetscBool frtol_changed;
148:     PetscBool gatol_changed;
149:     PetscBool grtol_changed;
150:     PetscBool gttol_changed;
151:     PetscBool fmin_changed;
152:     PetscBool catol_changed;
153:     PetscBool crtol_changed;
154:     PetscBool steptol_changed;
155:     PetscBool trust0_changed;
156:     PetscBool printreason;
157:     PetscBool viewsolution;
158:     PetscBool viewgradient;
159:     PetscBool viewconstraints;
160:     PetscBool viewhessian;
161:     PetscBool viewjacobian;

163:     TaoSubsetType subset_type;
164:     PetscInt      hist_max;/* Number of iteration histories to keep */
165:     PetscReal     *hist_obj; /* obj value at each iteration */
166:     PetscReal     *hist_resid; /* residual at each iteration */
167:     PetscReal     *hist_cnorm; /* constraint norm at each iteration */
168:     PetscInt      *hist_lits; /* number of ksp its at each TAO iteration */
169:     PetscInt      hist_len;
170:     PetscBool     hist_reset;
171:     PetscBool     hist_malloc;
172: };

174: extern PetscLogEvent Tao_Solve, Tao_ObjectiveEval, Tao_ObjGradientEval, Tao_GradientEval, Tao_HessianEval, Tao_ConstraintsEval, Tao_JacobianEval;

178: PETSC_STATIC_INLINE PetscErrorCode TaoLogConvergenceHistory(Tao tao, PetscReal obj, PetscReal resid, PetscReal cnorm, PetscInt totits)
179: {
181:   if (tao->hist_max > tao->hist_len) {
182:     if (tao->hist_obj) tao->hist_obj[tao->hist_len]=obj;
183:     if (tao->hist_resid) tao->hist_resid[tao->hist_len]=resid;
184:     if (tao->hist_cnorm) tao->hist_cnorm[tao->hist_len]=cnorm;
185:     if (tao->hist_lits) {
186:       if (tao->hist_len <= 0) {
187:         tao->hist_lits[0] = totits;
188:       } else {
189:         tao->hist_lits[tao->hist_len]=totits - tao->hist_lits[tao->hist_len-1];
190:       }
191:     }
192:     tao->hist_len++;
193:   }
194:   return(0);
195: }

197: PETSC_INTERN PetscErrorCode TaoVecGetSubVec(Vec, IS, TaoSubsetType, PetscReal, Vec*);
198: PETSC_INTERN PetscErrorCode TaoMatGetSubMat(Mat, IS, Vec, TaoSubsetType, Mat*);

200: #endif