LMqrsolv.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
5 // Copyright (C) 2012 Desire Nuentsa <desire.nuentsa_wakam@inria.fr>
6 //
7 // This code initially comes from MINPACK whose original authors are:
8 // Copyright Jorge More - Argonne National Laboratory
9 // Copyright Burt Garbow - Argonne National Laboratory
10 // Copyright Ken Hillstrom - Argonne National Laboratory
11 //
12 // This Source Code Form is subject to the terms of the Minpack license
13 // (a BSD-like license) described in the campaigned CopyrightMINPACK.txt file.
14 
15 #ifndef EIGEN_LMQRSOLV_H
16 #define EIGEN_LMQRSOLV_H
17 
18 // IWYU pragma: private
19 #include "./InternalHeaderCheck.h"
20 
21 namespace Eigen {
22 
23 namespace internal {
24 
25 template <typename Scalar, int Rows, int Cols, typename PermIndex>
29  /* Local variables */
30  Index i, j, k;
31  Scalar temp;
32  Index n = s.cols();
35 
36  /* Function Body */
37  // the following will only change the lower triangular part of s, including
38  // the diagonal, though the diagonal is restored afterward
39 
40  /* copy r and (q transpose)*b to preserve input and initialize s. */
41  /* in particular, save the diagonal elements of r in x. */
42  x = s.diagonal();
43  wa = qtb;
44 
45  s.topLeftCorner(n, n).template triangularView<StrictlyLower>() = s.topLeftCorner(n, n).transpose();
46  /* eliminate the diagonal matrix d using a givens rotation. */
47  for (j = 0; j < n; ++j) {
48  /* prepare the row of d to be eliminated, locating the */
49  /* diagonal element using p from the qr factorization. */
50  const PermIndex l = iPerm.indices()(j);
51  if (diag[l] == 0.) break;
52  sdiag.tail(n - j).setZero();
53  sdiag[j] = diag[l];
54 
55  /* the transformations to eliminate the row of d */
56  /* modify only a single element of (q transpose)*b */
57  /* beyond the first n, which is initially zero. */
58  Scalar qtbpj = 0.;
59  for (k = j; k < n; ++k) {
60  /* determine a givens rotation which eliminates the */
61  /* appropriate element in the current row of d. */
62  givens.makeGivens(-s(k, k), sdiag[k]);
63 
64  /* compute the modified diagonal element of r and */
65  /* the modified element of ((q transpose)*b,0). */
66  s(k, k) = givens.c() * s(k, k) + givens.s() * sdiag[k];
67  temp = givens.c() * wa[k] + givens.s() * qtbpj;
68  qtbpj = -givens.s() * wa[k] + givens.c() * qtbpj;
69  wa[k] = temp;
70 
71  /* accumulate the transformation in the row of s. */
72  for (i = k + 1; i < n; ++i) {
73  temp = givens.c() * s(i, k) + givens.s() * sdiag[i];
74  sdiag[i] = -givens.s() * s(i, k) + givens.c() * sdiag[i];
75  s(i, k) = temp;
76  }
77  }
78  }
79 
80  /* solve the triangular system for z. if the system is */
81  /* singular, then obtain a least squares solution. */
82  Index nsing;
83  for (nsing = 0; nsing < n && sdiag[nsing] != 0; nsing++) {
84  }
85 
86  wa.tail(n - nsing).setZero();
87  s.topLeftCorner(nsing, nsing).transpose().template triangularView<Upper>().solveInPlace(wa.head(nsing));
88 
89  // restore
90  sdiag = s.diagonal();
91  s.diagonal() = x;
92 
93  /* permute the components of z back to components of x. */
94  x = iPerm * wa;
95 }
96 
97 template <typename Scalar, int Options_, typename Index>
101  /* Local variables */
102  typedef SparseMatrix<Scalar, RowMajor, Index> FactorType;
103  Index i, j, k, l;
104  Scalar temp;
105  Index n = s.cols();
107  JacobiRotation<Scalar> givens;
108 
109  /* Function Body */
110  // the following will only change the lower triangular part of s, including
111  // the diagonal, though the diagonal is restored afterward
112 
113  /* copy r and (q transpose)*b to preserve input and initialize R. */
114  wa = qtb;
115  FactorType R(s);
116  // Eliminate the diagonal matrix d using a givens rotation
117  for (j = 0; j < n; ++j) {
118  // Prepare the row of d to be eliminated, locating the
119  // diagonal element using p from the qr factorization
120  l = iPerm.indices()(j);
121  if (diag(l) == Scalar(0)) break;
122  sdiag.tail(n - j).setZero();
123  sdiag[j] = diag[l];
124  // the transformations to eliminate the row of d
125  // modify only a single element of (q transpose)*b
126  // beyond the first n, which is initially zero.
127 
128  Scalar qtbpj = 0;
129  // Browse the nonzero elements of row j of the upper triangular s
130  for (k = j; k < n; ++k) {
131  typename FactorType::InnerIterator itk(R, k);
132  for (; itk; ++itk) {
133  if (itk.index() < k)
134  continue;
135  else
136  break;
137  }
138  // At this point, we have the diagonal element R(k,k)
139  // Determine a givens rotation which eliminates
140  // the appropriate element in the current row of d
141  givens.makeGivens(-itk.value(), sdiag(k));
142 
143  // Compute the modified diagonal element of r and
144  // the modified element of ((q transpose)*b,0).
145  itk.valueRef() = givens.c() * itk.value() + givens.s() * sdiag(k);
146  temp = givens.c() * wa(k) + givens.s() * qtbpj;
147  qtbpj = -givens.s() * wa(k) + givens.c() * qtbpj;
148  wa(k) = temp;
149 
150  // Accumulate the transformation in the remaining k row/column of R
151  for (++itk; itk; ++itk) {
152  i = itk.index();
153  temp = givens.c() * itk.value() + givens.s() * sdiag(i);
154  sdiag(i) = -givens.s() * itk.value() + givens.c() * sdiag(i);
155  itk.valueRef() = temp;
156  }
157  }
158  }
159 
160  // Solve the triangular system for z. If the system is
161  // singular, then obtain a least squares solution
162  Index nsing;
163  for (nsing = 0; nsing < n && sdiag(nsing) != 0; nsing++) {
164  }
165 
166  wa.tail(n - nsing).setZero();
167  // x = wa;
168  wa.head(nsing) = R.topLeftCorner(nsing, nsing).template triangularView<Upper>().solve /*InPlace*/ (wa.head(nsing));
169 
170  sdiag = R.diagonal();
171  // Permute the components of z back to components of x
172  x = iPerm * wa;
173 }
174 } // end namespace internal
175 
176 } // end namespace Eigen
177 
178 #endif // EIGEN_LMQRSOLV_H
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
@ R
Definition: StatisticsVector.h:21
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
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Permutation matrix.
Definition: PermutationMatrix.h:280
const IndicesType & indices() const
Definition: PermutationMatrix.h:334
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
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 lmqrsolv(Matrix< Scalar, Rows, Cols > &s, const PermutationMatrix< Dynamic, Dynamic, PermIndex > &iPerm, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &sdiag)
Definition: LMqrsolv.h:26
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
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2