Optimist  0.0.0
A C++ library for optimization
Loading...
Searching...
No Matches
NelderMead.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_OPTIMIZER_NELDER_MEAD_HH
14#define OPTIMIST_OPTIMIZER_NELDER_MEAD_HH
15
16#include "Optimist/Optimizer.hh"
17
18namespace Optimist
19{
20 namespace Optimizer
21 {
22
23 /*\
24 | _ _ _ _ __ __ _
25 | | \ | | ___| | __| | ___ _ __| \/ | ___ __ _ __| |
26 | | \| |/ _ \ |/ _` |/ _ \ '__| |\/| |/ _ \/ _` |/ _` |
27 | | |\ | __/ | (_| | __/ | | | | | __/ (_| | (_| |
28 | |_| \_|\___|_|\__,_|\___|_| |_| |_|\___|\__,_|\__,_|
29 |
30 \*/
31
32
41 template <typename Real, Integer N>
42 class NelderMead : public Optimizer<Real, N, NelderMead<Real, N>>
43 {
44 public:
45 static constexpr bool requires_function{true};
46 static constexpr bool requires_first_derivative{false};
47 static constexpr bool requires_second_derivative{false};
48
53
54 using Simplex = std::vector<Vector>;
55
56 private:
59 Real m_alpha{1.0};
60 Real m_gamma{2.0};
61 Real m_rho{0.5};
62 Real m_sigma{0.5};
63
64 public:
69 this->m_s.resize(N);
70 }
71
76 std::string name_impl() const {return "NelderMead";}
77
83 void initialize(const Vector & x_ini, Real delta)
84 {
85 this->m_s.clear();
86 this->m_s.push_back(x_ini);
87 Vector x;
88 for (Integer i{0}; i < N; ++i) {
89 x = x_ini;
90 x(i) += delta;
91 this->m_s.push_back(x);
92 }
93 this->m_c = this->centroid();
94 }
95
101 {
102 return 1.0/Real(N) * std::accumulate(this->m_s.begin(), this->m_s.end(), this->m_c);
103 }
104
108 void sort(FunctionWrapper function)
109 {
110 std::sort(this->m_s.begin(), this->m_s.end(),
111 [function, this] (const Vector & a, const Vector & b) {
112 Real function_a, function_b;
113 this->evaluate_function(function, a, function_a);
114 this->evaluate_function(function, b, function_b);
115 return function_a < function_b;
116 }
117 );
118 this->m_c = this->centroid();
119 }
120
126 {
127 return this->m_c + this->m_alpha * (this->m_c - this->m_s.back());
128 }
129
135 Vector expansion(const Vector & x) const
136 {
137 return this->m_c + this->m_gamma * (x - this->m_c);
138 }
139
146 {
147 return this->m_c + this->m_rho * (x - this->m_c);
148 }
149
155 {
156 return this->m_c + this->m_rho * (this->m_s.back() - this->m_c);
157 }
158
162 void shrink()
163 {
164 for (Integer i{1}; i < static_cast<Integer>(this->m_s.size()); ++i) {
165 this->m_s.at(i) = this->m_s.at(1) + this->m_sigma * (this->m_s.at(i) - this->m_s.at(1));
166 }
167 }
168
177 bool solve_impl(FunctionWrapper function, Vector const & x_ini, Vector & x_sol)
178 {
179 // Setup internal variables
180 this->reset();
181
182 // Print header
183 if (this->m_verbose) {this->header();}
184
185 // Initialize variables
186 Vector x_r, x_c, x_e, simplex_size;
187 Real function_x_r, function_x_1, function_x_last, function_x_e, function_x_c;
188
189 // Set initial iteration
190 this->initialize(x_ini, 100.0);
191
192 // Algorithm iterations
193 for (this->m_iterations = static_cast<Integer>(1); this->m_iterations < this->m_max_iterations; ++this->m_iterations)
194 {
195 // Sort the simplex points based on the function values
196 this->sort(function);
197
198 // Store trace
199 this->store_trace(this->m_s.back());
200
201 // Compute the reflection point
202 x_r = this->reflection();
203 this->evaluate_function(function, x_r, function_x_r);
204 this->evaluate_function(function, this->m_s.at(1), function_x_1);
205 this->evaluate_function(function, this->m_s.back(), function_x_last);
206
207 if (function_x_1 <= function_x_r && function_x_r < function_x_last)
208 {
209 this->m_s.back() = x_r;
210 if (this->m_verbose) {this->info(function_x_1, "Reflection");}
211 }
212 else if (function_x_r < function_x_1)
213 {
214 // Expansion
215 x_e = this->expansion(x_r);
216 this->evaluate_function(function, x_e, function_x_e);
217 if (function_x_e < function_x_r) {
218 this->m_s.back() = x_e;
219 } else {
220 this->m_s.back() = x_r;
221 }
222 if (this->m_verbose) {this->info(function_x_1, "Expansion");}
223 }
224 else if (function_x_r < function_x_last)
225 {
226 // Outward contraction
227 x_c = this->outward_contraction(x_r);
228 this->evaluate_function(function, x_c, function_x_c);
229 if (function_x_c < function_x_r) {
230 this->m_s.back() = x_c;
231 } else {
232 this->shrink();
233 }
234 if (this->m_verbose) {this->info(function_x_1, "Outward contraction");}
235 }
236 else
237 {
238 // Inward contraction
239 x_c = this->inward_contraction();
240 this->evaluate_function(function, x_c, function_x_c);
241 if (function_x_c < function_x_last) {
242 this->m_s.back() = x_c;
243 } else {
244 this->shrink();
245 }
246 if (this->m_verbose) {this->info(function_x_1, "Inward contraction");}
247 }
248 simplex_size = this->m_s.back() - this->m_s.front();
249 if (simplex_size.norm() < this->m_tolerance) {
250 this->m_converged = true;
251 break;
252 }
253 }
254
255 // Print bottom
256 if (this->m_verbose) {this->bottom();}
257
258 // Convergence data
259 x_sol = this->m_s.at(0);
260 return this->m_converged;
261 }
262
263 }; // class NelderMead
264
265 } // namespace Optimizer
266
267} // namespace Optimist
268
269#endif // OPTIMIST_OPTIMIZER_NELDER_MEAD_HH
Vector inward_contraction() const
Definition NelderMead.hh:154
static constexpr bool requires_first_derivative
Definition NelderMead.hh:46
Vector outward_contraction(const Vector &x) const
Definition NelderMead.hh:145
Vector expansion(const Vector &x) const
Computes the expansion point.
Definition NelderMead.hh:135
Simplex m_s
Definition NelderMead.hh:57
Real m_gamma
Definition NelderMead.hh:60
typename Optimizer< Real, N, NelderMead< Real, N > >::FunctionWrapper FunctionWrapper
Definition NelderMead.hh:51
std::string name_impl() const
Definition NelderMead.hh:76
Real m_sigma
Definition NelderMead.hh:62
Vector centroid() const
Definition NelderMead.hh:100
NelderMead()
Definition NelderMead.hh:68
Real m_alpha
Definition NelderMead.hh:59
typename Optimizer< Real, N, NelderMead< Real, N > >::Matrix Matrix
Definition NelderMead.hh:50
void sort(FunctionWrapper function)
Sort the points of the simplex based on the function values and updates the centroid.
Definition NelderMead.hh:108
Vector m_c
Definition NelderMead.hh:58
bool solve_impl(FunctionWrapper function, Vector const &x_ini, Vector &x_sol)
Definition NelderMead.hh:177
void initialize(const Vector &x_ini, Real delta)
Definition NelderMead.hh:83
typename Optimizer< Real, N, NelderMead< Real, N > >::Vector Vector
Definition NelderMead.hh:49
static constexpr bool requires_second_derivative
Definition NelderMead.hh:47
Vector reflection() const
Computes the reflection point.
Definition NelderMead.hh:125
Real m_rho
Definition NelderMead.hh:61
std::vector< Vector > Simplex
Definition NelderMead.hh:54
static constexpr bool requires_function
Definition NelderMead.hh:45
void shrink()
Definition NelderMead.hh:162
bool solve(FunctionWrapper function, Vector const &x_ini, Vector &x_sol)
Definition Optimizer.hh:160
void store_trace(const InputType &x)
Definition Solver.hh:750
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
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