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