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  {
41  triangular_solve_matrix<
42  Scalar, Index, Side==OnTheLeft?OnTheRight:OnTheLeft,
43  (Mode&UnitDiag) | ((Mode&Upper) ? Lower : Upper),
44  NumTraits<Scalar>::IsComplex && Conjugate,
45  TriStorageOrder==RowMajor ? ColMajor : RowMajor, ColMajor>
46  ::run(size, cols, tri, triStride, _other, otherStride);
47  }
48 };
49 
50 /* Optimized triangular solver with multiple right hand side and the triangular matrix on the left
51  */
52 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
53 struct triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conjugate,TriStorageOrder,ColMajor>
54 {
55  static EIGEN_DONT_INLINE void run(
56  Index size, Index otherSize,
57  const Scalar* _tri, Index triStride,
58  Scalar* _other, Index otherStride)
59  {
60  Index cols = otherSize;
61  const_blas_data_mapper<Scalar, Index, TriStorageOrder> tri(_tri,triStride);
62  blas_data_mapper<Scalar, Index, ColMajor> other(_other,otherStride);
63 
64  typedef gebp_traits<Scalar,Scalar> Traits;
65  enum {
66  SmallPanelWidth = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
67  IsLower = (Mode&Lower) == Lower
68  };
69 
70  Index kc = size; // cache block size along the K direction
71  Index mc = size; // cache block size along the M direction
72  Index nc = cols; // cache block size along the N direction
73  computeProductBlockingSizes<Scalar,Scalar,4>(kc, mc, nc);
74 
75  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
76  std::size_t sizeB = sizeW + kc*cols;
77  ei_declare_aligned_stack_constructed_variable(Scalar, blockA, kc*mc, 0);
78  ei_declare_aligned_stack_constructed_variable(Scalar, allocatedBlockB, sizeB, 0);
79  Scalar* blockB = allocatedBlockB + sizeW;
80  Scalar* blockW = allocatedBlockB;
81 
82  conj_if<Conjugate> conj;
83  gebp_kernel<Scalar, Scalar, Index, Traits::mr, Traits::nr, Conjugate, false> gebp_kernel;
84  gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, TriStorageOrder> pack_lhs;
85  gemm_pack_rhs<Scalar, Index, Traits::nr, ColMajor, false, true> pack_rhs;
86 
87  // the goal here is to subdivise the Rhs panels such that we keep some cache
88  // coherence when accessing the rhs elements
89  std::ptrdiff_t l1, l2;
91  Index subcols = cols>0 ? l2/(4 * sizeof(Scalar) * otherStride) : 0;
92  subcols = std::max<Index>((subcols/Traits::nr)*Traits::nr, Traits::nr);
93 
94  for(Index k2=IsLower ? 0 : size;
95  IsLower ? k2<size : k2>0;
96  IsLower ? k2+=kc : k2-=kc)
97  {
98  const Index actual_kc = (std::min)(IsLower ? size-k2 : k2, kc);
99 
100  // We have selected and packed a big horizontal panel R1 of rhs. Let B be the packed copy of this panel,
101  // and R2 the remaining part of rhs. The corresponding vertical panel of lhs is split into
102  // A11 (the triangular part) and A21 the remaining rectangular part.
103  // Then the high level algorithm is:
104  // - B = R1 => general block copy (done during the next step)
105  // - R1 = A11^-1 B => tricky part
106  // - update B from the new R1 => actually this has to be performed continuously during the above step
107  // - R2 -= A21 * B => GEPP
108 
109  // The tricky part: compute R1 = A11^-1 B while updating B from R1
110  // The idea is to split A11 into multiple small vertical panels.
111  // Each panel can be split into a small triangular part T1k which is processed without optimization,
112  // and the remaining small part T2k which is processed using gebp with appropriate block strides
113  for(Index j2=0; j2<cols; j2+=subcols)
114  {
115  Index actual_cols = (std::min)(cols-j2,subcols);
116  // for each small vertical panels [T1k^T, T2k^T]^T of lhs
117  for (Index k1=0; k1<actual_kc; k1+=SmallPanelWidth)
118  {
119  Index actualPanelWidth = std::min<Index>(actual_kc-k1, SmallPanelWidth);
120  // tr solve
121  for (Index k=0; k<actualPanelWidth; ++k)
122  {
123  // TODO write a small kernel handling this (can be shared with trsv)
124  Index i = IsLower ? k2+k1+k : k2-k1-k-1;
125  Index s = IsLower ? k2+k1 : i+1;
126  Index rs = actualPanelWidth - k - 1; // remaining size
127 
128  Scalar a = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(tri(i,i));
129  for (Index j=j2; j<j2+actual_cols; ++j)
130  {
131  if (TriStorageOrder==RowMajor)
132  {
133  Scalar b(0);
134  const Scalar* l = &tri(i,s);
135  Scalar* r = &other(s,j);
136  for (Index i3=0; i3<k; ++i3)
137  b += conj(l[i3]) * r[i3];
138 
139  other(i,j) = (other(i,j) - b)*a;
140  }
141  else
142  {
143  Index s = IsLower ? i+1 : i-rs;
144  Scalar b = (other(i,j) *= a);
145  Scalar* r = &other(s,j);
146  const Scalar* l = &tri(s,i);
147  for (Index i3=0;i3<rs;++i3)
148  r[i3] -= b * conj(l[i3]);
149  }
150  }
151  }
152 
153  Index lengthTarget = actual_kc-k1-actualPanelWidth;
154  Index startBlock = IsLower ? k2+k1 : k2-k1-actualPanelWidth;
155  Index blockBOffset = IsLower ? k1 : lengthTarget;
156 
157  // update the respective rows of B from other
158  pack_rhs(blockB+actual_kc*j2, &other(startBlock,j2), otherStride, actualPanelWidth, actual_cols, actual_kc, blockBOffset);
159 
160  // GEBP
161  if (lengthTarget>0)
162  {
163  Index startTarget = IsLower ? k2+k1+actualPanelWidth : k2-actual_kc;
164 
165  pack_lhs(blockA, &tri(startTarget,startBlock), triStride, actualPanelWidth, lengthTarget);
166 
167  gebp_kernel(&other(startTarget,j2), otherStride, blockA, blockB+actual_kc*j2, lengthTarget, actualPanelWidth, actual_cols, Scalar(-1),
168  actualPanelWidth, actual_kc, 0, blockBOffset, blockW);
169  }
170  }
171  }
172 
173  // R2 -= A21 * B => GEPP
174  {
175  Index start = IsLower ? k2+kc : 0;
176  Index end = IsLower ? size : k2-kc;
177  for(Index i2=start; i2<end; i2+=mc)
178  {
179  const Index actual_mc = (std::min)(mc,end-i2);
180  if (actual_mc>0)
181  {
182  pack_lhs(blockA, &tri(i2, IsLower ? k2 : k2-kc), triStride, actual_kc, actual_mc);
183 
184  gebp_kernel(_other+i2, otherStride, blockA, blockB, actual_mc, actual_kc, cols, Scalar(-1));
185  }
186  }
187  }
188  }
189  }
190 };
191 
192 /* Optimized triangular solver with multiple left hand sides and the trinagular matrix on the right
193  */
194 template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
195 struct triangular_solve_matrix<Scalar,Index,OnTheRight,Mode,Conjugate,TriStorageOrder,ColMajor>
196 {
197  static EIGEN_DONT_INLINE void run(
198  Index size, Index otherSize,
199  const Scalar* _tri, Index triStride,
200  Scalar* _other, Index otherStride)
201  {
202  Index rows = otherSize;
203  const_blas_data_mapper<Scalar, Index, TriStorageOrder> rhs(_tri,triStride);
204  blas_data_mapper<Scalar, Index, ColMajor> lhs(_other,otherStride);
205 
206  typedef gebp_traits<Scalar,Scalar> Traits;
207  enum {
208  RhsStorageOrder = TriStorageOrder,
209  SmallPanelWidth = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
210  IsLower = (Mode&Lower) == Lower
211  };
212 
213 // Index kc = std::min<Index>(Traits::Max_kc/4,size); // cache block size along the K direction
214 // Index mc = std::min<Index>(Traits::Max_mc,size); // cache block size along the M direction
215  // check that !!!!
216  Index kc = size; // cache block size along the K direction
217  Index mc = size; // cache block size along the M direction
218  Index nc = rows; // cache block size along the N direction
219  computeProductBlockingSizes<Scalar,Scalar,4>(kc, mc, nc);
220 
221  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
222  std::size_t sizeB = sizeW + kc*size;
223  ei_declare_aligned_stack_constructed_variable(Scalar, blockA, kc*mc, 0);
224  ei_declare_aligned_stack_constructed_variable(Scalar, allocatedBlockB, sizeB, 0);
225  Scalar* blockB = allocatedBlockB + sizeW;
226 
227  conj_if<Conjugate> conj;
228  gebp_kernel<Scalar,Scalar, Index, Traits::mr, Traits::nr, false, Conjugate> gebp_kernel;
229  gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder> pack_rhs;
230  gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder,false,true> pack_rhs_panel;
231  gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, ColMajor, false, true> pack_lhs_panel;
232 
233  for(Index k2=IsLower ? size : 0;
234  IsLower ? k2>0 : k2<size;
235  IsLower ? k2-=kc : k2+=kc)
236  {
237  const Index actual_kc = (std::min)(IsLower ? k2 : size-k2, kc);
238  Index actual_k2 = IsLower ? k2-actual_kc : k2 ;
239 
240  Index startPanel = IsLower ? 0 : k2+actual_kc;
241  Index rs = IsLower ? actual_k2 : size - actual_k2 - actual_kc;
242  Scalar* geb = blockB+actual_kc*actual_kc;
243 
244  if (rs>0) pack_rhs(geb, &rhs(actual_k2,startPanel), triStride, actual_kc, rs);
245 
246  // triangular packing (we only pack the panels off the diagonal,
247  // neglecting the blocks overlapping the diagonal
248  {
249  for (Index j2=0; j2<actual_kc; j2+=SmallPanelWidth)
250  {
251  Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
252  Index actual_j2 = actual_k2 + j2;
253  Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
254  Index panelLength = IsLower ? actual_kc-j2-actualPanelWidth : j2;
255 
256  if (panelLength>0)
257  pack_rhs_panel(blockB+j2*actual_kc,
258  &rhs(actual_k2+panelOffset, actual_j2), triStride,
259  panelLength, actualPanelWidth,
260  actual_kc, panelOffset);
261  }
262  }
263 
264  for(Index i2=0; i2<rows; i2+=mc)
265  {
266  const Index actual_mc = (std::min)(mc,rows-i2);
267 
268  // triangular solver kernel
269  {
270  // for each small block of the diagonal (=> vertical panels of rhs)
271  for (Index j2 = IsLower
272  ? (actual_kc - ((actual_kc%SmallPanelWidth) ? Index(actual_kc%SmallPanelWidth)
273  : Index(SmallPanelWidth)))
274  : 0;
275  IsLower ? j2>=0 : j2<actual_kc;
276  IsLower ? j2-=SmallPanelWidth : j2+=SmallPanelWidth)
277  {
278  Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
279  Index absolute_j2 = actual_k2 + j2;
280  Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
281  Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;
282 
283  // GEBP
284  if(panelLength>0)
285  {
286  gebp_kernel(&lhs(i2,absolute_j2), otherStride,
287  blockA, blockB+j2*actual_kc,
288  actual_mc, panelLength, actualPanelWidth,
289  Scalar(-1),
290  actual_kc, actual_kc, // strides
291  panelOffset, panelOffset, // offsets
292  allocatedBlockB); // workspace
293  }
294 
295  // unblocked triangular solve
296  for (Index k=0; k<actualPanelWidth; ++k)
297  {
298  Index j = IsLower ? absolute_j2+actualPanelWidth-k-1 : absolute_j2+k;
299 
300  Scalar* r = &lhs(i2,j);
301  for (Index k3=0; k3<k; ++k3)
302  {
303  Scalar b = conj(rhs(IsLower ? j+1+k3 : absolute_j2+k3,j));
304  Scalar* a = &lhs(i2,IsLower ? j+1+k3 : absolute_j2+k3);
305  for (Index i=0; i<actual_mc; ++i)
306  r[i] -= a[i] * b;
307  }
308  Scalar b = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(rhs(j,j));
309  for (Index i=0; i<actual_mc; ++i)
310  r[i] *= b;
311  }
312 
313  // pack the just computed part of lhs to A
314  pack_lhs_panel(blockA, _other+absolute_j2*otherStride+i2, otherStride,
315  actualPanelWidth, actual_mc,
316  actual_kc, j2);
317  }
318  }
319 
320  if (rs>0)
321  gebp_kernel(_other+i2+startPanel*otherStride, otherStride, blockA, geb,
322  actual_mc, actual_kc, rs, Scalar(-1),
323  -1, -1, 0, 0, allocatedBlockB);
324  }
325  }
326  }
327 };
328 
329 } // end namespace internal
330 
331 } // end namespace Eigen
332 
333 #endif // EIGEN_TRIANGULAR_SOLVER_MATRIX_H