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