Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMath Namespace Reference

TMath. More...

Namespaces

namespace  ROOTDict
 

Classes

struct  Limits
 

Functions

Double_t Abs (Double_t d)
 Returns the absolute value of parameter Double_t d.
 
Float_t Abs (Float_t d)
 Returns the absolute value of parameter Float_t d.
 
Int_t Abs (Int_t d)
 Returns the absolute value of parameter Int_t d.
 
Long64_t Abs (Long64_t d)
 Returns the absolute value of parameter Long64_t d.
 
Long_t Abs (Long_t d)
 Returns the absolute value of parameter Long_t d.
 
LongDouble_t Abs (LongDouble_t d)
 Returns the absolute value of parameter LongDouble_t d.
 
Short_t Abs (Short_t d)
 Returns the absolute value of parameter Short_t d.
 
Double_t ACos (Double_t)
 Returns the principal value of the arc cosine of x, expressed in radians.
 
Double_t ACosH (Double_t)
 Returns the nonnegative area hyperbolic cosine of x.
 
Bool_t AreEqualAbs (Double_t af, Double_t bf, Double_t epsilon)
 Comparing floating points.
 
Bool_t AreEqualRel (Double_t af, Double_t bf, Double_t relPrec)
 Comparing floating points.
 
Double_t ASin (Double_t)
 Returns the principal value of the arc sine of x, expressed in radians.
 
Double_t ASinH (Double_t)
 Returns the area hyperbolic sine of x.
 
Double_t ATan (Double_t)
 Returns the principal value of the arc tangent of x, expressed in radians.
 
Double_t ATan2 (Double_t y, Double_t x)
 Returns the principal value of the arc tangent of y/x, expressed in radians.
 
Double_t ATanH (Double_t)
 Returns the area hyperbolic tangent of x.
 
Double_t BesselI (Int_t n, Double_t x)
 Computes the Integer Order Modified Bessel function I_n(x) for n=0,1,2,... and any real x.
 
Double_t BesselI0 (Double_t x)
 Integer order modified Bessel function K_n(x)
 
Double_t BesselI1 (Double_t x)
 Modified Bessel function K_0(x)
 
Double_t BesselJ0 (Double_t x)
 Modified Bessel function K_1(x)
 
Double_t BesselJ1 (Double_t x)
 Bessel function J0(x) for any real x.
 
Double_t BesselK (Int_t n, Double_t x)
 Integer order modified Bessel function I_n(x)
 
Double_t BesselK0 (Double_t x)
 Modified Bessel function I_0(x)
 
Double_t BesselK1 (Double_t x)
 Modified Bessel function I_1(x)
 
Double_t BesselY0 (Double_t x)
 Bessel function J1(x) for any real x.
 
Double_t BesselY1 (Double_t x)
 Bessel function Y0(x) for positive x.
 
Double_t Beta (Double_t p, Double_t q)
 Calculates Beta-function Gamma(p)*Gamma(q)/Gamma(p+q).
 
Double_t BetaCf (Double_t x, Double_t a, Double_t b)
 Continued fraction evaluation by modified Lentz's method used in calculation of incomplete Beta function.
 
Double_t BetaDist (Double_t x, Double_t p, Double_t q)
 Computes the probability density function of the Beta distribution (the cumulative distribution function is computed in BetaDistI).
 
Double_t BetaDistI (Double_t x, Double_t p, Double_t q)
 Computes the cumulative distribution function of the Beta distribution, i.e.
 
Double_t BetaIncomplete (Double_t x, Double_t a, Double_t b)
 Calculates the incomplete Beta-function.
 
template<typename Iterator , typename Element >
Iterator BinarySearch (Iterator first, Iterator last, Element value)
 Binary search in an array defined by its iterators.
 
template<typename T >
Long64_t BinarySearch (Long64_t n, const T **array, T value)
 Binary search in an array of n values to locate value.
 
template<typename T >
Long64_t BinarySearch (Long64_t n, const T *array, T value)
 Binary search in an array of n values to locate value.
 
Double_t Binomial (Int_t n, Int_t k)
 Calculates the binomial coefficient n over k.
 
Double_t BinomialI (Double_t p, Int_t n, Int_t k)
 Suppose an event occurs with probability p per trial Then the probability P of its occurring k or more times in n trials is termed a cumulative binomial probability the formula is:
 
Double_t BreitWigner (Double_t x, Double_t mean=0, Double_t gamma=1)
 Calculates a Breit Wigner function with mean and gamma.
 
Double_t BreitWignerRelativistic (Double_t x, Double_t median=0, Double_t gamma=1)
 Calculates a Relativistic Breit Wigner function with median and gamma.
 
void BubbleHigh (Int_t Narr, Double_t *arr1, Int_t *arr2)
 Bubble sort variant to obtain the order of an array's elements into an index in order to do more useful things than the standard built in functions.
 
void BubbleLow (Int_t Narr, Double_t *arr1, Int_t *arr2)
 Opposite ordering of the array arr2[] to that of BubbleHigh.
 
constexpr Double_t C ()
 Velocity of light in \( m s^{-1} \).
 
Double_t CauchyDist (Double_t x, Double_t t=0, Double_t s=1)
 Computes the density of Cauchy distribution at point x by default, standard Cauchy distribution is used (t=0, s=1)
 
constexpr Double_t Ccgs ()
 \( cm s^{-1} \)
 
Double_t Ceil (Double_t x)
 Rounds x upward, returning the smallest integral value that is not less than x.
 
Int_t CeilNint (Double_t x)
 Returns the nearest integer of TMath::Ceil(x).
 
Double_t ChisquareQuantile (Double_t p, Double_t ndf)
 Evaluate the quantiles of the chi-squared probability distribution function.
 
Double_t Cos (Double_t)
 Returns the cosine of an angle of x radians.
 
Double_t CosH (Double_t)
 Returns the hyperbolic cosine of x.
 
template<typename T >
T * Cross (const T v1[3], const T v2[3], T out[3])
 Calculates the Cross Product of two vectors: out = [v1 x v2].
 
constexpr Double_t CUncertainty ()
 Speed of light uncertainty.
 
constexpr Double_t DegToRad ()
 Conversion from degree to radian: \( \frac{\pi}{180} \).
 
Double_t DiLog (Double_t x)
 Modified Struve functions of order 1.
 
constexpr Double_t E ()
 Base of natural log: \( e \).
 
Double_t Erf (Double_t x)
 Computation of the error function erf(x).
 
Double_t Erfc (Double_t x)
 Computes the complementary error function erfc(x).
 
Double_t ErfcInverse (Double_t x)
 Returns the inverse of the complementary error function.
 
Double_t ErfInverse (Double_t x)
 Returns the inverse error function.
 
constexpr Double_t EulerGamma ()
 Euler-Mascheroni Constant.
 
Bool_t Even (Long_t a)
 Returns true if a is even.
 
Double_t Exp (Double_t x)
 Returns the base-e exponential function of x, which is e raised to the power x.
 
Double_t Factorial (Int_t i)
 Computes factorial(n).
 
Double_t FDist (Double_t F, Double_t N, Double_t M)
 Computes the density function of F-distribution (probability function, integral of density, is computed in FDistI).
 
Double_t FDistI (Double_t F, Double_t N, Double_t M)
 Calculates the cumulative distribution function of F-distribution (see ROOT::Math::fdistribution_cdf).
 
Int_t Finite (Double_t x)
 Check if it is finite with a mask in order to be consistent in presence of fast math.
 
Int_t Finite (Float_t x)
 Check if it is finite with a mask in order to be consistent in presence of fast math.
 
Double_t Floor (Double_t x)
 Rounds x downward, returning the largest integral value that is not greater than x.
 
Int_t FloorNint (Double_t x)
 Returns the nearest integer of TMath::Floor(x).
 
Double_t Freq (Double_t x)
 Computation of the normal frequency function freq(x).
 
constexpr Double_t G ()
 Gravitational constant in: \( m^{3} kg^{-1} s^{-2} \).
 
Double_t GamCf (Double_t a, Double_t x)
 Computation of the incomplete gamma function P(a,x) via its continued fraction representation.
 
Double_t Gamma (Double_t a, Double_t x)
 Computation of the normalized lower incomplete gamma function P(a,x) as defined in the Handbook of Mathematical Functions by Abramowitz and Stegun, formula 6.5.1 on page 260 .
 
Double_t Gamma (Double_t z)
 Computation of gamma(z) for all z.
 
Double_t GammaDist (Double_t x, Double_t gamma, Double_t mu=0, Double_t beta=1)
 Computes the density function of Gamma distribution at point x.
 
Double_t GamSer (Double_t a, Double_t x)
 Computation of the incomplete gamma function P(a,x) via its series representation.
 
