4 #ifndef ROOT_Math_CholeskyDecomp
5 #define ROOT_Math_CholeskyDecomp
35 namespace CholeskyDecompHelpers {
40 template<
class F,
unsigned N,
class M>
struct _inverter;
42 template<
class F,
unsigned N,
class V>
struct _solver;
98 fOk = _decomposer<F, N, M>()(
fL, m);
117 fOk = _decomposer<F, N, PackedArrayAdapter<G> >()(
118 fL, PackedArrayAdapter<G>(m));
126 operator bool()
const {
return fOk; }
136 template<
class V>
bool Solve(V& rhs)
const
139 if (
fOk) _solver<F,N,V>()(rhs,
fL);
return fOk;
151 if (
fOk) _inverter<F,N,M>()(m,
fL);
return fOk;
169 PackedArrayAdapter<G> adapted(m);
170 _inverter<F,N,PackedArrayAdapter<G> >()(adapted,
fL);
181 template<
class M>
bool getL(M&
m)
const
183 if (!
fOk)
return false;
184 for (
unsigned i = 0; i <
N; ++i) {
186 for (
unsigned j = i + 1; j <
N; ++j)
189 for (
unsigned j = 0; j <= i; ++j)
190 m(i, j) =
fL[i * (i + 1) / 2 + j];
193 m(i, i) =
F(1) /
m(i, i);
206 template<
typename G>
bool getL(
G*
m)
const
208 if (!
fOk)
return false;
210 for (
unsigned i = 0; i < (
N * (
N + 1)) / 2; ++i)
214 for (
unsigned i = 0; i <
N; ++i)
215 m[(i * (i + 1)) / 2 + i] =
F(1) /
fL[(i * (i + 1)) / 2 + i];
227 if (!
fOk)
return false;
228 for (
unsigned i = 0; i <
N; ++i) {
230 for (
unsigned j = i + 1; j <
N; ++j)
233 for (
unsigned j = 0; j <= i; ++j)
234 m(j, i) =
fL[i * (i + 1) / 2 + j];
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);
258 if (!
fOk)
return false;
260 for (
unsigned i = 0; i < (
N * (
N + 1)) / 2; ++i)
264 for (
unsigned i = 1; i <
N; base1 += ++i) {
265 for (
unsigned j = 0; j < i; ++j) {
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];
333 fOk = _decomposerGenDim<F, M>()(
fL, m,
fN);
352 fOk = _decomposerGenDim<F, PackedArrayAdapter<G> >()(
353 fL, PackedArrayAdapter<G>(m),
fN);
364 operator bool()
const {
return fOk; }
374 template<
class V>
bool Solve(V& rhs)
const
377 if (
fOk) _solverGenDim<F,V>()(rhs,
fL,
fN);
return fOk;
389 if (
fOk) _inverterGenDim<F,M>()(m,
fL,
fN);
return fOk;
407 PackedArrayAdapter<G> adapted(m);
408 _inverterGenDim<F,PackedArrayAdapter<G> >()(adapted,
fL,
fN);
419 template<
class M>
bool getL(M&
m)
const
421 if (!
fOk)
return false;
422 for (
unsigned i = 0; i <
fN; ++i) {
424 for (
unsigned j = i + 1; j <
fN; ++j)
427 for (
unsigned j = 0; j <= i; ++j)
428 m(i, j) =
fL[i * (i + 1) / 2 + j];
431 m(i, i) =
F(1) /
m(i, i);
444 template<
typename G>
bool getL(
G*
m)
const
446 if (!
fOk)
return false;
448 for (
unsigned i = 0; i < (
fN * (
fN + 1)) / 2; ++i)
452 for (
unsigned i = 0; i <
fN; ++i)
453 m[(i * (i + 1)) / 2 + i] =
F(1) /
fL[(i * (i + 1)) / 2 + i];
465 if (!
fOk)
return false;
466 for (
unsigned i = 0; i <
fN; ++i) {
468 for (
unsigned j = i + 1; j <
fN; ++j)
471 for (
unsigned j = 0; j <= i; ++j)
472 m(j, i) =
fL[i * (i + 1) / 2 + j];
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);
496 if (!
fOk)
return false;
498 for (
unsigned i = 0; i < (
fN * (
fN + 1)) / 2; ++i)
502 for (
unsigned i = 1; i <
fN; base1 += ++i) {
503 for (
unsigned j = 0; j < i; ++j) {
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];
515 namespace CholeskyDecompHelpers {
517 template<
typename G>
class PackedArrayAdapter
526 {
return fArr[((i * (i + 1)) / 2) + j]; }
529 {
return fArr[((i * (i + 1)) / 2) + j]; }
532 template<
class F,
class M>
struct _decomposerGenDim
551 for (
unsigned i = 0; i <
N; base1 += ++i) {
555 for (
unsigned j = 0; j < i; base2 += ++j) {
557 for (
unsigned k = j; k--; )
558 tmp -= base1[k] * base2[k];
559 base1[j] = tmp *= base2[j];
561 tmpdiag += tmp * tmp;
564 tmpdiag = src(i, i) - tmpdiag;
566 if (tmpdiag <=
F(0.0))
return false;
574 template<
class F,
unsigned N,
class M>
struct _decomposer
583 template<
class F,
class M>
struct _inverterGenDim
589 F *
l =
new F[N * (N + 1) / 2];
590 std::copy(src, src + ((N * (N + 1)) / 2), l);
593 for (
unsigned i = 1; i <
N; base1 += ++i) {
594 for (
unsigned j = 0; j < i; ++j) {
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];
604 for (
unsigned i = N; i--; ) {
605 for (
unsigned j = i + 1; j--; ) {
607 base1 = &l[(N * (N - 1)) / 2];
608 for (
unsigned k = N; k-- > i; base1 -= k)
609 tmp += base1[i] * base1[j];
618 template<
class F,
unsigned N,
class M>
struct _inverter
626 template<
class F,
class V>
struct _solverGenDim
632 for (
unsigned k = 0; k <
N; ++k) {
633 const unsigned base = (k * (k + 1)) / 2;
635 for (
unsigned i = k; i--; )
636 sum += rhs[i] * l[base + i];
638 rhs[k] = (rhs[k] - sum) * l[base + k];
641 for (
unsigned k = N; k--; ) {
643 for (
unsigned i = N; --i > k; )
644 sum += rhs[i] * l[(i * (i + 1)) / 2 + k];
646 rhs[k] = (rhs[k] - sum) * l[(k * (k + 1)) / 2 + k];
652 template<
class F,
unsigned N,
class V>
struct _solver
665 if (src(0,0) <=
F(0.0))
return false;
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;
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;
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;
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;
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;
706 if (src(0,0) <=
F(0.0))
return false;
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;
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;
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;
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;
739 if (src(0,0) <=
F(0.0))
return false;
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;
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;
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;
765 if (src(0,0) <=
F(0.0))
return false;
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;
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;
785 if (src(0,0) <=
F(0.0))
return false;
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;
800 if (src(0,0) <=
F(0.0))
return false;
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]) *
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];
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];
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];
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];
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];
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];
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];
959 const F li21 = -src[1] * src[0] * src[2];
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];
972 dst(0,0) = src[0]*src[0];
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];
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];
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];
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];
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];
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];
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];
1055 rhs[1] = (y1-(l[4]*rhs[2]))*l[2];
1056 rhs[0] = (y0-(l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1066 const F y0 = rhs[0] * l[0];
1067 const F y1 = (rhs[1]-l[1]*y0)*l[2];
1070 rhs[0] = (y0-(l[1]*rhs[1]))*l[0];
1080 rhs[0] *= l[0] * l[0];
1097 #endif // ROOT_Math_CHOLESKYDECOMP
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()
destructor
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
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.
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
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)
G * fArr
pointer to first array element
bool operator()(F *dst, const M &src) const
method to do the decomposition
PackedArrayAdapter(G *arr)
constructor
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
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
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
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) ...