Warning, file /include/eigen3/Eigen/src/LU/FullPivLU.h was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #ifndef EIGEN_LU_H
0011 #define EIGEN_LU_H
0012
0013 namespace Eigen {
0014
0015 namespace internal {
0016 template<typename _MatrixType> struct traits<FullPivLU<_MatrixType> >
0017 : traits<_MatrixType>
0018 {
0019 typedef MatrixXpr XprKind;
0020 typedef SolverStorage StorageKind;
0021 typedef int StorageIndex;
0022 enum { Flags = 0 };
0023 };
0024
0025 }
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060 template<typename _MatrixType> class FullPivLU
0061 : public SolverBase<FullPivLU<_MatrixType> >
0062 {
0063 public:
0064 typedef _MatrixType MatrixType;
0065 typedef SolverBase<FullPivLU> Base;
0066 friend class SolverBase<FullPivLU>;
0067
0068 EIGEN_GENERIC_PUBLIC_INTERFACE(FullPivLU)
0069 enum {
0070 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0071 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
0072 };
0073 typedef typename internal::plain_row_type<MatrixType, StorageIndex>::type IntRowVectorType;
0074 typedef typename internal::plain_col_type<MatrixType, StorageIndex>::type IntColVectorType;
0075 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
0076 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
0077 typedef typename MatrixType::PlainObject PlainObject;
0078
0079
0080
0081
0082
0083
0084
0085 FullPivLU();
0086
0087
0088
0089
0090
0091
0092
0093 FullPivLU(Index rows, Index cols);
0094
0095
0096
0097
0098
0099
0100 template<typename InputType>
0101 explicit FullPivLU(const EigenBase<InputType>& matrix);
0102
0103
0104
0105
0106
0107
0108
0109 template<typename InputType>
0110 explicit FullPivLU(EigenBase<InputType>& matrix);
0111
0112
0113
0114
0115
0116
0117
0118
0119 template<typename InputType>
0120 FullPivLU& compute(const EigenBase<InputType>& matrix) {
0121 m_lu = matrix.derived();
0122 computeInPlace();
0123 return *this;
0124 }
0125
0126
0127
0128
0129
0130
0131
0132 inline const MatrixType& matrixLU() const
0133 {
0134 eigen_assert(m_isInitialized && "LU is not initialized.");
0135 return m_lu;
0136 }
0137
0138
0139
0140
0141
0142
0143
0144
0145 inline Index nonzeroPivots() const
0146 {
0147 eigen_assert(m_isInitialized && "LU is not initialized.");
0148 return m_nonzero_pivots;
0149 }
0150
0151
0152
0153
0154 RealScalar maxPivot() const { return m_maxpivot; }
0155
0156
0157
0158
0159
0160 EIGEN_DEVICE_FUNC inline const PermutationPType& permutationP() const
0161 {
0162 eigen_assert(m_isInitialized && "LU is not initialized.");
0163 return m_p;
0164 }
0165
0166
0167
0168
0169
0170 inline const PermutationQType& permutationQ() const
0171 {
0172 eigen_assert(m_isInitialized && "LU is not initialized.");
0173 return m_q;
0174 }
0175
0176
0177
0178
0179
0180
0181
0182
0183
0184
0185
0186
0187
0188
0189
0190 inline const internal::kernel_retval<FullPivLU> kernel() const
0191 {
0192 eigen_assert(m_isInitialized && "LU is not initialized.");
0193 return internal::kernel_retval<FullPivLU>(*this);
0194 }
0195
0196
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206
0207
0208
0209
0210
0211
0212
0213
0214
0215 inline const internal::image_retval<FullPivLU>
0216 image(const MatrixType& originalMatrix) const
0217 {
0218 eigen_assert(m_isInitialized && "LU is not initialized.");
0219 return internal::image_retval<FullPivLU>(*this, originalMatrix);
0220 }
0221
0222 #ifdef EIGEN_PARSED_BY_DOXYGEN
0223
0224
0225
0226
0227
0228
0229
0230
0231
0232
0233
0234
0235
0236
0237
0238
0239
0240
0241
0242 template<typename Rhs>
0243 inline const Solve<FullPivLU, Rhs>
0244 solve(const MatrixBase<Rhs>& b) const;
0245 #endif
0246
0247
0248
0249
0250 inline RealScalar rcond() const
0251 {
0252 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
0253 return internal::rcond_estimate_helper(m_l1_norm, *this);
0254 }
0255
0256
0257
0258
0259
0260
0261
0262
0263
0264
0265
0266
0267
0268
0269
0270
0271 typename internal::traits<MatrixType>::Scalar determinant() const;
0272
0273
0274
0275
0276
0277
0278
0279
0280
0281
0282
0283
0284
0285
0286
0287
0288
0289
0290 FullPivLU& setThreshold(const RealScalar& threshold)
0291 {
0292 m_usePrescribedThreshold = true;
0293 m_prescribedThreshold = threshold;
0294 return *this;
0295 }
0296
0297
0298
0299
0300
0301
0302
0303
0304
0305 FullPivLU& setThreshold(Default_t)
0306 {
0307 m_usePrescribedThreshold = false;
0308 return *this;
0309 }
0310
0311
0312
0313
0314
0315 RealScalar threshold() const
0316 {
0317 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
0318 return m_usePrescribedThreshold ? m_prescribedThreshold
0319
0320
0321 : NumTraits<Scalar>::epsilon() * RealScalar(m_lu.diagonalSize());
0322 }
0323
0324
0325
0326
0327
0328
0329
0330 inline Index rank() const
0331 {
0332 using std::abs;
0333 eigen_assert(m_isInitialized && "LU is not initialized.");
0334 RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
0335 Index result = 0;
0336 for(Index i = 0; i < m_nonzero_pivots; ++i)
0337 result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
0338 return result;
0339 }
0340
0341
0342
0343
0344
0345
0346
0347 inline Index dimensionOfKernel() const
0348 {
0349 eigen_assert(m_isInitialized && "LU is not initialized.");
0350 return cols() - rank();
0351 }
0352
0353
0354
0355
0356
0357
0358
0359
0360 inline bool isInjective() const
0361 {
0362 eigen_assert(m_isInitialized && "LU is not initialized.");
0363 return rank() == cols();
0364 }
0365
0366
0367
0368
0369
0370
0371
0372
0373 inline bool isSurjective() const
0374 {
0375 eigen_assert(m_isInitialized && "LU is not initialized.");
0376 return rank() == rows();
0377 }
0378
0379
0380
0381
0382
0383
0384
0385 inline bool isInvertible() const
0386 {
0387 eigen_assert(m_isInitialized && "LU is not initialized.");
0388 return isInjective() && (m_lu.rows() == m_lu.cols());
0389 }
0390
0391
0392
0393
0394
0395
0396
0397
0398 inline const Inverse<FullPivLU> inverse() const
0399 {
0400 eigen_assert(m_isInitialized && "LU is not initialized.");
0401 eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
0402 return Inverse<FullPivLU>(*this);
0403 }
0404
0405 MatrixType reconstructedMatrix() const;
0406
0407 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0408 inline Index rows() const EIGEN_NOEXCEPT { return m_lu.rows(); }
0409 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
0410 inline Index cols() const EIGEN_NOEXCEPT { return m_lu.cols(); }
0411
0412 #ifndef EIGEN_PARSED_BY_DOXYGEN
0413 template<typename RhsType, typename DstType>
0414 void _solve_impl(const RhsType &rhs, DstType &dst) const;
0415
0416 template<bool Conjugate, typename RhsType, typename DstType>
0417 void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
0418 #endif
0419
0420 protected:
0421
0422 static void check_template_parameters()
0423 {
0424 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
0425 }
0426
0427 void computeInPlace();
0428
0429 MatrixType m_lu;
0430 PermutationPType m_p;
0431 PermutationQType m_q;
0432 IntColVectorType m_rowsTranspositions;
0433 IntRowVectorType m_colsTranspositions;
0434 Index m_nonzero_pivots;
0435 RealScalar m_l1_norm;
0436 RealScalar m_maxpivot, m_prescribedThreshold;
0437 signed char m_det_pq;
0438 bool m_isInitialized, m_usePrescribedThreshold;
0439 };
0440
0441 template<typename MatrixType>
0442 FullPivLU<MatrixType>::FullPivLU()
0443 : m_isInitialized(false), m_usePrescribedThreshold(false)
0444 {
0445 }
0446
0447 template<typename MatrixType>
0448 FullPivLU<MatrixType>::FullPivLU(Index rows, Index cols)
0449 : m_lu(rows, cols),
0450 m_p(rows),
0451 m_q(cols),
0452 m_rowsTranspositions(rows),
0453 m_colsTranspositions(cols),
0454 m_isInitialized(false),
0455 m_usePrescribedThreshold(false)
0456 {
0457 }
0458
0459 template<typename MatrixType>
0460 template<typename InputType>
0461 FullPivLU<MatrixType>::FullPivLU(const EigenBase<InputType>& matrix)
0462 : m_lu(matrix.rows(), matrix.cols()),
0463 m_p(matrix.rows()),
0464 m_q(matrix.cols()),
0465 m_rowsTranspositions(matrix.rows()),
0466 m_colsTranspositions(matrix.cols()),
0467 m_isInitialized(false),
0468 m_usePrescribedThreshold(false)
0469 {
0470 compute(matrix.derived());
0471 }
0472
0473 template<typename MatrixType>
0474 template<typename InputType>
0475 FullPivLU<MatrixType>::FullPivLU(EigenBase<InputType>& matrix)
0476 : m_lu(matrix.derived()),
0477 m_p(matrix.rows()),
0478 m_q(matrix.cols()),
0479 m_rowsTranspositions(matrix.rows()),
0480 m_colsTranspositions(matrix.cols()),
0481 m_isInitialized(false),
0482 m_usePrescribedThreshold(false)
0483 {
0484 computeInPlace();
0485 }
0486
0487 template<typename MatrixType>
0488 void FullPivLU<MatrixType>::computeInPlace()
0489 {
0490 check_template_parameters();
0491
0492
0493 eigen_assert(m_lu.rows()<=NumTraits<int>::highest() && m_lu.cols()<=NumTraits<int>::highest());
0494
0495 m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
0496
0497 const Index size = m_lu.diagonalSize();
0498 const Index rows = m_lu.rows();
0499 const Index cols = m_lu.cols();
0500
0501
0502
0503 m_rowsTranspositions.resize(m_lu.rows());
0504 m_colsTranspositions.resize(m_lu.cols());
0505 Index number_of_transpositions = 0;
0506
0507 m_nonzero_pivots = size;
0508 m_maxpivot = RealScalar(0);
0509
0510 for(Index k = 0; k < size; ++k)
0511 {
0512
0513
0514
0515 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
0516 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
0517 typedef typename Scoring::result_type Score;
0518 Score biggest_in_corner;
0519 biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
0520 .unaryExpr(Scoring())
0521 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
0522 row_of_biggest_in_corner += k;
0523 col_of_biggest_in_corner += k;
0524
0525 if(biggest_in_corner==Score(0))
0526 {
0527
0528
0529 m_nonzero_pivots = k;
0530 for(Index i = k; i < size; ++i)
0531 {
0532 m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
0533 m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
0534 }
0535 break;
0536 }
0537
0538 RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
0539 if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
0540
0541
0542
0543
0544 m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
0545 m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
0546 if(k != row_of_biggest_in_corner) {
0547 m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
0548 ++number_of_transpositions;
0549 }
0550 if(k != col_of_biggest_in_corner) {
0551 m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
0552 ++number_of_transpositions;
0553 }
0554
0555
0556
0557
0558 if(k<rows-1)
0559 m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
0560 if(k<size-1)
0561 m_lu.block(k+1,k+1,rows-k-1,cols-k-1).noalias() -= m_lu.col(k).tail(rows-k-1) * m_lu.row(k).tail(cols-k-1);
0562 }
0563
0564
0565
0566
0567 m_p.setIdentity(rows);
0568 for(Index k = size-1; k >= 0; --k)
0569 m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
0570
0571 m_q.setIdentity(cols);
0572 for(Index k = 0; k < size; ++k)
0573 m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
0574
0575 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
0576
0577 m_isInitialized = true;
0578 }
0579
0580 template<typename MatrixType>
0581 typename internal::traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() const
0582 {
0583 eigen_assert(m_isInitialized && "LU is not initialized.");
0584 eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
0585 return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
0586 }
0587
0588
0589
0590
0591 template<typename MatrixType>
0592 MatrixType FullPivLU<MatrixType>::reconstructedMatrix() const
0593 {
0594 eigen_assert(m_isInitialized && "LU is not initialized.");
0595 const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
0596
0597 MatrixType res(m_lu.rows(),m_lu.cols());
0598
0599 res = m_lu.leftCols(smalldim)
0600 .template triangularView<UnitLower>().toDenseMatrix()
0601 * m_lu.topRows(smalldim)
0602 .template triangularView<Upper>().toDenseMatrix();
0603
0604
0605 res = m_p.inverse() * res;
0606
0607
0608 res = res * m_q.inverse();
0609
0610 return res;
0611 }
0612
0613
0614
0615 namespace internal {
0616 template<typename _MatrixType>
0617 struct kernel_retval<FullPivLU<_MatrixType> >
0618 : kernel_retval_base<FullPivLU<_MatrixType> >
0619 {
0620 EIGEN_MAKE_KERNEL_HELPERS(FullPivLU<_MatrixType>)
0621
0622 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
0623 MatrixType::MaxColsAtCompileTime,
0624 MatrixType::MaxRowsAtCompileTime)
0625 };
0626
0627 template<typename Dest> void evalTo(Dest& dst) const
0628 {
0629 using std::abs;
0630 const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
0631 if(dimker == 0)
0632 {
0633
0634
0635
0636 dst.setZero();
0637 return;
0638 }
0639
0640
0641
0642
0643
0644
0645
0646
0647
0648
0649
0650
0651
0652
0653
0654
0655
0656 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
0657 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
0658 Index p = 0;
0659 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
0660 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
0661 pivots.coeffRef(p++) = i;
0662 eigen_internal_assert(p == rank());
0663
0664
0665
0666
0667
0668 Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
0669 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
0670 m(dec().matrixLU().block(0, 0, rank(), cols));
0671 for(Index i = 0; i < rank(); ++i)
0672 {
0673 if(i) m.row(i).head(i).setZero();
0674 m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
0675 }
0676 m.block(0, 0, rank(), rank());
0677 m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
0678 for(Index i = 0; i < rank(); ++i)
0679 m.col(i).swap(m.col(pivots.coeff(i)));
0680
0681
0682
0683
0684 m.topLeftCorner(rank(), rank())
0685 .template triangularView<Upper>().solveInPlace(
0686 m.topRightCorner(rank(), dimker)
0687 );
0688
0689
0690 for(Index i = rank()-1; i >= 0; --i)
0691 m.col(i).swap(m.col(pivots.coeff(i)));
0692
0693
0694 for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
0695 for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
0696 for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
0697 }
0698 };
0699
0700
0701
0702 template<typename _MatrixType>
0703 struct image_retval<FullPivLU<_MatrixType> >
0704 : image_retval_base<FullPivLU<_MatrixType> >
0705 {
0706 EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
0707
0708 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
0709 MatrixType::MaxColsAtCompileTime,
0710 MatrixType::MaxRowsAtCompileTime)
0711 };
0712
0713 template<typename Dest> void evalTo(Dest& dst) const
0714 {
0715 using std::abs;
0716 if(rank() == 0)
0717 {
0718
0719
0720
0721 dst.setZero();
0722 return;
0723 }
0724
0725 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
0726 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
0727 Index p = 0;
0728 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
0729 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
0730 pivots.coeffRef(p++) = i;
0731 eigen_internal_assert(p == rank());
0732
0733 for(Index i = 0; i < rank(); ++i)
0734 dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
0735 }
0736 };
0737
0738
0739
0740 }
0741
0742 #ifndef EIGEN_PARSED_BY_DOXYGEN
0743 template<typename _MatrixType>
0744 template<typename RhsType, typename DstType>
0745 void FullPivLU<_MatrixType>::_solve_impl(const RhsType &rhs, DstType &dst) const
0746 {
0747
0748
0749
0750
0751
0752
0753
0754
0755 const Index rows = this->rows(),
0756 cols = this->cols(),
0757 nonzero_pivots = this->rank();
0758 const Index smalldim = (std::min)(rows, cols);
0759
0760 if(nonzero_pivots == 0)
0761 {
0762 dst.setZero();
0763 return;
0764 }
0765
0766 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
0767
0768
0769 c = permutationP() * rhs;
0770
0771
0772 m_lu.topLeftCorner(smalldim,smalldim)
0773 .template triangularView<UnitLower>()
0774 .solveInPlace(c.topRows(smalldim));
0775 if(rows>cols)
0776 c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
0777
0778
0779 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
0780 .template triangularView<Upper>()
0781 .solveInPlace(c.topRows(nonzero_pivots));
0782
0783
0784 for(Index i = 0; i < nonzero_pivots; ++i)
0785 dst.row(permutationQ().indices().coeff(i)) = c.row(i);
0786 for(Index i = nonzero_pivots; i < m_lu.cols(); ++i)
0787 dst.row(permutationQ().indices().coeff(i)).setZero();
0788 }
0789
0790 template<typename _MatrixType>
0791 template<bool Conjugate, typename RhsType, typename DstType>
0792 void FullPivLU<_MatrixType>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
0793 {
0794
0795
0796
0797
0798
0799
0800
0801
0802
0803
0804
0805 const Index rows = this->rows(), cols = this->cols(),
0806 nonzero_pivots = this->rank();
0807 const Index smalldim = (std::min)(rows, cols);
0808
0809 if(nonzero_pivots == 0)
0810 {
0811 dst.setZero();
0812 return;
0813 }
0814
0815 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
0816
0817
0818 c = permutationQ().inverse() * rhs;
0819
0820
0821 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
0822 .template triangularView<Upper>()
0823 .transpose()
0824 .template conjugateIf<Conjugate>()
0825 .solveInPlace(c.topRows(nonzero_pivots));
0826
0827
0828 m_lu.topLeftCorner(smalldim, smalldim)
0829 .template triangularView<UnitLower>()
0830 .transpose()
0831 .template conjugateIf<Conjugate>()
0832 .solveInPlace(c.topRows(smalldim));
0833
0834
0835 PermutationPType invp = permutationP().inverse().eval();
0836 for(Index i = 0; i < smalldim; ++i)
0837 dst.row(invp.indices().coeff(i)) = c.row(i);
0838 for(Index i = smalldim; i < rows; ++i)
0839 dst.row(invp.indices().coeff(i)).setZero();
0840 }
0841
0842 #endif
0843
0844 namespace internal {
0845
0846
0847
0848 template<typename DstXprType, typename MatrixType>
0849 struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename FullPivLU<MatrixType>::Scalar>, Dense2Dense>
0850 {
0851 typedef FullPivLU<MatrixType> LuType;
0852 typedef Inverse<LuType> SrcXprType;
0853 static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename MatrixType::Scalar> &)
0854 {
0855 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
0856 }
0857 };
0858 }
0859
0860
0861
0862
0863
0864
0865
0866
0867
0868 template<typename Derived>
0869 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
0870 MatrixBase<Derived>::fullPivLu() const
0871 {
0872 return FullPivLU<PlainObject>(eval());
0873 }
0874
0875 }
0876
0877 #endif