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>
28namespace py = pybind11;
30#define STRINGIFY(x) #x
31#define MACRO_STRINGIFY(x) STRINGIFY(x)
34 using py::operator
""_a;
37 options.enable_function_signatures();
38 options.enable_user_defined_docstrings();
40 m.doc() =
"Quala Quasi-Newton algorithms";
45 m.attr(
"__version__") =
"dev";
48 py::class_<quala::LBFGSParams::CBFGSParams>(
49 m,
"LBFGSParamsCBFGS",
50 "C++ documentation: :cpp:member:`quala::LBFGSParams::CBFGSParams `")
52 .def(py::init(&kwargs_to_struct<quala::LBFGSParams::CBFGSParams>))
53 .def(
"to_dict", &struct_to_dict<quala::LBFGSParams::CBFGSParams>)
56 .def(
"__bool__", &quala::LBFGSParams::CBFGSParams::operator
bool);
58 py::class_<quala::LBFGSParams>(
59 m,
"LBFGSParams",
"C++ documentation: :cpp:class:`quala::LBFGSParams`")
61 .def(py::init(&kwargs_to_struct<quala::LBFGSParams>))
62 .def(
"to_dict", &struct_to_dict<quala::LBFGSParams>)
69 auto lbfgs = py::class_<quala::LBFGS>(
70 m,
"LBFGS",
"C++ documentation: :cpp:class:`quala::LBFGS`");
71 auto lbfgssign = py::enum_<quala::LBFGS::Sign>(
72 lbfgs,
"Sign",
"C++ documentation :cpp:enum:`quala::LBFGS::Sign`");
78 .def(py::init<quala::LBFGS::Params>(),
"params"_a)
80 return {kwargs_to_struct<quala::LBFGS::Params>(params)};
83 .def(py::init<quala::LBFGS::Params, quala::length_t>(),
"params"_a,
86 return {kwargs_to_struct<quala::LBFGS::Params>(params),
n};
90 "yᵀs"_a,
"sᵀs"_a,
"pᵀp"_a)
96 if (xk.size() !=
self.n())
97 throw std::invalid_argument(
"xk dimension mismatch");
98 if (xkp1.size() !=
self.n())
99 throw std::invalid_argument(
"xkp1 dimension mismatch");
100 if (pk.size() !=
self.n())
101 throw std::invalid_argument(
"pk dimension mismatch");
102 if (pkp1.size() !=
self.n())
103 throw std::invalid_argument(
"pkp1 dimension mismatch");
104 return self.update(xk, xkp1, pk, pkp1, sign, forced);
106 "xk"_a,
"xkp1"_a,
"pk"_a,
"pkp1"_a,
112 if (sk.size() !=
self.n())
113 throw std::invalid_argument(
"sk dimension mismatch");
114 if (yk.size() !=
self.n())
115 throw std::invalid_argument(
"yk dimension mismatch");
116 return self.update_sy(sk, yk, pkp1Tpkp1, forced);
118 "sk"_a,
"yk"_a,
"pkp1Tpkp1"_a,
"forced"_a =
false)
122 if (q.size() !=
self.n())
123 throw std::invalid_argument(
"q dimension mismatch");
124 return self.apply(q, γ);
130 const std::vector<quala::vec::Index> &J) {
131 return self.apply(q, γ, J);
157 py::class_<quala::AndersonAccelParams>(
158 m,
"AndersonAccelParams",
159 "C++ documentation: :cpp:class:`quala::AndersonAccelParams`")
161 .def(py::init(&kwargs_to_struct<quala::AndersonAccelParams>))
162 .def(
"to_dict", &struct_to_dict<quala::AndersonAccelParams>)
165 py::class_<quala::AndersonAccel>(
167 "C++ documentation: :cpp:class:`quala::AndersonAccel`")
168 .def(py::init<quala::AndersonAccel::Params>(),
"params"_a)
171 kwargs_to_struct<quala::AndersonAccel::Params>(params)};
174 .def(py::init<quala::AndersonAccel::Params, quala::length_t>(),
176 .def(py::init([](py::dict params,
178 return {kwargs_to_struct<quala::AndersonAccel::Params>(params),
186 if (g_0.size() !=
self.n())
187 throw std::invalid_argument(
"g_0 dimension mismatch");
188 if (r_0.size() !=
self.n())
189 throw std::invalid_argument(
"r_0 dimension mismatch");
190 self.initialize(g_0, std::move(r_0));
197 if (g_k.size() !=
self.n())
198 throw std::invalid_argument(
"g_k dimension mismatch");
199 if (r_k.size() !=
self.n())
200 throw std::invalid_argument(
"r_k dimension mismatch");
201 if (x_k_aa.size() !=
self.n())
202 throw std::invalid_argument(
"x_k_aa dimension mismatch");
203 self.compute(g_k, std::move(r_k), x_k_aa);
205 "g_k"_a,
"r_k"_a,
"x_k_aa"_a)
209 if (g_k.size() !=
self.n())
210 throw std::invalid_argument(
"g_k dimension mismatch");
211 if (r_k.size() !=
self.n())
212 throw std::invalid_argument(
"r_k dimension mismatch");
214 self.compute(g_k, std::move(r_k), x_k_aa);
222 py::class_<quala::BroydenGoodParams>(
223 m,
"BroydenGoodParams",
224 "C++ documentation: :cpp:class:`quala::BroydenGoodParams`")
226 .def(py::init(&kwargs_to_struct<quala::BroydenGoodParams>))
227 .def(
"to_dict", &struct_to_dict<quala::BroydenGoodParams>)
230 .def_readwrite(
"force_pos_def",
234 py::class_<quala::BroydenGood>(
235 m,
"BroydenGood",
"C++ documentation: :cpp:class:`quala::BroydenGood`")
236 .def(py::init<quala::BroydenGood::Params>(),
"params"_a)
238 return {kwargs_to_struct<quala::BroydenGood::Params>(params)};
241 .def(py::init<quala::BroydenGood::Params, quala::length_t>(),
243 .def(py::init([](py::dict params,
245 return {kwargs_to_struct<quala::BroydenGood::Params>(params),
254 if (xk.size() !=
self.n())
255 throw std::invalid_argument(
"xk dimension mismatch");
256 if (xkp1.size() !=
self.n())
257 throw std::invalid_argument(
"xkp1 dimension mismatch");
258 if (pk.size() !=
self.n())
259 throw std::invalid_argument(
"pk dimension mismatch");
260 if (pkp1.size() !=
self.n())
261 throw std::invalid_argument(
"pkp1 dimension mismatch");
262 return self.update(xk, xkp1, pk, pkp1, forced);
264 "xk"_a,
"xkp1"_a,
"pk"_a,
"pkp1"_a,
"forced"_a =
false)
269 if (sk.size() !=
self.n())
270 throw std::invalid_argument(
"sk dimension mismatch");
271 if (yk.size() !=
self.n())
272 throw std::invalid_argument(
"yk dimension mismatch");
273 return self.update_sy(sk, yk, forced);
275 "sk"_a,
"yk"_a,
"forced"_a =
false)
279 if (q.size() !=
self.n())
280 throw std::invalid_argument(
"q dimension mismatch");
281 return self.apply(q);
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.
static bool update_valid(const LBFGSParams ¶ms, 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...
length_t current_history() const
Get the number of previous s and y vectors currently stored in the buffer.
length_t n() const
Get the size of the s and y vectors 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.
Sign
The sign of the vectors passed to the LBFGS::update method.
void scale_y(real_t factor)
Scale the stored y vectors by the given factor.
This file defines mappings from Python dicts (kwargs) to simple parameter structs.
Eigen::Ref< const vec > crvec
Default type for immutable references to vectors.
real_t min_abs_s
Reject update if .
real_t min_div_abs
Reject update if .
CBFGSParams cbfgs
Parameters in the cautious BFGS update condition.
bool restarted
If set to true, the buffer is cleared after memory iterations.
length_t memory
Length of the history to keep.
index_t length_t
Default type for vector sizes.
bool force_pos_def
If set to true, the inverse Hessian estimate should remain definite, i.e.
realvec vec
Default type for vectors.
Eigen::Index index_t
Default type for vector indices.
double real_t
Default floating point type.
real_t min_div_fac
Reject update if .
Eigen::Ref< vec > rvec
Default type for mutable references to vectors.
#define MACRO_STRINGIFY(x)
PYBIND11_MODULE(QUALA_MODULE_NAME, m)
real_t ϵ
Set to zero to disable CBFGS check.