TriangularSolverMatrix.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // Eigen is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 3 of the License, or (at your option) any later version.
10 //
11 // Alternatively, you can redistribute it and/or
12 // modify it under the terms of the GNU General Public License as
13 // published by the Free Software Foundation; either version 2 of
14 // the License, or (at your option) any later version.
15 //
16 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU Lesser General Public
22 // License and a copy of the GNU General Public License along with
23 // Eigen. If not, see <http://www.gnu.org/licenses/>.
24 
25 #ifndef EIGEN_TRIANGULAR_SOLVER_MATRIX_H
26 #define EIGEN_TRIANGULAR_SOLVER_MATRIX_H
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
32 // if the rhs is row major, let's transpose the product
33 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder>
34 struct triangular_solve_matrix<Scalar,Index,Side,Mode,Conjugate,TriStorageOrder,RowMajor>
35 {
36  static EIGEN_DONT_INLINE void run(
37  Index size, Index cols,
38  const Scalar* tri, Index triStride,
39  Scalar* _other, Index otherStride,
40  level3_blocking<Scalar,Scalar>& blocking)
41  {
42  triangular_solve_matrix<
43  Scalar, Index, Side==OnTheLeft?OnTheRight:OnTheLeft,
44  (Mode&UnitDiag) | ((Mode&Upper) ? Lower : Upper),
45  NumTraits<Scalar>::IsComplex && Conjugate,
46  TriStorageOrder==RowMajor ? ColMajor : RowMajor, ColMajor>
47  ::run(size, cols, tri, triStride, _other, otherStride, blocking);
48  }
49 };
50 
51 /* Optimized triangular solver with multiple right hand side and the triangular matrix on the left
52  */
53 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
54 struct triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conjugate,TriStorageOrder,ColMajor>
55 {
56  static EIGEN_DONT_INLINE void run(
57  Index size, Index otherSize,
58  const Scalar* _tri, Index triStride,
59  Scalar* _other, Index otherStride,
60  level3_blocking<Scalar,Scalar>& blocking)
61  {
62  Index cols = otherSize;
63  const_blas_data_mapper<Scalar, Index, TriStorageOrder> tri(_tri,triStride);
64  blas_data_mapper<Scalar, Index, ColMajor> other(_other,otherStride);
65 
66  typedef gebp_traits<Scalar,Scalar> Traits;
67  enum {
68  SmallPanelWidth = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
69  IsLower = (Mode&Lower) == Lower
70  };
71 
72  Index kc = blocking.kc(); // cache block size along the K direction
73  Index mc = (std::min)(size,blocking.mc()); // cache block size along the M direction
74 
75  std::size_t sizeA = kc*mc;
76  std::size_t sizeB = kc*cols;
77  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
78 
79  ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
80  ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
81  ei_declare_aligned_stack_constructed_variable(Scalar, blockW, sizeW, blocking.blockW());
82 
83  conj_if<Conjugate> conj;
84  gebp_kernel<Scalar, Scalar, Index, Traits::mr, Traits::nr, Conjugate, false> gebp_kernel;
85  gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, TriStorageOrder> pack_lhs;
86  gemm_pack_rhs<Scalar, Index, Traits::nr, ColMajor, false, true> pack_rhs;
87 
88  // the goal here is to subdivise the Rhs panels such that we keep some cache
89  // coherence when accessing the rhs elements
90  std::ptrdiff_t l1, l2;
92  Index subcols = cols>0 ? l2/(4 * sizeof(Scalar) * otherStride) : 0;
93  subcols = std::max<Index>((subcols/Traits::nr)*Traits::nr, Traits::nr);
94 
95  for(Index k2=IsLower ? 0 : size;
96  IsLower ? k2<size : k2>0;
97  IsLower ? k2+=kc : k2-=kc)
98  {
99  const Index actual_kc = (std::min)(IsLower ? size-k2 : k2, kc);
100 
101  // We have selected and packed a big horizontal panel R1 of rhs. Let B be the packed copy of this panel,
102  // and R2 the remaining part of rhs. The corresponding vertical panel of lhs is split into
103  // A11 (the triangular part) and A21 the remaining rectangular part.
104  // Then the high level algorithm is:
105  // - B = R1 => general block copy (done during the next step)
106  // - R1 = A11^-1 B => tricky part
107  // - update B from the new R1 => actually this has to be performed continuously during the above step
108  // - R2 -= A21 * B => GEPP
109 
110  // The tricky part: compute R1 = A11^-1 B while updating B from R1
111  // The idea is to split A11 into multiple small vertical panels.
112  // Each panel can be split into a small triangular part T1k which is processed without optimization,
113  // and the remaining small part T2k which is processed using gebp with appropriate block strides
114  for(Index j2=0; j2<cols; j2+=subcols)
115  {
116  Index actual_cols = (std::min)(cols-j2,subcols);
117  // for each small vertical panels [T1k^T, T2k^T]^T of lhs
118  for (Index k1=0; k1<actual_kc; k1+=SmallPanelWidth)
119  {
120  Index actualPanelWidth = std::min<Index>(actual_kc-k1, SmallPanelWidth);
121  // tr solve
122  for (Index k=0; k<actualPanelWidth; ++k)
123  {
124  // TODO write a small kernel handling this (can be shared with trsv)
125  Index i = IsLower ? k2+k1+k : k2-k1-k-1;
126  Index s = IsLower ? k2+k1 : i+1;
127  Index rs = actualPanelWidth - k - 1; // remaining size
128 
129  Scalar a = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(tri(i,i));
130  for (Index j=j2; j<j2+actual_cols; ++j)
131  {
132  if (TriStorageOrder==RowMajor)
133  {
134  Scalar b(0);
135  const Scalar* l = &tri(i,s);
136  Scalar* r = &other(s,j);
137  for (Index i3=0; i3<k; ++i3)
138  b += conj(l[i3]) * r[i3];
139 
140  other(i,j) = (other(i,j) - b)*a;
141  }
142  else
143  {
144  Index s = IsLower ? i+1 : i-rs;
145  Scalar b = (other(i,j) *= a);
146  Scalar* r = &other(s,j);
147  const Scalar* l = &tri(s,i);
148  for (Index i3=0;i3<rs;++i3)
149  r[i3] -= b * conj(l[i3]);
150  }
151  }
152  }
153 
154  Index lengthTarget = actual_kc-k1-actualPanelWidth;
155  Index startBlock = IsLower ? k2+k1 : k2-k1-actualPanelWidth;
156  Index blockBOffset = IsLower ? k1 : lengthTarget;
157 
158  // update the respective rows of B from other
159  pack_rhs(blockB+actual_kc*j2, &other(startBlock,j2), otherStride, actualPanelWidth, actual_cols, actual_kc, blockBOffset);
160 
161  // GEBP
162  if (lengthTarget>0)
163  {
164  Index startTarget = IsLower ? k2+k1+actualPanelWidth : k2-actual_kc;
165 
166  pack_lhs(blockA, &tri(startTarget,startBlock), triStride, actualPanelWidth, lengthTarget);
167 
168  gebp_kernel(&other(startTarget,j2), otherStride, blockA, blockB+actual_kc*j2, lengthTarget, actualPanelWidth, actual_cols, Scalar(-1),
169  actualPanelWidth, actual_kc, 0, blockBOffset, blockW);
170  }
171  }
172  }
173 
174  // R2 -= A21 * B => GEPP
175  {
176  Index start = IsLower ? k2+kc : 0;
177  Index end = IsLower ? size : k2-kc;
178  for(Index i2=start; i2<end; i2+=mc)
179  {
180  const Index actual_mc = (std::min)(mc,end-i2);
181  if (actual_mc>0)
182  {
183  pack_lhs(blockA, &tri(i2, IsLower ? k2 : k2-kc), triStride, actual_kc, actual_mc);
184 
185  gebp_kernel(_other+i2, otherStride, blockA, blockB, actual_mc, actual_kc, cols, Scalar(-1), -1, -1, 0, 0, blockW);
186  }
187  }
188  }
189  }
190  }
191 };
192 
193 /* Optimized triangular solver with multiple left hand sides and the trinagular matrix on the right
194  */
195 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
196 struct triangular_solve_matrix<Scalar,Index,OnTheRight,Mode,Conjugate,TriStorageOrder,ColMajor>
197 {
198  static EIGEN_DONT_INLINE void run(
199  Index size, Index otherSize,
200  const Scalar* _tri, Index triStride,
201  Scalar* _other, Index otherStride,
202  level3_blocking<Scalar,Scalar>& blocking)
203  {
204  Index rows = otherSize;
205  const_blas_data_mapper<Scalar, Index, TriStorageOrder> rhs(_tri,triStride);
206  blas_data_mapper<Scalar, Index, ColMajor> lhs(_other,otherStride);
207 
208  typedef gebp_traits<Scalar,Scalar> Traits;
209  enum {
210  RhsStorageOrder = TriStorageOrder,
211  SmallPanelWidth = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
212  IsLower = (Mode&Lower) == Lower
213  };
214 
215  Index kc = blocking.kc(); // cache block size along the K direction
216  Index mc = (std::min)(rows,blocking.mc()); // cache block size along the M direction
217 
218  std::size_t sizeA = kc*mc;
219  std::size_t sizeB = kc*size;
220  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
221 
222  ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
223  ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
224  ei_declare_aligned_stack_constructed_variable(Scalar, blockW, sizeW, blocking.blockW());
225 
226  conj_if<Conjugate> conj;
227  gebp_kernel<Scalar,Scalar, Index, Traits::mr, Traits::nr, false, Conjugate> gebp_kernel;
228  gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder> pack_rhs;
229  gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder,false,true> pack_rhs_panel;
230  gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, ColMajor, false, true> pack_lhs_panel;
231 
232  for(Index k2=IsLower ? size : 0;
233  IsLower ? k2>0 : k2<size;
234  IsLower ? k2-=kc : k2+=kc)
235  {
236  const Index actual_kc = (std::min)(IsLower ? k2 : size-k2, kc);
237  Index actual_k2 = IsLower ? k2-actual_kc : k2 ;
238 
239  Index startPanel = IsLower ? 0 : k2+actual_kc;
240  Index rs = IsLower ? actual_k2 : size - actual_k2 - actual_kc;
241  Scalar* geb = blockB+actual_kc*actual_kc;
242 
243  if (rs>0) pack_rhs(geb, &rhs(actual_k2,startPanel), triStride, actual_kc, rs);
244 
245  // triangular packing (we only pack the panels off the diagonal,
246  // neglecting the blocks overlapping the diagonal
247  {
248  for (Index j2=0; j2<actual_kc; j2+=SmallPanelWidth)
249  {
250  Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
251  Index actual_j2 = actual_k2 + j2;
252  Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
253  Index panelLength = IsLower ? actual_kc-j2-actualPanelWidth : j2;
254 
255  if (panelLength>0)
256  pack_rhs_panel(blockB+j2*actual_kc,
257  &rhs(actual_k2+panelOffset, actual_j2), triStride,
258  panelLength, actualPanelWidth,
259  actual_kc, panelOffset);
260  }
261  }
262 
263  for(Index i2=0; i2<rows; i2+=mc)
264  {
265  const Index actual_mc = (std::min)(mc,rows-i2);
266 
267  // triangular solver kernel
268  {
269  // for each small block of the diagonal (=> vertical panels of rhs)
270  for (Index j2 = IsLower
271  ? (actual_kc - ((actual_kc%SmallPanelWidth) ? Index(actual_kc%SmallPanelWidth)
272  : Index(SmallPanelWidth)))
273  : 0;
274  IsLower ? j2>=0 : j2<actual_kc;
275  IsLower ? j2-=SmallPanelWidth : j2+=SmallPanelWidth)
276  {
277  Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
278  Index absolute_j2 = actual_k2 + j2;
279  Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
280  Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;
281 
282  // GEBP
283  if(panelLength>0)
284  {
285  gebp_kernel(&lhs(i2,absolute_j2), otherStride,
286  blockA, blockB+j2*actual_kc,
287  actual_mc, panelLength, actualPanelWidth,
288  Scalar(-1),
289  actual_kc, actual_kc, // strides
290  panelOffset, panelOffset, // offsets
291  blockW); // workspace
292  }
293 
294  // unblocked triangular solve
295  for (Index k=0; k<actualPanelWidth; ++k)
296  {
297  Index j = IsLower ? absolute_j2+actualPanelWidth-k-1 : absolute_j2+k;
298 
299  Scalar* r = &lhs(i2,j);
300  for (Index k3=0; k3<k; ++k3)
301  {
302  Scalar b = conj(rhs(IsLower ? j+1+k3 : absolute_j2+k3,j));
303  Scalar* a = &lhs(i2,IsLower ? j+1+k3 : absolute_j2+k3);
304  for (Index i=0; i<actual_mc; ++i)
305  r[i] -= a[i] * b;
306  }
307  Scalar b = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(rhs(j,j));
308  for (Index i=0; i<actual_mc; ++i)
309  r[i] *= b;
310  }
311 
312  // pack the just computed part of lhs to A
313  pack_lhs_panel(blockA, _other+absolute_j2*otherStride+i2, otherStride,
314  actualPanelWidth, actual_mc,
315  actual_kc, j2);
316  }
317  }
318 
319  if (rs>0)
320  gebp_kernel(_other+i2+startPanel*otherStride, otherStride, blockA, geb,
321  actual_mc, actual_kc, rs, Scalar(-1),
322  -1, -1, 0, 0, blockW);
323  }
324  }
325  }
326 };
327 
328 } // end namespace internal
329 
330 } // end namespace Eigen
331 
332 #endif // EIGEN_TRIANGULAR_SOLVER_MATRIX_H