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