LMcovar.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 // This code initially comes from MINPACK whose original authors are:
5 // Copyright Jorge More - Argonne National Laboratory
6 // Copyright Burt Garbow - Argonne National Laboratory
7 // Copyright Ken Hillstrom - Argonne National Laboratory
8 //
9 // This Source Code Form is subject to the terms of the Minpack license
10 // (a BSD-like license) described in the campaigned CopyrightMINPACK.txt file.
11 
12 #ifndef EIGEN_LMCOVAR_H
13 #define EIGEN_LMCOVAR_H
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21 
22 template <typename Scalar>
23 void covar(Matrix<Scalar, Dynamic, Dynamic>& r, const VectorXi& ipvt,
25  using std::abs;
26  /* Local variables */
27  Index i, j, k, l, ii, jj;
28  bool sing;
29  Scalar temp;
30 
31  /* Function Body */
32  const Index n = r.cols();
33  const Scalar tolr = tol * abs(r(0, 0));
35  eigen_assert(ipvt.size() == n);
36 
37  /* form the inverse of r in the full upper triangle of r. */
38  l = -1;
39  for (k = 0; k < n; ++k)
40  if (abs(r(k, k)) > tolr) {
41  r(k, k) = 1. / r(k, k);
42  for (j = 0; j <= k - 1; ++j) {
43  temp = r(k, k) * r(j, k);
44  r(j, k) = 0.;
45  r.col(k).head(j + 1) -= r.col(j).head(j + 1) * temp;
46  }
47  l = k;
48  }
49 
50  /* form the full upper triangle of the inverse of (r transpose)*r */
51  /* in the full upper triangle of r. */
52  for (k = 0; k <= l; ++k) {
53  for (j = 0; j <= k - 1; ++j) r.col(j).head(j + 1) += r.col(k).head(j + 1) * r(j, k);
54  r.col(k).head(k + 1) *= r(k, k);
55  }
56 
57  /* form the full lower triangle of the covariance matrix */
58  /* in the strict lower triangle of r and in wa. */
59  for (j = 0; j < n; ++j) {
60  jj = ipvt[j];
61  sing = j > l;
62  for (i = 0; i <= j; ++i) {
63  if (sing) r(i, j) = 0.;
64  ii = ipvt[i];
65  if (ii > jj) r(ii, jj) = r(i, j);
66  if (ii < jj) r(jj, ii) = r(i, j);
67  }
68  wa[jj] = r(j, j);
69  }
70 
71  /* symmetrize the covariance matrix in r. */
72  r.topLeftCorner(n, n).template triangularView<StrictlyUpper>() = r.topLeftCorner(n, n).transpose();
73  r.diagonal() = wa;
74 }
75 
76 } // end namespace internal
77 
78 } // end namespace Eigen
79 
80 #endif // EIGEN_LMCOVAR_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define eigen_assert(x)
Definition: Macros.h:910
SCALAR Scalar
Definition: bench_gemm.cpp:45
char char char int int * k
Definition: level2_impl.h:374
void covar(Matrix< Scalar, Dynamic, Dynamic > &r, const VectorXi &ipvt, Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
Definition: LMcovar.h:23
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
r
Definition: UniformPSDSelfTest.py:20
Definition: Eigen_Colamd.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2