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;
82 F fL[
N * (
N + 1) / 2];
98 fOk = _decomposer<F, N, M>()(fL, m);
117 fOk = _decomposer<F, N, PackedArrayAdapter<G> >()(
118 fL, PackedArrayAdapter<G>(m));
123 bool ok()
const {
return fOk; }
126 operator bool()
const {
return fOk; }
136 template<
class V>
bool Solve(V& rhs)
const 139 if (fOk) _solver<F,N,V>()(rhs, fL);
152 if (fOk) _inverter<F,N,M>()(m, fL);
171 PackedArrayAdapter<G> adapted(m);
172 _inverter<F,N,PackedArrayAdapter<G> >()(adapted, fL);
183 template<
class M>
bool getL(M&
m)
const 185 if (!fOk)
return false;
186 for (
unsigned i = 0; i <
N; ++i) {
188 for (
unsigned j = i + 1; j <
N; ++j)
191 for (
unsigned j = 0; j <= i; ++j)
192 m(i, j) = fL[i * (i + 1) / 2 + j];
195 m(i, i) =
F(1) /
m(i, i);
208 template<
typename G>
bool getL(
G*
m)
const 210 if (!fOk)
return false;
212 for (
unsigned i = 0; i < (
N * (
N + 1)) / 2; ++i)
216 for (
unsigned i = 0; i <
N; ++i)
217 m[(i * (i + 1)) / 2 + i] =
F(1) / fL[(i * (i + 1)) / 2 + i];
229 if (!fOk)
return false;
230 for (
unsigned i = 0; i <
N; ++i) {
232 for (
unsigned j = i + 1; j <
N; ++j)
235 for (
unsigned j = 0; j <= i; ++j)
236 m(j, i) = fL[i * (i + 1) / 2 + j];
239 for (
unsigned i = 1; i <
N; ++i) {
240 for (
unsigned j = 0; j < i; ++j) {
241 typename M::value_type tmp =
F(0);
242 for (
unsigned k = i; k-- > j;)
243 tmp -=
m(k, i) *
m(j, k);
244 m(j, i) = tmp *
m(i, i);
260 if (!fOk)
return false;
262 for (
unsigned i = 0; i < (
N * (
N + 1)) / 2; ++i)
266 for (
unsigned i = 1; i <
N; base1 += ++i) {
267 for (
unsigned j = 0; j < i; ++j) {
269 const G *base2 = &m[(i * (i - 1)) / 2];
270 for (
unsigned k = i; k-- > j; base2 -= k)
271 tmp -= base1[k] * base2[j];
272 base1[j] = tmp * base1[i];
332 fN(N), fL(new
F[(fN * (fN + 1)) / 2]), fOk(false)
335 fOk = _decomposerGenDim<F, M>()(fL, m, fN);
350 fN(N), fL(new
F[(fN * (fN + 1)) / 2]), fOk(false)
354 fOk = _decomposerGenDim<F, PackedArrayAdapter<G> >()(
355 fL, PackedArrayAdapter<G>(m), fN);
363 bool ok()
const {
return fOk; }
366 operator bool()
const {
return fOk; }
376 template<
class V>
bool Solve(V& rhs)
const 379 if (fOk) _solverGenDim<F,V>()(rhs, fL, fN);
392 if (fOk) _inverterGenDim<F,M>()(m, fL, fN);
411 PackedArrayAdapter<G> adapted(m);
412 _inverterGenDim<F,PackedArrayAdapter<G> >()(adapted, fL, fN);
423 template<
class M>
bool getL(M&
m)
const 425 if (!fOk)
return false;
426 for (
unsigned i = 0; i < fN; ++i) {
428 for (
unsigned j = i + 1; j < fN; ++j)
431 for (
unsigned j = 0; j <= i; ++j)
432 m(i, j) = fL[i * (i + 1) / 2 + j];
435 m(i, i) =
F(1) /
m(i, i);
448 template<
typename G>
bool getL(
G*
m)
const 450 if (!fOk)
return false;
452 for (
unsigned i = 0; i < (fN * (fN + 1)) / 2; ++i)
456 for (
unsigned i = 0; i < fN; ++i)
457 m[(i * (i + 1)) / 2 + i] =
F(1) / fL[(i * (i + 1)) / 2 + i];
469 if (!fOk)
return false;
470 for (
unsigned i = 0; i < fN; ++i) {
472 for (
unsigned j = i + 1; j < fN; ++j)
475 for (
unsigned j = 0; j <= i; ++j)
476 m(j, i) = fL[i * (i + 1) / 2 + j];
479 for (
unsigned i = 1; i < fN; ++i) {
480 for (
unsigned j = 0; j < i; ++j) {
481 typename M::value_type tmp =
F(0);
482 for (
unsigned k = i; k-- > j;)
483 tmp -=
m(k, i) *
m(j, k);
484 m(j, i) = tmp *
m(i, i);
500 if (!fOk)
return false;
502 for (
unsigned i = 0; i < (fN * (fN + 1)) / 2; ++i)
506 for (
unsigned i = 1; i < fN; base1 += ++i) {
507 for (
unsigned j = 0; j < i; ++j) {
509 const G *base2 = &m[(i * (i - 1)) / 2];
510 for (
unsigned k = i; k-- > j; base2 -= k)
511 tmp -= base1[k] * base2[j];
512 base1[j] = tmp * base1[i];
519 namespace CholeskyDecompHelpers {
521 template<
typename G>
class PackedArrayAdapter
530 {
return fArr[((i * (i + 1)) / 2) + j]; }
533 {
return fArr[((i * (i + 1)) / 2) + j]; }
555 for (
unsigned i = 0; i <
N; base1 += ++i) {
559 for (
unsigned j = 0; j < i; base2 += ++j) {
561 for (
unsigned k = j; k--; )
562 tmp -= base1[k] * base2[k];
563 base1[j] = tmp *= base2[j];
565 tmpdiag += tmp * tmp;
568 tmpdiag = src(i, i) - tmpdiag;
570 if (tmpdiag <=
F(0.0))
return false;
578 template<
class F,
unsigned N,
class M>
struct _decomposer 593 F *
l =
new F[N * (N + 1) / 2];
594 std::copy(src, src + ((N * (N + 1)) / 2), l);
597 for (
unsigned i = 1; i <
N; base1 += ++i) {
598 for (
unsigned j = 0; j < i; ++j) {
600 const F *base2 = &l[(i * (i - 1)) / 2];
601 for (
unsigned k = i; k-- > j; base2 -= k)
602 tmp -= base1[k] * base2[j];
603 base1[j] = tmp * base1[i];
608 for (
unsigned i = N; i--; ) {
609 for (
unsigned j = i + 1; j--; ) {
611 base1 = &l[(N * (N - 1)) / 2];
612 for (
unsigned k = N; k-- > i; base1 -= k)
613 tmp += base1[i] * base1[j];
622 template<
class F,
unsigned N,
class M>
struct _inverter 636 for (
unsigned k = 0; k <
N; ++k) {
637 const unsigned base = (k * (k + 1)) / 2;
639 for (
unsigned i = k; i--; )
640 sum += rhs[i] * l[base + i];
642 rhs[k] = (rhs[k] -
sum) * l[base + k];
645 for (
unsigned k = N; k--; ) {
647 for (
unsigned i = N; --i > k; )
648 sum += rhs[i] * l[(i * (i + 1)) / 2 + k];
650 rhs[k] = (rhs[k] -
sum) * l[(k * (k + 1)) / 2 + k];
656 template<
class F,
unsigned N,
class V>
struct _solver 669 if (src(0,0) <=
F(0.0))
return false;
671 dst[1] = src(1,0) * dst[0];
672 dst[2] = src(1,1) - dst[1] * dst[1];
673 if (dst[2] <=
F(0.0))
return false;
675 dst[3] = src(2,0) * dst[0];
676 dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
677 dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
678 if (dst[5] <=
F(0.0))
return false;
680 dst[6] = src(3,0) * dst[0];
681 dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
682 dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
683 dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
684 if (dst[9] <=
F(0.0))
return false;
686 dst[10] = src(4,0) * dst[0];
687 dst[11] = (src(4,1) - dst[1] * dst[10]) * dst[2];
688 dst[12] = (src(4,2) - dst[3] * dst[10] - dst[4] * dst[11]) * dst[5];
689 dst[13] = (src(4,3) - dst[6] * dst[10] - dst[7] * dst[11] - dst[8] * dst[12]) * dst[9];
690 dst[14] = src(4,4) - (dst[10]*dst[10]+dst[11]*dst[11]+dst[12]*dst[12]+dst[13]*dst[13]);
691 if (dst[14] <=
F(0.0))
return false;
693 dst[15] = src(5,0) * dst[0];
694 dst[16] = (src(5,1) - dst[1] * dst[15]) * dst[2];
695 dst[17] = (src(5,2) - dst[3] * dst[15] - dst[4] * dst[16]) * dst[5];
696 dst[18] = (src(5,3) - dst[6] * dst[15] - dst[7] * dst[16] - dst[8] * dst[17]) * dst[9];
697 dst[19] = (src(5,4) - dst[10] * dst[15] - dst[11] * dst[16] - dst[12] * dst[17] - dst[13] * dst[18]) * dst[14];
698 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]);
699 if (dst[20] <=
F(0.0))
return false;
710 if (src(0,0) <=
F(0.0))
return false;
712 dst[1] = src(1,0) * dst[0];
713 dst[2] = src(1,1) - dst[1] * dst[1];
714 if (dst[2] <=
F(0.0))
return false;
716 dst[3] = src(2,0) * dst[0];
717 dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
718 dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
719 if (dst[5] <=
F(0.0))
return false;
721 dst[6] = src(3,0) * dst[0];
722 dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
723 dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
724 dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
725 if (dst[9] <=
F(0.0))
return false;
727 dst[10] = src(4,0) * dst[0];
728 dst[11] = (src(4,1) - dst[1] * dst[10]) * dst[2];
729 dst[12] = (src(4,2) - dst[3] * dst[10] - dst[4] * dst[11]) * dst[5];
730 dst[13] = (src(4,3) - dst[6] * dst[10] - dst[7] * dst[11] - dst[8] * dst[12]) * dst[9];
731 dst[14] = src(4,4) - (dst[10]*dst[10]+dst[11]*dst[11]+dst[12]*dst[12]+dst[13]*dst[13]);
732 if (dst[14] <=
F(0.0))
return false;
743 if (src(0,0) <=
F(0.0))
return false;
745 dst[1] = src(1,0) * dst[0];
746 dst[2] = src(1,1) - dst[1] * dst[1];
747 if (dst[2] <=
F(0.0))
return false;
749 dst[3] = src(2,0) * dst[0];
750 dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
751 dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
752 if (dst[5] <=
F(0.0))
return false;
754 dst[6] = src(3,0) * dst[0];
755 dst[7] = (src(3,1) - dst[1] * dst[6]) * dst[2];
756 dst[8] = (src(3,2) - dst[3] * dst[6] - dst[4] * dst[7]) * dst[5];
757 dst[9] = src(3,3) - (dst[6] * dst[6] + dst[7] * dst[7] + dst[8] * dst[8]);
758 if (dst[9] <=
F(0.0))
return false;
769 if (src(0,0) <=
F(0.0))
return false;
771 dst[1] = src(1,0) * dst[0];
772 dst[2] = src(1,1) - dst[1] * dst[1];
773 if (dst[2] <=
F(0.0))
return false;
775 dst[3] = src(2,0) * dst[0];
776 dst[4] = (src(2,1) - dst[1] * dst[3]) * dst[2];
777 dst[5] = src(2,2) - (dst[3] * dst[3] + dst[4] * dst[4]);
778 if (dst[5] <=
F(0.0))
return false;
789 if (src(0,0) <=
F(0.0))
return false;
791 dst[1] = src(1,0) * dst[0];
792 dst[2] = src(1,1) - dst[1] * dst[1];
793 if (dst[2] <=
F(0.0))
return false;
804 if (src(0,0) <=
F(0.0))
return false;
814 bool operator()(
F* dst,
const M& src)
const;
823 const F li21 = -src[1] * src[0] * src[2];
824 const F li32 = -src[4] * src[2] * src[5];
825 const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
826 const F li43 = -src[8] * src[9] * src[5];
827 const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
828 const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
829 src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
830 const F li54 = -src[13] * src[14] * src[9];
831 const F li53 = (src[13] * src[8] * src[9] - src[12]) * src[5] * src[14];
832 const F li52 = (-src[4] * src[8] * src[13] * src[5] * src[9] +
833 src[4] * src[12] * src[5] + src[7] * src[13] * src[9] - src[11]) * src[2] * src[14];
834 const F li51 = (src[1]*src[4]*src[8]*src[13]*src[2]*src[5]*src[9] -
835 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] +
836 src[11]*src[1]*src[2] + src[12]*src[3]*src[5] + src[13]*src[6]*src[9] -src[10]) * src[0] * src[14];
837 const F li65 = -src[19] * src[20] * src[14];
838 const F li64 = (src[19] * src[13] * src[14] - src[18]) * src[9] * src[20];
839 const F li63 = (-src[8] * src[13] * src[19] * src[9] * src[14] +
840 src[8] * src[18] * src[9] + src[12] * src[19] * src[14] - src[17]) * src[5] * src[20];
841 const F li62 = (src[4]*src[8]*src[13]*src[19]*src[5]*src[9]*src[14] -
842 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] +
843 src[17]*src[4]*src[5] + src[18]*src[7]*src[9] + src[19]*src[11]*src[14] - src[16]) * src[2] * src[20];
844 const F li61 = (-src[19]*src[13]*src[8]*src[4]*src[1]*src[2]*src[5]*src[9]*src[14] +
845 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] +
846 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] -
847 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] -
848 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] +
849 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]) *
852 dst(0,0) = li61*li61 + li51*li51 + li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
853 dst(1,0) = li61*li62 + li51*li52 + li41*li42 + li31*li32 + li21*src[2];
854 dst(1,1) = li62*li62 + li52*li52 + li42*li42 + li32*li32 + src[2]*src[2];
855 dst(2,0) = li61*li63 + li51*li53 + li41*li43 + li31*src[5];
856 dst(2,1) = li62*li63 + li52*li53 + li42*li43 + li32*src[5];
857 dst(2,2) = li63*li63 + li53*li53 + li43*li43 + src[5]*src[5];
858 dst(3,0) = li61*li64 + li51*li54 + li41*src[9];
859 dst(3,1) = li62*li64 + li52*li54 + li42*src[9];
860 dst(3,2) = li63*li64 + li53*li54 + li43*src[9];
861 dst(3,3) = li64*li64 + li54*li54 + src[9]*src[9];
862 dst(4,0) = li61*li65 + li51*src[14];
863 dst(4,1) = li62*li65 + li52*src[14];
864 dst(4,2) = li63*li65 + li53*src[14];
865 dst(4,3) = li64*li65 + li54*src[14];
866 dst(4,4) = li65*li65 + src[14]*src[14];
867 dst(5,0) = li61*src[20];
868 dst(5,1) = li62*src[20];
869 dst(5,2) = li63*src[20];
870 dst(5,3) = li64*src[20];
871 dst(5,4) = li65*src[20];
872 dst(5,5) = src[20]*src[20];
881 const F li21 = -src[1] * src[0] * src[2];
882 const F li32 = -src[4] * src[2] * src[5];
883 const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
884 const F li43 = -src[8] * src[9] * src[5];
885 const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
886 const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
887 src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
888 const F li54 = -src[13] * src[14] * src[9];
889 const F li53 = (src[13] * src[8] * src[9] - src[12]) * src[5] * src[14];
890 const F li52 = (-src[4] * src[8] * src[13] * src[5] * src[9] +
891 src[4] * src[12] * src[5] + src[7] * src[13] * src[9] - src[11]) * src[2] * src[14];
892 const F li51 = (src[1]*src[4]*src[8]*src[13]*src[2]*src[5]*src[9] -
893 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] +
894 src[11]*src[1]*src[2] + src[12]*src[3]*src[5] + src[13]*src[6]*src[9] -src[10]) * src[0] * src[14];
896 dst(0,0) = li51*li51 + li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
897 dst(1,0) = li51*li52 + li41*li42 + li31*li32 + li21*src[2];
898 dst(1,1) = li52*li52 + li42*li42 + li32*li32 + src[2]*src[2];
899 dst(2,0) = li51*li53 + li41*li43 + li31*src[5];
900 dst(2,1) = li52*li53 + li42*li43 + li32*src[5];
901 dst(2,2) = li53*li53 + li43*li43 + src[5]*src[5];
902 dst(3,0) = li51*li54 + li41*src[9];
903 dst(3,1) = li52*li54 + li42*src[9];
904 dst(3,2) = li53*li54 + li43*src[9];
905 dst(3,3) = li54*li54 + src[9]*src[9];
906 dst(4,0) = li51*src[14];
907 dst(4,1) = li52*src[14];
908 dst(4,2) = li53*src[14];
909 dst(4,3) = li54*src[14];
910 dst(4,4) = src[14]*src[14];
919 const F li21 = -src[1] * src[0] * src[2];
920 const F li32 = -src[4] * src[2] * src[5];
921 const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
922 const F li43 = -src[8] * src[9] * src[5];
923 const F li42 = (src[4] * src[8] * src[5] - src[7]) * src[2] * src[9];
924 const F li41 = (-src[1] * src[4] * src[8] * src[2] * src[5] +
925 src[1] * src[7] * src[2] + src[3] * src[8] * src[5] - src[6]) * src[0] * src[9];
927 dst(0,0) = li41*li41 + li31*li31 + li21*li21 + src[0]*src[0];
928 dst(1,0) = li41*li42 + li31*li32 + li21*src[2];
929 dst(1,1) = li42*li42 + li32*li32 + src[2]*src[2];
930 dst(2,0) = li41*li43 + li31*src[5];
931 dst(2,1) = li42*li43 + li32*src[5];
932 dst(2,2) = li43*li43 + src[5]*src[5];
933 dst(3,0) = li41*src[9];
934 dst(3,1) = li42*src[9];
935 dst(3,2) = li43*src[9];
936 dst(3,3) = src[9]*src[9];
945 const F li21 = -src[1] * src[0] * src[2];
946 const F li32 = -src[4] * src[2] * src[5];
947 const F li31 = (src[1] * src[4] * src[2] - src[3]) * src[0] * src[5];
949 dst(0,0) = li31*li31 + li21*li21 + src[0]*src[0];
950 dst(1,0) = li31*li32 + li21*src[2];
951 dst(1,1) = li32*li32 + src[2]*src[2];
952 dst(2,0) = li31*src[5];
953 dst(2,1) = li32*src[5];
954 dst(2,2) = src[5]*src[5];
963 const F li21 = -src[1] * src[0] * src[2];
965 dst(0,0) = li21*li21 + src[0]*src[0];
966 dst(1,0) = li21*src[2];
967 dst(1,1) = src[2]*src[2];
976 dst(0,0) = src[0]*src[0];
984 void operator()(M& dst,
const F* src)
const;
994 const F y0 = rhs[0] * l[0];
995 const F y1 = (rhs[1]-l[1]*y0)*l[2];
996 const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
997 const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
998 const F y4 = (rhs[4]-(l[10]*y0+l[11]*y1+l[12]*y2+l[13]*y3))*l[14];
999 const F y5 = (rhs[5]-(l[15]*y0+l[16]*y1+l[17]*y2+l[18]*y3+l[19]*y4))*l[20];
1001 rhs[5] = y5 * l[20];
1002 rhs[4] = (y4-l[19]*rhs[5])*l[14];
1003 rhs[3] = (y3-(l[18]*rhs[5]+l[13]*rhs[4]))*l[9];
1004 rhs[2] = (y2-(l[17]*rhs[5]+l[12]*rhs[4]+l[8]*rhs[3]))*l[5];
1005 rhs[1] = (y1-(l[16]*rhs[5]+l[11]*rhs[4]+l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1006 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];
1016 const F y0 = rhs[0] * l[0];
1017 const F y1 = (rhs[1]-l[1]*y0)*l[2];
1018 const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1019 const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
1020 const F y4 = (rhs[4]-(l[10]*y0+l[11]*y1+l[12]*y2+l[13]*y3))*l[14];
1022 rhs[4] = (y4)*l[14];
1023 rhs[3] = (y3-(l[13]*rhs[4]))*l[9];
1024 rhs[2] = (y2-(l[12]*rhs[4]+l[8]*rhs[3]))*l[5];
1025 rhs[1] = (y1-(l[11]*rhs[4]+l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1026 rhs[0] = (y0-(l[10]*rhs[4]+l[6]*rhs[3]+l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1036 const F y0 = rhs[0] * l[0];
1037 const F y1 = (rhs[1]-l[1]*y0)*l[2];
1038 const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1039 const F y3 = (rhs[3]-(l[6]*y0+l[7]*y1+l[8]*y2))*l[9];
1042 rhs[2] = (y2-(l[8]*rhs[3]))*l[5];
1043 rhs[1] = (y1-(l[7]*rhs[3]+l[4]*rhs[2]))*l[2];
1044 rhs[0] = (y0-(l[6]*rhs[3]+l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1054 const F y0 = rhs[0] * l[0];
1055 const F y1 = (rhs[1]-l[1]*y0)*l[2];
1056 const F y2 = (rhs[2]-(l[3]*y0+l[4]*y1))*l[5];
1059 rhs[1] = (y1-(l[4]*rhs[2]))*l[2];
1060 rhs[0] = (y0-(l[3]*rhs[2]+l[1]*rhs[1]))*l[0];
1070 const F y0 = rhs[0] * l[0];
1071 const F y1 = (rhs[1]-l[1]*y0)*l[2];
1074 rhs[0] = (y0-(l[1]*rhs[1]))*l[0];
1084 rhs[0] *= l[0] * l[0];
1092 void operator()(V& rhs,
const F*
l)
const;
1101 #endif // ROOT_Math_CHOLESKYDECOMP bool ok() const
returns true if decomposition was successful
bool getL(G *m) const
obtain the decomposed matrix L
static long int sum(long int i)
const G operator()(unsigned i, unsigned j) const
read access to elements (make sure that j <= i)
void operator()(V &rhs, const F *l) const
method to solve the linear system
This namespace contains pre-defined functions to be used in conjuction with TExecutor::Map and TExecu...
~CholeskyDecompGenDim()
destructor
bool Invert(G *m) const
place the inverse into m
CholeskyDecompGenDim(unsigned N, G *m)
perform a Cholesky decomposition
bool Invert(M &m) const
place the inverse into m
bool fOk
flag indicating a successful decomposition
bool getL(M &m) const
obtain the decomposed matrix L
unsigned fN
dimensionality dimensionality of the problem
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
CholeskyDecomp(G *m)
perform a Cholesky decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
adapter for packed arrays (to SMatrix indexing conventions)
void operator()(M &dst, const F *src) const
method to do the inversion
void operator()(M &dst, const F *src) const
method to do the inversion
bool getLi(M &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
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
void operator()(V &rhs, const F *l) const
method to solve the linear system
bool fOk
flag indicating a successful decomposition
void operator()(M &dst, const F *src) const
method to do the inversion
class to compute the Cholesky decomposition of a matrix
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool getL(G *m) const
obtain the decomposed matrix L
bool getLi(G *m) const
obtain the inverse of 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 getLi(G *m) const
obtain the inverse of the decomposed matrix L
void operator()(V &rhs, const F *l) const
method to solve the linear system
struct to solve a linear system using its Cholesky decomposition (generalised dimensionality) ...
bool operator()(F *dst, const M &src) const
method to do the decomposition
struct to do a Cholesky decomposition (general dimensionality)
G * fArr
pointer to first array element
PackedArrayAdapter(G *arr)
constructor
struct to obtain the inverse from a Cholesky decomposition
bool getLi(M &m) const
obtain the inverse of the decomposed matrix L
G & operator()(unsigned i, unsigned j)
write access to elements (make sure that j <= i)
bool Invert(M &m) const
place the inverse into m
struct to do a Cholesky decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool Solve(V &rhs) const
solves a linear system for the given right hand side
void operator()(M &dst, const F *src) const
method to do the inversion
class to compute the Cholesky decomposition of a matrix
CholeskyDecompGenDim(unsigned N, const M &m)
perform a Cholesky decomposition
bool getL(M &m) const
obtain the decomposed matrix L
CholeskyDecomp(const M &m)
perform a Cholesky decomposition
bool operator()(F *dst, const M &src) const
method to do the decomposition
F * fL
lower triangular matrix L
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
bool Invert(G *m) const
place the inverse into m
Namespace for new Math classes and functions.
bool Solve(V &rhs) const
solves a linear system for the given right hand side
bool operator()(F *dst, const M &src) const
method to do the decomposition
bool ok() const
returns true if decomposition was successful
struct to solve a linear system using its Cholesky decomposition
void operator()(V &rhs, const F *l, unsigned N) const
method to solve the linear system
struct to obtain the inverse from a Cholesky decomposition (general dimensionality) ...