MOAB
4.9.3pre
|
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2006-2009 Benoit Jacob <[email protected]> 00005 // Copyright (C) 2009 Gael Guennebaud <[email protected]> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00010 00011 #ifndef EIGEN_PARTIALLU_H 00012 #define EIGEN_PARTIALLU_H 00013 00014 namespace Eigen { 00015 00016 namespace internal { 00017 template<typename _MatrixType> struct traits<PartialPivLU<_MatrixType> > 00018 : traits<_MatrixType> 00019 { 00020 typedef MatrixXpr XprKind; 00021 typedef SolverStorage StorageKind; 00022 typedef traits<_MatrixType> BaseTraits; 00023 enum { 00024 Flags = BaseTraits::Flags & RowMajorBit, 00025 CoeffReadCost = Dynamic 00026 }; 00027 }; 00028 00029 } // end namespace internal 00030 00062 template<typename _MatrixType> class PartialPivLU 00063 : public SolverBase<PartialPivLU<_MatrixType> > 00064 { 00065 public: 00066 00067 typedef _MatrixType MatrixType; 00068 typedef SolverBase<PartialPivLU> Base; 00069 EIGEN_GENERIC_PUBLIC_INTERFACE(PartialPivLU) 00070 // FIXME StorageIndex defined in EIGEN_GENERIC_PUBLIC_INTERFACE should be int 00071 enum { 00072 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, 00073 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime 00074 }; 00075 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType; 00076 typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType; 00077 typedef typename MatrixType::PlainObject PlainObject; 00078 00079 00086 PartialPivLU(); 00087 00094 explicit PartialPivLU(Index size); 00095 00103 template<typename InputType> 00104 explicit PartialPivLU(const EigenBase<InputType>& matrix); 00105 00106 template<typename InputType> 00107 PartialPivLU& compute(const EigenBase<InputType>& matrix); 00108 00115 inline const MatrixType& matrixLU() const 00116 { 00117 eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); 00118 return m_lu; 00119 } 00120 00123 inline const PermutationType& permutationP() const 00124 { 00125 eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); 00126 return m_p; 00127 } 00128 00146 // FIXME this is a copy-paste of the base-class member to add the isInitialized assertion. 00147 template<typename Rhs> 00148 inline const Solve<PartialPivLU, Rhs> 00149 solve(const MatrixBase<Rhs>& b) const 00150 { 00151 eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); 00152 return Solve<PartialPivLU, Rhs>(*this, b.derived()); 00153 } 00154 00162 inline const Inverse<PartialPivLU> inverse() const 00163 { 00164 eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); 00165 return Inverse<PartialPivLU>(*this); 00166 } 00167 00181 typename internal::traits<MatrixType>::Scalar determinant() const; 00182 00183 MatrixType reconstructedMatrix() const; 00184 00185 inline Index rows() const { return m_lu.rows(); } 00186 inline Index cols() const { return m_lu.cols(); } 00187 00188 #ifndef EIGEN_PARSED_BY_DOXYGEN 00189 template<typename RhsType, typename DstType> 00190 EIGEN_DEVICE_FUNC 00191 void _solve_impl(const RhsType &rhs, DstType &dst) const { 00192 /* The decomposition PA = LU can be rewritten as A = P^{-1} L U. 00193 * So we proceed as follows: 00194 * Step 1: compute c = Pb. 00195 * Step 2: replace c by the solution x to Lx = c. 00196 * Step 3: replace c by the solution x to Ux = c. 00197 */ 00198 00199 eigen_assert(rhs.rows() == m_lu.rows()); 00200 00201 // Step 1 00202 dst = permutationP() * rhs; 00203 00204 // Step 2 00205 m_lu.template triangularView<UnitLower>().solveInPlace(dst); 00206 00207 // Step 3 00208 m_lu.template triangularView<Upper>().solveInPlace(dst); 00209 } 00210 00211 template<bool Conjugate, typename RhsType, typename DstType> 00212 EIGEN_DEVICE_FUNC 00213 void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const { 00214 /* The decomposition PA = LU can be rewritten as A = P^{-1} L U. 00215 * So we proceed as follows: 00216 * Step 1: compute c = Pb. 00217 * Step 2: replace c by the solution x to Lx = c. 00218 * Step 3: replace c by the solution x to Ux = c. 00219 */ 00220 00221 eigen_assert(rhs.rows() == m_lu.cols()); 00222 00223 if (Conjugate) { 00224 // Step 1 00225 dst = m_lu.template triangularView<Upper>().adjoint().solve(rhs); 00226 // Step 2 00227 m_lu.template triangularView<UnitLower>().adjoint().solveInPlace(dst); 00228 } else { 00229 // Step 1 00230 dst = m_lu.template triangularView<Upper>().transpose().solve(rhs); 00231 // Step 2 00232 m_lu.template triangularView<UnitLower>().transpose().solveInPlace(dst); 00233 } 00234 // Step 3 00235 dst = permutationP().transpose() * dst; 00236 } 00237 #endif 00238 00239 protected: 00240 00241 static void check_template_parameters() 00242 { 00243 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar); 00244 } 00245 00246 MatrixType m_lu; 00247 PermutationType m_p; 00248 TranspositionType m_rowsTranspositions; 00249 Index m_det_p; 00250 bool m_isInitialized; 00251 }; 00252 00253 template<typename MatrixType> 00254 PartialPivLU<MatrixType>::PartialPivLU() 00255 : m_lu(), 00256 m_p(), 00257 m_rowsTranspositions(), 00258 m_det_p(0), 00259 m_isInitialized(false) 00260 { 00261 } 00262 00263 template<typename MatrixType> 00264 PartialPivLU<MatrixType>::PartialPivLU(Index size) 00265 : m_lu(size, size), 00266 m_p(size), 00267 m_rowsTranspositions(size), 00268 m_det_p(0), 00269 m_isInitialized(false) 00270 { 00271 } 00272 00273 template<typename MatrixType> 00274 template<typename InputType> 00275 PartialPivLU<MatrixType>::PartialPivLU(const EigenBase<InputType>& matrix) 00276 : m_lu(matrix.rows(), matrix.rows()), 00277 m_p(matrix.rows()), 00278 m_rowsTranspositions(matrix.rows()), 00279 m_det_p(0), 00280 m_isInitialized(false) 00281 { 00282 compute(matrix.derived()); 00283 } 00284 00285 namespace internal { 00286 00288 template<typename Scalar, int StorageOrder, typename PivIndex> 00289 struct partial_lu_impl 00290 { 00291 // FIXME add a stride to Map, so that the following mapping becomes easier, 00292 // another option would be to create an expression being able to automatically 00293 // warp any Map, Matrix, and Block expressions as a unique type, but since that's exactly 00294 // a Map + stride, why not adding a stride to Map, and convenient ctors from a Matrix, 00295 // and Block. 00296 typedef Map<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > MapLU; 00297 typedef Block<MapLU, Dynamic, Dynamic> MatrixType; 00298 typedef Block<MatrixType,Dynamic,Dynamic> BlockType; 00299 typedef typename MatrixType::RealScalar RealScalar; 00300 00311 static Index unblocked_lu(MatrixType& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions) 00312 { 00313 typedef scalar_score_coeff_op<Scalar> Scoring; 00314 typedef typename Scoring::result_type Score; 00315 const Index rows = lu.rows(); 00316 const Index cols = lu.cols(); 00317 const Index size = (std::min)(rows,cols); 00318 nb_transpositions = 0; 00319 Index first_zero_pivot = -1; 00320 for(Index k = 0; k < size; ++k) 00321 { 00322 Index rrows = rows-k-1; 00323 Index rcols = cols-k-1; 00324 00325 Index row_of_biggest_in_col; 00326 Score biggest_in_corner 00327 = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col); 00328 row_of_biggest_in_col += k; 00329 00330 row_transpositions[k] = PivIndex(row_of_biggest_in_col); 00331 00332 if(biggest_in_corner != Score(0)) 00333 { 00334 if(k != row_of_biggest_in_col) 00335 { 00336 lu.row(k).swap(lu.row(row_of_biggest_in_col)); 00337 ++nb_transpositions; 00338 } 00339 00340 // FIXME shall we introduce a safe quotient expression in cas 1/lu.coeff(k,k) 00341 // overflow but not the actual quotient? 00342 lu.col(k).tail(rrows) /= lu.coeff(k,k); 00343 } 00344 else if(first_zero_pivot==-1) 00345 { 00346 // the pivot is exactly zero, we record the index of the first pivot which is exactly 0, 00347 // and continue the factorization such we still have A = PLU 00348 first_zero_pivot = k; 00349 } 00350 00351 if(k<rows-1) 00352 lu.bottomRightCorner(rrows,rcols).noalias() -= lu.col(k).tail(rrows) * lu.row(k).tail(rcols); 00353 } 00354 return first_zero_pivot; 00355 } 00356 00372 static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256) 00373 { 00374 MapLU lu1(lu_data,StorageOrder==RowMajor?rows:luStride,StorageOrder==RowMajor?luStride:cols); 00375 MatrixType lu(lu1,0,0,rows,cols); 00376 00377 const Index size = (std::min)(rows,cols); 00378 00379 // if the matrix is too small, no blocking: 00380 if(size<=16) 00381 { 00382 return unblocked_lu(lu, row_transpositions, nb_transpositions); 00383 } 00384 00385 // automatically adjust the number of subdivisions to the size 00386 // of the matrix so that there is enough sub blocks: 00387 Index blockSize; 00388 { 00389 blockSize = size/8; 00390 blockSize = (blockSize/16)*16; 00391 blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize); 00392 } 00393 00394 nb_transpositions = 0; 00395 Index first_zero_pivot = -1; 00396 for(Index k = 0; k < size; k+=blockSize) 00397 { 00398 Index bs = (std::min)(size-k,blockSize); // actual size of the block 00399 Index trows = rows - k - bs; // trailing rows 00400 Index tsize = size - k - bs; // trailing size 00401 00402 // partition the matrix: 00403 // A00 | A01 | A02 00404 // lu = A_0 | A_1 | A_2 = A10 | A11 | A12 00405 // A20 | A21 | A22 00406 BlockType A_0(lu,0,0,rows,k); 00407 BlockType A_2(lu,0,k+bs,rows,tsize); 00408 BlockType A11(lu,k,k,bs,bs); 00409 BlockType A12(lu,k,k+bs,bs,tsize); 00410 BlockType A21(lu,k+bs,k,trows,bs); 00411 BlockType A22(lu,k+bs,k+bs,trows,tsize); 00412 00413 PivIndex nb_transpositions_in_panel; 00414 // recursively call the blocked LU algorithm on [A11^T A21^T]^T 00415 // with a very small blocking size: 00416 Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride, 00417 row_transpositions+k, nb_transpositions_in_panel, 16); 00418 if(ret>=0 && first_zero_pivot==-1) 00419 first_zero_pivot = k+ret; 00420 00421 nb_transpositions += nb_transpositions_in_panel; 00422 // update permutations and apply them to A_0 00423 for(Index i=k; i<k+bs; ++i) 00424 { 00425 Index piv = (row_transpositions[i] += k); 00426 A_0.row(i).swap(A_0.row(piv)); 00427 } 00428 00429 if(trows) 00430 { 00431 // apply permutations to A_2 00432 for(Index i=k;i<k+bs; ++i) 00433 A_2.row(i).swap(A_2.row(row_transpositions[i])); 00434 00435 // A12 = A11^-1 A12 00436 A11.template triangularView<UnitLower>().solveInPlace(A12); 00437 00438 A22.noalias() -= A21 * A12; 00439 } 00440 } 00441 return first_zero_pivot; 00442 } 00443 }; 00444 00447 template<typename MatrixType, typename TranspositionType> 00448 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename TranspositionType::StorageIndex& nb_transpositions) 00449 { 00450 eigen_assert(lu.cols() == row_transpositions.size()); 00451 eigen_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1); 00452 00453 partial_lu_impl 00454 <typename MatrixType::Scalar, MatrixType::Flags&RowMajorBit?RowMajor:ColMajor, typename TranspositionType::StorageIndex> 00455 ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions); 00456 } 00457 00458 } // end namespace internal 00459 00460 template<typename MatrixType> 00461 template<typename InputType> 00462 PartialPivLU<MatrixType>& PartialPivLU<MatrixType>::compute(const EigenBase<InputType>& matrix) 00463 { 00464 check_template_parameters(); 00465 00466 // the row permutation is stored as int indices, so just to be sure: 00467 eigen_assert(matrix.rows()<NumTraits<int>::highest()); 00468 00469 m_lu = matrix.derived(); 00470 00471 eigen_assert(matrix.rows() == matrix.cols() && "PartialPivLU is only for square (and moreover invertible) matrices"); 00472 const Index size = matrix.rows(); 00473 00474 m_rowsTranspositions.resize(size); 00475 00476 typename TranspositionType::StorageIndex nb_transpositions; 00477 internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions); 00478 m_det_p = (nb_transpositions%2) ? -1 : 1; 00479 00480 m_p = m_rowsTranspositions; 00481 00482 m_isInitialized = true; 00483 return *this; 00484 } 00485 00486 template<typename MatrixType> 00487 typename internal::traits<MatrixType>::Scalar PartialPivLU<MatrixType>::determinant() const 00488 { 00489 eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); 00490 return Scalar(m_det_p) * m_lu.diagonal().prod(); 00491 } 00492 00496 template<typename MatrixType> 00497 MatrixType PartialPivLU<MatrixType>::reconstructedMatrix() const 00498 { 00499 eigen_assert(m_isInitialized && "LU is not initialized."); 00500 // LU 00501 MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix() 00502 * m_lu.template triangularView<Upper>(); 00503 00504 // P^{-1}(LU) 00505 res = m_p.inverse() * res; 00506 00507 return res; 00508 } 00509 00510 /***** Implementation details *****************************************************/ 00511 00512 namespace internal { 00513 00514 /***** Implementation of inverse() *****************************************************/ 00515 template<typename DstXprType, typename MatrixType, typename Scalar> 00516 struct Assignment<DstXprType, Inverse<PartialPivLU<MatrixType> >, internal::assign_op<Scalar>, Dense2Dense, Scalar> 00517 { 00518 typedef PartialPivLU<MatrixType> LuType; 00519 typedef Inverse<LuType> SrcXprType; 00520 static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar> &) 00521 { 00522 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols())); 00523 } 00524 }; 00525 } // end namespace internal 00526 00527 /******** MatrixBase methods *******/ 00528 00535 #ifndef __CUDACC__ 00536 template<typename Derived> 00537 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject> 00538 MatrixBase<Derived>::partialPivLu() const 00539 { 00540 return PartialPivLU<PlainObject>(eval()); 00541 } 00542 #endif 00543 00552 #ifndef __CUDACC__ 00553 template<typename Derived> 00554 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject> 00555 MatrixBase<Derived>::lu() const 00556 { 00557 return PartialPivLU<PlainObject>(eval()); 00558 } 00559 #endif 00560 00561 } // end namespace Eigen 00562 00563 #endif // EIGEN_PARTIALLU_H