Actual source code: taoimpl.h

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: #ifndef __TAO_IMPL_H

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

  8: PETSC_EXTERN PetscBool TaoRegisterAllCalled;
  9: PETSC_EXTERN PetscErrorCode TaoRegisterAll(void);

 11: typedef struct _TaoOps *TaoOps;

 13: struct _TaoOps {
 14:   /* Methods set by application */
 15:     PetscErrorCode (*computeobjective)(Tao, Vec, PetscReal*, void*);
 16:     PetscErrorCode (*computeobjectiveandgradient)(Tao, Vec, PetscReal*, Vec, void*);
 17:     PetscErrorCode (*computegradient)(Tao, Vec, Vec, void*);
 18:     PetscErrorCode (*computehessian)(Tao, Vec, Mat, Mat,  void*);
 19:     PetscErrorCode (*computeresidual)(Tao, Vec, Vec, void*);
 20:     PetscErrorCode (*computeresidualjacobian)(Tao, Vec, Mat, Mat, void*);
 21:     PetscErrorCode (*computeconstraints)(Tao, Vec, Vec, void*);
 22:     PetscErrorCode (*computeinequalityconstraints)(Tao, Vec, Vec, void*);
 23:     PetscErrorCode (*computeequalityconstraints)(Tao, Vec, Vec, void*);
 24:     PetscErrorCode (*computejacobian)(Tao, Vec, Mat, Mat,  void*);
 25:     PetscErrorCode (*computejacobianstate)(Tao, Vec, Mat, Mat, Mat,  void*);
 26:     PetscErrorCode (*computejacobiandesign)(Tao, Vec, Mat, void*);
 27:     PetscErrorCode (*computejacobianinequality)(Tao, Vec, Mat, Mat,  void*);
 28:     PetscErrorCode (*computejacobianequality)(Tao, Vec, Mat, Mat,  void*);
 29:     PetscErrorCode (*computebounds)(Tao, Vec, Vec, void*);
 30:     PetscErrorCode (*update)(Tao,PetscInt,void*);
 31:     PetscErrorCode (*convergencetest)(Tao,void*);
 32:     PetscErrorCode (*convergencedestroy)(void*);

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

 43: #define MAXTAOMONITORS 10

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

 65:     PetscErrorCode (*monitor[MAXTAOMONITORS])(Tao,void*);
 66:     PetscErrorCode (*monitordestroy[MAXTAOMONITORS])(void**);
 67:     void *monitorcontext[MAXTAOMONITORS];
 68:     PetscInt numbermonitors;
 69:     void *cnvP;
 70:     TaoConvergedReason reason;

 72:     PetscBool setupcalled;
 73:     void *data;

 75:     Vec solution;
 76:     Vec gradient;
 77:     Vec stepdirection;
 78:     Vec XL;
 79:     Vec XU;
 80:     Vec IL;
 81:     Vec IU;
 82:     Vec DI;
 83:     Vec DE;
 84:     Mat hessian;
 85:     Mat hessian_pre;
 86:     Mat gradient_norm;
 87:     Vec gradient_norm_tmp;
 88:     Vec ls_res;
 89:     Mat ls_jac;
 90:     Mat ls_jac_pre;
 91:     Vec res_weights_v;
 92:     PetscInt res_weights_n;
 93:     PetscInt *res_weights_rows;
 94:     PetscInt *res_weights_cols;
 95:     PetscReal *res_weights_w;
 96:     Vec constraints;
 97:     Vec constraints_equality;
 98:     Vec constraints_inequality;
 99:     Mat jacobian;
100:     Mat jacobian_pre;
101:     Mat jacobian_inequality;
102:     Mat jacobian_inequality_pre;
103:     Mat jacobian_equality;
104:     Mat jacobian_equality_pre;
105:     Mat jacobian_state;
106:     Mat jacobian_state_inv;
107:     Mat jacobian_design;
108:     Mat jacobian_state_pre;
109:     Mat jacobian_design_pre;
110:     IS state_is;
111:     IS design_is;
112:     PetscReal step;
113:     PetscReal residual;
114:     PetscReal gnorm0;
115:     PetscReal cnorm;
116:     PetscReal cnorm0;
117:     PetscReal fc;


120:     PetscInt  max_it;
121:     PetscInt  max_funcs;
122:     PetscInt  max_constraints;
123:     PetscInt  nfuncs;
124:     PetscInt  ngrads;
125:     PetscInt  nfuncgrads;
126:     PetscInt  nhess;
127:     PetscInt  niter;
128:     PetscInt  ntotalits;
129:     PetscInt  nconstraints;
130:     PetscInt  niconstraints;
131:     PetscInt  neconstraints;
132:     PetscInt  njac;
133:     PetscInt  njac_equality;
134:     PetscInt  njac_inequality;
135:     PetscInt  njac_state;
136:     PetscInt  njac_design;

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


142:     TaoLineSearch linesearch;
143:     PetscBool lsflag; /* goes up when line search fails */
144:     KSP ksp;
145:     PetscReal trust0; /* initial trust region radius */
146:     PetscReal trust;  /* Current trust region */

148:     PetscReal gatol;
149:     PetscReal grtol;
150:     PetscReal gttol;
151:     PetscReal catol;
152:     PetscReal crtol;
153:     PetscReal steptol;
154:     PetscReal fmin;
155:     PetscBool max_funcs_changed;
156:     PetscBool max_it_changed;
157:     PetscBool gatol_changed;
158:     PetscBool grtol_changed;
159:     PetscBool gttol_changed;
160:     PetscBool fmin_changed;
161:     PetscBool catol_changed;
162:     PetscBool crtol_changed;
163:     PetscBool steptol_changed;
164:     PetscBool trust0_changed;
165:     PetscBool printreason;
166:     PetscBool viewsolution;
167:     PetscBool viewgradient;
168:     PetscBool viewconstraints;
169:     PetscBool viewhessian;
170:     PetscBool viewjacobian;
171:     PetscBool bounded;
172:     PetscBool header_printed;

174:     TaoSubsetType subset_type;
175:     PetscInt      hist_max;/* Number of iteration histories to keep */
176:     PetscReal     *hist_obj; /* obj value at each iteration */
177:     PetscReal     *hist_resid; /* residual at each iteration */
178:     PetscReal     *hist_cnorm; /* constraint norm at each iteration */
179:     PetscInt      *hist_lits; /* number of ksp its at each TAO iteration */
180:     PetscInt      hist_len;
181:     PetscBool     hist_reset;
182:     PetscBool     hist_malloc;
183: };

