Optimist  0.0.0
A C++ library for optimization
Loading...
Searching...
No Matches
QuasiNewton.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_QUASINEWTON_HH
14#define OPTIMIST_ROOTFINDER_QUASINEWTON_HH
15
17
18namespace Optimist
19{
20 namespace RootFinder
21 {
22
23 /*\
24 | ___ _ _ _ _
25 | / _ \ _ _ __ _ ___(_) \ | | _____ _| |_ ___ _ __
26 | | | | | | | |/ _` / __| | \| |/ _ \ \ /\ / / __/ _ \| '_ \
27 | | |_| | |_| | (_| \__ \ | |\ | __/\ V V /| || (_) | | | |
28 | \__\_\\__,_|\__,_|___/_|_| \_|\___| \_/\_/ \__\___/|_| |_|
29 |
30 \*/
31
40 template <typename Real, Integer N, typename DerivedSolver>
41 class QuasiNewton : public RootFinder<Real, N, DerivedSolver>
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
50 using Method = enum class Method : Integer {GOOD = 0, BAD = 1, COMBINED = 2};
55 using RootFinder<Real, N, DerivedSolver>::solve;
56
57 private:
58 Method m_method{Method::COMBINED};
59
60 public:
65
70 std::string name_impl() const
71 {
72 return static_cast<const DerivedSolver *>(this)->name_impl();
73 }
74
84 bool solve_impl(FunctionWrapper function, JacobianWrapper jacobian, Vector const & x_ini,
85 Vector & x_sol)
86 {
87 // Setup internal variables
88 this->reset();
89
90 // Print header
91 if (this->m_verbose) {this->header();}
92
93 // Initialize variables
94 bool damped;
95 Real residuals, step_norm;
96 Vector x_old, x_new, function_old, function_new, step_old, step_new, delta_x_old, delta_x_new,
97 delta_function_old, delta_function_new;
98 Matrix jacobian_old, jacobian_new;
99
100 // Set initial iteration
101 x_old = x_ini;
102 this->evaluate_function(function, x_old, function_old);
103 this->evaluate_jacobian(jacobian, x_old, jacobian_old);
104
105 // Algorithm iterations
106 Real tolerance_residuals{this->m_tolerance};
107 Real tolerance_step_norm{this->m_tolerance * this->m_tolerance};
108 for (this->m_iterations = static_cast<Integer>(1); this->m_iterations < this->m_max_iterations; ++this->m_iterations)
109 {
110 // Store trace
111 this->store_trace(x_old);
112
113 // Calculate step
114 step_old = -jacobian_old * function_old;
115
116 // Check convergence
117 residuals = function_old.norm();
118 step_norm = step_old.norm();
119 if (this->m_verbose) {this->info(residuals);}
120 if (residuals < tolerance_residuals || step_norm < tolerance_step_norm) {
121 this->m_converged = true;
122 break;
123 }
124
125 if (this->m_damped)
126 {
127 // Relax the iteration process
128 damped = this->damp(function, x_old, function_old, step_old, x_new, function_new, step_new);
130 "Optimist::RootFinder::QuasiNewton::solve(...): damping failed.");
131 } else {
132 // Update point
133 x_new = x_old + step_old;
134 this->evaluate_function(function, x_new, function_new);
135 }
136
137 // Update jacobian approximation
138 delta_x_new = x_new - x_old;
139 delta_function_new = function_new - function_old;
140 this->update(
141 delta_x_old, delta_function_old, jacobian_old, // Old step data
142 delta_x_new, delta_function_new, function_new, jacobian_new // New step data
143 );
144
145 // Update internal variables
146 x_old = x_new;
147 function_old = function_new;
148 delta_x_old = delta_x_new;
149 delta_function_old = delta_function_new;
150 step_old = step_new;
151 jacobian_old = jacobian_new;
152 }
153
154 // Print bottom
155 if (this->m_verbose) {this->bottom();}
156
157 // Convergence data
158 x_sol = x_old;
159 return this->m_converged;
160 }
161
172 void update(
173 Vector const & delta_x_old, Vector const & delta_function_old, Matrix const & jacobian_old,
174 Vector const & delta_x_new, Vector const & delta_function_new, Vector const & function_new,
175 Matrix & jacobian_new
176 ) {
177 static_cast<DerivedSolver *>(this)->update_impl(
178 delta_x_old, delta_function_old, jacobian_old, // Old step data
179 delta_x_new, delta_function_new, function_new, jacobian_new // New step data
180 );
181 }
182
183 }; // class QuasiNewton
184
185 } // namespace RootFinder
186
187} // namespace Optimist
188
189#endif // OPTIMIST_ROOTFINDER_QUASINEWTON_HH
#define OPTIMIST_ASSERT_WARNING(COND, MSG)
Definition Optimist.hh:60
#define OPTIMIST_BASIC_CONSTANTS(Real)
Definition Optimist.hh:70
void update(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 &function_new, Matrix &jacobian_new)
Definition QuasiNewton.hh:172
QuasiNewton()
Definition QuasiNewton.hh:64
bool solve_impl(FunctionWrapper function, JacobianWrapper jacobian, Vector const &x_ini, Vector &x_sol)
Definition QuasiNewton.hh:84
Method m_method
Definition QuasiNewton.hh:58
typename RootFinder< Real, N, DerivedSolver >::FunctionWrapper FunctionWrapper
Definition QuasiNewton.hh:53
static constexpr bool requires_first_derivative
Definition QuasiNewton.hh:45
enum class Method :Integer {GOOD=0, BAD=1, COMBINED=2} Method
Definition QuasiNewton.hh:50
std::string name_impl() const
Definition QuasiNewton.hh:70
typename RootFinder< Real, N, DerivedSolver >::JacobianWrapper JacobianWrapper
Definition QuasiNewton.hh:54
typename RootFinder< Real, N, DerivedSolver >::Matrix Matrix
Definition QuasiNewton.hh:52
static constexpr bool requires_second_derivative
Definition QuasiNewton.hh:46
typename RootFinder< Real, N, DerivedSolver >::Vector Vector
Definition QuasiNewton.hh:51
static constexpr bool requires_function
Definition QuasiNewton.hh:44
RootFinder()
Definition RootFinder.hh:83
typename Solver< Real, N, N, DerivedSolver >::FirstDerivativeWrapper JacobianWrapper
Definition RootFinder.hh:77
typename Solver< Real, N, N, DerivedSolver >::InputType Vector
Definition RootFinder.hh:69
void evaluate_jacobian(JacobianWrapper jacobian, const Vector &x, Matrix &out)
Definition RootFinder.hh:140
typename Solver< Real, N, N, DerivedSolver >::FunctionWrapper FunctionWrapper
Definition RootFinder.hh:76
bool solve(FunctionWrapper function, Vector const &x_ini, Vector &x_sol)
Definition RootFinder.hh:164
typename Solver< Real, N, N, DerivedSolver >::FirstDerivativeType Matrix
Definition RootFinder.hh:72
Real m_tolerance
Definition Solver.hh:92
bool m_converged
Definition Solver.hh:99
void store_trace(const InputType &x)
Definition Solver.hh:750
void header()
Definition Solver.hh:799
void reset()
Definition Solver.hh:693
bool damp(FunctionWrapper function, InputType const &x_old, InputType const &function_old, InputType const &step_old, InputType &x_new, InputType &function_new, InputType &step_new)
Definition Solver.hh:763
void bottom()
Definition Solver.hh:829
void evaluate_function(FunctionWrapper function, const InputType &x, OutputType &out)
Definition Solver.hh:710
void info(Real residuals, std::string const &notes="-")
Definition Solver.hh:853
Integer m_iterations
Definition Solver.hh:85
Integer m_max_iterations
Definition Solver.hh:86
Namespace for the Optimist library.
Definition Optimist.hh:87
OPTIMIST_DEFAULT_INTEGER_TYPE Integer
The Integer type as used for the API.
Definition Optimist.hh:95