13#ifndef OPTIMIST_OPTIMIZER_NELDER_MEAD_HH
14#define OPTIMIST_OPTIMIZER_NELDER_MEAD_HH
41 template <
typename Real, Integer N>
76 std::string
name_impl()
const {
return "NelderMead";}
86 this->m_s.push_back(x_ini);
88 for (
Integer i{0}; i < N; ++i) {
91 this->m_s.push_back(x);
102 return 1.0/Real(N) * std::accumulate(this->m_s.begin(), this->m_s.end(), this->m_c);
110 std::sort(this->m_s.begin(), this->m_s.end(),
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;
127 return this->m_c + this->
m_alpha * (this->m_c - this->m_s.back());
137 return this->m_c + this->
m_gamma * (x - this->m_c);
147 return this->m_c + this->
m_rho * (x - this->m_c);
156 return this->m_c + this->
m_rho * (this->m_s.back() - this->m_c);
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));
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;
196 this->
sort(function);
207 if (function_x_1 <= function_x_r && function_x_r < function_x_last)
209 this->m_s.back() = x_r;
210 if (this->
m_verbose) {this->
info(function_x_1,
"Reflection");}
212 else if (function_x_r < function_x_1)
217 if (function_x_e < function_x_r) {
218 this->m_s.back() = x_e;
220 this->m_s.back() = x_r;
222 if (this->
m_verbose) {this->
info(function_x_1,
"Expansion");}
224 else if (function_x_r < function_x_last)
229 if (function_x_c < function_x_r) {
230 this->m_s.back() = x_c;
234 if (this->
m_verbose) {this->
info(function_x_1,
"Outward contraction");}
241 if (function_x_c < function_x_last) {
242 this->m_s.back() = x_c;
246 if (this->
m_verbose) {this->
info(function_x_1,
"Inward contraction");}
248 simplex_size = this->m_s.back() - this->m_s.front();
249 if (simplex_size.norm() < this->m_tolerance) {
259 x_sol = this->m_s.at(0);
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
Optimizer()
Definition Optimizer.hh:80
bool m_verbose
Definition Solver.hh:93
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
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 ¬es="-")
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