// @(#)root/smatrix:$Id$ // Authors: T. Glebe, L. Moneta 2005 #ifndef ROOT_Math_Functions #define ROOT_Math_Functions // ******************************************************************** // // source: // // type: source code // // created: 16. Mar 2001 // // author: Thorsten Glebe // HERA-B Collaboration // Max-Planck-Institut fuer Kernphysik // Saupfercheckweg 1 // 69117 Heidelberg // Germany // E-mail: T.Glebe@mpi-hd.mpg.de // // Description: Functions which are not applied like a unary operator // // changes: // 16 Mar 2001 (TG) creation // 03 Apr 2001 (TG) minimum added, doc++ comments added // 07 Apr 2001 (TG) Lmag2, Lmag added // 19 Apr 2001 (TG) added #include // 24 Apr 2001 (TG) added sign() // 26 Jul 2001 (TG) round() added // 27 Sep 2001 (TG) added Expr declaration // // ******************************************************************** #include #include "Math/Expression.h" /** @defgroup TempFunction Generic Template Functions @ingroup SMatrixGroup These functions apply for any type T, such as a scalar, a vector or a matrix. */ /** @defgroup VectFunction Vector Template Functions @ingroup SMatrixGroup These functions apply to SVector types (and also to Vector expressions) and can return a vector expression or a scalar, like in the Dot product, or a matrix, like in the Tensor product */ namespace ROOT { namespace Math { template class SVector; /** square Template function to compute \f$x\cdot x \f$, for any type T returning a type T @ingroup TempFunction @author T. Glebe */ //============================================================================== // square: x*x //============================================================================== template inline const T Square(const T& x) { return x*x; } /** maximum. Template to find max(a,b) where a,b are of type T @ingroup TempFunction @author T. Glebe */ //============================================================================== // maximum //============================================================================== template inline const T Maximum(const T& lhs, const T& rhs) { return (lhs > rhs) ? lhs : rhs; } /** minimum. Template to find min(a,b) where a,b are of type T @ingroup TempFunction @author T. Glebe */ //============================================================================== // minimum //============================================================================== template inline const T Minimum(const T& lhs, const T& rhs) { return (lhs < rhs) ? lhs : rhs; } /** round. Template to compute nearest integer value for any type T @ingroup TempFunction @author T. Glebe */ //============================================================================== // round //============================================================================== template inline int Round(const T& x) { return (x-static_cast(x) < 0.5) ? static_cast(x) : static_cast(x+1); } /** sign. Template to compute the sign of a number @ingroup TempFunction @author T. Glebe */ //============================================================================== // sign //============================================================================== template inline int Sign(const T& x) { return (x==0)? 0 : (x<0)? -1 : 1; } //============================================================================== // meta_dot //============================================================================== template struct meta_dot { template static inline T f(const A& lhs, const B& rhs, const T& x) { return lhs.apply(I) * rhs.apply(I) + meta_dot::f(lhs,rhs,x); } }; //============================================================================== // meta_dot<0> //============================================================================== template <> struct meta_dot<0> { template static inline T f(const A& lhs, const B& rhs, const T& /*x */) { return lhs.apply(0) * rhs.apply(0); } }; /** Vector dot product. Template to compute \f$\vec{a}\cdot\vec{b} = \sum_i a_i\cdot b_i \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // dot //============================================================================== template inline T Dot(const SVector& lhs, const SVector& rhs) { return meta_dot::f(lhs,rhs, T()); } //============================================================================== // dot //============================================================================== template inline T Dot(const SVector& lhs, const VecExpr& rhs) { return meta_dot::f(lhs,rhs, T()); } //============================================================================== // dot //============================================================================== template inline T Dot(const VecExpr& lhs, const SVector& rhs) { return meta_dot::f(lhs,rhs, T()); } //============================================================================== // dot //============================================================================== template inline T Dot(const VecExpr& lhs, const VecExpr& rhs) { return meta_dot::f(rhs,lhs, T()); } //============================================================================== // meta_mag //============================================================================== template struct meta_mag { template static inline T f(const A& rhs, const T& x) { return Square(rhs.apply(I)) + meta_mag::f(rhs, x); } }; //============================================================================== // meta_mag<0> //============================================================================== template <> struct meta_mag<0> { template static inline T f(const A& rhs, const T& ) { return Square(rhs.apply(0)); } }; /** Vector magnitude square Template to compute \f$|\vec{v}|^2 = \sum_iv_i^2 \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // mag2 //============================================================================== template inline T Mag2(const SVector& rhs) { return meta_mag::f(rhs, T()); } //============================================================================== // mag2 //============================================================================== template inline T Mag2(const VecExpr& rhs) { return meta_mag::f(rhs, T()); } /** Vector magnitude (Euclidian norm) Compute : \f$ |\vec{v}| = \sqrt{\sum_iv_i^2} \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // mag //============================================================================== template inline T Mag(const SVector& rhs) { return std::sqrt(Mag2(rhs)); } //============================================================================== // mag //============================================================================== template inline T Mag(const VecExpr& rhs) { return std::sqrt(Mag2(rhs)); } /** Lmag2: Square of Minkowski Lorentz-Vector norm (only for 4D Vectors) Template to compute \f$ |\vec{v}|^2 = v_0^2 - v_1^2 - v_2^2 -v_3^2 \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // Lmag2 //============================================================================== template inline T Lmag2(const SVector& rhs) { return Square(rhs[0]) - Square(rhs[1]) - Square(rhs[2]) - Square(rhs[3]); } //============================================================================== // Lmag2 //============================================================================== template inline T Lmag2(const VecExpr& rhs) { return Square(rhs.apply(0)) - Square(rhs.apply(1)) - Square(rhs.apply(2)) - Square(rhs.apply(3)); } /** Lmag: Minkowski Lorentz-Vector norm (only for 4-dim vectors) Length of a vector Lorentz-Vector: \f$ |\vec{v}| = \sqrt{v_0^2 - v_1^2 - v_2^2 -v_3^2} \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // Lmag //============================================================================== template inline T Lmag(const SVector& rhs) { return std::sqrt(Lmag2(rhs)); } //============================================================================== // Lmag //============================================================================== template inline T Lmag(const VecExpr& rhs) { return std::sqrt(Lmag2(rhs)); } /** Vector Cross Product (only for 3-dim vectors) \f$ \vec{c} = \vec{a}\times\vec{b} \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // cross product //============================================================================== template inline SVector Cross(const SVector& lhs, const SVector& rhs) { return SVector(lhs.apply(1)*rhs.apply(2) - lhs.apply(2)*rhs.apply(1), lhs.apply(2)*rhs.apply(0) - lhs.apply(0)*rhs.apply(2), lhs.apply(0)*rhs.apply(1) - lhs.apply(1)*rhs.apply(0)); } //============================================================================== // cross product //============================================================================== template inline SVector Cross(const VecExpr& lhs, const SVector& rhs) { return SVector(lhs.apply(1)*rhs.apply(2) - lhs.apply(2)*rhs.apply(1), lhs.apply(2)*rhs.apply(0) - lhs.apply(0)*rhs.apply(2), lhs.apply(0)*rhs.apply(1) - lhs.apply(1)*rhs.apply(0)); } //============================================================================== // cross product //============================================================================== template inline SVector Cross(const SVector& lhs, const VecExpr& rhs) { return SVector(lhs.apply(1)*rhs.apply(2) - lhs.apply(2)*rhs.apply(1), lhs.apply(2)*rhs.apply(0) - lhs.apply(0)*rhs.apply(2), lhs.apply(0)*rhs.apply(1) - lhs.apply(1)*rhs.apply(0)); } //============================================================================== // cross product //============================================================================== template inline SVector Cross(const VecExpr& lhs, const VecExpr& rhs) { return SVector(lhs.apply(1)*rhs.apply(2) - lhs.apply(2)*rhs.apply(1), lhs.apply(2)*rhs.apply(0) - lhs.apply(0)*rhs.apply(2), lhs.apply(0)*rhs.apply(1) - lhs.apply(1)*rhs.apply(0)); } /** Unit. Return a vector of unit length: \f$ \vec{e}_v = \vec{v}/|\vec{v}| \f$. @ingroup VectFunction @author T. Glebe */ //============================================================================== // unit: returns a unit vector //============================================================================== template inline SVector Unit(const SVector& rhs) { return SVector(rhs).Unit(); } //============================================================================== // unit: returns a unit vector //============================================================================== template inline SVector Unit(const VecExpr& rhs) { return SVector(rhs).Unit(); } #ifdef XXX //============================================================================== // unit: returns a unit vector (worse performance) //============================================================================== template inline VecExpr, SVector, Constant, T>, T, D> unit(const SVector& lhs) { typedef BinaryOp, SVector, Constant, T> DivOpBinOp; return VecExpr(DivOpBinOp(DivOp(),lhs,Constant(mag(lhs)))); } //============================================================================== // unit: returns a unit vector (worse performance) //============================================================================== template inline VecExpr, VecExpr, Constant, T>, T, D> unit(const VecExpr& lhs) { typedef BinaryOp, VecExpr, Constant, T> DivOpBinOp; return VecExpr(DivOpBinOp(DivOp(),lhs,Constant(mag(lhs)))); } #endif } // namespace Math } // namespace ROOT #endif /* ROOT_Math_Functions */