Actual source code: petsc-taoimpl.h

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #ifndef __TAO_IMPL_H

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

  8: typedef struct _TaoOps *TaoOps;

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

 27:     PetscErrorCode (*convergencetest)(Tao,void*);
 28:     PetscErrorCode (*convergencedestroy)(void*);

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

 39: #define MAXTAOMONITORS 10

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

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

 66:     PetscBool setupcalled;
 67:     void *data;

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


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

122:     PetscInt  ksp_its;


125:     TaoLineSearch linesearch;
126:     PetscBool lsflag; /* goes up when line search fails */
127:     KSP ksp;
128:     PetscReal trust0; /* initial trust region radius */
129:     PetscReal trust;  /* Current trust region */

131:     PetscReal fatol;
132:     PetscReal frtol;
133:     PetscReal gatol;
134:     PetscReal grtol;
135:     PetscReal gttol;
136:     PetscReal catol;
137:     PetscReal crtol;
138:     PetscReal xtol;
139:     PetscReal steptol;
140:     PetscReal fmin;

142:     PetscBool printreason;
143:     PetscBool viewsolution;
144:     PetscBool viewgradient;
145:     PetscBool viewconstraints;
146:     PetscBool viewhessian;
147:     PetscBool viewjacobian;

149:     TaoSubsetType subset_type;
150:     PetscInt      hist_max;/* Number of iteration histories to keep */
151:     PetscReal     *hist_obj; /* obj value at each iteration */
152:     PetscReal     *hist_resid; /* residual at each iteration */
153:     PetscReal     *hist_cnorm; /* constraint norm at each iteration */
154:     PetscInt      hist_len;
155:     PetscBool     hist_reset;
156: };

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

160: #define TaoLogHistory(tao,obj,resid,cnorm) \
161:   { if (tao->hist_max > tao->hist_len) \
162:       { if (tao->hist_obj) tao->hist_obj[tao->hist_len]=obj;\
163:         if (tao->hist_resid) tao->hist_resid[tao->hist_len]=resid;\
164:         if (tao->hist_cnorm) tao->hist_cnorm[tao->hist_len]=cnorm;} \
165:     tao->hist_len++;\
166:   }

168: #endif

170: PETSC_INTERN PetscErrorCode TaoVecGetSubVec(Vec, IS, TaoSubsetType, PetscReal, Vec*);
171: PETSC_INTERN PetscErrorCode TaoMatGetSubMat(Mat, IS, Vec, TaoSubsetType, Mat*);