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