GlobiPack  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
GlobiPack_BrentsLineSearch_def.hpp
1 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // GlobiPack: Collection of Scalar 1D globalizaton utilities
6 // Copyright (2009) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
39 //
40 // ***********************************************************************
41 // @HEADER
42 */
43 
44 #ifndef GLOBIPACK_BRENTS_LINE_SEARCH_DEF_HPP
45 #define GLOBIPACK_BRENTS_LINE_SEARCH_DEF_HPP
46 
47 
48 #include "GlobiPack_BrentsLineSearch_decl.hpp"
49 #include "Teuchos_TabularOutputter.hpp"
50 
51 
52 namespace GlobiPack {
53 
54 
55 // Constructor/Initializers/Accessors
56 
57 
58 template<typename Scalar>
60 {}
61 
62 
63 template<typename Scalar>
66 {
67  return bracket_;
68 }
69 
70 
71 template<typename Scalar>
74 {
75  return brentsMin_;
76 }
77 
78 
79 // Overridden from ParameterListAcceptor
80 
81 
82 template<class Scalar>
84  RCP<ParameterList> const& paramList
85  )
86 {
87  //typedef ScalarTraits<Scalar> ST; // unused
88  namespace BLSU = BrentsLineSearchUtils;
89  using Teuchos::sublist;
90  paramList->validateParametersAndSetDefaults(*this->getValidParameters());
91  bracket_.setParameterList(sublist(paramList, BLSU::bracket_name, true));
92  brentsMin_.setParameterList(sublist(paramList, BLSU::minimize_name, true));
93  setMyParamList(paramList);
94 }
95 
96 
97 template<class Scalar>
100 {
101  namespace BLSU = BrentsLineSearchUtils;
102  static RCP<const ParameterList> validPL;
103  if (is_null(validPL)) {
106  pl->sublist(BLSU::bracket_name).setParameters(
107  *bracket_.getValidParameters()
108  ).disableRecursiveValidation();
109  pl->sublist(BLSU::minimize_name).setParameters(
110  *brentsMin_.getValidParameters()
111  ).disableRecursiveValidation();
112  validPL = pl;
113  }
114  return validPL;
115 }
116 
117 
118 // Overrridden from LineSearchBase
119 
120 
121 template<typename Scalar>
123 {
124  return false;
125 }
126 
127 
128 template<typename Scalar>
130 {
131  return false;
132 }
133 
134 
135 template<typename Scalar>
137  const MeritFunc1DBase<Scalar> &phi,
138  const PointEval1D<Scalar> &point_k,
139  const Ptr<PointEval1D<Scalar> > &point_kp1,
140  const Ptr<int> &numIters
141  ) const
142 {
143 
144  using Teuchos::as;
145  using Teuchos::OSTab;
146  using Teuchos::outArg;
147  using Teuchos::inOutArg;
148 
149 #ifdef TEUCHOS_DEBUG
150  typedef ScalarTraits<Scalar> ST;
151  typedef PointEval1D<Scalar> PE1D;
152 
153  TEUCHOS_ASSERT_EQUALITY(point_k.alpha, ST::zero());
154  TEUCHOS_ASSERT_INEQUALITY(point_k.phi, !=, PE1D::valNotGiven());
155  TEUCHOS_ASSERT_EQUALITY(point_k.Dphi, PE1D::valNotGiven());
156  TEUCHOS_ASSERT(!is_null(point_kp1));
157  TEUCHOS_ASSERT_INEQUALITY(point_kp1->alpha, >, ST::zero());
158  TEUCHOS_ASSERT_INEQUALITY(point_kp1->phi, !=, PE1D::valNotGiven());
159  TEUCHOS_ASSERT_EQUALITY(point_kp1->Dphi, PE1D::valNotGiven());
160 #endif
161 
162  const RCP<Teuchos::FancyOStream> out = this->getOStream();
163  bracket_.setOStream(out);
164  brentsMin_.setOStream(out);
165 
166  *out << "\nStarting bracketing and brents 1D minimization linesearch ...\n";
167 
168  OSTab tab(out);
169 
170  int totalNumIters = 0;
171 
172  PointEval1D<Scalar> p_l = point_k;
173  PointEval1D<Scalar> &p_m = *point_kp1; // Update in place!
175 
176  bool success = true;
177 
178  // A) Bracket the minimum
179 
180  int numBracketIters = -1;
181 
182  const bool bracketSuccess = bracket_.bracketMinimum(
183  phi, inOutArg(p_l), inOutArg(p_m), outArg(p_u), outArg(numBracketIters) );
184 
185  if (!bracketSuccess) success = false;
186 
187  totalNumIters += numBracketIters;
188 
189  // B) Do approximate mimimization in the bracket
190 
191  if (bracketSuccess) {
192 
193  int numBrentsIters = -1;
194 
195  const bool brentsSuccess = brentsMin_.approxMinimize(
196  phi, p_l, inOutArg(p_m), p_u, outArg(numBrentsIters) );
197 
198  if (!brentsSuccess) success = false;
199 
200  totalNumIters += numBrentsIters;
201 
202  }
203 
204  // C) Overall success?
205 
206  if (!is_null(numIters))
207  *numIters = totalNumIters;
208 
209  return success;
210 
211 }
212 
213 
214 } // namespace GlobiPack
215 
216 
217 #endif // GLOBIPACK_BRENTS_LINE_SEARCH_DEF_HPP
virtual bool requiresBaseDeriv() const
Returns true.
Scalar phi
The value of the merit function phi(alpha).
void setParameterList(RCP< ParameterList > const &paramList)
bool is_null(const boost::shared_ptr< T > &p)
Scalar Dphi
The value of the derivative of the merit function Dphi(alpha).
basic_OSTab< char > OSTab
Simple concrete class that implements a 1D algorithm to mimimize a 1D function.
Simple concrete class that implements a 1D algorithm to bracket the minimum of a 1D merit function...
virtual bool doLineSearch(const MeritFunc1DBase< Scalar > &phi, const PointEval1D< Scalar > &point_k, const Ptr< PointEval1D< Scalar > > &point_kp1, const Ptr< int > &numIters) const
const Brents1DMinimization< Scalar > & brentsMin() const
For unit testing only .
#define TEUCHOS_ASSERT_INEQUALITY(val1, comp, val2)
Scalar alpha
The value of the unknown alpha.
Represents the evaluation point of the merit function phi(alpha) and/or is derivative Dphi(alpha)...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
ParameterList & setParameters(const ParameterList &source)
const GoldenQuadInterpBracket< Scalar > & bracket() const
For unit testing only .
TypeTo as(const TypeFrom &t)
Base class for 1D merit fucntions used in globalization methods.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
RCP< const ParameterList > getValidParameters() const
virtual bool requiresDerivEvals() const
Returns false.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
BrentsLineSearch()
Construct with default parameters.