Home | History | Annotate | Download | only in QR

Lines Matching refs:Matrix

22   * \brief Householder QR decomposition of a matrix
24 * \tparam _MatrixType the type of the matrix of which we are computing the QR decomposition
26 * This class performs a QR decomposition of a matrix \b A into matrices \b Q and \b R
31 * by using Householder transformations. Here, \b Q a unitary matrix and \b R an upper triangular matrix.
59 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, (MatrixType::Flags&RowMajorBit) ? RowMajor : ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
84 /** \brief Constructs a QR factorization from a given matrix
86 * This constructor computes the QR factorization of the matrix \a matrix by calling
90 * HouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols());
91 * qr.compute(matrix);
97 explicit HouseholderQR(const EigenBase<InputType>& matrix)
98 : m_qr(matrix.rows(), matrix.cols()),
99 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
100 m_temp(matrix.cols()),
103 compute(matrix.derived());
107 /** \brief Constructs a QR factorization from a given matrix
115 explicit HouseholderQR(EigenBase<InputType>& matrix)
116 : m_qr(matrix.derived()),
117 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
118 m_temp(matrix.cols()),
124 /** This method finds a solution x to the equation Ax=b, where A is the matrix of which
146 /** This method returns an expression of the unitary matrix Q as a sequence of Householder transformations.
148 * The returned expression can directly be used to perform matrix products. It can also be assigned to a dense Matrix object.
149 * Here is an example showing how to recover the full or thin matrix Q, as well as how to perform matrix products using operator*:
160 /** \returns a reference to the matrix where the Householder QR decomposition is stored
170 HouseholderQR& compute(const EigenBase<InputType>& matrix) {
171 m_qr = matrix.derived();
176 /** \returns the absolute value of the determinant of the matrix of which
178 * (that is, O(n) where n is the dimension of the square matrix)
191 /** \returns the natural log of the absolute value of the determinant of the matrix of which
193 * (that is, O(n) where n is the dimension of the square matrix)
240 eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
248 eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
266 typedef Matrix<Scalar,MatrixQR::ColsAtCompileTime,1> TempType;
306 typedef Matrix<Scalar,Dynamic,1,ColMajor,MatrixQR::MaxColsAtCompileTime,1> TempType;
323 // partition the matrix:
357 // Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
372 /** Performs the QR factorization of the given matrix \a matrix. The result of