Logo ROOT   6.10/09
Reference Guide
TMath.h
Go to the documentation of this file.
1 // @(#)root/mathcore:$Id$
2 // Authors: Rene Brun, Anna Kreshuk, Eddy Offermann, Fons Rademakers 29/07/95
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #ifndef ROOT_TMath
13 #define ROOT_TMath
14 
15 
16 //////////////////////////////////////////////////////////////////////////
17 // //
18 // TMath //
19 // //
20 // Encapsulate most frequently used Math functions. //
21 // NB. The basic functions Min, Max, Abs and Sign are defined //
22 // in TMathBase. //
23 // //
24 //////////////////////////////////////////////////////////////////////////
25 
26 #include "Rtypes.h"
27 #include "TMathBase.h"
28 
29 #include "TError.h"
30 #include <algorithm>
31 #include <limits>
32 #include <cmath>
33 
34 namespace TMath {
35 
36  /* ************************* */
37  /* * Fundamental constants * */
38  /* ************************* */
39 
40 constexpr Double_t Pi()
41 {
42  return 3.14159265358979323846;
43 }
44 constexpr Double_t TwoPi()
45 {
46  return 2.0 * Pi();
47 }
48 constexpr Double_t PiOver2()
49 {
50  return Pi() / 2.0;
51 }
52 constexpr Double_t PiOver4()
53 {
54  return Pi() / 4.0;
55 }
56 constexpr Double_t InvPi()
57 {
58  return 1.0 / Pi();
59 }
60 constexpr Double_t RadToDeg()
61 {
62  return 180.0 / Pi();
63 }
64 constexpr Double_t DegToRad()
65 {
66  return Pi() / 180.0;
67 }
68 constexpr Double_t Sqrt2()
69 {
70  return 1.4142135623730950488016887242097;
71 }
72 
73 // e (base of natural log)
74 constexpr Double_t E()
75 {
76  return 2.71828182845904523536;
77 }
78 
79 // natural log of 10 (to convert log to ln)
80 constexpr Double_t Ln10()
81 {
82  return 2.30258509299404568402;
83 }
84 
85 // base-10 log of e (to convert ln to log)
86 constexpr Double_t LogE()
87 {
88  return 0.43429448190325182765;
89 }
90 
91 // velocity of light
92 constexpr Double_t C()
93 {
94  return 2.99792458e8;
95 } // m s^-1
96 constexpr Double_t Ccgs()
97 {
98  return 100.0 * C();
99 } // cm s^-1
101 {
102  return 0.0;
103 } // exact
104 
105 // gravitational constant
106 constexpr Double_t G()
107 {
108  return 6.673e-11;
109 } // m^3 kg^-1 s^-2
110 constexpr Double_t Gcgs()
111 {
112  return G() / 1000.0;
113 } // cm^3 g^-1 s^-2
115 {
116  return 0.010e-11;
117 }
118 
119 // G over h-bar C
120 constexpr Double_t GhbarC()
121 {
122  return 6.707e-39;
123 } // (GeV/c^2)^-2
125 {
126  return 0.010e-39;
127 }
128 
129 // standard acceleration of gravity
130 constexpr Double_t Gn()
131 {
132  return 9.80665;
133 } // m s^-2
135 {
136  return 0.0;
137 } // exact
138 
139 // Planck's constant
140 constexpr Double_t H()
141 {
142  return 6.62606876e-34;
143 } // J s
144 constexpr Double_t Hcgs()
145 {
146  return 1.0e7 * H();
147 } // erg s
149 {
150  return 0.00000052e-34;
151 }
152 
153 // h-bar (h over 2 pi)
154 constexpr Double_t Hbar()
155 {
156  return 1.054571596e-34;
157 } // J s
158 constexpr Double_t Hbarcgs()
159 {
160  return 1.0e7 * Hbar();
161 } // erg s
163 {
164  return 0.000000082e-34;
165 }
166 
167 // hc (h * c)
168 constexpr Double_t HC()
169 {
170  return H() * C();
171 } // J m
172 constexpr Double_t HCcgs()
173 {
174  return Hcgs() * Ccgs();
175 } // erg cm
176 
177 // Boltzmann's constant
178 constexpr Double_t K()
179 {
180  return 1.3806503e-23;
181 } // J K^-1
182 constexpr Double_t Kcgs()
183 {
184  return 1.0e7 * K();
185 } // erg K^-1
187 {
188  return 0.0000024e-23;
189 }
190 
191 // Stefan-Boltzmann constant
192 constexpr Double_t Sigma()
193 {
194  return 5.6704e-8;
195 } // W m^-2 K^-4
197 {
198  return 0.000040e-8;
199 }
200 
201 // Avogadro constant (Avogadro's Number)
202 constexpr Double_t Na()
203 {
204  return 6.02214199e+23;
205 } // mol^-1
207 {
208  return 0.00000047e+23;
209 }
210 
211 // universal gas constant (Na * K)
212 // http://scienceworld.wolfram.com/physics/UniversalGasConstant.html
213 constexpr Double_t R()
214 {
215  return K() * Na();
216 } // J K^-1 mol^-1
218 {
219  return R() * ((KUncertainty() / K()) + (NaUncertainty() / Na()));
220 }
221 
222 // Molecular weight of dry air
223 // 1976 US Standard Atmosphere,
224 // also see http://atmos.nmsu.edu/jsdap/encyclopediawork.html
225 constexpr Double_t MWair()
226 {
227  return 28.9644;
228 } // kg kmol^-1 (or gm mol^-1)
229 
230 // Dry Air Gas Constant (R / MWair)
231 // http://atmos.nmsu.edu/education_and_outreach/encyclopedia/gas_constant.htm
232 constexpr Double_t Rgair()
233 {
234  return (1000.0 * R()) / MWair();
235 } // J kg^-1 K^-1
236 
237 // Euler-Mascheroni Constant
238 constexpr Double_t EulerGamma()
239 {
240  return 0.577215664901532860606512090082402431042;
241 }
242 
243 // Elementary charge
244 constexpr Double_t Qe()
245 {
246  return 1.602176462e-19;
247 } // C
249 {
250  return 0.000000063e-19;
251 }
252 
253 /* ************************** */
254 /* * Mathematical Functions * */
255 /* ************************** */
256 
257 /* ***************************** */
258 /* * Trigonometrical Functions * */
259 /* ***************************** */
260 inline Double_t Sin(Double_t);
261 inline Double_t Cos(Double_t);
262 inline Double_t Tan(Double_t);
263 inline Double_t SinH(Double_t);
264 inline Double_t CosH(Double_t);
265 inline Double_t TanH(Double_t);
266 inline Double_t ASin(Double_t);
267 inline Double_t ACos(Double_t);
268 inline Double_t ATan(Double_t);
274 
275 /* ************************ */
276 /* * Elementary Functions * */
277 /* ************************ */
278 inline Double_t Ceil(Double_t x);
279 inline Int_t CeilNint(Double_t x);
280 inline Double_t Floor(Double_t x);
281 inline Int_t FloorNint(Double_t x);
282 template <typename T>
283 inline Int_t Nint(T x);
284 
285 inline Double_t Sq(Double_t x);
286 inline Double_t Sqrt(Double_t x);
287 inline Double_t Exp(Double_t x);
288 inline Double_t Ldexp(Double_t x, Int_t exp);
292 inline LongDouble_t Power(Long64_t x, Long64_t y);
293 inline Double_t Power(Double_t x, Double_t y);
294 inline Double_t Power(Double_t x, Int_t y);
295 inline Double_t Log(Double_t x);
297 inline Double_t Log10(Double_t x);
298 inline Int_t Finite(Double_t x);
299 inline Int_t Finite(Float_t x);
300 inline Int_t IsNaN(Double_t x);
301 inline Int_t IsNaN(Float_t x);
302 
303 inline Double_t QuietNaN();
304 inline Double_t SignalingNaN();
305 inline Double_t Infinity();
306 
307 template <typename T>
308 struct Limits {
309  inline static T Min();
310  inline static T Max();
311  inline static T Epsilon();
312  };
313 
314  // Some integer math
315  Long_t Hypot(Long_t x, Long_t y); // sqrt(px*px + py*py)
316 
317  // Comparing floating points
319  //return kTRUE if absolute difference between af and bf is less than epsilon
320  return TMath::Abs(af-bf) < epsilon;
321  }
322  inline Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec) {
323  //return kTRUE if relative difference between af and bf is less than relPrec
324  return TMath::Abs(af - bf) <= 0.5 * relPrec * (TMath::Abs(af) + TMath::Abs(bf)) ||
325  TMath::Abs(af - bf) < Limits<Double_t>::Min(); // handle denormals
326  }
327 
328  /* ******************** */
329  /* * Array Algorithms * */
330  /* ******************** */
331 
332  // Min, Max of an array
333  template <typename T> T MinElement(Long64_t n, const T *a);
334  template <typename T> T MaxElement(Long64_t n, const T *a);
335 
336  // Locate Min, Max element number in an array
337  template <typename T> Long64_t LocMin(Long64_t n, const T *a);
338  template <typename Iterator> Iterator LocMin(Iterator first, Iterator last);
339  template <typename T> Long64_t LocMax(Long64_t n, const T *a);
340  template <typename Iterator> Iterator LocMax(Iterator first, Iterator last);
341 
342  // Binary search
343  template <typename T> Long64_t BinarySearch(Long64_t n, const T *array, T value);
344  template <typename T> Long64_t BinarySearch(Long64_t n, const T **array, T value);
345  template <typename Iterator, typename Element> Iterator BinarySearch(Iterator first, Iterator last, Element value);
346 
347  // Hashing
348  ULong_t Hash(const void *txt, Int_t ntxt);
349  ULong_t Hash(const char *str);
350 
351  // Sorting
352  template <typename Element, typename Index>
353  void Sort(Index n, const Element* a, Index* index, Bool_t down=kTRUE);
354  template <typename Iterator, typename IndexIterator>
355  void SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE);
356 
357  void BubbleHigh(Int_t Narr, Double_t *arr1, Int_t *arr2);
358  void BubbleLow (Int_t Narr, Double_t *arr1, Int_t *arr2);
359 
360  Bool_t Permute(Int_t n, Int_t *a); // Find permutations
361 
362  /* ************************* */
363  /* * Geometrical Functions * */
364  /* ************************* */
365 
366  //Sample quantiles
367  void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob,
368  Bool_t isSorted=kTRUE, Int_t *index = 0, Int_t type=7);
369 
370  // IsInside
371  template <typename T> Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y);
372 
373  // Calculate the Cross Product of two vectors
374  template <typename T> T *Cross(const T v1[3],const T v2[3], T out[3]);
375 
376  Float_t Normalize(Float_t v[3]); // Normalize a vector
377  Double_t Normalize(Double_t v[3]); // Normalize a vector
378 
379  //Calculate the Normalized Cross Product of two vectors
380  template <typename T> inline T NormCross(const T v1[3],const T v2[3],T out[3]);
381 
382  // Calculate a normal vector of a plane
383  template <typename T> T *Normal2Plane(const T v1[3],const T v2[3],const T v3[3], T normal[3]);
384 
385  /* ************************ */
386  /* * Polynomial Functions * */
387  /* ************************ */
388 
389  Bool_t RootsCubic(const Double_t coef[4],Double_t &a, Double_t &b, Double_t &c);
390 
391  /* *********************** */
392  /* * Statistic Functions * */
393  /* *********************** */
394 
395  Double_t Binomial(Int_t n,Int_t k); // Calculate the binomial coefficient n over k
404  Double_t KolmogorovTest(Int_t na, const Double_t *a, Int_t nb, const Double_t *b, Option_t *option);
413  Double_t Prob(Double_t chi2,Int_t ndf);
417  Double_t Vavilov(Double_t x, Double_t kappa, Double_t beta2);
418  Double_t VavilovI(Double_t x, Double_t kappa, Double_t beta2);
419  Double_t Voigt(Double_t x, Double_t sigma, Double_t lg, Int_t r = 4);
420 
421  /* ************************** */
422  /* * Statistics over arrays * */
423  /* ************************** */
424 
425  //Mean, Geometric Mean, Median, RMS(sigma)
426 
427  template <typename T> Double_t Mean(Long64_t n, const T *a, const Double_t *w=0);
428  template <typename Iterator> Double_t Mean(Iterator first, Iterator last);
429  template <typename Iterator, typename WeightIterator> Double_t Mean(Iterator first, Iterator last, WeightIterator wfirst);
430 
431  template <typename T> Double_t GeomMean(Long64_t n, const T *a);
432  template <typename Iterator> Double_t GeomMean(Iterator first, Iterator last);
433 
434  template <typename T> Double_t RMS(Long64_t n, const T *a, const Double_t *w=0);
435  template <typename Iterator> Double_t RMS(Iterator first, Iterator last);
436  template <typename Iterator, typename WeightIterator> Double_t RMS(Iterator first, Iterator last, WeightIterator wfirst);
437 
438  template <typename T> Double_t StdDev(Long64_t n, const T *a, const Double_t * w = 0) { return RMS<T>(n,a,w); }
439  template <typename Iterator> Double_t StdDev(Iterator first, Iterator last) { return RMS<Iterator>(first,last); }
440  template <typename Iterator, typename WeightIterator> Double_t StdDev(Iterator first, Iterator last, WeightIterator wfirst) { return RMS<Iterator,WeightIterator>(first,last,wfirst); }
441 
442  template <typename T> Double_t Median(Long64_t n, const T *a, const Double_t *w=0, Long64_t *work=0);
443 
444  //k-th order statistic
445  template <class Element, typename Size> Element KOrdStat(Size n, const Element *a, Size k, Size *work = 0);
446 
447  /* ******************* */
448  /* * Special Functions */
449  /* ******************* */
450 
456 
457  // Bessel functions
458  Double_t BesselI(Int_t n,Double_t x); // integer order modified Bessel function I_n(x)
459  Double_t BesselK(Int_t n,Double_t x); // integer order modified Bessel function K_n(x)
460  Double_t BesselI0(Double_t x); // modified Bessel function I_0(x)
461  Double_t BesselK0(Double_t x); // modified Bessel function K_0(x)
462  Double_t BesselI1(Double_t x); // modified Bessel function I_1(x)
463  Double_t BesselK1(Double_t x); // modified Bessel function K_1(x)
464  Double_t BesselJ0(Double_t x); // Bessel function J0(x) for any real x
465  Double_t BesselJ1(Double_t x); // Bessel function J1(x) for any real x
466  Double_t BesselY0(Double_t x); // Bessel function Y0(x) for positive x
467  Double_t BesselY1(Double_t x); // Bessel function Y1(x) for positive x
468  Double_t StruveH0(Double_t x); // Struve functions of order 0
469  Double_t StruveH1(Double_t x); // Struve functions of order 1
470  Double_t StruveL0(Double_t x); // Modified Struve functions of order 0
471  Double_t StruveL1(Double_t x); // Modified Struve functions of order 1
472 
474  Double_t Erf(Double_t x);
483 }
484 
485 
486 //---- Trig and other functions ------------------------------------------------
487 
488 #include <float.h>
489 
490 #if defined(R__WIN32) && !defined(__CINT__)
491 # ifndef finite
492 # define finite _finite
493 # define isnan _isnan
494 # endif
495 #endif
496 #if defined(R__AIX) || defined(R__SOLARIS_CC50) || \
497  defined(R__HPUX11) || defined(R__GLIBC) || \
498  (defined(R__MACOSX) )
499 // math functions are defined inline so we have to include them here
500 # include <math.h>
501 # ifdef R__SOLARIS_CC50
502  extern "C" { int finite(double); }
503 # endif
504 // # if defined(R__GLIBC) && defined(__STRICT_ANSI__)
505 // # ifndef finite
506 // # define finite __finite
507 // # endif
508 // # ifndef isnan
509 // # define isnan __isnan
510 // # endif
511 // # endif
512 #else
513 // don't want to include complete <math.h>
514 extern "C" {
515  extern double sin(double);
516  extern double cos(double);
517  extern double tan(double);
518  extern double sinh(double);
519  extern double cosh(double);
520  extern double tanh(double);
521  extern double asin(double);
522  extern double acos(double);
523  extern double atan(double);
524  extern double atan2(double, double);
525  extern double sqrt(double);
526  extern double exp(double);
527  extern double pow(double, double);
528  extern double log(double);
529  extern double log10(double);
530 #ifndef R__WIN32
531 # if !defined(finite)
532  extern int finite(double);
533 # endif
534 # if !defined(isnan)
535  extern int isnan(double);
536 # endif
537  extern double ldexp(double, int);
538  extern double ceil(double);
539  extern double floor(double);
540 #else
541  _CRTIMP double ldexp(double, int);
542  _CRTIMP double ceil(double);
543  _CRTIMP double floor(double);
544 #endif
545 }
546 #endif
547 
549  { return sin(x); }
550 
552  { return cos(x); }
553 
555  { return tan(x); }
556 
558  { return sinh(x); }
559 
561  { return cosh(x); }
562 
564  { return tanh(x); }
565 
567  { if (x < -1.) return -TMath::Pi()/2;
568  if (x > 1.) return TMath::Pi()/2;
569  return asin(x);
570  }
571 
573  { if (x < -1.) return TMath::Pi();
574  if (x > 1.) return 0;
575  return acos(x);
576  }
577 
579  { return atan(x); }
580 
582  { if (x != 0) return atan2(y, x);
583  if (y == 0) return 0;
584  if (y > 0) return Pi()/2;
585  else return -Pi()/2;
586  }
587 
589  { return x*x; }
590 
592  { return sqrt(x); }
593 
595  { return ceil(x); }
596 
598  { return TMath::Nint(ceil(x)); }
599 
601  { return floor(x); }
602 
604  { return TMath::Nint(floor(x)); }
605 
606 template<typename T>
608 {
609  // Round to nearest integer. Rounds half integers to the nearest
610  // even integer.
611  int i;
612  if (x >= 0) {
613  i = int(x + 0.5);
614  if ( i & 1 && x + 0.5 == T(i) ) i--;
615  } else {
616  i = int(x - 0.5);
617  if ( i & 1 && x - 0.5 == T(i) ) i++;
618  }
619  return i;
620 }
621 
623  { return exp(x); }
624 
626  { return ldexp(x, exp); }
627 
629  { return std::pow(x,y); }
630 
632  { return std::pow(x,(LongDouble_t)y); }
633 
635  { return std::pow(x,y); }
636 
638  { return pow(x, y); }
639 
641 #ifdef R__ANSISTREAM
642  return std::pow(x, y);
643 #else
644  return pow(x, (Double_t) y);
645 #endif
646 }
647 
648 
650  { return log(x); }
651 
653  { return log10(x); }
654 
656 #if defined(R__FAST_MATH)
657 /* Check if it is finite with a mask in order to be consistent in presence of
658  * fast math.
659  * Inspired from the CMSSW FWCore/Utilities package
660  */
661 {
662  const unsigned long long mask = 0x7FF0000000000000LL;
663  union { unsigned long long l; double d;} v;
664  v.d =x;
665  return (v.l&mask)!=mask;
666 }
667 #else
668 # if defined(R__HPUX11)
669  { return isfinite(x); }
670 # elif defined(R__MACOSX)
671 # ifdef isfinite
672  // from math.h
673  { return isfinite(x); }
674 # else
675  // from cmath
676  { return std::isfinite(x); }
677 # endif
678 # else
679  { return finite(x); }
680 # endif
681 #endif
682 
684 #if defined(R__FAST_MATH)
685 /* Check if it is finite with a mask in order to be consistent in presence of
686  * fast math.
687  * Inspired from the CMSSW FWCore/Utilities package
688  */
689 {
690  const unsigned int mask = 0x7f800000;
691  union { unsigned int l; float d;} v;
692  v.d =x;
693  return (v.l&mask)!=mask;
694 }
695 #else
696 { return std::isfinite(x); }
697 #endif
698 
699 /* This namespace provides all the routines necessary for checking if a number
700  * is a NaN also in presence of optimisations affecting the behaviour of the
701  * floating point calculations.
702  * Inspired from the CMSSW FWCore/Utilities package
703  */
704 #if defined (R__FAST_MATH)
705 namespace ROOT {
706 namespace Internal {
707 namespace Math {
708 // abridged from GNU libc 2.6.1 - in detail from
709 // math/math_private.h
710 // sysdeps/ieee754/ldbl-96/math_ldbl.h
711 
712 // part of ths file:
713  /*
714  * ====================================================
715  * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
716  *
717  * Developed at SunPro, a Sun Microsystems, Inc. business.
718  * Permission to use, copy, modify, and distribute this
719  * software is freely granted, provided that this notice
720  * is preserved.
721  * ====================================================
722  */
723 
724  // A union which permits us to convert between a double and two 32 bit ints.
725  typedef union {
726  Double_t value;
727  struct {
728  UInt_t lsw;
729  UInt_t msw;
730  } parts;
731  } ieee_double_shape_type;
732 
733 #define EXTRACT_WORDS(ix0,ix1,d) \
734  do { \
735  ieee_double_shape_type ew_u; \
736  ew_u.value = (d); \
737  (ix0) = ew_u.parts.msw; \
738  (ix1) = ew_u.parts.lsw; \
739  } while (0)
740 
741  inline Int_t IsNaN(Double_t x)
742  {
743  UInt_t hx, lx;
744 
745  EXTRACT_WORDS(hx, lx, x);
746 
747  lx |= hx & 0xfffff;
748  hx &= 0x7ff00000;
749  return (hx == 0x7ff00000) && (lx != 0);
750  }
751 
752  typedef union {
753  Float_t value;
754  UInt_t word;
755  } ieee_float_shape_type;
756 
757 #define GET_FLOAT_WORD(i,d) \
758  do { \
759  ieee_float_shape_type gf_u; \
760  gf_u.value = (d); \
761  (i) = gf_u.word; \
762  } while (0)
763 
764  inline Int_t IsNaN(Float_t x)
765  {
766  UInt_t wx;
767  GET_FLOAT_WORD (wx, x);
768  wx &= 0x7fffffff;
769  return (Bool_t)(wx > 0x7f800000);
770  }
771 } } } // end NS ROOT::Internal::Math
772 #endif // End R__FAST_MATH
773 
774 #if defined(R__FAST_MATH)
775  inline Int_t TMath::IsNaN(Double_t x) { return ROOT::Internal::Math::IsNaN(x); }
776  inline Int_t TMath::IsNaN(Float_t x) { return ROOT::Internal::Math::IsNaN(x); }
777 #else
778  inline Int_t TMath::IsNaN(Double_t x) { return std::isnan(x); }
779  inline Int_t TMath::IsNaN(Float_t x) { return std::isnan(x); }
780 #endif
781 
782 //--------wrapper to numeric_limits
783 //____________________________________________________________________________
785  // returns a quiet NaN as defined by IEEE 754
786  // see http://en.wikipedia.org/wiki/NaN#Quiet_NaN
787  return std::numeric_limits<Double_t>::quiet_NaN();
788 }
789 
790 //____________________________________________________________________________
792  // returns a signaling NaN as defined by IEEE 754
793  // see http://en.wikipedia.org/wiki/NaN#Signaling_NaN
794  return std::numeric_limits<Double_t>::signaling_NaN();
795 }
796 
798  // returns an infinity as defined by the IEEE standard
799  return std::numeric_limits<Double_t>::infinity();
800 }
801 
802 template<typename T>
804  // returns maximum representation for type T
805  return (std::numeric_limits<T>::min)(); //N.B. use this signature to avoid class with macro min() on Windows
806 }
807 
808 template<typename T>
810  // returns minimum double representation
811  return (std::numeric_limits<T>::max)(); //N.B. use this signature to avoid class with macro max() on Windows
812 }
813 
814 template<typename T>
816  // returns minimum double representation
818 }
819 
820 
821 //-------- Advanced -------------
822 
823 template <typename T> inline T TMath::NormCross(const T v1[3],const T v2[3],T out[3])
824 {
825  // Calculate the Normalized Cross Product of two vectors
826  return Normalize(Cross(v1,v2,out));
827 }
828 
829 template <typename T>
831  // Return minimum of array a of length n.
832 
833  return *std::min_element(a,a+n);
834 }
835 
836 template <typename T>
838  // Return maximum of array a of length n.
839 
840  return *std::max_element(a,a+n);
841 }
842 
843 template <typename T>
845  // Return index of array with the minimum element.
846  // If more than one element is minimum returns first found.
847 
848  // Implement here since this one is found to be faster (mainly on 64 bit machines)
849  // than stl generic implementation.
850  // When performing the comparison, the STL implementation needs to de-reference both the array iterator
851  // and the iterator pointing to the resulting minimum location
852 
853  if (n <= 0 || !a) return -1;
854  T xmin = a[0];
855  Long64_t loc = 0;
856  for (Long64_t i = 1; i < n; i++) {
857  if (xmin > a[i]) {
858  xmin = a[i];
859  loc = i;
860  }
861  }
862  return loc;
863 }
864 
865 template <typename Iterator>
866 Iterator TMath::LocMin(Iterator first, Iterator last) {
867  // Return index of array with the minimum element.
868  // If more than one element is minimum returns first found.
869  return std::min_element(first, last);
870 }
871 
872 template <typename T>
874  // Return index of array with the maximum element.
875  // If more than one element is maximum returns first found.
876 
877  // Implement here since it is faster (see comment in LocMin function)
878 
879  if (n <= 0 || !a) return -1;
880  T xmax = a[0];
881  Long64_t loc = 0;
882  for (Long64_t i = 1; i < n; i++) {
883  if (xmax < a[i]) {
884  xmax = a[i];
885  loc = i;
886  }
887  }
888  return loc;
889 }
890 
891 template <typename Iterator>
892 Iterator TMath::LocMax(Iterator first, Iterator last)
893 {
894  // Return index of array with the maximum element.
895  // If more than one element is maximum returns first found.
896 
897  return std::max_element(first, last);
898 }
899 
900 template<typename T>
901 struct CompareDesc {
902 
903  CompareDesc(T d) : fData(d) {}
904 
905  template<typename Index>
906  bool operator()(Index i1, Index i2) {
907  return *(fData + i1) > *(fData + i2);
908  }
909 
911 };
912 
913 template<typename T>
914 struct CompareAsc {
915 
916  CompareAsc(T d) : fData(d) {}
917 
918  template<typename Index>
919  bool operator()(Index i1, Index i2) {
920  return *(fData + i1) < *(fData + i2);
921  }
922 
924 };
925 
926 template <typename Iterator>
927 Double_t TMath::Mean(Iterator first, Iterator last)
928 {
929  // Return the weighted mean of an array defined by the iterators.
930 
931  Double_t sum = 0;
932  Double_t sumw = 0;
933  while ( first != last )
934  {
935  sum += *first;
936  sumw += 1;
937  first++;
938  }
939 
940  return sum/sumw;
941 }
942 
943 template <typename Iterator, typename WeightIterator>
944 Double_t TMath::Mean(Iterator first, Iterator last, WeightIterator w)
945 {
946  // Return the weighted mean of an array defined by the first and
947  // last iterators. The w iterator should point to the first element
948  // of a vector of weights of the same size as the main array.
949 
950  Double_t sum = 0;
951  Double_t sumw = 0;
952  int i = 0;
953  while ( first != last ) {
954  if ( *w < 0) {
955  ::Error("TMath::Mean","w[%d] = %.4e < 0 ?!",i,*w);
956  return 0;
957  }
958  sum += (*w) * (*first);
959  sumw += (*w) ;
960  ++w;
961  ++first;
962  ++i;
963  }
964  if (sumw <= 0) {
965  ::Error("TMath::Mean","sum of weights == 0 ?!");
966  return 0;
967  }
968 
969  return sum/sumw;
970 }
971 
972 template <typename T>
974 {
975  // Return the weighted mean of an array a with length n.
976 
977  if (w) {
978  return TMath::Mean(a, a+n, w);
979  } else {
980  return TMath::Mean(a, a+n);
981  }
982 }
983 
984 template <typename Iterator>
985 Double_t TMath::GeomMean(Iterator first, Iterator last)
986 {
987  // Return the geometric mean of an array defined by the iterators.
988  // geometric_mean = (Prod_i=0,n-1 |a[i]|)^1/n
989 
990  Double_t logsum = 0.;
991  Long64_t n = 0;
992  while ( first != last ) {
993  if (*first == 0) return 0.;
994  Double_t absa = (Double_t) TMath::Abs(*first);
995  logsum += TMath::Log(absa);
996  ++first;
997  ++n;
998  }
999 
1000  return TMath::Exp(logsum/n);
1001 }
1002 
1003 template <typename T>
1005 {
1006  // Return the geometric mean of an array a of size n.
1007  // geometric_mean = (Prod_i=0,n-1 |a[i]|)^1/n
1008 
1009  return TMath::GeomMean(a, a+n);
1010 }
1011 
1012 template <typename Iterator>
1013 Double_t TMath::RMS(Iterator first, Iterator last)
1014 {
1015  // Return the Standard Deviation of an array defined by the iterators.
1016  // Note that this function returns the sigma(standard deviation) and
1017  // not the root mean square of the array.
1018 
1019  // Use the two pass algorithm, which is slower (! a factor of 2) but much more
1020  // precise. Since we have a vector the 2 pass algorithm is still faster than the
1021  // Welford algorithm. (See also ROOT-5545)
1022 
1023  Double_t n = 0;
1024 
1025  Double_t tot = 0;
1026  Double_t mean = TMath::Mean(first,last);
1027  while ( first != last ) {
1028  Double_t x = Double_t(*first);
1029  tot += (x - mean)*(x - mean);
1030  ++first;
1031  ++n;
1032  }
1033  Double_t rms = (n > 1) ? TMath::Sqrt(tot/(n-1)) : 0.0;
1034  return rms;
1035 }
1036 
1037 template <typename Iterator, typename WeightIterator>
1038 Double_t TMath::RMS(Iterator first, Iterator last, WeightIterator w)
1039 {
1040  // Return the weighted Standard Deviation of an array defined by the iterators.
1041  // Note that this function returns the sigma(standard deviation) and
1042  // not the root mean square of the array.
1043 
1044  // As in the unweighted case use the two pass algorithm
1045 
1046  Double_t tot = 0;
1047  Double_t sumw = 0;
1048  Double_t sumw2 = 0;
1049  Double_t mean = TMath::Mean(first,last,w);
1050  while ( first != last ) {
1051  Double_t x = Double_t(*first);
1052  sumw += *w;
1053  sumw2 += (*w) * (*w);
1054  tot += (*w) * (x - mean)*(x - mean);
1055  ++first;
1056  ++w;
1057  }
1058  // use the correction neff/(neff -1) for the unbiased formula
1059  Double_t rms = TMath::Sqrt(tot * sumw/ (sumw*sumw - sumw2) );
1060  return rms;
1061 }
1062 
1063 
1064 template <typename T>
1065 Double_t TMath::RMS(Long64_t n, const T *a, const Double_t * w)
1066 {
1067  // Return the Standard Deviation of an array a with length n.
1068  // Note that this function returns the sigma(standard deviation) and
1069  // not the root mean square of the array.
1070 
1071  return (w) ? TMath::RMS(a, a+n, w) : TMath::RMS(a, a+n);
1072 }
1073 
1074 template <typename Iterator, typename Element>
1075 Iterator TMath::BinarySearch(Iterator first, Iterator last, Element value)
1076 {
1077  // Binary search in an array defined by its iterators.
1078  //
1079  // The values in the iterators range are supposed to be sorted
1080  // prior to this call. If match is found, function returns
1081  // position of element. If no match found, function gives nearest
1082  // element smaller than value.
1083 
1084  Iterator pind;
1085  pind = std::lower_bound(first, last, value);
1086  if ( (pind != last) && (*pind == value) )
1087  return pind;
1088  else
1089  return ( pind - 1);
1090 }
1091 
1092 
1093 template <typename T> Long64_t TMath::BinarySearch(Long64_t n, const T *array, T value)
1094 {
1095  // Binary search in an array of n values to locate value.
1096  //
1097  // Array is supposed to be sorted prior to this call.
1098  // If match is found, function returns position of element.
1099  // If no match found, function gives nearest element smaller than value.
1100 
1101  const T* pind;
1102  pind = std::lower_bound(array, array + n, value);
1103  if ( (pind != array + n) && (*pind == value) )
1104  return (pind - array);
1105  else
1106  return ( pind - array - 1);
1107 }
1108 
1109 template <typename T> Long64_t TMath::BinarySearch(Long64_t n, const T **array, T value)
1110 {
1111  // Binary search in an array of n values to locate value.
1112  //
1113  // Array is supposed to be sorted prior to this call.
1114  // If match is found, function returns position of element.
1115  // If no match found, function gives nearest element smaller than value.
1116 
1117  const T* pind;
1118  pind = std::lower_bound(*array, *array + n, value);
1119  if ( (pind != *array + n) && (*pind == value) )
1120  return (pind - *array);
1121  else
1122  return ( pind - *array - 1);
1123 }
1124 
1125 template <typename Iterator, typename IndexIterator>
1126 void TMath::SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down)
1127 {
1128  // Sort the n1 elements of the Short_t array defined by its
1129  // iterators. In output the array index contains the indices of
1130  // the sorted array. If down is false sort in increasing order
1131  // (default is decreasing order).
1132 
1133  // NOTE that the array index must be created with a length bigger
1134  // or equal than the main array before calling this function.
1135 
1136  int i = 0;
1137 
1138  IndexIterator cindex = index;
1139  for ( Iterator cfirst = first; cfirst != last; ++cfirst )
1140  {
1141  *cindex = i++;
1142  ++cindex;
1143  }
1144 
1145  if ( down )
1146  std::sort(index, cindex, CompareDesc<Iterator>(first) );
1147  else
1148  std::sort(index, cindex, CompareAsc<Iterator>(first) );
1149 }
1150 
1151 template <typename Element, typename Index> void TMath::Sort(Index n, const Element* a, Index* index, Bool_t down)
1152 {
1153  // Sort the n elements of the array a of generic templated type Element.
1154  // In output the array index of type Index contains the indices of the sorted array.
1155  // If down is false sort in increasing order (default is decreasing order).
1156 
1157  // NOTE that the array index must be created with a length >= n
1158  // before calling this function.
1159  // NOTE also that the size type for n must be the same type used for the index array
1160  // (templated type Index)
1161 
1162  for(Index i = 0; i < n; i++) { index[i] = i; }
1163  if ( down )
1164  std::sort(index, index + n, CompareDesc<const Element*>(a) );
1165  else
1166  std::sort(index, index + n, CompareAsc<const Element*>(a) );
1167 }
1168 
1169 template <typename T> T *TMath::Cross(const T v1[3],const T v2[3], T out[3])
1170 {
1171  // Calculate the Cross Product of two vectors:
1172  // out = [v1 x v2]
1173 
1174  out[0] = v1[1] * v2[2] - v1[2] * v2[1];
1175  out[1] = v1[2] * v2[0] - v1[0] * v2[2];
1176  out[2] = v1[0] * v2[1] - v1[1] * v2[0];
1177 
1178  return out;
1179 }
1180 
1181 template <typename T> T * TMath::Normal2Plane(const T p1[3],const T p2[3],const T p3[3], T normal[3])
1182 {
1183  // Calculate a normal vector of a plane.
1184  //
1185  // Input:
1186  // Float_t *p1,*p2,*p3 - 3 3D points belonged the plane to define it.
1187  //
1188  // Return:
1189  // Pointer to 3D normal vector (normalized)
1190 
1191  T v1[3], v2[3];
1192 
1193  v1[0] = p2[0] - p1[0];
1194  v1[1] = p2[1] - p1[1];
1195  v1[2] = p2[2] - p1[2];
1196 
1197  v2[0] = p3[0] - p1[0];
1198  v2[1] = p3[1] - p1[1];
1199  v2[2] = p3[2] - p1[2];
1200 
1201  NormCross(v1,v2,normal);
1202  return normal;
1203 }
1204 
1205 template <typename T> Bool_t TMath::IsInside(T xp, T yp, Int_t np, T *x, T *y)
1206 {
1207  // Function which returns kTRUE if point xp,yp lies inside the
1208  // polygon defined by the np points in arrays x and y, kFALSE otherwise.
1209  // Note that the polygon may be open or closed.
1210 
1211  Int_t i, j = np-1 ;
1212  Bool_t oddNodes = kFALSE;
1213 
1214  for (i=0; i<np; i++) {
1215  if ((y[i]<yp && y[j]>=yp) || (y[j]<yp && y[i]>=yp)) {
1216  if (x[i]+(yp-y[i])/(y[j]-y[i])*(x[j]-x[i])<xp) {
1217  oddNodes = !oddNodes;
1218  }
1219  }
1220  j=i;
1221  }
1222 
1223  return oddNodes;
1224 }
1225 
1226 template <typename T> Double_t TMath::Median(Long64_t n, const T *a, const Double_t *w, Long64_t *work)
1227 {
1228  // Return the median of the array a where each entry i has weight w[i] .
1229  // Both arrays have a length of at least n . The median is a number obtained
1230  // from the sorted array a through
1231  //
1232  // median = (a[jl]+a[jh])/2. where (using also the sorted index on the array w)
1233  //
1234  // sum_i=0,jl w[i] <= sumTot/2
1235  // sum_i=0,jh w[i] >= sumTot/2
1236  // sumTot = sum_i=0,n w[i]
1237  //
1238  // If w=0, the algorithm defaults to the median definition where it is
1239  // a number that divides the sorted sequence into 2 halves.
1240  // When n is odd or n > 1000, the median is kth element k = (n + 1) / 2.
1241  // when n is even and n < 1000the median is a mean of the elements k = n/2 and k = n/2 + 1.
1242  //
1243  // If the weights are supplied (w not 0) all weights must be >= 0
1244  //
1245  // If work is supplied, it is used to store the sorting index and assumed to be
1246  // >= n . If work=0, local storage is used, either on the stack if n < kWorkMax
1247  // or on the heap for n >= kWorkMax .
1248 
1249  const Int_t kWorkMax = 100;
1250 
1251  if (n <= 0 || !a) return 0;
1252  Bool_t isAllocated = kFALSE;
1253  Double_t median;
1254  Long64_t *ind;
1255  Long64_t workLocal[kWorkMax];
1256 
1257  if (work) {
1258  ind = work;
1259  } else {
1260  ind = workLocal;
1261  if (n > kWorkMax) {
1262  isAllocated = kTRUE;
1263  ind = new Long64_t[n];
1264  }
1265  }
1266 
1267  if (w) {
1268  Double_t sumTot2 = 0;
1269  for (Int_t j = 0; j < n; j++) {
1270  if (w[j] < 0) {
1271  ::Error("TMath::Median","w[%d] = %.4e < 0 ?!",j,w[j]);
1272  if (isAllocated) delete [] ind;
1273  return 0;
1274  }
1275  sumTot2 += w[j];
1276  }
1277 
1278  sumTot2 /= 2.;
1279 
1280  Sort(n, a, ind, kFALSE);
1281 
1282  Double_t sum = 0.;
1283  Int_t jl;
1284  for (jl = 0; jl < n; jl++) {
1285  sum += w[ind[jl]];
1286  if (sum >= sumTot2) break;
1287  }
1288 
1289  Int_t jh;
1290  sum = 2.*sumTot2;
1291  for (jh = n-1; jh >= 0; jh--) {
1292  sum -= w[ind[jh]];
1293  if (sum <= sumTot2) break;
1294  }
1295 
1296  median = 0.5*(a[ind[jl]]+a[ind[jh]]);
1297 
1298  } else {
1299 
1300  if (n%2 == 1)
1301  median = KOrdStat(n, a,n/2, ind);
1302  else {
1303  median = 0.5*(KOrdStat(n, a, n/2 -1, ind)+KOrdStat(n, a, n/2, ind));
1304  }
1305  }
1306 
1307  if (isAllocated)
1308  delete [] ind;
1309  return median;
1310 }
1311 
1312 
1313 
1314 
1315 template <class Element, typename Size>
1316 Element TMath::KOrdStat(Size n, const Element *a, Size k, Size *work)
1317 {
1318  // Returns k_th order statistic of the array a of size n
1319  // (k_th smallest element out of n elements).
1320  //
1321  // C-convention is used for array indexing, so if you want
1322  // the second smallest element, call KOrdStat(n, a, 1).
1323  //
1324  // If work is supplied, it is used to store the sorting index and
1325  // assumed to be >= n. If work=0, local storage is used, either on
1326  // the stack if n < kWorkMax or on the heap for n >= kWorkMax.
1327  // Note that the work index array will not contain the sorted indices but
1328  // all indeces of the smaller element in arbitrary order in work[0,...,k-1] and
1329  // all indeces of the larger element in arbitrary order in work[k+1,..,n-1]
1330  // work[k] will contain instead the index of the returned element.
1331  //
1332  // Taken from "Numerical Recipes in C++" without the index array
1333  // implemented by Anna Khreshuk.
1334  //
1335  // See also the declarations at the top of this file
1336 
1337  const Int_t kWorkMax = 100;
1338 
1339  typedef Size Index;
1340 
1341  Bool_t isAllocated = kFALSE;
1342  Size i, ir, j, l, mid;
1343  Index arr;
1344  Index *ind;
1345  Index workLocal[kWorkMax];
1346  Index temp;
1347 
1348  if (work) {
1349  ind = work;
1350  } else {
1351  ind = workLocal;
1352  if (n > kWorkMax) {
1353  isAllocated = kTRUE;
1354  ind = new Index[n];
1355  }
1356  }
1357 
1358  for (Size ii=0; ii<n; ii++) {
1359  ind[ii]=ii;
1360  }
1361  Size rk = k;
1362  l=0;
1363  ir = n-1;
1364  for(;;) {
1365  if (ir<=l+1) { //active partition contains 1 or 2 elements
1366  if (ir == l+1 && a[ind[ir]]<a[ind[l]])
1367  {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1368  Element tmp = a[ind[rk]];
1369  if (isAllocated)
1370  delete [] ind;
1371  return tmp;
1372  } else {
1373  mid = (l+ir) >> 1; //choose median of left, center and right
1374  {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}//elements as partitioning element arr.
1375  if (a[ind[l]]>a[ind[ir]]) //also rearrange so that a[l]<=a[l+1]
1376  {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1377 
1378  if (a[ind[l+1]]>a[ind[ir]])
1379  {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}
1380 
1381  if (a[ind[l]]>a[ind[l+1]])
1382  {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}
1383 
1384  i=l+1; //initialize pointers for partitioning
1385  j=ir;
1386  arr = ind[l+1];
1387  for (;;){
1388  do i++; while (a[ind[i]]<a[arr]);
1389  do j--; while (a[ind[j]]>a[arr]);
1390  if (j<i) break; //pointers crossed, partitioning complete
1391  {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
1392  }
1393  ind[l+1]=ind[j];
1394  ind[j]=arr;
1395  if (j>=rk) ir = j-1; //keep active the partition that
1396  if (j<=rk) l=i; //contains the k_th element
1397  }
1398  }
1399 }
1400 
1401 #endif
constexpr Double_t H()
Definition: TMath.h:140
Double_t ACosH(Double_t)
Definition: TMath.cxx:80
Double_t BesselI(Int_t n, Double_t x)
Compute the Integer Order Modified Bessel function I_n(x) for n=0,1,2,...
Definition: TMath.cxx:1559
constexpr Double_t QeUncertainty()
Definition: TMath.h:248
double par[1]
Definition: unuranDistr.cxx:38
constexpr Double_t HC()
Definition: TMath.h:168
Double_t FDist(Double_t F, Double_t N, Double_t M)
Computes the density function of F-distribution (probability function, integral of density...
Definition: TMath.cxx:2230
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition: TMath.cxx:472
Double_t ErfInverse(Double_t x)
returns the inverse error function x must be <-1<x<1
Definition: TMath.cxx:206
static long int sum(long int i)
Definition: Factory.cxx:2162
Double_t BreitWigner(Double_t x, Double_t mean=0, Double_t gamma=1)
Calculate a Breit Wigner function with mean and gamma.
Definition: TMath.cxx:441
float xmin
Definition: THbookFile.cxx:93
constexpr Double_t K()
Definition: TMath.h:178
Double_t FDistI(Double_t F, Double_t N, Double_t M)
Calculates the cumulative distribution function of F-distribution, this function occurs in the statis...
Definition: TMath.cxx:2248
Double_t PoissonI(Double_t x, Double_t par)
compute the Poisson distribution function for (x,par) This is a non-smooth function.
Definition: TMath.cxx:602
Double_t Floor(Double_t x)
Definition: TMath.h:600
double tanh(double)
Double_t TanH(Double_t)
Definition: TMath.h:563
Double_t Sq(Double_t x)
Definition: TMath.h:588
long long Long64_t
Definition: RtypesCore.h:69
Float_t Normalize(Float_t v[3])
Normalize a vector v in place.
Definition: TMath.cxx:498
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t mid
Definition: TRolke.cxx:630
Long64_t LocMax(Long64_t n, const T *a)
Definition: TMath.h:873
bool operator()(Index i1, Index i2)
Definition: TMath.h:906
static double p3(double t, double a, double b, double c, double d)
Namespace for new ROOT classes and functions.
Definition: StringConv.hxx:21
Double_t Log(Double_t x)
Definition: TMath.h:649
float Float_t
Definition: RtypesCore.h:53
constexpr Double_t Sqrt2()
Definition: TMath.h:68
const char Option_t
Definition: RtypesCore.h:62
Double_t KolmogorovProb(Double_t z)
Calculates the Kolmogorov distribution function,.
Definition: TMath.cxx:665
constexpr Double_t Hcgs()
Definition: TMath.h:144
constexpr Double_t TwoPi()
Definition: TMath.h:44
const char * Size
Definition: TXMLSetup.cxx:56
Double_t QuietNaN()
Definition: TMath.h:784
double T(double x)
Definition: ChebyshevPol.h:34
Double_t LaplaceDistI(Double_t x, Double_t alpha=0, Double_t beta=1)
Computes the distribution function of Laplace distribution at point x, with location parameter alpha ...
Definition: TMath.cxx:2329
Double_t StudentQuantile(Double_t p, Double_t ndf, Bool_t lower_tail=kTRUE)
Computes quantiles of the Student&#39;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&#39;s t-quantiles" "Communications of the ACM", 13(10), October 1970.
Definition: TMath.cxx:2613
Double_t NormQuantile(Double_t p)
Computes quantiles for standard normal distribution N(0, 1) at probability p ALGORITHM AS241 APPL...
Definition: TMath.cxx:2402
Double_t BetaDist(Double_t x, Double_t p, Double_t q)
Computes the probability density function of the Beta distribution (the distribution function is comp...
Definition: TMath.cxx:2039
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...
Definition: TMath.cxx:2313
constexpr Double_t Rgair()
Definition: TMath.h:232
Double_t Log2(Double_t x)
Definition: TMath.cxx:104
#define N
Double_t StruveH1(Double_t x)
Struve Functions of Order 1.
Definition: TMath.cxx:1815
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
CompareAsc(T d)
Definition: TMath.h:916
TArc * a
Definition: textangle.C:12
constexpr Double_t SigmaUncertainty()
Definition: TMath.h:196
void BubbleHigh(Int_t Narr, Double_t *arr1, Int_t *arr2)
Bubble sort variant to obtain the order of an array&#39;s elements into an index in order to do more usef...
Definition: TMath.cxx:1283
Int_t FloorNint(Double_t x)
Definition: TMath.h:603
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:352
constexpr Double_t InvPi()
Definition: TMath.h:56
Double_t GeomMean(Long64_t n, const T *a)
Definition: TMath.h:1004
Double_t BesselJ1(Double_t x)
Returns the Bessel function J1(x) for any real x.
Definition: TMath.cxx:1638
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
double cos(double)
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...
Definition: TMath.cxx:624
Double_t StudentI(Double_t T, Double_t ndf)
Calculates the cumulative distribution function of Student&#39;s t-distribution second parameter stands f...
Definition: TMath.cxx:2589
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:628
double beta(double x, double y)
Calculates the beta function.
Double_t Ldexp(Double_t x, Int_t exp)
Definition: TMath.h:625
T * Normal2Plane(const T v1[3], const T v2[3], const T v3[3], T normal[3])
Definition: TMath.h:1181
constexpr Double_t Gn()
Definition: TMath.h:130
constexpr Double_t Ln10()
Definition: TMath.h:80
constexpr Double_t Hbarcgs()
Definition: TMath.h:158
Double_t RMS(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:1065
T NormCross(const T v1[3], const T v2[3], T out[3])
Definition: TMath.h:823
double sqrt(double)
Double_t BesselY1(Double_t x)
Returns the Bessel function Y1(x) for positive x.
Definition: TMath.cxx:1708
T fData
Definition: TMath.h:923
Double_t ChisquareQuantile(Double_t p, Double_t ndf)
Evaluate the quantiles of the chi-squared probability distribution function.
Definition: TMath.cxx:2147
Double_t x[n]
Definition: legend1.C:17
double acos(double)
double sinh(double)
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=0, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities Parameters: x -the data sample n ...
Definition: TMath.cxx:1177
double log10(double)
Double_t ASinH(Double_t)
Definition: TMath.cxx:67
Double_t Log10(Double_t x)
Definition: TMath.h:652
Double_t BesselI1(Double_t x)
Compute the modified Bessel function I_1(x) for any real x.
Definition: TMath.cxx:1463
static double p2(double t, double a, double b, double c)
Double_t Freq(Double_t x)
Computation of the normal frequency function freq(x).
Definition: TMath.cxx:268
constexpr Double_t GhbarCUncertainty()
Definition: TMath.h:124
constexpr Double_t HCcgs()
Definition: TMath.h:172
double pow(double, double)
bool operator()(Index i1, Index i2)
Definition: TMath.h:919
constexpr Double_t DegToRad()
Definition: TMath.h:64
constexpr Double_t MWair()
Definition: TMath.h:225
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Definition: TMath.h:1151
constexpr Double_t Ccgs()
Definition: TMath.h:96
Int_t Finite(Double_t x)
Definition: TMath.h:655
const Double_t sigma
constexpr Double_t PiOver4()
Definition: TMath.h:52
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:581
Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y)
Definition: TMath.h:1205
constexpr Double_t Pi()
Definition: TMath.h:40
double sin(double)
Double_t Infinity()
Definition: TMath.h:797
void BubbleLow(Int_t Narr, Double_t *arr1, Int_t *arr2)
Opposite ordering of the array arr2[] to that of BubbleHigh.
Definition: TMath.cxx:1322
void Error(const char *location, const char *msgfmt,...)
Double_t Erfc(Double_t x)
Compute the complementary error function erfc(x).
Definition: TMath.cxx:197
constexpr Double_t GnUncertainty()
Definition: TMath.h:134
Double_t DiLog(Double_t x)
The DiLogarithm function Code translated by R.Brun from CERNLIB DILOG function C332.
Definition: TMath.cxx:113
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Definition: TMath.h:318
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Definition: TMath.h:322
#define F(x, y, z)
int isnan(double)
double gamma(double x)
constexpr Double_t Qe()
Definition: TMath.h:244
constexpr Double_t G()
Definition: TMath.h:106
TRandom2 r(17)
Double_t VavilovI(Double_t x, Double_t kappa, Double_t beta2)
Returns the value of the Vavilov distribution function Parameters: 1st - the point were the density f...
Definition: TMath.cxx:2738
SVector< double, 2 > v
Definition: Dict.h:5
Double_t BesselY0(Double_t x)
Returns the Bessel function Y0(x) for positive x.
Definition: TMath.cxx:1674
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition: TMath.cxx:187
Double_t BetaDistI(Double_t x, Double_t p, Double_t q)
Computes the distribution function of the Beta distribution.
Definition: TMath.cxx:2057
constexpr Double_t NaUncertainty()
Definition: TMath.h:206
constexpr Double_t RUncertainty()
Definition: TMath.h:217
Double_t StruveL1(Double_t x)
Modified Struve Function of Order 1.
Definition: TMath.cxx:1938
double cosh(double)
Double_t Voigt(Double_t x, Double_t sigma, Double_t lg, Int_t r=4)
Computation of Voigt function (normalised).
Definition: TMath.cxx:879
Bool_t Permute(Int_t n, Int_t *a)
Simple recursive algorithm to find the permutations of n natural numbers, not necessarily all distinc...
Definition: TMath.cxx:2503
Double_t StruveL0(Double_t x)
Modified Struve Function of Order 0.
Definition: TMath.cxx:1892
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:973
unsigned int UInt_t
Definition: RtypesCore.h:42
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 ...
Definition: TMath.cxx:788
TMarker * m
Definition: textangle.C:8
Double_t Binomial(Int_t n, Int_t k)
Calculate the binomial coefficient n over k.
Definition: TMath.cxx:2078
double floor(double)
Double_t LandauI(Double_t x)
Returns the value of the Landau distribution function at point x.
Definition: TMath.cxx:2767
TLine * l
Definition: textangle.C:4
long double LongDouble_t
Definition: RtypesCore.h:57
Double_t ACos(Double_t)
Definition: TMath.h:572
Double_t BesselJ0(Double_t x)
Returns the Bessel function J0(x) for any real x.
Definition: TMath.cxx:1603
static double p1(double t, double a, double b)
float xmax
Definition: THbookFile.cxx:93
Double_t ErfcInverse(Double_t x)
Definition: TMath.cxx:233
Double_t SignalingNaN()
Definition: TMath.h:791
T * Cross(const T v1[3], const T v2[3], T out[3])
Definition: TMath.h:1169
double asin(double)
ULong_t Hash(const void *txt, Int_t ntxt)
Calculates hash index from any char string.
Definition: TMath.cxx:1377
constexpr Double_t Na()
Definition: TMath.h:202
constexpr Double_t HUncertainty()
Definition: TMath.h:148
constexpr Double_t GUncertainty()
Definition: TMath.h:114
REAL epsilon
Definition: triangle.c:617
static T Min()
Definition: TMath.h:803
constexpr Double_t E()
Definition: TMath.h:74
Double_t Poisson(Double_t x, Double_t par)
Compute the Poisson distribution function for (x,par) The Poisson PDF is implemented by means of Eule...
Definition: TMath.cxx:571
Double_t Student(Double_t T, Double_t ndf)
Computes density function for Student&#39;s t- distribution (the probability function (integral of densit...
Definition: TMath.cxx:2568
Double_t Beta(Double_t p, Double_t q)
Calculates Beta-function Gamma(p)*Gamma(q)/Gamma(p+q).
Definition: TMath.cxx:1979
Double_t Cos(Double_t)
Definition: TMath.h:551
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:452
const Bool_t kFALSE
Definition: RtypesCore.h:92
Double_t BesselI0(Double_t x)
Compute the modified Bessel function I_0(x) for any real x.
Definition: TMath.cxx:1395
Double_t BesselK(Int_t n, Double_t x)
Compute the Integer Order Modified Bessel function K_n(x) for n=0,1,2,...
Definition: TMath.cxx:1530
constexpr Double_t HbarUncertainty()
Definition: TMath.h:162
Double_t StdDev(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:438
long Long_t
Definition: RtypesCore.h:50
RooCmdArg Index(RooCategory &icat)
Double_t Exp(Double_t x)
Definition: TMath.h:622
double Double_t
Definition: RtypesCore.h:55
Double_t StruveH0(Double_t x)
Struve Functions of Order 0.
Definition: TMath.cxx:1746
double atan2(double, double)
int type
Definition: TGX11.cxx:120
Double_t Median(Long64_t n, const T *a, const Double_t *w=0, Long64_t *work=0)
Definition: TMath.h:1226
unsigned long ULong_t
Definition: RtypesCore.h:51
T MaxElement(Long64_t n, const T *a)
Definition: TMath.h:837
Double_t y[n]
Definition: legend1.C:17
double atan(double)
Double_t Hypot(Double_t x, Double_t y)
Definition: TMath.cxx:60
constexpr Double_t C()
Definition: TMath.h:92
constexpr Double_t Hbar()
Definition: TMath.h:154
Namespace for new Math classes and functions.
Double_t BesselK0(Double_t x)
Compute the modified Bessel function K_0(x) for positive real x.
Definition: TMath.cxx:1429
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 mor...
Definition: TMath.cxx:2103
constexpr Double_t EulerGamma()
Definition: TMath.h:238
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
constexpr Double_t Kcgs()
Definition: TMath.h:182
constexpr Double_t RadToDeg()
Definition: TMath.h:60
constexpr Double_t CUncertainty()
Definition: TMath.h:100
Double_t BetaCf(Double_t x, Double_t a, Double_t b)
Continued fraction evaluation by modified Lentz&#39;s method used in calculation of incomplete Beta funct...
Definition: TMath.cxx:1988
Element KOrdStat(Size n, const Element *a, Size k, Size *work=0)
Definition: TMath.h:1316
CompareDesc(T d)
Definition: TMath.h:903
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.
Definition: TMath.cxx:2296
Int_t IsNaN(Float_t x)
Definition: TMath.h:779
double ceil(double)
Int_t IsNaN(Double_t x)
Definition: TMath.h:778
int finite(double)
double tan(double)
Double_t Factorial(Int_t i)
Compute factorial(n).
Definition: TMath.cxx:250
Double_t Sin(Double_t)
Definition: TMath.h:548
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
Double_t Ceil(Double_t x)
Definition: TMath.h:594
constexpr Double_t LogE()
Definition: TMath.h:86
Double_t ASin(Double_t)
Definition: TMath.h:566
Double_t BesselK1(Double_t x)
Compute the modified Bessel function K_1(x) for positive real x.
Definition: TMath.cxx:1498
Double_t LnGamma(Double_t z)
Computation of ln[gamma(z)] for all z.
Definition: TMath.cxx:489
constexpr Double_t Gcgs()
Definition: TMath.h:110
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.
Definition: TMath.cxx:2384
constexpr Double_t Sigma()
Definition: TMath.h:192
Double_t SinH(Double_t)
Definition: TMath.h:557
constexpr Double_t GhbarC()
Definition: TMath.h:120
Definition: first.py:1
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 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: TMath.cxx:1081
def normal(shape, name=None)
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
double exp(double)
Double_t BetaIncomplete(Double_t x, Double_t a, Double_t b)
Calculates the incomplete Beta-function.
Definition: TMath.cxx:2070
Long64_t LocMin(Long64_t n, const T *a)
Definition: TMath.h:844
static T Max()
Definition: TMath.h:809
float * q
Definition: THbookFile.cxx:87
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 us...
Definition: TMath.cxx:2131
constexpr Double_t R()
Definition: TMath.h:213
static T Epsilon()
Definition: TMath.h:815
const Bool_t kTRUE
Definition: RtypesCore.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
Int_t Nint(T x)
Definition: TMath.h:607
const Int_t n
Definition: legend1.C:16
Double_t ATanH(Double_t)
Definition: TMath.cxx:93
Double_t CosH(Double_t)
Definition: TMath.h:560
Int_t CeilNint(Double_t x)
Definition: TMath.h:597
Double_t Tan(Double_t)
Definition: TMath.h:554
constexpr Double_t PiOver2()
Definition: TMath.h:48
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Definition: TMath.h:1093
double log(double)
double ldexp(double, int)
void SortItr(Iterator first, Iterator last, IndexIterator index, Bool_t down=kTRUE)
Definition: TMath.h:1126
Double_t Vavilov(Double_t x, Double_t kappa, Double_t beta2)
Returns the value of the Vavilov density function Parameters: 1st - the point were the density functi...
Definition: TMath.cxx:2710
constexpr Double_t KUncertainty()
Definition: TMath.h:186
T MinElement(Long64_t n, const T *a)
Definition: TMath.h:830
Double_t ATan(Double_t)
Definition: TMath.h:578