11 #ifndef EIGEN_PARTIALLU_H    12 #define EIGEN_PARTIALLU_H    29 template<
typename T,
typename Derived>
    35 template<
typename T,
typename Derived>
    76   : 
public SolverBase<PartialPivLU<_MatrixType> >
    85       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
    86       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
   115     template<
typename InputType>
   125     template<
typename InputType>
   128     template<
typename InputType>
   143       eigen_assert(m_isInitialized && 
"PartialPivLU is not initialized.");
   151       eigen_assert(m_isInitialized && 
"PartialPivLU is not initialized.");
   173     template<
typename Rhs>
   177       eigen_assert(m_isInitialized && 
"PartialPivLU is not initialized.");
   186       eigen_assert(m_isInitialized && 
"PartialPivLU is not initialized.");
   199       eigen_assert(m_isInitialized && 
"PartialPivLU is not initialized.");
   216     Scalar determinant() 
const;
   218     MatrixType reconstructedMatrix() 
const;
   223     #ifndef EIGEN_PARSED_BY_DOXYGEN   224     template<
typename RhsType, 
typename DstType>
   237       dst = permutationP() * rhs;
   240       m_lu.template triangularView<UnitLower>().solveInPlace(dst);
   243       m_lu.template triangularView<Upper>().solveInPlace(dst);
   246     template<
bool Conjugate, 
typename RhsType, 
typename DstType>
   260         dst = m_lu.template triangularView<Upper>().adjoint().solve(rhs);
   262         m_lu.template triangularView<UnitLower>().adjoint().solveInPlace(dst);
   265         dst = m_lu.template triangularView<Upper>().transpose().solve(rhs);
   267         m_lu.template triangularView<UnitLower>().transpose().solveInPlace(dst);
   270       dst = permutationP().transpose() * dst;
   291 template<
typename MatrixType>
   295     m_rowsTranspositions(),
   298     m_isInitialized(false)
   302 template<
typename MatrixType>
   313 template<
typename MatrixType>
   314 template<
typename InputType>
   326 template<
typename MatrixType>
   327 template<
typename InputType>
   342 template<
typename Scalar, 
int StorageOrder, 
typename PivIndex>
   365   static Index unblocked_lu(MatrixType& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
   368     typedef typename Scoring::result_type Score;
   372     nb_transpositions = 0;
   373     Index first_zero_pivot = -1;
   376       Index rrows = rows-k-1;
   377       Index rcols = cols-k-1;
   379       Index row_of_biggest_in_col;
   380       Score biggest_in_corner
   381         = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
   382       row_of_biggest_in_col += k;
   384       row_transpositions[k] = PivIndex(row_of_biggest_in_col);
   386       if(biggest_in_corner != Score(0))
   388         if(k != row_of_biggest_in_col)
   390           lu.row(k).swap(lu.row(row_of_biggest_in_col));
   396         lu.col(k).tail(rrows) /= lu.coeff(k,k);
   398       else if(first_zero_pivot==-1)
   402         first_zero_pivot = k;
   406         lu.bottomRightCorner(rrows,rcols).noalias() -= lu.col(k).tail(rrows) * lu.row(k).tail(rcols);
   408     return first_zero_pivot;
   428     MapLU lu1(lu_data,StorageOrder==
RowMajor?rows:luStride,StorageOrder==
RowMajor?luStride:cols);
   429     MatrixType lu(lu1,0,0,rows,cols);
   436       return unblocked_lu(lu, row_transpositions, nb_transpositions);
   444       blockSize = (blockSize/16)*16;
   448     nb_transpositions = 0;
   449     Index first_zero_pivot = -1;
   453       Index trows = rows - k - bs; 
   454       Index tsize = size - k - bs; 
   460       BlockType A_0(lu,0,0,rows,k);
   461       BlockType A_2(lu,0,k+bs,rows,tsize);
   462       BlockType A11(lu,k,k,bs,bs);
   463       BlockType A12(lu,k,k+bs,bs,tsize);
   464       BlockType A21(lu,k+bs,k,trows,bs);
   465       BlockType A22(lu,k+bs,k+bs,trows,tsize);
   467       PivIndex nb_transpositions_in_panel;
   470       Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
   471                    row_transpositions+k, nb_transpositions_in_panel, 16);
   472       if(ret>=0 && first_zero_pivot==-1)
   473         first_zero_pivot = k+ret;
   475       nb_transpositions += nb_transpositions_in_panel;
   477       for(
Index i=k; i<k+bs; ++i)
   479         Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
   480         A_0.row(i).swap(A_0.row(piv));
   486         for(
Index i=k;i<k+bs; ++i)
   487           A_2.row(i).swap(A_2.row(row_transpositions[i]));
   490         A11.template triangularView<UnitLower>().solveInPlace(A12);
   492         A22.noalias() -= A21 * A12;
   495     return first_zero_pivot;
   501 template<
