16#ifndef dealii_line_minimization_h
17#define dealii_line_minimization_h
51 template <
typename NumberType>
52 std_cxx17::optional<NumberType>
54 const NumberType f_low,
55 const NumberType g_low,
56 const NumberType x_hi,
57 const NumberType f_hi);
68 template <
typename NumberType>
69 std_cxx17::optional<NumberType>
71 const NumberType f_low,
72 const NumberType g_low,
73 const NumberType x_hi,
74 const NumberType f_hi,
75 const NumberType g_hi);
84 template <
typename NumberType>
85 std_cxx17::optional<NumberType>
87 const NumberType f_low,
88 const NumberType g_low,
89 const NumberType x_hi,
90 const NumberType f_hi,
91 const NumberType x_rec,
92 const NumberType f_rec);
105 template <
typename NumberType>
108 const NumberType f_low,
109 const NumberType g_low,
110 const NumberType x_hi,
111 const NumberType f_hi,
112 const NumberType g_hi,
116 const std::pair<NumberType, NumberType> bounds);
122 template <
typename NumberType>
125 const NumberType f_low,
126 const NumberType g_low,
127 const NumberType x_hi,
128 const NumberType f_hi,
129 const NumberType g_hi,
133 const std::pair<NumberType, NumberType> bounds);
322 template <
typename NumberType>
323 std::pair<NumberType, unsigned int>
325 const std::function<std::pair<NumberType, NumberType>(
const NumberType x)>
330 NumberType(
const NumberType x_low,
331 const NumberType f_low,
332 const NumberType g_low,
333 const NumberType x_hi,
334 const NumberType f_hi,
335 const NumberType g_hi,
339 const std::pair<NumberType, NumberType> bounds)> &interpolate,
341 const NumberType eta = 0.9,
342 const NumberType mu = 0.01,
343 const NumberType a_max = std::numeric_limits<NumberType>::max(),
344 const unsigned int max_evaluations = 20,
345 const bool debug_output =
false);
354 template <
typename NumberType>
355 std_cxx17::optional<NumberType>
363 const NumberType denom = (2. * g1 * x2 - 2. * g1 * x1 - 2. * f2 + 2. * f1);
367 return (g1 * (x2 * x2 - x1 * x1) + 2. * (f1 - f2) * x1) / denom;
372 template <
typename NumberType>
373 std_cxx17::optional<NumberType>
382 const NumberType beta1 = g1 + g2 - 3. * (f1 - f2) / (x1 - x2);
383 const NumberType s = beta1 * beta1 - g1 * g2;
388 const NumberType denom =
389 x1 < x2 ? g2 - g1 + 2. * beta2 : g1 - g2 + 2. * beta2;
393 return x1 < x2 ? x2 - (x2 - x1) * (g2 + beta2 - beta1) / denom :
394 x1 - (x1 - x2) * (g1 + beta2 - beta1) / denom;
399 template <
typename NumberType>
400 std_cxx17::optional<NumberType>
418 const NumberType x2_shift = x2 - x1;
419 const NumberType x3_shift = x3 - x1;
420 const NumberType r1 = f2 - f1 - g1 * x2_shift;
421 const NumberType r2 = f3 - f1 - g1 * x3_shift;
422 const NumberType denom =
423 Utilities::fixed_power<2>(x2_shift * x3_shift) * (x2_shift - x3_shift);
427 const NumberType A = (r1 * Utilities::fixed_power<2>(x3_shift) -
428 r2 * Utilities::fixed_power<2>(x2_shift)) /
430 const NumberType B = (r2 * Utilities::fixed_power<3>(x2_shift) -
431 r1 * Utilities::fixed_power<3>(x3_shift)) /
433 const NumberType &
C = g1;
436 const NumberType radical = B * B - A *
C * 3;
440 return x1 + (-B +
std::sqrt(radical)) / (A * 3);
445 template <
typename NumberType>
456 const std::pair<NumberType, NumberType> bounds)
465 std_cxx17::optional<NumberType> res =
cubic_fit(x1, f1, g1, x2, f2, g2);
466 if (res && *res >= bounds.first && *res <= bounds.second)
471 if (res && *res >= bounds.first && *res <= bounds.second)
476 return (bounds.first + bounds.second) * 0.5;
481 template <
typename NumberType>
492 const std::pair<NumberType, NumberType> bounds)
501 std_cxx17::optional<NumberType> res =
505 if (res && *res >= bounds.first && *res <= bounds.second)
510 if (res && *res >= bounds.first && *res <= bounds.second)
515 return (bounds.first + bounds.second) * 0.5;
520 template <
typename NumberType>
521 std::pair<NumberType, unsigned int>
523 const std::function<std::pair<NumberType, NumberType>(
const NumberType x)>
528 NumberType(
const NumberType x_low,
529 const NumberType f_low,
530 const NumberType g_low,
531 const NumberType x_hi,
532 const NumberType f_hi,
533 const NumberType g_hi,
537 const std::pair<NumberType, NumberType> bounds)> &choose,
539 const NumberType eta,
541 const NumberType a_max,
542 const unsigned int max_evaluations,
543 const bool debug_output)
554 const NumberType tau1 = 9.;
559 const NumberType tau2 = 0.1;
560 const NumberType tau3 = 0.5;
562 const NumberType g0_abs =
std::abs(g0);
563 const NumberType f_min = f0 + a_max * mu * g0;
567 const auto w1 = [&](
const NumberType a,
const NumberType f) {
568 return f <= f0 + a * mu * g0;
573 const auto w2 = [&](
const NumberType g) {
580 const NumberType x = std::numeric_limits<NumberType>::signaling_NaN();
581 NumberType a_lo = x, f_lo = x, g_lo = x;
582 NumberType a_hi = x, f_hi = x, g_hi = x;
583 NumberType ai = x, fi = x, gi = x;
588 NumberType f_prev, g_prev, a_prev;
594 while (i < max_evaluations)
596 const auto fgi = func(ai);
602 deallog <<
"Bracketing phase: " << i << std::endl
603 << ai <<
' ' << fi <<
' ' << gi <<
' ' << w1(ai, fi) <<
' '
604 << w2(gi) <<
' ' << f_min << std::endl;
607 if (fi <= f_min || ai == a_max)
610 deallog <<
"Reached the maximum step size." << std::endl;
611 return std::make_pair(ai, i);
615 (fi >= f_prev && i > 1))
630 deallog <<
"Satisfied both Wolfe conditions during Bracketing."
634 return std::make_pair(ai, i);
651 std::make_pair(2. * ai - a_prev,
652 std::min(a_max, ai + tau1 * (ai - a_prev)));
669 "Could not find the initial bracket within the given number of iterations."));
686 if (
std::abs(a_lo) > std::numeric_limits<NumberType>::epsilon() &&
687 std::abs(a_hi) > std::numeric_limits<NumberType>::epsilon())
695 while (i < max_evaluations)
697 const NumberType a_lo_safe = a_lo + tau2 * (a_hi - a_lo);
698 const NumberType a_hi_safe = a_hi - tau3 * (a_hi - a_lo);
699 const auto bounds = std::minmax(a_lo_safe, a_hi_safe);
702 a_lo, f_lo, g_lo, a_hi, f_hi, g_hi, a_rec, f_rec, g_rec, bounds);
704 const std::pair<NumberType, NumberType> fgi = func(ai);
710 deallog <<
"Sectioning phase: " << i << std::endl
711 << a_lo <<
' ' << f_lo <<
' ' << g_lo <<
' ' << w1(a_lo, f_lo)
712 <<
' ' << w2(g_lo) << std::endl
713 << a_hi <<
' ' << f_hi <<
' ' << g_hi <<
' ' << w1(a_hi, f_hi)
714 <<
' ' << w2(g_hi) << std::endl
715 << ai <<
' ' << fi <<
' ' << gi <<
' ' << w1(ai, fi) <<
' '
716 << w2(gi) << std::endl;
718 if (!w1(ai, fi) || fi >= f_lo)
734 deallog <<
"Satisfied both Wolfe conditions." << std::endl;
736 return std::make_pair(ai, i);
739 if (gi * (a_hi - a_lo) >= 0)
768 "Could not could complete the sectioning phase within the given number of iterations."));
769 return std::make_pair(std::numeric_limits<NumberType>::signaling_NaN(), i);
void add(const T &element)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::pair< NumberType, unsigned int > line_search(const std::function< std::pair< NumberType, NumberType >(const NumberType x)> &func, const NumberType f0, const NumberType g0, const std::function< NumberType(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi, const NumberType g_hi, const FiniteSizeHistory< NumberType > &x_rec, const FiniteSizeHistory< NumberType > &f_rec, const FiniteSizeHistory< NumberType > &g_rec, const std::pair< NumberType, NumberType > bounds)> &interpolate, const NumberType a1, const NumberType eta=0.9, const NumberType mu=0.01, const NumberType a_max=std::numeric_limits< NumberType >::max(), const unsigned int max_evaluations=20, const bool debug_output=false)
std_cxx17::optional< NumberType > cubic_fit(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi, const NumberType g_hi)
NumberType poly_fit(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi, const NumberType g_hi, const FiniteSizeHistory< NumberType > &x_rec, const FiniteSizeHistory< NumberType > &f_rec, const FiniteSizeHistory< NumberType > &g_rec, const std::pair< NumberType, NumberType > bounds)
NumberType poly_fit_three_points(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi, const NumberType g_hi, const FiniteSizeHistory< NumberType > &x_rec, const FiniteSizeHistory< NumberType > &f_rec, const FiniteSizeHistory< NumberType > &g_rec, const std::pair< NumberType, NumberType > bounds)
std_cxx17::optional< NumberType > cubic_fit_three_points(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi, const NumberType x_rec, const NumberType f_rec)
std_cxx17::optional< NumberType > quadratic_fit(const NumberType x_low, const NumberType f_low, const NumberType g_low, const NumberType x_hi, const NumberType f_hi)
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)