openGPMP
Open Source Mathematics Package
quasi.hpp
Go to the documentation of this file.
1 /*************************************************************************
2  *
3  * Project
4  * _____ _____ __ __ _____
5  * / ____| __ \| \/ | __ \
6  * ___ _ __ ___ _ __ | | __| |__) | \ / | |__) |
7  * / _ \| '_ \ / _ \ '_ \| | |_ | ___/| |\/| | ___/
8  *| (_) | |_) | __/ | | | |__| | | | | | | |
9  * \___/| .__/ \___|_| |_|\_____|_| |_| |_|_|
10  * | |
11  * |_|
12  *
13  * Copyright (C) Akiel Aries, <akiel@akiel.org>, et al.
14  *
15  * This software is licensed as described in the file LICENSE, which
16  * you should have received as part of this distribution. The terms
17  * among other details are referenced in the official documentation
18  * seen here : https://akielaries.github.io/openGPMP/ along with
19  * important files seen in this project.
20  *
21  * You may opt to use, copy, modify, merge, publish, distribute
22  * and/or sell copies of the Software, and permit persons to whom
23  * the Software is furnished to do so, under the terms of the
24  * LICENSE file. As this is an Open Source effort, all implementations
25  * must be of the same methodology.
26  *
27  *
28  *
29  * This software is distributed on an AS IS basis, WITHOUT
30  * WARRANTY OF ANY KIND, either express or implied.
31  *
32  ************************************************************************/
33 
34 #ifndef QUASI_HPP
35 #define QUASI_HPP
36 
37 #include <functional>
38 #include <vector>
39 
40 namespace gpmp {
41 
42 namespace optim {
43 
48 class QuasiNewton {
49  public:
60  std::vector<double> bhhh_optimize(
61  const std::function<double(const std::vector<double> &)> &func,
62  const std::vector<double> &initial_point,
63  double tolerance,
64  size_t max_iterations);
65 
74  std::vector<double> calculate_gradient(
75  const std::function<double(const std::vector<double> &)> &func,
76  const std::vector<double> &point,
77  double epsilon);
78 
85  std::vector<std::vector<double>>
86  calculate_bhhh_matrix(const std::vector<double> &gradient);
87 
96  std::vector<double>
97  update_point(const std::vector<double> &current_point,
98  const std::vector<double> &gradient,
99  const std::vector<std::vector<double>> &bhhh_matrix);
100 
111  std::vector<double> bfgs_optimize(
112  const std::function<double(const std::vector<double> &)> &func,
113  const std::vector<double> &initial_point,
114  double tolerance,
115  size_t max_iterations);
116 
124  std::vector<double> calculate_search_direction(
125  const std::vector<double> &gradient,
126  const std::vector<std::vector<double>> &hessian_inverse);
127 
136  double
137  line_search(const std::function<double(const std::vector<double> &)> &func,
138  const std::vector<double> &current_point,
139  const std::vector<double> &search_direction);
140 
149  std::vector<double>
150  update_point(const std::vector<double> &current_point,
151  const std::vector<double> &search_direction,
152  double step_size);
153 
162  std::vector<double>
163  calculate_gradient_difference(const std::vector<double> &next_point,
164  const std::vector<double> &current_point,
165  const std::vector<double> &gradient);
166 
175  std::vector<std::vector<double>> update_hessian_inverse(
176  const std::vector<std::vector<double>> &hessian_inverse,
177  const std::vector<double> &gradient_difference,
178  const std::vector<double> &search_direction);
179 
187  double dot_product(const std::vector<double> &a,
188  const std::vector<double> &b);
189 
197  std::vector<double> vector_subtraction(const std::vector<double> &a,
198  const std::vector<double> &b) const;
199 
212  std::tuple<std::vector<double>, double>
213  lbfgs_optimize(const std::function<double(const std::vector<double> &)> &f,
214  const std::vector<double> &initial_point,
215  double tolerance = 1e-4,
216  size_t max_iterations = 100,
217  size_t memory_size = 5);
218 };
219 
220 } // namespace optim
221 
222 } // namespace gpmp
223 #endif
Class implementing Quasi-Newton optimization methods.
Definition: quasi.hpp:48
double line_search(const std::function< double(const std::vector< double > &)> &func, const std::vector< double > &current_point, const std::vector< double > &search_direction)
Perform line search to find an appropriate step size.
Definition: quasi.cpp:221
std::vector< double > bhhh_optimize(const std::function< double(const std::vector< double > &)> &func, const std::vector< double > &initial_point, double tolerance, size_t max_iterations)
Optimize a function using the Berndt–Hall–Hall–Hausman (BHHH) algorithm.
Definition: quasi.cpp:58
std::vector< double > calculate_search_direction(const std::vector< double > &gradient, const std::vector< std::vector< double >> &hessian_inverse)
Calculate the search direction using the BFGS method.
Definition: quasi.cpp:205
std::vector< std::vector< double > > calculate_bhhh_matrix(const std::vector< double > &gradient)
Calculate the BHHH matrix from the gradient.
Definition: quasi.cpp:117
std::vector< double > calculate_gradient_difference(const std::vector< double > &next_point, const std::vector< double > &current_point, const std::vector< double > &gradient)
Calculate the gradient difference between two points.
Definition: quasi.cpp:271
double dot_product(const std::vector< double > &a, const std::vector< double > &b)
Calculate the dot product of two vectors.
Definition: quasi.cpp:311
std::vector< double > calculate_gradient(const std::function< double(const std::vector< double > &)> &func, const std::vector< double > &point, double epsilon)
Calculate the gradient of a function at a given point.
Definition: quasi.cpp:96
std::vector< double > update_point(const std::vector< double > &current_point, const std::vector< double > &gradient, const std::vector< std::vector< double >> &bhhh_matrix)
Update the current point using the BHHH matrix.
Definition: quasi.cpp:131
std::tuple< std::vector< double >, double > lbfgs_optimize(const std::function< double(const std::vector< double > &)> &f, const std::vector< double > &initial_point, double tolerance=1e-4, size_t max_iterations=100, size_t memory_size=5)
Optimize a function using the Limited-Memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS) algorithm.
Definition: quasi.cpp:327
std::vector< double > bfgs_optimize(const std::function< double(const std::vector< double > &)> &func, const std::vector< double > &initial_point, double tolerance, size_t max_iterations)
Optimize a function using the Broyden–Fletcher–Goldfarb–Shanno (BFGS) algorithm.
Definition: quasi.cpp:145
std::vector< std::vector< double > > update_hessian_inverse(const std::vector< std::vector< double >> &hessian_inverse, const std::vector< double > &gradient_difference, const std::vector< double > &search_direction)
Update the inverse of the Hessian matrix using the BFGS method.
Definition: quasi.cpp:287
std::vector< double > vector_subtraction(const std::vector< double > &a, const std::vector< double > &b) const
Subtract two vectors element-wise.
Definition: quasi.cpp:40
The source C++ openGPMP namespace.