typename MatrixType, 
typename TranspositionType>
   505   eigen_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
   509     ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
   514 template<
typename MatrixType>
   524   eigen_assert(
m_lu.rows() == 
m_lu.cols() && 
"PartialPivLU is only for square (and moreover invertible) matrices");
   531   m_det_p = (nb_transpositions%2) ? -1 : 1;
   538 template<
typename MatrixType>
   548 template<
typename MatrixType>
   553   MatrixType res = 
m_lu.template triangularView<UnitLower>().toDenseMatrix()
   554                  * 
m_lu.template triangularView<Upper>();
   567 template<
typename DstXprType, 
typename MatrixType>
   587 template<
typename Derived>
   602 template<
typename Derived>
   611 #endif // EIGEN_PARTIALLU_H 
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
MatrixType::PlainObject PlainObject
PartialPivLU & compute(const EigenBase< InputType > &matrix)
PartialPivLU()
Default Constructor. 
const MatrixType & matrixLU() const
internal::traits< PartialPivLU< _MatrixType > >::Scalar Scalar
Block< MatrixType, Dynamic, Dynamic > BlockType
PartialPivLU< MatrixType > LuType
Map< Matrix< Scalar, Dynamic, Dynamic, StorageOrder > > MapLU
InverseReturnType inverse() const
A matrix or vector expression mapping an existing array of data. 
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< typename DstXprType::Scalar, typename LuType::Scalar > &)
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixType reconstructedMatrix() const
LU decomposition of a matrix with partial pivoting, and related features. 
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
MatrixType::RealScalar RealScalar
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator. 
Eigen::Index Index
The interface type of indices. 
const unsigned int RowMajorBit
EIGEN_DEVICE_FUNC Index cols() const
Expression of the inverse of another expression. 
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
void partial_lu_inplace(MatrixType &lu, TranspositionType &row_transpositions, typename TranspositionType::StorageIndex &nb_transpositions)
SolverBase< PartialPivLU > Base
static void check_template_parameters()
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API. 
Scalar determinant() const
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
const Inverse< PartialPivLU > inverse() const
void resize(Index newSize)
const PartialPivLU< PlainObject > partialPivLu() const
A matrix or vector expression mapping an existing expression. 
EIGEN_DEVICE_FUNC Index rows() const
IndicesType::Scalar StorageIndex
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationType
traits< _MatrixType > BaseTraits
Expression of a fixed-size or dynamic-size block. 
SolverStorage StorageKind
int64_t max(int64_t a, const int b)
static Index unblocked_lu(MatrixType &lu, PivIndex *row_transpositions, PivIndex &nb_transpositions)
const PermutationType & permutationP() const
Block< MapLU, Dynamic, Dynamic > MatrixType
TranspositionType m_rowsTranspositions
const PartialPivLU< PlainObject > lu() const
Inverse< LuType > SrcXprType
Pseudo expression representing a solving operation. 
EIGEN_DEVICE_FUNC Index size() const
A base class for matrix decomposition and solvers. 
static Index blocked_lu(Index rows, Index cols, Scalar *lu_data, Index luStride, PivIndex *row_transpositions, PivIndex &nb_transpositions, Index maxBlockSize=256)
EIGEN_DEVICE_FUNC Derived & derived()
Base class for all dense matrices, vectors, and expressions. 
EIGEN_DEVICE_FUNC const XprTypeNestedCleaned & nestedExpression() const
EIGEN_DEVICE_FUNC void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const