TriangularSolverMatrix.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_TRIANGULAR_SOLVER_MATRIX_H
00026 #define EIGEN_TRIANGULAR_SOLVER_MATRIX_H
00027 
00028 namespace Eigen { 
00029 
00030 namespace internal {
00031 
00032 // if the rhs is row major, let's transpose the product
00033 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder>
00034 struct triangular_solve_matrix<Scalar,Index,Side,Mode,Conjugate,TriStorageOrder,RowMajor>
00035 {
00036   static EIGEN_DONT_INLINE void run(
00037     Index size, Index cols,
00038     const Scalar*  tri, Index triStride,
00039     Scalar* _other, Index otherStride)
00040   {
00041     triangular_solve_matrix<
00042       Scalar, Index, Side==OnTheLeft?OnTheRight:OnTheLeft,
00043       (Mode&UnitDiag) | ((Mode&Upper) ? Lower : Upper),
00044       NumTraits<Scalar>::IsComplex && Conjugate,
00045       TriStorageOrder==RowMajor ? ColMajor : RowMajor, ColMajor>
00046       ::run(size, cols, tri, triStride, _other, otherStride);
00047   }
00048 };
00049 
00050 /* Optimized triangular solver with multiple right hand side and the triangular matrix on the left
00051  */
00052 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
00053 struct triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conjugate,TriStorageOrder,ColMajor>
00054 {
00055   static EIGEN_DONT_INLINE void run(
00056     Index size, Index otherSize,
00057     const Scalar* _tri, Index triStride,
00058     Scalar* _other, Index otherStride)
00059   {
00060     Index cols = otherSize;
00061     const_blas_data_mapper<Scalar, Index, TriStorageOrder> tri(_tri,triStride);
00062     blas_data_mapper<Scalar, Index, ColMajor> other(_other,otherStride);
00063 
00064     typedef gebp_traits<Scalar,Scalar> Traits;
00065     enum {
00066       SmallPanelWidth   = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
00067       IsLower = (Mode&Lower) == Lower
00068     };
00069 
00070     Index kc = size; // cache block size along the K direction
00071     Index mc = size;  // cache block size along the M direction
00072     Index nc = cols;  // cache block size along the N direction
00073     computeProductBlockingSizes<Scalar,Scalar,4>(kc, mc, nc);
00074 
00075     std::size_t sizeW = kc*Traits::WorkSpaceFactor;
00076     std::size_t sizeB = sizeW + kc*cols;
00077     ei_declare_aligned_stack_constructed_variable(Scalar, blockA, kc*mc, 0);
00078     ei_declare_aligned_stack_constructed_variable(Scalar, allocatedBlockB, sizeB, 0);
00079     Scalar* blockB = allocatedBlockB + sizeW;
00080     Scalar* blockW = allocatedBlockB;
00081 
00082     conj_if<Conjugate> conj;
00083     gebp_kernel<Scalar, Scalar, Index, Traits::mr, Traits::nr, Conjugate, false> gebp_kernel;
00084     gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, TriStorageOrder> pack_lhs;
00085     gemm_pack_rhs<Scalar, Index, Traits::nr, ColMajor, false, true> pack_rhs;
00086 
00087     // the goal here is to subdivise the Rhs panels such that we keep some cache
00088     // coherence when accessing the rhs elements
00089     std::ptrdiff_t l1, l2;
00090     manage_caching_sizes(GetAction, &l1, &l2);
00091     Index subcols = cols>0 ? l2/(4 * sizeof(Scalar) * otherStride) : 0;
00092     subcols = std::max<Index>((subcols/Traits::nr)*Traits::nr, Traits::nr);
00093 
00094     for(Index k2=IsLower ? 0 : size;
00095         IsLower ? k2<size : k2>0;
00096         IsLower ? k2+=kc : k2-=kc)
00097     {
00098       const Index actual_kc = (std::min)(IsLower ? size-k2 : k2, kc);
00099 
00100       // We have selected and packed a big horizontal panel R1 of rhs. Let B be the packed copy of this panel,
00101       // and R2 the remaining part of rhs. The corresponding vertical panel of lhs is split into
00102       // A11 (the triangular part) and A21 the remaining rectangular part.
00103       // Then the high level algorithm is:
00104       //  - B = R1                    => general block copy (done during the next step)
00105       //  - R1 = A11^-1 B             => tricky part
00106       //  - update B from the new R1  => actually this has to be performed continuously during the above step
00107       //  - R2 -= A21 * B             => GEPP
00108 
00109       // The tricky part: compute R1 = A11^-1 B while updating B from R1
00110       // The idea is to split A11 into multiple small vertical panels.
00111       // Each panel can be split into a small triangular part T1k which is processed without optimization,
00112       // and the remaining small part T2k which is processed using gebp with appropriate block strides
00113       for(Index j2=0; j2<cols; j2+=subcols)
00114       {
00115         Index actual_cols = (std::min)(cols-j2,subcols);
00116         // for each small vertical panels [T1k^T, T2k^T]^T of lhs
00117         for (Index k1=0; k1<actual_kc; k1+=SmallPanelWidth)
00118         {
00119           Index actualPanelWidth = std::min<Index>(actual_kc-k1, SmallPanelWidth);
00120           // tr solve
00121           for (Index k=0; k<actualPanelWidth; ++k)
00122           {
00123             // TODO write a small kernel handling this (can be shared with trsv)
00124             Index i  = IsLower ? k2+k1+k : k2-k1-k-1;
00125             Index s  = IsLower ? k2+k1 : i+1;
00126             Index rs = actualPanelWidth - k - 1; // remaining size
00127 
00128             Scalar a = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(tri(i,i));
00129             for (Index j=j2; j<j2+actual_cols; ++j)
00130             {
00131               if (TriStorageOrder==RowMajor)
00132               {
00133                 Scalar b(0);
00134                 const Scalar* l = &tri(i,s);
00135                 Scalar* r = &other(s,j);
00136                 for (Index i3=0; i3<k; ++i3)
00137                   b += conj(l[i3]) * r[i3];
00138 
00139                 other(i,j) = (other(i,j) - b)*a;
00140               }
00141               else
00142               {
00143                 Index s = IsLower ? i+1 : i-rs;
00144                 Scalar b = (other(i,j) *= a);
00145                 Scalar* r = &other(s,j);
00146                 const Scalar* l = &tri(s,i);
00147                 for (Index i3=0;i3<rs;++i3)
00148                   r[i3] -= b * conj(l[i3]);
00149               }
00150             }
00151           }
00152 
00153           Index lengthTarget = actual_kc-k1-actualPanelWidth;
00154           Index startBlock   = IsLower ? k2+k1 : k2-k1-actualPanelWidth;
00155           Index blockBOffset = IsLower ? k1 : lengthTarget;
00156 
00157           // update the respective rows of B from other
00158           pack_rhs(blockB+actual_kc*j2, &other(startBlock,j2), otherStride, actualPanelWidth, actual_cols, actual_kc, blockBOffset);
00159 
00160           // GEBP
00161           if (lengthTarget>0)
00162           {
00163             Index startTarget  = IsLower ? k2+k1+actualPanelWidth : k2-actual_kc;
00164 
00165             pack_lhs(blockA, &tri(startTarget,startBlock), triStride, actualPanelWidth, lengthTarget);
00166 
00167             gebp_kernel(&other(startTarget,j2), otherStride, blockA, blockB+actual_kc*j2, lengthTarget, actualPanelWidth, actual_cols, Scalar(-1),
00168                         actualPanelWidth, actual_kc, 0, blockBOffset, blockW);
00169           }
00170         }
00171       }
00172       
00173       // R2 -= A21 * B => GEPP
00174       {
00175         Index start = IsLower ? k2+kc : 0;
00176         Index end   = IsLower ? size : k2-kc;
00177         for(Index i2=start; i2<end; i2+=mc)
00178         {
00179           const Index actual_mc = (std::min)(mc,end-i2);
00180           if (actual_mc>0)
00181           {
00182             pack_lhs(blockA, &tri(i2, IsLower ? k2 : k2-kc), triStride, actual_kc, actual_mc);
00183 
00184             gebp_kernel(_other+i2, otherStride, blockA, blockB, actual_mc, actual_kc, cols, Scalar(-1));
00185           }
00186         }
00187       }
00188     }
00189   }
00190 };
00191 
00192 /* Optimized triangular solver with multiple left hand sides and the trinagular matrix on the right
00193  */
00194 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
00195 struct triangular_solve_matrix<Scalar,Index,OnTheRight,Mode,Conjugate,TriStorageOrder,ColMajor>
00196 {
00197   static EIGEN_DONT_INLINE void run(
00198     Index size, Index otherSize,
00199     const Scalar* _tri, Index triStride,
00200     Scalar* _other, Index otherStride)
00201   {
00202     Index rows = otherSize;
00203     const_blas_data_mapper<Scalar, Index, TriStorageOrder> rhs(_tri,triStride);
00204     blas_data_mapper<Scalar, Index, ColMajor> lhs(_other,otherStride);
00205 
00206     typedef gebp_traits<Scalar,Scalar> Traits;
00207     enum {
00208       RhsStorageOrder   = TriStorageOrder,
00209       SmallPanelWidth   = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
00210       IsLower = (Mode&Lower) == Lower
00211     };
00212 
00213 //     Index kc = std::min<Index>(Traits::Max_kc/4,size); // cache block size along the K direction
00214 //     Index mc = std::min<Index>(Traits::Max_mc,size);   // cache block size along the M direction
00215     // check that !!!!
00216     Index kc = size; // cache block size along the K direction
00217     Index mc = size;  // cache block size along the M direction
00218     Index nc = rows;  // cache block size along the N direction
00219     computeProductBlockingSizes<Scalar,Scalar,4>(kc, mc, nc);
00220 
00221     std::size_t sizeW = kc*Traits::WorkSpaceFactor;
00222     std::size_t sizeB = sizeW + kc*size;
00223     ei_declare_aligned_stack_constructed_variable(Scalar, blockA, kc*mc, 0);
00224     ei_declare_aligned_stack_constructed_variable(Scalar, allocatedBlockB, sizeB, 0);
00225     Scalar* blockB = allocatedBlockB + sizeW;
00226 
00227     conj_if<Conjugate> conj;
00228     gebp_kernel<Scalar,Scalar, Index, Traits::mr, Traits::nr, false, Conjugate> gebp_kernel;
00229     gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder> pack_rhs;
00230     gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder,false,true> pack_rhs_panel;
00231     gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, ColMajor, false, true> pack_lhs_panel;
00232 
00233     for(Index k2=IsLower ? size : 0;
00234         IsLower ? k2>0 : k2<size;
00235         IsLower ? k2-=kc : k2+=kc)
00236     {
00237       const Index actual_kc = (std::min)(IsLower ? k2 : size-k2, kc);
00238       Index actual_k2 = IsLower ? k2-actual_kc : k2 ;
00239 
00240       Index startPanel = IsLower ? 0 : k2+actual_kc;
00241       Index rs = IsLower ? actual_k2 : size - actual_k2 - actual_kc;
00242       Scalar* geb = blockB+actual_kc*actual_kc;
00243 
00244       if (rs>0) pack_rhs(geb, &rhs(actual_k2,startPanel), triStride, actual_kc, rs);
00245 
00246       // triangular packing (we only pack the panels off the diagonal,
00247       // neglecting the blocks overlapping the diagonal
00248       {
00249         for (Index j2=0; j2<actual_kc; j2+=SmallPanelWidth)
00250         {
00251           Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
00252           Index actual_j2 = actual_k2 + j2;
00253           Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
00254           Index panelLength = IsLower ? actual_kc-j2-actualPanelWidth : j2;
00255 
00256           if (panelLength>0)
00257           pack_rhs_panel(blockB+j2*actual_kc,
00258                          &rhs(actual_k2+panelOffset, actual_j2), triStride,
00259                          panelLength, actualPanelWidth,
00260                          actual_kc, panelOffset);
00261         }
00262       }
00263 
00264       for(Index i2=0; i2<rows; i2+=mc)
00265       {
00266         const Index actual_mc = (std::min)(mc,rows-i2);
00267 
00268         // triangular solver kernel
00269         {
00270           // for each small block of the diagonal (=> vertical panels of rhs)
00271           for (Index j2 = IsLower
00272                       ? (actual_kc - ((actual_kc%SmallPanelWidth) ? Index(actual_kc%SmallPanelWidth)
00273                                                                   : Index(SmallPanelWidth)))
00274                       : 0;
00275                IsLower ? j2>=0 : j2<actual_kc;
00276                IsLower ? j2-=SmallPanelWidth : j2+=SmallPanelWidth)
00277           {
00278             Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
00279             Index absolute_j2 = actual_k2 + j2;
00280             Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
00281             Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;
00282 
00283             // GEBP
00284             if(panelLength>0)
00285             {
00286               gebp_kernel(&lhs(i2,absolute_j2), otherStride,
00287                           blockA, blockB+j2*actual_kc,
00288                           actual_mc, panelLength, actualPanelWidth,
00289                           Scalar(-1),
00290                           actual_kc, actual_kc, // strides
00291                           panelOffset, panelOffset, // offsets
00292                           allocatedBlockB);  // workspace
00293             }
00294 
00295             // unblocked triangular solve
00296             for (Index k=0; k<actualPanelWidth; ++k)
00297             {
00298               Index j = IsLower ? absolute_j2+actualPanelWidth-k-1 : absolute_j2+k;
00299 
00300               Scalar* r = &lhs(i2,j);
00301               for (Index k3=0; k3<k; ++k3)
00302               {
00303                 Scalar b = conj(rhs(IsLower ? j+1+k3 : absolute_j2+k3,j));
00304                 Scalar* a = &lhs(i2,IsLower ? j+1+k3 : absolute_j2+k3);
00305                 for (Index i=0; i<actual_mc; ++i)
00306                   r[i] -= a[i] * b;
00307               }
00308               Scalar b = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(rhs(j,j));
00309               for (Index i=0; i<actual_mc; ++i)
00310                 r[i] *= b;
00311             }
00312 
00313             // pack the just computed part of lhs to A
00314             pack_lhs_panel(blockA, _other+absolute_j2*otherStride+i2, otherStride,
00315                            actualPanelWidth, actual_mc,
00316                            actual_kc, j2);
00317           }
00318         }
00319 
00320         if (rs>0)
00321           gebp_kernel(_other+i2+startPanel*otherStride, otherStride, blockA, geb,
00322                       actual_mc, actual_kc, rs, Scalar(-1),
00323                       -1, -1, 0, 0, allocatedBlockB);
00324       }
00325     }
00326   }
00327 };
00328 
00329 } // end namespace internal
00330 
00331 } // end namespace Eigen
00332 
00333 #endif // EIGEN_TRIANGULAR_SOLVER_MATRIX_H