16#ifndef dealii_line_minimization_h
17#define dealii_line_minimization_h
48 template <
typename NumberType>
49 std_cxx17::optional<NumberType>
51 const NumberType f_low,
52 const NumberType g_low,
53 const NumberType x_hi,
54 const NumberType f_hi);
65 template <
typename NumberType>
66 std_cxx17::optional<NumberType>
68 const NumberType f_low,
69 const NumberType g_low,
70 const NumberType x_hi,
71 const NumberType f_hi,
72 const NumberType g_hi);
81 template <
typename NumberType>
82 std_cxx17::optional<NumberType>
84 const NumberType f_low,
85 const NumberType g_low,
86 const NumberType x_hi,
87 const NumberType f_hi,
88 const NumberType x_rec,
89 const NumberType f_rec);
102 template <
typename NumberType>
105 const NumberType f_low,
106 const NumberType g_low,
107 const NumberType x_hi,
108 const NumberType f_hi,
109 const NumberType g_hi,
113 const std::pair<NumberType, NumberType> bounds);
119 template <
typename NumberType>
122 const NumberType f_low,
123 const NumberType g_low,
124 const NumberType x_hi,
125 const NumberType f_hi,
126 const NumberType g_hi,
130 const std::pair<NumberType, NumberType> bounds);
319 template <
typename NumberType>
320 std::pair<NumberType, unsigned int>
322 const std::function<std::pair<NumberType, NumberType>(
const NumberType x)>
327 NumberType(
const NumberType x_low,
328 const NumberType f_low,
329 const NumberType g_low,
330 const NumberType x_hi,
331 const NumberType f_hi,
332 const NumberType g_hi,
336 const std::pair<NumberType, NumberType> bounds)> &interpolate,
338 const NumberType eta = 0.9,
339 const NumberType mu = 0.01,
340 const NumberType a_max = std::numeric_limits<NumberType>::max(),
341 const unsigned int max_evaluations = 20,
342 const bool debug_output =
false);
351 template <
typename NumberType>
352 std_cxx17::optional<NumberType>
360 const NumberType denom = (2. * g1 * x2 - 2. * g1 * x1 - 2. * f2 + 2. * f1);
364 return (g1 * (x2 * x2 - x1 * x1) + 2. * (f1 - f2) * x1) / denom;
369 template <
typename NumberType>
370 std_cxx17::optional<NumberType>
379 const NumberType beta1 = g1 + g2 - 3. * (f1 - f2) / (x1 - x2);
380 const NumberType s = beta1 * beta1 - g1 * g2;
385 const NumberType denom =
386 x1 < x2 ? g2 - g1 + 2. * beta2 : g1 - g2 + 2. * beta2;
390 return x1 < x2 ? x2 - (x2 - x1) * (g2 + beta2 - beta1) / denom :
391 x1 - (x1 - x2) * (g1 + beta2 - beta1) / denom;
396 template <
typename NumberType>
397 std_cxx17::optional<NumberType>
415 const NumberType x2_shift = x2 - x1;
416 const NumberType x3_shift = x3 - x1;
417 const NumberType r1 = f2 - f1 - g1 * x2_shift;
418 const NumberType r2 = f3 - f1 - g1 * x3_shift;
419 const NumberType denom =
420 std::pow(x2_shift * x3_shift, 2) * (x2_shift - x3_shift);
428 const NumberType &
C = g1;
431 const NumberType radical = B * B - A *
C * 3;
435 return x1 + (-B +
std::sqrt(radical)) / (A * 3);
440 template <
typename NumberType>
451 const std::pair<NumberType, NumberType> bounds)
460 std_cxx17::optional<NumberType> res =
cubic_fit(x1, f1, g1, x2, f2, g2);
461 if (res && *res >= bounds.first && *res <= bounds.second)
466 if (res && *res >= bounds.first && *res <= bounds.second)
471 return (bounds.first + bounds.second) * 0.5;
476 template <
typename NumberType>
487 const std::pair<NumberType, NumberType> bounds)
496 std_cxx17::optional<NumberType> res =
500 if (res && *res >= bounds.first && *res <= bounds.second)
505 if (res && *res >= bounds.first && *res <= bounds.second)
510 return (bounds.first + bounds.second) * 0.5;
515 template <
typename NumberType>
516 std::pair<NumberType, unsigned int>
518 const std::function<std::pair<NumberType, NumberType>(
const NumberType x)>
523 NumberType(
const NumberType x_low,
524 const NumberType f_low,
525 const NumberType g_low,
526 const NumberType x_hi,
527 const NumberType f_hi,
528 const NumberType g_hi,
532 const std::pair<NumberType, NumberType> bounds)> &choose,
534 const NumberType eta,
536 const NumberType a_max,
537 const unsigned int max_evaluations,
538 const bool debug_output)
549 const NumberType tau1 = 9.;
554 const NumberType tau2 = 0.1;
555 const NumberType tau3 = 0.5;
557 const NumberType g0_abs =
std::abs(g0);
558 const NumberType f_min = f0 + a_max * mu * g0;
562 const auto w1 = [&](
const NumberType a,
const NumberType f) {
563 return f <= f0 + a * mu * g0;
568 const auto w2 = [&](
const NumberType g) {
575 const NumberType x = std::numeric_limits<NumberType>::signaling_NaN();
576 NumberType a_lo = x, f_lo = x, g_lo = x;
577 NumberType a_hi = x, f_hi = x, g_hi = x;
578 NumberType ai = x, fi = x, gi = x;
583 NumberType f_prev, g_prev, a_prev;
589 while (i < max_evaluations)
591 const auto fgi = func(ai);
597 deallog <<
"Bracketing phase: " << i << std::endl
598 << ai <<
' ' << fi <<
' ' << gi <<
' ' << w1(ai, fi) <<
' '
599 << w2(gi) <<
' ' << f_min << std::endl;
602 if (fi <= f_min || ai == a_max)
605 deallog <<
"Reached the maximum step size." << std::endl;
606 return std::make_pair(ai, i);
610 (fi >= f_prev && i > 1))
625 deallog <<
"Satisfied both Wolfe conditions during Bracketing."
629 return std::make_pair(ai, i);
646 std::make_pair(2. * ai - a_prev,
647 std::min(a_max, ai + tau1 * (ai - a_prev)));
664 "Could not find the initial bracket within the given number of iterations."));
681 if (
std::abs(a_lo) > std::numeric_limits<NumberType>::epsilon() &&
682 std::abs(a_hi) > std::numeric_limits<NumberType>::epsilon())
690 while (i < max_evaluations)
692 const NumberType a_lo_safe = a_lo + tau2 * (a_hi - a_lo);
693 const NumberType a_hi_safe = a_hi - tau3 * (a_hi - a_lo);
694 const auto bounds = std::minmax(a_lo_safe, a_hi_safe);
697 a_lo, f_lo, g_lo, a_hi, f_hi, g_hi, a_rec, f_rec, g_rec, bounds);
699 const std::pair<NumberType, NumberType> fgi = func(ai);
705 deallog <<
"Sectioning phase: " << i << std::endl
706 << a_lo <<
' ' << f_lo <<
' ' << g_lo <<
' ' << w1(a_lo, f_lo)
707 <<
' ' << w2(g_lo) << std::endl
708 << a_hi <<
' ' << f_hi <<
' ' << g_hi <<
' ' << w1(a_hi, f_hi)
709 <<
' ' << w2(g_hi) << std::endl
710 << ai <<
' ' << fi <<
' ' << gi <<
' ' << w1(ai, fi) <<
' '
711 << w2(gi) << std::endl;
713 if (!w1(ai, fi) || fi >= f_lo)
729 deallog <<
"Satisfied both Wolfe conditions." << std::endl;
731 return std::make_pair(ai, i);
734 if (gi * (a_hi - a_lo) >= 0)
763 "Could not could complete the sectioning phase within the given number of iterations."));
764 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 > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)