Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "TMathBase.h"
16
17#include "TError.h"
18#include <algorithm>
19#include <limits>
20#include <cmath>
21#include <vector>
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
31namespace TMath {
32
33////////////////////////////////////////////////////////////////////////////////
34// Fundamental constants
35
36////////////////////////////////////////////////////////////////////////////////
37/// \f$ \pi\f$
38constexpr Double_t Pi()
39{
40 return 3.14159265358979323846;
41}
42
43////////////////////////////////////////////////////////////////////////////////
44/// \f$ 2\pi\f$
45constexpr Double_t TwoPi()
46{
47 return 2.0 * Pi();
48}
49
50////////////////////////////////////////////////////////////////////////////////
51/// \f$ \frac{\pi}{2} \f$
52constexpr Double_t PiOver2()
53{
54 return Pi() / 2.0;
55}
56
57////////////////////////////////////////////////////////////////////////////////
58/// \f$ \frac{\pi}{4} \f$
59constexpr Double_t PiOver4()
60{
61 return Pi() / 4.0;
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// \f$ \frac{1.}{\pi}\f$
66constexpr Double_t InvPi()
67{
68 return 1.0 / Pi();
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// Conversion from radian to degree: \f$ \frac{180}{\pi} \f$
73constexpr Double_t RadToDeg()
74{
75 return 180.0 / Pi();
76}
77
78////////////////////////////////////////////////////////////////////////////////
79/// Conversion from degree to radian: \f$ \frac{\pi}{180} \f$
80constexpr Double_t DegToRad()
81{
82 return Pi() / 180.0;
83}
84
85////////////////////////////////////////////////////////////////////////////////
86/// \f$ \sqrt{2} \f$
87constexpr Double_t Sqrt2()
88{
89 return 1.4142135623730950488016887242097;
90}
91
92////////////////////////////////////////////////////////////////////////////////
93/// Base of natural log: \f$ e \f$
94constexpr Double_t E()
95{
96 return 2.71828182845904523536;
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Natural log of 10 (to convert log to ln)
101constexpr Double_t Ln10()
102{
103 return 2.30258509299404568402;
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Base-10 log of e (to convert ln to log)
108constexpr Double_t LogE()
109{
110 return 0.43429448190325182765;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Velocity of light in \f$ m s^{-1} \f$
115constexpr Double_t C()
116{
117 return 2.99792458e8;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// \f$ cm s^{-1} \f$
122constexpr Double_t Ccgs()
123{
124 return 100.0 * C();
125}
126
127////////////////////////////////////////////////////////////////////////////////
128/// Speed of light uncertainty.
130{
131 return 0.0;
132}
133
134////////////////////////////////////////////////////////////////////////////////
135/// Gravitational constant in: \f$ m^{3} kg^{-1} s^{-2} \f$
136constexpr Double_t G()
137{
138 // use 2018 value from NIST (https://physics.nist.gov/cgi-bin/cuu/Value?bg|search_for=G)
139 return 6.67430e-11;
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// \f$ cm^{3} g^{-1} s^{-2} \f$
144constexpr Double_t Gcgs()
145{
146 return G() * 1000.0;
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Gravitational constant uncertainty.
152{
153 // use 2018 value from NIST
154 return 0.00015e-11;
155}
156
157////////////////////////////////////////////////////////////////////////////////
158/// \f$ \frac{G}{\hbar C} \f$ in \f$ (GeV/c^{2})^{-2} \f$
159constexpr Double_t GhbarC()
160{
161 // use new value from NIST (2018)
162 return 6.70883e-39;
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// \f$ \frac{G}{\hbar C} \f$ uncertainty.
168{
169 // use new value from NIST (2018)
170 return 0.00015e-39;
171}
172
173////////////////////////////////////////////////////////////////////////////////
174/// Standard acceleration of gravity in \f$ m s^{-2} \f$
175constexpr Double_t Gn()
176{
177 return 9.80665;
178}
179
180////////////////////////////////////////////////////////////////////////////////
181/// Standard acceleration of gravity uncertainty.
183{
184 return 0.0;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188/// Planck's constant in \f$ J s \f$: \f$ h \f$
189constexpr Double_t H()
190{
191 return 6.62607015e-34;
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// \f$ erg s \f$
196constexpr Double_t Hcgs()
197{
198 return 1.0e7 * H();
199}
200
201////////////////////////////////////////////////////////////////////////////////
202/// Planck's constant uncertainty.
204{
205 // Planck constant is exact according to 2019 redefinition
206 // (https://en.wikipedia.org/wiki/2019_redefinition_of_the_SI_base_units)
207 return 0.0;
208}
209
210////////////////////////////////////////////////////////////////////////////////
211/// \f$ \hbar \f$ in \f$ J s \f$: \f$ \hbar = \frac{h}{2\pi} \f$
212constexpr Double_t Hbar()
213{
214 return 1.054571817e-34;
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// \f$ erg s \f$
219constexpr Double_t Hbarcgs()
220{
221 return 1.0e7 * Hbar();
222}
223
224////////////////////////////////////////////////////////////////////////////////
225/// \f$ \hbar \f$ uncertainty.
227{
228 // hbar is an exact constant
229 return 0.0;
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// \f$ hc \f$ in \f$ J m \f$
234constexpr Double_t HC()
235{
236 return H() * C();
237}
238
239////////////////////////////////////////////////////////////////////////////////
240/// \f$ erg cm \f$
241constexpr Double_t HCcgs()
242{
243 return Hcgs() * Ccgs();
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Boltzmann's constant in \f$ J K^{-1} \f$: \f$ k \f$
248constexpr Double_t K()
249{
250 return 1.380649e-23;
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// \f$ erg K^{-1} \f$
255constexpr Double_t Kcgs()
256{
257 return 1.0e7 * K();
258}
259
260////////////////////////////////////////////////////////////////////////////////
261/// Boltzmann's constant uncertainty.
263{
264 // constant is exact according to 2019 redefinition
265 // (https://en.wikipedia.org/wiki/2019_redefinition_of_the_SI_base_units)
266 return 0.0;
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Stefan-Boltzmann constant in \f$ W m^{-2} K^{-4}\f$: \f$ \sigma \f$
271constexpr Double_t Sigma()
272{
273 return 5.670373e-8;
274}
275
276////////////////////////////////////////////////////////////////////////////////
277/// Stefan-Boltzmann constant uncertainty.
279{
280 return 0.0;
281}
282
283////////////////////////////////////////////////////////////////////////////////
284/// Avogadro constant (Avogadro's Number) in \f$ mol^{-1} \f$
285constexpr Double_t Na()
286{
287 return 6.02214076e+23;
288}
289
290////////////////////////////////////////////////////////////////////////////////
291/// Avogadro constant (Avogadro's Number) uncertainty.
293{
294 // constant is exact according to 2019 redefinition
295 // (https://en.wikipedia.org/wiki/2019_redefinition_of_the_SI_base_units)
296 return 0.0;
297}
298
299////////////////////////////////////////////////////////////////////////////////
300/// [Universal gas constant](http://scienceworld.wolfram.com/physics/UniversalGasConstant.html)
301/// (\f$ Na K \f$) in \f$ J K^{-1} mol^{-1} \f$
302//
303constexpr Double_t R()
304{
305 return K() * Na();
306}
307
308////////////////////////////////////////////////////////////////////////////////
309/// Universal gas constant uncertainty.
311{
312 return R() * ((KUncertainty() / K()) + (NaUncertainty() / Na()));
313}
314
315////////////////////////////////////////////////////////////////////////////////
316/// [Molecular weight of dry air 1976 US Standard Atmosphere](http://atmos.nmsu.edu/jsdap/encyclopediawork.html)
317/// in \f$ kg kmol^{-1} \f$ or \f$ gm mol^{-1} \f$
318constexpr Double_t MWair()
319{
320 return 28.9644;
321}
322
323////////////////////////////////////////////////////////////////////////////////
324/// [Dry Air Gas Constant (R / MWair)](http://atmos.nmsu.edu/education_and_outreach/encyclopedia/gas_constant.htm)
325/// in \f$ J kg^{-1} K^{-1} \f$
326constexpr Double_t Rgair()
327{
328 return (1000.0 * R()) / MWair();
329}
330
331////////////////////////////////////////////////////////////////////////////////
332/// Euler-Mascheroni Constant.
334{
335 return 0.577215664901532860606512090082402431042;
336}
337
338////////////////////////////////////////////////////////////////////////////////
339/// Elementary charge in \f$ C \f$ .
340constexpr Double_t Qe()
341{
342 return 1.602176634e-19;
343}
344
345////////////////////////////////////////////////////////////////////////////////
346/// Elementary charge uncertainty.
348{
349 // constant is exact according to 2019 redefinition
350 // (https://en.wikipedia.org/wiki/2019_redefinition_of_the_SI_base_units)
351 return 0.0;
352}
353
354////////////////////////////////////////////////////////////////////////////////
355// Mathematical Functions
356
357////////////////////////////////////////////////////////////////////////////////
358// Trigonometrical Functions
359
360inline Double_t Sin(Double_t);
361inline Double_t Cos(Double_t);
362inline Double_t Tan(Double_t);
363inline Double_t SinH(Double_t);
364inline Double_t CosH(Double_t);
365inline Double_t TanH(Double_t);
366inline Double_t ASin(Double_t);
367inline Double_t ACos(Double_t);
368inline Double_t ATan(Double_t);
374
375////////////////////////////////////////////////////////////////////////////////
376// Elementary Functions
377
378inline Double_t Ceil(Double_t x);
379inline Int_t CeilNint(Double_t x);
380inline Double_t Floor(Double_t x);
381inline Int_t FloorNint(Double_t x);
382template <typename T>
383inline Int_t Nint(T x);
384
385inline Double_t Sq(Double_t x);
386inline Double_t Sqrt(Double_t x);
387inline Double_t Exp(Double_t x);
388inline Double_t Ldexp(Double_t x, Int_t exp);
394inline Double_t Power(Double_t x, Int_t y);
395inline Double_t Log(Double_t x);
397inline Double_t Log10(Double_t x);
398inline Int_t Finite(Double_t x);
399inline Int_t Finite(Float_t x);
400inline Bool_t IsNaN(Double_t x);
401inline Bool_t IsNaN(Float_t x);
402
403inline Double_t QuietNaN();
404inline Double_t SignalingNaN();
405inline Double_t Infinity();
406
407template <typename T>
408struct Limits {
409 inline static T Min();
410 inline static T Max();
411 inline static T Epsilon();
412 };
413
414 // Some integer math
416
417 /// Comparing floating points.
418 /// Returns `kTRUE` if the absolute difference between `af` and `bf` is less than `epsilon`.
420 return af == bf || // shortcut for exact equality (necessary because otherwise (inf != inf))
421 TMath::Abs(af-bf) < epsilon ||
422 TMath::Abs(af - bf) < Limits<Double_t>::Min(); // handle 0 < 0 case
423
424 }
425 /// Comparing floating points.
426 /// Returns `kTRUE` if the relative difference between `af` and `bf` is less than `relPrec`.
428 return af == bf || // shortcut for exact equality (necessary because otherwise (inf != inf))
429 TMath::Abs(af - bf) <= 0.5 * relPrec * (TMath::Abs(af) + TMath::Abs(bf)) ||
430 TMath::Abs(af - bf) < Limits<Double_t>::Min(); // handle denormals
431 }
432
433 /////////////////////////////////////////////////////////////////////////////
434 // Array Algorithms
435
436 // Min, Max of an array
437 template <typename T> T MinElement(Long64_t n, const T *a);
438 template <typename T> T MaxElement(Long64_t n, const T *a);
439
440 // Locate Min, Max element number in an array
441 template <typename T> Long64_t LocMin(Long64_t n, const T *a);
442 template <typename Iterator> Iterator LocMin(Iterator first, Iterator last);
443 template <typename T> Long64_t LocMax(Long64_t n, const T *a);
444 template <typename Iterator> Iterator LocMax(Iterator first, Iterator last);
445
446 // Derivatives of an array
447 template <typename T> T *Gradient(Long64_t n, T *f, double h = 1);
448 template <typename T> T *Laplacian(Long64_t n, T *f, double h = 1);
449
450 // Hashing
451 ULong_t Hash(const void *txt, Int_t ntxt);
452 ULong_t Hash(const char *str);
453
456
457 Bool_t Permute(Int_t n, Int_t *a); // Find permutations
458
459 /////////////////////////////////////////////////////////////////////////////
460 // Geometrical Functions
461
462 //Sample quantiles
464 Bool_t isSorted=kTRUE, Int_t *index = nullptr, Int_t type=7);
465
466 // IsInside
467 template <typename T> Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y);
468
469 // Calculate the Cross Product of two vectors
470 template <typename T> T *Cross(const T v1[3],const T v2[3], T out[3]);
471
472 Float_t Normalize(Float_t v[3]); // Normalize a vector
473 Double_t Normalize(Double_t v[3]); // Normalize a vector
474
475 //Calculate the Normalized Cross Product of two vectors
476 template <typename T> inline T NormCross(const T v1[3],const T v2[3],T out[3]);
477
478 // Calculate a normal vector of a plane
479 template <typename T> T *Normal2Plane(const T v1[3],const T v2[3],const T v3[3], T normal[3]);
480
481 /////////////////////////////////////////////////////////////////////////////
482 // Polynomial Functions
483
484 Bool_t RootsCubic(const Double_t coef[4],Double_t &a, Double_t &b, Double_t &c);
485
486 /////////////////////////////////////////////////////////////////////////////
487 // Statistic Functions
488
489 Double_t Binomial(Int_t n,Int_t k); // Calculate the binomial coefficient n over k
515
516 /////////////////////////////////////////////////////////////////////////////
517 // Statistics over arrays
518
519 // Mean, Geometric Mean, Median, RMS(sigma), ModeHalfSample
520
521 template <typename T> Double_t Mean(Long64_t n, const T *a, const Double_t *w=nullptr);
522 template <typename Iterator> Double_t Mean(Iterator first, Iterator last);
523 template <typename Iterator, typename WeightIterator> Double_t Mean(Iterator first, Iterator last, WeightIterator wfirst);
524
525 template <typename T> Double_t GeomMean(Long64_t n, const T *a);
526 template <typename Iterator> Double_t GeomMean(Iterator first, Iterator last);
527
528 template <typename T> Double_t RMS(Long64_t n, const T *a, const Double_t *w=nullptr);
529 template <typename Iterator> Double_t RMS(Iterator first, Iterator last);
530 template <typename Iterator, typename WeightIterator> Double_t RMS(Iterator first, Iterator last, WeightIterator wfirst);
531
532 template <typename T> Double_t StdDev(Long64_t n, const T *a, const Double_t * w = nullptr) { return RMS<T>(n,a,w); } /// Same as RMS
533 template <typename Iterator> Double_t StdDev(Iterator first, Iterator last) { return RMS<Iterator>(first,last); } /// Same as RMS
534 template <typename Iterator, typename WeightIterator> Double_t StdDev(Iterator first, Iterator last, WeightIterator wfirst) { return RMS<Iterator,WeightIterator>(first,last,wfirst); } /// Same as RMS
535
536 template <typename T> Double_t Median(Long64_t n, const T *a, const Double_t *w=nullptr, Long64_t *work=nullptr);
537 template <typename T> Double_t ModeHalfSample(Long64_t n, const T *a, const Double_t *w = nullptr);
538
539 //k-th order statistic
540 template <class Element, typename Size> Element KOrdStat(Size n, const Element *a, Size k, Size *work = 0);
541
542 /////////////////////////////////////////////////////////////////////////////
543 // Special Functions
544
550
551 // Bessel functions
552 Double_t BesselI(Int_t n,Double_t x); /// Integer order modified Bessel function I_n(x)
553 Double_t BesselK(Int_t n,Double_t x); /// Integer order modified Bessel function K_n(x)
554 Double_t BesselI0(Double_t x); /// Modified Bessel function I_0(x)
555 Double_t BesselK0(Double_t x); /// Modified Bessel function K_0(x)
556 Double_t BesselI1(Double_t x); /// Modified Bessel function I_1(x)
557 Double_t BesselK1(Double_t x); /// Modified Bessel function K_1(x)
558 Double_t BesselJ0(Double_t x); /// Bessel function J0(x) for any real x
559 Double_t BesselJ1(Double_t x); /// Bessel function J1(x) for any real x
560 Double_t BesselY0(Double_t x); /// Bessel function Y0(x) for positive x
561 Double_t BesselY1(Double_t x); /// Bessel function Y1(x) for positive x
562 Double_t StruveH0(Double_t x); /// Struve functions of order 0
563 Double_t StruveH1(Double_t x); /// Struve functions of order 1
564 Double_t StruveL0(Double_t x); /// Modified Struve functions of order 0
565 Double_t StruveL1(Double_t x); /// Modified Struve functions of order 1
566
575 Double_t GammaDist(Double_t x, Double_t gamma, Double_t mu=0, Double_t beta=1);
577}
578
579////////////////////////////////////////////////////////////////////////////////
580// Trig and other functions
581
582#include <cfloat>
583#include <cmath>
584
585#if defined(R__WIN32)
586# ifndef finite
587# define finite _finite
588# endif
589#endif
590
591////////////////////////////////////////////////////////////////////////////////
592/// Returns the sine of an angle of `x` radians.
593
595 { return sin(x); }
596
597////////////////////////////////////////////////////////////////////////////////
598/// Returns the cosine of an angle of `x` radians.
599
601 { return cos(x); }
602
603////////////////////////////////////////////////////////////////////////////////
604/// Returns the tangent of an angle of `x` radians.
605
607 { return tan(x); }
608
609////////////////////////////////////////////////////////////////////////////////
610/// Returns the hyperbolic sine of `x.
611
613 { return sinh(x); }
614
615////////////////////////////////////////////////////////////////////////////////
616/// Returns the hyperbolic cosine of `x`.
617
619 { return cosh(x); }
620
621////////////////////////////////////////////////////////////////////////////////
622/// Returns the hyperbolic tangent of `x`.
623
625 { return tanh(x); }
626
627////////////////////////////////////////////////////////////////////////////////
628/// Returns the principal value of the arc sine of `x`, expressed in radians.
629
631 {
632 return asin(x);
633 }
634
635////////////////////////////////////////////////////////////////////////////////
636/// Returns the principal value of the arc cosine of `x`, expressed in radians.
637
639 {
640 return acos(x);
641 }
642
643////////////////////////////////////////////////////////////////////////////////
644/// Returns the principal value of the arc tangent of `x`, expressed in radians.
645
647 { return atan(x); }
648
649////////////////////////////////////////////////////////////////////////////////
650/// Returns the principal value of the arc tangent of `y/x`, expressed in radians.
651
653 { if (x != 0) return atan2(y, x);
654 if (y == 0) return 0;
655 if (y > 0) return Pi()/2;
656 else return -Pi()/2;
657 }
658
659////////////////////////////////////////////////////////////////////////////////
660/// Returns `x*x`.
661
663 { return x*x; }
664
665////////////////////////////////////////////////////////////////////////////////
666/// Returns the square root of x.
667
669 { return sqrt(x); }
670
671////////////////////////////////////////////////////////////////////////////////
672/// Rounds `x` upward, returning the smallest integral value that is not less than `x`.
673
675 { return ceil(x); }
676
677////////////////////////////////////////////////////////////////////////////////
678/// Returns the nearest integer of `TMath::Ceil(x)`.
679
681 { return TMath::Nint(ceil(x)); }
682
683////////////////////////////////////////////////////////////////////////////////
684/// Rounds `x` downward, returning the largest integral value that is not greater than `x`.
685
687 { return floor(x); }
688
689////////////////////////////////////////////////////////////////////////////////
690/// Returns the nearest integer of `TMath::Floor(x)`.
691
693 { return TMath::Nint(floor(x)); }
694
695////////////////////////////////////////////////////////////////////////////////
696/// Round to nearest integer. Rounds half integers to the nearest even integer.
697
698template<typename T>
700{
701 int i;
702 if (x >= 0) {
703 i = int(x + 0.5);
704 if ( i & 1 && x + 0.5 == T(i) ) i--;
705 } else {
706 i = int(x - 0.5);
707 if ( i & 1 && x - 0.5 == T(i) ) i++;
708 }
709 return i;
710}
711
712////////////////////////////////////////////////////////////////////////////////
713/// Returns the base-e exponential function of x, which is e raised to the power `x`.
714
716 { return exp(x); }
717
718////////////////////////////////////////////////////////////////////////////////
719/// Returns the result of multiplying `x` (the significant) by 2 raised to the power of `exp` (the exponent).
720
722 { return ldexp(x, exp); }
723
724////////////////////////////////////////////////////////////////////////////////
725/// Returns `x` raised to the power `y`.
726
728 { return std::pow(x,y); }
729
730////////////////////////////////////////////////////////////////////////////////
731/// Returns `x` raised to the power `y`.
732
734 { return std::pow(x,(LongDouble_t)y); }
735
736////////////////////////////////////////////////////////////////////////////////
737/// Returns `x` raised to the power `y`.
738
740 { return std::pow(x,y); }
741
742////////////////////////////////////////////////////////////////////////////////
743/// Returns `x` raised to the power `y`.
744
746 { return pow(x, y); }
747
748////////////////////////////////////////////////////////////////////////////////
749/// Returns `x` raised to the power `y`.
750
752#ifdef R__ANSISTREAM
753 return std::pow(x, y);
754#else
755 return pow(x, (Double_t) y);
756#endif
757}
758
759////////////////////////////////////////////////////////////////////////////////
760/// Returns the natural logarithm of `x`.
761
763 { return log(x); }
764
765////////////////////////////////////////////////////////////////////////////////
766/// Returns the common (base-10) logarithm of `x`.
767
769 { return log10(x); }
770
771////////////////////////////////////////////////////////////////////////////////
772/// Check if it is finite with a mask in order to be consistent in presence of
773/// fast math.
774/// Inspired from the CMSSW FWCore/Utilities package
775
777#if defined(R__FAST_MATH)
778
779{
780 const unsigned long long mask = 0x7FF0000000000000LL;
781 union { unsigned long long l; double d;} v;
782 v.d =x;
783 return (v.l&mask)!=mask;
784}
785#else
786# if defined(R__HPUX11)
787 { return isfinite(x); }
788# elif defined(R__MACOSX)
789# ifdef isfinite
790 // from math.h
791 { return isfinite(x); }
792# else
793 // from cmath
794 { return std::isfinite(x); }
795# endif
796# else
797 { return finite(x); }
798# endif
799#endif
800
801////////////////////////////////////////////////////////////////////////////////
802/// Check if it is finite with a mask in order to be consistent in presence of
803/// fast math.
804/// Inspired from the CMSSW FWCore/Utilities package
805
807#if defined(R__FAST_MATH)
808
809{
810 const unsigned int mask = 0x7f800000;
811 union { unsigned int l; float d;} v;
812 v.d =x;
813 return (v.l&mask)!=mask;
814}
815#else
816{ return std::isfinite(x); }
817#endif
818
819// This namespace provides all the routines necessary for checking if a number
820// is a NaN also in presence of optimisations affecting the behaviour of the
821// floating point calculations.
822// Inspired from the CMSSW FWCore/Utilities package
823
824#if defined (R__FAST_MATH)
825namespace ROOT {
826namespace Internal {
827namespace Math {
828// abridged from GNU libc 2.6.1 - in detail from
829// math/math_private.h
830// sysdeps/ieee754/ldbl-96/math_ldbl.h
831
832// part of this file:
833 /*
834 * ====================================================
835 * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
836 *
837 * Developed at SunPro, a Sun Microsystems, Inc. business.
838 * Permission to use, copy, modify, and distribute this
839 * software is freely granted, provided that this notice
840 * is preserved.
841 * ====================================================
842 */
843
844 // A union which permits us to convert between a double and two 32 bit ints.
845 typedef union {
847 struct {
848 UInt_t lsw;
849 UInt_t msw;
850 } parts;
852
853#define EXTRACT_WORDS(ix0,ix1,d) \
854 do { \
855 ieee_double_shape_type ew_u; \
856 ew_u.value = (d); \
857 (ix0) = ew_u.parts.msw; \
858 (ix1) = ew_u.parts.lsw; \
859 } while (0)
860
861 inline Bool_t IsNaN(Double_t x)
862 {
863 UInt_t hx, lx;
864
865 EXTRACT_WORDS(hx, lx, x);
866
867 lx |= hx & 0xfffff;
868 hx &= 0x7ff00000;
869 return (hx == 0x7ff00000) && (lx != 0);
870 }
871
872 typedef union {
874 UInt_t word;
876
877#define GET_FLOAT_WORD(i,d) \
878 do { \
879 ieee_float_shape_type gf_u; \
880 gf_u.value = (d); \
881 (i) = gf_u.word; \
882 } while (0)
883
884 inline Bool_t IsNaN(Float_t x)
885 {
886 UInt_t wx;
888 wx &= 0x7fffffff;
889 return (Bool_t)(wx > 0x7f800000);
890 }
891} } } // end NS ROOT::Internal::Math
892#endif // End R__FAST_MATH
893
894#if defined(R__FAST_MATH)
895 inline Bool_t TMath::IsNaN(Double_t x) { return ROOT::Internal::Math::IsNaN(x); }
896 inline Bool_t TMath::IsNaN(Float_t x) { return ROOT::Internal::Math::IsNaN(x); }
897#else
898 inline Bool_t TMath::IsNaN(Double_t x) { return std::isnan(x); }
899 inline Bool_t TMath::IsNaN(Float_t x) { return std::isnan(x); }
900#endif
901
902////////////////////////////////////////////////////////////////////////////////
903// Wrapper to numeric_limits
904
905////////////////////////////////////////////////////////////////////////////////
906/// Returns a quiet NaN as [defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Quiet_NaN).
907
909
910 return std::numeric_limits<Double_t>::quiet_NaN();
911}
912
913////////////////////////////////////////////////////////////////////////////////
914/// Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
915
917 return std::numeric_limits<Double_t>::signaling_NaN();
918}
919
920////////////////////////////////////////////////////////////////////////////////
921/// Returns an infinity as defined by the IEEE standard.
922
924 return std::numeric_limits<Double_t>::infinity();
925}
926
927////////////////////////////////////////////////////////////////////////////////
928/// Returns maximum representation for type T.
929
930template<typename T>
932 return (std::numeric_limits<T>::min)(); //N.B. use this signature to avoid class with macro min() on Windows
933}
934
935////////////////////////////////////////////////////////////////////////////////
936/// Returns minimum double representation.
937
938template<typename T>
940 return (std::numeric_limits<T>::max)(); //N.B. use this signature to avoid class with macro max() on Windows
941}
942
943////////////////////////////////////////////////////////////////////////////////
944/// Returns minimum double representation.
945
946template<typename T>
948 return std::numeric_limits<T>::epsilon();
949}
950
951////////////////////////////////////////////////////////////////////////////////
952// Advanced.
953
954////////////////////////////////////////////////////////////////////////////////
955/// Calculates the Normalized Cross Product of two vectors.
956
957template <typename T> inline T TMath::NormCross(const T v1[3],const T v2[3],T out[3])
958{
959 return Normalize(Cross(v1,v2,out));
960}
961
962////////////////////////////////////////////////////////////////////////////////
963/// Returns minimum of array a of length n.
964
965template <typename T>
967 return *std::min_element(a,a+n);
968}
969
970////////////////////////////////////////////////////////////////////////////////
971/// Returns maximum of array a of length n.
972
973template <typename T>
975 return *std::max_element(a,a+n);
976}
977
978////////////////////////////////////////////////////////////////////////////////
979/// Returns index of array with the minimum element.
980/// If more than one element is minimum returns first found.
981///
982/// Implement here since this one is found to be faster (mainly on 64 bit machines)
983/// than stl generic implementation.
984/// When performing the comparison, the STL implementation needs to de-reference both the array iterator
985/// and the iterator pointing to the resulting minimum location
986
987template <typename T>
989 if (n <= 0 || !a) return -1;
990 T xmin = a[0];
991 Long64_t loc = 0;
992 for (Long64_t i = 1; i < n; i++) {
993 if (xmin > a[i]) {
994 xmin = a[i];
995 loc = i;
996 }
997 }
998 return loc;
999}
1000
1001////////////////////////////////////////////////////////////////////////////////
1002/// Returns index of array with the minimum element.
1003/// If more than one element is minimum returns first found.
1004
1005template <typename Iterator>
1006Iterator TMath::LocMin(Iterator first, Iterator last) {
1007
1008 return std::min_element(first, last);
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// \brief Calculate the one-dimensional gradient of an array with length n.
1013/// The first value in the returned array is a forward difference,
1014/// the next n-2 values are central differences, and the last is a backward difference.
1015///
1016/// \note Function leads to undefined behavior if n does not match the length of f
1017/// \param n the number of points in the array
1018/// \param f the array of points.
1019/// \param h the step size. The default step size is 1.
1020/// \return an array of size n with the gradient. Returns nullptr if n < 2 or f empty. Ownership is transferred to the
1021/// caller.
1022
1023template <typename T>
1024T *TMath::Gradient(const Long64_t n, T *f, const double h)
1025{
1026 if (!f) {
1027 ::Error("TMath::Gradient", "Input parameter f is empty.");
1028 return nullptr;
1029 } else if (n < 2) {
1030 ::Error("TMath::Gradient", "Input parameter n=%lld is smaller than 2.", n);
1031 return nullptr;
1032 }
1033 Long64_t i = 1;
1034 T *result = new T[n];
1035
1036 // Forward difference
1037 result[0] = (f[1] - f[0]) / h;
1038
1039 // Central difference
1040 while (i < n - 1) {
1041 result[i] = (f[i + 1] - f[i - 1]) / (2 * h);
1042 i++;
1043 }
1044 // Backward difference
1045 result[i] = (f[i] - f[i - 1]) / h;
1046 return result;
1047}
1048
1049////////////////////////////////////////////////////////////////////////////////
1050/// \brief Calculate the Laplacian of an array with length n.
1051/// The first value in the returned array is a forward difference,
1052/// the next n-2 values are central differences, and the last is a backward difference.
1053///
1054/// \note Function leads to undefined behavior if n does not match the length of f
1055/// \param n the number of points in the array
1056/// \param f the array of points.
1057/// \param h the step size. The default step size is 1.
1058/// \return an array of size n with the laplacian. Returns nullptr if n < 4 or f empty. Ownership is transferred to the
1059/// caller.
1060
1061template <typename T>
1062T *TMath::Laplacian(const Long64_t n, T *f, const double h)
1063{
1064 if (!f) {
1065 ::Error("TMath::Laplacian", "Input parameter f is empty.");
1066 return nullptr;
1067 } else if (n < 4) {
1068 ::Error("TMath::Laplacian", "Input parameter n=%lld is smaller than 4.", n);
1069 return nullptr;
1070 }
1071 Long64_t i = 1;
1072 T *result = new T[n];
1073
1074 // Forward difference
1075 result[0] = (4 * f[2] + 2 * f[0] - 5 * f[1] - f[3]) / (4 * h * h);
1076
1077 // Central difference
1078 while (i < n - 1) {
1079 result[i] = (f[i + 1] + f[i - 1] - 2 * f[i]) / (4 * h * h);
1080 i++;
1081 }
1082 // Backward difference
1083 result[i] = (2 * f[i] - 5 * f[i - 1] + 4 * f[i - 2] - f[i - 3]) / (4 * h * h);
1084 return result;
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// Returns index of array with the maximum element.
1089/// If more than one element is maximum returns first found.
1090///
1091/// Implement here since it is faster (see comment in LocMin function)
1092
1093template <typename T>
1095 if (n <= 0 || !a) return -1;
1096 T xmax = a[0];
1097 Long64_t loc = 0;
1098 for (Long64_t i = 1; i < n; i++) {
1099 if (xmax < a[i]) {
1100 xmax = a[i];
1101 loc = i;
1102 }
1103 }
1104 return loc;
1105}
1106
1107////////////////////////////////////////////////////////////////////////////////
1108/// Returns index of array with the maximum element.
1109/// If more than one element is maximum returns first found.
1110
1111template <typename Iterator>
1112Iterator TMath::LocMax(Iterator first, Iterator last)
1113{
1114
1115 return std::max_element(first, last);
1116}
1117
1118////////////////////////////////////////////////////////////////////////////////
1119/// Returns the weighted mean of an array defined by the iterators.
1120
1121template <typename Iterator>
1122Double_t TMath::Mean(Iterator first, Iterator last)
1123{
1124 Double_t sum = 0;
1125 Double_t sumw = 0;
1126 while ( first != last )
1127 {
1128 sum += *first;
1129 sumw += 1;
1130 first++;
1131 }
1132
1133 return sum/sumw;
1134}
1135
1136////////////////////////////////////////////////////////////////////////////////
1137/// Returns the weighted mean of an array defined by the first and
1138/// last iterators. The w iterator should point to the first element
1139/// of a vector of weights of the same size as the main array.
1140
1141template <typename Iterator, typename WeightIterator>
1142Double_t TMath::Mean(Iterator first, Iterator last, WeightIterator w)
1143{
1144
1145 Double_t sum = 0;
1146 Double_t sumw = 0;
1147 int i = 0;
1148 while ( first != last ) {
1149 if ( *w < 0) {
1150 ::Error("TMath::Mean","w[%d] = %.4e < 0 ?!",i,*w);
1151 return 0;
1152 }
1153 sum += (*w) * (*first);
1154 sumw += (*w) ;
1155 ++w;
1156 ++first;
1157 ++i;
1158 }
1159 if (sumw <= 0) {
1160 ::Error("TMath::Mean","sum of weights == 0 ?!");
1161 return 0;
1162 }
1163
1164 return sum/sumw;
1165}
1166
1167////////////////////////////////////////////////////////////////////////////////
1168/// Returns the weighted mean of an array a with length n.
1169
1170template <typename T>
1172{
1173 if (w) {
1174 return TMath::Mean(a, a+n, w);
1175 } else {
1176 return TMath::Mean(a, a+n);
1177 }
1178}
1179
1180////////////////////////////////////////////////////////////////////////////////
1181/// Returns the geometric mean of an array defined by the iterators.
1182/// \f[ GeomMean = (\prod_{i=0}^{n-1} |a[i]|)^{1/n} \f]
1183
1184template <typename Iterator>
1185Double_t TMath::GeomMean(Iterator first, Iterator last)
1186{
1187 Double_t logsum = 0.;
1188 Long64_t n = 0;
1189 while ( first != last ) {
1190 if (*first == 0) return 0.;
1191 Double_t absa = (Double_t) TMath::Abs(*first);
1193 ++first;
1194 ++n;
1195 }
1196
1197 return TMath::Exp(logsum/n);
1198}
1199
1200////////////////////////////////////////////////////////////////////////////////
1201/// Returns the geometric mean of an array a of size n.
1202/// \f[ GeomMean = (\prod_{i=0}^{n-1} |a[i]|)^{1/n} \f]
1203
1204template <typename T>
1206{
1207 return TMath::GeomMean(a, a+n);
1208}
1209
1210////////////////////////////////////////////////////////////////////////////////
1211/// Returns the Standard Deviation of an array defined by the iterators.
1212/// Note that this function returns the sigma(standard deviation) and
1213/// not the root mean square of the array.
1214///
1215/// Use the two pass algorithm, which is slower (! a factor of 2) but much more
1216/// precise. Since we have a vector the 2 pass algorithm is still faster than the
1217/// Welford algorithm. (See also ROOT-5545)
1218
1219template <typename Iterator>
1220Double_t TMath::RMS(Iterator first, Iterator last)
1221{
1222
1223 Double_t n = 0;
1224
1225 Double_t tot = 0;
1226 Double_t mean = TMath::Mean(first,last);
1227 while ( first != last ) {
1228 Double_t x = Double_t(*first);
1229 tot += (x - mean)*(x - mean);
1230 ++first;
1231 ++n;
1232 }
1233 Double_t rms = (n > 1) ? TMath::Sqrt(tot/(n-1)) : 0.0;
1234 return rms;
1235}
1236
1237////////////////////////////////////////////////////////////////////////////////
1238/// Returns the weighted Standard Deviation of an array defined by the iterators.
1239/// Note that this function returns the sigma(standard deviation) and
1240/// not the root mean square of the array.
1241///
1242/// As in the unweighted case use the two pass algorithm
1243
1244template <typename Iterator, typename WeightIterator>
1245Double_t TMath::RMS(Iterator first, Iterator last, WeightIterator w)
1246{
1247 Double_t tot = 0;
1248 Double_t sumw = 0;
1249 Double_t sumw2 = 0;
1250 Double_t mean = TMath::Mean(first,last,w);
1251 while ( first != last ) {
1252 Double_t x = Double_t(*first);
1253 sumw += *w;
1254 sumw2 += (*w) * (*w);
1255 tot += (*w) * (x - mean)*(x - mean);
1256 ++first;
1257 ++w;
1258 }
1259 // use the correction neff/(neff -1) for the unbiased formula
1261 return rms;
1262}
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// Returns the Standard Deviation of an array a with length n.
1266/// Note that this function returns the sigma(standard deviation) and
1267/// not the root mean square of the array.
1268
1269template <typename T>
1271{
1272 return (w) ? TMath::RMS(a, a+n, w) : TMath::RMS(a, a+n);
1273}
1274
1275////////////////////////////////////////////////////////////////////////////////
1276/// Calculates the Cross Product of two vectors:
1277/// out = [v1 x v2]
1278
1279template <typename T> T *TMath::Cross(const T v1[3],const T v2[3], T out[3])
1280{
1281 out[0] = v1[1] * v2[2] - v1[2] * v2[1];
1282 out[1] = v1[2] * v2[0] - v1[0] * v2[2];
1283 out[2] = v1[0] * v2[1] - v1[1] * v2[0];
1284
1285 return out;
1286}
1287
1288////////////////////////////////////////////////////////////////////////////////
1289/// Calculates a normal vector of a plane.
1290///
1291/// \param[in] p1, p2,p3 3 3D points belonged the plane to define it.
1292/// \param[out] normal Pointer to 3D normal vector (normalized)
1293
1294template <typename T> T * TMath::Normal2Plane(const T p1[3],const T p2[3],const T p3[3], T normal[3])
1295{
1296 T v1[3], v2[3];
1297
1298 v1[0] = p2[0] - p1[0];
1299 v1[1] = p2[1] - p1[1];
1300 v1[2] = p2[2] - p1[2];
1301
1302 v2[0] = p3[0] - p1[0];
1303 v2[1] = p3[1] - p1[1];
1304 v2[2] = p3[2] - p1[2];
1305
1307 return normal;
1308}
1309
1310////////////////////////////////////////////////////////////////////////////////
1311/// Function which returns kTRUE if point xp,yp lies inside the
1312/// polygon defined by the np points in arrays x and y, kFALSE otherwise.
1313/// Note that the polygon may be open or closed.
1314
1315template <typename T> Bool_t TMath::IsInside(T xp, T yp, Int_t np, T *x, T *y)
1316{
1317 Int_t i, j = np-1 ;
1319
1320 for (i=0; i<np; i++) {
1321 if ((y[i]<yp && y[j]>=yp) || (y[j]<yp && y[i]>=yp)) {
1322 if (x[i]+(yp-y[i])/(y[j]-y[i])*(x[j]-x[i])<xp) {
1323 oddNodes = !oddNodes;
1324 }
1325 }
1326 j=i;
1327 }
1328
1329 return oddNodes;
1330}
1331
1332////////////////////////////////////////////////////////////////////////////////
1333/// Returns the median of the array a where each entry i has weight w[i] .
1334/// Both arrays have a length of at least n . The median is a number obtained
1335/// from the sorted array a through
1336///
1337/// median = (a[jl]+a[jh])/2. where (using also the sorted index on the array w)
1338///
1339/// sum_i=0,jl w[i] <= sumTot/2
1340/// sum_i=0,jh w[i] >= sumTot/2
1341/// sumTot = sum_i=0,n w[i]
1342///
1343/// If w=0, the algorithm defaults to the median definition where it is
1344/// a number that divides the sorted sequence into 2 halves.
1345/// When n is odd or n > 1000, the median is kth element k = (n + 1) / 2.
1346/// when n is even and n < 1000the median is a mean of the elements k = n/2 and k = n/2 + 1.
1347///
1348/// If the weights are supplied (w not 0) all weights must be >= 0
1349///
1350/// If work is supplied, it is used to store the sorting index and assumed to be
1351/// >= n . If work=0, local storage is used, either on the stack if n < kWorkMax
1352/// or on the heap for n >= kWorkMax .
1353
1354template <typename T> Double_t TMath::Median(Long64_t n, const T *a, const Double_t *w, Long64_t *work)
1355{
1356
1357 const Int_t kWorkMax = 100;
1358
1359 if (n <= 0 || !a) return 0;
1362 Long64_t *ind;
1363 Long64_t workLocal[kWorkMax];
1364
1365 if (work) {
1366 ind = work;
1367 } else {
1368 ind = workLocal;
1369 if (n > kWorkMax) {
1371 ind = new Long64_t[n];
1372 }
1373 }
1374
1375 if (w) {
1376 Double_t sumTot2 = 0;
1377 for (Int_t j = 0; j < n; j++) {
1378 if (w[j] < 0) {
1379 ::Error("TMath::Median","w[%d] = %.4e < 0 ?!",j,w[j]);
1380 if (isAllocated) delete [] ind;
1381 return 0;
1382 }
1383 sumTot2 += w[j];
1384 }
1385
1386 sumTot2 /= 2.;
1387
1388 Sort(n, a, ind, kFALSE);
1389
1390 Double_t sum = 0.;
1391 Int_t jl;
1392 for (jl = 0; jl < n; jl++) {
1393 sum += w[ind[jl]];
1394 if (sum >= sumTot2) break;
1395 }
1396
1397 Int_t jh;
1398 sum = 2.*sumTot2;
1399 for (jh = n-1; jh >= 0; jh--) {
1400 sum -= w[ind[jh]];
1401 if (sum <= sumTot2) break;
1402 }
1403
1404 median = 0.5*(a[ind[jl]]+a[ind[jh]]);
1405
1406 } else {
1407
1408 if (n%2 == 1)
1409 median = KOrdStat(n, a,n/2, ind);
1410 else {
1411 median = 0.5*(KOrdStat(n, a, n/2 -1, ind)+KOrdStat(n, a, n/2, ind));
1412 }
1413 }
1414
1415 if (isAllocated)
1416 delete [] ind;
1417 return median;
1418}
1419
1420////////////////////////////////////////////////////////////////////////////////
1421/// Returns the half-sample mode of the array a where each entry i has weight w[i].
1422/// Both arrays must have a length of n. The mode is a number obtained
1423/// according to the paper http://arxiv.org/ftp/math/papers/0505/0505419.pdf (page 19).
1424/// The algorithm searches for the smallest range of sorted a values that
1425/// contains ~n/2 values. Then it re-applies the same algorithm to
1426/// that range until the range has <= 3 elements in it.
1427///
1428/// \note See David R. Bickel: "On a Fast, Robust Estimator of the Mode"
1429/// http://arxiv.org/ftp/math/papers/0505/0505419.pdf (page 19)
1430/// Initial implementation suggestion by Jean-Fran\c{c}ois Caron, 2014 on JIRA-6849
1431///
1432/// \tparam T data type of the values array `a`
1433/// \param n Number of elements in array `a`
1434/// \param a Pointer to array of values (owned by user, preallocated), must have length of n, and none of them should be NaN. The array must not be necessarily sorted nor contain unique values.
1435/// \param w Pointer to array of weights (owned by user, preallocated), must have length of n. If the weights are supplied (w not nullptr) all weights must be >= 0 to properly work.
1436/// \return NaN if n <= 0, a[0] if n == 1, TMath::Mean(n,a,w) if n == 2,
1437/// the (weighted) mean of the closest two if n == 3; for the rest of the cases,
1438/// if all values of `a` are unique or have equal weights (or `w` is nullptr), then the Bickel half-sample-mode algorithm is applied; otherwise the
1439/// value associated with the highest sum of weights (if `w` is not nullptr)
1440/// is returned.
1441/// \note If in a given iteration step, there are several equally small subranges, then the first occurrence is chosen.
1442template <typename T> Double_t TMath::ModeHalfSample(Long64_t n, const T *a, const Double_t *w)
1443{
1444 if (n <= 0)
1445 return TMath::QuietNaN();
1446 else if (n == 1)
1447 return a[0];
1448 else if (n == 2)
1449 return TMath::Mean(n, a, w);
1450
1451 if (w && std::adjacent_find(w, w + n, std::not_equal_to<>()) == w + n) // If all weights are equal, ignore those
1452 w = nullptr;
1453
1454 // Sort the array and remove duplicates (if w != nullptr)
1455 std::vector<T> values;
1456 std::vector<Double_t> weights;
1457 values.reserve(n);
1458 weights.reserve(n);
1459 for (Long64_t i = 0; i < n; ++i) {
1460 const T value = a[i];
1461 const Double_t weight = w ? w[i] : 1;
1462 const auto vstart = values.begin();
1463 const auto vstop = values.end();
1464 const auto viter = std::lower_bound(vstart, vstop, value);
1465 const auto vidx = std::distance(vstart, viter);
1466 if (viter != vstop && w) {
1467 weights[vidx] += weight;
1468 } else {
1469 values.insert(viter, value);
1470 weights.insert(weights.begin() + vidx, weight);
1471 }
1472 }
1473
1474 const size_t sn = values.size();
1475 if (sn <= 0)
1476 return TMath::QuietNaN();
1477 else if (sn == 1)
1478 return values[0];
1479 else if (sn == 2)
1480 return TMath::Mean(sn, values.data(), weights.data());
1481
1482 if (sn == 3) {
1483 const double d1_0 = values[1] - values[0];
1484 const double d2_1 = values[2] - values[1];
1485 if (d2_1 < d1_0)
1486 return TMath::Mean(2, values.data() + 1, weights.data() + 1);
1487 else if (d2_1 > d1_0)
1488 return TMath::Mean(2, values.data(), weights.data());
1489 else
1490 return values[1];
1491 }
1492
1493 const auto wstart = weights.begin();
1494 const auto wstop = weights.end();
1495 if (std::adjacent_find(wstart, wstop, std::not_equal_to<>()) != wstop) {
1496 // Get highest (sum of) weight element
1497 const auto wmaxiter = std::max_element(wstart, wstop);
1498 const auto wmaxidx = std::distance(wstart, wmaxiter);
1499 return values[wmaxidx];
1500 }
1501
1502 // All elements are unique and have equal weights
1503
1504 // Initialize search
1505 n = sn;
1506 size_t jMin = 0;
1507
1508 // Do recursive calls dividing each time the interval by two
1509 while (n > 3) {
1510 Double_t min_v_range = values[jMin + n - 1] - values[jMin];
1511 const size_t N = std::ceil(n * 0.5);
1512 const size_t start = jMin;
1513 const size_t stop = start + n - N + 1; // +1 since we use < and not <=
1514 // Find sequentally what v_range is smallest by sliding the half-window
1515 for (size_t i = start; i < stop; i++) {
1516 Double_t range = values[i + N - 1] - values[i];
1517 if (range < min_v_range) {
1519 jMin = i;
1520 }
1521 }
1522 // assert(min_v_range == values[N - 1 + start] - values[start]);
1523 n = N;
1524 }
1525 if (n == 3) {
1526 const double d1_0 = values[jMin + 1] - values[jMin + 0];
1527 const double d2_1 = values[jMin + 2] - values[jMin + 1];
1528 if (d2_1 < d1_0)
1529 return (values[jMin + 1] + values[jMin + 2]) * 0.5;
1530 else if (d2_1 > d1_0)
1531 return (values[jMin] + values[jMin + 1]) * 0.5;
1532 else
1533 return values[jMin + 1];
1534 } else if (n == 2) {
1535 return (values[jMin] + values[jMin + 1]) * 0.5;
1536 } else {
1537 Error("ModeHalfSample", "Error in recursive algorithm, returning NaN"); // this should not happen
1538 return TMath::QuietNaN();
1539 }
1540}
1541
1542////////////////////////////////////////////////////////////////////////////////
1543/// Returns k_th order statistic of the array a of size n
1544/// (k_th smallest element out of n elements).
1545///
1546/// C-convention is used for array indexing, so if you want
1547/// the second smallest element, call KOrdStat(n, a, 1).
1548///
1549/// If work is supplied, it is used to store the sorting index and
1550/// assumed to be >= n. If work=0, local storage is used, either on
1551/// the stack if n < kWorkMax or on the heap for n >= kWorkMax.
1552/// Note that the work index array will not contain the sorted indices but
1553/// all indices of the smaller element in arbitrary order in work[0,...,k-1] and
1554/// all indices of the larger element in arbitrary order in work[k+1,..,n-1]
1555/// work[k] will contain instead the index of the returned element.
1556///
1557/// Taken from "Numerical Recipes in C++" without the index array
1558/// implemented by Anna Khreshuk.
1559///
1560/// See also the declarations at the top of this file
1561
1562template <class Element, typename Size>
1563Element TMath::KOrdStat(Size n, const Element *a, Size k, Size *work)
1564{
1565
1566 const Int_t kWorkMax = 100;
1567
1568 typedef Size Index;
1569
1571 Size i, ir, j, l, mid;
1572 Index arr;
1573 Index *ind;
1574 Index workLocal[kWorkMax];
1575 Index temp;
1576
1577 if (work) {
1578 ind = work;
1579 } else {
1580 ind = workLocal;
1581 if (n > kWorkMax) {
1583 ind = new Index[n];
1584 }
1585 }
1586
1587 for (Size ii=0; ii<n; ii++) {
1588 ind[ii]=ii;
1589 }
1590 Size rk = k;
1591 l=0;
1592 ir = n-1;
1593 for(;;) {
1594 if (ir<=l+1) { //active partition contains 1 or 2 elements
1595 if (ir == l+1 && a[ind[ir]]<a[ind[l]])
1596 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1597 Element tmp = a[ind[rk]];
1598 if (isAllocated)
1599 delete [] ind;
1600 return tmp;
1601 } else {
1602 mid = (l+ir) >> 1; //choose median of left, center and right
1603 {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}//elements as partitioning element arr.
1604 if (a[ind[l]]>a[ind[ir]]) //also rearrange so that a[l]<=a[l+1]
1605 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1606
1607 if (a[ind[l+1]]>a[ind[ir]])
1608 {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}
1609
1610 if (a[ind[l]]>a[ind[l+1]])
1611 {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}
1612
1613 i=l+1; //initialize pointers for partitioning
1614 j=ir;
1615 arr = ind[l+1];
1616 for (;;){
1617 do i++; while (a[ind[i]]<a[arr]);
1618 do j--; while (a[ind[j]]>a[arr]);
1619 if (j<i) break; //pointers crossed, partitioning complete
1620 {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
1621 }
1622 ind[l+1]=ind[j];
1623 ind[j]=arr;
1624 if (j>=rk) ir = j-1; //keep active the partition that
1625 if (j<=rk) l=i; //contains the k_th element
1626 }
1627 }
1628}
1629
1630#endif
#define IsNaN(a)
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
unsigned long ULong_t
Unsigned long integer 4 bytes (unsigned long). Size depends on architecture.
Definition RtypesCore.h:69
long Long_t
Signed long integer 4 bytes (long). Size depends on architecture.
Definition RtypesCore.h:68
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int)
Definition RtypesCore.h:60
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
long double LongDouble_t
Long Double (not portable)
Definition RtypesCore.h:75
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:208
#define N
TGLVector3 Cross(const TGLVector3 &v1, const TGLVector3 &v2)
Definition TGLUtil.h:323
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t mask
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
float xmin
float * q
float xmax
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Namespace for new Math classes and functions.
Namespace for new ROOT classes and functions.
TMath.
Definition TMathBase.h:35
Double_t FDistI(Double_t F, Double_t N, Double_t M)
Calculates the cumulative distribution function of F-distribution (see ROOT::Math::fdistribution_cdf)...
Definition TMath.cxx:2300
Double_t GeomMean(Long64_t n, const T *a)
Returns the geometric mean of an array a of size n.
Definition TMath.h:1205
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:2440
T * Laplacian(Long64_t n, T *f, double h=1)
Calculate the Laplacian of an array with length n.
Definition TMath.h:1062
constexpr Double_t G()
Gravitational constant in: .
Definition TMath.h:136
Double_t CosH(Double_t)
Returns the hyperbolic cosine of x.
Definition TMath.h:618
T * Normal2Plane(const T v1[3], const T v2[3], const T v3[3], T normal[3])
Calculates a normal vector of a plane.
Definition TMath.h:1294
Double_t DiLog(Double_t x)
Modified Struve functions of order 1.
Definition TMath.cxx:116
Double_t BetaDist(Double_t x, Double_t p, Double_t q)
Computes the probability density function of the Beta distribution (the cumulative distribution funct...
Definition TMath.cxx:2071
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:638
Double_t BesselI(Int_t n, Double_t x)
Computes the Integer Order Modified Bessel function I_n(x) for n=0,1,2,... and any real x.
Definition TMath.cxx:1590
Element KOrdStat(Size n, const Element *a, Size k, Size *work=0)
Returns k_th order statistic of the array a of size n (k_th smallest element out of n elements).
Definition TMath.h:1563
T * Gradient(Long64_t n, T *f, double h=1)
Calculate the one-dimensional gradient of an array with length n.
Definition TMath.h:1024
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculates a gaussian function with mean and sigma.
Definition TMath.cxx:471
Bool_t IsNaN(Double_t x)
Definition TMath.h:898
constexpr Double_t GUncertainty()
Gravitational constant uncertainty.
Definition TMath.h:151
constexpr Double_t C()
Velocity of light in .
Definition TMath.h:115
Double_t Factorial(Int_t i)
Computes factorial(n).
Definition TMath.cxx:252
Double_t RMS(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the Standard Deviation of an array a with length n.
Definition TMath.h:1270
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:805
constexpr Double_t GhbarCUncertainty()
uncertainty.
Definition TMath.h:167
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:988
constexpr Double_t Ccgs()
Definition TMath.h:122
constexpr Double_t SigmaUncertainty()
Stefan-Boltzmann constant uncertainty.
Definition TMath.h:278
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:699
Double_t ModeHalfSample(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the half-sample mode of the array a where each entry i has weight w[i].
Definition TMath.h:1442
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:2144
Double_t Vavilov(Double_t x, Double_t kappa, Double_t beta2)
Returns the value of the Vavilov probability density function.
Definition TMath.cxx:2780
Double_t Binomial(Int_t n, Int_t k)
Calculates the binomial coefficient n over k.
Definition TMath.cxx:2112
Float_t Normalize(Float_t v[3])
Normalize a vector v in place.
Definition TMath.cxx:518
constexpr Double_t NaUncertainty()
Avogadro constant (Avogadro's Number) uncertainty.
Definition TMath.h:292
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:637
Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y)
Function which returns kTRUE if point xp,yp lies inside the polygon defined by the np points in array...
Definition TMath.h:1315
Double_t ASin(Double_t)
Returns the principal value of the arc sine of x, expressed in radians.
Definition TMath.h:630
Double_t Log2(Double_t x)
Returns the binary (base-2) logarithm of x.
Definition TMath.cxx:107
Double_t BesselK1(Double_t x)
Modified Bessel function I_1(x)
Definition TMath.cxx:1529
Double_t Median(Long64_t n, const T *a, const Double_t *w=nullptr, Long64_t *work=nullptr)
Same as RMS.
Definition TMath.h:1354
void BubbleHigh(Int_t Narr, Double_t *arr1, Int_t *arr2)
Bubble sort variant to obtain the order of an array's elements into an index in order to do more usef...
Definition TMath.cxx:1314
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:715
Double_t BesselI1(Double_t x)
Modified Bessel function K_0(x)
Definition TMath.cxx:1494
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition TMath.cxx:190
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:2560
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:908
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:686
Double_t PoissonI(Double_t x, Double_t par)
Computes the Discrete Poisson distribution function for (x,par).
Definition TMath.cxx:615
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:2178
Double_t ATan(Double_t)
Returns the principal value of the arc tangent of x, expressed in radians.
Definition TMath.h:646
Double_t StruveL1(Double_t x)
Modified Struve functions of order 0.
Definition TMath.cxx:1970
constexpr Double_t Gn()
Standard acceleration of gravity in .
Definition TMath.h:175
Double_t ASinH(Double_t)
Returns the area hyperbolic sine of x.
Definition TMath.cxx:67
Double_t LaplaceDistI(Double_t x, Double_t alpha=0, Double_t beta=1)
Computes the cumulative distribution function (lower tail integral) of Laplace distribution at point ...
Definition TMath.cxx:2383
ULong_t Hash(const void *txt, Int_t ntxt)
Calculates hash index from any char string.
Definition TMath.cxx:1408
constexpr Double_t QeUncertainty()
Elementary charge uncertainty.
Definition TMath.h:347
constexpr Double_t K()
Boltzmann's constant in : .
Definition TMath.h:248
Double_t BreitWigner(Double_t x, Double_t mean=0, Double_t gamma=1)
Calculates a Breit Wigner function with mean and gamma.
Definition TMath.cxx:442
constexpr Double_t Sqrt2()
Definition TMath.h:87
constexpr Double_t KUncertainty()
Boltzmann's constant uncertainty.
Definition TMath.h:262
constexpr Double_t Hbarcgs()
Definition TMath.h:219
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition TMath.cxx:492
Double_t Voigt(Double_t x, Double_t sigma, Double_t lg, Int_t r=4)
Computation of Voigt function (normalised).
Definition TMath.cxx:898
Double_t Student(Double_t T, Double_t ndf)
Computes density function for Student's t- distribution.
Definition TMath.cxx:2625
constexpr Double_t CUncertainty()
Speed of light uncertainty.
Definition TMath.h:129
constexpr Double_t Qe()
Elementary charge in .
Definition TMath.h:340
Double_t Ceil(Double_t x)
Rounds x upward, returning the smallest integral value that is not less than x.
Definition TMath.h:674
constexpr Double_t PiOver2()
Definition TMath.h:52
constexpr Double_t HCcgs()
Definition TMath.h:241
T MinElement(Long64_t n, const T *a)
Returns minimum of array a of length n.
Definition TMath.h:966
Double_t BetaDistI(Double_t x, Double_t p, Double_t q)
Computes the cumulative distribution function of the Beta distribution, i.e.
Definition TMath.cxx:2090
T NormCross(const T v1[3], const T v2[3], T out[3])
Calculates the Normalized Cross Product of two vectors.
Definition TMath.h:957
Int_t Finite(Double_t x)
Check if it is finite with a mask in order to be consistent in presence of fast math.
Definition TMath.h:776
Double_t TanH(Double_t)
Returns the hyperbolic tangent of x.
Definition TMath.h:624
Int_t FloorNint(Double_t x)
Returns the nearest integer of TMath::Floor(x).
Definition TMath.h:692
Double_t ACosH(Double_t)
Returns the nonnegative area hyperbolic cosine of x.
Definition TMath.cxx:81
Double_t BesselK0(Double_t x)
Modified Bessel function I_0(x)
Definition TMath.cxx:1460
Double_t BesselY0(Double_t x)
Bessel function J1(x) for any real x.
Definition TMath.cxx:1705
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:652
constexpr Double_t RUncertainty()
Universal gas constant uncertainty.
Definition TMath.h:310
Double_t BetaCf(Double_t x, Double_t a, Double_t b)
Continued fraction evaluation by modified Lentz's method used in calculation of incomplete Beta funct...
Definition TMath.cxx:2020
Long64_t LocMax(Long64_t n, const T *a)
Returns index of array with the maximum element.
Definition TMath.h:1094
Double_t ErfInverse(Double_t x)
Returns the inverse error function.
Definition TMath.cxx:208
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:2367
constexpr Double_t E()
Base of natural log: .
Definition TMath.h:94
constexpr Double_t GnUncertainty()
Standard acceleration of gravity uncertainty.
Definition TMath.h:182
constexpr Double_t Hcgs()
Definition TMath.h:196
constexpr Double_t HUncertainty()
Planck's constant uncertainty.
Definition TMath.h:203
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:762
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:80
Double_t Erfc(Double_t x)
Computes the complementary error function erfc(x).
Definition TMath.cxx:199
Double_t VavilovI(Double_t x, Double_t kappa, Double_t beta2)
Returns the value of the Vavilov cumulative distribution function (lower tail integral of the probabi...
Definition TMath.cxx:2817
constexpr Double_t Sigma()
Stefan-Boltzmann constant in : .
Definition TMath.h:271
Double_t Beta(Double_t p, Double_t q)
Calculates Beta-function Gamma(p)*Gamma(q)/Gamma(p+q).
Definition TMath.cxx:2011
constexpr Double_t Kcgs()
Definition TMath.h:255
Double_t Sq(Double_t x)
Returns x*x.
Definition TMath.h:662
Double_t Poisson(Double_t x, Double_t par)
Computes the Poisson distribution function for (x,par).
Definition TMath.cxx:587
constexpr Double_t H()
Planck's constant in : .
Definition TMath.h:189
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:668
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:727
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:680
Double_t Ldexp(Double_t x, Int_t exp)
Returns the result of multiplying x (the significant) by 2 raised to the power of exp (the exponent).
Definition TMath.h:721
Double_t BesselJ0(Double_t x)
Modified Bessel function K_1(x)
Definition TMath.cxx:1634
constexpr Double_t LogE()
Base-10 log of e (to convert ln to log)
Definition TMath.h:108
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:353
constexpr Double_t MWair()
Molecular weight of dry air 1976 US Standard Atmosphere in or
Definition TMath.h:318
constexpr Double_t Gcgs()
Definition TMath.h:144
Double_t StruveL0(Double_t x)
Struve functions of order 1.
Definition TMath.cxx:1923
Double_t NormQuantile(Double_t p)
Computes quantiles for standard normal distribution N(0, 1) at probability p.
Definition TMath.cxx:2459
constexpr Double_t Ln10()
Natural log of 10 (to convert log to ln)
Definition TMath.h:101
Double_t Hypot(Double_t x, Double_t y)
Returns sqrt(x*x + y*y)
Definition TMath.cxx:59
constexpr Double_t EulerGamma()
Euler-Mascheroni Constant.
Definition TMath.h:333
constexpr Double_t PiOver4()
Definition TMath.h:59
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:600
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=nullptr, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.
Definition TMath.cxx:1207
constexpr Double_t Pi()
Definition TMath.h:38
Double_t StruveH0(Double_t x)
Bessel function Y1(x) for positive x.
Definition TMath.cxx:1777
constexpr Double_t R()
Universal gas constant ( ) in
Definition TMath.h:303
Double_t LnGamma(Double_t z)
Computation of ln[gamma(z)] for all z.
Definition TMath.cxx:509
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Comparing floating points.
Definition TMath.h:427
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Comparing floating points.
Definition TMath.h:419
Double_t Mean(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the weighted mean of an array a with length n.
Definition TMath.h:1171
Double_t KolmogorovProb(Double_t z)
Calculates the Kolmogorov distribution function,.
Definition TMath.cxx:679
constexpr Double_t InvPi()
Definition TMath.h:66
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.
Definition TMath.cxx:1107
Double_t ChisquareQuantile(Double_t p, Double_t ndf)
Evaluate the quantiles of the chi-squared probability distribution function.
Definition TMath.cxx:2196
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:594
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:2280
Double_t SignalingNaN()
Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
Definition TMath.h:916
Double_t BreitWignerRelativistic(Double_t x, Double_t median=0, Double_t gamma=1)
Calculates a Relativistic Breit Wigner function with median and gamma.
Definition TMath.cxx:452
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:432
T * Cross(const T v1[3], const T v2[3], T out[3])
Calculates the Cross Product of two vectors: out = [v1 x v2].
Definition TMath.h:1279
void BubbleLow(Int_t Narr, Double_t *arr1, Int_t *arr2)
Opposite ordering of the array arr2[] to that of BubbleHigh.
Definition TMath.cxx:1353
Double_t BesselK(Int_t n, Double_t x)
Integer order modified Bessel function I_n(x)
Definition TMath.cxx:1561
constexpr Double_t Na()
Avogadro constant (Avogadro's Number) in .
Definition TMath.h:285
T MaxElement(Long64_t n, const T *a)
Returns maximum of array a of length n.
Definition TMath.h:974
Double_t BesselJ1(Double_t x)
Bessel function J0(x) for any real x.
Definition TMath.cxx:1669
Double_t BetaIncomplete(Double_t x, Double_t a, Double_t b)
Calculates the incomplete Beta-function.
Definition TMath.cxx:2103
constexpr Double_t Rgair()
Dry Air Gas Constant (R / MWair) in
Definition TMath.h:326
constexpr Double_t Hbar()
in :
Definition TMath.h:212
Double_t StruveH1(Double_t x)
Struve functions of order 0.
Definition TMath.cxx:1846
Double_t Freq(Double_t x)
Computation of the normal frequency function freq(x).
Definition TMath.cxx:270
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:606
Double_t LandauI(Double_t x)
Returns the cumulative (lower tail integral) of the Landau distribution function at point x.
Definition TMath.cxx:2847
Double_t StdDev(Long64_t n, const T *a, const Double_t *w=nullptr)
Definition TMath.h:532
Double_t ATanH(Double_t)
Returns the area hyperbolic tangent of x.
Definition TMath.cxx:95
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:73
Double_t BesselI0(Double_t x)
Integer order modified Bessel function K_n(x)
Definition TMath.cxx:1426
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:768
Double_t StudentI(Double_t T, Double_t ndf)
Calculates the cumulative distribution function of Student's t-distribution second parameter stands f...
Definition TMath.cxx:2648
Double_t StudentQuantile(Double_t p, Double_t ndf, Bool_t lower_tail=kTRUE)
Computes quantiles of the Student's t-distribution 1st argument is the probability,...
Definition TMath.cxx:2676
Double_t BesselY1(Double_t x)
Bessel function Y0(x) for positive x.
Definition TMath.cxx:1739
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
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:2350
constexpr Double_t GhbarC()
in
Definition TMath.h:159
constexpr Double_t HC()
in
Definition TMath.h:234
constexpr Double_t TwoPi()
Definition TMath.h:45
constexpr Double_t HbarUncertainty()
uncertainty.
Definition TMath.h:226
Double_t Infinity()
Returns an infinity as defined by the IEEE standard.
Definition TMath.h:923
Double_t ErfcInverse(Double_t x)
Returns the inverse of the complementary error function.
Definition TMath.cxx:242
Double_t SinH(Double_t)
Returns the hyperbolic sine of `x.
Definition TMath.h:612
static T Min()
Returns maximum representation for type T.
Definition TMath.h:931
static T Epsilon()
Returns minimum double representation.
Definition TMath.h:947
static T Max()
Returns minimum double representation.
Definition TMath.h:939
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2340