42 #ifndef LINE_SEARCH_FILTER_STEP_H
43 #define LINE_SEARCH_FILTER_STEP_H
55 namespace MoochoPack {
62 :
f(new_f),
theta(new_theta),
iter(new_iter) {}
168 ,
const std::string obj_iq_name =
"f"
169 ,
const std::string grad_obj_iq_name =
"Gf"
201 const std::string& leading_str )
const;
226 const IterQuantityAccess<VectorMutable>& x,
227 const IterQuantityAccess<value_type>& f,
228 const IterQuantityAccess<VectorMutable>* c,
229 const IterQuantityAccess<VectorMutable>* h,
230 const bool throw_excpt
237 const AlgorithmState& s
244 const IterQuantityAccess<value_type>& f_iq
249 const IterQuantityAccess<value_type>& f_iq,
257 const VectorMutable& d,
259 IterQuantityAccess<VectorMutable>& x,
260 IterQuantityAccess<value_type>& f,
261 IterQuantityAccess<VectorMutable>* c,
262 IterQuantityAccess<VectorMutable>* h,
276 const IterQuantityAccess<VectorMutable>* c,
277 const IterQuantityAccess<VectorMutable>* h,
283 const IterQuantityAccess<value_type> &f,
315 #endif // LINE_SEARCH_FILTER_STEP_H
Base type for all objects that perform steps in an Algorithm.
void UpdateFilter(IterationPack::AlgorithmState &s) const
const std::string FILTER_IQ_STRING
bool CheckFractionalReduction(const IterQuantityAccess< value_type > &f_iq, const value_type gamma_f_used, const value_type theta_kp1, const value_type theta_k) const
value_type CalculateAlphaMin(const value_type gamma_f_used, const value_type Gf_t_dk, const value_type theta_k, const value_type theta_small) const
bool CheckArmijo(const value_type Gf_t_dk, const value_type alpha_k, const IterQuantityAccess< value_type > &f_iq) const
~LineSearchFilter_Step()
Destructor.
CastIQMember< Filter_T > filter_
value_type CalculateGammaFUsed(const IterQuantityAccess< value_type > &f, const value_type theta_k, const EJournalOutputLevel olevel, std::ostream &out) const
LineSearchFilter_Step(Teuchos::RCP< NLPInterfacePack::NLP > nlp, const std::string obj_iq_name="f", const std::string grad_obj_iq_name="Gf", const value_type &gamma_theta=1e-5, const value_type &gamma_f=1e-5, const value_type &f_min=F_MIN_UNBOUNDED, const value_type &gamma_alpha=5e-2, const value_type &delta=1e-4, const value_type &s_theta=1.1, const value_type &s_f=2.3, const value_type &theta_small_fact=1e-4, const value_type &theta_max=1e10, const value_type &eta_f=1e-4, const value_type &back_track_frac=0.5)
Constructor.
bool ValidatePoint(const IterQuantityAccess< VectorMutable > &x, const IterQuantityAccess< value_type > &f, const IterQuantityAccess< VectorMutable > *c, const IterQuantityAccess< VectorMutable > *h, const bool throw_excpt) const
CastIQMember< VectorMutable > grad_obj_f_
ITeration quantity access for objective gradient.
bool do_step(Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
std::list< FilterEntry > Filter_T
value_type CalculateTheta_k(const IterQuantityAccess< VectorMutable > *c, const IterQuantityAccess< VectorMutable > *h, int k) const
FilterEntry(value_type new_f, value_type new_theta, int new_iter)
EJournalOutputLevel
enum for journal output.
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
Teuchos::RCP< NLPInterfacePack::NLP > nlp_
Abstacts a set of iteration quantities for an iterative algorithm.
void AugmentFilter(const value_type gamma_f_used, const value_type f_with_boundary, const value_type theta_with_boundary, IterationPack::AlgorithmState &s, const EJournalOutputLevel olevel, std::ostream &out) const
CastIQMember< value_type > obj_f_
Iteration quantity access for objective value.
bool ShouldSwitchToArmijo(const value_type Gf_t_dk, const value_type alpha_k, const value_type theta_k, const value_type theta_small) const
static value_type F_MIN_UNBOUNDED
STANDARD_MEMBER_COMPOSITION_MEMBERS(value_type, gamma_theta)
Feasibility decrease fraction.
bool CheckFilterAcceptability(const value_type f, const value_type theta, const AlgorithmState &s) const
AbstractLinAlgPack::value_type value_type
Acts as the central hub for an iterative algorithm.
void UpdatePoint(const VectorMutable &d, value_type alpha, IterQuantityAccess< VectorMutable > &x, IterQuantityAccess< value_type > &f, IterQuantityAccess< VectorMutable > *c, IterQuantityAccess< VectorMutable > *h, NLP &nlp) const
Filter line-search step class.