15 #ifndef EIGEN_LMQRSOLV_H
16 #define EIGEN_LMQRSOLV_H
25 template <
typename Scalar,
int Rows,
int Cols,
typename PermIndex>
45 s.topLeftCorner(
n,
n).template triangularView<StrictlyLower>() =
s.topLeftCorner(
n,
n).transpose();
47 for (
j = 0;
j <
n; ++
j) {
50 const PermIndex l = iPerm.
indices()(
j);
51 if (
diag[l] == 0.)
break;
59 for (
k =
j;
k <
n; ++
k) {
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;
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];
83 for (nsing = 0; nsing <
n && sdiag[nsing] != 0; nsing++) {
87 s.topLeftCorner(nsing, nsing).transpose().template triangularView<Upper>().solveInPlace(wa.head(nsing));
97 template <
typename Scalar,
int Options_,
typename Index>
117 for (
j = 0;
j <
n; ++
j) {
130 for (
k =
j;
k <
n; ++
k) {
131 typename FactorType::InnerIterator itk(
R,
k);
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;
151 for (++itk; itk; ++itk) {
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;
163 for (nsing = 0; nsing <
n && sdiag(nsing) != 0; nsing++) {
168 wa.head(nsing) =
R.topLeftCorner(nsing, nsing).template triangularView<Upper>().solve (wa.head(nsing));
170 sdiag =
R.diagonal();
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