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


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:33:18