quala 0.0.1a1
Quasi-Newton and other accelerators
quala.py.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
30PYBIND11_MODULE(QUALA_MODULE_NAME, m) {
31 using py::operator""_a;
32
33 py::options options;
34 options.enable_function_signatures();
35 options.enable_user_defined_docstrings();
36
37 m.doc() = "Quala Quasi-Newton algorithms";
38
39#ifdef QUALA_VERSION_INFO
40 m.attr("__version__") = QUALA_VERSION_INFO;
41#else
42 m.attr("__version__") = "dev";
43#endif
44
45 py::class_<quala::LBFGSParams::CBFGSParams>(
46 m, "LBFGSParamsCBFGS", "C++ documentation: :cpp:member:`quala::LBFGSParams::CBFGSParams `")
47 .def(py::init())
48 .def(py::init(&kwargs_to_struct<quala::LBFGSParams::CBFGSParams>))
49 .def("to_dict", &struct_to_dict<quala::LBFGSParams::CBFGSParams>)
50 .def_readwrite("α", &quala::LBFGSParams::CBFGSParams::α)
51 .def_readwrite("ϵ", &quala::LBFGSParams::CBFGSParams::ϵ)
52 .def("__bool__", &quala::LBFGSParams::CBFGSParams::operator bool);
53
54 py::class_<quala::LBFGSParams>(m, "LBFGSParams",
55 "C++ documentation: :cpp:class:`quala::LBFGSParams`")
56 .def(py::init())
57 .def(py::init(&kwargs_to_struct<quala::LBFGSParams>))
58 .def("to_dict", &struct_to_dict<quala::LBFGSParams>)
59 .def_readwrite("memory", &quala::LBFGSParams::memory)
60 .def_readwrite("min_div_fac", &quala::LBFGSParams::min_div_fac)
61 .def_readwrite("min_abs_s", &quala::LBFGSParams::min_abs_s)
62 .def_readwrite("force_pos_def", &quala::LBFGSParams::force_pos_def)
63 .def_readwrite("cbfgs", &quala::LBFGSParams::cbfgs);
64
65 auto lbfgs =
66 py::class_<quala::LBFGS>(m, "LBFGS", "C++ documentation: :cpp:class:`quala::LBFGS`");
67 auto lbfgssign = py::enum_<quala::LBFGS::Sign>(
68 lbfgs, "Sign", "C++ documentation :cpp:enum:`quala::LBFGS::Sign`");
69 lbfgssign //
70 .value("Positive", quala::LBFGS::Sign::Positive)
71 .value("Negative", quala::LBFGS::Sign::Negative)
72 .export_values();
73 lbfgs //
74 .def(py::init<quala::LBFGS::Params>(), "params"_a)
75 .def(py::init([](py::dict params) -> quala::LBFGS {
76 return {kwargs_to_struct<quala::LBFGS::Params>(params)};
77 }),
78 "params"_a)
79 .def(py::init<quala::LBFGS::Params, quala::length_t>(), "params"_a, "n"_a)
80 .def(py::init([](py::dict params, quala::length_t n) -> quala::LBFGS {
81 return {kwargs_to_struct<quala::LBFGS::Params>(params), n};
82 }),
83 "params"_a, "n"_a)
84 .def_static("update_valid", quala::LBFGS::update_valid, "params"_a, "yᵀs"_a, "sᵀs"_a,
85 "pᵀp"_a)
86 .def(
87 "update",
89 quala::crvec pkp1, quala::LBFGS::Sign sign, bool forced) {
90 if (xk.size() != self.n())
91 throw std::invalid_argument("xk dimension mismatch");
92 if (xkp1.size() != self.n())
93 throw std::invalid_argument("xkp1 dimension mismatch");
94 if (pk.size() != self.n())
95 throw std::invalid_argument("pk dimension mismatch");
96 if (pkp1.size() != self.n())
97 throw std::invalid_argument("pkp1 dimension mismatch");
98 return self.update(xk, xkp1, pk, pkp1, sign, forced);
99 },
100 "xk"_a, "xkp1"_a, "pk"_a, "pkp1"_a, "sign"_a = quala::LBFGS::Sign::Positive,
101 "forced"_a = false)
102 .def(
103 "update_sy",
104 [](quala::LBFGS &self, quala::crvec sk, quala::crvec yk, quala::real_t pkp1Tpkp1,
105 bool forced) {
106 if (sk.size() != self.n())
107 throw std::invalid_argument("sk dimension mismatch");
108 if (yk.size() != self.n())
109 throw std::invalid_argument("yk dimension mismatch");
110 return self.update_sy(sk, yk, pkp1Tpkp1, forced);
111 },
112 "sk"_a, "yk"_a, "pkp1Tpkp1"_a, "forced"_a = false)
113 .def(
114 "apply",
115 [](quala::LBFGS &self, quala::rvec q, quala::real_t γ) {
116 if (q.size() != self.n())
117 throw std::invalid_argument("q dimension mismatch");
118 return self.apply(q, γ);
119 },
120 "q"_a, "γ"_a)
121 .def(
122 "apply",
123 [](quala::LBFGS &self, quala::rvec q, quala::real_t γ,
124 const std::vector<quala::vec::Index> &J) { return self.apply(q, γ, J); },
125 "q"_a, "γ"_a, "J"_a)
126 .def("reset", &quala::LBFGS::reset)
127 .def("current_history", &quala::LBFGS::current_history)
128 .def("resize", &quala::LBFGS::resize, "n"_a)
129 .def("scale_y", &quala::LBFGS::scale_y, "factor"_a)
130 .def_property_readonly("n", &quala::LBFGS::n)
131 .def("s", [](quala::LBFGS &self, quala::index_t i) -> quala::rvec { return self.s(i); })
132 .def("y", [](quala::LBFGS &self, quala::index_t i) -> quala::rvec { return self.y(i); })
133 .def("ρ", [](quala::LBFGS &self, quala::index_t i) -> quala::real_t & { return self.ρ(i); })
134 .def("α", [](quala::LBFGS &self, quala::index_t i) -> quala::real_t & { return self.α(i); })
135 .def_property_readonly("params", &quala::LBFGS::get_params);
136
137 py::class_<quala::AndersonAccelParams>(
138 m, "AndersonAccelParams", "C++ documentation: :cpp:class:`quala::AndersonAccelParams`")
139 .def(py::init())
140 .def(py::init(&kwargs_to_struct<quala::AndersonAccelParams>))
141 .def("to_dict", &struct_to_dict<quala::AndersonAccelParams>)
142 .def_readwrite("memory", &quala::AndersonAccelParams::memory);
143
144 py::class_<quala::AndersonAccel>(m, "AndersonAccel",
145 "C++ documentation: :cpp:class:`quala::AndersonAccel`")
146 .def(py::init<quala::AndersonAccel::Params>(), "params"_a)
147 .def(py::init([](py::dict params) -> quala::AndersonAccel {
148 return {kwargs_to_struct<quala::AndersonAccel::Params>(params)};
149 }),
150 "params"_a)
151 .def(py::init<quala::AndersonAccel::Params, quala::length_t>(), "params"_a, "n"_a)
152 .def(py::init([](py::dict params, quala::length_t n) -> quala::AndersonAccel {
153 return {kwargs_to_struct<quala::AndersonAccel::Params>(params), n};
154 }),
155 "params"_a, "n"_a)
156 .def("resize", &quala::AndersonAccel::resize, "n"_a)
157 .def(
158 "initialize",
159 [](quala::AndersonAccel &self, quala::crvec g_0, quala::vec r_0) {
160 if (g_0.size() != self.n())
161 throw std::invalid_argument("g_0 dimension mismatch");
162 if (r_0.size() != self.n())
163 throw std::invalid_argument("r_0 dimension mismatch");
164 self.initialize(g_0, std::move(r_0));
165 },
166 "g_0"_a, "r_0"_a)
167 .def(
168 "compute",
169 [](quala::AndersonAccel &self, quala::crvec g_k, quala::vec r_k, quala::rvec x_k_aa) {
170 if (g_k.size() != self.n())
171 throw std::invalid_argument("g_k dimension mismatch");
172 if (r_k.size() != self.n())
173 throw std::invalid_argument("r_k dimension mismatch");
174 if (x_k_aa.size() != self.n())
175 throw std::invalid_argument("x_k_aa dimension mismatch");
176 self.compute(g_k, std::move(r_k), x_k_aa);
177 },
178 "g_k"_a, "r_k"_a, "x_k_aa"_a)
179 .def(
180 "compute",
181 [](quala::AndersonAccel &self, quala::crvec g_k, quala::vec r_k) {
182 if (g_k.size() != self.n())
183 throw std::invalid_argument("g_k dimension mismatch");
184 if (r_k.size() != self.n())
185 throw std::invalid_argument("r_k dimension mismatch");
186 quala::vec x_k_aa(self.n());
187 self.compute(g_k, std::move(r_k), x_k_aa);
188 return x_k_aa;
189 },
190 "g_k"_a, "r_k"_a)
191 .def("reset", &quala::AndersonAccel::reset)
192 .def("current_history", &quala::AndersonAccel::current_history)
193 .def_property_readonly("params", &quala::AndersonAccel::get_params);
194
195 py::class_<quala::BroydenGoodParams>(m, "BroydenGoodParams",
196 "C++ documentation: :cpp:class:`quala::BroydenGoodParams`")
197 .def(py::init())
198 .def(py::init(&kwargs_to_struct<quala::BroydenGoodParams>))
199 .def("to_dict", &struct_to_dict<quala::BroydenGoodParams>)
200 .def_readwrite("memory", &quala::BroydenGoodParams::memory)
201 .def_readwrite("min_div_abs", &quala::BroydenGoodParams::min_div_abs)
202 .def_readwrite("force_pos_def", &quala::BroydenGoodParams::force_pos_def)
203 .def_readwrite("restarted", &quala::BroydenGoodParams::restarted)
204 .def_readwrite("powell_damping_factor", &quala::BroydenGoodParams::powell_damping_factor)
205 .def_readwrite("min_stepsize", &quala::BroydenGoodParams::min_stepsize);
206
207 py::class_<quala::BroydenGood>(m, "BroydenGood",
208 "C++ documentation: :cpp:class:`quala::BroydenGood`")
209 .def(py::init<quala::BroydenGood::Params>(), "params"_a)
210 .def(py::init([](py::dict params) -> quala::BroydenGood {
211 return {kwargs_to_struct<quala::BroydenGood::Params>(params)};
212 }),
213 "params"_a)
214 .def(py::init<quala::BroydenGood::Params, quala::length_t>(), "params"_a, "n"_a)
215 .def(py::init([](py::dict params, quala::length_t n) -> quala::BroydenGood {
216 return {kwargs_to_struct<quala::BroydenGood::Params>(params), n};
217 }),
218 "params"_a, "n"_a)
219 .def("resize", &quala::BroydenGood::resize, "n"_a)
220 .def(
221 "update",
223 quala::crvec pkp1, bool forced) {
224 if (xk.size() != self.n())
225 throw std::invalid_argument("xk dimension mismatch");
226 if (xkp1.size() != self.n())
227 throw std::invalid_argument("xkp1 dimension mismatch");
228 if (pk.size() != self.n())
229 throw std::invalid_argument("pk dimension mismatch");
230 if (pkp1.size() != self.n())
231 throw std::invalid_argument("pkp1 dimension mismatch");
232 return self.update(xk, xkp1, pk, pkp1, forced);
233 },
234 "xk"_a, "xkp1"_a, "pk"_a, "pkp1"_a, "forced"_a = false)
235 .def(
236 "update_sy",
237 [](quala::BroydenGood &self, quala::crvec sk, quala::crvec yk, bool forced) {
238 if (sk.size() != self.n())
239 throw std::invalid_argument("sk dimension mismatch");
240 if (yk.size() != self.n())
241 throw std::invalid_argument("yk dimension mismatch");
242 return self.update_sy(sk, yk, forced);
243 },
244 "sk"_a, "yk"_a, "forced"_a = false)
245 .def(
246 "apply",
248 if (q.size() != self.n())
249 throw std::invalid_argument("q dimension mismatch");
250 return self.apply(q, γ);
251 },
252 "q"_a, "γ"_a = -1)
253 .def("reset", &quala::BroydenGood::reset)
254 .def("current_history", &quala::BroydenGood::current_history)
255 .def_property_readonly("params", &quala::BroydenGood::get_params);
256}
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
auto s(index_t i)
Definition: decl/lbfgs.hpp:136
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
bool apply(rvec q, real_t γ=-1)
Apply the inverse Hessian approximation to the given vector q.
Definition: lbfgs.hpp:64
void resize(length_t n)
Re-allocate storage for a problem with a different size.
Definition: lbfgs.hpp:185
real_t & α(index_t i)
Definition: decl/lbfgs.hpp:142
real_t & ρ(index_t i)
Definition: decl/lbfgs.hpp:140
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
bool update(crvec xₖ, crvec xₙₑₓₜ, crvec pₖ, crvec pₙₑₓₜ, Sign sign=Sign::Positive, bool forced=false)
Update the inverse Hessian approximation using the new vectors xₖ₊₁ and pₖ₊₁.
Definition: lbfgs.hpp:56
auto y(index_t i)
Definition: decl/lbfgs.hpp:138
bool update_sy(const anymat< VecS > &s, const anymat< VecY > &y, real_t pₙₑₓₜᵀpₙₑₓₜ, bool forced=false)
Update the inverse Hessian approximation using the new vectors sₖ = xₖ₊₁ - xₖ and yₖ = pₖ₊₁ - pₖ.
Definition: lbfgs.hpp:34
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
real_t min_stepsize
Minimum automatic step size.
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
PYBIND11_MODULE(QUALA_MODULE_NAME, m)
Definition: quala.py.cpp:30
real_t ϵ
Set to zero to disable CBFGS check.
Definition: decl/lbfgs.hpp:20