• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@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 template <int SegSizeAtCompileTime> struct LU_kernel_bmod
18 {
19   /** \internal
20     * \brief Performs numeric block updates from a given supernode to a single column
21     *
22     * \param segsize Size of the segment (and blocks ) to use for updates
23     * \param[in,out] dense Packed values of the original matrix
24     * \param tempv temporary vector to use for updates
25     * \param lusup array containing the supernodes
26     * \param lda Leading dimension in the supernode
27     * \param nrow Number of rows in the rectangular part of the supernode
28     * \param lsub compressed row subscripts of supernodes
29     * \param lptr pointer to the first column of the current supernode in lsub
30     * \param no_zeros Number of nonzeros elements before the diagonal part of the supernode
31     */
32   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
33   static EIGEN_DONT_INLINE void run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
34                                     const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
35 };
36 
37 template <int SegSizeAtCompileTime>
38 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
run(const Index segsize,BlockScalarVector & dense,ScalarVector & tempv,ScalarVector & lusup,Index & luptr,const Index lda,const Index nrow,IndexVector & lsub,const Index lptr,const Index no_zeros)39 EIGEN_DONT_INLINE void LU_kernel_bmod<SegSizeAtCompileTime>::run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
40                                                                   const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
41 {
42   typedef typename ScalarVector::Scalar Scalar;
43   // First, copy U[*,j] segment from dense(*) to tempv(*)
44   // The result of triangular solve is in tempv[*];
45     // The result of matric-vector update is in dense[*]
46   Index isub = lptr + no_zeros;
47   Index i;
48   Index irow;
49   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
50   {
51     irow = lsub(isub);
52     tempv(i) = dense(irow);
53     ++isub;
54   }
55   // Dense triangular solve -- start effective triangle
56   luptr += lda * no_zeros + no_zeros;
57   // Form Eigen matrix and vector
58   Map<Matrix<Scalar,SegSizeAtCompileTime,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(lda) );
59   Map<Matrix<Scalar,SegSizeAtCompileTime,1> > u(tempv.data(), segsize);
60 
61   u = A.template triangularView<UnitLower>().solve(u);
62 
63   // Dense matrix-vector product y <-- B*x
64   luptr += segsize;
65   const Index PacketSize = internal::packet_traits<Scalar>::size;
66   Index ldl = internal::first_multiple(nrow, PacketSize);
67   Map<Matrix<Scalar,Dynamic,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > B( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(lda) );
68   Index aligned_offset = internal::first_default_aligned(tempv.data()+segsize, PacketSize);
69   Index aligned_with_B_offset = (PacketSize-internal::first_default_aligned(B.data(), PacketSize))%PacketSize;
70   Map<Matrix<Scalar,Dynamic,1>, 0, OuterStride<> > l(tempv.data()+segsize+aligned_offset+aligned_with_B_offset, nrow, OuterStride<>(ldl) );
71 
72   l.setZero();
73   internal::sparselu_gemm<Scalar>(l.rows(), l.cols(), B.cols(), B.data(), B.outerStride(), u.data(), u.outerStride(), l.data(), l.outerStride());
74 
75   // Scatter tempv[] into SPA dense[] as a temporary storage
76   isub = lptr + no_zeros;
77   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
78   {
79     irow = lsub(isub++);
80     dense(irow) = tempv(i);
81   }
82 
83   // Scatter l into SPA dense[]
84   for (i = 0; i < nrow; i++)
85   {
86     irow = lsub(isub++);
87     dense(irow) -= l(i);
88   }
89 }
90 
91 template <> struct LU_kernel_bmod<1>
92 {
93   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
94   static EIGEN_DONT_INLINE void run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
95                                     const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
96 };
97 
98 
99 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
100 EIGEN_DONT_INLINE void LU_kernel_bmod<1>::run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
101                                               const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
102 {
103   typedef typename ScalarVector::Scalar Scalar;
104   typedef typename IndexVector::Scalar StorageIndex;
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 StorageIndex*  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