16#ifndef dealii_solver_bfgs_h
17#define dealii_solver_bfgs_h
55template <
typename VectorType>
62 using Number =
typename VectorType::value_type;
109 const std::function<
Number(
const VectorType &x, VectorType &g)> &compute,
120 boost::signals2::connection
123 Number(
Number &f, VectorType &x, VectorType &g,
const VectorType &p)>
152 boost::signals2::connection
154 const std::function<
void(VectorType & g,
168 boost::signals2::signal<
169 Number(
Number &f, VectorType &x, VectorType &g,
const VectorType &p)>
175 boost::signals2::signal<void(VectorType & g,
185template <
typename VectorType>
187 const unsigned int max_history_size_,
188 const bool debug_output_)
189 : max_history_size(max_history_size_)
190 , debug_output(debug_output_)
195template <
typename VectorType>
197 const AdditionalData &data)
199 , additional_data(data)
204template <
class VectorType>
205boost::signals2::connection
208 Number(Number &f, VectorType &x, VectorType &g,
const VectorType &p)> &slot)
210 Assert(line_search_signal.empty(),
211 ExcMessage(
"One should not attach more than one line search signal."));
212 return line_search_signal.connect(slot);
217template <
class VectorType>
218boost::signals2::connection
220 const std::function<
void(VectorType & g,
224 Assert(preconditioner_signal.empty(),
226 "One should not attach more than one preconditioner signal."));
227 return preconditioner_signal.connect(slot);
232template <
typename VectorType>
235 const std::function<
typename VectorType::value_type(
const VectorType &x,
236 VectorType &f)> &compute,
246 bool first_step =
true;
250 if (line_search_signal.empty())
253 const auto default_line_min =
254 [&](Number &f, VectorType &x, VectorType &g,
const VectorType &p) {
256 const Number g0 = g * p;
259 "Function does not decrease along the current direction"));
267 Number df = f_prev - f;
268 Assert(first_step || df >= 0.,
269 ExcMessage(
"Function value is not decreasing"));
273 (first_step ? 1. :
std::min(1., -1.01 * 2. * df / g0));
278 const auto line_func =
279 [&](
const Number &x_line) -> std::pair<Number, Number> {
283 const Number g_line = g * p;
284 return std::make_pair(f, g_line);
288 const auto res = LineMinimization::line_search<Number>(
292 LineMinimization::poly_fit<Number>,
302 this->connect_line_search_slot(default_line_min);
311 VectorType g(x), p(x), y_k(x), s_k(x);
313 std::vector<Number> c1;
314 c1.reserve(additional_data.max_history_size);
329 conv = this->iteration_status(k, g.l2_norm(), x);
335 if (additional_data.debug_output)
336 deallog <<
"Iteration " << k <<
" history " << m << std::endl
337 <<
"f=" << f << std::endl;
343 for (
unsigned int i = 0; i < m; ++i)
345 c1[i] = rho[i] * (s[i] * p);
349 if (!preconditioner_signal.empty())
350 preconditioner_signal(p, s, y);
353 for (
int i = m - 1; i >= 0; --i)
356 const Number c2 = rho[i] * (y[i] * p);
357 p.add(c1[i] - c2, s[i]);
364 const Number alpha = line_search_signal(f, x, g, p)
369 if (additional_data.debug_output)
370 deallog <<
"Line search a=" << alpha <<
" f=" << f << std::endl;
374 const Number g_l2 = g.l2_norm();
375 conv = this->iteration_status(k, g_l2, x);
380 const Number curvature = s_k * y_k;
381 if (additional_data.debug_output)
382 deallog <<
"Curvature " << curvature << std::endl;
384 if (curvature > 0. && additional_data.max_history_size > 0)
388 rho.add(1. / curvature);
void add(const T &element)
typename VectorType::value_type Number
boost::signals2::connection connect_line_search_slot(const std::function< Number(Number &f, VectorType &x, VectorType &g, const VectorType &p)> &slot)
boost::signals2::connection connect_preconditioner_slot(const std::function< void(VectorType &g, const FiniteSizeHistory< VectorType > &s, const FiniteSizeHistory< VectorType > &y)> &slot)
void solve(const std::function< Number(const VectorType &x, VectorType &g)> &compute, VectorType &x)
boost::signals2::signal< void(VectorType &g, const FiniteSizeHistory< VectorType > &s, const FiniteSizeHistory< VectorType > &y)> preconditioner_signal
const AdditionalData additional_data
boost::signals2::signal< Number(Number &f, VectorType &x, VectorType &g, const VectorType &p)> line_search_signal
SolverBFGS(SolverControl &residual_control, const AdditionalData &data=AdditionalData())
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
AdditionalData(const unsigned int max_history_size=5, const bool debug_output=false)
unsigned int max_history_size