qrsolv.h
Go to the documentation of this file.
1 // IWYU pragma: private
3 
4 namespace Eigen {
5 
6 namespace internal {
7 
8 // TODO : once qrsolv2 is removed, use ColPivHouseholderQR or PermutationMatrix instead of ipvt
9 template <typename Scalar>
11  // TODO : use a PermutationMatrix once lmpar is no more:
12  const VectorXi &ipvt, const Matrix<Scalar, Dynamic, 1> &diag, const Matrix<Scalar, Dynamic, 1> &qtb,
14 
15 {
16  typedef DenseIndex Index;
17 
18  /* Local variables */
19  Index i, j, k, l;
20  Scalar temp;
21  Index n = s.cols();
24 
25  /* Function Body */
26  // the following will only change the lower triangular part of s, including
27  // the diagonal, though the diagonal is restored afterward
28 
29  /* copy r and (q transpose)*b to preserve input and initialize s. */
30  /* in particular, save the diagonal elements of r in x. */
31  x = s.diagonal();
32  wa = qtb;
33 
34  s.topLeftCorner(n, n).template triangularView<StrictlyLower>() = s.topLeftCorner(n, n).transpose();
35 
36  /* eliminate the diagonal matrix d using a givens rotation. */
37  for (j = 0; j < n; ++j) {
38  /* prepare the row of d to be eliminated, locating the */
39  /* diagonal element using p from the qr factorization. */
40  l = ipvt[j];
41  if (diag[l] == 0.) break;
42  sdiag.tail(n - j).setZero();
43  sdiag[j] = diag[l];
44 
45  /* the transformations to eliminate the row of d */
46  /* modify only a single element of (q transpose)*b */
47  /* beyond the first n, which is initially zero. */
48  Scalar qtbpj = 0.;
49  for (k = j; k < n; ++k) {
50  /* determine a givens rotation which eliminates the */
51  /* appropriate element in the current row of d. */
52  givens.makeGivens(-s(k, k), sdiag[k]);
53 
54  /* compute the modified diagonal element of r and */
55  /* the modified element of ((q transpose)*b,0). */
56  s(k, k) = givens.c() * s(k, k) + givens.s() * sdiag[k];
57  temp = givens.c() * wa[k] + givens.s() * qtbpj;
58  qtbpj = -givens.s() * wa[k] + givens.c() * qtbpj;
59  wa[k] = temp;
60 
61  /* accumulate the transformation in the row of s. */
62  for (i = k + 1; i < n; ++i) {
63  temp = givens.c() * s(i, k) + givens.s() * sdiag[i];
64  sdiag[i] = -givens.s() * s(i, k) + givens.c() * sdiag[i];
65  s(i, k) = temp;
66  }
67  }
68  }
69 
70  /* solve the triangular system for z. if the system is */
71  /* singular, then obtain a least squares solution. */
72  Index nsing;
73  for (nsing = 0; nsing < n && sdiag[nsing] != 0; nsing++) {
74  }
75 
76  wa.tail(n - nsing).setZero();
77  s.topLeftCorner(nsing, nsing).transpose().template triangularView<Upper>().solveInPlace(wa.head(nsing));
78 
79  // restore
80  sdiag = s.diagonal();
81  s.diagonal() = x;
82 
83  /* permute the components of z back to components of x. */
84  for (j = 0; j < n; ++j) x[ipvt[j]] = wa[j];
85 }
86 
87 } // end namespace internal
88 
89 } // end namespace Eigen
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
SCALAR Scalar
Definition: bench_gemm.cpp:45
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:38
EIGEN_DEVICE_FUNC Scalar & s()
Definition: Jacobi.h:50
EIGEN_DEVICE_FUNC void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:152
EIGEN_DEVICE_FUNC Scalar & c()
Definition: Jacobi.h:48
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
RealScalar s
Definition: level1_cplx_impl.h:130
const char const char const char * diag
Definition: level2_impl.h:86
char char char int int * k
Definition: level2_impl.h:374
void qrsolv(Matrix< Scalar, Dynamic, Dynamic > &s, const VectorXi &ipvt, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &sdiag)
Definition: qrsolv.h:10
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:75
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2