Optimist  0.0.0
A C++ library for optimization
Loading...
Searching...
No Matches
Broyden.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *\
2 * Copyright (c) 2025, Davide Stocco, Mattia Piazza and Enrico Bertolazzi. *
3 * *
4 * The Optimist project is distributed under the BSD 2-Clause License. *
5 * *
6 * Davide Stocco Mattia Piazza Enrico Bertolazzi *
7 * University of Trento University of Trento University of Trento *
8 * davide.stocco@unitn.it mattia.piazza@unitn.it enrico.bertolazzi@unitn.it *
9\* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
10
11#pragma once
12
13#ifndef OPTIMIST_ROOTFINDER_BROYDEN_HH
14#define OPTIMIST_ROOTFINDER_BROYDEN_HH
15
17
18namespace Optimist
19{
20 namespace RootFinder
21 {
22
23 /*\
24 | ____ _
25 | | __ ) _ __ ___ _ _ __| | ___ _ __
26 | | _ \| '__/ _ \| | | |/ _` |/ _ \ '_ \
27 | | |_) | | | (_) | |_| | (_| | __/ | | |
28 | |____/|_| \___/ \__, |\__,_|\___|_| |_|
29 | |___/
30 \*/
31
40 template <typename Real, Integer N>
41 class Broyden : public QuasiNewton<Real, N, Broyden<Real, N>>
42 {
43 public:
44 static constexpr bool requires_function{true};
45 static constexpr bool requires_first_derivative{true};
46 static constexpr bool requires_second_derivative{false};
47
49
52
56 using Method = enum class Method : Integer {GOOD = 0, BAD = 1, COMBINED = 2};
57
58 private:
59 Method m_method{Method::COMBINED};
60
61 public:
66
71 std::string name_impl() const
72 {
73 std::ostringstream os;
74 os << "Broyden";
75 if (this->m_method == Method::GOOD) {
76 os << "Good";
77 } else if (this->m_method == Method::BAD) {
78 os << "Bad";
79 } else if (this->m_method == Method::COMBINED) {
80 os << "Combined";
81 }
82 return os.str();
83 }
84
89 Method method() const {return this->m_method;}
90
95 void method(Method t_method) {this->m_method = t_method;}
96
100 void enable_good_method() {this->m_method = Method::GOOD;}
101
105 void enable_bad_method() {this->m_method = Method::BAD;}
106
110 void enable_combined_method() {this->m_method = Method::COMBINED;}
111
116 void set_method(Method t_method) {this->m_method = t_method;}
117
129 Vector const & delta_x_old, Vector const & delta_function_old, Matrix const & jacobian_old,
130 Vector const & delta_x_new, Vector const & delta_function_new, Vector const & /*function_new*/,
131 Matrix & jacobian_new
132 ) {
133 Vector tmp_1(jacobian_old * delta_function_new);
134 Real tmp_2{delta_function_new.squaredNorm()};
135 // Selection criteria: |(dx_new'*dx_old) / (dx_new'*J_old*dF_new)| < |(dF_new'*dF_old) / (dF_new'*dF_new)|
136 if (this->m_method == Method::COMBINED || this->m_method == Method::GOOD || this->iterations() < Integer(2) ||
137 std::abs(delta_x_new.transpose() * delta_x_old) / std::abs(delta_x_new.transpose() * tmp_1)
138 < std::abs(delta_function_new.transpose() * delta_function_old) / tmp_2) {
139 // Broyden's Good solver: J_new = J_old - (J_old*dF_new-dx_new) / (C'*dF_new)*C', where C = J_old'*dx_new;
140 Vector C_g(jacobian_old.transpose() * delta_x_new);
141 jacobian_new = jacobian_old - (tmp_1 - delta_x_new) / (C_g.transpose() * delta_function_new) * C_g.transpose();
142 } else {
143 // Broyden's Bad solver: J_new = J_old - (J_old*dF_old-dx_new) / (C'*dF_old)*C', where C = dF_old;
144 jacobian_new = jacobian_old - (tmp_1 - delta_x_new) / tmp_2 * delta_function_old.transpose();
145 }
146 }
147
148 }; // class Broyden
149
150 } // namespace RootFinder
151
152} // namespace Optimist
153
154#endif // OPTIMIST_ROOTFINDER_BROYDEN_HH
#define OPTIMIST_BASIC_CONSTANTS(Real)
Definition Optimist.hh:71
static constexpr bool requires_function
Definition Broyden.hh:44
void method(Method t_method)
Definition Broyden.hh:95
Broyden()
Definition Broyden.hh:65
static constexpr bool requires_first_derivative
Definition Broyden.hh:45
void enable_good_method()
Definition Broyden.hh:100
void enable_bad_method()
Definition Broyden.hh:105
static constexpr bool requires_second_derivative
Definition Broyden.hh:46
Method m_method
Definition Broyden.hh:59
std::string name_impl() const
Definition Broyden.hh:71
void enable_combined_method()
Definition Broyden.hh:110
Method method() const
Definition Broyden.hh:89
void update_impl(Vector const &delta_x_old, Vector const &delta_function_old, Matrix const &jacobian_old, Vector const &delta_x_new, Vector const &delta_function_new, Vector const &, Matrix &jacobian_new)
Definition Broyden.hh:128
enum class Method :Integer {GOOD=0, BAD=1, COMBINED=2} Method
Definition Broyden.hh:56
void set_method(Method t_method)
Definition Broyden.hh:116
Class container for the QuasiNewton's method.
Definition QuasiNewton.hh:42
typename SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >::InputType Vector
Definition RootFinder.hh:61
typename SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >::FirstDerivativeType Matrix
Definition RootFinder.hh:64
Integer iterations() const
Definition SolverBase.hh:329
Namespace for the Optimist library.
Definition Optimist.hh:88
OPTIMIST_DEFAULT_INTEGER_TYPE Integer
The Integer type as used for the API.
Definition Optimist.hh:96