Double_t Gaus (Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
 Calculates a gaussian function with mean and sigma.
 
constexpr Double_t Gcgs ()
 \( cm^{3} g^{-1} s^{-2} \)
 
template<typename Iterator >
Double_t GeomMean (Iterator first, Iterator last)
 Returns the geometric mean of an array defined by the iterators.
 
template<typename T >
Double_t GeomMean (Long64_t n, const T *a)
 Returns the geometric mean of an array a of size n.
 
constexpr Double_t GhbarC ()
 \( \frac{G}{\hbar C} \) in \( (GeV/c^{2})^{-2} \)
 
constexpr Double_t GhbarCUncertainty ()
 \( \frac{G}{\hbar C} \) uncertainty.
 
constexpr Double_t Gn ()
 Standard acceleration of gravity in \( m s^{-2} \).
 
constexpr Double_t GnUncertainty ()
 Standard acceleration of gravity uncertainty.
 
constexpr Double_t GUncertainty ()
 Gravitational constant uncertainty.
 
constexpr Double_t H ()
 Planck's constant in \( J s \): \( h \).
 
ULong_t Hash (const char *str)
 
ULong_t Hash (const void *txt, Int_t ntxt)
 Calculates hash index from any char string.
 
constexpr Double_t Hbar ()
 \( \hbar \) in \( J s \): \( \hbar = \frac{h}{2\pi} \)
 
constexpr Double_t Hbarcgs ()
 \( erg s \)
 
constexpr Double_t HbarUncertainty ()
 \( \hbar \) uncertainty.
 
constexpr Double_t HC ()
 \( hc \) in \( J m \)
 
constexpr Double_t HCcgs ()
 \( erg cm \)
 
constexpr Double_t Hcgs ()
 \( erg s \)
 
constexpr Double_t HUncertainty ()
 Planck's constant uncertainty.
 
Double_t Hypot (Double_t x, Double_t y)
 Returns sqrt(x*x + y*y)
 
Long_t Hypot (Long_t x, Long_t y)
 Returns sqrt(x*x + y*y)
 
Double_t Infinity ()
 Returns an infinity as defined by the IEEE standard.
 
constexpr Double_t InvPi ()
 \( \frac{1.}{\pi}\)
 
template<typename T >
Bool_t IsInside (T xp, T yp, Int_t np, T *x, T *y)
 Function which returns kTRUE if point xp,yp lies inside the polygon defined by the np points in arrays x and y, kFALSE otherwise.
 
Bool_t IsNaN (Double_t x)
 
Bool_t IsNaN (Float_t x)
 
constexpr Double_t K ()
 Boltzmann's constant in \( J K^{-1} \): \( k \).
 
constexpr Double_t Kcgs ()
 \( erg K^{-1} \)
 
Double_t KolmogorovProb (Double_t z)
 Calculates the Kolmogorov distribution function,.
 
Double_t KolmogorovTest (Int_t na, const Double_t *a, Int_t nb, const Double_t *b, Option_t *option)
 Statistical test whether two one-dimensional sets of points are compatible with coming from the same parent distribution, using the Kolmogorov test.
 
template<class Element , typename Size >
Element KOrdStat (Size n, const Element *a, Size k, Size *work=0)
 Returns k_th order statistic of the array a of size n (k_th smallest element out of n elements).
 
constexpr Double_t KUncertainty ()
 Boltzmann's constant uncertainty.
 
Double_t Landau (Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
 The LANDAU function.
 
Double_t LandauI (Double_t x)
 Returns the cumulative (lower tail integral) of the Landau distribution function at point x.
 
Double_t LaplaceDist (Double_t x, Double_t alpha=0, Double_t beta=1)
 Computes the probability density function of Laplace distribution at point x, with location parameter alpha and shape parameter beta.
 
Double_t LaplaceDistI (Double_t x, Double_t alpha=0, Double_t beta=1)
 Computes the cumulative distribution function (lower tail integral) of Laplace distribution at point x, with location parameter alpha and shape parameter beta.
 
Double_t Ldexp (Double_t x, Int_t exp)
 Returns the result of multiplying x (the significant) by 2 raised to the power of exp (the exponent).
 
constexpr Double_t Ln10 ()
 Natural log of 10 (to convert log to ln)
 
Double_t LnGamma (Double_t z)
 Computation of ln[gamma(z)] for all z.
 
template<typename Iterator >
Iterator LocMax (Iterator first, Iterator last)
 Returns index of array with the maximum element.
 
template<typename T >
Long64_t LocMax (Long64_t n, const T *a)
 Returns index of array with the maximum element.
 
template<typename Iterator >
Iterator LocMin (Iterator first, Iterator last)
 Returns index of array with the minimum element.
 
template<typename T >
Long64_t LocMin (Long64_t n, const T *a)
 Returns index of array with the minimum element.
 
Double_t Log (Double_t x)
 Returns the natural logarithm of x.
 
Double_t Log10 (Double_t x)
 Returns the common (base-10) logarithm of x.
 
Double_t Log2 (Double_t x)
 Returns the binary (base-2) logarithm of x.
 
constexpr Double_t LogE ()
 Base-10 log of e (to convert ln to log)
 
Double_t LogNormal (Double_t x, Double_t sigma, Double_t theta=0, Double_t m=1)
 Computes the density of LogNormal distribution at point x.
 
Double_t Max (Double_t a, Double_t b)
 Returns the largest of a and b.
 
Float_t Max (Float_t a, Float_t b)
 Returns the largest of a and b.
 
Int_t Max (Int_t a, Int_t b)
 Returns the largest of a and b.
 
Long64_t Max (Long64_t a, Long64_t b)
 Returns the largest of a and b.
 
Long_t Max (Long_t a, Long_t b)
 Returns the largest of a and b.
 
Short_t Max (Short_t a, Short_t b)
 Returns the largest of a and b.
 
UInt_t Max (UInt_t a, UInt_t b)
 Returns the largest of a and b.
 
ULong64_t Max (ULong64_t a, ULong64_t b)
 Returns the largest of a and b.
 
ULong_t Max (ULong_t a, ULong_t b)
 Returns the largest of a and b.
 
UShort_t Max (UShort_t a, UShort_t b)
 Returns the largest of a and b.
 
template<typename T >
MaxElement (Long64_t n, const T *a)
 Returns maximum of array a of length n.
 
template<typename Iterator >
Double_t Mean (Iterator first, Iterator last)
 Returns the weighted mean of an array defined by the iterators.
 
template<typename Iterator , typename WeightIterator >
Double_t Mean (Iterator first, Iterator last, WeightIterator wfirst)
 Returns the weighted mean of an array defined by the first and last iterators.
 
template<typename T >
Double_t Mean (Long64_t n, const T *a, const Double_t *w=nullptr)
 Returns the weighted mean of an array a with length n.
 
template<typename T >
Double_t Median (Long64_t n, const T *a, const Double_t *w=nullptr, Long64_t *work=nullptr)
 Same as RMS.
 
Double_t Min (Double_t a, Double_t b)
 Returns the smallest of a and b.
 
Float_t Min (Float_t a, Float_t b)
 Returns the smallest of a and b.
 
Int_t Min (Int_t a, Int_t b)
 Returns the smallest of a and b.
 
Long64_t Min (Long64_t a, Long64_t b)
 Returns the smallest of a and b.
 
Long_t Min (Long_t a, Long_t b)
 Returns the smallest of a and b.
 
Short_t Min (Short_t a, Short_t b)
 Returns the smallest of a and b.
 
UInt_t Min (UInt_t a, UInt_t b)
 Returns the smallest of a and b.
 
ULong64_t Min (ULong64_t a, ULong64_t b)
 Returns the smallest of a and b.
 
ULong_t Min (ULong_t a, ULong_t b)
 Returns the smallest of a and b.
 
UShort_t Min (UShort_t a, UShort_t b)
 Returns the smallest of a and b.
 
template<typename T >
MinElement (Long64_t n, const T *a)
 Returns minimum of array a of length n.
 
constexpr Double_t MWair ()
 Molecular weight of dry air 1976 US Standard Atmosphere in \( kg kmol^{-1} \) or \( gm mol^{-1} \)
 
constexpr Double_t Na ()
 Avogadro constant (Avogadro's Number) in \( mol^{-1} \).
 
constexpr Double_t NaUncertainty ()
 Avogadro constant (Avogadro's Number) uncertainty.
 
Long_t NextPrime (Long_t x)
 
template<typename T >
Int_t Nint (T x)
 Round to nearest integer. Rounds half integers to the nearest even integer.
 
template<typename T >
T * Normal2Plane (const T v1[3], const T v2[3], const T v3[3], T normal[3])
 Calculates a normal vector of a plane.
 
Double_t Normalize (Double_t v[3])
 Normalize a vector v in place.
 
Float_t Normalize (Float_t v[3])
 Normalize a vector v in place.
 
template<typename T >
NormCross (const T v1[3], const T v2[3], T out[3])
 Calculates the Normalized Cross Product of two vectors.
 
Double_t NormQuantile (Double_t p)
 Computes quantiles for standard normal distribution N(0, 1) at probability p.
 
Bool_t Odd (Long_t a)
 Returns true if a is odd.
 
Bool_t Permute (Int_t n, Int_t *a)
 Simple recursive algorithm to find the permutations of n natural numbers, not necessarily all distinct adapted from CERNLIB routine PERMU.
 
constexpr Double_t Pi ()
 \( \pi\)
 
constexpr Double_t PiOver2 ()
 \( \frac{\pi}{2} \)
 
constexpr Double_t PiOver4 ()
 \( \frac{\pi}{4} \)
 
Double_t Poisson (Double_t x, Double_t par)
 Computes the Poisson distribution function for (x,par).
 
Double_t PoissonI (Double_t x, Double_t par)
 Computes the Discrete Poisson distribution function for (x,par).
 
Double_t Power (Double_t x, Double_t y)
 Returns x raised to the power y.
 
Double_t Power (Double_t x, Int_t y)
 Returns x raised to the power y.
 
LongDouble_t Power (Long64_t x, Long64_t y)
 Returns x raised to the power y.
 
LongDouble_t Power (LongDouble_t x, Long64_t y)
 Returns x raised to the power y.
 
LongDouble_t Power (LongDouble_t x, LongDouble_t y)
 Returns x raised to the power y.
 
Double_t Prob (Double_t chi2, Int_t ndf)
 Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf).
 
constexpr Double_t Qe ()
 Elementary charge in \( C \) .
 
constexpr Double_t QeUncertainty ()
 Elementary charge uncertainty.
 
void Quantiles (Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=nullptr, Int_t type=7)
 Computes sample quantiles, corresponding to the given probabilities.
 
Double_t QuietNaN ()
 Returns a quiet NaN as defined by IEEE 754.
 
constexpr Double_t R ()
 Universal gas constant ( \( Na K \)) in \( J K^{-1} mol^{-1} \)
 
constexpr Double_t RadToDeg ()
 Conversion from radian to degree: \( \frac{180}{\pi} \).
 
Double_t Range (Double_t lb, Double_t ub, Double_t x)
 Returns x if lb < x < up, lb if x < lb and ub if x > ub.
 
Int_t Range (Int_t lb, Int_t ub, Int_t x)
 Returns x if lb < x < up, lb if x < lb and ub if x > ub.
 
Long_t Range (Long_t lb, Long_t ub, Long_t x)
 Returns x if lb < x < up, lb if x < lb and ub if x > ub.
 
Short_t Range (Short_t lb, Short_t ub, Short_t x)
 Returns x if lb < x < up, lb if x < lb and ub if x > ub.
 
ULong_t Range (ULong_t lb, ULong_t ub, ULong_t x)
 Returns x if lb < x < up, lb if x < lb and ub if x > ub.
 
constexpr Double_t Rgair ()
 Dry Air Gas Constant (R / MWair) in \( J kg^{-1} K^{-1} \)
 
template<typename Iterator >
Double_t RMS (Iterator first, Iterator last)
 Returns the Standard Deviation of an array defined by the iterators.
 
template<typename Iterator , typename WeightIterator >
Double_t RMS (Iterator first, Iterator last, WeightIterator wfirst)
 Returns the weighted Standard Deviation of an array defined by the iterators.
 
template<typename T >
Double_t RMS (Long64_t n, const T *a, const Double_t *w=nullptr)
 Returns the Standard Deviation of an array a with length n.
 
Bool_t RootsCubic (const Double_t coef[4], Double_t &a, Double_t &b, Double_t &c)
 Calculates roots of polynomial of 3rd order a*x^3 + b*x^2 + c*x + d, where.
 
constexpr Double_t RUncertainty ()
 Universal gas constant uncertainty.
 
constexpr Double_t Sigma ()
 Stefan-Boltzmann constant in \( W m^{-2} K^{-4}\): \( \sigma \).
 
constexpr Double_t SigmaUncertainty ()
 Stefan-Boltzmann constant uncertainty.
 
Double_t Sign (Double_t a, Double_t b)
 Returns a value with the magnitude of a and the sign of b.
 
Float_t Sign (Float_t a, Float_t b)
 Returns a value with the magnitude of a and the sign of b.
 
LongDouble_t Sign (LongDouble_t a, LongDouble_t b)
 Returns a value with the magnitude of a and the sign of b.
 
template<typename T1 , typename T2 >
T1 Sign (T1 a, T2 b)
 Returns a value with the magnitude of a and the sign of b.
 
Double_t SignalingNaN ()
 Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
 
Bool_t SignBit (Double_t a)
 Returns whether the sign of Double_t a is negative.
 
Bool_t SignBit (Float_t a)
 Returns whether the sign of Float_t a is negative.
 
template<typename Integer >
Bool_t SignBit (Integer a)
 Returns whether the sign of Integer a is negative.
 
Bool_t SignBit (LongDouble_t a)
 Returns whether the sign of LongDouble_t a is negative.
 
Double_t Sin (Double_t)
 Returns the sine of an angle of x radians.
 
Double_t SinH (Double_t)
 Returns the hyperbolic sine of `x.
 
template<typename Element , typename Index >
void Sort (Index n, const Element *a, Index *index, Bool_t down=kTRUE)
 Sort the n elements of the array a of generic templated type Element.
 
template<typename Iterator , typename IndexIterator >
void SortItr (Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE)
 Sort the n1 elements of the Short_t array defined by its iterators.
 
Double_t Sq (Double_t x)
 Returns x*x.
 
Double_t Sqrt (Double_t x)
 Returns the square root of x.
 
constexpr Double_t Sqrt2 ()
 \( \sqrt{2} \)
 
template<typename Iterator >
Double_t StdDev (Iterator first, Iterator last)
 Same as RMS.
 
template<typename Iterator , typename WeightIterator >
Double_t StdDev (Iterator first, Iterator last, WeightIterator wfirst)
 Same as RMS.
 
template<typename T >
Double_t StdDev (Long64_t n, const T *a, const Double_t *w=nullptr)
 
Double_t StruveH0 (Double_t x)
 Bessel function Y1(x) for positive x.
 
Double_t StruveH1 (Double_t x)
 Struve functions of order 0.
 
Double_t StruveL0 (Double_t x)
 Struve functions of order 1.
 
Double_t StruveL1 (Double_t x)
 Modified Struve functions of order 0.
 
Double_t Student (Double_t T, Double_t ndf)
 Computes density function for Student's t- distribution (the probability function (integral of density) is computed in StudentI).
 
Double_t StudentI (Double_t T, Double_t ndf)
 Calculates the cumulative distribution function of Student's t-distribution second parameter stands for number of degrees of freedom, not for the number of samples if x has Student's t-distribution, the function returns the probability of x being less than T.
 
Double_t StudentQuantile (Double_t p, Double_t ndf, Bool_t lower_tail=kTRUE)
 Computes quantiles of the Student's t-distribution 1st argument is the probability, at which the quantile is computed 2nd argument - the number of degrees of freedom of the Student distribution When the 3rd argument lower_tail is kTRUE (default)- the algorithm returns such x0, that.
 
Double_t Tan (Double_t)
 Returns the tangent of an angle of x radians.
 
Double_t TanH (Double_t)
 Returns the hyperbolic tangent of x.
 
constexpr Double_t TwoPi ()
 \( 2\pi\)
 
Double_t Vavilov (Double_t x, Double_t kappa, Double_t beta2)
 Returns the value of the Vavilov probability density function.
 
Double_t VavilovDenEval (Double_t rlam, Double_t *AC, Double_t *HC, Int_t itype)
 Internal function, called by Vavilov and VavilovSet.
 
Double_t VavilovI (Double_t x, Double_t kappa, Double_t beta2)
 Returns the value of the Vavilov cumulative distribution function (lower tail integral of the probability distribution function)
 
void VavilovSet (Double_t rkappa, Double_t beta2, Bool_t mode, Double_t *WCM, Double_t *AC, Double_t *HC, Int_t &itype, Int_t &npt)
 Internal function, called by Vavilov and VavilovI.
 
Double_t Voigt (Double_t x, Double_t sigma, Double_t lg, Int_t r=4)
 Computation of Voigt function (normalised).
 

Detailed Description

TMath.

Encapsulate most frequently used Math functions. NB. The basic functions Min, Max, Abs and Sign are defined in TMathBase.

Function Documentation

◆ Abs() [1/7]

Double_t TMath::Abs ( Double_t  d)
inline

Returns the absolute value of parameter Double_t d.

Definition at line 143 of file TMathBase.h.

◆ Abs() [2/7]

Float_t TMath::Abs ( Float_t  d)
inline

Returns the absolute value of parameter Float_t d.

Definition at line 139 of file TMathBase.h.

◆ Abs() [3/7]

Int_t TMath::Abs ( Int_t  d)
inline

Returns the absolute value of parameter Int_t d.

Definition at line 127 of file TMathBase.h.

◆ Abs() [4/7]

Long64_t TMath::Abs ( Long64_t  d)
inline

Returns the absolute value of parameter Long64_t d.

Definition at line 135 of file TMathBase.h.

◆ Abs() [5/7]

Long_t TMath::Abs ( Long_t  d)
inline

Returns the absolute value of parameter Long_t d.

Definition at line 131 of file TMathBase.h.

◆ Abs() [6/7]

LongDouble_t TMath::Abs ( LongDouble_t  d)
inline

Returns the absolute value of parameter LongDouble_t d.

Definition at line 147 of file TMathBase.h.

◆ Abs() [7/7]

Short_t TMath::Abs ( Short_t  d)
inline

Returns the absolute value of parameter Short_t d.

Definition at line 123 of file TMathBase.h.

◆ ACos()

Double_t TMath::ACos ( Double_t  x)
inline

Returns the principal value of the arc cosine of x, expressed in radians.

Definition at line 632 of file TMath.h.

◆ ACosH()

Double_t TMath::ACosH ( Double_t  x)

Returns the nonnegative area hyperbolic cosine of x.

Definition at line 81 of file TMath.cxx.

◆ AreEqualAbs()

Bool_t TMath::AreEqualAbs ( Double_t  af,
Double_t  bf,
Double_t  epsilon 
)
inline

Comparing floating points.

Returns kTRUE if the absolute difference between af and bf is less than epsilon.

Definition at line 418 of file TMath.h.

◆ AreEqualRel()

Bool_t TMath::AreEqualRel ( Double_t  af,
Double_t  bf,
Double_t  relPrec 
)
inline

Comparing floating points.

Returns kTRUE if the relative difference between af and bf is less than relPrec.

Definition at line 426 of file TMath.h.

◆ ASin()

Double_t TMath::ASin ( Double_t  x)
inline

Returns the principal value of the arc sine of x, expressed in radians.

Definition at line 624 of file TMath.h.

◆ ASinH()

Double_t TMath::ASinH ( Double_t  x)

Returns the area hyperbolic sine of x.

Definition at line 67 of file TMath.cxx.

◆ ATan()

Double_t TMath::ATan ( Double_t  x)
inline

Returns the principal value of the arc tangent of x, expressed in radians.

Definition at line 640 of file TMath.h.

◆ ATan2()

Double_t TMath::ATan2 ( Double_t  y,
Double_t  x 
)
inline

Returns the principal value of the arc tangent of y/x, expressed in radians.

Definition at line 646 of file TMath.h.

◆ ATanH()

Double_t TMath::ATanH ( Double_t  x)

Returns the area hyperbolic tangent of x.

Definition at line 95 of file TMath.cxx.

◆ BesselI()

Double_t TMath::BesselI ( Int_t  n,
Double_t  x 
)

Computes the Integer Order Modified Bessel function I_n(x) for n=0,1,2,... and any real x.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1590 of file TMath.cxx.

◆ BesselI0()

Double_t TMath::BesselI0 ( Double_t  x)

Integer order modified Bessel function K_n(x)

Computes the modified Bessel function I_0(x) for any real x.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1426 of file TMath.cxx.

◆ BesselI1()

Double_t TMath::BesselI1 ( Double_t  x)

Modified Bessel function K_0(x)

Computes the modified Bessel function I_1(x) for any real x.

M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions, Applied Mathematics Series vol. 55 (1964), Washington.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1494 of file TMath.cxx.

◆ BesselJ0()

Double_t TMath::BesselJ0 ( Double_t  x)

Modified Bessel function K_1(x)

Returns the Bessel function J0(x) for any real x.

Definition at line 1634 of file TMath.cxx.

◆ BesselJ1()

Double_t TMath::BesselJ1 ( Double_t  x)

Bessel function J0(x) for any real x.

Returns the Bessel function J1(x) for any real x.

Definition at line 1669 of file TMath.cxx.

◆ BesselK()

Double_t TMath::BesselK ( Int_t  n,
Double_t  x 
)

Integer order modified Bessel function I_n(x)

Computes the Integer Order Modified Bessel function K_n(x) for n=0,1,2,... and positive real x.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1561 of file TMath.cxx.

◆ BesselK0()

Double_t TMath::BesselK0 ( Double_t  x)

Modified Bessel function I_0(x)

Computes the modified Bessel function K_0(x) for positive real x.

M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions, Applied Mathematics Series vol. 55 (1964), Washington.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1460 of file TMath.cxx.

◆ BesselK1()

Double_t TMath::BesselK1 ( Double_t  x)

Modified Bessel function I_1(x)

Computes the modified Bessel function K_1(x) for positive real x.

M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions, Applied Mathematics Series vol. 55 (1964), Washington.

Author
NvE 12-mar-2000 UU-SAP Utrecht

Definition at line 1529 of file TMath.cxx.

◆ BesselY0()

Double_t TMath::BesselY0 ( Double_t  x)

Bessel function J1(x) for any real x.

Returns the Bessel function Y0(x) for positive x.

Definition at line 1705 of file TMath.cxx.

◆ BesselY1()

Double_t TMath::BesselY1 ( Double_t  x)

Bessel function Y0(x) for positive x.

Returns the Bessel function Y1(x) for positive x.

Definition at line 1739 of file TMath.cxx.

◆ Beta()

Double_t TMath::Beta ( Double_t  p,
Double_t  q 
)

Calculates Beta-function Gamma(p)*Gamma(q)/Gamma(p+q).

Definition at line 2011 of file TMath.cxx.

◆ BetaCf()

Double_t TMath::BetaCf ( Double_t  x,
Double_t  a,
Double_t  b 
)

Continued fraction evaluation by modified Lentz's method used in calculation of incomplete Beta function.

Definition at line 2020 of file TMath.cxx.

◆ BetaDist()

Double_t TMath::BetaDist ( Double_t  x,
Double_t  p,
Double_t  q 
)

Computes the probability density function of the Beta distribution (the cumulative distribution function is computed in BetaDistI).

The first argument is the point, where the function will be computed, second and third are the function parameters. Since the Beta distribution is bounded on both sides, it's often used to represent processes with natural lower and upper limits.

Definition at line 2071 of file TMath.cxx.

◆ BetaDistI()

Double_t TMath::BetaDistI ( Double_t  x,
Double_t  p,
Double_t  q 
)

Computes the cumulative distribution function of the Beta distribution, i.e.

the lower tail integral of TMath::BetaDist The first argument is the point, where the function will be computed, second and third are the function parameters. Since the Beta distribution is bounded on both sides, it's often used to represent processes with natural lower and upper limits.

Definition at line 2090 of file TMath.cxx.

◆ BetaIncomplete()

Double_t TMath::BetaIncomplete ( Double_t  x,
Double_t  a,
Double_t  b 
)

Calculates the incomplete Beta-function.

Definition at line 2103 of file TMath.cxx.

◆ BinarySearch() [1/3]

template<typename Iterator , typename Element >
Iterator TMath::BinarySearch ( Iterator  first,
Iterator  last,
Element  value 
)

Binary search in an array defined by its iterators.

The values in the iterators range are supposed to be sorted prior to this call. If match is found, function returns position of element. If no match found, function gives nearest element smaller than value.

Definition at line 332 of file TMathBase.h.

◆ BinarySearch() [2/3]

template<typename T >
Long64_t TMath::BinarySearch ( Long64_t  n,
const T **  array,
value 
)

Binary search in an array of n values to locate value.

Array is supposed to be sorted prior to this call. If match is found, function returns position of element. If no match found, function gives nearest element smaller than value.

Definition at line 362 of file TMathBase.h.

◆ BinarySearch() [3/3]

template<typename T >
Long64_t TMath::BinarySearch ( Long64_t  n,
const T *  array,
value 
)

Binary search in an array of n values to locate value.

Array is supposed to be sorted prior to this call. If match is found, function returns position of element. If no match found, function gives nearest element smaller than value.

Definition at line 347 of file TMathBase.h.

◆ Binomial()

Double_t TMath::Binomial ( Int_t  n,
Int_t  k 
)

Calculates the binomial coefficient n over k.

Definition at line 2111 of file TMath.cxx.

◆ BinomialI()

Double_t TMath::BinomialI ( Double_t  p,
Int_t  n,
Int_t  k 
)

Suppose an event occurs with probability p per trial Then the probability P of its occurring k or more times in n trials is termed a cumulative binomial probability the formula is:

P = sum_from_j=k_to_n(TMath::Binomial (n, j)**TMath::Power (p, j)*TMath::Power (1-p, n-j)
winID h TVirtualViewer3D TVirtualGLPainter p
const Int_t n
Definition legend1.C:16
Double_t Binomial(Int_t n, Int_t k)
Calculates the binomial coefficient n over k.
Definition TMath.cxx:2111
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721

For n larger than 12 BetaIncomplete is a much better way to evaluate the sum than would be the straightforward sum calculation for n smaller than 12 either method is acceptable ("Numerical Recipes")

Note this function is not exactly implementing the cumulative or the complement of the cumulative of the Binomial distrinution. It is equivalent to ROOT::Math::binomial_cdf_c(k-1,p,n)

Author
Anna Kreshuk

Definition at line 2141 of file TMath.cxx.

◆ BreitWigner()

Double_t TMath::BreitWigner ( Double_t  x,
Double_t  mean = 0,
Double_t  gamma = 1 
)

Calculates a Breit Wigner function with mean and gamma.

Definition at line 442 of file TMath.cxx.

◆ BreitWignerRelativistic()

Double_t TMath::BreitWignerRelativistic ( Double_t  x,
Double_t  median = 0,
Double_t  gamma = 1 
)

Calculates a Relativistic Breit Wigner function with median and gamma.

Definition at line 452 of file TMath.cxx.

◆ BubbleHigh()

void TMath::BubbleHigh ( Int_t  Narr,
Double_t arr1,
Int_t arr2 
)

Bubble sort variant to obtain the order of an array's elements into an index in order to do more useful things than the standard built in functions.

Parameters
[in]Narrnumber of array elements
[in]*arr1is unchanged;
[in]*arr2is the array of indices corresponding to the descending value of arr1 with arr2[0] corresponding to the largest arr1 value and arr2[Narr] the smallest.
Author
Adrian Bevan (bevan.nosp@m.@sla.nosp@m.c.sta.nosp@m.nfor.nosp@m.d.edu)

Definition at line 1314 of file TMath.cxx.

◆ BubbleLow()

void TMath::BubbleLow ( Int_t  Narr,
Double_t arr1,
Int_t arr2 
)

Opposite ordering of the array arr2[] to that of BubbleHigh.

Author
Adrian Bevan (bevan.nosp@m.@sla.nosp@m.c.sta.nosp@m.nfor.nosp@m.d.edu)

Definition at line 1353 of file TMath.cxx.

◆ C()

constexpr Double_t TMath::C ( )
constexpr

Velocity of light in \( m s^{-1} \).

Definition at line 114 of file TMath.h.

◆ CauchyDist()

Double_t TMath::CauchyDist ( Double_t  x,
Double_t  t = 0,
Double_t  s = 1 
)

Computes the density of Cauchy distribution at point x by default, standard Cauchy distribution is used (t=0, s=1)

  • t is the location parameter
  • s is the scale parameter

The Cauchy distribution, also called Lorentzian distribution, is a continuous distribution describing resonance behavior The mean and standard deviation of the Cauchy distribution are undefined. The practical meaning of this is that collecting 1,000 data points gives no more accurate an estimate of the mean and standard deviation than does a single point. The formula was taken from "Engineering Statistics Handbook" on site http://www.itl.nist.gov/div898/handbook/eda/section3/eda3663.htm

Example:

TF1* fc = new TF1("fc", "TMath::CauchyDist(x, [0], [1])", -5, 5);
fc->SetParameters(0, 1);
fc->Draw();
1-Dim function class
Definition TF1.h:233
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1335
virtual void SetParameters(const Double_t *params)
Definition TF1.h:670
Author
Anna Kreshuk

Definition at line 2175 of file TMath.cxx.

◆ Ccgs()

constexpr Double_t TMath::Ccgs ( )
constexpr

\( cm s^{-1} \)

Definition at line 121 of file TMath.h.

◆ Ceil()

Double_t TMath::Ceil ( Double_t  x)
inline

Rounds x upward, returning the smallest integral value that is not less than x.

Definition at line 668 of file TMath.h.

◆ CeilNint()

Int_t TMath::CeilNint ( Double_t  x)
inline

Returns the nearest integer of TMath::Ceil(x).

Definition at line 674 of file TMath.h.

◆ ChisquareQuantile()

Double_t TMath::ChisquareQuantile ( Double_t  p,
Double_t  ndf 
)

Evaluate the quantiles of the chi-squared probability distribution function.

Algorithm AS 91 Appl. Statist. (1975) Vol.24, P.35Incorporates the suggested changes in AS R85 (vol.40(1), pp.233-5, 1991)

Parameters
[in]pthe probability value, at which the quantile is computed
[in]ndfnumber of degrees of freedom
Author
Anna Kreshuk

Definition at line 2193 of file TMath.cxx.

◆ Cos()

Double_t TMath::Cos ( Double_t  x)
inline

Returns the cosine of an angle of x radians.

Definition at line 594 of file TMath.h.

◆ CosH()

Double_t TMath::CosH ( Double_t  x)
inline

Returns the hyperbolic cosine of x.

Definition at line 612 of file TMath.h.

◆ Cross()

template<typename T >
T * TMath::Cross ( const T  v1[3],
const T  v2[3],
out[3] 
)

Calculates the Cross Product of two vectors: out = [v1 x v2].

Definition at line 1197 of file TMath.h.

◆ CUncertainty()

constexpr Double_t TMath::CUncertainty ( )
constexpr

Speed of light uncertainty.

Definition at line 128 of file TMath.h.

◆ DegToRad()

constexpr Double_t TMath::DegToRad ( )
constexpr

Conversion from degree to radian: \( \frac{\pi}{180} \).

Definition at line 79 of file TMath.h.

◆ DiLog()

Double_t TMath::DiLog ( Double_t  x)

Modified Struve functions of order 1.

The DiLogarithm function Code translated by R.Brun from CERNLIB DILOG function C332.

Definition at line 116 of file TMath.cxx.

◆ E()

constexpr Double_t TMath::E ( )
constexpr

Base of natural log: \( e \).

Definition at line 93 of file TMath.h.

◆ Erf()

Double_t TMath::Erf ( Double_t  x)

Computation of the error function erf(x).

Erf(x) = (2/sqrt(pi)) Integral(exp(-t^2))dt between 0 and x

Definition at line 190 of file TMath.cxx.

◆ Erfc()

Double_t TMath::Erfc ( Double_t  x)

Computes the complementary error function erfc(x).

Erfc(x) = (2/sqrt(pi)) Integral(exp(-t^2))dt between x and infinity

Definition at line 199 of file TMath.cxx.

◆ ErfcInverse()

Double_t TMath::ErfcInverse ( Double_t  x)

Returns the inverse of the complementary error function.

x must be 0<x<2 implement using the quantile of the normal distribution instead of ErfInverse for better numerical precision for large x

Definition at line 242 of file TMath.cxx.

◆ ErfInverse()

Double_t TMath::ErfInverse ( Double_t  x)

Returns the inverse error function.

x must be <-1<x<1

Definition at line 208 of file TMath.cxx.

◆ EulerGamma()

constexpr Double_t TMath::EulerGamma ( )
constexpr

Euler-Mascheroni Constant.

Definition at line 332 of file TMath.h.

◆ Even()

Bool_t TMath::Even ( Long_t  a)
inline

Returns true if a is even.

Definition at line 113 of file TMathBase.h.

◆ Exp()

Double_t TMath::Exp ( Double_t  x)
inline

Returns the base-e exponential function of x, which is e raised to the power x.

Definition at line 709 of file TMath.h.

◆ Factorial()

Double_t TMath::Factorial ( Int_t  i)

Computes factorial(n).

Definition at line 252 of file TMath.cxx.

◆ FDist()

Double_t TMath::FDist ( Double_t  F,
Double_t  N,
Double_t  M 
)

Computes the density function of F-distribution (probability function, integral of density, is computed in FDistI).

Parameters N and M stand for degrees of freedom of chi-squares mentioned above parameter F is the actual variable x of the density function p(x) and the point at which the density function is calculated.

About F distribution:

F-distribution arises in testing whether two random samples have the same variance. It is the ratio of two chi-square distributions, with N and M degrees of freedom respectively, where each chi-square is first divided by it's number of degrees of freedom.

Author
Anna Kreshuk

Definition at line 2277 of file TMath.cxx.

◆ FDistI()

Double_t TMath::FDistI ( Double_t  F,
Double_t  N,
Double_t  M 
)

Calculates the cumulative distribution function of F-distribution (see ROOT::Math::fdistribution_cdf).

This function occurs in the statistical test of whether two observed samples have the same variance. For this test a certain statistic F, the ratio of observed dispersion of the first sample to that of the second sample, is calculated. N and M stand for numbers of degrees of freedom in the samples 1-FDistI() is the significance level at which the hypothesis "1 has smaller variance than 2" can be rejected. A small numerical value of 1 - FDistI() implies a very significant rejection, in turn implying high confidence in the hypothesis "1 has variance greater than 2".

Author
Anna Kreshuk

Definition at line 2297 of file TMath.cxx.

◆ Finite() [1/2]

Int_t TMath::Finite ( Double_t  x)
inline

Check if it is finite with a mask in order to be consistent in presence of fast math.

Inspired from the CMSSW FWCore/Utilities package

Definition at line 770 of file TMath.h.

◆ Finite() [2/2]

Int_t TMath::Finite ( Float_t  x)
inline

Check if it is finite with a mask in order to be consistent in presence of fast math.

Inspired from the CMSSW FWCore/Utilities package

Definition at line 800 of file TMath.h.

◆ Floor()

Double_t TMath::Floor ( Double_t  x)
inline

Rounds x downward, returning the largest integral value that is not greater than x.

Definition at line 680 of file TMath.h.

◆ FloorNint()

Int_t TMath::FloorNint ( Double_t  x)
inline

Returns the nearest integer of TMath::Floor(x).

Definition at line 686 of file TMath.h.

◆ Freq()

Double_t TMath::Freq ( Double_t  x)

Computation of the normal frequency function freq(x).

Freq(x) = (1/sqrt(2pi)) Integral(exp(-t^2/2))dt between -infinity and x.

Translated from CERNLIB C300 by Rene Brun.

Definition at line 270 of file TMath.cxx.

◆ G()

constexpr Double_t TMath::G ( )
constexpr

Gravitational constant in: \( m^{3} kg^{-1} s^{-2} \).

Definition at line 135 of file TMath.h.

◆ GamCf()

Double_t TMath::GamCf ( Double_t  a,
Double_t  x 
)

Computation of the incomplete gamma function P(a,x) via its continued fraction representation.

Author
NvE 14-nov-1998 UU-SAP Utrecht

Definition at line 380 of file TMath.cxx.

◆ Gamma() [1/2]

Double_t TMath::Gamma ( Double_t  a,
Double_t  x 
)

Computation of the normalized lower incomplete gamma function P(a,x) as defined in the Handbook of Mathematical Functions by Abramowitz and Stegun, formula 6.5.1 on page 260 .

Its normalization is such that TMath::Gamma(a,+infinity) = 1 .

\[ P(a, x) = \frac{1}{\Gamma(a)} \int_{0}^{x} t^{a-1} e^{-t} dt \]

Author
NvE 14-nov-1998 UU-SAP Utrecht

Definition at line 369 of file TMath.cxx.

◆ Gamma() [2/2]

Double_t TMath::Gamma ( Double_t  z)

Computation of gamma(z) for all z.

C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.

Definition at line 353 of file TMath.cxx.

◆ GammaDist()

Double_t TMath::GammaDist ( Double_t  x,
Double_t  gamma,
Double_t  mu = 0,
Double_t  beta = 1 
)

Computes the density function of Gamma distribution at point x.

Parameters
[in]xevaluation point
[in]gammashape parameter
[in]mulocation parameter
[in]betascale parameter

The definition can be found in "Engineering Statistics Handbook" on site http://www.itl.nist.gov/div898/handbook/eda/section3/eda366b.htm use now implementation in ROOT::Math::gamma_pdf

Definition at line 2347 of file TMath.cxx.

◆ GamSer()

Double_t TMath::GamSer ( Double_t  a,
Double_t  x 
)

Computation of the incomplete gamma function P(a,x) via its series representation.

Author
NvE 14-nov-1998 UU-SAP Utrecht

Definition at line 417 of file TMath.cxx.

◆ Gaus()

Double_t TMath::Gaus ( Double_t  x,
Double_t  mean = 0,
Double_t  sigma = 1,
Bool_t  norm = kFALSE 
)

Calculates a gaussian function with mean and sigma.

If norm=kTRUE (default is kFALSE) the result is divided by sqrt(2*Pi)*sigma.

Definition at line 471 of file TMath.cxx.

◆ Gcgs()

constexpr Double_t TMath::Gcgs ( )
constexpr

\( cm^{3} g^{-1} s^{-2} \)

Definition at line 143 of file TMath.h.

◆ GeomMean() [1/2]

template<typename Iterator >
Double_t TMath::GeomMean ( Iterator  first,
Iterator  last 
)

Returns the geometric mean of an array defined by the iterators.

\[ GeomMean = (\prod_{i=0}^{n-1} |a[i]|)^{1/n} \]

Definition at line 1103 of file TMath.h.

◆ GeomMean() [2/2]

template<typename T >
Double_t TMath::GeomMean ( Long64_t  n,
const T *  a 
)

Returns the geometric mean of an array a of size n.

\[ GeomMean = (\prod_{i=0}^{n-1} |a[i]|)^{1/n} \]

Definition at line 1123 of file TMath.h.

◆ GhbarC()

constexpr Double_t TMath::GhbarC ( )
constexpr

\( \frac{G}{\hbar C} \) in \( (GeV/c^{2})^{-2} \)

Definition at line 158 of file TMath.h.

◆ GhbarCUncertainty()

constexpr Double_t TMath::GhbarCUncertainty ( )
constexpr

\( \frac{G}{\hbar C} \) uncertainty.

Definition at line 166 of file TMath.h.

◆ Gn()

constexpr Double_t TMath::Gn ( )
constexpr

Standard acceleration of gravity in \( m s^{-2} \).

Definition at line 174 of file TMath.h.

◆ GnUncertainty()

constexpr Double_t TMath::GnUncertainty ( )
constexpr

Standard acceleration of gravity uncertainty.

Definition at line 181 of file TMath.h.

◆ GUncertainty()

constexpr Double_t TMath::GUncertainty ( )
constexpr

Gravitational constant uncertainty.

Definition at line 150 of file TMath.h.

◆ H()

constexpr Double_t TMath::H ( )
constexpr

Planck's constant in \( J s \): \( h \).

Definition at line 188 of file TMath.h.

◆ Hash() [1/2]

ULong_t TMath::Hash ( const char *  str)

Definition at line 1416 of file TMath.cxx.

◆ Hash() [2/2]

ULong_t TMath::Hash ( const void *  txt,
Int_t  ntxt 
)

Calculates hash index from any char string.

Based on pre-calculated table of 256 specially selected numbers. These numbers are selected in such a way, that for string length == 4 (integer number) the hash is unambiguous, i.e. from hash value we can recalculate input (no degeneration).

The quality of hash method is good enough, that "random" numbers made as R = Hash(1), Hash(2), ...Hash(N) tested by <R>, <R*R>, <Ri*Ri+1> gives the same result as for libc rand().

  • For string: i = TMath::Hash(string,nstring);
  • For int: i = TMath::Hash(&intword,sizeof(int));
  • For pointer: i = TMath::Hash(&pointer,sizeof(void*));
          V.Perev
    
    This function is kept for back compatibility. The code previously in this function has been moved to the static function TString::Hash

Definition at line 1408 of file TMath.cxx.

◆ Hbar()

constexpr Double_t TMath::Hbar ( )
constexpr

\( \hbar \) in \( J s \): \( \hbar = \frac{h}{2\pi} \)

Definition at line 211 of file TMath.h.

◆ Hbarcgs()

constexpr Double_t TMath::Hbarcgs ( )
constexpr

\( erg s \)

Definition at line 218 of file TMath.h.

◆ HbarUncertainty()

constexpr Double_t TMath::HbarUncertainty ( )
constexpr

\( \hbar \) uncertainty.

Definition at line 225 of file TMath.h.

◆ HC()

constexpr Double_t TMath::HC ( )
constexpr

\( hc \) in \( J m \)

Definition at line 233 of file TMath.h.

◆ HCcgs()

constexpr Double_t TMath::HCcgs ( )
constexpr

\( erg cm \)

Definition at line 240 of file TMath.h.

◆ Hcgs()

constexpr Double_t TMath::Hcgs ( )
constexpr

\( erg s \)

Definition at line 195 of file TMath.h.

◆ HUncertainty()

constexpr Double_t TMath::HUncertainty ( )
constexpr

Planck's constant uncertainty.

Definition at line 202 of file TMath.h.

◆ Hypot() [1/2]

Double_t TMath::Hypot ( Double_t  x,
Double_t  y 
)

Returns sqrt(x*x + y*y)

Definition at line 59 of file TMath.cxx.

◆ Hypot() [2/2]

Long_t TMath::Hypot ( Long_t  x,
Long_t  y 
)

Returns sqrt(x*x + y*y)

Definition at line 51 of file TMath.cxx.

◆ Infinity()

Double_t TMath::Infinity ( )
inline

Returns an infinity as defined by the IEEE standard.

Definition at line 917 of file TMath.h.

◆ InvPi()

constexpr Double_t TMath::InvPi ( )
constexpr

\( \frac{1.}{\pi}\)

Definition at line 65 of file TMath.h.

◆ IsInside()

template<typename T >
Bool_t TMath::IsInside ( xp,
yp,
Int_t  np,
T *  x,
T *  y 
)

Function which returns kTRUE if point xp,yp lies inside the polygon defined by the np points in arrays x and y, kFALSE otherwise.

Note that the polygon may be open or closed.

Definition at line 1233 of file TMath.h.

◆ IsNaN() [1/2]

Bool_t TMath::IsNaN ( Double_t  x)
inline

Definition at line 892 of file TMath.h.

◆ IsNaN() [2/2]

Bool_t TMath::IsNaN ( Float_t  x)
inline

Definition at line 893 of file TMath.h.

◆ K()

constexpr Double_t TMath::K ( )
constexpr

Boltzmann's constant in \( J K^{-1} \): \( k \).

Definition at line 247 of file TMath.h.

◆ Kcgs()

constexpr Double_t TMath::Kcgs ( )
constexpr

\( erg K^{-1} \)

Definition at line 254 of file TMath.h.

◆ KolmogorovProb()

Double_t TMath::KolmogorovProb ( Double_t  z)

Calculates the Kolmogorov distribution function,.

\[ P(z) = 2 \sum_{j=1}^{\infty} (-1)^{j-1} e^{-2 j^2 z^2} \]

which gives the probability that Kolmogorov's test statistic will exceed the value z assuming the null hypothesis. This gives a very powerful test for comparing two one-dimensional distributions. see, for example, Eadie et al, "statistical Methods in Experimental Physics', pp 269-270). This function returns the confidence level for the null hypothesis, where: - \_form#647, and - \_form#648 is the maximum deviation between a hypothetical distribution function and an experimental distribution with - \_form#354 events NOTE: To compare two experimental distributions with m and n events, use \_form#649 Accuracy: The function is far too accurate for any imaginable application. Probabilities less than \_form#650 are returned as zero. However, remember that the formula is only valid for "large" n.

Theta function inversion formula is used for z <= 1

This function was translated by Rene Brun from PROBKL in CERNLIB.

Definition at line 679 of file TMath.cxx.

◆ KolmogorovTest()

Double_t TMath::KolmogorovTest ( Int_t  na,
const Double_t a,
Int_t  nb,
const Double_t b,
Option_t option 
)

Statistical test whether two one-dimensional sets of points are compatible with coming from the same parent distribution, using the Kolmogorov test.

That is, it is used to compare two experimental distributions of unbinned data.

Input:

a,b: One-dimensional arrays of length na, nb, respectively. The elements of a and b must be given in ascending order. option is a character string to specify options "D" Put out a line of "Debug" printout "M" Return the Maximum Kolmogorov distance instead of prob

Output:

The returned value prob is a calculated confidence level which gives a statistical test for compatibility of a and b. Values of prob close to zero are taken as indicating a small probability of compatibility. For two point sets drawn randomly from the same parent distribution, the value of prob should be uniformly distributed between zero and one. in case of error the function return -1 If the 2 sets have a different number of points, the minimum of the two sets is used.

Method:

The Kolmogorov test is used. The test statistic is the maximum deviation between the two integrated distribution functions, multiplied by the normalizing factor (rdmax*sqrt(na*nb/(na+nb)).

Code adapted by Rene Brun from CERNLIB routine TKOLMO (Fred James) (W.T. Eadie, D. Drijard, F.E. James, M. Roos and B. Sadoulet, Statistical Methods in Experimental Physics, (North-Holland, Amsterdam 1971) 269-271)

Method Improvement by Jason A Detwiler (JADetwiler@lbl.gov)

The nuts-and-bolts of the TMath::KolmogorovTest() algorithm is a for-loop over the two sorted arrays a and b representing empirical distribution functions. The for-loop handles 3 cases: when the next points to be evaluated satisfy a>b, a<b, or a=b:

for (Int_t i=0;i<na+nb;i++) {
if (a[ia-1] < b[ib-1]) {
rdiff -= sa;
ia++;
if (ia > na) {ok = kTRUE; break;}
} else if (a[ia-1] > b[ib-1]) {
rdiff += sb;
ib++;
if (ib > nb) {ok = kTRUE; break;}
} else {
rdiff += sb - sa;
ia++;
ib++;
if (ia > na) {ok = kTRUE; break;}
if (ib > nb) {ok = kTRUE; break;}
}
rdmax = TMath::Max(rdmax,TMath::Abs(rdiff));
}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123

For the last case, a=b, the algorithm advances each array by one index in an attempt to move through the equality. However, this is incorrect when one or the other of a or b (or both) have a repeated value, call it x. For the KS statistic to be computed properly, rdiff needs to be calculated after all of the a and b at x have been tallied (this is due to the definition of the empirical distribution function; another way to convince yourself that the old CERNLIB method is wrong is that it implies that the function defined as the difference between a and b is multi-valued at x – besides being ugly, this would invalidate Kolmogorov's theorem).

The solution is to just add while-loops into the equality-case handling to perform the tally:

} else {
double x = a[ia-1];
while(a[ia-1] == x && ia <= na) {
rdiff -= sa;
ia++;
}
while(b[ib-1] == x && ib <= nb) {
rdiff += sb;
ib++;
}
if (ia > na) {ok = kTRUE; break;}
if (ib > nb) {ok = kTRUE; break;}
}
Double_t x[n]
Definition legend1.C:17

Note:

A good description of the Kolmogorov test can be seen at: http://www.itl.nist.gov/div898/handbook/eda/section3/eda35g.htm

Definition at line 805 of file TMath.cxx.

◆ KOrdStat()

template<class Element , typename Size >
Element TMath::KOrdStat ( Size  n,
const Element *  a,
Size  k,
Size *  work = 0 
)

Returns k_th order statistic of the array a of size n (k_th smallest element out of n elements).

C-convention is used for array indexing, so if you want the second smallest element, call KOrdStat(n, a, 1).

If work is supplied, it is used to store the sorting index and assumed to be >= n. If work=0, local storage is used, either on the stack if n < kWorkMax or on the heap for n >= kWorkMax. Note that the work index array will not contain the sorted indices but all indices of the smaller element in arbitrary order in work[0,...,k-1] and all indices of the larger element in arbitrary order in work[k+1,..,n-1] work[k] will contain instead the index of the returned element.

Taken from "Numerical Recipes in C++" without the index array implemented by Anna Khreshuk.

See also the declarations at the top of this file

Definition at line 1359 of file TMath.h.

◆ KUncertainty()

constexpr Double_t TMath::KUncertainty ( )
constexpr

Boltzmann's constant uncertainty.

Definition at line 261 of file TMath.h.

◆ Landau()

Double_t TMath::Landau ( Double_t  x,
Double_t  mu = 0,
Double_t  sigma = 1,
Bool_t  norm = kFALSE 
)

The LANDAU function.

mu is a location parameter and correspond approximately to the most probable value and sigma is a scale parameter (not the sigma of the full distribution which is not defined) Note that for mu=0 and sigma=1 (default values) the exact location of the maximum of the distribution (most proper value) is at x = -0.22278 This function has been adapted from the CERNLIB routine G110 denlan. If norm=kTRUE (default is kFALSE) the result is divided by sigma

Definition at line 492 of file TMath.cxx.

◆ LandauI()

Double_t TMath::LandauI ( Double_t  x)

Returns the cumulative (lower tail integral) of the Landau distribution function at point x.

(see ROOT::Math::landau_cdf) The algorithm was taken from the Cernlib function dislan(G110) Reference: K.S.Kolbig and B.Schorr, "A program package for the Landau distribution", Computer Phys.Comm., 31(1984), 97-111

Definition at line 2845 of file TMath.cxx.

◆ LaplaceDist()

Double_t TMath::LaplaceDist ( Double_t  x,
Double_t  alpha = 0,
Double_t  beta = 1 
)

Computes the probability density function of Laplace distribution at point x, with location parameter alpha and shape parameter beta.

By default, alpha=0, beta=1 This distribution is known under different names, most common is double exponential distribution, but it also appears as the two-tailed exponential or the bilateral exponential distribution

Definition at line 2364 of file TMath.cxx.

◆ LaplaceDistI()

Double_t TMath::LaplaceDistI ( Double_t  x,
Double_t  alpha = 0,
Double_t  beta = 1 
)

Computes the cumulative distribution function (lower tail integral) of Laplace distribution at point x, with location parameter alpha and shape parameter beta.

By default, alpha=0, beta=1 This distribution is known under different names, most common is double exponential distribution, but it also appears as the two-tailed exponential or the bilateral exponential distribution

Definition at line 2380 of file TMath.cxx.

◆ Ldexp()

Double_t TMath::Ldexp ( Double_t  x,
Int_t  exp 
)
inline

Returns the result of multiplying x (the significant) by 2 raised to the power of exp (the exponent).

Definition at line 715 of file TMath.h.

◆ Ln10()

constexpr Double_t TMath::Ln10 ( )
constexpr

Natural log of 10 (to convert log to ln)

Definition at line 100 of file TMath.h.

◆ LnGamma()

Double_t TMath::LnGamma ( Double_t  z)

Computation of ln[gamma(z)] for all z.

C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.

The accuracy of the result is better than 2e-10.

Author
NvE 14-nov-1998 UU-SAP Utrecht

Definition at line 509 of file TMath.cxx.

◆ LocMax() [1/2]

template<typename Iterator >
Iterator TMath::LocMax ( Iterator  first,
Iterator  last 
)

Returns index of array with the maximum element.

If more than one element is maximum returns first found.

Definition at line 1030 of file TMath.h.

◆ LocMax() [2/2]

template<typename T >
Long64_t TMath::LocMax ( Long64_t  n,
const T *  a 
)

Returns index of array with the maximum element.

If more than one element is maximum returns first found.

Implement here since it is faster (see comment in LocMin function)

Definition at line 1012 of file TMath.h.

◆ LocMin() [1/2]

template<typename Iterator >
Iterator TMath::LocMin ( Iterator  first,
Iterator  last 
)

Returns index of array with the minimum element.

If more than one element is minimum returns first found.

Definition at line 1000 of file TMath.h.

◆ LocMin() [2/2]

template<typename T >
Long64_t TMath::LocMin ( Long64_t  n,
const T *  a 
)

Returns index of array with the minimum element.

If more than one element is minimum returns first found.

Implement here since this one is found to be faster (mainly on 64 bit machines) than stl generic implementation. When performing the comparison, the STL implementation needs to de-reference both the array iterator and the iterator pointing to the resulting minimum location

Definition at line 982 of file TMath.h.

◆ Log()

Double_t TMath::Log ( Double_t  x)
inline

Returns the natural logarithm of x.

Definition at line 756 of file TMath.h.

◆ Log10()

Double_t TMath::Log10 ( Double_t  x)
inline

Returns the common (base-10) logarithm of x.

Definition at line 762 of file TMath.h.

◆ Log2()

Double_t TMath::Log2 ( Double_t  x)

Returns the binary (base-2) logarithm of x.

Definition at line 107 of file TMath.cxx.

◆ LogE()

constexpr Double_t TMath::LogE ( )
constexpr

Base-10 log of e (to convert ln to log)

Definition at line 107 of file TMath.h.

◆ LogNormal()

Double_t TMath::LogNormal ( Double_t  x,
Double_t  sigma,
Double_t  theta = 0,
Double_t  m = 1 
)

Computes the density of LogNormal distribution at point x.

Variable X has lognormal distribution if Y=Ln(X) has normal distribution

Parameters
[in]xis the evaluation point
[in]sigmais the shape parameter
[in]thetais the location parameter
[in]mis the scale parameter

The formula was taken from "Engineering Statistics Handbook" on site http://www.itl.nist.gov/div898/handbook/eda/section3/eda3669.htm Implementation using ROOT::Math::lognormal_pdf

Definition at line 2437 of file TMath.cxx.

◆ Max() [1/10]

Double_t TMath::Max ( Double_t  a,
Double_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Double_t.

Definition at line 295 of file TMathBase.h.

◆ Max() [2/10]

Float_t TMath::Max ( Float_t  a,
Float_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Float_t.

Definition at line 290 of file TMathBase.h.

◆ Max() [3/10]

Int_t TMath::Max ( Int_t  a,
Int_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Int_t.

Definition at line 260 of file TMathBase.h.

◆ Max() [4/10]

Long64_t TMath::Max ( Long64_t  a,
Long64_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Long64_t.

Definition at line 280 of file TMathBase.h.

◆ Max() [5/10]

Long_t TMath::Max ( Long_t  a,
Long_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Long_t.

Definition at line 270 of file TMathBase.h.

◆ Max() [6/10]

Short_t TMath::Max ( Short_t  a,
Short_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are Short_t.

Definition at line 250 of file TMathBase.h.

◆ Max() [7/10]

UInt_t TMath::Max ( UInt_t  a,
UInt_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are UInt_t.

Definition at line 265 of file TMathBase.h.

◆ Max() [8/10]

ULong64_t TMath::Max ( ULong64_t  a,
ULong64_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are ULong64_t.

Definition at line 285 of file TMathBase.h.

◆ Max() [9/10]

ULong_t TMath::Max ( ULong_t  a,
ULong_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are ULong_t.

Definition at line 275 of file TMathBase.h.

◆ Max() [10/10]

UShort_t TMath::Max ( UShort_t  a,
UShort_t  b 
)
inline

Returns the largest of a and b.

If both are equivalent, a is returned. aand b are UShort_t.

Definition at line 255 of file TMathBase.h.

◆ MaxElement()

template<typename T >
T TMath::MaxElement ( Long64_t  n,
const T *  a 
)

Returns maximum of array a of length n.

Definition at line 968 of file TMath.h.

◆ Mean() [1/3]

template<typename Iterator >
Double_t TMath::Mean ( Iterator  first,
Iterator  last 
)

Returns the weighted mean of an array defined by the iterators.

Definition at line 1040 of file TMath.h.

◆ Mean() [2/3]

template<typename Iterator , typename WeightIterator >
Double_t TMath::Mean ( Iterator  first,
Iterator  last,
WeightIterator  w 
)

Returns the weighted mean of an array defined by the first and last iterators.

The w iterator should point to the first element of a vector of weights of the same size as the main array.

Definition at line 1060 of file TMath.h.

◆ Mean() [3/3]

template<typename T >
Double_t TMath::Mean ( Long64_t  n,
const T *  a,
const Double_t w = nullptr 
)

Returns the weighted mean of an array a with length n.

Definition at line 1089 of file TMath.h.

◆ Median()

template<typename T >
Double_t TMath::Median ( Long64_t  n,
const T *  a,
const Double_t w = nullptr,
Long64_t work = nullptr 
)

Same as RMS.

Returns the median of the array a where each entry i has weight w[i] .

Both arrays have a length of at least n . The median is a number obtained from the sorted array a through

median = (a[jl]+a[jh])/2. where (using also the sorted index on the array w)

sum_i=0,jl w[i] <= sumTot/2 sum_i=0,jh w[i] >= sumTot/2 sumTot = sum_i=0,n w[i]

If w=0, the algorithm defaults to the median definition where it is a number that divides the sorted sequence into 2 halves. When n is odd or n > 1000, the median is kth element k = (n + 1) / 2. when n is even and n < 1000the median is a mean of the elements k = n/2 and k = n/2 + 1.

If the weights are supplied (w not 0) all weights must be >= 0

If work is supplied, it is used to store the sorting index and assumed to be >= n . If work=0, local storage is used, either on the stack if n < kWorkMax or on the heap for n >= kWorkMax .

Definition at line 1272 of file TMath.h.

◆ Min() [1/10]

Double_t TMath::Min ( Double_t  a,
Double_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Double_t.

Definition at line 243 of file TMathBase.h.

◆ Min() [2/10]

Float_t TMath::Min ( Float_t  a,
Float_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Float_t.

Definition at line 238 of file TMathBase.h.

◆ Min() [3/10]

Int_t TMath::Min ( Int_t  a,
Int_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Int_t.

Definition at line 208 of file TMathBase.h.

◆ Min() [4/10]

Long64_t TMath::Min ( Long64_t  a,
Long64_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Long64_t.

Definition at line 228 of file TMathBase.h.

◆ Min() [5/10]

Long_t TMath::Min ( Long_t  a,
Long_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Long_t.

Definition at line 218 of file TMathBase.h.

◆ Min() [6/10]

Short_t TMath::Min ( Short_t  a,
Short_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Short_t.

Definition at line 198 of file TMathBase.h.

◆ Min() [7/10]

UInt_t TMath::Min ( UInt_t  a,
UInt_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are Short_t.

Definition at line 213 of file TMathBase.h.

◆ Min() [8/10]

ULong64_t TMath::Min ( ULong64_t  a,
ULong64_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are ULong64_t.

Definition at line 233 of file TMathBase.h.

◆ Min() [9/10]

ULong_t TMath::Min ( ULong_t  a,
ULong_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are ULong_t.

Definition at line 223 of file TMathBase.h.

◆ Min() [10/10]

UShort_t TMath::Min ( UShort_t  a,
UShort_t  b 
)
inline

Returns the smallest of a and b.

If both are equivalent, a is returned. aand b are UShort_t.

Definition at line 203 of file TMathBase.h.

◆ MinElement()

template<typename T >
T TMath::MinElement ( Long64_t  n,
const T *  a 
)

Returns minimum of array a of length n.

Definition at line 960 of file TMath.h.

◆ MWair()

constexpr Double_t TMath::MWair ( )
constexpr

Molecular weight of dry air 1976 US Standard Atmosphere in \( kg kmol^{-1} \) or \( gm mol^{-1} \)

Definition at line 317 of file TMath.h.

◆ Na()

constexpr Double_t TMath::Na ( )
constexpr

Avogadro constant (Avogadro's Number) in \( mol^{-1} \).

Definition at line 284 of file TMath.h.

◆ NaUncertainty()

constexpr Double_t TMath::NaUncertainty ( )
constexpr

Avogadro constant (Avogadro's Number) uncertainty.

Definition at line 291 of file TMath.h.

◆ NextPrime()

Long_t TMath::NextPrime ( Long_t  x)

◆ Nint()

template<typename T >
Int_t TMath::Nint ( x)
inline

Round to nearest integer. Rounds half integers to the nearest even integer.

Definition at line 693 of file TMath.h.

◆ Normal2Plane()

template<typename T >
T * TMath::Normal2Plane ( const T  p1[3],
const T  p2[3],
const T  p3[3],
normal[3] 
)

Calculates a normal vector of a plane.

Parameters
[in]p1,p2,p33 3D points belonged the plane to define it.
[out]normalPointer to 3D normal vector (normalized)

Definition at line 1212 of file TMath.h.

◆ Normalize() [1/2]

Double_t TMath::Normalize ( Double_t  v[3])

Normalize a vector v in place.

Returns the norm of the original vector. This implementation (thanks Kevin Lynch krlyn.nosp@m.ch@b.nosp@m.u.edu) is protected against possible overflows.

Definition at line 535 of file TMath.cxx.

◆ Normalize() [2/2]

Float_t TMath::Normalize ( Float_t  v[3])

Normalize a vector v in place.

Returns the norm of the original vector.

Definition at line 518 of file TMath.cxx.

◆ NormCross()

template<typename T >
T TMath::NormCross ( const T  v1[3],
const T  v2[3],
out[3] 
)
inline

Calculates the Normalized Cross Product of two vectors.

Definition at line 951 of file TMath.h.

◆ NormQuantile()

Double_t TMath::NormQuantile ( Double_t  p)

Computes quantiles for standard normal distribution N(0, 1) at probability p.

ALGORITHM AS241 APPL. STATIST. (1988) VOL. 37, NO. 3, 477-484.

Definition at line 2456 of file TMath.cxx.

◆ Odd()

Bool_t TMath::Odd ( Long_t  a)
inline

Returns true if a is odd.

Definition at line 117 of file TMathBase.h.

◆ Permute()

Bool_t TMath::Permute ( Int_t  n,
Int_t a 
)

Simple recursive algorithm to find the permutations of n natural numbers, not necessarily all distinct adapted from CERNLIB routine PERMU.

The input array has to be initialised with a non descending sequence. The method returns kFALSE when all combinations are exhausted.

Definition at line 2557 of file TMath.cxx.

◆ Pi()

constexpr Double_t TMath::Pi ( )
constexpr

\( \pi\)

Definition at line 37 of file TMath.h.

◆ PiOver2()

constexpr Double_t TMath::PiOver2 ( )
constexpr

\( \frac{\pi}{2} \)

Definition at line 51 of file TMath.h.

◆ PiOver4()

constexpr Double_t TMath::PiOver4 ( )
constexpr

\( \frac{\pi}{4} \)

Definition at line 58 of file TMath.h.

◆ Poisson()

Double_t TMath::Poisson ( Double_t  x,
Double_t  par 
)

Computes the Poisson distribution function for (x,par).

The Poisson PDF is implemented by means of Euler's Gamma-function (for the factorial), so for any x integer argument it is the correct Poisson distribution. BUT for non-integer x values, it IS NOT equal to the Poisson distribution !

Definition at line 587 of file TMath.cxx.

◆ PoissonI()

Double_t TMath::PoissonI ( Double_t  x,
Double_t  par 
)

Computes the Discrete Poisson distribution function for (x,par).

This is a discrete and a non-smooth function. This function is equivalent to ROOT::Math::poisson_pdf

Definition at line 615 of file TMath.cxx.

◆ Power() [1/5]

Double_t TMath::Power ( Double_t  x,
Double_t  y 
)
inline

Returns x raised to the power y.

Definition at line 739 of file TMath.h.

◆ Power() [2/5]

Double_t TMath::Power ( Double_t  x,
Int_t  y 
)
inline

Returns x raised to the power y.

Definition at line 745 of file TMath.h.

◆ Power() [3/5]

LongDouble_t TMath::Power ( Long64_t  x,
Long64_t  y 
)
inline

Returns x raised to the power y.

Definition at line 733 of file TMath.h.

◆ Power() [4/5]

LongDouble_t TMath::Power ( LongDouble_t  x,
Long64_t  y 
)
inline

Returns x raised to the power y.

Definition at line 727 of file TMath.h.

◆ Power() [5/5]

LongDouble_t TMath::Power ( LongDouble_t  x,
LongDouble_t  y 
)
inline

Returns x raised to the power y.

Definition at line 721 of file TMath.h.

◆ Prob()

Double_t TMath::Prob ( Double_t  chi2,
Int_t  ndf 
)

Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf).

Calculations are based on the incomplete gamma function P(a,x), where a=ndf/2 and x=chi2/2.

P(a,x) represents the probability that the observed Chi-squared for a correct model should be less than the value chi2.

The returned probability corresponds to 1-P(a,x), which denotes the probability that an observed Chi-squared exceeds the value chi2 by chance, even for a correct model.

Author
NvE 14-nov-1998 UU-SAP Utrecht

Definition at line 637 of file TMath.cxx.

◆ Qe()

constexpr Double_t TMath::Qe ( )
constexpr

Elementary charge in \( C \) .

Definition at line 339 of file TMath.h.

◆ QeUncertainty()

constexpr Double_t TMath::QeUncertainty ( )
constexpr

Elementary charge uncertainty.

Definition at line 346 of file TMath.h.

◆ Quantiles()

void TMath::Quantiles ( Int_t  n,
Int_t  nprob,
Double_t x,
Double_t quantiles,
Double_t prob,
Bool_t  isSorted = kTRUE,
Int_t index = nullptr,
Int_t  type = 7 
)

Computes sample quantiles, corresponding to the given probabilities.

Parameters
[in]xthe data sample
[in]nits size
[out]quantilescomputed quantiles are returned in there
[in]probprobabilities where to compute quantiles
[in]nprobsize of prob array
[in]isSortedis the input array x sorted ?
[in]indexparameter index
[in]typemethod to compute (from 1 to 9).

NOTE:

When the input is not sorted, an array of integers of size n needs to be allocated. It can be passed by the user in parameter index, or, if not passed, it will be allocated inside the function

Following types are provided:

  • Discontinuous:
    • type=1 - inverse of the empirical distribution function
    • type=2 - like type 1, but with averaging at discontinuities
    • type=3 - SAS definition: nearest even order statistic
  • Piecewise linear continuous:

    • In this case, sample quantiles can be obtained by linear interpolation between the k-th order statistic and p(k). -type=4 - linear interpolation of empirical cdf, p(k)=k/n;
    • type=5 - a very popular definition, p(k) = (k-0.5)/n;
    • type=6 - used by Minitab and SPSS, p(k) = k/(n+1);
    • type=7 - used by S-Plus and R, p(k) = (k-1)/(n-1);
    • type=8 - resulting sample quantiles are approximately median unbiased regardless of the distribution of x. p(k) = (k-1/3)/(n+1/3);
    • type=9 - resulting sample quantiles are approximately unbiased, when the sample comes from Normal distribution. p(k)=(k-3/8)/(n+1/4);

    default type = 7

References:

  1. Hyndman, R.J and Fan, Y, (1996) "Sample quantiles in statistical packages" American Statistician, 50, 361-365
  2. R Project documentation for the function quantile of package {stats}

Definition at line 1207 of file TMath.cxx.

◆ QuietNaN()

Double_t TMath::QuietNaN ( )
inline

Returns a quiet NaN as defined by IEEE 754.

Definition at line 902 of file TMath.h.

◆ R()

constexpr Double_t TMath::R ( )
constexpr

Universal gas constant ( \( Na K \)) in \( J K^{-1} mol^{-1} \)

Definition at line 302 of file TMath.h.

◆ RadToDeg()

constexpr Double_t TMath::RadToDeg ( )
constexpr

Conversion from radian to degree: \( \frac{180}{\pi} \).

Definition at line 72 of file TMath.h.

◆ Range() [1/5]

Double_t TMath::Range ( Double_t  lb,
Double_t  ub,
Double_t  x 
)
inline

Returns x if lb < x < up, lb if x < lb and ub if x > ub.

lb, ub and x are Double_t.

Definition at line 322 of file TMathBase.h.

◆ Range() [2/5]

Int_t TMath::Range ( Int_t  lb,
Int_t  ub,
Int_t  x 
)
inline

Returns x if lb < x < up, lb if x < lb and ub if x > ub.

lb, ub and x are Int_t.

Definition at line 307 of file TMathBase.h.

◆ Range() [3/5]

Long_t TMath::Range ( Long_t  lb,
Long_t  ub,
Long_t  x 
)
inline

Returns x if lb < x < up, lb if x < lb and ub if x > ub.

lb, ub and x are Long_t.

Definition at line 312 of file TMathBase.h.

◆ Range() [4/5]

Short_t TMath::Range ( Short_t  lb,
Short_t  ub,
Short_t  x 
)
inline

Returns x if lb < x < up, lb if x < lb and ub if x > ub.

lb, ub and x are Short_t.

Definition at line 302 of file TMathBase.h.

◆ Range() [5/5]

ULong_t TMath::Range ( ULong_t  lb,
ULong_t  ub,
ULong_t  x 
)
inline

Returns x if lb < x < up, lb if x < lb and ub if x > ub.

lb, ub and x are ULong_t.

Definition at line 317 of file TMathBase.h.

◆ Rgair()

constexpr Double_t TMath::Rgair ( )
constexpr

Dry Air Gas Constant (R / MWair) in \( J kg^{-1} K^{-1} \)

Definition at line 325 of file TMath.h.

◆ RMS() [1/3]

template<typename Iterator >
Double_t TMath::RMS ( Iterator  first,
Iterator  last 
)

Returns the Standard Deviation of an array defined by the iterators.

Note that this function returns the sigma(standard deviation) and not the root mean square of the array.

Use the two pass algorithm, which is slower (! a factor of 2) but much more precise. Since we have a vector the 2 pass algorithm is still faster than the Welford algorithm. (See also ROOT-5545)

Definition at line 1138 of file TMath.h.

◆ RMS() [2/3]

template<typename Iterator , typename WeightIterator >
Double_t TMath::RMS ( Iterator  first,
Iterator  last,
WeightIterator  w 
)

Returns the weighted Standard Deviation of an array defined by the iterators.

Note that this function returns the sigma(standard deviation) and not the root mean square of the array.

As in the unweighted case use the two pass algorithm

Definition at line 1163 of file TMath.h.

◆ RMS() [3/3]

template<typename T >
Double_t TMath::RMS ( Long64_t  n,
const T *  a,
const Double_t w = nullptr 
)

Returns the Standard Deviation of an array a with length n.

Note that this function returns the sigma(standard deviation) and not the root mean square of the array.

Definition at line 1188 of file TMath.h.

◆ RootsCubic()

Bool_t TMath::RootsCubic ( const Double_t  coef[4],
Double_t a,
Double_t b,
Double_t c 
)

Calculates roots of polynomial of 3rd order a*x^3 + b*x^2 + c*x + d, where.

  • a == coef[3],
  • b == coef[2],
  • c == coef[1],
  • d == coef[0]

coef[3] must be different from 0

If the boolean returned by the method is false: ==> there are 3 real roots a,b,c

If the boolean returned by the method is true: ==> there is one real root a and 2 complex conjugates roots (b+i*c,b-i*c)

Author
Francois-Xavier Gentit

Definition at line 1107 of file TMath.cxx.

◆ RUncertainty()

constexpr Double_t TMath::RUncertainty ( )
constexpr

Universal gas constant uncertainty.

Definition at line 309 of file TMath.h.

◆ Sigma()

constexpr Double_t TMath::Sigma ( )
constexpr

Stefan-Boltzmann constant in \( W m^{-2} K^{-4}\): \( \sigma \).

Definition at line 270 of file TMath.h.

◆ SigmaUncertainty()

constexpr Double_t TMath::SigmaUncertainty ( )
constexpr

Stefan-Boltzmann constant uncertainty.

Definition at line 277 of file TMath.h.

◆ Sign() [1/4]

Double_t TMath::Sign ( Double_t  a,
Double_t  b 
)
inline

Returns a value with the magnitude of a and the sign of b.

aand b are Double_t.

Definition at line 185 of file TMathBase.h.

◆ Sign() [2/4]

Float_t TMath::Sign ( Float_t  a,
Float_t  b 
)
inline

Returns a value with the magnitude of a and the sign of b.

aand b are Short_t.

Definition at line 180 of file TMathBase.h.

◆ Sign() [3/4]

LongDouble_t TMath::Sign ( LongDouble_t  a,
LongDouble_t  b 
)
inline

Returns a value with the magnitude of a and the sign of b.

aand b are LongDouble_t.

Definition at line 190 of file TMathBase.h.

◆ Sign() [4/4]

template<typename T1 , typename T2 >
T1 TMath::Sign ( T1  a,
T2  b 
)
inline

Returns a value with the magnitude of a and the sign of b.

Definition at line 175 of file TMathBase.h.

◆ SignalingNaN()

Double_t TMath::SignalingNaN ( )
inline

Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).

Definition at line 910 of file TMath.h.

◆ SignBit() [1/4]

Bool_t TMath::SignBit ( Double_t  a)
inline

Returns whether the sign of Double_t a is negative.

Definition at line 163 of file TMathBase.h.

◆ SignBit() [2/4]

Bool_t TMath::SignBit ( Float_t  a)
inline

Returns whether the sign of Float_t a is negative.

Definition at line 159 of file TMathBase.h.

◆ SignBit() [3/4]

template<typename Integer >
Bool_t TMath::SignBit ( Integer  a)
inline

Returns whether the sign of Integer a is negative.

Definition at line 155 of file TMathBase.h.

◆ SignBit() [4/4]

Bool_t TMath::SignBit ( LongDouble_t  a)
inline

Returns whether the sign of LongDouble_t a is negative.

Definition at line 167 of file TMathBase.h.

◆ Sin()

Double_t TMath::Sin ( Double_t  x)
inline

Returns the sine of an angle of x radians.

Definition at line 588 of file TMath.h.

◆ SinH()

Double_t TMath::SinH ( Double_t  x)
inline

Returns the hyperbolic sine of `x.

Definition at line 606 of file TMath.h.

◆ Sort()

template<typename Element , typename Index >
void TMath::Sort ( Index  n,
const Element *  a,
Index *  index,
Bool_t  down = kTRUE 
)

Sort the n elements of the array a of generic templated type Element.

In output the array index of type Index contains the indices of the sorted array. If down is false sort in increasing order (default is decreasing order).

NOTE that the array index must be created with a length >= n before calling this function. NOTE also that the size type for n must be the same type used for the index array (templated type Index)

Definition at line 431 of file TMathBase.h.

◆ SortItr()

template<typename Iterator , typename IndexIterator >
void TMath::SortItr ( Iterator  first,
Iterator  last,
IndexIterator  index,
Bool_t  down = kTRUE 
)

Sort the n1 elements of the Short_t array defined by its iterators.

In output the array index contains the indices of the sorted array. If down is false sort in increasing order (default is decreasing order).

NOTE that the array index must be created with a length bigger or equal than the main array before calling this function.

Definition at line 406 of file TMathBase.h.

◆ Sq()

Double_t TMath::Sq ( Double_t  x)
inline

Returns x*x.

Definition at line 656 of file TMath.h.

◆ Sqrt()

Double_t TMath::Sqrt ( Double_t  x)
inline

Returns the square root of x.

Definition at line 662 of file TMath.h.

◆ Sqrt2()

constexpr Double_t TMath::Sqrt2 ( )
constexpr

\( \sqrt{2} \)

Definition at line 86 of file TMath.h.

◆ StdDev() [1/3]

template<typename Iterator >
Double_t TMath::StdDev ( Iterator  first,
Iterator  last 
)

Same as RMS.

Definition at line 528 of file TMath.h.

◆ StdDev() [2/3]

template<typename Iterator , typename WeightIterator >
Double_t TMath::StdDev ( Iterator  first,
Iterator  last,
WeightIterator  wfirst 
)

Same as RMS.

Definition at line 529 of file TMath.h.

◆ StdDev() [3/3]

template<typename T >
Double_t TMath::StdDev ( Long64_t  n,
const T *  a,
const Double_t w = nullptr 
)

Definition at line 527 of file TMath.h.

◆ StruveH0()

Double_t TMath::StruveH0 ( Double_t  x)

Bessel function Y1(x) for positive x.

Struve Functions of Order 0.

Converted from CERNLIB M342 by Rene Brun.

Definition at line 1777 of file TMath.cxx.

◆ StruveH1()

Double_t TMath::StruveH1 ( Double_t  x)

Struve functions of order 0.

Struve Functions of Order 1.

Converted from CERNLIB M342 by Rene Brun.

Definition at line 1846 of file TMath.cxx.

◆ StruveL0()

Double_t TMath::StruveL0 ( Double_t  x)

Struve functions of order 1.

Modified Struve Function of Order 0.

Author
Kirill Filimonov.

Definition at line 1923 of file TMath.cxx.

◆ StruveL1()

Double_t TMath::StruveL1 ( Double_t  x)

Modified Struve functions of order 0.

Modified Struve Function of Order 1.

Author
Kirill Filimonov.

Definition at line 1970 of file TMath.cxx.

◆ Student()

Double_t TMath::Student ( Double_t  T,
Double_t  ndf 
)

Computes density function for Student's t- distribution (the probability function (integral of density) is computed in StudentI).

First parameter stands for x - the actual variable of the density function p(x) and the point at which the density is calculated. Second parameter stands for number of degrees of freedom.

About Student distribution: Student's t-distribution is used for many significance tests, for example, for the Student's t-tests for the statistical significance of difference between two sample means and for confidence intervals for the difference between two population means.

Example: suppose we have a random sample of size n drawn from normal distribution with mean Mu and st.deviation Sigma. Then the variable

t = (sample_mean - Mu)/(sample_deviation / sqrt(n))

has Student's t-distribution with n-1 degrees of freedom.

NOTE that this function's second argument is number of degrees of freedom, not the sample size.

As the number of degrees of freedom grows, t-distribution approaches Normal(0,1) distribution.

Author
Anna Kreshuk

Definition at line 2623 of file TMath.cxx.

◆ StudentI()

Double_t TMath::StudentI ( Double_t  T,
Double_t  ndf 
)

Calculates the cumulative distribution function of Student's t-distribution second parameter stands for number of degrees of freedom, not for the number of samples if x has Student's t-distribution, the function returns the probability of x being less than T.

This is equivalent to ROOT::Math::tdistribution_cdf(T,ndf)

Author
Anna Kreshuk

Definition at line 2646 of file TMath.cxx.

◆ StudentQuantile()

Double_t TMath::StudentQuantile ( Double_t  p,
Double_t  ndf,
Bool_t  lower_tail = kTRUE 
)

Computes quantiles of the Student's t-distribution 1st argument is the probability, at which the quantile is computed 2nd argument - the number of degrees of freedom of the Student distribution When the 3rd argument lower_tail is kTRUE (default)- the algorithm returns such x0, that.

P(x < x0)=p

upper tail (lower_tail is kFALSE)- the algorithm returns such x0, that

P(x > x0)=p

the algorithm was taken from: G.W.Hill, "Algorithm 396, Student's t-quantiles" "Communications of the ACM", 13(10), October 1970

Definition at line 2674 of file TMath.cxx.

◆ Tan()

Double_t TMath::Tan ( Double_t  x)
inline

Returns the tangent of an angle of x radians.

Definition at line 600 of file TMath.h.

◆ TanH()

Double_t TMath::TanH ( Double_t  x)
inline

Returns the hyperbolic tangent of x.

Definition at line 618 of file TMath.h.

◆ TwoPi()

constexpr Double_t TMath::TwoPi ( )
constexpr

\( 2\pi\)

Definition at line 44 of file TMath.h.

◆ Vavilov()

Double_t TMath::Vavilov ( Double_t  x,
Double_t  kappa,
Double_t  beta2 
)

Returns the value of the Vavilov probability density function.

Parameters
[in]xthe point were the density function is evaluated
[in]kappavalue of kappa (distribution parameter)
[in]beta2value of beta2 (distribution parameter)

The algorithm was taken from the CernLib function vavden(G115) Reference: A.Rotondi and P.Montagna, Fast Calculation of Vavilov distribution Nucl.Instr. and Meth. B47(1990), 215-224

Accuracy: quote from the reference above:

"The results of our code have been compared with the values of the Vavilov density function computed numerically in an accurate way: our approximation shows a difference of less than 3% around the peak of the density function, slowly increasing going towards the extreme tails to the right and to the left"

For a more accurate implementation see the documentation in the Vavilov class and ROOT::Math::vavilov_accurate_pdf

Definition at line 2778 of file TMath.cxx.

◆ VavilovDenEval()

Double_t TMath::VavilovDenEval ( Double_t  rlam,
Double_t AC,
Double_t HC,
Int_t  itype 
)

Internal function, called by Vavilov and VavilovSet.

Definition at line 3155 of file TMath.cxx.

◆ VavilovI()

Double_t TMath::VavilovI ( Double_t  x,
Double_t  kappa,
Double_t  beta2 
)

Returns the value of the Vavilov cumulative distribution function (lower tail integral of the probability distribution function)

Parameters
[in]xthe point were the density function is evaluated
[in]kappavalue of kappa (distribution parameter)
[in]beta2value of beta2 (distribution parameter)

The algorithm was taken from the CernLib function vavden(G115)

Reference: A.Rotondi and P.Montagna, Fast Calculation of Vavilov distribution Nucl.Instr. and Meth. B47(1990), 215-224

Accuracy: quote from the reference above:

"The results of our code have been compared with the values of the Vavilov density function computed numerically in an accurate way: our approximation shows a difference of less than 3% around the peak of the density function, slowly increasing going towards the extreme tails to the right and to the left"

For a more accurate implementation see the documentation of the Vavilov class and the cumulative ROOT::Math::vavilov_accurate_cdf

Definition at line 2815 of file TMath.cxx.

◆ VavilovSet()

void TMath::VavilovSet ( Double_t  rkappa,
Double_t  beta2,
Bool_t  mode,
Double_t WCM,
Double_t AC,
Double_t HC,
Int_t itype,
Int_t npt 
)

Internal function, called by Vavilov and VavilovI.

Definition at line 2854 of file TMath.cxx.

◆ Voigt()

Double_t TMath::Voigt ( Double_t  xx,
Double_t  sigma,
Double_t  lg,
Int_t  r = 4 
)

Computation of Voigt function (normalised).

Voigt is a convolution of the two functions:

\[ gauss(xx) = \frac{1}{(\sqrt{2\pi} sigma)} e^{\frac{xx^{2}}{(2 sigma{^2})}} \]

and

\[ lorentz(xx) = \frac{ \frac{1}{\pi} \frac{lg}{2} }{ (xx^{2} + \frac{lg^{2}}{4}) } \]

.

The Voigt function is known to be the real part of Faddeeva function also called complex error function [2].

The algorithm was developed by J. Humlicek [1]. This code is based on fortran code presented by R. J. Wells [2]. Translated and adapted by Miha D. Puc

To calculate the Faddeeva function with relative error less than 10^(-r). r can be set by the user subject to the constraints 2 <= r <= 5.

Definition at line 898 of file TMath.cxx.