Eigen  3.3.3
SparseLU_kernel_bmod.h
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
00005 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
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 SPARSELU_KERNEL_BMOD_H
00012 #define SPARSELU_KERNEL_BMOD_H
00013 
00014 namespace Eigen {
00015 namespace internal {
00016   
00017 template <int SegSizeAtCompileTime> struct LU_kernel_bmod
00018 {
00032   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
00033   static EIGEN_DONT_INLINE void run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
00034                                     const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
00035 };
00036 
00037 template <int SegSizeAtCompileTime>
00038 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
00039 EIGEN_DONT_INLINE void LU_kernel_bmod<SegSizeAtCompileTime>::run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
00040                                                                   const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
00041 {
00042   typedef typename ScalarVector::Scalar Scalar;
00043   // First, copy U[*,j] segment from dense(*) to tempv(*)
00044   // The result of triangular solve is in tempv[*]; 
00045     // The result of matric-vector update is in dense[*]
00046   Index isub = lptr + no_zeros; 
00047   Index i;
00048   Index irow;
00049   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
00050   {
00051     irow = lsub(isub); 
00052     tempv(i) = dense(irow); 
00053     ++isub; 
00054   }
00055   // Dense triangular solve -- start effective triangle
00056   luptr += lda * no_zeros + no_zeros; 
00057   // Form Eigen matrix and vector 
00058   Map<Matrix<Scalar,SegSizeAtCompileTime,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(lda) );
00059   Map<Matrix<Scalar,SegSizeAtCompileTime,1> > u(tempv.data(), segsize);
00060   
00061   u = A.template triangularView<UnitLower>().solve(u); 
00062   
00063   // Dense matrix-vector product y <-- B*x 
00064   luptr += segsize;
00065   const Index PacketSize = internal::packet_traits<Scalar>::size;
00066   Index ldl = internal::first_multiple(nrow, PacketSize);
00067   Map<Matrix<Scalar,Dynamic,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > B( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(lda) );
00068   Index aligned_offset = internal::first_default_aligned(tempv.data()+segsize, PacketSize);
00069   Index aligned_with_B_offset = (PacketSize-internal::first_default_aligned(B.data(), PacketSize))%PacketSize;
00070   Map<Matrix<Scalar,Dynamic,1>, 0, OuterStride<> > l(tempv.data()+segsize+aligned_offset+aligned_with_B_offset, nrow, OuterStride<>(ldl) );
00071   
00072   l.setZero();
00073   internal::sparselu_gemm<Scalar>(l.rows(), l.cols(), B.cols(), B.data(), B.outerStride(), u.data(), u.outerStride(), l.data(), l.outerStride());
00074   
00075   // Scatter tempv[] into SPA dense[] as a temporary storage 
00076   isub = lptr + no_zeros;
00077   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
00078   {
00079     irow = lsub(isub++); 
00080     dense(irow) = tempv(i);
00081   }
00082   
00083   // Scatter l into SPA dense[]
00084   for (i = 0; i < nrow; i++)
00085   {
00086     irow = lsub(isub++); 
00087     dense(irow) -= l(i);
00088   } 
00089 }
00090 
00091 template <> struct LU_kernel_bmod<1>
00092 {
00093   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
00094   static EIGEN_DONT_INLINE void run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
00095                                     const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
00096 };
00097 
00098 
00099 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
00100 EIGEN_DONT_INLINE void LU_kernel_bmod<1>::run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
00101                                               const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
00102 {
00103   typedef typename ScalarVector::Scalar Scalar;
00104   typedef typename IndexVector::Scalar StorageIndex;
00105   Scalar f = dense(lsub(lptr + no_zeros));
00106   luptr += lda * no_zeros + no_zeros + 1;
00107   const Scalar* a(lusup.data() + luptr);
00108   const StorageIndex*  irow(lsub.data()+lptr + no_zeros + 1);
00109   Index i = 0;
00110   for (; i+1 < nrow; i+=2)
00111   {
00112     Index i0 = *(irow++);
00113     Index i1 = *(irow++);
00114     Scalar a0 = *(a++);
00115     Scalar a1 = *(a++);
00116     Scalar d0 = dense.coeff(i0);
00117     Scalar d1 = dense.coeff(i1);
00118     d0 -= f*a0;
00119     d1 -= f*a1;
00120     dense.coeffRef(i0) = d0;
00121     dense.coeffRef(i1) = d1;
00122   }
00123   if(i<nrow)
00124     dense.coeffRef(*(irow++)) -= f * *(a++);
00125 }
00126 
00127 } // end namespace internal
00128 
00129 } // end namespace Eigen
00130 #endif // SPARSELU_KERNEL_BMOD_H
 All Classes Functions Variables Typedefs Enumerations Enumerator Friends