4#ifndef LBFGSPP_BFGS_MAT_H
5#define LBFGSPP_BFGS_MAT_H
25template <
typename Scalar,
bool LBFGSB = false>
29 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
30 using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>;
31 using RefConstVec = Eigen::Ref<const Vector>;
32 using IndexSet = std::vector<int>;
50 BKLDLT<Scalar> m_permMsolver;
59 inline void reset(
int n,
int m)
72 m_permMinv.resize(2 * m, 2 * m);
74 m_permMinv.diagonal().setOnes();
79 inline void add_correction(
const RefConstVec& s,
const RefConstVec& y)
81 const int loc = m_ptr % m_m;
83 m_s.col(loc).noalias() = s;
84 m_y.col(loc).noalias() = y;
87 const Scalar ys = m_s.col(loc).dot(m_y.col(loc));
90 m_theta = m_y.col(loc).squaredNorm() / ys;
105 m_permMinv(loc, loc) = -ys;
109 Vector Ss = m_s.leftCols(m_ncorr).transpose() * m_s.col(loc);
110 m_permMinv.block(m_m + loc, m_m, 1, m_ncorr).noalias() = Ss.transpose();
111 m_permMinv.block(m_m, m_m + loc, m_ncorr, 1).noalias() = Ss;
125 const int len = m_ncorr - 1;
128 m_permMinv.block(m_m, loc, m_m, 1).setZero();
133 int yloc = (loc + m_m - 1) % m_m;
134 for (
int i = 0; i < len; i++)
136 m_permMinv(m_m + loc, yloc) = m_s.col(loc).dot(m_y.col(yloc));
137 yloc = (yloc + m_m - 1) % m_m;
141 m_permMinv.block(m_m, m_m, m_m, m_m) *= m_theta;
142 m_permMsolver.compute(m_permMinv);
143 m_permMinv.block(m_m, m_m, m_m, m_m) /= m_theta;
150 inline void apply_Hv(
const Vector& v,
const Scalar& a, Vector& res)
152 res.resize(v.size());
157 res.noalias() = a * v;
159 for (
int i = 0; i < m_ncorr; i++)
161 j = (j + m_m - 1) % m_m;
162 m_alpha[j] = m_s.col(j).dot(res) / m_ys[j];
163 res.noalias() -= m_alpha[j] * m_y.col(j);
170 for (
int i = 0; i < m_ncorr; i++)
172 const Scalar beta = m_y.col(j).dot(res) / m_ys[j];
173 res.noalias() += (m_alpha[j] - beta) * m_s.col(j);
181 inline Scalar theta()
const {
return m_theta; }
184 inline int num_corrections()
const {
return m_ncorr; }
189 inline void apply_Wtv(
const Vector& v, Vector& res)
const
191 res.resize(2 * m_ncorr);
192 res.head(m_ncorr).noalias() = m_y.leftCols(m_ncorr).transpose() * v;
193 res.tail(m_ncorr).noalias() = m_theta * m_s.leftCols(m_ncorr).transpose() * v;
199 inline Vector Wb(
int b)
const
201 Vector res(2 * m_ncorr);
202 for (
int j = 0; j < m_ncorr; j++)
205 res[m_ncorr + j] = m_s(b, j);
207 res.tail(m_ncorr) *= m_theta;
212 inline Matrix Wb(
const IndexSet& b)
const
214 const int nb = b.size();
215 const int* bptr = b.data();
216 Matrix res(nb, 2 * m_ncorr);
218 for (
int j = 0; j < m_ncorr; j++)
220 const Scalar* Yptr = &m_y(0, j);
221 const Scalar* Sptr = &m_s(0, j);
222 Scalar* resYptr = res.data() + j * nb;
223 Scalar* resSptr = resYptr + m_ncorr * nb;
224 for (
int i = 0; i < nb; i++)
226 const int row = bptr[i];
227 resYptr[i] = Yptr[row];
228 resSptr[i] = Sptr[row];
235 inline void apply_Mv(
const Vector& v, Vector& res)
const
237 res.resize(2 * m_ncorr);
241 Vector vpadding = Vector::Zero(2 * m_m);
242 vpadding.head(m_ncorr).noalias() = v.head(m_ncorr);
243 vpadding.segment(m_m, m_ncorr).noalias() = v.tail(m_ncorr);
246 m_permMsolver.solve_inplace(vpadding);
248 res.head(m_ncorr).noalias() = vpadding.head(m_ncorr);
249 res.tail(m_ncorr).noalias() = vpadding.segment(m_m, m_ncorr);
256 inline bool apply_WtPv(
const IndexSet& P_set,
const Vector& v, Vector& res,
bool test_zero =
false)
const
258 const int* Pptr = P_set.data();
259 const Scalar* vptr = v.data();
260 int nP = P_set.size();
264 std::vector<Scalar> v_reduced;
267 P_reduced.reserve(nP);
268 for (
int i = 0; i < nP; i++)
270 if (vptr[i] != Scalar(0))
272 P_reduced.push_back(Pptr[i]);
273 v_reduced.push_back(vptr[i]);
276 Pptr = P_reduced.data();
277 vptr = v_reduced.data();
278 nP = P_reduced.size();
281 res.resize(2 * m_ncorr);
282 if (m_ncorr < 1 || nP < 1)
288 for (
int j = 0; j < m_ncorr; j++)
290 Scalar resy = Scalar(0), ress = Scalar(0);
291 const Scalar* yptr = &m_y(0, j);
292 const Scalar* sptr = &m_s(0, j);
293 for (
int i = 0; i < nP; i++)
295 const int row = Pptr[i];
296 resy += yptr[row] * vptr[i];
297 ress += sptr[row] * vptr[i];
300 res[m_ncorr + j] = ress;
302 res.tail(m_ncorr) *= m_theta;
309 inline bool apply_PtWMv(
const IndexSet& P_set,
const Vector& v, Vector& res,
const Scalar& scale)
const
311 const int nP = P_set.size();
314 if (m_ncorr < 1 || nP < 1)
320 Mv.tail(m_ncorr) *= m_theta;
321 for (
int j = 0; j < m_ncorr; j++)
323 const Scalar* yptr = &m_y(0, j);
324 const Scalar* sptr = &m_s(0, j);
325 const Scalar Mvy = Mv[j], Mvs = Mv[m_ncorr + j];
326 for (
int i = 0; i < nP; i++)
328 const int row = P_set[i];
329 res[i] += Mvy * yptr[row] + Mvs * sptr[row];
336 inline bool apply_PtWMv(
const Matrix& WP,
const Vector& v, Vector& res,
const Scalar& scale)
const
338 const int nP = WP.rows();
340 if (m_ncorr < 1 || nP < 1)
349 Mv.tail(m_ncorr) *= m_theta;
350 res.noalias() = scale * (WP * Mv);
360 inline void compute_FtBAb(
361 const Matrix& WF,
const IndexSet& fv_set,
const IndexSet& newact_set,
const Vector& Wd,
const Vector& drt,
364 const int nact = newact_set.size();
365 const int nfree = WF.rows();
367 if (m_ncorr < 1 || nact < 1 || nfree < 1)
374 Vector rhs(2 * m_ncorr);
379 for (
int i = 0; i < nact; i++)
380 Ad[i] = drt[newact_set[i]];
381 apply_WtPv(newact_set, Ad, rhs);
387 for (
int i = 0; i < nfree; i++)
388 Fd[i] = drt[fv_set[i]];
390 rhs.noalias() = WF.transpose() * Fd;
391 rhs.tail(m_ncorr) *= m_theta;
392 rhs.noalias() = Wd - rhs;
395 apply_PtWMv(WF, rhs, res, Scalar(-1));
403 inline void solve_PtBP(
const Matrix& WP,
const Vector& v, Vector& res)
const
405 const int nP = WP.rows();
407 if (m_ncorr < 1 || nP < 1)
409 res.noalias() = v / m_theta;
415 Matrix mid(2 * m_ncorr, 2 * m_ncorr);
417 for (
int j = 0; j < m_ncorr; j++)
419 mid.col(j).segment(j, m_ncorr - j).noalias() = m_permMinv.col(j).segment(j, m_ncorr - j) -
420 WP.block(0, j, nP, m_ncorr - j).transpose() * WP.col(j) / m_theta;
423 mid.block(m_ncorr, 0, m_ncorr, m_ncorr).noalias() = m_permMinv.block(m_m, 0, m_ncorr, m_ncorr) -
424 WP.rightCols(m_ncorr).transpose() * WP.leftCols(m_ncorr);
426 for (
int j = 0; j < m_ncorr; j++)
428 mid.col(m_ncorr + j).segment(m_ncorr + j, m_ncorr - j).noalias() = m_theta *
429 (m_permMinv.col(m_m + j).segment(m_m + j, m_ncorr - j) - WP.rightCols(m_ncorr - j).transpose() * WP.col(m_ncorr + j));
432 BKLDLT<Scalar> midsolver(mid);
434 Vector WPv = WP.transpose() * v;
435 WPv.tail(m_ncorr) *= m_theta;
436 midsolver.solve_inplace(WPv);
437 WPv.tail(m_ncorr) *= m_theta;
438 res.noalias() = v / m_theta + (WP * WPv) / (m_theta * m_theta);
444 inline bool apply_PtBQv(
const Matrix& WP,
const IndexSet& Q_set,
const Vector& v, Vector& res,
bool test_zero =
false)
const
446 const int nP = WP.rows();
447 const int nQ = Q_set.size();
449 if (m_ncorr < 1 || nP < 1 || nQ < 1)
456 bool nonzero = apply_WtPv(Q_set, v, WQtv, test_zero);
464 apply_Mv(WQtv, MWQtv);
465 MWQtv.tail(m_ncorr) *= m_theta;
466 res.noalias() = -WP * MWQtv;
470 inline bool apply_PtBQv(
const Matrix& WP,
const Matrix& WQ,
const Vector& v, Vector& res)
const
472 const int nP = WP.rows();
473 const int nQ = WQ.rows();
475 if (m_ncorr < 1 || nP < 1 || nQ < 1)
482 Vector WQtv = WQ.transpose() * v;
483 WQtv.tail(m_ncorr) *= m_theta;
485 apply_Mv(WQtv, MWQtv);
486 MWQtv.tail(m_ncorr) *= m_theta;
487 res.noalias() = -WP * MWQtv;