aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Eigen2Support/QR.h
blob: 60fc21f56b9eebd380e43ca7738f804fcd46b0bb (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2011 Benoit Jacob <jacob.benoit.1@gmail.com>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.

#ifndef EIGEN2_QR_H
#define EIGEN2_QR_H

namespace Eigen { 

template<typename MatrixType>
class QR : public HouseholderQR<MatrixType>
{
  public:

    typedef HouseholderQR<MatrixType> Base;
    typedef Block<const MatrixType, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> MatrixRBlockType;

    QR() : Base() {}

    template<typename T>
    explicit QR(const T& t) : Base(t) {}

    template<typename OtherDerived, typename ResultType>
    bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
    {
      *result = static_cast<const Base*>(this)->solve(b);
      return true;
    }

    MatrixType matrixQ(void) const {
      MatrixType ret = MatrixType::Identity(this->rows(), this->cols());
      ret = this->householderQ() * ret;
      return ret;
    }

    bool isFullRank() const {
      return true;
    }
    
    const TriangularView<MatrixRBlockType, UpperTriangular>
    matrixR(void) const
    {
      int cols = this->cols();
      return MatrixRBlockType(this->matrixQR(), 0, 0, cols, cols).template triangularView<UpperTriangular>();
    }
};

/** \return the QR decomposition of \c *this.
  *
  * \sa class QR
  */
template<typename Derived>
const QR<typename MatrixBase<Derived>::PlainObject>
MatrixBase<Derived>::qr() const
{
  return QR<PlainObject>(eval());
}

} // end namespace Eigen

#endif // EIGEN2_QR_H