185: PETSC_EXTERN PetscLogEvent TAO_Solve;
186: PETSC_EXTERN PetscLogEvent TAO_ObjectiveEval;
187: PETSC_EXTERN PetscLogEvent TAO_GradientEval;
188: PETSC_EXTERN PetscLogEvent TAO_ObjGradEval;
189: PETSC_EXTERN PetscLogEvent TAO_HessianEval;
190: PETSC_EXTERN PetscLogEvent TAO_ConstraintsEval;
191: PETSC_EXTERN PetscLogEvent TAO_JacobianEval;

193: PETSC_STATIC_INLINE PetscErrorCode TaoLogConvergenceHistory(Tao tao, PetscReal obj, PetscReal resid, PetscReal cnorm, PetscInt totits)
194: {
196:   if (tao->hist_max > tao->hist_len) {
197:     if (tao->hist_obj) tao->hist_obj[tao->hist_len]=obj;
198:     if (tao->hist_resid) tao->hist_resid[tao->hist_len]=resid;
199:     if (tao->hist_cnorm) tao->hist_cnorm[tao->hist_len]=cnorm;
200:     if (tao->hist_lits) {
201:       if (tao->hist_len <= 0) {
202:         tao->hist_lits[0] = totits;
203:       } else {
204:         tao->hist_lits[tao->hist_len]=totits - tao->hist_lits[tao->hist_len-1];
205:       }
206:     }
207:     tao->hist_len++;
208:   }
209:   return(0);
210: }

212: PETSC_INTERN PetscErrorCode TaoVecGetSubVec(Vec, IS, TaoSubsetType, PetscReal, Vec*);
213: PETSC_INTERN PetscErrorCode TaoMatGetSubMat(Mat, IS, Vec, TaoSubsetType, Mat*);
214: PETSC_INTERN PetscErrorCode TaoGradientNorm(Tao, Vec, NormType, PetscReal*);
215: PETSC_INTERN PetscErrorCode TaoEstimateActiveBounds(Vec, Vec, Vec, Vec, Vec, Vec, PetscReal, PetscReal*, IS*, IS*, IS*, IS*, IS*);
216: PETSC_INTERN PetscErrorCode TaoBoundStep(Vec, Vec, Vec, IS, IS, IS, PetscReal, Vec);
217: PETSC_INTERN PetscErrorCode TaoBoundSolution(Vec, Vec, Vec, PetscReal, PetscInt*, Vec);

219: #endif