44 #ifndef GLOBIPACK_GOLDEN_BRACKET_QUAD_INTERP_DEF_HPP
45 #define GLOBIPACK_GOLDEN_BRACKET_QUAD_INTERP_DEF_HPP
58 template<
typename Scalar>
66 template<
typename Scalar>
72 setMyParamList(paramList);
76 template<
typename Scalar>
93 template<
typename Scalar>
122 const Scalar GOLDEN_RATIO = 1.618033988749895;
123 const Scalar SMALL_DIV = 1e-20;
124 const Scalar MAX_EXTRAP_FACTOR = 100.0;
125 const int MAX_TOTAL_ITERS = 30;
127 *out <<
"\nStarting golden quadratic interpolating bracketing of the minimum ...\n\n";
132 Scalar &alpha_l = pointLower->alpha, &phi_l = pointLower->phi;
133 Scalar &alpha_m = pointMiddle->alpha, &phi_m = pointMiddle->phi;
134 Scalar &alpha_u = pointUpper->alpha = ST::nan(), &phi_u = pointUpper->phi = ST::nan();
136 Scalar tmp = ST::nan(), q = ST::nan(), r = ST::nan();
138 const Scalar zero = ST::zero();
142 const Scalar goldinv = 1.0/(1.0+GOLDEN_RATIO);
144 TabularOutputter tblout(out);
146 tblout.pushFieldSpec(
"itr", TO::INT);
147 tblout.pushFieldSpec(
"alpha_l", TO::DOUBLE);
148 tblout.pushFieldSpec(
"alpha_m", TO::DOUBLE);
149 tblout.pushFieldSpec(
"alpha_u", TO::DOUBLE);
150 tblout.pushFieldSpec(
"phi_l", TO::DOUBLE);
151 tblout.pushFieldSpec(
"phi_m", TO::DOUBLE);
152 tblout.pushFieldSpec(
"phi_u", TO::DOUBLE);
153 tblout.pushFieldSpec(
"step type ", TO::STRING);
155 tblout.outputHeader();
159 std::string stepType =
"";
165 tblout.outputField(
"-");
166 tblout.outputField(alpha_l);
167 tblout.outputField(alpha_m);
168 tblout.outputField(
"-");
169 tblout.outputField(phi_l);
170 tblout.outputField(phi_m);
171 tblout.outputField(
"-");
172 tblout.outputField(
"start");
175 for (; icount < MAX_TOTAL_ITERS; ++icount) {
183 stepType =
"golden back";
186 alpha_m = goldinv * (alpha_u + GOLDEN_RATIO*alpha_l);
187 phi_m = computeValue<Scalar>(phi, alpha_m);
189 tblout.outputField(icount);
190 tblout.outputField(alpha_l);
191 tblout.outputField(alpha_m);
192 tblout.outputField(alpha_u);
193 tblout.outputField(phi_l);
194 tblout.outputField(phi_m);
195 tblout.outputField(phi_u);
196 tblout.outputField(stepType);
201 if (alpha_u == zero) {
204 alpha_u = alpha_m + (GOLDEN_RATIO-1.0) * (alpha_m-alpha_l);
205 phi_u = computeValue<Scalar>(phi, alpha_u);
212 bool bracketedMin =
false;
214 for (; icount < MAX_TOTAL_ITERS; ++icount) {
223 q = (phi_m-phi_l)*(alpha_m-alpha_u);
224 r = (phi_m-phi_u)*(alpha_m-alpha_l);
227 tmp = ST::magnitude(q-r);
228 tmp = (tmp > SMALL_DIV ? tmp : SMALL_DIV);
229 tmp = (q-r >= 0 ? tmp : -tmp);
232 alpha_m - (q*(alpha_m-alpha_u) - r*(alpha_m-alpha_l))/(2.0*tmp);
235 const Scalar alpha_lim = alpha_m + MAX_EXTRAP_FACTOR * (alpha_u-alpha_m);
238 bool skipToNextIter =
false;
239 Scalar phi_quad = ST::nan();
240 if ( (alpha_m-alpha_quad)*(alpha_quad-alpha_u) > zero ) {
241 phi_quad = computeValue<Scalar>(phi, alpha_quad);
242 if (phi_quad < phi_u) {
245 alpha_m = alpha_quad;
247 skipToNextIter =
true;
248 stepType =
"alpha_quad middle";
250 else if (phi_quad > phi_m) {
251 alpha_u = alpha_quad;
253 skipToNextIter =
true;
254 stepType =
"alpha_quad upper";
257 alpha_quad = alpha_u + GOLDEN_RATIO*(alpha_u-alpha_m);
258 phi_quad = computeValue<Scalar>(phi, alpha_quad);
262 if (!skipToNextIter) {
264 if ((alpha_u-alpha_quad)*(alpha_quad-alpha_lim) > zero) {
265 phi_quad = computeValue<Scalar>(phi, alpha_quad);
266 stepType =
"[alpha_u, alpha_lim]";
267 if (phi_quad < phi_u) {
269 alpha_u = alpha_quad;
270 alpha_quad = alpha_u + GOLDEN_RATIO*(alpha_u-alpha_m);
273 phi_quad = computeValue<Scalar>(phi, alpha_quad);
274 stepType =
"phi_quad < phi_u";
277 else if ((alpha_quad-alpha_lim)*(alpha_lim-alpha_u) >= zero ) {
278 alpha_quad = alpha_lim;
279 phi_quad = computeValue<Scalar>(phi, alpha_quad);
280 stepType =
"[alpha_lim, inf]";
283 alpha_quad = alpha_u + GOLDEN_RATIO*(alpha_u-alpha_m);
284 phi_quad = computeValue<Scalar>(phi, alpha_quad);
285 stepType =
"[0, alpha_m]";
293 alpha_u = alpha_quad;
298 tblout.outputField(icount);
299 tblout.outputField(alpha_l);
300 tblout.outputField(alpha_m);
301 tblout.outputField(alpha_u);
302 tblout.outputField(phi_l);
303 tblout.outputField(phi_m);
304 tblout.outputField(phi_u);
305 tblout.outputField(stepType);
310 if (icount >= MAX_TOTAL_ITERS) {
311 *out <<
"\nExceeded maximum number of iterations!.\n";
328 #endif // GLOBIPACK_GOLDEN_BRACKET_QUAD_INTERP_DEF_HPP
bool is_null(const boost::shared_ptr< T > &p)
GoldenQuadInterpBracket()
Construct with default parameters.
#define TEUCHOS_ASSERT_INEQUALITY(val1, comp, val2)
void setParameterList(RCP< ParameterList > const ¶mList)
bool bracketMinimum(const MeritFunc1DBase< Scalar > &phi, const Ptr< PointEval1D< Scalar > > &pointLower, const Ptr< PointEval1D< Scalar > > &pointMiddle, const Ptr< PointEval1D< Scalar > > &pointUpper, const Ptr< int > &numIters=Teuchos::null) const
Bracket the minimum of a 1D function.
basic_OSTab< char > OSTab
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)
TypeTo as(const TypeFrom &t)
RCP< const ParameterList > getValidParameters() const
Base class for 1D merit fucntions used in globalization methods.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)