47 #include "MoochoPack_TangentialStepWithInequStd_Step.hpp"
48 #include "MoochoPack_moocho_algo_conversion.hpp"
49 #include "MoochoPack_Exceptions.hpp"
50 #include "IterationPack_print_algorithm_step.hpp"
51 #include "ConstrainedOptPack_MatrixIdentConcat.hpp"
52 #include "AbstractLinAlgPack_MatrixOpOut.hpp"
53 #include "AbstractLinAlgPack_VectorOut.hpp"
54 #include "AbstractLinAlgPack_VectorAuxiliaryOps.hpp"
55 #include "AbstractLinAlgPack_LinAlgOpPack.hpp"
56 #include "Teuchos_dyn_cast.hpp"
61 T my_max(
const T& v1,
const T& v2 ) {
return v1 > v2 ? v1 : v2; }
64 T my_min(
const T& v1,
const T& v2 ) {
return v1 < v2 ? v1 : v2; }
67 namespace MoochoPack {
69 TangentialStepWithInequStd_Step::TangentialStepWithInequStd_Step(
70 const qp_solver_ptr_t &qp_solver
71 ,
const qp_tester_ptr_t &qp_tester
72 ,value_type warm_start_frac
74 ,
bool primal_feasible_point_error
75 ,
bool dual_feasible_point_error
77 :qp_solver_(qp_solver)
78 ,qp_tester_(qp_tester)
79 ,warm_start_frac_(warm_start_frac)
80 ,qp_testing_(qp_testing)
81 ,primal_feasible_point_error_(primal_feasible_point_error)
82 ,dual_feasible_point_error_(dual_feasible_point_error)
102 using LinAlgOpPack::V_MtV;
105 typedef VectorMutable::vec_mut_ptr_t vec_mut_ptr_t;
107 NLPAlgo &algo = rsqp_algo(_algo);
109 EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
110 EJournalOutputLevel ns_olevel = algo.algo_cntr().null_space_journal_output_level();
115 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
116 using IterationPack::print_algorithm_step;
117 print_algorithm_step( algo, step_poss, type, assoc_step_poss, out );
124 r = s.equ_decomp().size();
127 IterQuantityAccess<value_type>
128 &alpha_iq = s.alpha(),
131 IterQuantityAccess<VectorMutable>
135 *c_iq = m > 0 ? &s.c() : NULL,
136 *lambda_iq = m > 0 ? &s.lambda() : NULL,
139 &qp_grad_iq = s.qp_grad(),
144 IterQuantityAccess<MatrixOp>
147 *Uy_iq = (m > r) ? &s.Uy() : NULL;
148 IterQuantityAccess<MatrixSymOp>
150 IterQuantityAccess<ActSetStats>
151 &act_set_stats_iq = act_set_stats_(s);
154 VectorMutable *Ypy_k = (m ? &Ypy_iq.get_k(0) : NULL);
155 const MatrixOp &Z_k = Z_iq.get_k(0);
156 VectorMutable &pz_k = pz_iq.set_k(0);
157 VectorMutable &Zpz_k = Zpz_iq.set_k(0);
163 &qp_grad_k = ( qp_grad_iq.set_k(0) = rGf_iq.get_k(0) );
166 if( w_iq.updated_k(0) ) {
167 if(zeta_iq.updated_k(0))
168 Vp_StV( &qp_grad_k, zeta_iq.get_k(0), w_iq.get_k(0) );
170 Vp_V( &qp_grad_k, w_iq.get_k(0) );
179 bl = s.space_x().create_member(),
180 bu = s.space_x().create_member();
184 V_VmV( bl.get(), dl_iq.get_k(0), *Ypy_k );
186 V_VmV( bu.get(), du_iq.get_k(0), *Ypy_k );
189 *bl = dl_iq.get_k(0);
190 *bu = du_iq.get_k(0);
194 if( static_cast<int>(ns_olevel) >= static_cast<int>(PRINT_VECTORS) ) {
195 out <<
"\nqp_grad_k = \n" << qp_grad_k;
197 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_VECTORS) ) {
198 out <<
"\nbl = \n" << *bl;
199 out <<
"\nbu = \n" << *bu;
205 bool do_warm_start =
false;
206 if( act_set_stats_iq.updated_k(-1) ) {
207 if( static_cast<int>(olevel) >=
static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
208 out <<
"\nDetermining if the QP should use a warm start ...\n";
218 = ( num_adds == ActSetStats::NOT_KNOWN || num_active == 0
220 : my_max(((
double)(num_active)-num_adds-num_drops) / num_active, 0.0 ) );
221 do_warm_start = ( num_active > 0 && frac_same >= warm_start_frac() );
222 if( static_cast<int>(olevel) >=
static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
223 out <<
"\nnum_active = " << num_active;
225 out <<
"\nmax(num_active-num_adds-num_drops,0)/(num_active) = "
226 <<
"max("<<num_active<<
"-"<<num_adds<<
"-"<<num_drops<<
",0)/("<<num_active<<
") = "
232 out <<
"warm_start_frac = " << warm_start_frac();
235 out <<
"\nUse a warm start this time!\n";
237 out <<
"\nDon't use a warm start this time!\n";
246 if( do_warm_start ) {
250 nu_iq.set_k(0) = 0.0;
252 VectorMutable &nu_k = nu_iq.get_k(0);
268 const value_type qp_bnd_inf = NLP::infinite_bound();
270 const Vector &qp_g = qp_grad_k;
271 const MatrixSymOp &qp_G = rHL_iq.get_k(0);
272 const value_type qp_etaL = 0.0;
273 vec_mut_ptr_t qp_dL = Teuchos::null;
274 vec_mut_ptr_t qp_dU = Teuchos::null;
276 qp_E = Teuchos::null;
278 vec_mut_ptr_t qp_b = Teuchos::null;
279 vec_mut_ptr_t qp_eL = Teuchos::null;
280 vec_mut_ptr_t qp_eU = Teuchos::null;
282 qp_F = Teuchos::null;
284 vec_mut_ptr_t qp_f = Teuchos::null;
285 value_type qp_eta = 0.0;
286 VectorMutable &qp_d = pz_k;
287 vec_mut_ptr_t qp_nu = Teuchos::null;
288 vec_mut_ptr_t qp_mu = Teuchos::null;
289 vec_mut_ptr_t qp_Ed = Teuchos::null;
290 vec_mut_ptr_t qp_lambda = Teuchos::null;
314 const MatrixIdentConcat
315 *Zvr =
dynamic_cast<const MatrixIdentConcat*
>( &Z_k );
317 var_dep = Zvr ? Zvr->D_rng() : Range1D::Invalid,
318 var_indep = Zvr ? Zvr->I_rng() : Range1D();
323 = ( m ? (Ypy_indep=Ypy_k->sub_view(var_indep))->norm_inf() : 0.0);
325 if( (
int)olevel >= (int)PRINT_ALGORITHM_STEPS )
327 <<
"\nDetermine if we can use simple bounds on pz ...\n"
328 <<
" m = " << m << std::endl
329 <<
" dynamic_cast<const MatrixIdentConcat*>(&Z_k) = " << Zvr << std::endl
330 <<
" ||Ypy_k(var_indep)||inf = " << Ypy_indep_norm_inf << std::endl;
337 num_bounded( *bl->sub_view(var_dep), *bu->sub_view(var_dep), qp_bnd_inf )
345 ( Zvr != NULL && ( Ypy_indep_norm_inf == 0.0 || bounded_var_dep == 0 ) )
348 if( (
int)olevel >= (int)PRINT_ALGORITHM_STEPS )
350 << (use_simple_pz_bounds
351 ?
"\nUsing simple bounds on pz ...\n"
352 :
"\nUsing bounds on full Z*pz ...\n")
354 ?
"\nThere are finite bounds on dependent variables. Adding extra inequality constrints for D*pz ...\n"
355 :
"\nThere are no finite bounds on dependent variables. There will be no extra inequality constraints added on D*pz ...\n" ) ;
357 if( use_simple_pz_bounds ) {
359 qp_dL = bl->sub_view(var_indep);
360 qp_dU = bu->sub_view(var_indep);
361 qp_nu = nu_k.sub_view(var_indep);
362 if( m && bounded_var_dep ) {
365 qp_b = Ypy_k->sub_view(var_dep);
366 qp_eL = bl->sub_view(var_dep);
367 qp_eU = bu->sub_view(var_dep);
368 qp_mu = nu_k.sub_view(var_dep);
369 qp_Ed = Zpz_k.sub_view(var_dep);
375 else if( !use_simple_pz_bounds ) {
390 Range1D equ_undecomp = s.equ_undecomp();
391 if( m > r && m > 0 ) {
393 qp_f = s.space_c().sub_space(equ_undecomp)->create_member();
395 Vp_V( qp_f.get(), *c_iq->get_k(0).sub_view(equ_undecomp) );
398 qp_lambda = lambda_iq->set_k(0).sub_view(equ_undecomp);
403 QPSolverRelaxed::EOutputLevel
407 qp_olevel = QPSolverRelaxed::PRINT_NONE;
409 case PRINT_BASIC_ALGORITHM_INFO:
410 qp_olevel = QPSolverRelaxed::PRINT_NONE;
412 case PRINT_ALGORITHM_STEPS:
413 qp_olevel = QPSolverRelaxed::PRINT_BASIC_INFO;
415 case PRINT_ACTIVE_SET:
416 qp_olevel = QPSolverRelaxed::PRINT_ITER_SUMMARY;
419 qp_olevel = QPSolverRelaxed::PRINT_ITER_VECTORS;
421 case PRINT_ITERATION_QUANTITIES:
422 qp_olevel = QPSolverRelaxed::PRINT_EVERY_THING;
433 qp_solver().infinite_bound(qp_bnd_inf);
436 qp_solver().solve_qp(
437 int(olevel) ==
int(PRINT_NOTHING) ? NULL : &out
439 ,( algo.algo_cntr().check_results()
440 ? QPSolverRelaxed::RUN_TESTS : QPSolverRelaxed::NO_TESTS )
441 ,qp_g, qp_G, qp_etaL, qp_dL.get(), qp_dU.get()
442 ,qp_E.
get(), qp_trans_E, qp_E.
get() ? qp_b.get() : NULL
443 ,qp_E.
get() ? qp_eL.get() : NULL, qp_E.
get() ? qp_eU.get() : NULL
444 ,qp_F.
get(), qp_trans_F, qp_F.
get() ? qp_f.get() : NULL
448 ,qp_mu.get(), qp_E.
get() ? qp_Ed.get() : NULL
449 ,qp_F.
get() ? qp_lambda.get() : NULL
456 std::ostringstream omsg;
457 bool throw_qp_failure =
false;
458 if( qp_testing() == QP_TEST
459 || ( qp_testing() == QP_TEST_DEFAULT && algo.algo_cntr().check_results() ) )
461 if(
int(olevel) >= int(PRINT_ALGORITHM_STEPS) ) {
462 out <<
"\nChecking the optimality conditions of the reduced QP subproblem ...\n";
464 if(!qp_tester().check_optimality_conditions(
465 solution_type,qp_solver().infinite_bound()
466 ,
int(olevel) ==
int(PRINT_NOTHING) ? NULL : &out
467 ,
int(olevel) >=
int(PRINT_VECTORS) ?
true :
false
468 ,
int(olevel) >=
int(PRINT_ITERATION_QUANTITIES) ?
true :
false
469 ,qp_g, qp_G, qp_etaL, qp_dL.get(), qp_dU.get()
470 ,qp_E.
get(), qp_trans_E, qp_E.
get() ? qp_b.get() : NULL
471 ,qp_E.
get() ? qp_eL.get() : NULL, qp_E.
get() ? qp_eU.get() : NULL
472 ,qp_F.
get(), qp_trans_F, qp_F.
get() ? qp_f.get() : NULL
476 ,qp_mu.get(), qp_E.
get() ? qp_Ed.get() : NULL
477 ,qp_F.
get() ? qp_lambda.get() : NULL
481 omsg <<
"\n*** Alert! at least one of the QP optimality conditions did not check out.\n";
482 if( static_cast<int>(olevel) >=
static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
485 throw_qp_failure =
true;
492 if( !use_simple_pz_bounds ) {
495 else if( use_simple_pz_bounds ) {
497 *Zpz_k.sub_view(var_indep) = pz_k;
498 if( m && !bounded_var_dep ) {
509 qp_solver_stats_(s).set_k(0) = qp_solver().get_qp_stats();
512 const value_type eps = std::numeric_limits<value_type>::epsilon();
513 if( fabs(qp_eta - 0.0) > eps ) {
514 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
516 <<
"\n*** Alert! the QP was infeasible (eta = "<<qp_eta<<
"). Cutting back Ypy_k = (1.0 - eta)*Ypy_k ...\n";
518 Vt_S( Ypy_k, 1.0 - qp_eta );
522 eta_iq.set_k(0) = qp_eta;
527 switch(solution_type) {
528 case QPSolverStats::OPTIMAL_SOLUTION:
530 case QPSolverStats::PRIMAL_FEASIBLE_POINT:
533 <<
"\n*** Alert! the returned QP solution is PRIMAL_FEASIBLE_POINT but not optimal!\n";
534 if( primal_feasible_point_error() )
536 <<
"\n*** primal_feasible_point_error == true, this is an error!\n";
537 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
540 throw_qp_failure = primal_feasible_point_error();
543 case QPSolverStats::DUAL_FEASIBLE_POINT:
546 <<
"\n*** Alert! the returned QP solution is DUAL_FEASIBLE_POINT"
547 <<
"\n*** but not optimal so we cut back the step ...\n";
548 if( dual_feasible_point_error() )
550 <<
"\n*** dual_feasible_point_error == true, this is an error!\n";
551 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
559 zero = s.space_x().create_member(0.0),
560 d_tmp = s.space_x().create_member();
561 V_VpV( d_tmp.get(), *Ypy_k, Zpz_k );
562 const std::pair<value_type,value_type>
565 u = my_min( u_steps.first, 1.0 );
566 alpha_iq.set_k(0) = u;
567 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
568 out <<
"\nFinding u s.t. dl <= u*(Ypy_k+Zpz_k) <= du\n"
569 <<
"max step length u = " << u << std::endl
570 <<
"alpha_k = u = " << alpha_iq.get_k(0) << std::endl;
572 throw_qp_failure = dual_feasible_point_error();
575 case QPSolverStats::SUBOPTIMAL_POINT:
578 <<
"\n*** Alert!, the returned QP solution is SUBOPTIMAL_POINT!\n";
579 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
582 throw_qp_failure =
true;
592 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
593 out <<
"\n||pz_k||inf = " << s.pz().get_k(0).norm_inf()
594 <<
"\nnu_k.nz() = " << s.nu().get_k(0).nz()
595 <<
"\nmax(|nu_k(i)|) = " << s.nu().get_k(0).norm_inf()
598 if( m > r ) out <<
"\n||lambda_k(undecomp)||inf = " << s.lambda().get_k(0).norm_inf();
599 out <<
"\n||Zpz_k||2 = " << s.Zpz().get_k(0).norm_2()
601 if(qp_eta > 0.0) out <<
"\n||Ypy||2 = " << s.Ypy().get_k(0).norm_2();
605 if( static_cast<int>(ns_olevel) >= static_cast<int>(PRINT_VECTORS) ) {
606 out <<
"\npz_k = \n" << s.pz().get_k(0);
608 out <<
"\nnu_k(var_indep) = \n" << *s.nu().get_k(0).sub_view(var_indep);
611 if( static_cast<int>(ns_olevel) >= static_cast<int>(PRINT_VECTORS) ) {
613 out <<
"\nZpz(var_indep)_k = \n" << *s.Zpz().get_k(0).sub_view(var_indep);
617 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_VECTORS) ) {
619 out <<
"\nZpz(var_dep)_k = \n" << *s.Zpz().get_k(0).sub_view(var_dep);
620 out <<
"\nZpz_k = \n" << s.Zpz().get_k(0);
624 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_VECTORS) ) {
625 out <<
"\nnu_k = \n" << s.nu().get_k(0);
627 out <<
"\nnu_k(var_dep) = \n" << *s.nu().get_k(0).sub_view(var_dep);
629 out <<
"\nlambda_k(equ_undecomp) = \n" << *s.lambda().get_k(0).sub_view(equ_undecomp);
630 if(qp_eta > 0.0) out <<
"\nYpy = \n" << s.Ypy().get_k(0);
633 if( qp_eta == 1.0 ) {
635 <<
"TangentialStepWithInequStd_Step::do_step(...) : Error, a QP relaxation parameter\n"
636 <<
"of eta = " << qp_eta <<
" was calculated and therefore it must be assumed\n"
637 <<
"that the NLP's constraints are infeasible\n"
638 <<
"Throwing an InfeasibleConstraints exception!\n";
639 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
642 throw InfeasibleConstraints(omsg.str());
645 if( throw_qp_failure )
646 throw QPFailure( omsg.str(), qp_solver().get_qp_stats() );
653 ,
poss_type assoc_step_poss, std::ostream& out,
const std::string& L
657 << L <<
"*** Calculate the null-space step by solving a constrained QP\n"
658 << L <<
"qp_grad_k = rGf_k\n"
659 << L <<
"if w_k is updated then set qp_grad_k = qp_grad_k + zeta_k * w_k\n"
660 << L <<
"bl = dl_k - Ypy_k\n"
661 << L <<
"bu = du_k - Ypy_k\n"
662 << L <<
"etaL = 0.0\n"
663 << L <<
"*** Determine if we can use simple bounds on pz or not\n"
664 << L <<
"if num_bounded(bl_k(var_dep),bu_k(var_dep)) > 0 then\n"
665 << L <<
" bounded_var_dep = true\n"
667 << L <<
" bounded_var_dep = false\n"
671 << L <<
" ( Z_k is a variable reduction null space matrix\n"
673 << L <<
" ( ||Ypy_k(var_indep)||inf == 0 or bounded_var_dep==false ) )\n"
675 << L <<
" use_simple_pz_bounds = true\n"
677 << L <<
" use_simple_pz_bounds = false\n"
679 << L <<
"*** Setup QP arguments\n"
680 << L <<
"qp_g = qp_grad_k\n"
681 << L <<
"qp_G = rHL_k\n"
682 << L <<
"if (use_simple_pz_bounds == true) then\n"
683 << L <<
" qp_dL = bl(var_indep), qp_dU = bu(var_indep))\n"
684 << L <<
" if (m > 0) then\n"
685 << L <<
" qp_E = Z_k.D, qp_b = Ypy_k(var_dep)\n"
686 << L <<
" qp_eL = bl(var_dep), qp_eU = bu(var_dep)\n"
688 << L <<
"elseif (use_simple_pz_bounds == false) then\n"
689 << L <<
" qp_dL = -inf, qp_dU = +inf\n"
690 << L <<
" qp_E = Z_k, qp_b = Ypy_k\n"
691 << L <<
" qp_eL = bl, qp_eU = bu\n"
693 << L <<
"if (m > r) then\n"
694 << L <<
" qp_F = V_k, qp_f = Uy_k*py_k + c_k(equ_undecomp)\n"
696 << L <<
" qp_F = empty, qp_f = empty\n"
698 << L <<
"Given active-set statistics (act_set_stats_km1)\n"
699 << L <<
" frac_same = max(num_active-num_adds-num_drops,0)/(num_active)\n"
700 << L <<
"Use a warm start when frac_same >= warm_start_frac\n"
701 << L <<
"Solve the following QP to compute qp_d, qp_eta, qp_Ed = qp_E * qp_d\n"
702 << L <<
",qp_nu, qp_mu and qp_lambda (" <<
typeName(qp_solver()) <<
"):\n"
703 << L <<
" min qp_g' * qp_d + 1/2 * qp_d' * qp_G * qp_d + M(eta)\n"
704 << L <<
" qp_d <: R^(n-r)\n"
706 << L <<
" etaL <= eta\n"
707 << L <<
" qp_dL <= qp_d <= qp_dU [qp_nu]\n"
708 << L <<
" qp_eL <= qp_E * qp_d + (1-eta)*qp_b <= qp_eU [qp_mu]\n"
709 << L <<
" qp_F * d_qp + (1-eta) * qp_f = 0 [qp_lambda]\n"
710 << L <<
"if (qp_testing==QP_TEST) or (fd_deriv_testing==QP_TEST_DEFAULT\n"
711 << L <<
"and check_results==true) then\n"
712 << L <<
" Check the optimality conditions of the above QP\n"
713 << L <<
" if the optimality conditions do not check out then\n"
714 << L <<
" set throw_qp_failure = true\n"
717 << L <<
"*** Set the solution to the QP subproblem\n"
718 << L <<
"pz_k = qp_d\n"
719 << L <<
"eta_k = qp_eta\n"
720 << L <<
"if (use_simple_pz_bounds == true) then\n"
721 << L <<
" nu_k(var_dep) = qp_mu, nu_k(var_indep) = qp_nu\n"
722 << L <<
" Zpz_k(var_dep) = qp_Ed, Zpz_k(var_indep) = pz_k\n"
723 << L <<
"elseif (use_simple_pz_bounds == false) then\n"
724 << L <<
" nu_k = qp_mu\n"
725 << L <<
" Zpz_k = qp_Ed\n"
727 << L <<
"if m > r then\n"
728 << L <<
" lambda_k(equ_undecomp) = qp_lambda\n"
730 << L <<
"if (eta_k > 0) then set Ypy_k = (1-eta_k) * Ypy_k\n"
731 << L <<
"if QP solution is suboptimal then\n"
732 << L <<
" throw_qp_failure = true\n"
733 << L <<
"elseif QP solution is primal feasible (not optimal) then\n"
734 << L <<
" throw_qp_failure = primal_feasible_point_error\n"
735 << L <<
"elseif QP solution is dual feasible (not optimal) then\n"
736 << L <<
" find max u s.t.\n"
737 << L <<
" dl_k <= u*(Ypy_k+Zpz_k) <= du_k\n"
738 << L <<
" alpha_k = u\n"
739 << L <<
" throw_qp_failure = dual_feasible_point_error\n"
741 << L <<
"if (eta_k == 1.0) then\n"
742 << L <<
" The constraints are infeasible!\n"
743 << L <<
" throw InfeasibleConstraints(...)\n"
745 << L <<
"if (throw_qp_failure == true) then\n"
746 << L <<
" throw QPFailure(...)\n"
void Vt_S(VectorMutable *v_lhs, const value_type &alpha)
void Vp_StV(VectorMutable *v_lhs, const value_type &alpha, const Vector &v_rhs)
T_To & dyn_cast(T_From &from)
Thrown if a the QP failed and was not corrected.
rSQP Algorithm control class.
void V_StV(VectorMutable *v_lhs, value_type alpha, const V &V_rhs)
void print_step(const Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss, std::ostream &out, const std::string &leading_str) const
Class for storing statistics about the changes in the active set of an SQP algorithm.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual std::ostream & journal_out() const
void V_VpV(VectorMutable *v_lhs, const V1 &V1_rhs1, const V2 &V2_rhs2)
Reduced space SQP state encapsulation interface.
AlgorithmTracker & track()
size_type num_bounded(const Vector &xl, const Vector &xu, value_type inf_bound)
bool do_step(Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
std::pair< value_type, value_type > max_near_feas_step(const Vector &x, const Vector &d, const Vector &xl, const Vector &xu, value_type max_bnd_viol)
NLPAlgoState & rsqp_state()
<<std aggr>="">> members for algo_cntr
void V_VmV(VectorMutable *v_lhs, const V1 &V1_rhs1, const V2 &V2_rhs2)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
std::string typeName(const T &t)
void Vp_V(VectorMutable *v_lhs, const V &V_rhs)