10#ifndef ROL_TYPEE_COMPOSITESTEPALGORITHM_DEF_H
11#define ROL_TYPEE_COMPOSITESTEPALGORITHM_DEF_H
16template<
typename Real>
21 status_->add(makePtr<ConstraintStatusTest<Real>>(list_));
27 Real one(1), two(2), p8(0.8),
zero(0), oem8(1.e-8);
28 ParameterList& cslist = list_.sublist(
"Step").sublist(
"Composite Step");
31 tolOSS_ = cslist.sublist(
"Optimality System Solver").get(
"Nominal Relative Tolerance", 1e-8);
32 tolOSSfixed_ = cslist.sublist(
"Optimality System Solver").get(
"Fix Tolerance",
true);
33 maxiterCG_ = cslist.sublist(
"Tangential Subproblem Solver").get(
"Iteration Limit", 20);
34 tolCG_ = cslist.sublist(
"Tangential Subproblem Solver").get(
"Relative Tolerance", 1e-2);
35 Delta_ = cslist.get(
"Initial Radius", 1e2);
36 useConHess_ = cslist.get(
"Use Constraint Hessian",
true);
37 verbosity_ = list_.sublist(
"General").get(
"Output Level", 0);
38 printHeader_ = (verbosity_ > 2);
64 infoQN_ = infoQN_ || infoALL_;
65 infoLM_ = infoLM_ || infoALL_;
66 infoTS_ = infoTS_ || infoALL_;
67 infoAC_ = infoAC_ || infoALL_;
68 infoLS_ = infoLS_ || infoALL_;
81template<
typename Real>
82void CompositeStepAlgorithm<Real>::computeTrial(Vector<Real> &s,
83 const Vector<Real> &x,
84 const Vector<Real> &l,
86 Constraint<Real> &con,
89 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
91 ROL::Ptr<Vector<Real> > n = xvec_->clone();
92 ROL::Ptr<Vector<Real> > c = cvec_->clone();
93 ROL::Ptr<Vector<Real> > t = xvec_->clone();
94 ROL::Ptr<Vector<Real> > tCP = xvec_->clone();
95 ROL::Ptr<Vector<Real> > g = gvec_->clone();
96 ROL::Ptr<Vector<Real> > gf = gvec_->clone();
97 ROL::Ptr<Vector<Real> > Wg = xvec_->clone();
98 ROL::Ptr<Vector<Real> > ajl = gvec_->clone();
101 ROL::Ptr<Vector<Real> > l_new = lvec_->clone();
102 ROL::Ptr<Vector<Real> > c_new = cvec_->clone();
103 ROL::Ptr<Vector<Real> > g_new = gvec_->clone();
104 ROL::Ptr<Vector<Real> > gf_new = gvec_->clone();
107 f = obj.value(x, zerotol);
110 obj.gradient(*gf, x, zerotol);
112 con.value(*c, x, zerotol);
115 computeQuasinormalStep(*n, *c, x, zeta_*Delta_, con, os);
118 con.applyAdjointJacobian(*ajl, l, x, zerotol);
124 solveTangentialSubproblem(*t, *tCP, *Wg, x, *g, *n, l, Delta_, obj, con, os);
125 totalIterCG_ += iterCG_;
128 accept(s, *n, *t, f_new, *c_new, *gf_new, *l_new, *g_new, x, l, f, *gf, *c, *g, *tCP, *Wg, obj, con, os);
132template<
typename Real>
133void CompositeStepAlgorithm<Real>::computeLagrangeMultiplier(Vector<Real> &l,
134 const Vector<Real> &x,
135 const Vector<Real> &gf,
136 Constraint<Real> &con,
140 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
141 std::vector<Real> augiters;
145 os <<
"\n Lagrange multiplier step\n";
150 Ptr<Vector<Real> > ajl = gvec_->clone();
151 con.applyAdjointJacobian(*ajl, l, x, zerotol);
154 Ptr<Vector<Real> > b1 = gvec_->clone();
155 Ptr<Vector<Real> > b2 = cvec_->clone();
157 b1->set(gf); b1->plus(*ajl); b1->scale(-one);
162 Ptr<Vector<Real> > v1 = xvec_->clone();
163 Ptr<Vector<Real> > v2 = lvec_->clone();
166 Real b1norm = b1->norm();
167 Real tol = setTolOSS(lmhtol_*b1norm);
170 augiters = con.solveAugmentedSystem(*v1, *v2, *b1, *b2, x, tol);
172 totalIterLS_ = totalIterLS_ + augiters.size();
173 printInfoLS(augiters, os);
182template<
typename Real>
183void CompositeStepAlgorithm<Real>::computeQuasinormalStep(Vector<Real> &n,
184 const Vector<Real> &c,
185 const Vector<Real> &x,
187 Constraint<Real> &con,
192 os <<
"\n Quasi-normal step\n";
198 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
199 std::vector<Real> augiters;
202 Ptr<Vector<Real> > nCP = xvec_->clone();
203 Ptr<Vector<Real> > nCPdual = gvec_->clone();
204 Ptr<Vector<Real> > nN = xvec_->clone();
205 Ptr<Vector<Real> > ctemp = cvec_->clone();
206 Ptr<Vector<Real> > dualc0 = lvec_->clone();
207 dualc0->set(c.dual());
208 con.applyAdjointJacobian(*nCPdual, *dualc0, x, zerotol);
209 nCP->set(nCPdual->dual());
210 con.applyJacobian(*ctemp, *nCP, x, zerotol);
212 Real normsquare_ctemp = ctemp->dot(*ctemp);
213 if (normsquare_ctemp !=
zero) {
214 nCP->scale( -(nCP->dot(*nCP))/normsquare_ctemp );
219 Real norm_nCP = nCP->norm();
220 if (norm_nCP >= delta) {
222 n.scale( delta/norm_nCP );
225 os <<
" taking partial Cauchy step\n";
234 con.applyJacobian(*ctemp, *nCP, x, zerotol);
236 Real tol = setTolOSS(qntol_*ctemp->norm());
239 nCPdual->set(nCP->dual());
240 nCPdual->scale(-one);
242 Ptr<Vector<Real> > dn = xvec_->clone();
243 Ptr<Vector<Real> > y = lvec_->clone();
245 augiters = con.solveAugmentedSystem(*dn, *y, *nCPdual, *ctemp, x, tol);
247 totalIterLS_ = totalIterLS_ + augiters.size();
248 printInfoLS(augiters, os);
255 Real norm_nN = nN->norm();
256 if (norm_nN <= delta) {
261 os <<
" taking full Newton step\n";
269 Real aa = dn->dot(*dn);
270 Real bb = dn->dot(*nCP);
271 Real cc = norm_nCP*norm_nCP - delta*delta;
272 Real tau = (-bb+sqrt(bb*bb-aa*cc))/aa;
277 os <<
" taking dogleg step\n";
285template<
typename Real>
286void CompositeStepAlgorithm<Real>::solveTangentialSubproblem(Vector<Real> &t,
289 const Vector<Real> &x,
290 const Vector<Real> &g,
291 const Vector<Real> &n,
292 const Vector<Real> &l,
294 Objective<Real> &obj,
295 Constraint<Real> &con,
299 bool orthocheck =
true;
301 Real tol_ortho = 0.5;
308 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
309 std::vector<Real> augiters;
314 Ptr<Vector<Real> > r = gvec_->clone();
315 Ptr<Vector<Real> > pdesc = xvec_->clone();
316 Ptr<Vector<Real> > tprev = xvec_->clone();
317 Ptr<Vector<Real> > Wr = xvec_->clone();
318 Ptr<Vector<Real> > Hp = gvec_->clone();
319 Ptr<Vector<Real> > xtemp = xvec_->clone();
320 Ptr<Vector<Real> > gtemp = gvec_->clone();
321 Ptr<Vector<Real> > ltemp = lvec_->clone();
322 Ptr<Vector<Real> > czero = cvec_->clone();
325 obj.hessVec(*gtemp, n, x, zerotol);
328 con.applyAdjointHessian(*gtemp, l, n, x, zerotol);
331 Real normg = r->norm();
339 std::vector<Real> normWr(maxiterCG_+1,
zero);
341 std::vector<Ptr<Vector<Real > > > p;
342 std::vector<Ptr<Vector<Real > > > Hps;
343 std::vector<Ptr<Vector<Real > > > rs;
344 std::vector<Ptr<Vector<Real > > > Wrs;
349 std::ios_base::fmtflags osFlags(os.flags());
350 os <<
"\n Tangential subproblem\n";
351 os << std::setw(6) << std::right <<
"iter" << std::setw(18) <<
"||Wr||/||Wr0||" << std::setw(15) <<
"||s||";
352 os << std::setw(15) <<
"delta" << std::setw(15) <<
"||c'(x)s||" <<
"\n";
359 os <<
" >>> Tangential subproblem: Initial gradient is zero! \n";
362 iterCG_ = 0; Wg.zero(); flagCG_ = 0;
367 while (iterCG_ < maxiterCG_) {
377 Real tol = setTolOSS(pgtol_);
378 augiters = con.solveAugmentedSystem(*Wr, *ltemp, *r, *czero, x, tol);
380 totalIterLS_ = totalIterLS_ + augiters.size();
381 printInfoLS(augiters, os);
386 Wrs.push_back(xvec_->clone());
387 (Wrs[iterCG_-1])->set(*Wr);
390 if (normWg ==
zero) {
395 os <<
" Initial projected residual is close to zero! \n";
404 rs.push_back(xvec_->clone());
406 (rs[0])->set(r->dual());
411 Real tol = setTolOSS(projtol_);
412 augiters = con.solveAugmentedSystem(*Wr, *ltemp, *r, *czero, x, tol);
414 totalIterLS_ = totalIterLS_ + augiters.size();
415 printInfoLS(augiters, os);
418 Wrs.push_back(xvec_->clone());
419 (Wrs[iterCG_-1])->set(*Wr);
423 normWr[iterCG_-1] = Wr->norm();
426 Ptr<Vector<Real> > ct = cvec_->clone();
427 con.applyJacobian(*ct, t, x, zerotol);
428 Real linc = ct->norm();
429 std::ios_base::fmtflags osFlags(os.flags());
430 os << std::scientific << std::setprecision(6);
431 os << std::setw(6) << std::right << iterCG_-1 << std::setw(18) << normWr[iterCG_-1]/normWg << std::setw(15) << t.norm();
432 os << std::setw(15) << delta << std::setw(15) << linc <<
"\n";
437 if (normWr[iterCG_-1]/normWg < tolCG_) {
442 os <<
" || W(g + H*(n+s)) || <= cgtol*|| W(g + H*n)|| \n";
450 LA::Matrix<Real> Wrr(iterCG_,iterCG_);
451 LA::Matrix<Real> T(iterCG_,iterCG_);
452 LA::Matrix<Real> Tm1(iterCG_,iterCG_);
453 for (
int i=0; i<iterCG_; i++) {
454 for (
int j=0; j<iterCG_; j++) {
455 Wrr(i,j) = (Wrs[i])->dot(*rs[j]);
456 T(i,j) = Wrr(i,j)/(normWr[i]*normWr[j]);
459 Tm1(i,j) = Tm1(i,j) - one;
463 if (Tm1.normOne() >= tol_ortho) {
464 LAPACK<int,Real> lapack;
465 std::vector<int> ipiv(iterCG_);
467 std::vector<Real> work(3*iterCG_);
469 lapack.GETRF(iterCG_, iterCG_, T.values(), T.stride(), &ipiv[0], &info);
470 lapack.GETRI(iterCG_, T.values(), T.stride(), &ipiv[0], &work[0], 3*iterCG_, &info);
472 for (
int i=0; i<iterCG_; i++) {
473 Tm1(i,i) = Tm1(i,i) - one;
475 if (Tm1.normOne() > S_max) {
479 os <<
" large nonorthogonality in W(R)'*R detected \n";
488 p.push_back(xvec_->clone());
489 (p[iterCG_-1])->set(*Wr);
490 (p[iterCG_-1])->scale(-one);
491 for (
int j=1; j<iterCG_; j++) {
492 Real scal = (p[iterCG_-1])->dot(*(Hps[j-1])) / (p[j-1])->dot(*(Hps[j-1]));
493 Ptr<Vector<Real> > pj = xvec_->clone();
496 (p[iterCG_-1])->plus(*pj);
500 Hps.push_back(xvec_->clone());
502 obj.hessVec(*Hp, *(p[iterCG_-1]), x, zerotol);
504 con.applyAdjointHessian(*gtemp, l, *(p[iterCG_-1]), x, zerotol);
509 (Hps[iterCG_-1])->set(Hp->dual());
511 pHp = (p[iterCG_-1])->dot(*(Hps[iterCG_-1]));
513 rp = (p[iterCG_-1])->dot(*(rs[iterCG_-1]));
515 normp = (p[iterCG_-1])->norm();
520 pdesc->set(*(p[iterCG_-1]));
521 if ((std::abs(rp) >= rptol*normp*normr) && (sgn(rp) == 1)) {
525 Real a = pdesc->dot(*pdesc);
526 Real b = pdesc->dot(t);
527 Real c = t.dot(t) - delta*delta;
529 Real theta = (-b + std::sqrt(b*b - a*c)) / a;
530 xtemp->set(*(p[iterCG_-1]));
539 os <<
" negative curvature detected \n";
546 if (std::abs(rp) < rptol*normp*normr) {
550 os <<
" Zero alpha due to inexactness. \n";
560 xtemp->set(*(p[iterCG_-1]));
566 if (normt >= delta) {
567 pdesc->set(*(p[iterCG_-1]));
571 Real a = pdesc->dot(*pdesc);
572 Real b = pdesc->dot(*tprev);
573 Real c = tprev->dot(*tprev) - delta*delta;
575 Real theta = (-b + std::sqrt(b*b - a*c)) / a;
576 xtemp->set(*(p[iterCG_-1]));
587 os <<
" trust-region condition active \n";
594 xtemp->set(*(Hps[iterCG_-1]));
597 r->plus(xtemp->dual());
600 rs.push_back(xvec_->clone());
602 (rs[iterCG_])->set(r->dual());
612 os <<
" maximum number of iterations reached \n";
619template<
typename Real>
620void CompositeStepAlgorithm<Real>::accept(Vector<Real> &s, Vector<Real> &n, Vector<Real> &t, Real f_new, Vector<Real> &c_new,
621 Vector<Real> &gf_new, Vector<Real> &l_new, Vector<Real> &g_new,
622 const Vector<Real> &x,
const Vector<Real> &l, Real f,
const Vector<Real> &gf,
const Vector<Real> &c,
623 const Vector<Real> &g, Vector<Real> &tCP, Vector<Real> &Wg,
624 Objective<Real> &obj, Constraint<Real> &con, std::ostream &os) {
627 Real tol_red_tang = 1e-3;
628 Real tol_red_all = 1e-1;
631 Real tol_fdiff = 1e-12;
636 Real rpred_over_pred = 0.5*(1-eta_);
640 os <<
"\n Composite step acceptance\n";
648 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
649 std::vector<Real> augiters;
654 Real part_pred =
zero;
655 Real linc_preproj =
zero;
656 Real linc_postproj =
zero;
657 Real tangtol_start =
zero;
658 Real tangtol = tangtol_;
662 bool try_tCP =
false;
665 Ptr<Vector<Real> > xtrial = xvec_->clone();
666 Ptr<Vector<Real> > Jl = gvec_->clone();
667 Ptr<Vector<Real> > gfJl = gvec_->clone();
668 Ptr<Vector<Real> > Jnc = cvec_->clone();
669 Ptr<Vector<Real> > t_orig = xvec_->clone();
670 Ptr<Vector<Real> > t_dual = gvec_->clone();
671 Ptr<Vector<Real> > Jt_orig = cvec_->clone();
672 Ptr<Vector<Real> > t_m_tCP = xvec_->clone();
673 Ptr<Vector<Real> > ltemp = lvec_->clone();
674 Ptr<Vector<Real> > xtemp = xvec_->clone();
675 Ptr<Vector<Real> > rt = cvec_->clone();
676 Ptr<Vector<Real> > Hn = gvec_->clone();
677 Ptr<Vector<Real> > Hto = gvec_->clone();
678 Ptr<Vector<Real> > cxxvec = gvec_->clone();
679 Ptr<Vector<Real> > czero = cvec_->clone();
681 Real Jnc_normsquared =
zero;
682 Real c_normsquared =
zero;
686 con.applyAdjointJacobian(*Jl, l, x, zerotol);
687 con.applyJacobian(*Jnc, n, x, zerotol);
689 Jnc_normsquared = Jnc->dot(*Jnc);
690 c_normsquared = c.dot(c);
692 for (
int ct=0; ct<ct_max; ct++) {
696 t_m_tCP->scale(-one);
698 if (t_m_tCP->norm() ==
zero) {
703 con.applyJacobian(*Jt_orig, *t_orig, x, zerotol);
704 linc_preproj = Jt_orig->norm();
706 rpred = two*rpred_over_pred*pred;
709 tangtol_start = tangtol;
711 while (std::abs(rpred)/pred > rpred_over_pred) {
714 tangtol = tol_red_tang*tangtol;
716 if (tangtol < mintol) {
719 os <<
"\n The projection of the tangential step cannot be done with sufficient precision.\n";
720 os <<
" Is the quasi-normal step very small? Continuing with no global convergence guarantees.\n";
727 Real tol = setTolOSS(tangtol);
729 t_dual->set(t_orig->dual());
730 augiters = con.solveAugmentedSystem(t, *ltemp, *t_dual, *czero, x, tol);
732 totalIterLS_ = totalIterLS_ + augiters.size();
733 printInfoLS(augiters, os);
735 con.applyJacobian(*rt, t, x, zerotol);
736 linc_postproj = rt->norm();
743 obj.hessVec(*Hn, n, x, zerotol);
745 con.applyAdjointHessian(*cxxvec, l, n, x, zerotol);
748 obj.hessVec(*Hto, *t_orig, x, zerotol);
750 con.applyAdjointHessian(*cxxvec, l, *t_orig, x, zerotol);
757 obj.update(*xtrial,UpdateType::Trial,state_->iter);
758 con.update(*xtrial,UpdateType::Trial,state_->iter);
759 f_new = obj.value(*xtrial, zerotol);
760 obj.gradient(gf_new, *xtrial, zerotol);
761 con.value(c_new, *xtrial, zerotol);
763 computeLagrangeMultiplier(l_new, *xtrial, gf_new, con, os);
766 part_pred = - Wg.dot(*t_orig);
771 part_pred -= n.apply(*gfJl);
774 part_pred -= half*n.apply(*Hn);
777 part_pred -= half*t_orig->apply(*Hto);
779 ltemp->axpy(-one, l);
782 part_pred -= Jnc->apply(*ltemp);
784 if ( part_pred < -half*penalty_*(c_normsquared-Jnc_normsquared) ) {
785 penalty_ = ( -two * part_pred / (c_normsquared-Jnc_normsquared) ) + beta;
788 pred = part_pred + penalty_*(c_normsquared-Jnc_normsquared);
793 rpred = - rt->apply(*ltemp) - penalty_ * rt->dot(*rt) - two * penalty_ * rt->dot(*Jnc);
801 tangtol = tangtol_start;
807 if ( t_orig->norm()/xtemp->norm() < tntmax_ ) {
811 t_m_tCP->set(*t_orig);
812 t_m_tCP->scale(-one);
814 if ((t_m_tCP->norm() > 0) && try_tCP) {
817 os <<
" ---> now trying tangential Cauchy point\n";
825 os <<
" ---> recomputing quasi-normal step and re-solving tangential subproblem\n";
830 obj.update(x, UpdateType::Trial, state_->iter);
831 con.update(x, UpdateType::Trial, state_->iter);
845 lmhtol_ *= tol_red_all;
846 qntol_ *= tol_red_all;
847 pgtol_ *= tol_red_all;
848 projtol_ *= tol_red_all;
849 tangtol_ *= tol_red_all;
852 computeQuasinormalStep(n, c, x, zeta_*Delta_, con, os);
854 solveTangentialSubproblem(t, tCP, Wg, x, g, n, l, Delta_, obj, con, os);
855 totalIterCG_ += iterCG_;
870 if (std::abs(fdiff / (f+em24)) < tol_fdiff) {
876 ared = fdiff + (c.apply(l) - c_new.apply(l_new)) + penalty_*(c.dot(c) - c_new.dot(c_new));
889 std::ios_base::fmtflags osFlags(os.flags());
890 os << std::scientific << std::setprecision(6);
891 os <<
"\n Trial step info ...\n";
892 os <<
" n_norm = " << nnorm_ <<
"\n";
893 os <<
" t_norm = " << tnorm_ <<
"\n";
894 os <<
" s_norm = " << snorm_ <<
"\n";
895 os <<
" xtrial_norm = " << xtrial->norm() <<
"\n";
896 os <<
" f_old = " << f <<
"\n";
897 os <<
" f_trial = " << f_new <<
"\n";
898 os <<
" f_old-f_trial = " << f-f_new <<
"\n";
899 os <<
" ||c_old|| = " << c.norm() <<
"\n";
900 os <<
" ||c_trial|| = " << c_new.norm() <<
"\n";
901 os <<
" ||Jac*t_preproj|| = " << linc_preproj <<
"\n";
902 os <<
" ||Jac*t_postproj|| = " << linc_postproj <<
"\n";
903 os <<
" ||t_tilde||/||t|| = " << t_orig->norm() / t.norm() <<
"\n";
904 os <<
" ||t_tilde||/||n+t|| = " << t_orig->norm() / snorm_ <<
"\n";
905 os <<
" # projections = " << num_proj <<
"\n";
906 os <<
" penalty param = " << penalty_ <<
"\n";
907 os <<
" ared = " << ared_ <<
"\n";
908 os <<
" pred = " << pred_ <<
"\n";
909 os <<
" ared/pred = " << ared_/pred_ <<
"\n";
915template<
typename Real>
916void CompositeStepAlgorithm<Real>::updateRadius(Vector<Real> &x,
918 const Vector<Real> &s,
919 Objective<Real> &obj,
920 Constraint<Real> &con,
930 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
933 Ptr<Vector<Real> > g = gvec_->clone();
934 Ptr<Vector<Real> > ajl = gvec_->clone();
935 Ptr<Vector<Real> > gl = gvec_->clone();
936 Ptr<Vector<Real> > c = cvec_->clone();
941 if ((std::abs(ared_) < em12) && std::abs(pred_) < em12) {
947 Delta_ = std::max(seven*snorm_, Delta_);
949 else if (ratio >= zp8) {
950 Delta_ = std::max(two*snorm_, Delta_);
952 obj.update(x,UpdateType::Accept,state_->iter);
953 con.update(x,UpdateType::Accept,state_->iter);
957 Delta_ = half*std::max(nnorm_, tnorm_);
958 obj.update(x,UpdateType::Revert,state_->iter);
959 con.update(x,UpdateType::Revert,state_->iter);
963 Real val = obj.value(x, zerotol);
965 obj.gradient(*g, x, zerotol);
966 computeLagrangeMultiplier(l, x, *g, con, os);
967 con.applyAdjointJacobian(*ajl, l, x, zerotol);
968 gl->set(*g); gl->plus(*ajl);
970 con.value(*c, x, zerotol);
972 state_->gradientVec->set(*gl);
973 state_->constraintVec->set(*c);
976 state_->gnorm = gl->norm();
977 state_->cnorm = c->norm();
979 state_->snorm = snorm_;
986template<
typename Real>
987void CompositeStepAlgorithm<Real>::initialize(Vector<Real> &x,
988 const Vector<Real> &g,
990 const Vector<Real> &c,
991 Objective<Real> &obj,
992 Constraint<Real> &con,
994 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
995 TypeE::Algorithm<Real>::initialize(x,g,l,c);
1007 Ptr<Vector<Real> > ajl = gvec_->clone();
1008 Ptr<Vector<Real> > gl = gvec_->clone();
1011 obj.update(x,UpdateType::Initial,state_->iter);
1012 state_->value = obj.value(x, zerotol);
1014 con.update(x,UpdateType::Initial,state_->iter);
1015 con.value(*cvec_, x, zerotol);
1016 state_->cnorm = cvec_->norm();
1018 obj.gradient(*gvec_, x, zerotol);
1021 computeLagrangeMultiplier(l, x, *gvec_, con, os);
1022 con.applyAdjointJacobian(*ajl, l, x, zerotol);
1023 gl->set(*gvec_); gl->plus(*ajl);
1025 state_->gnorm = gl->norm();
1029template<
typename Real>
1030void CompositeStepAlgorithm<Real>::run(Vector<Real> &x,
1031 const Vector<Real> &g,
1032 Objective<Real> &obj,
1033 Constraint<Real> &econ,
1035 const Vector<Real> &eres,
1036 std::ostream &outStream) {
1038 initialize(x, g, emul, eres, obj, econ, outStream);
1041 if (verbosity_ > 0) writeOutput(outStream,
true);
1044 Ptr<Vector<Real> > s = x.clone();
1046 while (status_->check(*state_)) {
1047 computeTrial(*s, x, emul, obj, econ, outStream);
1048 updateRadius(x, emul, *s, obj, econ, outStream);
1051 if (verbosity_ > 0) writeOutput(outStream, printHeader_);
1054 if (verbosity_ > 0) TypeE::Algorithm<Real>::writeExitStatus(outStream);
1058template<
typename Real>
1059void CompositeStepAlgorithm<Real>::writeHeader(std::ostream& os)
const {
1060 std::ios_base::fmtflags osFlags(os.flags());
1062 os << std::string(144,
'-') << std::endl;
1063 os <<
"Composite Step status output definitions" << std::endl << std::endl;
1064 os <<
" iter - Number of iterates (steps taken)" << std::endl;
1065 os <<
" fval - Objective function value" << std::endl;
1066 os <<
" cnorm - Norm of the constraint violation" << std::endl;
1067 os <<
" gLnorm - Norm of the gradient of the Lagrangian" << std::endl;
1068 os <<
" snorm - Norm of the step" << std::endl;
1069 os <<
" delta - Trust-region radius" << std::endl;
1070 os <<
" nnorm - Norm of the quasinormal step" << std::endl;
1071 os <<
" tnorm - Norm of the tangential step" << std::endl;
1072 os <<
" #fval - Number of times the objective was computed" << std::endl;
1073 os <<
" #grad - Number of times the gradient was computed" << std::endl;
1074 os <<
" iterCG - Number of projected CG iterations" << std::endl;
1075 os <<
" flagCG - Flag returned by projected CG" << std::endl;
1076 os <<
" accept - Acceptance flag for the trial step" << std::endl;
1077 os <<
" linsys - Number of augmented solver calls/iterations" << std::endl;
1078 os << std::string(144,
'-') << std::endl;
1081 os << std::setw(6) << std::left <<
"iter";
1082 os << std::setw(15) << std::left <<
"fval";
1083 os << std::setw(15) << std::left <<
"cnorm";
1084 os << std::setw(15) << std::left <<
"gLnorm";
1085 os << std::setw(15) << std::left <<
"snorm";
1086 os << std::setw(10) << std::left <<
"delta";
1087 os << std::setw(10) << std::left <<
"nnorm";
1088 os << std::setw(10) << std::left <<
"tnorm";
1089 os << std::setw(8) << std::left <<
"#fval";
1090 os << std::setw(8) << std::left <<
"#grad";
1091 os << std::setw(8) << std::left <<
"iterCG";
1092 os << std::setw(8) << std::left <<
"flagCG";
1093 os << std::setw(8) << std::left <<
"accept";
1094 os << std::setw(8) << std::left <<
"linsys";
1100template<
typename Real>
1101void CompositeStepAlgorithm<Real>::writeName(std::ostream& os)
const {
1102 std::ios_base::fmtflags osFlags(os.flags());
1103 os << std::endl <<
"Composite-Step Trust-Region Solver (Type E, Equality Constraints)";
1109template<
typename Real>
1110void CompositeStepAlgorithm<Real>::writeOutput(std::ostream& os,
const bool print_header)
const {
1111 std::ios_base::fmtflags osFlags(os.flags());
1112 os << std::scientific << std::setprecision(6);
1113 if (state_->iter == 0) writeName(os);
1114 if (print_header) writeHeader(os);
1115 if (state_->iter == 0 ) {
1117 os << std::setw(6) << std::left << state_->iter;
1118 os << std::setw(15) << std::left << state_->value;
1119 os << std::setw(15) << std::left << state_->cnorm;
1120 os << std::setw(15) << std::left << state_->gnorm;
1121 os << std::setw(15) << std::left <<
"---";
1122 os << std::setw(10) << std::left <<
"---";
1123 os << std::setw(10) << std::left <<
"---";
1124 os << std::setw(10) << std::left <<
"---";
1125 os << std::setw(8) << std::left <<
"---";
1126 os << std::setw(8) << std::left <<
"---";
1127 os << std::setw(8) << std::left <<
"---";
1128 os << std::setw(8) << std::left <<
"---";
1129 os << std::setw(8) << std::left <<
"---";
1130 os << std::setw(8) << std::left <<
"---";
1135 os << std::setw(6) << std::left << state_->iter;
1136 os << std::setw(15) << std::left << state_->value;
1137 os << std::setw(15) << std::left << state_->cnorm;
1138 os << std::setw(15) << std::left << state_->gnorm;
1139 os << std::setw(15) << std::left << state_->snorm;
1140 os << std::scientific << std::setprecision(2);
1141 os << std::setw(10) << std::left << Delta_;
1142 os << std::setw(10) << std::left << nnorm_;
1143 os << std::setw(10) << std::left << tnorm_;
1144 os << std::scientific << std::setprecision(6);
1145 os << std::setw(8) << std::left << state_->nfval;
1146 os << std::setw(8) << std::left << state_->ngrad;
1147 os << std::setw(8) << std::left << iterCG_;
1148 os << std::setw(8) << std::left << flagCG_;
1149 os << std::setw(8) << std::left << flagAC_;
1150 os << std::left << totalCallLS_ <<
"/" << totalIterLS_;
1157template<
typename Real>
1159int CompositeStepAlgorithm<Real>::sgn(T val)
const {
1160 return (T(0) < val) - (val < T(0));
1164template<
typename Real>
1165void CompositeStepAlgorithm<Real>::printInfoLS(
const std::vector<Real> &res, std::ostream &os)
const {
1167 std::ios_base::fmtflags osFlags(os.flags());
1168 os << std::scientific << std::setprecision(8);
1169 os << std::endl <<
" Augmented System Solver:" << std::endl;
1170 os <<
" True Residual" << std::endl;
1171 for (
unsigned j=0; j<res.size(); j++) {
1172 os <<
" " << std::left << std::setw(14) << res[j] << std::endl;
1180template<
typename Real>
1181Real CompositeStepAlgorithm<Real>::setTolOSS(
const Real intol)
const {
1182 return tolOSSfixed_ ? tolOSS_ : intol;
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
CompositeStepAlgorithm(ParameterList &list)