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