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


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:12:29