ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CholeskyDecomp.h
Go to the documentation of this file.
1 // @(#)root/smatrix:$Id$
2 // Author: M. Schiller 2009
3 
4 #ifndef ROOT_Math_CholeskyDecomp
5 #define ROOT_Math_CholeskyDecomp
6 
7 /** @file
8  * header file containing the templated implementation of matrix inversion
9  * routines for use with ROOT's SMatrix classes (symmetric positive
10  * definite case)
11  *
12  * @author Manuel Schiller
13  * @date Aug 29 2008
14  * initial release inside LHCb
15  * @date May 7 2009
16  * factored code to provide a nice Cholesky decomposition class, along
17  * with separate methods for solving a single linear system and to
18  * obtain the inverse matrix from the decomposition
19  * @date July 15th 2013
20  * provide a version of that class which works if the dimension of the
21  * problem is only known at run time
22  * @date September 30th 2013
23  * provide routines to access the result of the decomposition L and its
24  * inverse
25  */
26 
27 #include <cmath>
28 #include <algorithm>
29 
30 namespace ROOT {
31 
32  namespace Math {
33 
34 /// helpers for CholeskyDecomp
35 namespace CholeskyDecompHelpers {
36  // forward decls
37  template<class F, class M> struct _decomposerGenDim;
38  template<class F, unsigned N, class M> struct _decomposer;
39  template<class F, class M> struct _inverterGenDim;
40  template<class F, unsigned N, class M> struct _inverter;
41  template<class F, class V> struct _solverGenDim;
42  template<class F, unsigned N, class V> struct _solver;
43  template<typename G> class PackedArrayAdapter;
44 }
45 
46 /// class to compute the Cholesky decomposition of a matrix
47 /** class to compute the Cholesky decomposition of a symmetric
48  * positive definite matrix
49  *
50  * provides routines to check if the decomposition succeeded (i.e. if
51  * matrix is positive definite and non-singular), to solve a linear
52  * system for the given matrix and to obtain its inverse
53  *
54  * the actual functionality is implemented in templated helper
55  * classes which have specializations for dimensions N = 1 to 6
56  * to achieve a gain in speed for common matrix sizes
57  *
58  * usage example:
59  * @code
60  * // let m be a symmetric positive definite SMatrix (use type float
61  * // for internal computations, matrix size is 4x4)
62  * CholeskyDecomp<float, 4> decomp(m);
63  * // check if the decomposition succeeded
64  * if (!decomp) {
65  * std::cerr << "decomposition failed!" << std::endl;
66  * } else {
67  * // let rhs be a vector; we seek a vector x such that m * x = rhs
68  * decomp.Solve(rhs);
69  * // rhs now contains the solution we are looking for
70  *
71  * // obtain the inverse of m, put it into m itself
72  * decomp.Invert(m);
73  * }
74  * @endcode
75  */
76 template<class F, unsigned N> class CholeskyDecomp
77 {
78 private:
79  /// lower triangular matrix L
80  /** lower triangular matrix L, packed storage, with diagonal
81  * elements pre-inverted */
82  F fL[N * (N + 1) / 2];
83  /// flag indicating a successful decomposition
84  bool fOk;
85 public:
86  /// perform a Cholesky decomposition
87  /** perfrom a Cholesky decomposition of a symmetric positive
88  * definite matrix m
89  *
90  * this is the constructor to uses with an SMatrix (and objects
91  * that behave like an SMatrix in terms of using
92  * operator()(int i, int j) for access to elements)
93  */
94  template<class M> CholeskyDecomp(const M& m) :
95  fL( ), fOk(false)
96  {
98  fOk = _decomposer<F, N, M>()(fL, m);
99  }
100 
101  /// perform a Cholesky decomposition
102  /** perfrom a Cholesky decomposition of a symmetric positive
103  * definite matrix m
104  *
105  * this is the constructor to use in special applications where
106  * plain arrays are used
107  *
108  * NOTE: the matrix is given in packed representation, matrix
109  * element m(i,j) (j <= i) is supposed to be in array element
110  * (i * (i + 1)) / 2 + j
111  */
112  template<typename G> CholeskyDecomp(G* m) :
113  fL(), fOk(false)
114  {
117  fOk = _decomposer<F, N, PackedArrayAdapter<G> >()(
118  fL, PackedArrayAdapter<G>(m));
119  }
120 
121  /// returns true if decomposition was successful
122  /** @returns true if decomposition was successful */
123  bool ok() const { return fOk; }
124  /// returns true if decomposition was successful
125  /** @returns true if decomposition was successful */
126  operator bool() const { return fOk; }
127 
128  /** @brief solves a linear system for the given right hand side
129  *
130  * Note that you can use both SVector classes and plain arrays for
131  * rhs. (Make sure that the sizes match!). It will work with any vector
132  * implementing the operator [i]
133  *
134  * @returns if the decomposition was successful
135  */
136  template<class V> bool Solve(V& rhs) const
137  {
139  if (fOk) _solver<F,N,V>()(rhs, fL); return fOk;
140  }
141 
142  /** @brief place the inverse into m
143  *
144  * This is the method to use with an SMatrix.
145  *
146  * @returns if the decomposition was successful
147  */
148  template<class M> bool Invert(M& m) const
149  {
151  if (fOk) _inverter<F,N,M>()(m, fL); return fOk;
152  }
153 
154  /** @brief place the inverse into m
155  *
156  * This is the method to use with a plain array.
157  *
158  * @returns if the decomposition was successful
159  *
160  * NOTE: the matrix is given in packed representation, matrix
161  * element m(i,j) (j <= i) is supposed to be in array element
162  * (i * (i + 1)) / 2 + j
163  */
164  template<typename G> bool Invert(G* m) const
165  {
168  if (fOk) {
169  PackedArrayAdapter<G> adapted(m);
170  _inverter<F,N,PackedArrayAdapter<G> >()(adapted, fL);
171  }
172  return fOk;
173  }
174 
175  /** @brief obtain the decomposed matrix L
176  *
177  * This is the method to use with a plain array.
178  *
179  * @returns if the decomposition was successful
180  */
181  template<class M> bool getL(M& m) const
182  {
183  if (!fOk) return false;
184  for (unsigned i = 0; i < N; ++i) {
185  // zero upper half of matrix
186  for (unsigned j = i + 1; j < N; ++j)
187  m(i, j) = F(0);
188  // copy the rest
189  for (unsigned j = 0; j <= i; ++j)
190  m(i, j) = fL[i * (i + 1) / 2 + j];
191  // adjust the diagonal - we save 1/L(i, i) in that position, so
192  // convert to what caller expects
193  m(i, i) = F(1) / m(i, i);
194  }
195  return true;
196  }
197 
198  /** @brief obtain the decomposed matrix L
199  *
200  * @returns if the decomposition was successful
201  *
202  * NOTE: the matrix is given in packed representation, matrix
203  * element m(i,j) (j <= i) is supposed to be in array element
204  * (i * (i + 1)) / 2 + j
205  */
206  template<typename G> bool getL(G* m) const
207  {
208  if (!fOk) return false;
209  // copy L
210  for (unsigned i = 0; i < (N * (N + 1)) / 2; ++i)
211  m[i] = fL[i];
212  // adjust diagonal - we save 1/L(i, i) in that position, so convert to
213  // what caller expects
214  for (unsigned i = 0; i < N; ++i)
215  m[(i * (i + 1)) / 2 + i] = F(1) / fL[(i * (i + 1)) / 2 + i];
216  return true;
217  }
218 
219  /** @brief obtain the inverse of the decomposed matrix L
220  *
221  * This is the method to use with a plain array.
222  *
223  * @returns if the decomposition was successful
224  */
225  template<class M> bool getLi(M& m) const
226  {
227  if (!fOk) return false;
228  for (unsigned i = 0; i < N; ++i) {
229  // zero lower half of matrix
230  for (unsigned j = i + 1; j < N; ++j)
231  m(j, i) = F(0);
232  // copy the rest
233  for (unsigned j = 0; j <= i; ++j)
234  m(j, i) = fL[i * (i + 1) / 2 + j];
235  }
236  // invert the off-diagonal part of what we just copied
237  for (unsigned i = 1; i < N; ++i) {
238  for (unsigned j = 0; j < i; ++j) {
239  typename M::value_type tmp = F(0);
240  for (unsigned k = i; k-- > j;)
241  tmp -= m(k, i) * m(j, k);
242  m(j, i) = tmp * m(i, i);
243  }
244  }
245  return true;
246  }
247 
248  /** @brief obtain the inverse of the decomposed matrix L
249  *
250  * @returns if the decomposition was successful
251  *
252  * NOTE: the matrix is given in packed representation, matrix
253  * element m(j,i) (j <= i) is supposed to be in array element
254  * (i * (i + 1)) / 2 + j
255  */
256  template<typename G> bool getLi(G* m) const
257  {
258  if (!fOk) return false;
259  // copy L
260  for (unsigned i = 0; i < (N * (N + 1)) / 2; ++i)
261  m[i] = fL[i];
262  // invert the off-diagonal part of what we just copied
263  G* base1 = &m[1];
264  for (unsigned i = 1; i < N; base1 += ++i) {
265  for (unsigned j = 0; j < i; ++j) {
266  G tmp = F(0);
267  const G *base2 = &m[(i * (i - 1)) / 2];
268  for (unsigned k = i; k-- > j; base2 -= k)
269  tmp -= base1[k] * base2[j];
270  base1[j] = tmp * base1[i];
271  }
272  }
273  return true;
274  }
275 };
276 
277 /// class to compute the Cholesky decomposition of a matrix
278 /** class to compute the Cholesky decomposition of a symmetric
279  * positive definite matrix when the dimensionality of the problem is not known
280  * at compile time
281  *
282  * provides routines to check if the decomposition succeeded (i.e. if
283  * matrix is positive definite and non-singular), to solve a linear
284  * system for the given matrix and to obtain its inverse
285  *
286  * the actual functionality is implemented in templated helper
287  * classes which have specializations for dimensions N = 1 to 6
288  * to achieve a gain in speed for common matrix sizes
289  *
290  * usage example:
291  * @code
292  * // let m be a symmetric positive definite SMatrix (use type float
293  * // for internal computations, matrix size is 4x4)
294  * CholeskyDecompGenDim<float> decomp(4, m);
295  * // check if the decomposition succeeded
296  * if (!decomp) {
297  * std::cerr << "decomposition failed!" << std::endl;
298  * } else {
299  * // let rhs be a vector; we seek a vector x such that m * x = rhs
300  * decomp.Solve(rhs);
301  * // rhs now contains the solution we are looking for
302  *
303  * // obtain the inverse of m, put it into m itself
304  * decomp.Invert(m);
305  * }
306  * @endcode
307  */
308 template<class F> class CholeskyDecompGenDim
309 {
310 private:
311  /** @brief dimensionality
312  * dimensionality of the problem */
313  unsigned fN;
314  /// lower triangular matrix L
315  /** lower triangular matrix L, packed storage, with diagonal
316  * elements pre-inverted */
317  F *fL;
318  /// flag indicating a successful decomposition
319  bool fOk;
320 public:
321  /// perform a Cholesky decomposition
322  /** perfrom a Cholesky decomposition of a symmetric positive
323  * definite matrix m
324  *
325  * this is the constructor to uses with an SMatrix (and objects
326  * that behave like an SMatrix in terms of using
327  * operator()(int i, int j) for access to elements)
328  */
329  template<class M> CholeskyDecompGenDim(unsigned N, const M& m) :
330  fN(N), fL(new F[(fN * (fN + 1)) / 2]), fOk(false)
331  {
333  fOk = _decomposerGenDim<F, M>()(fL, m, fN);
334  }
335 
336  /// perform a Cholesky decomposition
337  /** perfrom a Cholesky decomposition of a symmetric positive
338  * definite matrix m
339  *
340  * this is the constructor to use in special applications where
341  * plain arrays are used
342  *
343  * NOTE: the matrix is given in packed representation, matrix
344  * element m(i,j) (j <= i) is supposed to be in array element
345  * (i * (i + 1)) / 2 + j
346  */
347  template<typename G> CholeskyDecompGenDim(unsigned N, G* m) :
348  fN(N), fL(new F[(fN * (fN + 1)) / 2]), fOk(false)
349  {
352  fOk = _decomposerGenDim<F, PackedArrayAdapter<G> >()(
353  fL, PackedArrayAdapter<G>(m), fN);
354  }
355 
356  /// destructor
357  ~CholeskyDecompGenDim() { delete[] fL; }
358 
359  /// returns true if decomposition was successful
360  /** @returns true if decomposition was successful */
361  bool ok() const { return fOk; }
362  /// returns true if decomposition was successful
363  /** @returns true if decomposition was successful */
364  operator bool() const { return fOk; }
365 
366  /** @brief solves a linear system for the given right hand side
367  *
368  * Note that you can use both SVector classes and plain arrays for
369  * rhs. (Make sure that the sizes match!). It will work with any vector
370  * implementing the operator [i]
371  *
372  * @returns if the decomposition was successful
373  */
374  template<class V> bool Solve(V& rhs) const
375  {
377  if (fOk) _solverGenDim<F,V>()(rhs, fL, fN); return fOk;
378  }
379 
380  /** @brief place the inverse into m
381  *
382  * This is the method to use with an SMatrix.
383  *
384  * @returns if the decomposition was successful
385  */
386  template<class M> bool Invert(M& m) const
387  {
389  if (fOk) _inverterGenDim<F,M>()(m, fL, fN); return fOk;
390  }
391 
392  /** @brief place the inverse into m
393  *
394  * This is the method to use with a plain array.
395  *
396  * @returns if the decomposition was successful
397  *
398  * NOTE: the matrix is given in packed representation, matrix
399  * element m(i,j) (j <= i) is supposed to be in array element
400  * (i * (i + 1)) / 2 + j
401  */
402  template<typename G> bool Invert(G* m) const
403  {
406  if (fOk) {
407  PackedArrayAdapter<G> adapted(m);
408  _inverterGenDim<F,PackedArrayAdapter<G> >()(adapted, fL, fN);
409  }
410  return fOk;
411  }
412 
413  /** @brief obtain the decomposed matrix L
414  *
415  * This is the method to use with a plain array.
416  *
417  * @returns if the decomposition was successful
418  */
419  template<class M> bool getL(M& m) const
420  {
421  if (!fOk) return false;
422  for (unsigned i = 0; i < fN; ++i) {
423  // zero upper half of matrix
424  for (unsigned j = i + 1; j < fN; ++j)
425  m(i, j) = F(0);
426  // copy the rest
427  for (unsigned j = 0; j <= i; ++j)
428  m(i, j) = fL[i * (i + 1) / 2 + j];
429  // adjust the diagonal - we save 1/L(i, i) in that position, so
430  // convert to what caller expects
431  m(i, i) = F(1) / m(i, i);
432  }
433  return true;
434  }
435 
436  /** @brief obtain the decomposed matrix L
437  *
438  * @returns if the decomposition was successful
439  *
440  * NOTE: the matrix is given in packed representation, matrix
441  * element m(i,j) (j <= i) is supposed to be in array element
442  * (i * (i + 1)) / 2 + j
443  */
444  template<typename G> bool getL(G* m) const
445  {
446  if (!fOk) return false;
447  // copy L
448  for (unsigned i = 0; i < (fN * (fN + 1)) / 2; ++i)
449  m[i] = fL[i];
450  // adjust diagonal - we save 1/L(i, i) in that position, so convert to
451  // what caller expects
452  for (unsigned i = 0; i < fN; ++i)
453  m[(i * (i + 1)) / 2 + i] = F(1) / fL[(i * (i + 1)) / 2 + i];
454  return true;
455  }
456 
457  /** @brief obtain the inverse of the decomposed matrix L
458  *
459  * This is the method to use with a plain array.
460  *
461  * @returns if the decomposition was successful
462  */
463  template<class M> bool getLi(M& m) const
464  {
465  if (!fOk) return false;
466  for (unsigned i = 0; i < fN; ++i) {
467  // zero lower half of matrix
468  for (unsigned j = i + 1; j < fN; ++j)
469  m(j, i) = F(0);
470  // copy the rest
471  for (unsigned j = 0; j <= i; ++j)
472  m(j, i) = fL[i * (i + 1) / 2 + j];
473  }
474  // invert the off-diagonal part of what we just copied
475  for (unsigned i = 1; i < fN; ++i) {
476  for (unsigned j = 0; j < i; ++j) {
477  typename M::value_type tmp = F(0);
478  for (unsigned k = i; k-- > j;)
479  tmp -= m(k, i) * m(j, k);
480  m(j, i) = tmp * m(i, i);
481  }
482  }
483  return true;
484  }
485 
486  /** @brief obtain the inverse of the decomposed matrix L
487  *
488  * @returns if the decomposition was successful
489  *
490  * NOTE: the matrix is given in packed representation, matrix
491  * element m(j,i) (j <= i) is supposed to be in array element
492  * (i * (i + 1)) / 2 + j
493  */
494  template<typename G> bool getLi(G* m) const
495  {
496  if (!fOk) return false;
497  // copy L
498  for (unsigned i = 0; i < (fN * (fN + 1)) / 2; ++i)
499  m[i] = fL[i];
500  // invert the off-diagonal part of what we just copied
501  G* base1 = &m[1];
502  for (unsigned i = 1; i < fN; base1 += ++i) {
503  for (unsigned j = 0; j < i; ++j) {
504  G tmp = F(0);
505  const G *base2 = &m[(i * (i - 1)) / 2];
506  for (unsigned k = i; k-- > j; base2 -= k)
507  tmp -= base1[k] * base2[j];
508  base1[j] = tmp * base1[i];
509  }
510  }
511  return true;
512  }
513 };
514 
515 namespace CholeskyDecompHelpers {
516  /// adapter for packed arrays (to SMatrix indexing conventions)
517  template<typename G> class PackedArrayAdapter
518  {
519  private:
520  G* fArr; ///< pointer to first array element
521  public:
522  /// constructor
523  PackedArrayAdapter(G* arr) : fArr(arr) {}
524  /// read access to elements (make sure that j <= i)
525  const G operator()(unsigned i, unsigned j) const
526  { return fArr[((i * (i + 1)) / 2) + j]; }
527  /// write access to elements (make sure that j <= i)
528  G& operator()(unsigned i, unsigned j)
529  { return fArr[((i * (i + 1)) / 2) + j]; }
530  };
531  /// struct to do a Cholesky decomposition (general dimensionality)
532  template<class F, class M> struct _decomposerGenDim
533  {
534  /// method to do the decomposition
535  /** @returns if the decomposition was successful */
536  bool operator()(F* dst, const M& src, unsigned N) const
537  {
538  // perform Cholesky decomposition of matrix: M = L L^T
539  // only thing that can go wrong: trying to take square
540  // root of negative number or zero (matrix is
541  // ill-conditioned or singular in these cases)
542 
543  // element L(i,j) is at array position (i * (i+1)) / 2 + j
544 
545  // quirk: we may need to invert L later anyway, so we can
546  // invert elements on diagonale straight away (we only
547  // ever need their reciprocals!)
548 
549  // cache starting address of rows of L for speed reasons
550  F *base1 = &dst[0];
551  for (unsigned i = 0; i < N; base1 += ++i) {
552  F tmpdiag = F(0.0); // for element on diagonale
553  // calculate off-diagonal elements
554  F *base2 = &dst[0];
555  for (unsigned j = 0; j < i; base2 += ++j) {
556  F tmp = src(i, j);
557  for (unsigned k = j; k--; )
558  tmp -= base1[k] * base2[k];
559  base1[j] = tmp *= base2[j];
560  // keep track of contribution to element on diagonale
561  tmpdiag += tmp * tmp;
562  }
563  // keep truncation error small
564  tmpdiag = src(i, i) - tmpdiag;
565  // check if positive definite
566  if (tmpdiag <= F(0.0)) return false;
567  else base1[i] = std::sqrt(F(1.0) / tmpdiag);
568  }
569  return true;
570  }
571  };
572 
573  /// struct to do a Cholesky decomposition
574  template<class F, unsigned N, class M> struct _decomposer
575  {
576  /// method to do the decomposition
577  /** @returns if the decomposition was successful */
578  bool operator()(F* dst, const M& src) const
579  { return _decomposerGenDim<F, M>()(dst, src, N); }
580  };
581 
582  /// struct to obtain the inverse from a Cholesky decomposition (general dimensionality)
583  template<class F, class M> struct _inverterGenDim
584  {
585  /// method to do the inversion
586  void operator()(M& dst, const F* src, unsigned N) const
587  {
588  // make working copy
589  F * l = new F[N * (N + 1) / 2];
590  std::copy(src, src + ((N * (N + 1)) / 2), l);
591  // ok, next step: invert off-diagonal part of matrix
592  F* base1 = &l[1];
593  for (unsigned i = 1; i < N; base1 += ++i) {
594  for (unsigned j = 0; j < i; ++j) {
595  F tmp = F(0.0);
596  const F *base2 = &l[(i * (i - 1)) / 2];
597  for (unsigned k = i; k-- > j; base2 -= k)
598  tmp -= base1[k] * base2[j];
599  base1[j] = tmp * base1[i];
600  }
601  }
602 
603  // Li = L^(-1) formed, now calculate M^(-1) = Li^T Li
604  for (unsigned i = N; i--; ) {
605  for (unsigned j = i + 1; j--; ) {
606  F tmp = F(0.0);
607  base1 = &l[(N * (N - 1)) / 2];
608  for (unsigned k = N; k-- > i; base1 -= k)
609  tmp += base1[i] * base1[j];
610  dst(i, j) = tmp;
611  }
612  }
613  delete [] l;
614  }
615  };
616 
617  /// struct to obtain the inverse from a Cholesky decomposition
618  template<class F, unsigned N, class M> struct _inverter
619  {
620  /// method to do the inversion
621  void operator()(M& dst, const F* src) const
622  { return _inverterGenDim<F, M>()(dst, src, N); }
623  };
624 
625  /// struct to solve a linear system using its Cholesky decomposition (generalised dimensionality)
626  template<class F, class V> struct _solverGenDim
627  {
628  /// method to solve the linear system
629  void operator()(V& rhs, const F* l, unsigned N) const
630  {
631  // solve Ly = rhs
632  for (unsigned k = 0; k < N; ++k) {
633  const unsigned base = (k * (k + 1)) / 2;
634  F sum = F(0.0);
635  for (unsigned i = k; i--; )
636  sum += rhs[i] * l[base + i];
637  // elements on diagonale are pre-inverted!
638  rhs[k] = (rhs[k] - sum) * l[base + k];
639  }
640  // solve L^Tx = y
641  for (unsigned k = N; k--; ) {
642  F sum = F(0.0);
643  for (unsigned i = N; --i > k; )
644  sum += rhs[i] * l[(i * (i + 1)) / 2 + k];
645  // elements on diagonale are pre-inverted!
646  rhs[k] = (rhs[k] - sum) * l[(k * (k + 1)) / 2 + k];
647  }
648  }
649  };
650 
651  /// struct to solve a linear system using its Cholesky decomposition
652  template<class F, unsigned N, class V> struct _solver
653  {
654  /// method to solve the linear system
655  void operator()(V& rhs, const F* l) const
656  { _solverGenDim<F, V>()(rhs, l, N); }
657  };
658 
659  /// struct to do a Cholesky decomposition (specialized, N = 6)
660  template<class F, class M> struct _decomposer<F, 6, M>
661  {
662  /// method to do the decomposition
663  bool operator()(F* dst, const M& src) const
664  {
665  if (src(0,0) <= F(0.0)) return false;
666  dst[0] = std::sqrt(F(1.0) / src(0,0));
667  dst[1] = src(1,0) * dst[0];
668  dst[2] = src(1,1) - dst[1] * dst[1];
669  if (dst[2] <= F(0.0)) return false;
670  else dst[2] = std::sqrt(F(1.0) / dst[2]);
671  dst[3] = src(2,0) * dst[0];
672  dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
673  dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
674  if (dst[5] <= F(0.0)) return false;
675  else dst[5] = std::sqrt(F(1.0) / dst[5]);
676  dst[6] = src(3,0) * dst[0];
677  dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
678  dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
679  dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
680  if (dst[9] <= F(0.0)) return false;
681  else dst[9] = std::sqrt(F(1.0) / dst[9]);
682  dst[10] = src(4,0) * dst[0];
683  dst[11] = (src(4,1) - dst[1] * dst[10]) * dst[2];
684  dst[12] = (src(4,2) - dst[3] * dst[10] - dst[4] * dst[11]) * dst[5];
685  dst[13] = (src(4,3) - dst[6] * dst[10] - dst[7] * dst[11] - dst[8] * dst[12]) * dst[9];
686  dst[14] = src(4,4) - (dst[10]*dst[10]+dst[11]*dst[11]+dst[12]*dst[12]+dst[13]*dst[13]);
687  if (dst[14] <= F(0.0)) return false;
688  else dst[14] = std::sqrt(F(1.0) / dst[14]);
689  dst[15] = src(5,0) * dst[0];
690  dst[16] = (src(5,1) - dst[1] * dst[15]) * dst[2];
691  dst[17] = (src(5,2) - dst[3] * dst[15] - dst[4] * dst[16]) * dst[5];
692  dst[18] = (src(5,3) - dst[6] * dst[15] - dst[7] * dst[16] - dst[8] * dst[17]) * dst[9];
693  dst[19] = (src(5,4) - dst[10] * dst[15] - dst[11] * dst[16] - dst[12] * dst[17] - dst[13] * dst[18]) * dst[14];
694  dst[20] = src(5,5) - (dst[15]*dst[15]+dst[16]*dst[16]+dst[17]*dst[17]+dst[18]*dst[18]+dst[19]*dst[19]);
695  if (dst[20] <= F(0.0)) return false;
696  else dst[20] = std::sqrt(F(1.0) / dst[20]);
697  return true;
698  }
699  };
700  /// struct to do a Cholesky decomposition (specialized, N = 5)
701  template<class F, class M> struct _decomposer<F, 5, M>
702  {
703  /// method to do the decomposition
704  bool operator()(F* dst, const M& src) const
705  {
706  if (src(0,0) <= F(0.0)) return false;
707  dst[0] = std::sqrt(F(1.0) / src(0,0));
708  dst[1] = src(1,0) * dst[0];
709  dst[2] = src(1,1) - dst[1] * dst[1];
710  if (dst[2] <= F(0.0)) return false;
711  else dst[2] = std::sqrt(F(1.0) / dst[2]);
712  dst[3] = src(2,0) * dst[0];
713  dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
714  dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
715  if (dst[5] <= F(0.0)) return false;
716  else dst[5] = std::sqrt(F(1.0) / dst[5]);
717  dst[6] = src(3,0) * dst[0];
718  dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
719  dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
720  dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
721  if (dst[9] <= F(0.0)) return false;
722  else dst[9] = std::sqrt(F(1.0) / dst[9]);
723  dst[10] = src(4,0) * dst[0];
724  dst[11] = (src(4,1) - dst[1] * dst[10]) * dst[2];
725  dst[12] = (src(4,2) - dst[3] * dst[10] - dst[4] * dst[11]) * dst[5];
726  dst[13] = (src(4,3) - dst[6] * dst[10] - dst[7] * dst[11] - dst[8] * dst[12]) * dst[9];
727  dst[14] = src(4,4) - (dst[10]*dst[10]+dst[11]*dst[11]+dst[12]*dst[12]+dst[13]*dst[13]);
728  if (dst[14] <= F(0.0)) return false;
729  else dst[14] = std::sqrt(F(1.0) / dst[14]);
730  return true;
731  }
732  };
733  /// struct to do a Cholesky decomposition (specialized, N = 4)
734  template<class F, class M> struct _decomposer<F, 4, M>
735  {
736  /// method to do the decomposition
737  bool operator()(F* dst, const M& src) const
738  {
739  if (src(0,0) <= F(0.0)) return false;
740  dst[0] = std::sqrt(F(1.0) / src(0,0));
741  dst[1] = src(1,0) * dst[0];
742  dst[2] = src(1,1) - dst[1] * dst[1];
743  if (dst[2] <= F(0.0)) return false;
744  else dst[2] = std::sqrt(F(1.0) / dst[2]);
745  dst[3] = src(2,0) * dst[0];
746  dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
747  dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
748  if (dst[5] <= F(0.0)) return false;
749  else dst[5] = std::sqrt(F(1.0) / dst[5]);
750  dst[6] = src(3,0) * dst[0];
751  dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
752  dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
753  dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
754  if (dst[9] <= F(0.0)) return false;
755  else dst[9] = std::sqrt(F(1.0) / dst[9]);
756  return true;
757  }
758  };
759  /// struct to do a Cholesky decomposition (specialized, N = 3)
760  template<class F, class M> struct _decomposer<F, 3, M>
761  {
762  /// method to do the decomposition
763  bool operator()(F* dst, const M& src) const
764  {
765  if (src(0,0) <= F(0.0)) return false;
766  dst[0] = std::sqrt(F(1.0) / src(0,0));
767  dst[1] = src(1,0) * dst[0];
768  dst[2] = src(1,1) - dst[1] * dst[1];
769  if (dst[2] <= F(0.0)) return false;
770  else dst[2] = std::sqrt(F(1.0) / dst[2]);
771  dst[3] = src(2,0) * dst[0];
772  dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
773  dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
774  if (dst[5] <= F(0.0)) return false;
775  else dst[5] = std::sqrt(F(1.0) / dst[5]);
776  return true;
777  }
778  };
779  /// struct to do a Cholesky decomposition (specialized, N = 2)
780  template<class F, class M> struct _decomposer<F, 2, M>
781  {
782  /// method to do the decomposition
783  bool operator()(F* dst, const M& src) const
784  {
785  if (src(0,0) <= F(0.0)) return false;
786  dst[0] = std::sqrt(F(1.0) / src(0,0));
787  dst[1] = src(1,0) * dst[0];
788  dst[2] = src(1,1) - dst[1] * dst[1];
789  if (dst[2] <= F(0.0)) return false;
790  else dst[2] = std::sqrt(F(1.0) / dst[2]);
791  return true;
792  }
793  };
794  /// struct to do a Cholesky decomposition (specialized, N = 1)
795  template<class F, class M> struct _decomposer<F, 1, M>
796  {
797  /// method to do the decomposition
798  bool operator()(F* dst, const M& src) const
799  {
800  if (src(0,0) <= F(0.0)) return false;
801  dst[0] = std::sqrt(F(1.0) / src(0,0));
802  return true;
803  }
804  };
805  /// struct to do a Cholesky decomposition (specialized, N = 0)
806  template<class F, class M> struct _decomposer<F, 0, M>
807  {
808  private:
809  _decomposer() { };
810  bool operator()(F* dst, const M& src) const;
811  };
812 
813  /// struct to obtain the inverse from a Cholesky decomposition (N = 6)
814  template<class F, class M> struct _inverter<F,6,M>
815  {
816  /// method to do the inversion
817  inline void operator()(M& dst, const F* src) const
818  {
819  const F li21 = -src[1] * src[0] * src[2];
820  const F li32 = -src[4] * src[2] * src[5];
821  const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
822  const F li43 = -src[8] * src[9] * src[5];
823  const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
824  const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
825  src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
826  const F li54 = -src[13] * src[14] * src[9];
827  const F li53 = (src[13] * src[8] * src[9] - src[12]) * src[5] * src[14];
828  const F li52 = (-src[4] * src[8] * src[13] * src[5] * src[9] +
829  src[4] * src[12] * src[5] + src[7] * src[13] * src[9] - src[11]) * src[2] * src[14];
830  const F li51 = (src[1]*src[4]*src[8]*src[13]*src[2]*src[5]*src[9] -
831  src[13]*src[8]*src[3]*src[9]*src[5] - src[12]*src[4]*src[1]*src[2]*src[5] - src[13]*src[7]*src[1]*src[9]*src[2] +
832  src[11]*src[1]*src[2] + src[12]*src[3]*src[5] + src[13]*src[6]*src[9] -src[10]) * src[0] * src[14];
833  const F li65 = -src[19] * src[20] * src[14];
834  const F li64 = (src[19] * src[13] * src[14] - src[18]) * src[9] * src[20];
835  const F li63 = (-src[8] * src[13] * src[19] * src[9] * src[14] +
836  src[8] * src[18] * src[9] + src[12] * src[19] * src[14] - src[17]) * src[5] * src[20];
837  const F li62 = (src[4]*src[8]*src[13]*src[19]*src[5]*src[9]*src[14] -
838  src[18]*src[8]*src[4]*src[9]*src[5] - src[19]*src[12]*src[4]*src[14]*src[5] -src[19]*src[13]*src[7]*src[14]*src[9] +
839  src[17]*src[4]*src[5] + src[18]*src[7]*src[9] + src[19]*src[11]*src[14] - src[16]) * src[2] * src[20];
840  const F li61 = (-src[19]*src[13]*src[8]*src[4]*src[1]*src[2]*src[5]*src[9]*src[14] +
841  src[18]*src[8]*src[4]*src[1]*src[2]*src[5]*src[9] + src[19]*src[12]*src[4]*src[1]*src[2]*src[5]*src[14] +
842  src[19]*src[13]*src[7]*src[1]*src[2]*src[9]*src[14] + src[19]*src[13]*src[8]*src[3]*src[5]*src[9]*src[14] -
843  src[17]*src[4]*src[1]*src[2]*src[5] - src[18]*src[7]*src[1]*src[2]*src[9] - src[19]*src[11]*src[1]*src[2]*src[14] -
844  src[18]*src[8]*src[3]*src[5]*src[9] - src[19]*src[12]*src[3]*src[5]*src[14] - src[19]*src[13]*src[6]*src[9]*src[14] +
845  src[16]*src[1]*src[2] + src[17]*src[3]*src[5] + src[18]*src[6]*src[9] + src[19]*src[10]*src[14] - src[15]) *
846  src[0] * src[20];
847 
848  dst(0,0) = li61*li61 + li51*li51 + li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
849  dst(1,0) = li61*li62 + li51*li52 + li41*li42 + li31*li32 + li21*src[2];
850  dst(1,1) = li62*li62 + li52*li52 + li42*li42 + li32*li32 + src[2]*src[2];
851  dst(2,0) = li61*li63 + li51*li53 + li41*li43 + li31*src[5];
852  dst(2,1) = li62*li63 + li52*li53 + li42*li43 + li32*src[5];
853  dst(2,2) = li63*li63 + li53*li53 + li43*li43 + src[5]*src[5];
854  dst(3,0) = li61*li64 + li51*li54 + li41*src[9];
855  dst(3,1) = li62*li64 + li52*li54 + li42*src[9];
856  dst(3,2) = li63*li64 + li53*li54 + li43*src[9];
857  dst(3,3) = li64*li64 + li54*li54 + src[9]*src[9];
858  dst(4,0) = li61*li65 + li51*src[14];
859  dst(4,1) = li62*li65 + li52*src[14];
860  dst(4,2) = li63*li65 + li53*src[14];
861  dst(4,3) = li64*li65 + li54*src[14];
862  dst(4,4) = li65*li65 + src[14]*src[14];
863  dst(5,0) = li61*src[20];
864  dst(5,1) = li62*src[20];
865  dst(5,2) = li63*src[20];
866  dst(5,3) = li64*src[20];
867  dst(5,4) = li65*src[20];
868  dst(5,5) = src[20]*src[20];
869  }
870  };
871  /// struct to obtain the inverse from a Cholesky decomposition (N = 5)
872  template<class F, class M> struct _inverter<F,5,M>
873  {
874  /// method to do the inversion
875  inline void operator()(M& dst, const F* src) const
876  {
877  const F li21 = -src[1] * src[0] * src[2];
878  const F li32 = -src[4] * src[2] * src[5];
879  const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
880  const F li43 = -src[8] * src[9] * src[5];
881  const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
882  const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
883  src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
884  const F li54 = -src[13] * src[14] * src[9];
885  const F li53 = (src[13] * src[8] * src[9] - src[12]) * src[5] * src[14];
886  const F li52 = (-src[4] * src[8] * src[13] * src[5] * src[9] +
887  src[4] * src[12] * src[5] + src[7] * src[13] * src[9] - src[11]) * src[2] * src[14];
888  const F li51 = (src[1]*src[4]*src[8]*src[13]*src[2]*src[5]*src[9] -
889  src[13]*src[8]*src[3]*src[9]*src[5] - src[12]*src[4]*src[1]*src[2]*src[5] - src[13]*src[7]*src[1]*src[9]*src[2] +
890  src[11]*src[1]*src[2] + src[12]*src[3]*src[5] + src[13]*src[6]*src[9] -src[10]) * src[0] * src[14];
891 
892  dst(0,0) = li51*li51 + li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
893  dst(1,0) = li51*li52 + li41*li42 + li31*li32 + li21*src[2];
894  dst(1,1) = li52*li52 + li42*li42 + li32*li32 + src[2]*src[2];
895  dst(2,0) = li51*li53 + li41*li43 + li31*src[5];
896  dst(2,1) = li52*li53 + li42*li43 + li32*src[5];
897  dst(2,2) = li53*li53 + li43*li43 + src[5]*src[5];
898  dst(3,0) = li51*li54 + li41*src[9];
899  dst(3,1) = li52*li54 + li42*src[9];
900  dst(3,2) = li53*li54 + li43*src[9];
901  dst(3,3) = li54*li54 + src[9]*src[9];
902  dst(4,0) = li51*src[14];
903  dst(4,1) = li52*src[14];
904  dst(4,2) = li53*src[14];
905  dst(4,3) = li54*src[14];
906  dst(4,4) = src[14]*src[14];
907  }
908  };
909  /// struct to obtain the inverse from a Cholesky decomposition (N = 4)
910  template<class F, class M> struct _inverter<F,4,M>
911  {
912  /// method to do the inversion
913  inline void operator()(M& dst, const F* src) const
914  {
915  const F li21 = -src[1] * src[0] * src[2];
916  const F li32 = -src[4] * src[2] * src[5];
917  const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
918  const F li43 = -src[8] * src[9] * src[5];
919  const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
920  const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
921  src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
922 
923  dst(0,0) = li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
924  dst(1,0) = li41*li42 + li31*li32 + li21*src[2];
925  dst(1,1) = li42*li42 + li32*li32 + src[2]*src[2];
926  dst(2,0) = li41*li43 + li31*src[5];
927  dst(2,1) = li42*li43 + li32*src[5];
928  dst(2,2) = li43*li43 + src[5]*src[5];
929  dst(3,0) = li41*src[9];
930  dst(3,1) = li42*src[9];
931  dst(3,2) = li43*src[9];
932  dst(3,3) = src[9]*src[9];
933  }
934  };
935  /// struct to obtain the inverse from a Cholesky decomposition (N = 3)
936  template<class F, class M> struct _inverter<F,3,M>
937  {
938  /// method to do the inversion
939  inline void operator()(M& dst, const F* src) const
940  {
941  const F li21 = -src[1] * src[0] * src[2];
942  const F li32 = -src[4] * src[2] * src[5];
943  const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
944 
945  dst(0,0) = li31*li31 + li21*li21 + src[0]*src[0];
946  dst(1,0) = li31*li32 + li21*src[2];
947  dst(1,1) = li32*li32 + src[2]*src[2];
948  dst(2,0) = li31*src[5];
949  dst(2,1) = li32*src[5];
950  dst(2,2) = src[5]*src[5];
951  }
952  };
953  /// struct to obtain the inverse from a Cholesky decomposition (N = 2)
954  template<class F, class M> struct _inverter<F,2,M>
955  {
956  /// method to do the inversion
957  inline void operator()(M& dst, const F* src) const
958  {
959  const F li21 = -src[1] * src[0] * src[2];
960 
961  dst(0,0) = li21*li21 + src[0]*src[0];
962  dst(1,0) = li21*src[2];
963  dst(1,1) = src[2]*src[2];
964  }
965  };
966  /// struct to obtain the inverse from a Cholesky decomposition (N = 1)
967  template<class F, class M> struct _inverter<F,1,M>
968  {
969  /// method to do the inversion
970  inline void operator()(M& dst, const F* src) const
971  {
972  dst(0,0) = src[0]*src[0];
973  }
974  };
975  /// struct to obtain the inverse from a Cholesky decomposition (N = 0)
976  template<class F, class M> struct _inverter<F,0,M>
977  {
978  private:
979  _inverter();
980  void operator()(M& dst, const F* src) const;
981  };
982 
983  /// struct to solve a linear system using its Cholesky decomposition (N=6)
984  template<class F, class V> struct _solver<F,6,V>
985  {
986  /// method to solve the linear system
987  void operator()(V& rhs, const F* l) const
988  {
989  // solve Ly = rhs
990  const F y0 = rhs[0] * l[0];
991  const F y1 = (rhs[1]-l[1]*y0)*l[2];
992  const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
993  const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
994  const F y4 = (rhs[4]-(l[10]*y0+l[11]*y1+l[12]*y2+l[13]*y3))*l[14];
995  const F y5 = (rhs[5]-(l[15]*y0+l[16]*y1+l[17]*y2+l[18]*y3+l[19]*y4))*l[20];
996  // solve L^Tx = y, and put x into rhs
997  rhs[5] = y5 * l[20];
998  rhs[4] = (y4-l[19]*rhs[5])*l[14];
999  rhs[3] = (y3-(l[18]*rhs[5]+l[13]*rhs[4]))*l[9];
1000  rhs[2] = (y2-(l[17]*rhs[5]+l[12]*rhs[4]+l[8]*rhs[3]))*l[5];
1001  rhs[1] = (y1-(l[16]*rhs[5]+l[11]*rhs[4]+l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1002  rhs[0] = (y0-(l[15]*rhs[5]+l[10]*rhs[4]+l[6]*rhs[3]+l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1003  }
1004  };
1005  /// struct to solve a linear system using its Cholesky decomposition (N=5)
1006  template<class F, class V> struct _solver<F,5,V>
1007  {
1008  /// method to solve the linear system
1009  void operator()(V& rhs, const F* l) const
1010  {
1011  // solve Ly = rhs
1012  const F y0 = rhs[0] * l[0];
1013  const F y1 = (rhs[1]-l[1]*y0)*l[2];
1014  const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1015  const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
1016  const F y4 = (rhs[4]-(l[10]*y0+l[11]*y1+l[12]*y2+l[13]*y3))*l[14];
1017  // solve L^Tx = y, and put x into rhs
1018  rhs[4] = (y4)*l[14];
1019  rhs[3] = (y3-(l[13]*rhs[4]))*l[9];
1020  rhs[2] = (y2-(l[12]*rhs[4]+l[8]*rhs[3]))*l[5];
1021  rhs[1] = (y1-(l[11]*rhs[4]+l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1022  rhs[0] = (y0-(l[10]*rhs[4]+l[6]*rhs[3]+l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1023  }
1024  };
1025  /// struct to solve a linear system using its Cholesky decomposition (N=4)
1026  template<class F, class V> struct _solver<F,4,V>
1027  {
1028  /// method to solve the linear system
1029  void operator()(V& rhs, const F* l) const
1030  {
1031  // solve Ly = rhs
1032  const F y0 = rhs[0] * l[0];
1033  const F y1 = (rhs[1]-l[1]*y0)*l[2];
1034  const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1035  const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
1036  // solve L^Tx = y, and put x into rhs
1037  rhs[3] = (y3)*l[9];
1038  rhs[2] = (y2-(l[8]*rhs[3]))*l[5];
1039  rhs[1] = (y1-(l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1040  rhs[0] = (y0-(l[6]*rhs[3]+l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1041  }
1042  };
1043  /// struct to solve a linear system using its Cholesky decomposition (N=3)
1044  template<class F, class V> struct _solver<F,3,V>
1045  {
1046  /// method to solve the linear system
1047  void operator()(V& rhs, const F* l) const
1048  {
1049  // solve Ly = rhs
1050  const F y0 = rhs[0] * l[0];
1051  const F y1 = (rhs[1]-l[1]*y0)*l[2];
1052  const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1053  // solve L^Tx = y, and put x into rhs
1054  rhs[2] = (y2)*l[5];
1055  rhs[1] = (y1-(l[4]*rhs[2]))*l[2];
1056  rhs[0] = (y0-(l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1057  }
1058  };
1059  /// struct to solve a linear system using its Cholesky decomposition (N=2)
1060  template<class F, class V> struct _solver<F,2,V>
1061  {
1062  /// method to solve the linear system
1063  void operator()(V& rhs, const F* l) const
1064  {
1065  // solve Ly = rhs
1066  const F y0 = rhs[0] * l[0];
1067  const F y1 = (rhs[1]-l[1]*y0)*l[2];
1068  // solve L^Tx = y, and put x into rhs
1069  rhs[1] = (y1)*l[2];
1070  rhs[0] = (y0-(l[1]*rhs[1]))*l[0];
1071  }
1072  };
1073  /// struct to solve a linear system using its Cholesky decomposition (N=1)
1074  template<class F, class V> struct _solver<F,1,V>
1075  {
1076  /// method to solve the linear system
1077  void operator()(V& rhs, const F* l) const
1078  {
1079  // solve LL^T x = rhs, put y into rhs
1080  rhs[0] *= l[0] * l[0];
1081  }
1082  };
1083  /// struct to solve a linear system using its Cholesky decomposition (N=0)
1084  template<class F, class V> struct _solver<F,0,V>
1085  {
1086  private:
1087  _solver();
1088  void operator()(V& rhs, const F* l) const;
1089  };
1090 }
1091 
1092 
1093  } // namespace Math
1094 
1095 } // namespace ROOT
1096 
1097 #endif // ROOT_Math_CHOLESKYDECOMP
1098 
bool Solve(V &rhs) const
solves a linear system for the given right hand side
bool getLi(M &m) const
obtain the inverse of the decomposed matrix L
bool getL(M &m) const
obtain the decomposed matrix L
F fL[N *(N+1)/2]
lower triangular matrix L
bool Invert(M &m) const
place the inverse into m
const G operator()(unsigned i, unsigned j) const
read access to elements (make sure that j <= i)
CholeskyDecompGenDim(unsigned N, G *m)
perform a Cholesky decomposition
bool fOk
flag indicating a successful decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
unsigned fN
dimensionality dimensionality of the problem
CholeskyDecomp(G *m)
perform a Cholesky decomposition
void operator()(V &rhs, const F *l, unsigned N) const
method to solve the linear system
#define N
void operator()(M &dst, const F *src) const
method to do the inversion
adapter for packed arrays (to SMatrix indexing conventions)
ClassImp(TIterator) Bool_t TIterator return false
Compare two iterator objects.
Definition: TIterator.cxx:21
#define G(x, y, z)
bool fOk
flag indicating a successful decomposition
bool Invert(G *m) const
place the inverse into m
class to compute the Cholesky decomposition of a matrix
double sqrt(double)
bool operator()(F *dst, const M &src) const
method to do the decomposition
void operator()(M &dst, const F *src) const
method to do the inversion
bool getL(G *m) const
obtain the decomposed matrix L
void operator()(M &dst, const F *src) const
method to do the inversion
struct to solve a linear system using its Cholesky decomposition (generalised dimensionality) ...
void operator()(V &rhs, const F *l) const
method to solve the linear system
void operator()(M &dst, const F *src) const
method to do the inversion
struct to do a Cholesky decomposition (general dimensionality)
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool Invert(M &m) const
place the inverse into m
struct to obtain the inverse from a Cholesky decomposition
bool getLi(G *m) const
obtain the inverse of the decomposed matrix L
void operator()(M &dst, const F *src, unsigned N) const
method to do the inversion
#define F(x, y, z)
G & operator()(unsigned i, unsigned j)
write access to elements (make sure that j <= i)
struct to do a Cholesky decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
void operator()(V &rhs, const F *l) const
method to solve the linear system
bool Solve(V &rhs) const
solves a linear system for the given right hand side
bool Invert(G *m) const
place the inverse into m
bool ok() const
returns true if decomposition was successful
TMarker * m
Definition: textangle.C:8
void operator()(M &dst, const F *src) const
method to do the inversion
TLine * l
Definition: textangle.C:4
bool operator()(F *dst, const M &src) const
method to do the decomposition
class to compute the Cholesky decomposition of a matrix
bool getLi(M &m) const
obtain the inverse of the decomposed matrix L
CholeskyDecompGenDim(unsigned N, const M &m)
perform a Cholesky decomposition
void operator()(V &rhs, const F *l) const
method to solve the linear system
CholeskyDecomp(const M &m)
perform a Cholesky decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
void operator()(V &rhs, const F *l) const
method to solve the linear system
void operator()(V &rhs, const F *l) const
method to solve the linear system
F * fL
lower triangular matrix L
void operator()(V &rhs, const F *l) const
method to solve the linear system
bool getLi(G *m) const
obtain the inverse of the decomposed matrix L
bool getL(M &m) const
obtain the decomposed matrix L
bool getL(G *m) const
obtain the decomposed matrix L
void operator()(V &rhs, const F *l) const
method to solve the linear system
void operator()(M &dst, const F *src) const
method to do the inversion
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool ok() const
returns true if decomposition was successful
bool operator()(F *dst, const M &src, unsigned N) const
method to do the decomposition
void operator()(M &dst, const F *src) const
method to do the inversion
struct to solve a linear system using its Cholesky decomposition
struct to obtain the inverse from a Cholesky decomposition (general dimensionality) ...