Home | History | Annotate | Download | only in LU

Lines Matching full:matrix

19   * \brief LU decomposition of a matrix with complete pivoting, and related features
21 * \param MatrixType the type of the matrix of which we are computing the LU decomposition
23 * This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A
32 * decomposition is inherently more stable and/or flexible. For example, when computing the kernel of a matrix,
33 * working with the SVD allows to select the smallest singular values of the matrix, something that
39 * As an exemple, here is how the original matrix can be retrieved:
83 * \param matrix the matrix of which to compute the LU decomposition.
86 FullPivLU(const MatrixType& matrix);
88 /** Computes the LU decomposition of the given matrix.
90 * \param matrix the matrix of which to compute the LU decomposition.
95 FullPivLU& compute(const MatrixType& matrix);
97 /** \returns the LU decomposition matrix: the upper-triangular part is U, the
127 /** \returns the permutation matrix P
137 /** \returns the permutation matrix Q
147 /** \returns the kernel of the matrix, also called its null-space. The columns of the returned matrix
150 * \note If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.
167 /** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
170 * \param originalMatrix the original matrix, of which *this is the LU decomposition.
175 * \note If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.
193 /** \return a solution x to the equation Ax=b, where A is the matrix of which
196 * \param b the right-hand-side of the equation to solve. Can be a vector or a matrix,
198 * b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.
220 /** \returns the determinant of the matrix of which
222 * (that is, O(n) where n is the dimension of the square matrix)
288 /** \returns the rank of the matrix of which *this is the LU decomposition.
304 /** \returns the dimension of the kernel of the matrix of which *this is the LU decomposition.
316 /** \returns true if the matrix of which *this is the LU decomposition represents an injective
329 /** \returns true if the matrix of which *this is the LU decomposition represents a surjective
342 /** \returns true if the matrix of which *this is the LU decomposition is invertible.
354 /** \returns the inverse of the matrix of which *this is the LU decomposition.
356 * \note If this matrix is not invertible, the returned matrix has undefined coefficients.
357 * Use isInvertible() to first determine whether this matrix is invertible.
364 eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
404 FullPivLU<MatrixType>::FullPivLU(const MatrixType& matrix)
405 : m_lu(matrix.rows(), matrix.cols()),
406 m_p(matrix.rows()),
407 m_q(matrix.cols()),
408 m_rowsTranspositions(matrix.rows()),
409 m_colsTranspositions(matrix.cols()),
413 compute(matrix);
417 FullPivLU<MatrixType>& FullPivLU<MatrixType>::compute(const MatrixType& matrix)
420 m_lu = matrix;
422 const Index size = matrix.diagonalSize();
423 const Index rows = matrix.rows();
424 const Index cols = matrix.cols();
428 m_rowsTranspositions.resize(matrix.rows());
429 m_colsTranspositions.resize(matrix.cols());
505 eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
509 /** \returns the matrix represented by the decomposition,
562 * Lemma: If the matrix A has the LU decomposition PAQ = LU,
576 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
584 // we construct a temporaty trapezoid matrix m, by taking the U matrix and
588 Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
601 // ok, we have our trapezoid matrix, we can apply the triangular solver.
644 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());