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 #include "./InternalHeaderCheck.h"
19 
20 namespace Eigen {
21 
22 namespace internal {
23 
24 template <typename Scalar,int Rows, int Cols, typename PermIndex>
25 void lmqrsolv(
28  const Matrix<Scalar,Dynamic,1> &diag,
29  const Matrix<Scalar,Dynamic,1> &qtb,
32 {
33  /* Local variables */
34  Index i, j, k;
35  Scalar temp;
36  Index n = s.cols();
39 
40  /* Function Body */
41  // the following will only change the lower triangular part of s, including
42  // the diagonal, though the diagonal is restored afterward
43 
44  /* copy r and (q transpose)*b to preserve input and initialize s. */
45  /* in particular, save the diagonal elements of r in x. */
46  x = s.diagonal();
47  wa = qtb;
48 
49 
50  s.topLeftCorner(n,n).template triangularView<StrictlyLower>() = s.topLeftCorner(n,n).transpose();
51  /* eliminate the diagonal matrix d using a givens rotation. */
52  for (j = 0; j < n; ++j) {
53 
54  /* prepare the row of d to be eliminated, locating the */
55  /* diagonal element using p from the qr factorization. */
56  const PermIndex l = iPerm.indices()(j);
57  if (diag[l] == 0.)
58  break;
59  sdiag.tail(n-j).setZero();
60  sdiag[j] = diag[l];
61 
62  /* the transformations to eliminate the row of d */
63  /* modify only a single element of (q transpose)*b */
64  /* beyond the first n, which is initially zero. */
65  Scalar qtbpj = 0.;
66  for (k = j; k < n; ++k) {
67  /* determine a givens rotation which eliminates the */
68  /* appropriate element in the current row of d. */
69  givens.makeGivens(-s(k,k), sdiag[k]);
70 
71  /* compute the modified diagonal element of r and */
72  /* the modified element of ((q transpose)*b,0). */
73  s(k,k) = givens.c() * s(k,k) + givens.s() * sdiag[k];
74  temp = givens.c() * wa[k] + givens.s() * qtbpj;
75  qtbpj = -givens.s() * wa[k] + givens.c() * qtbpj;
76  wa[k] = temp;
77 
78  /* accumulate the transformation in the row of s. */
79  for (i = k+1; i<n; ++i) {
80  temp = givens.c() * s(i,k) + givens.s() * sdiag[i];
81  sdiag[i] = -givens.s() * s(i,k) + givens.c() * sdiag[i];
82  s(i,k) = temp;
83  }
84  }
85  }
86 
87  /* solve the triangular system for z. if the system is */
88  /* singular, then obtain a least squares solution. */
89  Index nsing;
90  for(nsing=0; nsing<n && sdiag[nsing]!=0; nsing++) {}
91 
92  wa.tail(n-nsing).setZero();
93  s.topLeftCorner(nsing, nsing).transpose().template triangularView<Upper>().solveInPlace(wa.head(nsing));
94 
95  // restore
96  sdiag = s.diagonal();
97  s.diagonal() = x;
98 
99  /* permute the components of z back to components of x. */
100  x = iPerm * wa;
101 }
102 
103 template <typename Scalar, int Options_, typename Index>
104 void lmqrsolv(
107  const Matrix<Scalar,Dynamic,1> &diag,
108  const Matrix<Scalar,Dynamic,1> &qtb,
111 {
112  /* Local variables */
113  typedef SparseMatrix<Scalar,RowMajor,Index> FactorType;
114  Index i, j, k, l;
115  Scalar temp;
116  Index n = s.cols();
118  JacobiRotation<Scalar> givens;
119 
120  /* Function Body */
121  // the following will only change the lower triangular part of s, including
122  // the diagonal, though the diagonal is restored afterward
123 
124  /* copy r and (q transpose)*b to preserve input and initialize R. */
125  wa = qtb;
126  FactorType R(s);
127  // Eliminate the diagonal matrix d using a givens rotation
128  for (j = 0; j < n; ++j)
129  {
130  // Prepare the row of d to be eliminated, locating the
131  // diagonal element using p from the qr factorization
132  l = iPerm.indices()(j);
133  if (diag(l) == Scalar(0))
134  break;
135  sdiag.tail(n-j).setZero();
136  sdiag[j] = diag[l];
137  // the transformations to eliminate the row of d
138  // modify only a single element of (q transpose)*b
139  // beyond the first n, which is initially zero.
140 
141  Scalar qtbpj = 0;
142  // Browse the nonzero elements of row j of the upper triangular s
143  for (k = j; k < n; ++k)
144  {
145  typename FactorType::InnerIterator itk(R,k);
146  for (; itk; ++itk){
147  if (itk.index() < k) continue;
148  else break;
149  }
150  //At this point, we have the diagonal element R(k,k)
151  // Determine a givens rotation which eliminates
152  // the appropriate element in the current row of d
153  givens.makeGivens(-itk.value(), sdiag(k));
154 
155  // Compute the modified diagonal element of r and
156  // the modified element of ((q transpose)*b,0).
157  itk.valueRef() = givens.c() * itk.value() + givens.s() * sdiag(k);
158  temp = givens.c() * wa(k) + givens.s() * qtbpj;
159  qtbpj = -givens.s() * wa(k) + givens.c() * qtbpj;
160  wa(k) = temp;
161 
162  // Accumulate the transformation in the remaining k row/column of R
163  for (++itk; itk; ++itk)
164  {
165  i = itk.index();
166  temp = givens.c() * itk.value() + givens.s() * sdiag(i);
167  sdiag(i) = -givens.s() * itk.value() + givens.c() * sdiag(i);
168  itk.valueRef() = temp;
169  }
170  }
171  }
172 
173  // Solve the triangular system for z. If the system is
174  // singular, then obtain a least squares solution
175  Index nsing;
176  for(nsing = 0; nsing<n && sdiag(nsing) !=0; nsing++) {}
177 
178  wa.tail(n-nsing).setZero();
179 // x = wa;
180  wa.head(nsing) = R.topLeftCorner(nsing,nsing).template triangularView<Upper>().solve/*InPlace*/(wa.head(nsing));
181 
182  sdiag = R.diagonal();
183  // Permute the components of z back to components of x
184  x = iPerm * wa;
185 }
186 } // end namespace internal
187 
188 } // end namespace Eigen
189 
190 #endif // EIGEN_LMQRSOLV_H
int n
int i
TransposeReturnType transpose()
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
DiagonalReturnType diagonal()
IndicesType & indices()
Derived & setZero(Index rows, Index cols)
Index cols() const
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:25
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
std::ptrdiff_t j