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;
31 using py::operator
""_a;
34 options.enable_function_signatures();
35 options.enable_user_defined_docstrings();
37 m.doc() =
"Quala Quasi-Newton algorithms";
39#ifdef QUALA_VERSION_INFO
40 m.attr(
"__version__") = QUALA_VERSION_INFO;
42 m.attr(
"__version__") =
"dev";
45 py::class_<quala::LBFGSParams::CBFGSParams>(
46 m,
"LBFGSParamsCBFGS",
"C++ documentation: :cpp:member:`quala::LBFGSParams::CBFGSParams `")
48 .def(py::init(&kwargs_to_struct<quala::LBFGSParams::CBFGSParams>))
49 .def(
"to_dict", &struct_to_dict<quala::LBFGSParams::CBFGSParams>)
52 .def(
"__bool__", &quala::LBFGSParams::CBFGSParams::operator
bool);
54 py::class_<quala::LBFGSParams>(m,
"LBFGSParams",
55 "C++ documentation: :cpp:class:`quala::LBFGSParams`")
57 .def(py::init(&kwargs_to_struct<quala::LBFGSParams>))
58 .def(
"to_dict", &struct_to_dict<quala::LBFGSParams>)
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`");
74 .def(py::init<quala::LBFGS::Params>(),
"params"_a)
76 return {kwargs_to_struct<quala::LBFGS::Params>(params)};
79 .def(py::init<quala::LBFGS::Params, quala::length_t>(),
"params"_a,
"n"_a)
81 return {kwargs_to_struct<quala::LBFGS::Params>(params), n};
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);
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);
112 "sk"_a,
"yk"_a,
"pkp1Tpkp1"_a,
"forced"_a =
false)
116 if (q.size() != self.
n())
117 throw std::invalid_argument(
"q dimension mismatch");
118 return self.
apply(q, γ);
124 const std::vector<quala::vec::Index> &J) {
return self.
apply(q, γ, J); },
137 py::class_<quala::AndersonAccelParams>(
138 m,
"AndersonAccelParams",
"C++ documentation: :cpp:class:`quala::AndersonAccelParams`")
140 .def(py::init(&kwargs_to_struct<quala::AndersonAccelParams>))
141 .def(
"to_dict", &struct_to_dict<quala::AndersonAccelParams>)
144 py::class_<quala::AndersonAccel>(m,
"AndersonAccel",
145 "C++ documentation: :cpp:class:`quala::AndersonAccel`")
146 .def(py::init<quala::AndersonAccel::Params>(),
"params"_a)
148 return {kwargs_to_struct<quala::AndersonAccel::Params>(params)};
151 .def(py::init<quala::AndersonAccel::Params, quala::length_t>(),
"params"_a,
"n"_a)
153 return {kwargs_to_struct<quala::AndersonAccel::Params>(params), n};
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));
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);
178 "g_k"_a,
"r_k"_a,
"x_k_aa"_a)
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");
187 self.compute(g_k, std::move(r_k), x_k_aa);
195 py::class_<quala::BroydenGoodParams>(m,
"BroydenGoodParams",
196 "C++ documentation: :cpp:class:`quala::BroydenGoodParams`")
198 .def(py::init(&kwargs_to_struct<quala::BroydenGoodParams>))
199 .def(
"to_dict", &struct_to_dict<quala::BroydenGoodParams>)
207 py::class_<quala::BroydenGood>(m,
"BroydenGood",
208 "C++ documentation: :cpp:class:`quala::BroydenGood`")
209 .def(py::init<quala::BroydenGood::Params>(),
"params"_a)
211 return {kwargs_to_struct<quala::BroydenGood::Params>(params)};
214 .def(py::init<quala::BroydenGood::Params, quala::length_t>(),
"params"_a,
"n"_a)
216 return {kwargs_to_struct<quala::BroydenGood::Params>(params), n};
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);
234 "xk"_a,
"xkp1"_a,
"pk"_a,
"pkp1"_a,
"forced"_a =
false)
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");
244 "sk"_a,
"yk"_a,
"forced"_a =
false)
248 if (q.size() != self.
n())
249 throw std::invalid_argument(
"q dimension mismatch");
250 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.
bool apply(rvec q, real_t γ=-1)
Apply the inverse Hessian approximation to the given vector q.
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.
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ₖ₊₁.
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ₖ.
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.
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.
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.
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.
PYBIND11_MODULE(QUALA_MODULE_NAME, m)
real_t ϵ
Set to zero to disable CBFGS check.