MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MoochoPack_LineSearchFilter_Step.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef LINE_SEARCH_FILTER_STEP_H
43 #define LINE_SEARCH_FILTER_STEP_H
44 
45 #include <list>
46 
47 #include "MoochoPack_Types.hpp"
52 
54 
55 namespace MoochoPack {
56 
57 // structure for storing filter points
59 {
60 public:
61  FilterEntry( value_type new_f, value_type new_theta, int new_iter)
62  : f(new_f), theta(new_theta), iter(new_iter) {}
63 
66  int iter;
67 };
68 
69 // It is critical that an std::list be used because of the way iterators are
70 // used!
71 typedef std::list< FilterEntry > Filter_T;
72 
73 const std::string FILTER_IQ_STRING = "LS_FilterEntries";
74 
80  : public IterationPack::AlgorithmStep // doxygen needs full path
81 {
82 public:
83 
86 
89 
91 
94 
100 
106 
114 
120 
126 
132 
138 
144  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, theta_small_fact );
145 
151 
157 
163 
168  ,const std::string obj_iq_name = "f"
169  ,const std::string grad_obj_iq_name = "Gf"
170  ,const value_type &gamma_theta = 1e-5
171  ,const value_type &gamma_f = 1e-5
172  ,const value_type &f_min = F_MIN_UNBOUNDED
173  ,const value_type &gamma_alpha = 5e-2
174  ,const value_type &delta = 1e-4
175  ,const value_type &s_theta = 1.1
176  ,const value_type &s_f = 2.3
177  ,const value_type &theta_small_fact = 1e-4
178  ,const value_type &theta_max = 1e10
179  ,const value_type &eta_f = 1e-4
180  ,const value_type &back_track_frac = 0.5
181  );
182 
184 
188 
190 
194  bool do_step( Algorithm& algo, poss_type step_poss,
196  poss_type assoc_step_poss);
198  void print_step( const Algorithm& algo, poss_type step_poss,
200  poss_type assoc_step_poss, std::ostream& out,
201  const std::string& leading_str ) const;
203 
204 private:
205 
206  // /////////////////////////////
207  // Private data members
208 
209  // Private Data
210  CastIQMember<Filter_T> filter_;
211 
213  CastIQMember<value_type> obj_f_;
214 
216  CastIQMember<VectorMutable> grad_obj_f_;
217 
218  // nlp to use for calculations
220 
221  // /////////////////////////////
222  // Private member functions
223 
224  // Validate input parameters - fix if possible
225  bool ValidatePoint(
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
231  ) const;
232 
233  // Check that new point is not within taboo region of filter
235  const value_type f,
236  const value_type theta,
237  const AlgorithmState& s
238  ) const;
239 
240  // Check the Armijo condition on f
241  bool CheckArmijo(
242  const value_type Gf_t_dk,
243  const value_type alpha_k,
244  const IterQuantityAccess<value_type>& f_iq
245  ) const;
246 
247  // Check if f or c has sufficient reduction
249  const IterQuantityAccess<value_type>& f_iq,
250  const value_type gamma_f_used,
251  const value_type theta_kp1,
252  const value_type theta_k
253  ) const;
254 
255  // Calculate the new point given d and alpha
256  void UpdatePoint(
257  const VectorMutable& d,
258  value_type alpha,
259  IterQuantityAccess<VectorMutable>& x,
260  IterQuantityAccess<value_type>& f,
261  IterQuantityAccess<VectorMutable>* c,
262  IterQuantityAccess<VectorMutable>* h,
263  NLP& nlp
264  ) const;
265 
266  // Calculate the minimum alpha before hitting restoration phase
268  const value_type gamma_f_used,
269  const value_type Gf_t_dk,
270  const value_type theta_k,
271  const value_type theta_small
272  ) const;
273 
274  // Calculate the constraint norm
276  const IterQuantityAccess<VectorMutable>* c,
277  const IterQuantityAccess<VectorMutable>* h,
278  int k
279  ) const;
280 
281  // Calculate the value of gamma_k to use
283  const IterQuantityAccess<value_type> &f,
284  const value_type theta_k,
285  const EJournalOutputLevel olevel,
286  std::ostream &out
287  ) const;
288 
289  // decide if we should switch to Armijo for objective
291  const value_type Gf_t_dk,
292  const value_type alpha_k,
293  const value_type theta_k,
294  const value_type theta_small
295  ) const;
296 
297  // Update the filter from the last iteration
299 
300  // Update the filter from the last iteration and Augment it with
301  // the new point
302  void AugmentFilter(
303  const value_type gamma_f_used,
304  const value_type f_with_boundary,
305  const value_type theta_with_boundary,
307  const EJournalOutputLevel olevel,
308  std::ostream &out
309  ) const;
310 
311 }; // end class LineSearchFilter_Step
312 
313 } // end namespace MoochoPack
314 
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
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.
std::ostream * out
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
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