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