CMS 3D CMS Logo

choleskyInversion.h
Go to the documentation of this file.
1 #ifndef DataFormat_Math_choleskyInversion_h
2 #define DataFormat_Math_choleskyInversion_h
3 
4 #include <cmath>
5 
6 #include <Eigen/Core>
7 
19 namespace math {
20  namespace cholesky {
21 
22  template <typename M1, typename M2>
23  inline constexpr void __attribute__((always_inline)) invert11(M1 const& src, M2& dst) {
24  using F = decltype(src(0, 0));
25  dst(0, 0) = F(1.0) / src(0, 0);
26  }
27 
28  template <typename M1, typename M2>
29  inline constexpr void __attribute__((always_inline)) invert22(M1 const& src, M2& dst) {
30  using F = decltype(src(0, 0));
31  auto luc0 = F(1.0) / src(0, 0);
32  auto luc1 = src(1, 0) * src(1, 0) * luc0;
33  auto luc2 = F(1.0) / (src(1, 1) - luc1);
34 
35  auto li21 = luc1 * luc0 * luc2;
36 
37  dst(0, 0) = li21 + luc0;
38  dst(1, 0) = -src(1, 0) * luc0 * luc2;
39  dst(1, 1) = luc2;
40  }
41 
42  template <typename M1, typename M2>
43  inline constexpr void __attribute__((always_inline)) invert33(M1 const& src, M2& dst) {
44  using F = decltype(src(0, 0));
45  auto luc0 = F(1.0) / src(0, 0);
46  auto luc1 = src(1, 0);
47  auto luc2 = src(1, 1) - luc0 * luc1 * luc1;
48  luc2 = F(1.0) / luc2;
49  auto luc3 = src(2, 0);
50  auto luc4 = (src(2, 1) - luc0 * luc1 * luc3);
51  auto luc5 = src(2, 2) - (luc0 * luc3 * luc3 + (luc2 * luc4) * luc4);
52  luc5 = F(1.0) / luc5;
53 
54  auto li21 = -luc0 * luc1;
55  auto li32 = -(luc2 * luc4);
56  auto li31 = (luc1 * (luc2 * luc4) - luc3) * luc0;
57 
58  dst(0, 0) = luc5 * li31 * li31 + li21 * li21 * luc2 + luc0;
59  dst(1, 0) = luc5 * li31 * li32 + li21 * luc2;
60  dst(1, 1) = luc5 * li32 * li32 + luc2;
61  dst(2, 0) = luc5 * li31;
62  dst(2, 1) = luc5 * li32;
63  dst(2, 2) = luc5;
64  }
65 
66  template <typename M1, typename M2>
67  inline constexpr void __attribute__((always_inline)) invert44(M1 const& src, M2& dst) {
68  using F = decltype(src(0, 0));
69  auto luc0 = F(1.0) / src(0, 0);
70  auto luc1 = src(1, 0);
71  auto luc2 = src(1, 1) - luc0 * luc1 * luc1;
72  luc2 = F(1.0) / luc2;
73  auto luc3 = src(2, 0);
74  auto luc4 = (src(2, 1) - luc0 * luc1 * luc3);
75  auto luc5 = src(2, 2) - (luc0 * luc3 * luc3 + luc2 * luc4 * luc4);
76  luc5 = F(1.0) / luc5;
77  auto luc6 = src(3, 0);
78  auto luc7 = (src(3, 1) - luc0 * luc1 * luc6);
79  auto luc8 = (src(3, 2) - luc0 * luc3 * luc6 - luc2 * luc4 * luc7);
80  auto luc9 = src(3, 3) - (luc0 * luc6 * luc6 + luc2 * luc7 * luc7 + luc8 * (luc8 * luc5));
81  luc9 = F(1.0) / luc9;
82 
83  auto li21 = -luc1 * luc0;
84  auto li32 = -luc2 * luc4;
85  auto li31 = (luc1 * (luc2 * luc4) - luc3) * luc0;
86  auto li43 = -(luc8 * luc5);
87  auto li42 = (luc4 * luc8 * luc5 - luc7) * luc2;
88  auto li41 = (-luc1 * (luc2 * luc4) * (luc8 * luc5) + luc1 * (luc2 * luc7) + luc3 * (luc8 * luc5) - luc6) * luc0;
89 
90  dst(0, 0) = luc9 * li41 * li41 + luc5 * li31 * li31 + luc2 * li21 * li21 + luc0;
91  dst(1, 0) = luc9 * li41 * li42 + luc5 * li31 * li32 + luc2 * li21;
92  dst(1, 1) = luc9 * li42 * li42 + luc5 * li32 * li32 + luc2;
93  dst(2, 0) = luc9 * li41 * li43 + luc5 * li31;
94  dst(2, 1) = luc9 * li42 * li43 + luc5 * li32;
95  dst(2, 2) = luc9 * li43 * li43 + luc5;
96  dst(3, 0) = luc9 * li41;
97  dst(3, 1) = luc9 * li42;
98  dst(3, 2) = luc9 * li43;
99  dst(3, 3) = luc9;
100  }
101 
102  template <typename M1, typename M2>
103  inline constexpr void __attribute__((always_inline)) invert55(M1 const& src, M2& dst) {
104  using F = decltype(src(0, 0));
105  auto luc0 = F(1.0) / src(0, 0);
106  auto luc1 = src(1, 0);
107  auto luc2 = src(1, 1) - luc0 * luc1 * luc1;
108  luc2 = F(1.0) / luc2;
109  auto luc3 = src(2, 0);
110  auto luc4 = (src(2, 1) - luc0 * luc1 * luc3);
111  auto luc5 = src(2, 2) - (luc0 * luc3 * luc3 + luc2 * luc4 * luc4);
112  luc5 = F(1.0) / luc5;
113  auto luc6 = src(3, 0);
114  auto luc7 = (src(3, 1) - luc0 * luc1 * luc6);
115  auto luc8 = (src(3, 2) - luc0 * luc3 * luc6 - luc2 * luc4 * luc7);
116  auto luc9 = src(3, 3) - (luc0 * luc6 * luc6 + luc2 * luc7 * luc7 + luc8 * (luc8 * luc5));
117  luc9 = F(1.0) / luc9;
118  auto luc10 = src(4, 0);
119  auto luc11 = (src(4, 1) - luc0 * luc1 * luc10);
120  auto luc12 = (src(4, 2) - luc0 * luc3 * luc10 - luc2 * luc4 * luc11);
121  auto luc13 = (src(4, 3) - luc0 * luc6 * luc10 - luc2 * luc7 * luc11 - luc5 * luc8 * luc12);
122  auto luc14 =
123  src(4, 4) - (luc0 * luc10 * luc10 + luc2 * luc11 * luc11 + luc5 * luc12 * luc12 + luc9 * luc13 * luc13);
124  luc14 = F(1.0) / luc14;
125 
126  auto li21 = -luc1 * luc0;
127  auto li32 = -luc2 * luc4;
128  auto li31 = (luc1 * (luc2 * luc4) - luc3) * luc0;
129  auto li43 = -(luc8 * luc5);
130  auto li42 = (luc4 * luc8 * luc5 - luc7) * luc2;
131  auto li41 = (-luc1 * (luc2 * luc4) * (luc8 * luc5) + luc1 * (luc2 * luc7) + luc3 * (luc8 * luc5) - luc6) * luc0;
132  auto li54 = -luc13 * luc9;
133  auto li53 = (luc13 * luc8 * luc9 - luc12) * luc5;
134  auto li52 = (-luc4 * luc8 * luc13 * luc5 * luc9 + luc4 * luc12 * luc5 + luc7 * luc13 * luc9 - luc11) * luc2;
135  auto li51 = (luc1 * luc4 * luc8 * luc13 * luc2 * luc5 * luc9 - luc13 * luc8 * luc3 * luc9 * luc5 -
136  luc12 * luc4 * luc1 * luc2 * luc5 - luc13 * luc7 * luc1 * luc9 * luc2 + luc11 * luc1 * luc2 +
137  luc12 * luc3 * luc5 + luc13 * luc6 * luc9 - luc10) *
138  luc0;
139 
140  dst(0, 0) = luc14 * li51 * li51 + luc9 * li41 * li41 + luc5 * li31 * li31 + luc2 * li21 * li21 + luc0;
141  dst(1, 0) = luc14 * li51 * li52 + luc9 * li41 * li42 + luc5 * li31 * li32 + luc2 * li21;
142  dst(1, 1) = luc14 * li52 * li52 + luc9 * li42 * li42 + luc5 * li32 * li32 + luc2;
143  dst(2, 0) = luc14 * li51 * li53 + luc9 * li41 * li43 + luc5 * li31;
144  dst(2, 1) = luc14 * li52 * li53 + luc9 * li42 * li43 + luc5 * li32;
145  dst(2, 2) = luc14 * li53 * li53 + luc9 * li43 * li43 + luc5;
146  dst(3, 0) = luc14 * li51 * li54 + luc9 * li41;
147  dst(3, 1) = luc14 * li52 * li54 + luc9 * li42;
148  dst(3, 2) = luc14 * li53 * li54 + luc9 * li43;
149  dst(3, 3) = luc14 * li54 * li54 + luc9;
150  dst(4, 0) = luc14 * li51;
151  dst(4, 1) = luc14 * li52;
152  dst(4, 2) = luc14 * li53;
153  dst(4, 3) = luc14 * li54;
154  dst(4, 4) = luc14;
155  }
156 
157  template <typename M1, typename M2>
158  inline constexpr void __attribute__((always_inline)) invert66(M1 const& src, M2& dst) {
159  using F = decltype(src(0, 0));
160  auto luc0 = F(1.0) / src(0, 0);
161  auto luc1 = src(1, 0);
162  auto luc2 = src(1, 1) - luc0 * luc1 * luc1;
163  luc2 = F(1.0) / luc2;
164  auto luc3 = src(2, 0);
165  auto luc4 = (src(2, 1) - luc0 * luc1 * luc3);
166  auto luc5 = src(2, 2) - (luc0 * luc3 * luc3 + luc2 * luc4 * luc4);
167  luc5 = F(1.0) / luc5;
168  auto luc6 = src(3, 0);
169  auto luc7 = (src(3, 1) - luc0 * luc1 * luc6);
170  auto luc8 = (src(3, 2) - luc0 * luc3 * luc6 - luc2 * luc4 * luc7);
171  auto luc9 = src(3, 3) - (luc0 * luc6 * luc6 + luc2 * luc7 * luc7 + luc8 * (luc8 * luc5));
172  luc9 = F(1.0) / luc9;
173  auto luc10 = src(4, 0);
174  auto luc11 = (src(4, 1) - luc0 * luc1 * luc10);
175  auto luc12 = (src(4, 2) - luc0 * luc3 * luc10 - luc2 * luc4 * luc11);
176  auto luc13 = (src(4, 3) - luc0 * luc6 * luc10 - luc2 * luc7 * luc11 - luc5 * luc8 * luc12);
177  auto luc14 =
178  src(4, 4) - (luc0 * luc10 * luc10 + luc2 * luc11 * luc11 + luc5 * luc12 * luc12 + luc9 * luc13 * luc13);
179  luc14 = F(1.0) / luc14;
180  auto luc15 = src(5, 0);
181  auto luc16 = (src(5, 1) - luc0 * luc1 * luc15);
182  auto luc17 = (src(5, 2) - luc0 * luc3 * luc15 - luc2 * luc4 * luc16);
183  auto luc18 = (src(5, 3) - luc0 * luc6 * luc15 - luc2 * luc7 * luc16 - luc5 * luc8 * luc17);
184  auto luc19 =
185  (src(5, 4) - luc0 * luc10 * luc15 - luc2 * luc11 * luc16 - luc5 * luc12 * luc17 - luc9 * luc13 * luc18);
186  auto luc20 = src(5, 5) - (luc0 * luc15 * luc15 + luc2 * luc16 * luc16 + luc5 * luc17 * luc17 +
187  luc9 * luc18 * luc18 + luc14 * luc19 * luc19);
188  luc20 = F(1.0) / luc20;
189 
190  auto li21 = -luc1 * luc0;
191  auto li32 = -luc2 * luc4;
192  auto li31 = (luc1 * (luc2 * luc4) - luc3) * luc0;
193  auto li43 = -(luc8 * luc5);
194  auto li42 = (luc4 * luc8 * luc5 - luc7) * luc2;
195  auto li41 = (-luc1 * (luc2 * luc4) * (luc8 * luc5) + luc1 * (luc2 * luc7) + luc3 * (luc8 * luc5) - luc6) * luc0;
196  auto li54 = -luc13 * luc9;
197  auto li53 = (luc13 * luc8 * luc9 - luc12) * luc5;
198  auto li52 = (-luc4 * luc8 * luc13 * luc5 * luc9 + luc4 * luc12 * luc5 + luc7 * luc13 * luc9 - luc11) * luc2;
199  auto li51 = (luc1 * luc4 * luc8 * luc13 * luc2 * luc5 * luc9 - luc13 * luc8 * luc3 * luc9 * luc5 -
200  luc12 * luc4 * luc1 * luc2 * luc5 - luc13 * luc7 * luc1 * luc9 * luc2 + luc11 * luc1 * luc2 +
201  luc12 * luc3 * luc5 + luc13 * luc6 * luc9 - luc10) *
202  luc0;
203 
204  auto li65 = -luc19 * luc14;
205  auto li64 = (luc19 * luc14 * luc13 - luc18) * luc9;
206  auto li63 =
207  (-luc8 * luc13 * (luc19 * luc14) * luc9 + luc8 * luc9 * luc18 + luc12 * (luc19 * luc14) - luc17) * luc5;
208  auto li62 = (luc4 * (luc8 * luc9) * luc13 * luc5 * (luc19 * luc14) - luc18 * luc4 * (luc8 * luc9) * luc5 -
209  luc19 * luc12 * luc4 * luc14 * luc5 - luc19 * luc13 * luc7 * luc14 * luc9 + luc17 * luc4 * luc5 +
210  luc18 * luc7 * luc9 + luc19 * luc11 * luc14 - luc16) *
211  luc2;
212  auto li61 =
213  (-luc19 * luc13 * luc8 * luc4 * luc1 * luc2 * luc5 * luc9 * luc14 +
214  luc18 * luc8 * luc4 * luc1 * luc2 * luc5 * luc9 + luc19 * luc12 * luc4 * luc1 * luc2 * luc5 * luc14 +
215  luc19 * luc13 * luc7 * luc1 * luc2 * luc9 * luc14 + luc19 * luc13 * luc8 * luc3 * luc5 * luc9 * luc14 -
216  luc17 * luc4 * luc1 * luc2 * luc5 - luc18 * luc7 * luc1 * luc2 * luc9 - luc19 * luc11 * luc1 * luc2 * luc14 -
217  luc18 * luc8 * luc3 * luc5 * luc9 - luc19 * luc12 * luc3 * luc5 * luc14 -
218  luc19 * luc13 * luc6 * luc9 * luc14 + luc16 * luc1 * luc2 + luc17 * luc3 * luc5 + luc18 * luc6 * luc9 +
219  luc19 * luc10 * luc14 - luc15) *
220  luc0;
221 
222  dst(0, 0) = luc20 * li61 * li61 + luc14 * li51 * li51 + luc9 * li41 * li41 + luc5 * li31 * li31 +
223  luc2 * li21 * li21 + luc0;
224  dst(1, 0) = luc20 * li61 * li62 + luc14 * li51 * li52 + luc9 * li41 * li42 + luc5 * li31 * li32 + luc2 * li21;
225  dst(1, 1) = luc20 * li62 * li62 + luc14 * li52 * li52 + luc9 * li42 * li42 + luc5 * li32 * li32 + luc2;
226  dst(2, 0) = luc20 * li61 * li63 + luc14 * li51 * li53 + luc9 * li41 * li43 + luc5 * li31;
227  dst(2, 1) = luc20 * li62 * li63 + luc14 * li52 * li53 + luc9 * li42 * li43 + luc5 * li32;
228  dst(2, 2) = luc20 * li63 * li63 + luc14 * li53 * li53 + luc9 * li43 * li43 + luc5;
229  dst(3, 0) = luc20 * li61 * li64 + luc14 * li51 * li54 + luc9 * li41;
230  dst(3, 1) = luc20 * li62 * li64 + luc14 * li52 * li54 + luc9 * li42;
231  dst(3, 2) = luc20 * li63 * li64 + luc14 * li53 * li54 + luc9 * li43;
232  dst(3, 3) = luc20 * li64 * li64 + luc14 * li54 * li54 + luc9;
233  dst(4, 0) = luc20 * li61 * li65 + luc14 * li51;
234  dst(4, 1) = luc20 * li62 * li65 + luc14 * li52;
235  dst(4, 2) = luc20 * li63 * li65 + luc14 * li53;
236  dst(4, 3) = luc20 * li64 * li65 + luc14 * li54;
237  dst(4, 4) = luc20 * li65 * li65 + luc14;
238  dst(5, 0) = luc20 * li61;
239  dst(5, 1) = luc20 * li62;
240  dst(5, 2) = luc20 * li63;
241  dst(5, 3) = luc20 * li64;
242  dst(5, 4) = luc20 * li65;
243  dst(5, 5) = luc20;
244  }
245 
246  template <typename M>
247  inline constexpr void symmetrize11(M& dst) {}
248 
249  template <typename M>
250  inline constexpr void symmetrize22(M& dst) {
251  dst(0, 1) = dst(1, 0);
252  }
253 
254  template <typename M>
255  inline constexpr void symmetrize33(M& dst) {
256  symmetrize22(dst);
257  dst(0, 2) = dst(2, 0);
258  dst(1, 2) = dst(2, 1);
259  }
260 
261  template <typename M>
262  inline constexpr void symmetrize44(M& dst) {
263  symmetrize33(dst);
264  dst(0, 3) = dst(3, 0);
265  dst(1, 3) = dst(3, 1);
266  dst(2, 3) = dst(3, 2);
267  }
268 
269  template <typename M>
270  inline constexpr void symmetrize55(M& dst) {
271  symmetrize44(dst);
272  dst(0, 4) = dst(4, 0);
273  dst(1, 4) = dst(4, 1);
274  dst(2, 4) = dst(4, 2);
275  dst(3, 4) = dst(4, 3);
276  }
277 
278  template <typename M>
279  inline constexpr void symmetrize66(M& dst) {
280  symmetrize55(dst);
281  dst(0, 5) = dst(5, 0);
282  dst(1, 5) = dst(5, 1);
283  dst(2, 5) = dst(5, 2);
284  dst(3, 5) = dst(5, 3);
285  dst(4, 5) = dst(5, 4);
286  }
287 
288  template <typename M1, typename M2, int N>
289  struct Inverter {
290  static constexpr void eval(M1 const& src, M2& dst) { dst = src.inverse(); }
291  };
292 
293  template <typename M1, typename M2>
294  struct Inverter<M1, M2, 1> {
295  static constexpr void eval(M1 const& src, M2& dst) { invert11(src, dst); }
296  };
297 
298  template <typename M1, typename M2>
299  struct Inverter<M1, M2, 2> {
300  static constexpr void __attribute__((always_inline)) eval(M1 const& src, M2& dst) {
301  invert22(src, dst);
302  symmetrize22(dst);
303  }
304  };
305 
306  template <typename M1, typename M2>
307  struct Inverter<M1, M2, 3> {
308  static constexpr void __attribute__((always_inline)) eval(M1 const& src, M2& dst) {
309  invert33(src, dst);
310  symmetrize33(dst);
311  }
312  };
313 
314  template <typename M1, typename M2>
315  struct Inverter<M1, M2, 4> {
316  static constexpr void __attribute__((always_inline)) eval(M1 const& src, M2& dst) {
317  invert44(src, dst);
318  symmetrize44(dst);
319  }
320  };
321 
322  template <typename M1, typename M2>
323  struct Inverter<M1, M2, 5> {
324  static constexpr void __attribute__((always_inline)) eval(M1 const& src, M2& dst) {
325  invert55(src, dst);
326  symmetrize55(dst);
327  }
328  };
329 
330  template <typename M1, typename M2>
331  struct Inverter<M1, M2, 6> {
332  static constexpr void __attribute__((always_inline)) eval(M1 const& src, M2& dst) {
333  invert66(src, dst);
334  symmetrize66(dst);
335  }
336  };
337 
338  // Eigen interface
339  template <typename D1, typename D2>
340  inline constexpr void __attribute__((always_inline))
341  invert(Eigen::DenseBase<D1> const& src, Eigen::DenseBase<D2>& dst) {
342  using M1 = Eigen::DenseBase<D1>;
343  using M2 = Eigen::DenseBase<D2>;
345  }
346 
347  } // namespace cholesky
348 } // namespace math
349 
350 #endif // DataFormat_Math_choleskyInversion_h
math::cholesky::luc13
auto luc13
Definition: choleskyInversion.h:121
math::cholesky::luc15
auto luc15
Definition: choleskyInversion.h:180
math::cholesky::li51
auto li51
Definition: choleskyInversion.h:135
math::cholesky::luc2
auto luc2
Definition: choleskyInversion.h:33
makeMuonMisalignmentScenario.cholesky
def cholesky(A)
Definition: makeMuonMisalignmentScenario.py:112
math::cholesky::luc9
auto luc9
Definition: choleskyInversion.h:80
math::cholesky::luc12
auto luc12
Definition: choleskyInversion.h:120
math::cholesky::li41
auto li41
Definition: choleskyInversion.h:88
math::cholesky::li63
auto li63
Definition: choleskyInversion.h:206
math::cholesky::li52
auto li52
Definition: choleskyInversion.h:134
math::cholesky::luc3
auto luc3
Definition: choleskyInversion.h:49
math::cholesky::__attribute__
constexpr void __attribute__((always_inline)) invert11(M1 const &src
math::cholesky::luc10
auto luc10
Definition: choleskyInversion.h:118
F
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
math::cholesky::luc6
auto luc6
Definition: choleskyInversion.h:77
math::cholesky::symmetrize66
constexpr void symmetrize66(M &dst)
Definition: choleskyInversion.h:279
math::cholesky::li43
auto li43
Definition: choleskyInversion.h:86
math::cholesky::luc19
auto luc19
Definition: choleskyInversion.h:184
math::cholesky::li65
auto li65
Definition: choleskyInversion.h:204
math::cholesky::li42
auto li42
Definition: choleskyInversion.h:87
math::cholesky::symmetrize55
constexpr void symmetrize55(M &dst)
Definition: choleskyInversion.h:270
math::cholesky::luc20
auto luc20
Definition: choleskyInversion.h:186
math::cholesky::luc11
auto luc11
Definition: choleskyInversion.h:119
math::cholesky::Inverter
Definition: choleskyInversion.h:289
math::cholesky::dst
constexpr void M2 & dst
Definition: choleskyInversion.h:23
TrackRefitter_38T_cff.src
src
Definition: TrackRefitter_38T_cff.py:24
math::cholesky::luc8
auto luc8
Definition: choleskyInversion.h:79
math::cholesky::li53
auto li53
Definition: choleskyInversion.h:133
math::cholesky::li61
auto li61
Definition: choleskyInversion.h:212
math::cholesky::luc7
auto luc7
Definition: choleskyInversion.h:78
math::cholesky::li64
auto li64
Definition: choleskyInversion.h:205
math::cholesky::Inverter::eval
static constexpr void eval(M1 const &src, M2 &dst)
Definition: choleskyInversion.h:290
math::cholesky::li31
auto li31
Definition: choleskyInversion.h:56
math::cholesky::symmetrize11
constexpr void symmetrize11(M &dst)
Definition: choleskyInversion.h:247
math::cholesky::li54
auto li54
Definition: choleskyInversion.h:132
math
Definition: choleskyInversion.h:19
l1tstage2_dqm_sourceclient-live_cfg.invert
invert
Definition: l1tstage2_dqm_sourceclient-live_cfg.py:86
math::cholesky::luc1
auto luc1
Definition: choleskyInversion.h:32
math::cholesky::li21
auto li21
Definition: choleskyInversion.h:35
math::cholesky::luc0
auto luc0
Definition: choleskyInversion.h:31
math::cholesky::luc4
auto luc4
Definition: choleskyInversion.h:50
math::cholesky::li62
auto li62
Definition: choleskyInversion.h:208
math::cholesky::luc18
auto luc18
Definition: choleskyInversion.h:183
math::cholesky::symmetrize22
constexpr void symmetrize22(M &dst)
Definition: choleskyInversion.h:250
math::cholesky::Inverter< M1, M2, 1 >::eval
static constexpr void eval(M1 const &src, M2 &dst)
Definition: choleskyInversion.h:295
math::cholesky::luc17
auto luc17
Definition: choleskyInversion.h:182
math::cholesky::symmetrize33
constexpr void symmetrize33(M &dst)
Definition: choleskyInversion.h:255
math::cholesky::li32
auto li32
Definition: choleskyInversion.h:55
math::cholesky::luc16
auto luc16
Definition: choleskyInversion.h:181
math::cholesky::luc5
auto luc5
Definition: choleskyInversion.h:51
math::cholesky::symmetrize44
constexpr void symmetrize44(M &dst)
Definition: choleskyInversion.h:262
math::cholesky::luc14
auto luc14
Definition: choleskyInversion.h:122