00001 #ifndef OptConstrNewtonLike_h
00002 #define OptConstrNewtonLike_h
00003
00004
00005
00006
00007
00008
00009 #include "Opt.h"
00010
00011 using std::ostream;
00012
00013 namespace OPTPP {
00014
00031 class OptConstrNewtonLike: public OptimizeClass {
00032 protected:
00033 virtual NLP1* nlprob() const = 0;
00034 int me;
00035 int mi;
00036 int grad_evals;
00037 NEWMAT::ColumnVector gprev;
00038 NEWMAT::ColumnVector z;
00039 NEWMAT::ColumnVector y;
00040 NEWMAT::ColumnVector s;
00041 NEWMAT::ColumnVector constrType;
00042 NEWMAT::ColumnVector constraintResidual;
00043 NEWMAT::ColumnVector gradl;
00044 NEWMAT::ColumnVector gradlprev;
00045 NEWMAT::Matrix constraintGradient;
00046 NEWMAT::Matrix constraintGradientPrev;
00047 NEWMAT::SymmetricMatrix Hessian;
00048 NEWMAT::SymmetricMatrix hessl;
00049 SearchStrategy strategy;
00050 DerivOption finitediff;
00051 MeritFcn mfcn;
00052 real TR_size;
00053 real gradMult;
00054 int searchSize;
00055 real cost;
00056 void defaultAcceptStep(int, int);
00057 NEWMAT::ColumnVector defaultComputeSearch(NEWMAT::SymmetricMatrix& );
00058 bool WarmStart;
00059 bool feas_flag;
00060 int max_feas_iter;
00061
00062 public:
00070 OptConstrNewtonLike() {}
00076 OptConstrNewtonLike(int n):
00077 OptimizeClass(n), me(0), mi(0), grad_evals(0), gprev(n),
00078 z(n), y(n), s(n), constrType(n), constraintResidual(n),
00079 gradl(n), gradlprev(n),constraintGradient(n,n), constraintGradientPrev(n,n),
00080 Hessian(n), hessl(n), strategy(TrustRegion), finitediff(ForwardDiff),
00081 mfcn(ArgaezTapia), TR_size(0.0),
00082 gradMult(0.1), searchSize(64), cost(0.0), WarmStart(false),
00083 feas_flag(false), max_feas_iter(3)
00084 {z = 0; y = 0; s = 0;}
00091 OptConstrNewtonLike(int n, UPDATEFCN u):
00092 OptimizeClass(n), me(0), mi(0), grad_evals(0), gprev(n),
00093 z(n), y(n), s(n), constrType(n), constraintResidual(n),
00094 gradl(n), gradlprev(n),constraintGradient(n,n), constraintGradientPrev(n,n),
00095 Hessian(n), hessl(n), strategy(TrustRegion), finitediff(ForwardDiff),
00096 mfcn(ArgaezTapia), TR_size(0.0),
00097 gradMult(0.1), searchSize(64), cost(0.0), WarmStart(false),
00098 feas_flag(false), max_feas_iter(3)
00099 {update_fcn = u; z = 0; y = 0; s = 0;}
00106 OptConstrNewtonLike(int n, TOLS t):
00107 OptimizeClass(n,t), me(0), mi(0), grad_evals(0), gprev(n),
00108 z(n), y(n), s(n), constrType(n), constraintResidual(n),
00109 gradl(n), gradlprev(n),constraintGradient(n,n), constraintGradientPrev(n,n),
00110 Hessian(n), hessl(n), strategy(TrustRegion), finitediff(ForwardDiff),
00111 mfcn(ArgaezTapia), TR_size(0.0),
00112 gradMult(0.1), searchSize(64), cost(0.0), WarmStart(false),
00113 feas_flag(false), max_feas_iter(3)
00114 {z = 0; y = 0; s = 0;}
00115
00119 virtual ~OptConstrNewtonLike(){}
00120
00121
00122
00123
00124
00125
00127 virtual void acceptStep(int k, int step_type)
00128 {defaultAcceptStep(k, step_type);}
00129
00131 virtual NEWMAT::ColumnVector computeSearch(NEWMAT::SymmetricMatrix& H )
00132 {return defaultComputeSearch(H);}
00133
00134 virtual void updateModel(int k, int ndim, NEWMAT::ColumnVector x)
00135 {OptimizeClass::defaultUpdateModel(k, ndim, x);}
00136
00137 virtual void reset();
00138
00139
00140
00142 virtual int checkConvg();
00143
00145 virtual int checkDeriv();
00146
00149 virtual int computeStep(NEWMAT::ColumnVector sk);
00150
00151 virtual real computeMaxStep(NEWMAT::ColumnVector sk);
00152
00154 virtual void initOpt();
00155
00157 virtual void initHessian();
00158
00160 virtual double initTrustRegionSize() const;
00161
00163 virtual void optimize();
00164
00166 virtual void readOptInput();
00167
00168
00169
00171 int checkAnalyticFDGrad();
00172
00176 int getFevals() const {return fcn_evals;}
00177
00181 int getGevals() const {return grad_evals;}
00182
00186 real getTRSize() const {return TR_size;}
00188 void setTRSize(real delta) {TR_size = delta;}
00189
00193 real getGradMult() const {return gradMult;}
00195 void setGradMult(real tau) {gradMult = tau;}
00196
00201 int getSearchSize() const {return searchSize;}
00203 void setSearchSize(int sss) {searchSize = sss;}
00204
00205 bool getWarmStart() const {return WarmStart;}
00206 void UseWarmStart(NEWMAT::SymmetricMatrix& H) {Hessian = H; WarmStart = true;}
00207
00211 void setSearchStrategy(SearchStrategy search) {strategy = search;}
00213 SearchStrategy getSearchStrategy() const {return strategy;}
00214
00218 void setDerivOption(DerivOption d) {finitediff = d;}
00220 DerivOption getDerivOption() const {return finitediff;}
00221
00225 NEWMAT::ColumnVector getEqualityMultiplier() const { return y;}
00227 void setEqualityMultiplier(const NEWMAT::ColumnVector& ymult) { y = ymult;}
00228
00232 NEWMAT::ColumnVector getInequalityMultiplier() const { return z;}
00234 void setInequalityMultiplier(const NEWMAT::ColumnVector& zmult){ z = zmult;}
00235
00239 NEWMAT::ColumnVector getSlacks() const { return s;}
00241 void setSlacks(const NEWMAT::ColumnVector& slackVar) { s = slackVar;}
00242
00246 MeritFcn getMeritFcn() const { return mfcn;}
00248 virtual void setMeritFcn(MeritFcn option) { mfcn = option;}
00249
00253 NEWMAT::ColumnVector getGradL() const { return gradl;}
00255 virtual void setGradL(NEWMAT::ColumnVector gradl_value) { gradl = gradl_value;}
00256
00260 NEWMAT::ColumnVector getGradLPrev() const { return gradlprev;}
00262 virtual void setGradLPrev(NEWMAT::ColumnVector gradl_value) { gradlprev = gradl_value;}
00263
00267 NEWMAT::ColumnVector getConstraintResidual() const { return constraintResidual;}
00269 virtual void setConstraintResidual(const NEWMAT::ColumnVector& constraint_value)
00270 { constraintResidual = constraint_value;}
00271
00275 NEWMAT::Matrix getConstraintGradient() const { return constraintGradient;}
00277 virtual void setConstraintGradient(const NEWMAT::Matrix& constraint_grad)
00278 { constraintGradient = constraint_grad;}
00279
00283 real getCost() const { return cost;}
00285 void setCost(real value) { cost = value;}
00286
00290 bool getFeasibilityRecoveryFlag() const {return feas_flag;}
00292 void setFeasibilityRecoveryFlag(bool flag) {feas_flag = flag;}
00293
00297 int getMaxFeasIter() const {return max_feas_iter;}
00299 void setMaxFeasIter(int k) {max_feas_iter = k;}
00300
00304 NEWMAT::SymmetricMatrix getHessian() const {return Hessian;}
00306 void setHessian(NEWMAT::SymmetricMatrix& H) {Hessian = H;}
00307
00309 virtual NEWMAT::SymmetricMatrix updateH(NEWMAT::SymmetricMatrix& H, int k) = 0;
00310
00312 NEWMAT::ColumnVector computeFFK1Ind(const NEWMAT::ColumnVector& xc);
00314 NEWMAT::ColumnVector computeFFK2Ind(const NEWMAT::ColumnVector& xc);
00316 NEWMAT::ColumnVector computeTapiaInd(const NEWMAT::ColumnVector& step);
00317
00319 void printStatus(char *);
00321 void printMultipliers(char *);
00323 void fPrintMultipliers(ostream *nlpout, char *);
00325 void fPrintSecSuff(ostream *nlpout, NEWMAT::ColumnVector& info);
00326
00327 friend int trustregion(NLP1*, ostream*, NEWMAT::SymmetricMatrix&,
00328 NEWMAT::ColumnVector&, NEWMAT::ColumnVector&,
00329 real&, real&, real stpmax, real stpmin);
00330
00331 friend int trustpds(NLP1*, ostream*, NEWMAT::SymmetricMatrix&,
00332 NEWMAT::ColumnVector&, NEWMAT::ColumnVector&,
00333 real&, real&, real stpmax, real stpmin, int);
00334 };
00335
00340 class OptConstrNewton1Deriv: public OptConstrNewtonLike {
00341 public:
00342
00343 OptConstrNewton1Deriv() {}
00344
00345 OptConstrNewton1Deriv(NLP1* p):
00346 OptConstrNewtonLike(p->getDim()), mem_nlp(p) {;}
00347
00348 OptConstrNewton1Deriv(NLP1* p, UPDATEFCN u):
00349 OptConstrNewtonLike(p->getDim(),u), mem_nlp(p) {;}
00350
00351 OptConstrNewton1Deriv(NLP1* p, TOLS t):
00352 OptConstrNewtonLike(p->getDim(),t), mem_nlp(p) {;}
00353
00354 virtual ~OptConstrNewton1Deriv(){}
00355
00356 private:
00357 NLP1* mem_nlp;
00358
00359 protected:
00363 NLP1* nlprob() const { return mem_nlp; }
00364 };
00365
00370 class OptConstrNewton2Deriv: public OptConstrNewtonLike {
00371 public:
00372
00373 OptConstrNewton2Deriv() {}
00374
00375 OptConstrNewton2Deriv(NLP2* p):
00376 OptConstrNewtonLike(p->getDim()), mem_nlp(p){;}
00377
00378 OptConstrNewton2Deriv(NLP2* p, UPDATEFCN u):
00379 OptConstrNewtonLike(p->getDim(),u), mem_nlp(p){;}
00380
00381 OptConstrNewton2Deriv(NLP2* p, TOLS t):
00382 OptConstrNewtonLike(p->getDim(),t), mem_nlp(p){;}
00383
00384 virtual ~OptConstrNewton2Deriv(){}
00385
00386 private:
00387 NLP2* mem_nlp;
00388
00389 protected:
00393 NLP1* nlprob() const { return mem_nlp;}
00397 NLP2* nlprob2() const { return mem_nlp;}
00398 };
00399
00400 }
00401
00402 #endif