4#ifndef ROOT_Math_MatrixFunctions
5#define ROOT_Math_MatrixFunctions
54 template <
class T,
unsigned int D>
61template <
class T,
unsigned int D1,
unsigned int D2,
class R>
62SVector<T,D1>
operator*(
const SMatrix<T,D1,D2,R>& rhs,
const SVector<T,D2>& lhs)
65 for(
unsigned int i=0; i<D1; ++i) {
66 const unsigned int rpos = i*D2;
67 for(
unsigned int j=0; j<D2; ++j) {
68 tmp[i] += rhs.apply(rpos+j) * lhs.apply(j);
83template <
unsigned int I>
85 template <
class A,
class B>
86 static inline typename A::value_type
f(
const A& lhs,
const B& rhs,
87 const unsigned int offset) {
98 template <
class A,
class B>
99 static inline typename A::value_type
f(
const A& lhs,
const B& rhs,
100 const unsigned int offset) {
101 return lhs.apply(
offset) * rhs.apply(0);
108template <
class Matrix,
class Vector,
unsigned int D2>
112 typedef typename Vector::value_type
T;
122 inline typename Matrix::value_type
apply(
unsigned int i)
const {
129 return rhs_.IsInUse(
p);
142template <
unsigned int I>
144 template <
class Matrix,
class Vector>
145 static inline typename Matrix::value_type
f(
const Matrix& lhs,
const Vector& rhs,
146 const unsigned int offset) {
147 return lhs.apply(Matrix::kCols*
I+
offset) * rhs.apply(
I) +
158 template <
class Matrix,
class Vector>
159 static inline typename Matrix::value_type
f(
const Matrix& lhs,
const Vector& rhs,
160 const unsigned int offset) {
161 return lhs.apply(
offset) * rhs.apply(0);
173template <
class Vector,
class Matrix,
unsigned int D1>
177 typedef typename Vector::value_type
T;
186 inline typename Matrix::value_type
apply(
unsigned int i)
const {
193 return lhs_.IsInUse(
p);
211template <
class T,
unsigned int D1,
unsigned int D2,
class R>
221template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
222inline VecExpr<VectorMatrixRowOp<SMatrix<T,D1,D2,R>, VecExpr<A,T,D2>, D2>, T, D1>
231template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
232inline VecExpr<VectorMatrixRowOp<Expr<A,T,D1,D2,R>, SVector<T,D2>, D2>, T, D1>
241template <
class A,
class B,
class T,
unsigned int D1,
unsigned int D2,
class R>
242inline VecExpr<VectorMatrixRowOp<Expr<A,T,D1,D2,R>, VecExpr<B,T,D2>, D2>, T, D1>
251template <
class T,
unsigned int D1,
unsigned int D2,
class R>
252inline VecExpr<VectorMatrixColOp<SVector<T,D1>, SMatrix<T,D1,D2,R>, D1>, T, D2>
261template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
262inline VecExpr<VectorMatrixColOp<SVector<T,D1>, Expr<A,T,D1,D2,R>, D1>, T, D2>
271template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
272inline VecExpr<VectorMatrixColOp<VecExpr<A,T,D1>, SMatrix<T,D1,D2,R>, D1>, T, D2>
281template <
class A,
class B,
class T,
unsigned int D1,
unsigned int D2,
class R>
282inline VecExpr<VectorMatrixColOp<VecExpr<A,T,D1>, Expr<B,T,D1,D2,R>, D1>, T, D2>
291template <
unsigned int I>
294 template <
class MatrixA,
class MatrixB>
295 static inline typename MatrixA::value_type
f(
const MatrixA& lhs,
297 const unsigned int offset) {
298 return lhs.apply(
offset/MatrixB::kCols*MatrixA::kCols +
I) *
299 rhs.apply(MatrixB::kCols*
I +
offset%MatrixB::kCols) +
304 template <
class MatrixA,
class MatrixB>
305 static inline typename MatrixA::value_type
g(
const MatrixA& lhs,
309 return lhs(i,
I) * rhs(
I , j) +
321 template <
class MatrixA,
class MatrixB>
322 static inline typename MatrixA::value_type
f(
const MatrixA& lhs,
324 const unsigned int offset) {
325 return lhs.apply(
offset/MatrixB::kCols*MatrixA::kCols) *
326 rhs.apply(
offset%MatrixB::kCols);
330 template <
class MatrixA,
class MatrixB>
331 static inline typename MatrixA::value_type
g(
const MatrixA& lhs,
333 unsigned int i,
unsigned int j) {
334 return lhs(i,0) * rhs(0,j);
347template <
class MatrixA,
class MatrixB,
class T,
unsigned int D>
358 inline T
apply(
unsigned int i)
const {
386template <
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
387inline Expr<MatrixMulOp<SMatrix<T,D1,D,R1>,
SMatrix<T,D,D2,R2>,T,D>, T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
390 return Expr<MatMulOp,T,D1,D2,
397template <
class A,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
398inline Expr<MatrixMulOp<SMatrix<T,D1,D,R1>, Expr<A,T,D,D2,R2>,T,D>, T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
400 typedef MatrixMulOp<SMatrix<T,D1,D,R1>,
Expr<A,T,D,D2,R2>,T,D> MatMulOp;
401 return Expr<MatMulOp,T,D1,D2,
408template <
class A,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
409inline Expr<MatrixMulOp<Expr<A,T,D1,D,R1>, SMatrix<T,D,D2,R2>,T,D>, T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
411 typedef MatrixMulOp<Expr<A,T,D1,D,R1>,
SMatrix<T,D,D2,R2>,T,D> MatMulOp;
412 return Expr<MatMulOp,T,D1,D2,
419template <
class A,
class B,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
420inline Expr<MatrixMulOp<Expr<A,T,D1,D,R1>, Expr<B,T,D,D2,R2>,T,D>, T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
422 typedef MatrixMulOp<Expr<A,T,D1,D,R1>,
Expr<B,T,D,D2,R2>, T,D> MatMulOp;
432template <
class MatrixA,
class MatrixB,
unsigned int D>
436 MatrixMulOp(
const MatrixA& lhs,
const MatrixB& rhs) :
443 inline typename MatrixA::value_type
apply(
unsigned int i)
const {
456template <
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
457inline Expr<MatrixMulOp<SMatrix<T,D1,D,R1>, SMatrix<T,D,D2,R2>, D>,
T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
458 operator*(
const SMatrix<T,D1,D,R1>& lhs,
const SMatrix<T,D,D2,R2>& rhs) {
459 typedef MatrixMulOp<SMatrix<T,D1,D,R1>, SMatrix<T,D,D2,R2>, D> MatMulOp;
460 return Expr<MatMulOp,T,D1,D2,typename MultPolicy<T,R1,R2>::RepType>(MatMulOp(lhs,rhs));
466template <
class A,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
467inline Expr<MatrixMulOp<SMatrix<T,D1,D,R1>, Expr<A,T,D,D2,R2>, D>,
T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
468 operator*(
const SMatrix<T,D1,D,R1>& lhs,
const Expr<A,T,D,D2,R2>& rhs) {
469 typedef MatrixMulOp<SMatrix<T,D1,D,R1>, Expr<A,T,D,D2,R2>, D> MatMulOp;
470 return Expr<MatMulOp,T,D1,D2,typename MultPolicy<T,R1,R2>::RepType>(MatMulOp(lhs,rhs));
476template <
class A,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
477inline Expr<MatrixMulOp<Expr<A,T,D1,D,R1>, SMatrix<T,D,D2,R2>, D>,
T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
478 operator*(
const Expr<A,T,D1,D,R1>& lhs,
const SMatrix<T,D,D2,R2>& rhs) {
479 typedef MatrixMulOp<Expr<A,T,D1,D,R1>, SMatrix<T,D,D2,R2>, D> MatMulOp;
480 return Expr<MatMulOp,T,D1,D2,typename MultPolicy<T,R1,R2>::RepType>(MatMulOp(lhs,rhs));
486template <
class A,
class B,
class T,
unsigned int D1,
unsigned int D,
unsigned int D2,
class R1,
class R2>
487inline Expr<MatrixMulOp<Expr<A,T,D1,D,R1>, Expr<B,T,D,D2,R2>, D>,
T, D1, D2,
typename MultPolicy<T,R1,R2>::RepType>
488 operator*(
const Expr<A,T,D1,D,R1>& lhs,
const Expr<B,T,D,D2,R2>& rhs) {
489 typedef MatrixMulOp<Expr<A,T,D1,D,R1>, Expr<B,T,D,D2,R2>, D> MatMulOp;
490 return Expr<MatMulOp,T,D1,D2,typename MultPolicy<T,R1,R2>::RepType>(MatMulOp(lhs,rhs));
502template <
class Matrix,
class T,
unsigned int D1,
unsigned int D2=D1>
513 inline T
apply(
unsigned int i)
const {
514 return rhs_.apply( (i%D1)*D2 + i/D1);
521 return rhs_.IsInUse(
p);
538template <
class T,
unsigned int D1,
unsigned int D2,
class R>
539inline Expr<TransposeOp<SMatrix<T,D1,D2,R>,T,D1,D2>, T, D2, D1,
typename TranspPolicy<T,D1,D2,R>::RepType>
549template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
550inline Expr<TransposeOp<Expr<A,T,D1,D2,R>,T,D1,D2>, T, D2, D1,
typename TranspPolicy<T,D1,D2,R>::RepType>
558#ifdef ENABLE_TEMPORARIES_TRANSPOSE
564template <
class T,
unsigned int D1,
unsigned int D2,
class R>
565inline SMatrix< T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>
566 Transpose(
const SMatrix<T,D1,D2, R>& rhs) {
567 typedef TransposeOp<SMatrix<T,D1,D2,R>,T,D1,D2> MatTrOp;
569 return SMatrix< T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>
570 ( Expr<MatTrOp, T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>(MatTrOp(rhs)) );
576template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
577inline SMatrix< T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>
578 Transpose(
const Expr<A,T,D1,D2,R>& rhs) {
579 typedef TransposeOp<Expr<A,T,D1,D2,R>,T,D1,D2> MatTrOp;
581 return SMatrix< T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>
582 ( Expr<MatTrOp, T, D2, D1, typename TranspPolicy<T,D1,D2,R>::RepType>(MatTrOp(rhs)) );
592template <
class T,
unsigned int D,
class R>
593inline T Product(
const SMatrix<T,D,D,R>& lhs,
const SVector<T,D>& rhs) {
594 return Dot(rhs, lhs * rhs);
600template <
class T,
unsigned int D,
class R>
601inline T Product(
const SVector<T,D>& lhs,
const SMatrix<T,D,D,R>& rhs) {
602 return Dot(lhs, rhs * lhs);
608template <
class A,
class T,
unsigned int D,
class R>
609inline T Product(
const SMatrix<T,D,D,R>& lhs,
const VecExpr<A,T,D>& rhs) {
610 return Dot(rhs, lhs * rhs);
616template <
class A,
class T,
unsigned int D,
class R>
617inline T Product(
const VecExpr<A,T,D>& lhs,
const SMatrix<T,D,D,R>& rhs) {
618 return Dot(lhs, rhs * lhs);
624template <
class A,
class T,
unsigned int D,
class R>
625inline T Product(
const SVector<T,D>& lhs,
const Expr<A,T,D,D,R>& rhs) {
626 return Dot(lhs, rhs * lhs);
632template <
class A,
class T,
unsigned int D,
class R>
633inline T Product(
const Expr<A,T,D,D,R>& lhs,
const SVector<T,D>& rhs) {
634 return Dot(rhs, lhs * rhs);
640template <
class A,
class B,
class T,
unsigned int D,
class R>
641inline T Product(
const Expr<A,T,D,D,R>& lhs,
const VecExpr<B,T,D>& rhs) {
642 return Dot(rhs, lhs * rhs);
648template <
class A,
class B,
class T,
unsigned int D,
class R>
649inline T Product(
const VecExpr<A,T,D>& lhs,
const Expr<B,T,D,D,R>& rhs) {
650 return Dot(lhs, rhs * lhs);
664template <
class T,
unsigned int D,
class R>
666 return Dot(rhs, lhs * rhs);
672template <
class T,
unsigned int D,
class R>
674 return Dot(lhs, rhs * lhs);
680template <
class A,
class T,
unsigned int D,
class R>
682 return Dot(rhs, lhs * rhs);
688template <
class A,
class T,
unsigned int D,
class R>
690 return Dot(lhs, rhs * lhs);
696template <
class A,
class T,
unsigned int D,
class R>
698 return Dot(lhs, rhs * lhs);
704template <
class A,
class T,
unsigned int D,
class R>
706 return Dot(rhs, lhs * rhs);
712template <
class A,
class B,
class T,
unsigned int D,
class R>
714 return Dot(rhs, lhs * rhs);
720template <
class A,
class B,
class T,
unsigned int D,
class R>
722 return Dot(lhs, rhs * lhs);
737template <
class T,
unsigned int D1,
unsigned int D2,
class R>
738inline SMatrix<T,D1,D1,MatRepSym<T,D1> >
Similarity(
const SMatrix<T,D1,D2,R>& lhs,
const SMatrix<T,D2,D2,
MatRepSym<T,D2> >& rhs) {
751template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
752inline SMatrix<T,D1,D1,MatRepSym<T,D1> >
Similarity(
const Expr<A,T,D1,D2,R>& lhs,
const SMatrix<T,D2,D2,
MatRepSym<T,D2> >& rhs) {
766template <
class T,
unsigned int D1>
767inline SMatrix<T,D1,D1,MatRepSym<T,D1> >
Similarity(
const SMatrix<T,D1,D1,MatRepSym<T,D1> >& lhs,
const SMatrix<T,D1,D1,MatRepSym<T,D1> >& rhs) {
768 SMatrix<T,D1,D1, MatRepStd<T,D1,D1> > tmp = lhs * rhs;
769 typedef SMatrix<T,D1,D1,MatRepSym<T,D1> > SMatrixSym;
787template <
class T,
unsigned int D1,
unsigned int D2,
class R>
788inline SMatrix<T,D2,D2,MatRepSym<T,D2> >
SimilarityT(
const SMatrix<T,D1,D2,R>& lhs,
const SMatrix<T,D1,D1,
MatRepSym<T,D1> >& rhs) {
801template <
class A,
class T,
unsigned int D1,
unsigned int D2,
class R>
802inline SMatrix<T,D2,D2,MatRepSym<T,D2> >
SimilarityT(
const Expr<A,T,D1,D2,R>& lhs,
const SMatrix<T,D1,D1,
MatRepSym<T,D1> >& rhs) {
836template <
class Vector1,
class Vector2>
848 inline typename Vector1::value_type
apply(
unsigned int i)
const {
849 return lhs_.apply( i/ Vector2::kSize) *
rhs_.apply( i % Vector2::kSize );
851 inline typename Vector1::value_type
operator() (
unsigned int i,
unsigned j)
const {
852 return lhs_.apply(i) *
rhs_.apply(j);
855 inline bool IsInUse (
const typename Vector1::value_type * )
const {
883template <
class T,
unsigned int D1,
unsigned int D2>
893 template <
class T,
unsigned int D1,
unsigned int D2,
class A>
894inline Expr<TensorMulOp<VecExpr<A,T,D1>, SVector<T,D2> >, T, D1, D2 >
903 template <
class T,
unsigned int D1,
unsigned int D2,
class A>
904inline Expr<TensorMulOp<SVector<T,D1>, VecExpr<A,T,D2> >, T, D1, D2 >
914 template <
class T,
unsigned int D1,
unsigned int D2,
class A,
class B>
915inline Expr<TensorMulOp<VecExpr<A,T,D1>, VecExpr<B,T,D2> >, T, D1, D2 >
928template <
class T,
unsigned int D1,
unsigned int D2>
929inline SMatrix<T,D1,D2,MatRepStd<T, D1, D2>>
TensorProd(
const SVector<T,D1>& lhs,
const SVector<T,D2>& rhs) {
930 SMatrix<T,D1,D2,MatRepStd<T, D1, D2>> tmp;
931 for (
unsigned int i=0; i< D1; ++i)
932 for (
unsigned int j=0; j< D2; ++j) {
933 tmp(i,j) = lhs[i]*rhs[j];
941template <
class T,
unsigned int D1,
unsigned int D2,
class A>
942inline SMatrix<T,D1,D2,MatRepStd<T, D1, D2>>
TensorProd(
const VecExpr<A,T,D1>& lhs,
const SVector<T,D2>& rhs) {
943 SMatrix<T,D1,D2,MatRepStd<T, D1, D2>> tmp;
944 for (
unsigned int i=0; i< D1; ++i)
945 for (
unsigned int j=0; j< D2; ++j)
946 tmp(i,j) = lhs.apply(i) * rhs.apply(j);
953template <
class T,
unsigned int D1,
unsigned int D2,
class A>
954inline SMatrix<T,D1,D2, MatRepStd<T, D1, D2>>
TensorProd(
const SVector<T,D1>& lhs,
const VecExpr<A,T,D2>& rhs) {
955 SMatrix<T,D1,D2,MatRepStd<T, D1, D2>> tmp;
956 for (
unsigned int i=0; i< D1; ++i)
957 for (
unsigned int j=0; j< D2; ++j)
958 tmp(i,j) = lhs.apply(i) * rhs.apply(j);
967template <
class T,
unsigned int D1,
unsigned int D2,
class A,
class B>
968inline SMatrix<T,D1,D2,MatRepStd<T, D1, D2>>
TensorProd(
const VecExpr<A,T,D1>& lhs,
const VecExpr<B,T,D2>& rhs) {
969 SMatrix<T,D1,D2,MatRepStd<T, D1, D2>> tmp;
970 for (
unsigned int i=0; i< D1; ++i)
971 for (
unsigned int j=0; j< D2; ++j)
972 tmp(i,j) = lhs.apply(i) * rhs.apply(j);
984template <
class T,
unsigned int D>
992template <
class T,
unsigned int D>
997 ifail = (ok) ? 0 : -1;
header file containing the templated implementation of matrix inversion routines for use with ROOT's ...
static Double_t Product(const Double_t *x, const Float_t *y)
Product.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
class to compute the Cholesky decomposition of a matrix
bool Solve(V &rhs) const
solves a linear system for the given right hand side
Expression wrapper class for Matrix objects.
MatRepSym Matrix storage representation for a symmetric matrix of dimension NxN This class is a templ...
Class for Matrix-Matrix multiplication.
bool IsInUse(const T *p) const
MatrixMulOp(const MatrixA &lhs, const MatrixB &rhs)
T operator()(unsigned int i, unsigned j) const
T apply(unsigned int i) const
calc
SMatrix: a generic fixed size D1 x D2 Matrix class.
SVector: a generic fixed size Vector class.
Class for Tensor Multiplication (outer product) of two vectors giving a matrix.
bool IsInUse(const typename Vector1::value_type *) const
TensorMulOp(const Vector1 &lhs, const Vector2 &rhs)
Vector1::value_type apply(unsigned int i) const
Vector2::kSize is the number of columns in the resulting matrix.
Vector1::value_type operator()(unsigned int i, unsigned j) const
Class for Transpose Operations.
T operator()(unsigned int i, unsigned j) const
T apply(unsigned int i) const
bool IsInUse(const T *p) const
TransposeOp(const Matrix &rhs)
Expression wrapper class for Vector objects.
Class for Vector-Matrix multiplication.
VectorMatrixColOp(const Vector &lhs, const Matrix &rhs)
bool IsInUse(const T *p) const
Matrix::value_type apply(unsigned int i) const
calc
VectorMatrixRowOp(const Matrix &lhs, const Vector &rhs)
bool IsInUse(const T *p) const
Matrix::value_type apply(unsigned int i) const
calc
T Similarity(const SMatrix< T, D, D, R > &lhs, const SVector< T, D > &rhs)
Similarity Vector - Matrix Product: v^T * A * v returning a scalar value of type T .
Expr< TransposeOp< SMatrix< T, D1, D2, R >, T, D1, D2 >, T, D2, D1, typename TranspPolicy< T, D1, D2, R >::RepType > Transpose(const SMatrix< T, D1, D2, R > &rhs)
Matrix Transpose B(i,j) = A(j,i) returning a matrix expression.
SMatrix< T, D2, D2, MatRepSym< T, D2 > > SimilarityT(const SMatrix< T, D1, D2, R > &lhs, const SMatrix< T, D1, D1, MatRepSym< T, D1 > > &rhs)
Transpose Similarity Matrix Product : B = U^T * A * U for A symmetric returning a symmetric matrix ex...
Expr< TensorMulOp< SVector< T, D1 >, SVector< T, D2 > >, T, D1, D2 > TensorProd(const SVector< T, D1 > &lhs, const SVector< T, D2 > &rhs)
Tensor Vector Product : M(i,j) = v(i) * v(j) returning a matrix expression.
Namespace for new Math classes and functions.
bool SolveChol(SMatrix< T, D, D, MatRepSym< T, D > > &mat, SVector< T, D > &vec)
AxisAngle operator*(RotationX const &r1, AxisAngle const &r2)
Multiplication of an axial rotation by an AxisAngle.
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
static void Evaluate(SMatrix< T, D, D, MatRepSym< T, D > > &lhs, const Expr< A, T, D, D, R > &rhs)
assign a symmetric matrix from an expression
MatRepStd< T, N1, N2 > RepType
MatRepStd< T, N2, N1 > RepType