proxsuite-nlp  0.10.0
A primal-dual augmented Lagrangian-type solver for nonlinear programming on manifolds.
Loading...
Searching...
No Matches
dense.hpp
Go to the documentation of this file.
1
6#pragma once
7
9
10namespace proxsuite {
11namespace nlp {
12namespace linalg {
13
14using Eigen::internal::LDLT_Traits;
15
16namespace backend {
17
18template <typename Scalar>
19void update_sign_matrix(SignMatrix &sign, const Scalar &akk) {
20 switch (sign) {
21 case SignMatrix::PositiveSemiDef:
22 if (akk < 0)
23 sign = SignMatrix::Indefinite;
24 break;
25 case SignMatrix::NegativeSemiDef:
26 if (akk > 0)
27 sign = SignMatrix::Indefinite;
28 break;
29 case SignMatrix::ZeroSign: {
30 if (akk > 0)
31 sign = SignMatrix::PositiveSemiDef;
32 else if (akk < 0)
33 sign = SignMatrix::NegativeSemiDef;
34 }
35 default:
36 break;
37 }
38}
39
44template <typename Derived>
45inline bool ldlt_in_place_unblocked(Eigen::MatrixBase<Derived> &a,
46 SignMatrix &sign) {
47 using Scalar = typename Derived::Scalar;
48 const isize n = a.rows();
49 if (n <= 1) {
50 if (n == 0)
51 sign = SignMatrix::ZeroSign;
52 else if (a(0, 0) > 0)
53 sign = SignMatrix::PositiveSemiDef;
54 else if (a(0, 0) < 0)
55 sign = SignMatrix::NegativeSemiDef;
56 else
57 sign = SignMatrix::ZeroSign;
58 return true;
59 }
60
61 isize j = 0;
62 while (true) {
63 auto l10 = a.row(j).head(j);
64 auto d0 = a.diagonal().head(j);
65 auto work = a.col(n - 1).head(j);
66
67 work = l10.transpose().cwiseProduct(d0);
68
69 Scalar &akk = a.coeffRef(j, j);
70 akk -= work.dot(l10);
71
72 update_sign_matrix(sign, akk);
73
74 if (j + 1 == n) {
75 return true;
76 }
77
78 const isize rem = n - j - 1;
79
80 auto l20 = a.bottomLeftCorner(rem, j);
81 auto l21 = a.col(j).tail(rem);
82
83 l21.noalias() -= l20 * work;
84 l21 *= 1 / akk;
85 ++j;
86 }
87}
88
89static constexpr isize UNBLK_THRESHOLD = 128;
90
93template <typename Derived>
94inline bool dense_ldlt_in_place(Eigen::MatrixBase<Derived> &a,
95 SignMatrix &sign) {
96 using PlainObject = typename Derived::PlainObject;
97 using MatrixRef = Eigen::Ref<PlainObject>;
98 const isize n = a.rows();
99 if (n <= UNBLK_THRESHOLD) {
100 return backend::ldlt_in_place_unblocked(a, sign);
101 } else {
102 const isize bs = (n + 1) / 2;
103 const isize rem = n - bs;
104
105 MatrixRef l00 = a.block(0, 0, bs, bs);
106 Eigen::Block<Derived> l10 = a.block(bs, 0, rem, bs);
107 MatrixRef l11 = a.block(bs, bs, rem, rem);
108
110 auto d0 = l00.diagonal();
111
112 l00.transpose()
113 .template triangularView<Eigen::UnitUpper>()
114 .template solveInPlace<Eigen::OnTheRight>(l10);
115
116 auto work = a.block(0, bs, bs, rem).transpose();
117 work = l10;
118 l10 = l10 * d0.asDiagonal().inverse();
119
120 l11.template triangularView<Eigen::Lower>() -= l10 * work.transpose();
121
122 return backend::dense_ldlt_in_place(l11, sign);
123 }
124}
125
128template <typename MatDerived, typename Rhs>
129inline bool dense_ldlt_solve_in_place(MatDerived &mat, Rhs &b) {
130 using Scalar = typename MatDerived::Scalar;
131 using Traits = LDLT_Traits<MatDerived, Eigen::Lower>;
132 Traits::getL(mat).solveInPlace(b);
133
134 using std::abs;
135 auto vecD(mat.diagonal());
136 const Scalar tol = std::numeric_limits<Scalar>::min();
137 for (isize i = 0; i < vecD.size(); ++i) {
138 if (abs(vecD(i)) > tol)
139 b.row(i) /= vecD(i);
140 else
141 b.row(i).setZero();
142 }
143
144 Traits::getU(mat).solveInPlace(b);
145 return true;
146}
147
148template <typename Scalar>
149inline void
151 typename math_types<Scalar>::MatrixRef res) {
152 using ConstMatrixRef = typename math_types<Scalar>::ConstMatrixRef;
153 using Traits = LDLT_Traits<ConstMatrixRef, Eigen::Lower>;
154 res = Traits::getU(mat) * res;
155
156 auto vecD = mat.diagonal();
157 res = vecD.asDiagonal() * res;
158
159 res = Traits::getL(mat) * res;
160}
161} // namespace backend
162
164template <typename _Scalar> struct DenseLDLT : ldlt_base<_Scalar> {
165 using Scalar = _Scalar;
168 using DView = typename Base::DView;
169 using MatrixType = MatrixXs;
170
171 DenseLDLT() = default;
172 explicit DenseLDLT(isize size) : Base(), m_matrix(size, size) {
173 m_matrix.setZero();
174 }
175
176 explicit DenseLDLT(MatrixRef a) : Base(), m_matrix(a) {
178 ? Eigen::Success
179 : Eigen::NumericalIssue;
180 }
181
182 DenseLDLT &compute(const ConstMatrixRef &mat) {
183 m_matrix = mat;
185 ? Eigen::Success
186 : Eigen::NumericalIssue;
187 return *this;
188 }
189
190 const MatrixXs &matrixLDLT() const { return m_matrix; }
191
192 template <typename Derived>
193 bool solveInPlace(Eigen::MatrixBase<Derived> &b) const {
195 }
196
197 template <typename Rhs>
198 typename Rhs::PlainObject solve(const Eigen::MatrixBase<Rhs> &rhs) const {
199 typename Rhs::PlainObject out = rhs;
200 solveInPlace(out);
201 return out;
202 }
203
204 MatrixXs reconstructedMatrix() const {
205 MatrixXs res(m_matrix.rows(), m_matrix.cols());
206 res.setIdentity();
208 return res;
209 }
210
211 inline DView vectorD() const { return Base::diag_view_impl(m_matrix); }
212
213protected:
215 using Base::m_info;
216 using Base::m_sign;
217};
218
219} // namespace linalg
220} // namespace nlp
221} // namespace proxsuite
#define PROXSUITE_NLP_DYNAMIC_TYPEDEFS(Scalar)
Definition math.hpp:26
static constexpr isize UNBLK_THRESHOLD
Definition dense.hpp:89
bool dense_ldlt_solve_in_place(MatDerived &mat, Rhs &b)
Definition dense.hpp:129
bool ldlt_in_place_unblocked(Eigen::MatrixBase< Derived > &a, SignMatrix &sign)
Definition dense.hpp:45
bool dense_ldlt_in_place(Eigen::MatrixBase< Derived > &a, SignMatrix &sign)
Definition dense.hpp:94
void update_sign_matrix(SignMatrix &sign, const Scalar &akk)
Definition dense.hpp:19
void dense_ldlt_reconstruct(typename math_types< Scalar >::ConstMatrixRef const &mat, typename math_types< Scalar >::MatrixRef res)
Definition dense.hpp:150
Main package namespace.
Definition bcl-params.hpp:5
A fast, recursive divide-and-conquer LDLT algorithm.
Definition dense.hpp:164
Eigen::ComputationInfo m_info
Definition ldlt-base.hpp:48
MatrixXs reconstructedMatrix() const
Definition dense.hpp:204
bool solveInPlace(Eigen::MatrixBase< Derived > &b) const
Definition dense.hpp:193
Rhs::PlainObject solve(const Eigen::MatrixBase< Rhs > &rhs) const
Definition dense.hpp:198
const MatrixXs & matrixLDLT() const
Definition dense.hpp:190
DenseLDLT & compute(const ConstMatrixRef &mat)
Definition dense.hpp:182
typename Base::DView DView
Definition dense.hpp:168
Base interface for LDLT solvers.
Definition ldlt-base.hpp:24
Eigen::ComputationInfo m_info
Definition ldlt-base.hpp:48
Eigen::Map< const VectorXs, Eigen::Unaligned, Eigen::InnerStride< Eigen::Dynamic > > DView
Definition ldlt-base.hpp:26
static DView diag_view_impl(Mat &&mat)
Definition ldlt-base.hpp:29
Typedefs for math (Eigen vectors, matrices) depending on scalar type.
Definition math.hpp:43