1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2012 Dsir Nuentsa-Wakam <desire.nuentsa_wakam (at) inria.fr> 5 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud (at) inria.fr> 6 // 7 // This Source Code Form is subject to the terms of the Mozilla 8 // Public License v. 2.0. If a copy of the MPL was not distributed 9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 10 11 #ifndef SPARSELU_KERNEL_BMOD_H 12 #define SPARSELU_KERNEL_BMOD_H 13 14 namespace Eigen { 15 namespace internal { 16 17 /** 18 * \brief Performs numeric block updates from a given supernode to a single column 19 * 20 * \param segsize Size of the segment (and blocks ) to use for updates 21 * \param[in,out] dense Packed values of the original matrix 22 * \param tempv temporary vector to use for updates 23 * \param lusup array containing the supernodes 24 * \param lda Leading dimension in the supernode 25 * \param nrow Number of rows in the rectangular part of the supernode 26 * \param lsub compressed row subscripts of supernodes 27 * \param lptr pointer to the first column of the current supernode in lsub 28 * \param no_zeros Number of nonzeros elements before the diagonal part of the supernode 29 * \return 0 on success 30 */ 31 template <int SegSizeAtCompileTime> struct LU_kernel_bmod 32 { 33 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector, typename Index> 34 static EIGEN_DONT_INLINE void run(const int segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda, 35 const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros); 36 }; 37 38 template <int SegSizeAtCompileTime> 39 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector, typename Index> 40 EIGEN_DONT_INLINE void LU_kernel_bmod<SegSizeAtCompileTime>::run(const int segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda, 41 const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros) 42 { 43 typedef typename ScalarVector::Scalar Scalar; 44 // First, copy U[*,j] segment from dense(*) to tempv(*) 45 // The result of triangular solve is in tempv[*]; 46 // The result of matric-vector update is in dense[*] 47 Index isub = lptr + no_zeros; 48 int i; 49 Index irow; 50 for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++) 51 { 52 irow = lsub(isub); 53 tempv(i) = dense(irow); 54 ++isub; 55 } 56 // Dense triangular solve -- start effective triangle 57 luptr += lda * no_zeros + no_zeros; 58 // Form Eigen matrix and vector 59 Map<Matrix<Scalar,SegSizeAtCompileTime,SegSizeAtCompileTime>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(lda) ); 60 Map<Matrix<Scalar,SegSizeAtCompileTime,1> > u(tempv.data(), segsize); 61 62 u = A.template triangularView<UnitLower>().solve(u); 63 64 // Dense matrix-vector product y <-- B*x 65 luptr += segsize; 66 const Index PacketSize = internal::packet_traits<Scalar>::size; 67 Index ldl = internal::first_multiple(nrow, PacketSize); 68 Map<Matrix<Scalar,Dynamic,SegSizeAtCompileTime>, 0, OuterStride<> > B( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(lda) ); 69 Index aligned_offset = internal::first_aligned(tempv.data()+segsize, PacketSize); 70 Index aligned_with_B_offset = (PacketSize-internal::first_aligned(B.data(), PacketSize))%PacketSize; 71 Map<Matrix<Scalar,Dynamic,1>, 0, OuterStride<> > l(tempv.data()+segsize+aligned_offset+aligned_with_B_offset, nrow, OuterStride<>(ldl) ); 72 73 l.setZero(); 74 internal::sparselu_gemm<Scalar>(l.rows(), l.cols(), B.cols(), B.data(), B.outerStride(), u.data(), u.outerStride(), l.data(), l.outerStride()); 75 76 // Scatter tempv[] into SPA dense[] as a temporary storage 77 isub = lptr + no_zeros; 78 for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++) 79 { 80 irow = lsub(isub++); 81 dense(irow) = tempv(i); 82 } 83 84 // Scatter l into SPA dense[] 85 for (i = 0; i < nrow; i++) 86 { 87 irow = lsub(isub++); 88 dense(irow) -= l(i); 89 } 90 } 91 92 template <> struct LU_kernel_bmod<1> 93 { 94 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector, typename Index> 95 static EIGEN_DONT_INLINE void run(const int /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr, 96 const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros); 97 }; 98 99 100 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector, typename Index> 101 EIGEN_DONT_INLINE void LU_kernel_bmod<1>::run(const int /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr, 102 const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros) 103 { 104 typedef typename ScalarVector::Scalar Scalar; 105 Scalar f = dense(lsub(lptr + no_zeros)); 106 luptr += lda * no_zeros + no_zeros + 1; 107 const Scalar* a(lusup.data() + luptr); 108 const /*typename IndexVector::Scalar*/Index* irow(lsub.data()+lptr + no_zeros + 1); 109 Index i = 0; 110 for (; i+1 < nrow; i+=2) 111 { 112 Index i0 = *(irow++); 113 Index i1 = *(irow++); 114 Scalar a0 = *(a++); 115 Scalar a1 = *(a++); 116 Scalar d0 = dense.coeff(i0); 117 Scalar d1 = dense.coeff(i1); 118 d0 -= f*a0; 119 d1 -= f*a1; 120 dense.coeffRef(i0) = d0; 121 dense.coeffRef(i1) = d1; 122 } 123 if(i<nrow) 124 dense.coeffRef(*(irow++)) -= f * *(a++); 125 } 126 127 } // end namespace internal 128 129 } // end namespace Eigen 130 #endif // SPARSELU_KERNEL_BMOD_H 131