quala main
Quasi-Newton and other accelerators
quala.cpp
Go to the documentation of this file.
1/**
2 * @file
3 * This file defines all Python bindings.
4 */
5
6#include <pybind11/attr.h>
7#include <pybind11/cast.h>
8#include <pybind11/chrono.h>
9#include <pybind11/detail/common.h>
10#include <pybind11/eigen.h>
11#include <pybind11/functional.h>
12#include <pybind11/iostream.h>
13#include <pybind11/pybind11.h>
14#include <pybind11/pytypes.h>
15#include <pybind11/stl.h>
16
17#include <memory>
18#include <optional>
19#include <sstream>
20#include <stdexcept>
21#include <string>
22#include <tuple>
23#include <type_traits>
24#include <utility>
25
26#include "kwargs-to-struct.hpp"
27
28namespace py = pybind11;
29
30#define STRINGIFY(x) #x
31#define MACRO_STRINGIFY(x) STRINGIFY(x)
32
33PYBIND11_MODULE(QUALA_MODULE_NAME, m) {
34 using py::operator""_a;
35
36 py::options options;
37 options.enable_function_signatures();
38 options.enable_user_defined_docstrings();
39
40 m.doc() = "Quala Quasi-Newton algorithms";
41
42#ifdef VERSION_INFO
43 m.attr("__version__") = MACRO_STRINGIFY(VERSION_INFO);
44#else
45 m.attr("__version__") = "dev";
46#endif
47
48 py::class_<quala::LBFGSParams::CBFGSParams>(
49 m, "LBFGSParamsCBFGS", "C++ documentation: :cpp:member:`quala::LBFGSParams::CBFGSParams `")
50 .def(py::init())
51 .def(py::init(&kwargs_to_struct<quala::LBFGSParams::CBFGSParams>))
52 .def("to_dict", &struct_to_dict<quala::LBFGSParams::CBFGSParams>)
53 .def_readwrite("α", &quala::LBFGSParams::CBFGSParams::α)
54 .def_readwrite("ϵ", &quala::LBFGSParams::CBFGSParams::ϵ)
55 .def("__bool__", &quala::LBFGSParams::CBFGSParams::operator bool);
56
57 py::class_<quala::LBFGSParams>(m, "LBFGSParams",
58 "C++ documentation: :cpp:class:`quala::LBFGSParams`")
59 .def(py::init())
60 .def(py::init(&kwargs_to_struct<quala::LBFGSParams>))
61 .def("to_dict", &struct_to_dict<quala::LBFGSParams>)
62 .def_readwrite("memory", &quala::LBFGSParams::memory)
63 .def_readwrite("min_div_fac", &quala::LBFGSParams::min_div_fac)
64 .def_readwrite("min_abs_s", &quala::LBFGSParams::min_abs_s)
65 .def_readwrite("force_pos_def", &quala::LBFGSParams::force_pos_def)
66 .def_readwrite("cbfgs", &quala::LBFGSParams::cbfgs);
67
68 auto lbfgs =
69 py::class_<quala::LBFGS>(m, "LBFGS", "C++ documentation: :cpp:class:`quala::LBFGS`");
70 auto lbfgssign = py::enum_<quala::LBFGS::Sign>(
71 lbfgs, "Sign", "C++ documentation :cpp:enum:`quala::LBFGS::Sign`");
72 lbfgssign //
73 .value("Positive", quala::LBFGS::Sign::Positive)
74 .value("Negative", quala::LBFGS::Sign::Negative)
75 .export_values();
76 lbfgs //
77 .def(py::init<quala::LBFGS::Params>(), "params"_a)
78 .def(py::init([](py::dict params) -> quala::LBFGS {
79 return {kwargs_to_struct<quala::LBFGS::Params>(params)};
80 }),
81 "params"_a)
82 .def(py::init<quala::LBFGS::Params, quala::length_t>(), "params"_a, "n"_a)
83 .def(py::init([](py::dict params, quala::length_t n) -> quala::LBFGS {
84 return {kwargs_to_struct<quala::LBFGS::Params>(params), n};
85 }),
86 "params"_a, "n"_a)
87 .def_static("update_valid", quala::LBFGS::update_valid, "params"_a, "yᵀs"_a, "sᵀs"_a,
88 "pᵀp"_a)
89 .def(
90 "update",
92 quala::crvec pkp1, quala::LBFGS::Sign sign, bool forced) {
93 if (xk.size() != self.n())
94 throw std::invalid_argument("xk dimension mismatch");
95 if (xkp1.size() != self.n())
96 throw std::invalid_argument("xkp1 dimension mismatch");
97 if (pk.size() != self.n())
98 throw std::invalid_argument("pk dimension mismatch");
99 if (pkp1.size() != self.n())
100 throw std::invalid_argument("pkp1 dimension mismatch");
101 return self.update(xk, xkp1, pk, pkp1, sign, forced);
102 },
103 "xk"_a, "xkp1"_a, "pk"_a, "pkp1"_a, "sign"_a = quala::LBFGS::Sign::Positive,
104 "forced"_a = false)
105 .def(
106 "update_sy",
107 [](quala::LBFGS &self, quala::crvec sk, quala::crvec yk, quala::real_t pkp1Tpkp1,
108 bool forced) {
109 if (sk.size() != self.n())
110 throw std::invalid_argument("sk dimension mismatch");
111 if (yk.size() != self.n())
112 throw std::invalid_argument("yk dimension mismatch");
113 return self.update_sy(sk, yk, pkp1Tpkp1, forced);
114 },
115 "sk"_a, "yk"_a, "pkp1Tpkp1"_a, "forced"_a = false)
116 .def(
117 "apply",
118 [](quala::LBFGS &self, quala::rvec q, quala::real_t γ) {
119 if (q.size() != self.n())
120 throw std::invalid_argument("q dimension mismatch");
121 return self.apply(q, γ);
122 },
123 "q"_a, "γ"_a)
124 .def(
125 "apply",
126 [](quala::LBFGS &self, quala::rvec q, quala::real_t γ,
127 const std::vector<quala::vec::Index> &J) { return self.apply(q, γ, J); },
128 "q"_a, "γ"_a, "J"_a)
129 .def("reset", &quala::LBFGS::reset)
130 .def("current_history", &quala::LBFGS::current_history)
131 .def("resize", &quala::LBFGS::resize, "n"_a)
132 .def("scale_y", &quala::LBFGS::scale_y, "factor"_a)
133 .def_property_readonly("n", &quala::LBFGS::n)
134 .def("s", [](quala::LBFGS &self, quala::index_t i) -> quala::rvec { return self.s(i); })
135 .def("y", [](quala::LBFGS &self, quala::index_t i) -> quala::rvec { return self.y(i); })
136 .def("ρ", [](quala::LBFGS &self, quala::index_t i) -> quala::real_t & { return self.ρ(i); })
137 .def("α", [](quala::LBFGS &self, quala::index_t i) -> quala::real_t & { return self.α(i); })
138 .def_property_readonly("params", &quala::LBFGS::get_params);
139
140 py::class_<quala::AndersonAccelParams>(
141 m, "AndersonAccelParams", "C++ documentation: :cpp:class:`quala::AndersonAccelParams`")
142 .def(py::init())
143 .def(py::init(&kwargs_to_struct<quala::AndersonAccelParams>))
144 .def("to_dict", &struct_to_dict<quala::AndersonAccelParams>)
145 .def_readwrite("memory", &quala::AndersonAccelParams::memory);
146
147 py::class_<quala::AndersonAccel>(m, "AndersonAccel",
148 "C++ documentation: :cpp:class:`quala::AndersonAccel`")
149 .def(py::init<quala::AndersonAccel::Params>(), "params"_a)
150 .def(py::init([](py::dict params) -> quala::AndersonAccel {
151 return {kwargs_to_struct<quala::AndersonAccel::Params>(params)};
152 }),
153 "params"_a)
154 .def(py::init<quala::AndersonAccel::Params, quala::length_t>(), "params"_a, "n"_a)
155 .def(py::init([](py::dict params, quala::length_t n) -> quala::AndersonAccel {
156 return {kwargs_to_struct<quala::AndersonAccel::Params>(params), n};
157 }),
158 "params"_a, "n"_a)
159 .def("resize", &quala::AndersonAccel::resize, "n"_a)
160 .def(
161 "initialize",
162 [](quala::AndersonAccel &self, quala::crvec g_0, quala::vec r_0) {
163 if (g_0.size() != self.n())
164 throw std::invalid_argument("g_0 dimension mismatch");
165 if (r_0.size() != self.n())
166 throw std::invalid_argument("r_0 dimension mismatch");
167 self.initialize(g_0, std::move(r_0));
168 },
169 "g_0"_a, "r_0"_a)
170 .def(
171 "compute",
172 [](quala::AndersonAccel &self, quala::crvec g_k, quala::vec r_k, quala::rvec x_k_aa) {
173 if (g_k.size() != self.n())
174 throw std::invalid_argument("g_k dimension mismatch");
175 if (r_k.size() != self.n())
176 throw std::invalid_argument("r_k dimension mismatch");
177 if (x_k_aa.size() != self.n())
178 throw std::invalid_argument("x_k_aa dimension mismatch");
179 self.compute(g_k, std::move(r_k), x_k_aa);
180 },
181 "g_k"_a, "r_k"_a, "x_k_aa"_a)
182 .def(
183 "compute",
184 [](quala::AndersonAccel &self, quala::crvec g_k, quala::vec r_k) {
185 if (g_k.size() != self.n())
186 throw std::invalid_argument("g_k dimension mismatch");
187 if (r_k.size() != self.n())
188 throw std::invalid_argument("r_k dimension mismatch");
189 quala::vec x_k_aa(self.n());
190 self.compute(g_k, std::move(r_k), x_k_aa);
191 return x_k_aa;
192 },
193 "g_k"_a, "r_k"_a)
194 .def("reset", &quala::AndersonAccel::reset)
195 .def("current_history", &quala::AndersonAccel::current_history)
196 .def_property_readonly("params", &quala::AndersonAccel::get_params);
197
198 py::class_<quala::BroydenGoodParams>(m, "BroydenGoodParams",
199 "C++ documentation: :cpp:class:`quala::BroydenGoodParams`")
200 .def(py::init())
201 .def(py::init(&kwargs_to_struct<quala::BroydenGoodParams>))
202 .def("to_dict", &struct_to_dict<quala::BroydenGoodParams>)
203 .def_readwrite("memory", &quala::BroydenGoodParams::memory)
204 .def_readwrite("min_div_abs", &quala::BroydenGoodParams::min_div_abs)
205 .def_readwrite("force_pos_def", &quala::BroydenGoodParams::force_pos_def)
206 .def_readwrite("restarted", &quala::BroydenGoodParams::restarted)
207 .def_readwrite("powell_damping_factor", &quala::BroydenGoodParams::powell_damping_factor);
208
209 py::class_<quala::BroydenGood>(m, "BroydenGood",
210 "C++ documentation: :cpp:class:`quala::BroydenGood`")
211 .def(py::init<quala::BroydenGood::Params>(), "params"_a)
212 .def(py::init([](py::dict params) -> quala::BroydenGood {
213 return {kwargs_to_struct<quala::BroydenGood::Params>(params)};
214 }),
215 "params"_a)
216 .def(py::init<quala::BroydenGood::Params, quala::length_t>(), "params"_a, "n"_a)
217 .def(py::init([](py::dict params, quala::length_t n) -> quala::BroydenGood {
218 return {kwargs_to_struct<quala::BroydenGood::Params>(params), n};
219 }),
220 "params"_a, "n"_a)
221 .def("resize", &quala::BroydenGood::resize, "n"_a)
222 .def(
223 "update",
225 quala::crvec pkp1, bool forced) {
226 if (xk.size() != self.n())
227 throw std::invalid_argument("xk dimension mismatch");
228 if (xkp1.size() != self.n())
229 throw std::invalid_argument("xkp1 dimension mismatch");
230 if (pk.size() != self.n())
231 throw std::invalid_argument("pk dimension mismatch");
232 if (pkp1.size() != self.n())
233 throw std::invalid_argument("pkp1 dimension mismatch");
234 return self.update(xk, xkp1, pk, pkp1, forced);
235 },
236 "xk"_a, "xkp1"_a, "pk"_a, "pkp1"_a, "forced"_a = false)
237 .def(
238 "update_sy",
239 [](quala::BroydenGood &self, quala::crvec sk, quala::crvec yk, bool forced) {
240 if (sk.size() != self.n())
241 throw std::invalid_argument("sk dimension mismatch");
242 if (yk.size() != self.n())
243 throw std::invalid_argument("yk dimension mismatch");
244 return self.update_sy(sk, yk, forced);
245 },
246 "sk"_a, "yk"_a, "forced"_a = false)
247 .def(
248 "apply",
250 if (q.size() != self.n())
251 throw std::invalid_argument("q dimension mismatch");
252 return self.apply(q, γ);
253 },
254 "q"_a, "γ"_a = -1)
255 .def("reset", &quala::BroydenGood::reset)
256 .def("current_history", &quala::BroydenGood::current_history)
257 .def_property_readonly("params", &quala::BroydenGood::get_params);
258}
Anderson Acceleration.
length_t current_history() const
Get the number of columns currently stored in the buffer.
const Params & get_params() const
Get the parameters.
void resize(length_t n)
Change the problem dimension.
void reset()
Reset the accelerator (but keep the last function value and residual, so calling initialize is not ne...
Broyden's “Good” method for solving systems of nonlinear equations .
length_t current_history() const
Get the number of previous s and y vectors currently stored in the buffer.
const Params & get_params() const
Get the parameters.
void resize(length_t n)
Re-allocate storage for a problem with a different size.
void reset()
Throw away the approximation and all previous vectors s and y.
Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm.
Definition: decl/lbfgs.hpp:72
static bool update_valid(const LBFGSParams &params, real_t yᵀs, real_t sᵀs, real_t pᵀp)
Check if the new vectors s and y allow for a valid BFGS update that preserves the positive definitene...
Definition: lbfgs.hpp:9
length_t current_history() const
Get the number of previous s and y vectors currently stored in the buffer.
Definition: decl/lbfgs.hpp:134
length_t n() const
Get the size of the s and y vectors in the buffer.
Definition: decl/lbfgs.hpp:125
const Params & get_params() const
Get the parameters.
Definition: decl/lbfgs.hpp:122
void resize(length_t n)
Re-allocate storage for a problem with a different size.
Definition: lbfgs.hpp:185
void reset()
Throw away the approximation and all previous vectors s and y.
Definition: lbfgs.hpp:180
Sign
The sign of the vectors passed to the LBFGS::update method.
Definition: decl/lbfgs.hpp:78
void scale_y(real_t factor)
Scale the stored y vectors by the given factor.
Definition: lbfgs.hpp:196
This file defines mappings from Python dicts (kwargs) to simple parameter structs.
real_t powell_damping_factor
Powell's trick, damping, prevents nonsingularity.
Eigen::Ref< const vec > crvec
Default type for immutable references to vectors.
Definition: vec.hpp:18
real_t min_abs_s
Reject update if .
Definition: decl/lbfgs.hpp:15
real_t min_div_abs
Reject update if .
CBFGSParams cbfgs
Parameters in the cautious BFGS update condition.
Definition: decl/lbfgs.hpp:26
bool restarted
If set to true, the buffer is cleared after memory iterations.
length_t memory
Length of the history to keep.
Definition: decl/lbfgs.hpp:11
index_t length_t
Default type for vector sizes.
Definition: vec.hpp:29
bool force_pos_def
If set to true, the inverse Hessian estimate should remain definite, i.e.
Definition: decl/lbfgs.hpp:33
realvec vec
Default type for vectors.
Definition: vec.hpp:14
Eigen::Index index_t
Default type for vector indices.
Definition: vec.hpp:27
double real_t
Default floating point type.
Definition: vec.hpp:8
real_t min_div_fac
Reject update if .
Definition: decl/lbfgs.hpp:13
Eigen::Ref< vec > rvec
Default type for mutable references to vectors.
Definition: vec.hpp:16
#define MACRO_STRINGIFY(x)
Definition: quala.cpp:31
PYBIND11_MODULE(QUALA_MODULE_NAME, m)
Definition: quala.cpp:33
real_t ϵ
Set to zero to disable CBFGS check.
Definition: decl/lbfgs.hpp:20