Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
UpperBidiagonalization.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2010 Benoit Jacob <[email protected]>
5// Copyright (C) 2013-2014 Gael Guennebaud <[email protected]>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_BIDIAGONALIZATION_H
12#define EIGEN_BIDIAGONALIZATION_H
13
14// IWYU pragma: private
15#include "./InternalHeaderCheck.h"
16
17namespace Eigen {
18
19namespace internal {
20// UpperBidiagonalization will probably be replaced by a Bidiagonalization class, don't want to make it stable API.
21// At the same time, it's useful to keep for now as it's about the only thing that is testing the BandMatrix class.
22
23template <typename MatrixType_>
24class UpperBidiagonalization {
25 public:
26 typedef MatrixType_ MatrixType;
27 enum {
28 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
29 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
30 ColsAtCompileTimeMinusOne = internal::decrement_size<ColsAtCompileTime>::ret
31 };
32 typedef typename MatrixType::Scalar Scalar;
33 typedef typename MatrixType::RealScalar RealScalar;
34 typedef Eigen::Index Index;
35 typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
36 typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
37 typedef BandMatrix<RealScalar, ColsAtCompileTime, ColsAtCompileTime, 1, 0, RowMajor> BidiagonalType;
38 typedef Matrix<Scalar, ColsAtCompileTime, 1> DiagVectorType;
39 typedef Matrix<Scalar, ColsAtCompileTimeMinusOne, 1> SuperDiagVectorType;
40 typedef HouseholderSequence<
41 const MatrixType, const internal::remove_all_t<typename Diagonal<const MatrixType, 0>::ConjugateReturnType> >
42 HouseholderUSequenceType;
43 typedef HouseholderSequence<const internal::remove_all_t<typename MatrixType::ConjugateReturnType>,
44 Diagonal<const MatrixType, 1>, OnTheRight>
45 HouseholderVSequenceType;
46
53 UpperBidiagonalization() : m_householder(), m_bidiagonal(0, 0), m_isInitialized(false) {}
54
55 explicit UpperBidiagonalization(const MatrixType& matrix)
56 : m_householder(matrix.rows(), matrix.cols()),
57 m_bidiagonal(matrix.cols(), matrix.cols()),
58 m_isInitialized(false) {
59 compute(matrix);
60 }
61
62 UpperBidiagonalization(Index rows, Index cols)
63 : m_householder(rows, cols), m_bidiagonal(cols, cols), m_isInitialized(false) {}
64
65 UpperBidiagonalization& compute(const MatrixType& matrix);
66 UpperBidiagonalization& computeUnblocked(const MatrixType& matrix);
67
68 const MatrixType& householder() const { return m_householder; }
69 const BidiagonalType& bidiagonal() const { return m_bidiagonal; }
70
71 const HouseholderUSequenceType householderU() const {
72 eigen_assert(m_isInitialized && "UpperBidiagonalization is not initialized.");
73 return HouseholderUSequenceType(m_householder, m_householder.diagonal().conjugate());
74 }
75
76 const HouseholderVSequenceType householderV() // const here gives nasty errors and i'm lazy
77 {
78 eigen_assert(m_isInitialized && "UpperBidiagonalization is not initialized.");
79 return HouseholderVSequenceType(m_householder.conjugate(), m_householder.const_derived().template diagonal<1>())
80 .setLength(m_householder.cols() - 1)
81 .setShift(1);
82 }
83
84 protected:
85 MatrixType m_householder;
86 BidiagonalType m_bidiagonal;
87 bool m_isInitialized;
88};
89
90// Standard upper bidiagonalization without fancy optimizations
91// This version should be faster for small matrix size
92template <typename MatrixType>
93void upperbidiagonalization_inplace_unblocked(MatrixType& mat, typename MatrixType::RealScalar* diagonal,
94 typename MatrixType::RealScalar* upper_diagonal,
95 typename MatrixType::Scalar* tempData = 0) {
96 typedef typename MatrixType::Scalar Scalar;
97
98 Index rows = mat.rows();
99 Index cols = mat.cols();
100
101 typedef Matrix<Scalar, Dynamic, 1, ColMajor, MatrixType::MaxRowsAtCompileTime, 1> TempType;
102 TempType tempVector;
103 if (tempData == 0) {
104 tempVector.resize(rows);
105 tempData = tempVector.data();
106 }
107
108 for (Index k = 0; /* breaks at k==cols-1 below */; ++k) {
109 Index remainingRows = rows - k;
110 Index remainingCols = cols - k - 1;
111
112 // construct left householder transform in-place in A
113 mat.col(k).tail(remainingRows).makeHouseholderInPlace(mat.coeffRef(k, k), diagonal[k]);
114 // apply householder transform to remaining part of A on the left
115 mat.bottomRightCorner(remainingRows, remainingCols)
116 .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows - 1), mat.coeff(k, k), tempData);
117
118 if (k == cols - 1) break;
119
120 // construct right householder transform in-place in mat
121 mat.row(k).tail(remainingCols).makeHouseholderInPlace(mat.coeffRef(k, k + 1), upper_diagonal[k]);
122 // apply householder transform to remaining part of mat on the left
123 mat.bottomRightCorner(remainingRows - 1, remainingCols)
124 .applyHouseholderOnTheRight(mat.row(k).tail(remainingCols - 1).adjoint(), mat.coeff(k, k + 1), tempData);
125 }
126}
127
145template <typename MatrixType>
146void upperbidiagonalization_blocked_helper(
147 MatrixType& A, typename MatrixType::RealScalar* diagonal, typename MatrixType::RealScalar* upper_diagonal, Index bs,
148 Ref<Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, traits<MatrixType>::Flags & RowMajorBit> > X,
149 Ref<Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, traits<MatrixType>::Flags & RowMajorBit> > Y) {
150 typedef typename MatrixType::Scalar Scalar;
151 typedef typename MatrixType::RealScalar RealScalar;
152 typedef typename NumTraits<RealScalar>::Literal Literal;
153 static constexpr int StorageOrder = (traits<MatrixType>::Flags & RowMajorBit) ? RowMajor : ColMajor;
154 typedef InnerStride<StorageOrder == ColMajor ? 1 : Dynamic> ColInnerStride;
155 typedef InnerStride<StorageOrder == ColMajor ? Dynamic : 1> RowInnerStride;
156 typedef Ref<Matrix<Scalar, Dynamic, 1>, 0, ColInnerStride> SubColumnType;
157 typedef Ref<Matrix<Scalar, 1, Dynamic>, 0, RowInnerStride> SubRowType;
158 typedef Ref<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > SubMatType;
159
160 Index brows = A.rows();
161 Index bcols = A.cols();
162
163 Scalar tau_u, tau_u_prev(0), tau_v;
164
165 for (Index k = 0; k < bs; ++k) {
166 Index remainingRows = brows - k;
167 Index remainingCols = bcols - k - 1;
168
169 SubMatType X_k1(X.block(k, 0, remainingRows, k));
170 SubMatType V_k1(A.block(k, 0, remainingRows, k));
171
172 // 1 - update the k-th column of A
173 SubColumnType v_k = A.col(k).tail(remainingRows);
174 v_k -= V_k1 * Y.row(k).head(k).adjoint();
175 if (k) v_k -= X_k1 * A.col(k).head(k);
176
177 // 2 - construct left Householder transform in-place
178 v_k.makeHouseholderInPlace(tau_v, diagonal[k]);
179
180 if (k + 1 < bcols) {
181 SubMatType Y_k(Y.block(k + 1, 0, remainingCols, k + 1));
182 SubMatType U_k1(A.block(0, k + 1, k, remainingCols));
183
184 // this eases the application of Householder transforAions
185 // A(k,k) will store tau_v later
186 A(k, k) = Scalar(1);
187
188 // 3 - Compute y_k^T = tau_v * ( A^T*v_k - Y_k-1*V_k-1^T*v_k - U_k-1*X_k-1^T*v_k )
189 {
190 SubColumnType y_k(Y.col(k).tail(remainingCols));
191
192 // let's use the beginning of column k of Y as a temporary vector
193 SubColumnType tmp(Y.col(k).head(k));
194 y_k.noalias() = A.block(k, k + 1, remainingRows, remainingCols).adjoint() * v_k; // bottleneck
195 tmp.noalias() = V_k1.adjoint() * v_k;
196 y_k.noalias() -= Y_k.leftCols(k) * tmp;
197 tmp.noalias() = X_k1.adjoint() * v_k;
198 y_k.noalias() -= U_k1.adjoint() * tmp;
199 y_k *= numext::conj(tau_v);
200 }
201
202 // 4 - update k-th row of A (it will become u_k)
203 SubRowType u_k(A.row(k).tail(remainingCols));
204 u_k = u_k.conjugate();
205 {
206 u_k -= Y_k * A.row(k).head(k + 1).adjoint();
207 if (k) u_k -= U_k1.adjoint() * X.row(k).head(k).adjoint();
208 }
209
210 // 5 - construct right Householder transform in-place
211 u_k.makeHouseholderInPlace(tau_u, upper_diagonal[k]);
212
213 // this eases the application of Householder transformations
214 // A(k,k+1) will store tau_u later
215 A(k, k + 1) = Scalar(1);
216
217 // 6 - Compute x_k = tau_u * ( A*u_k - X_k-1*U_k-1^T*u_k - V_k*Y_k^T*u_k )
218 {
219 SubColumnType x_k(X.col(k).tail(remainingRows - 1));
220
221 // let's use the beginning of column k of X as a temporary vectors
222 // note that tmp0 and tmp1 overlaps
223 SubColumnType tmp0(X.col(k).head(k)), tmp1(X.col(k).head(k + 1));
224
225 x_k.noalias() = A.block(k + 1, k + 1, remainingRows - 1, remainingCols) * u_k.transpose(); // bottleneck
226 tmp0.noalias() = U_k1 * u_k.transpose();
227 x_k.noalias() -= X_k1.bottomRows(remainingRows - 1) * tmp0;
228 tmp1.noalias() = Y_k.adjoint() * u_k.transpose();
229 x_k.noalias() -= A.block(k + 1, 0, remainingRows - 1, k + 1) * tmp1;
230 x_k *= numext::conj(tau_u);
231 tau_u = numext::conj(tau_u);
232 u_k = u_k.conjugate();
233 }
234
235 if (k > 0) A.coeffRef(k - 1, k) = tau_u_prev;
236 tau_u_prev = tau_u;
237 } else
238 A.coeffRef(k - 1, k) = tau_u_prev;
239
240 A.coeffRef(k, k) = tau_v;
241 }
242
243 if (bs < bcols) A.coeffRef(bs - 1, bs) = tau_u_prev;
244
245 // update A22
246 if (bcols > bs && brows > bs) {
247 SubMatType A11(A.bottomRightCorner(brows - bs, bcols - bs));
248 SubMatType A10(A.block(bs, 0, brows - bs, bs));
249 SubMatType A01(A.block(0, bs, bs, bcols - bs));
250 Scalar tmp = A01(bs - 1, 0);
251 A01(bs - 1, 0) = Literal(1);
252 A11.noalias() -= A10 * Y.topLeftCorner(bcols, bs).bottomRows(bcols - bs).adjoint();
253 A11.noalias() -= X.topLeftCorner(brows, bs).bottomRows(brows - bs) * A01;
254 A01(bs - 1, 0) = tmp;
255 }
256}
257
265template <typename MatrixType, typename BidiagType>
266void upperbidiagonalization_inplace_blocked(MatrixType& A, BidiagType& bidiagonal, Index maxBlockSize = 32,
267 typename MatrixType::Scalar* /*tempData*/ = 0) {
268 typedef typename MatrixType::Scalar Scalar;
269 typedef Block<MatrixType, Dynamic, Dynamic> BlockType;
270
271 Index rows = A.rows();
272 Index cols = A.cols();
273 Index size = (std::min)(rows, cols);
274
275 // X and Y are work space
276 static constexpr int StorageOrder = (traits<MatrixType>::Flags & RowMajorBit) ? RowMajor : ColMajor;
277 Matrix<Scalar, MatrixType::RowsAtCompileTime, Dynamic, StorageOrder, MatrixType::MaxRowsAtCompileTime> X(
278 rows, maxBlockSize);
279 Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic, StorageOrder, MatrixType::MaxColsAtCompileTime> Y(
280 cols, maxBlockSize);
281 Index blockSize = (std::min)(maxBlockSize, size);
282
283 Index k = 0;
284 for (k = 0; k < size; k += blockSize) {
285 Index bs = (std::min)(size - k, blockSize); // actual size of the block
286 Index brows = rows - k; // rows of the block
287 Index bcols = cols - k; // columns of the block
288
289 // partition the matrix A:
290 //
291 // | A00 A01 A02 |
292 // | |
293 // A = | A10 A11 A12 |
294 // | |
295 // | A20 A21 A22 |
296 //
297 // where A11 is a bs x bs diagonal block,
298 // and let:
299 // | A11 A12 |
300 // B = | |
301 // | A21 A22 |
302
303 BlockType B = A.block(k, k, brows, bcols);
304
305 // This stage performs the bidiagonalization of A11, A21, A12, and updating of A22.
306 // Finally, the algorithm continue on the updated A22.
307 //
308 // However, if B is too small, or A22 empty, then let's use an unblocked strategy
309
310 auto upper_diagonal = bidiagonal.template diagonal<1>();
311 typename MatrixType::RealScalar* upper_diagonal_ptr =
312 upper_diagonal.size() > 0 ? &upper_diagonal.coeffRef(k) : nullptr;
313
314 if (k + bs == cols || bcols < 48) // somewhat arbitrary threshold
315 {
316 upperbidiagonalization_inplace_unblocked(B, &(bidiagonal.template diagonal<0>().coeffRef(k)), upper_diagonal_ptr,
317 X.data());
318 break; // We're done
319 } else {
320 upperbidiagonalization_blocked_helper<BlockType>(B, &(bidiagonal.template diagonal<0>().coeffRef(k)),
321 upper_diagonal_ptr, bs, X.topLeftCorner(brows, bs),
322 Y.topLeftCorner(bcols, bs));
323 }
324 }
325}
326
327template <typename MatrixType_>
328UpperBidiagonalization<MatrixType_>& UpperBidiagonalization<MatrixType_>::computeUnblocked(const MatrixType_& matrix) {
329 Index rows = matrix.rows();
330 Index cols = matrix.cols();
331 EIGEN_ONLY_USED_FOR_DEBUG(cols);
332
333 eigen_assert(rows >= cols && "UpperBidiagonalization is only for Arices satisfying rows>=cols.");
334
335 m_householder = matrix;
336
337 ColVectorType temp(rows);
338
339 upperbidiagonalization_inplace_unblocked(m_householder, &(m_bidiagonal.template diagonal<0>().coeffRef(0)),
340 &(m_bidiagonal.template diagonal<1>().coeffRef(0)), temp.data());
341
342 m_isInitialized = true;
343 return *this;
344}
345
346template <typename MatrixType_>
347UpperBidiagonalization<MatrixType_>& UpperBidiagonalization<MatrixType_>::compute(const MatrixType_& matrix) {
348 Index rows = matrix.rows();
349 Index cols = matrix.cols();
350 EIGEN_ONLY_USED_FOR_DEBUG(rows);
351 EIGEN_ONLY_USED_FOR_DEBUG(cols);
352
353 eigen_assert(rows >= cols && "UpperBidiagonalization is only for Arices satisfying rows>=cols.");
354
355 m_householder = matrix;
356 upperbidiagonalization_inplace_blocked(m_householder, m_bidiagonal);
357
358 m_isInitialized = true;
359 return *this;
360}
361
362#if 0
367template<typename Derived>
368const UpperBidiagonalization<typename MatrixBase<Derived>::PlainObject>
369MatrixBase<Derived>::bidiagonalization() const
370{
371 return UpperBidiagonalization<PlainObject>(eval());
372}
373#endif
374
375} // end namespace internal
376
377} // end namespace Eigen
378
379#endif // EIGEN_BIDIAGONALIZATION_H
@ ColMajor
Definition Constants.h:318
@ RowMajor
Definition Constants.h:320
@ OnTheRight
Definition Constants.h:333
const unsigned int RowMajorBit
Definition Constants.h:70
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25