Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-19 09:06:43

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
0005 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
0006 //
0007 // This Source Code Form is subject to the terms of the Mozilla
0008 // Public License v. 2.0. If a copy of the MPL was not distributed
0009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0010 
0011 #ifndef SPARSELU_KERNEL_BMOD_H
0012 #define SPARSELU_KERNEL_BMOD_H
0013 
0014 namespace RivetEigen {
0015 namespace internal {
0016   
0017 template <int SegSizeAtCompileTime> struct LU_kernel_bmod
0018 {
0019   /** \internal
0020     * \brief Performs numeric block updates from a given supernode to a single column
0021     *
0022     * \param segsize Size of the segment (and blocks ) to use for updates
0023     * \param[in,out] dense Packed values of the original matrix
0024     * \param tempv temporary vector to use for updates
0025     * \param lusup array containing the supernodes
0026     * \param lda Leading dimension in the supernode
0027     * \param nrow Number of rows in the rectangular part of the supernode
0028     * \param lsub compressed row subscripts of supernodes
0029     * \param lptr pointer to the first column of the current supernode in lsub
0030     * \param no_zeros Number of nonzeros elements before the diagonal part of the supernode
0031     */
0032   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
0033   static EIGEN_DONT_INLINE void run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
0034                                     const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
0035 };
0036 
0037 template <int SegSizeAtCompileTime>
0038 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
0039 EIGEN_DONT_INLINE void LU_kernel_bmod<SegSizeAtCompileTime>::run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda,
0040                                                                   const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
0041 {
0042   typedef typename ScalarVector::Scalar Scalar;
0043   // First, copy U[*,j] segment from dense(*) to tempv(*)
0044   // The result of triangular solve is in tempv[*]; 
0045     // The result of matric-vector update is in dense[*]
0046   Index isub = lptr + no_zeros; 
0047   Index i;
0048   Index irow;
0049   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
0050   {
0051     irow = lsub(isub); 
0052     tempv(i) = dense(irow); 
0053     ++isub; 
0054   }
0055   // Dense triangular solve -- start effective triangle
0056   luptr += lda * no_zeros + no_zeros; 
0057   // Form Eigen matrix and vector 
0058   Map<Matrix<Scalar,SegSizeAtCompileTime,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(lda) );
0059   Map<Matrix<Scalar,SegSizeAtCompileTime,1> > u(tempv.data(), segsize);
0060   
0061   u = A.template triangularView<UnitLower>().solve(u); 
0062   
0063   // Dense matrix-vector product y <-- B*x 
0064   luptr += segsize;
0065   const Index PacketSize = internal::packet_traits<Scalar>::size;
0066   Index ldl = internal::first_multiple(nrow, PacketSize);
0067   Map<Matrix<Scalar,Dynamic,SegSizeAtCompileTime, ColMajor>, 0, OuterStride<> > B( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(lda) );
0068   Index aligned_offset = internal::first_default_aligned(tempv.data()+segsize, PacketSize);
0069   Index aligned_with_B_offset = (PacketSize-internal::first_default_aligned(B.data(), PacketSize))%PacketSize;
0070   Map<Matrix<Scalar,Dynamic,1>, 0, OuterStride<> > l(tempv.data()+segsize+aligned_offset+aligned_with_B_offset, nrow, OuterStride<>(ldl) );
0071   
0072   l.setZero();
0073   internal::sparselu_gemm<Scalar>(l.rows(), l.cols(), B.cols(), B.data(), B.outerStride(), u.data(), u.outerStride(), l.data(), l.outerStride());
0074   
0075   // Scatter tempv[] into SPA dense[] as a temporary storage 
0076   isub = lptr + no_zeros;
0077   for (i = 0; i < ((SegSizeAtCompileTime==Dynamic)?segsize:SegSizeAtCompileTime); i++)
0078   {
0079     irow = lsub(isub++); 
0080     dense(irow) = tempv(i);
0081   }
0082   
0083   // Scatter l into SPA dense[]
0084   for (i = 0; i < nrow; i++)
0085   {
0086     irow = lsub(isub++); 
0087     dense(irow) -= l(i);
0088   } 
0089 }
0090 
0091 template <> struct LU_kernel_bmod<1>
0092 {
0093   template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
0094   static EIGEN_DONT_INLINE void run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
0095                                     const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros);
0096 };
0097 
0098 
0099 template <typename BlockScalarVector, typename ScalarVector, typename IndexVector>
0100 EIGEN_DONT_INLINE void LU_kernel_bmod<1>::run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr,
0101                                               const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros)
0102 {
0103   typedef typename ScalarVector::Scalar Scalar;
0104   typedef typename IndexVector::Scalar StorageIndex;
0105   Scalar f = dense(lsub(lptr + no_zeros));
0106   luptr += lda * no_zeros + no_zeros + 1;
0107   const Scalar* a(lusup.data() + luptr);
0108   const StorageIndex*  irow(lsub.data()+lptr + no_zeros + 1);
0109   Index i = 0;
0110   for (; i+1 < nrow; i+=2)
0111   {
0112     Index i0 = *(irow++);
0113     Index i1 = *(irow++);
0114     Scalar a0 = *(a++);
0115     Scalar a1 = *(a++);
0116     Scalar d0 = dense.coeff(i0);
0117     Scalar d1 = dense.coeff(i1);
0118     d0 -= f*a0;
0119     d1 -= f*a1;
0120     dense.coeffRef(i0) = d0;
0121     dense.coeffRef(i1) = d1;
0122   }
0123   if(i<nrow)
0124     dense.coeffRef(*(irow++)) -= f * *(a++);
0125 }
0126 
0127 } // end namespace internal
0128 
0129 } // end namespace RivetEigen
0130 #endif // SPARSELU_KERNEL_BMOD_H