|
|
Rivet 4.0.0
|
2#ifndef RIVET_MathUtils_HH
3#define RIVET_MathUtils_HH
5#include "Rivet/Math/MathConstants.hh"
22 template < typename NUM>
23 inline typename std::enable_if<std::is_floating_point<NUM>::value, bool>::type
24 isZero(NUM val, double tolerance=1e-8) {
25 return fabs(val) < tolerance;
32 template < typename NUM>
33 inline typename std::enable_if<std::is_integral<NUM>::value, bool>::type
39 template < typename NUM>
40 inline typename std::enable_if<std::is_floating_point<NUM>::value, bool>::type
41 isNaN(NUM val) { return std::isnan(val); }
44 template < typename NUM>
45 inline typename std::enable_if<std::is_floating_point<NUM>::value, bool>::type
46 notNaN(NUM val) { return !std::isnan(val); }
53 template < typename N1, typename N2>
54 inline typename std::enable_if<
55 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value &&
56 (std::is_floating_point<N1>::value || std::is_floating_point<N2>::value), bool>::type
58 const double absavg = (std::abs(a) + std::abs(b))/2.0;
59 const double absdiff = std::abs(a - b);
60 const bool rtn = ( isZero(a) && isZero(b)) || absdiff < tolerance*absavg;
68 template < typename N1, typename N2>
69 inline typename std::enable_if<
70 std::is_integral<N1>::value && std::is_integral<N2>::value, bool>::type
79 template < typename N1, typename N2>
80 inline typename std::enable_if<
81 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value, bool>::type
90 template < typename N1, typename N2>
91 inline typename std::enable_if<
92 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value, bool>::type
98 template < typename N1, typename N2>
99 inline typename std::enable_if<
100 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value,
101 typename std::common_type<N1,N2>::type >::type
103 return a > b ? b : a;
107 template < typename N1, typename N2>
108 inline typename std::enable_if<
109 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value,
110 typename std::common_type<N1,N2>::type >::type
112 return a > b ? a : b;
130 template < typename N1, typename N2, typename N3>
131 inline typename std::enable_if<
132 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
135 if (lowbound == OPEN && highbound == OPEN) {
136 return (value > low && value < high);
137 } else if (lowbound == OPEN && highbound == CLOSED) {
138 return (value > low && value <= high);
139 } else if (lowbound == CLOSED && highbound == OPEN) {
140 return (value >= low && value < high);
142 return (value >= low && value <= high);
150 template < typename N1, typename N2, typename N3>
151 inline typename std::enable_if<
152 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
155 if (lowbound == OPEN && highbound == OPEN) {
156 return (value > low && value < high);
157 } else if (lowbound == OPEN && highbound == CLOSED) {
159 } else if (lowbound == CLOSED && highbound == OPEN) {
167 template < typename N1, typename N2, typename N3>
168 inline typename std::enable_if<
169 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
172 return inRange(value, lowhigh.first, lowhigh.second, lowbound, highbound);
181 template < typename N1, typename N2, typename N3>
182 inline typename std::enable_if<
183 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
185 return inRange(val, low, high, CLOSED, OPEN);
191 template < typename N1, typename N2, typename N3>
192 inline typename std::enable_if<
193 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
195 return inRange(val, low, high, CLOSED, CLOSED);
201 template < typename N1, typename N2, typename N3>
202 inline typename std::enable_if<
203 std::is_arithmetic<N1>::value && std::is_arithmetic<N2>::value && std::is_arithmetic<N3>::value, bool>::type
205 return inRange(val, low, high, OPEN, OPEN);
217 template < typename NUM>
218 inline typename std::enable_if<std::is_arithmetic<NUM>::value, NUM>::type
228 template < typename NUM>
229 inline typename std::enable_if<std::is_arithmetic<NUM>::value, NUM>::type
232 return sqrt(a*a + b*b);
240 template < typename NUM>
241 inline typename std::enable_if<std::is_arithmetic<NUM>::value, NUM>::type
244 return sqrt(a*a + b*b + c*c);
249 inline double safediv( double num, double den, double fail=0.0) {
250 return (! isZero(den)) ? num/den : fail;
254 template < typename NUM>
255 constexpr inline typename std::enable_if<std::is_arithmetic<NUM>::value, NUM>::type
257 if (exp == 0) return (NUM) 1;
258 else if (exp == 1) return val;
259 return val * intpow(val, exp-1);
263 template < typename NUM>
264 constexpr inline typename std::enable_if<std::is_arithmetic<NUM>::value, int>::type
266 if ( isZero(val)) return ZERO;
267 const int valsign = (val > 0) ? PLUS : MINUS;
278 inline double cdfBW( double x, double mu, double gamma) {
280 const double xn = (x - mu)/gamma;
281 return std::atan(xn)/M_PI + 0.5;
285 inline double invcdfBW( double p, double mu, double gamma) {
286 const double xn = std::tan(M_PI*(p-0.5));
287 return gamma*xn + mu;
302 inline vector<double> linspace( size_t nbins, double start, double end, bool include_end= true) {
305 const double interval = (end-start)/ static_cast<double>(nbins);
306 for ( size_t i = 0; i < nbins; ++i) {
307 rtn.push_back(start + i*interval);
309 assert(rtn.size() == nbins);
310 if (include_end) rtn.push_back(end);
326 inline vector<double> aspace( double step, double start, double end, bool include_end= true, double tol=1e-2) {
327 assert( (end-start)*step > 0);
331 if (next > end) break;
336 if (end - rtn[rtn.size()-1] > tol*step) rtn.push_back(end);
345 inline vector<double> fnspace( size_t nbins, double start, double end,
346 const std::function< double( double)>& fn, const std::function< double( double)>& invfn,
347 bool include_end= true) {
350 const double pmin = fn(start);
351 const double pmax = fn(end);
352 const vector<double> edges = linspace(nbins, pmin, pmax, false);
353 assert(edges.size() == nbins);
354 vector<double> rtn; rtn.reserve(nbins+1);
355 rtn.push_back(start);
356 for ( size_t i = 1; i < edges.size(); ++i) {
357 rtn.push_back(invfn(edges[i]));
359 assert(rtn.size() == nbins);
360 if (include_end) rtn.push_back(end);
374 inline vector<double> logspace( size_t nbins, double start, double end, bool include_end= true) {
375 return fnspace(nbins, start, end,
376 []( double x){ return std::log(x); },
377 []( double x){ return std::exp(x); },
391 inline vector<double> powspace( size_t nbins, double start, double end, double npow, bool include_end= true) {
393 return fnspace(nbins, start, end,
394 [&]( double x){ return std::pow(x, npow); },
395 [&]( double x){ return std::pow(x, 1/npow); },
410 inline vector<double> powdbnspace( size_t nbins, double start, double end, double npow, bool include_end= true) {
412 return fnspace(nbins, start, end,
413 [&]( double x){ return std::pow(x, npow+1) / (npow+1); },
414 [&]( double x){ return std::pow((npow+1) * x, 1/(npow+1)); },
426 inline vector<double> bwdbnspace( size_t nbins, double start, double end, double mu, double gamma, bool include_end= true) {
427 return fnspace(nbins, start, end,
428 [&]( double x){ return cdfBW(x, mu, gamma); },
429 [&]( double x){ return invcdfBW(x, mu, gamma); },
435 template < typename NUM, typename CONTAINER>
436 inline typename std::enable_if<std::is_arithmetic<NUM>::value && std::is_arithmetic<typename CONTAINER::value_type>::value, int>::type
437 _binIndex(NUM val, const CONTAINER& binedges, bool allow_overflow= false) {
438 if (val < *begin(binedges)) return -1;
440 if (val >= *(end(binedges)-1)) return allow_overflow ? int(binedges.size())-1 : -1;
441 auto it = std::upper_bound(begin(binedges), end(binedges), val);
442 return std::distance(begin(binedges), --it);
453 template < typename NUM1, typename NUM2>
454 inline typename std::enable_if<std::is_arithmetic<NUM1>::value && std::is_arithmetic<NUM2>::value, int>::type
455 binIndex(NUM1 val, std::initializer_list<NUM2> binedges, bool allow_overflow= false) {
456 return _binIndex(val, binedges, allow_overflow);
467 template < typename NUM, typename CONTAINER>
468 inline typename std::enable_if<std::is_arithmetic<NUM>::value && std::is_arithmetic<typename CONTAINER::value_type>::value, int>::type
469 binIndex(NUM val, const CONTAINER& binedges, bool allow_overflow= false) {
470 return _binIndex(val, binedges, allow_overflow);
481 template < typename NUM>
482 inline typename std::enable_if<std::is_arithmetic<NUM>::value, NUM>::type
484 if (sample.empty()) throw RangeError( "Can't compute median of an empty set");
485 vector<NUM> tmp = sample;
486 std::sort(tmp.begin(), tmp.end());
487 const size_t imid = tmp.size()/2;
488 if (sample.size() % 2 == 0) return (tmp.at(imid-1) + tmp.at(imid)) / 2.0;
489 else return tmp.at(imid);
495 template < typename NUM>
496 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
497 mean( const vector<NUM>& sample) {
498 if (sample.empty()) throw RangeError( "Can't compute mean of an empty set");
500 for ( size_t i = 0; i < sample.size(); ++i) {
503 return mean/sample.size();
508 template < typename NUM>
509 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
511 if (sample.empty()) throw RangeError( "Can't compute mean_err of an empty set");
513 for ( size_t i = 0; i < sample.size(); ++i) {
514 mean_e += sqrt(sample[i]);
516 return mean_e/sample.size();
522 template < typename NUM>
523 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
524 covariance( const vector<NUM>& sample1, const vector<NUM>& sample2) {
525 if (sample1.empty() || sample2.empty()) throw RangeError( "Can't compute covariance of an empty set");
526 if (sample1.size() != sample2.size()) throw RangeError( "Sizes of samples must be equal for covariance calculation");
527 const double mean1 = mean(sample1);
528 const double mean2 = mean(sample2);
529 const size_t N = sample1.size();
531 for ( size_t i = 0; i < N; i++) {
532 const double cov_i = (sample1[i] - mean1)*(sample2[i] - mean2);
535 if (N > 1) return cov/(N-1);
541 template < typename NUM>
542 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
544 if (sample1.empty() || sample2.empty()) throw RangeError( "Can't compute covariance_err of an empty set");
545 if (sample1.size() != sample2.size()) throw RangeError( "Sizes of samples must be equal for covariance_err calculation");
546 const double mean1 = mean(sample1);
547 const double mean2 = mean(sample2);
548 const double mean1_e = mean_err(sample1);
549 const double mean2_e = mean_err(sample2);
550 const size_t N = sample1.size();
552 for ( size_t i = 0; i < N; i++) {
553 const double cov_i = (sqrt(sample1[i]) - mean1_e)*(sample2[i] - mean2) +
554 (sample1[i] - mean1)*(sqrt(sample2[i]) - mean2_e);
557 if (N > 1) return cov_e/(N-1);
564 template < typename NUM>
565 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
566 correlation( const vector<NUM>& sample1, const vector<NUM>& sample2) {
567 const double cov = covariance(sample1, sample2);
568 const double var1 = covariance(sample1, sample1);
569 const double var2 = covariance(sample2, sample2);
571 const double corr_strength = correlation*sqrt(var2/var1);
572 return corr_strength;
577 template < typename NUM>
578 inline typename std::enable_if<std::is_arithmetic<NUM>::value, double>::type
580 const double cov = covariance(sample1, sample2);
581 const double var1 = covariance(sample1, sample1);
582 const double var2 = covariance(sample2, sample2);
591 cov/(2*pow(3./2., var1*var2)) * (var1_e * var2 + var1 * var2_e);
595 correlation/(2*sqrt(var2/var1)) * (var2_e/var1 - var2*var1_e/pow(2, var2));
597 return corr_strength_err;
610 inline double _mapAngleM2PITo2Pi( double angle) {
612 if ( isZero(rtn)) return 0;
619 double rtn = _mapAngleM2PITo2Pi( angle);
620 if ( isZero(rtn)) return 0;
623 assert(rtn > - PI && rtn <= PI);
629 double rtn = _mapAngleM2PITo2Pi( angle);
630 if ( isZero(rtn)) return 0;
631 if (rtn < 0) rtn += TWOPI;
632 if (rtn == TWOPI) rtn = 0;
633 assert(rtn >= 0 && rtn < TWOPI);
640 if ( isZero(rtn)) return 0;
641 assert(rtn > 0 && rtn <= PI);
655 throw Rivet::UserError( "The specified phi mapping scheme is not implemented");
670 return sign ? x : fabs(x);
677 const double x = eta1 - eta2;
678 return sign ? x : fabs(x);
685 const double x = y1 - y2;
686 return sign? x : fabs(x);
691 inline double deltaR2( double rap1, double phi1, double rap2, double phi2) {
692 const double dphi = deltaPhi(phi1, phi2);
693 return sqr(rap1-rap2) + sqr(dphi);
698 inline double deltaR( double rap1, double phi1, double rap2, double phi2) {
699 return sqrt( deltaR2(rap1, phi1, rap2, phi2));
705 throw std::runtime_error( "Divergent positive rapidity");
709 throw std::runtime_error( "Divergent negative rapidity");
712 return 0.5*log((E+pz)/(E-pz));
720 inline double mT( double pT1, double pT2, double dphi) {
721 return sqrt(2*pT1*pT2 * (1 - cos(dphi)) );
Definition MC_CENT_PPB_Projections.hh:10
double deltaR(double rap1, double phi1, double rap2, double phi2) Definition MathUtils.hh:698
double deltaPhi(double phi1, double phi2, bool sign=false) Calculate the difference between two angles in radians. Definition MathUtils.hh:668
vector< double > aspace(double step, double start, double end, bool include_end=true, double tol=1e-2) Make a list of values equally spaced by step between start and end inclusive. Definition MathUtils.hh:326
double deltaEta(double eta1, double eta2, bool sign=false) Definition MathUtils.hh:676
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&std::is_arithmetic< N3 >::value, bool >::type in_range(N1 val, N2 low, N3 high) Boolean function to determine if value is within the given range. Definition MathUtils.hh:184
PhiMapping Enum for range of to be mapped into. Definition MathConstants.hh:49
vector< double > logspace(size_t nbins, double start, double end, bool include_end=true) Make a list of nbins + 1 values exponentially spaced between start and end inclusive. Definition MathUtils.hh:374
constexpr std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type intpow(NUM val, unsigned int exp) A more efficient version of pow for raising numbers to integer powers. Definition MathUtils.hh:256
std::enable_if< std::is_floating_point< NUM >::value, bool >::type notNaN(NUM val) Check if a number is non-NaN. Definition MathUtils.hh:46
double mapAngle0To2Pi(double angle) Map an angle into the range [0, 2PI). Definition MathUtils.hh:628
std::enable_if< std::is_arithmetic< NUM >::value, double >::type correlation_err(const vector< NUM > &sample1, const vector< NUM > &sample2) Definition MathUtils.hh:579
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&std::is_arithmetic< N3 >::value, bool >::type in_closed_range(N1 val, N2 low, N3 high) Boolean function to determine if value is within the given range. Definition MathUtils.hh:194
static const double TWOPI A pre-defined value of . Definition MathConstants.hh:16
std::enable_if< std::is_arithmetic< NUM >::value, double >::type covariance_err(const vector< NUM > &sample1, const vector< NUM > &sample2) Definition MathUtils.hh:543
double deltaR2(double rap1, double phi1, double rap2, double phi2) Definition MathUtils.hh:691
double mT(double pT1, double pT2, double dphi) Definition MathUtils.hh:720
std::enable_if< std::is_arithmetic< NUM >::value, double >::type mean_err(const vector< NUM > &sample) Definition MathUtils.hh:510
std::enable_if< std::is_arithmetic< NUM >::value, double >::type covariance(const vector< NUM > &sample1, const vector< NUM > &sample2) Definition MathUtils.hh:524
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value, typenamestd::common_type< N1, N2 >::type >::type max(N1 a, N2 b) Get the maximum of two numbers. Definition MathUtils.hh:111
std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type add_quad(NUM a, NUM b) Named number-type addition in quadrature operation. Definition MathUtils.hh:231
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&std::is_arithmetic< N3 >::value, bool >::type fuzzyInRange(N1 value, N2 low, N3 high, RangeBoundary lowbound=CLOSED, RangeBoundary highbound=OPEN) Determine if value is in the range low to high, for floating point numbers. Definition MathUtils.hh:153
std::enable_if< std::is_floating_point< NUM >::value, bool >::type isNaN(NUM val) Check if a number is NaN. Definition MathUtils.hh:41
vector< double > fnspace(size_t nbins, double start, double end, const std::function< double(double)> &fn, const std::function< double(double)> &invfn, bool include_end=true) Definition MathUtils.hh:345
constexpr std::enable_if< std::is_arithmetic< NUM >::value, int >::type sign(NUM val) Find the sign of a number. Definition MathUtils.hh:265
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value, typenamestd::common_type< N1, N2 >::type >::type min(N1 a, N2 b) Get the minimum of two numbers. Definition MathUtils.hh:102
double mapAngleMPiToPi(double angle) Map an angle into the range (-PI, PI]. Definition MathUtils.hh:618
std::enable_if< std::is_arithmetic< NUM1 >::value &&std::is_arithmetic< NUM2 >::value, int >::type binIndex(NUM1 val, std::initializer_list< NUM2 > binedges, bool allow_overflow=false) Return the bin index of the given value, val, given a vector of bin edges. Definition MathUtils.hh:455
RangeBoundary Definition MathUtils.hh:125
std::enable_if< std::is_arithmetic< NUM >::value, double >::type correlation(const vector< NUM > &sample1, const vector< NUM > &sample2) Definition MathUtils.hh:566
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&std::is_arithmetic< N3 >::value, bool >::type inRange(N1 value, N2 low, N3 high, RangeBoundary lowbound=CLOSED, RangeBoundary highbound=OPEN) Determine if value is in the range low to high, for floating point numbers. Definition MathUtils.hh:133
static const double PI Definition MathConstants.hh:13
std::enable_if< std::is_arithmetic< NUM >::value, double >::type mean(const vector< NUM > &sample) Definition MathUtils.hh:497
vector< double > powspace(size_t nbins, double start, double end, double npow, bool include_end=true) Make a list of nbins + 1 values power-law spaced between start and end inclusive. Definition MathUtils.hh:391
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value, bool >::type fuzzyLessEquals(N1 a, N2 b, double tolerance=1e-5) Compare two floating point numbers for <= with a degree of fuzziness. Definition MathUtils.hh:93
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&(std::is_floating_point< N1 >::value||std::is_floating_point< N2 >::value), bool >::type fuzzyEquals(N1 a, N2 b, double tolerance=1e-5) Compare two numbers for equality with a degree of fuzziness. Definition MathUtils.hh:57
double cdfBW(double x, double mu, double gamma) CDF for the Breit-Wigner distribution. Definition MathUtils.hh:278
double safediv(double num, double den, double fail=0.0) Definition MathUtils.hh:249
std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type median(const vector< NUM > &sample) Definition MathUtils.hh:483
double deltaRap(double y1, double y2, bool sign=false) Definition MathUtils.hh:684
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value &&std::is_arithmetic< N3 >::value, bool >::type in_open_range(N1 val, N2 low, N3 high) Boolean function to determine if value is within the given range. Definition MathUtils.hh:204
double mapAngle(double angle, PhiMapping mapping) Map an angle into the enum-specified range. Definition MathUtils.hh:646
vector< double > linspace(size_t nbins, double start, double end, bool include_end=true) Make a list of nbins + 1 values equally spaced between start and end inclusive. Definition MathUtils.hh:302
double mapAngle0ToPi(double angle) Map an angle into the range [0, PI]. Definition MathUtils.hh:638
std::enable_if< std::is_arithmetic< NUM >::value, NUM >::type sqr(NUM a) Named number-type squaring operation. Definition MathUtils.hh:219
double invcdfBW(double p, double mu, double gamma) Inverse CDF for the Breit-Wigner distribution. Definition MathUtils.hh:285
std::enable_if< std::is_floating_point< NUM >::value, bool >::type isZero(NUM val, double tolerance=1e-8) Compare a number to zero. Definition MathUtils.hh:24
double angle(const Vector2 &a, const Vector2 &b) Angle (in radians) between two 2-vectors. Definition Vector2.hh:177
std::enable_if< std::is_arithmetic< N1 >::value &&std::is_arithmetic< N2 >::value, bool >::type fuzzyGtrEquals(N1 a, N2 b, double tolerance=1e-5) Compare two numbers for >= with a degree of fuzziness. Definition MathUtils.hh:82
vector< double > powdbnspace(size_t nbins, double start, double end, double npow, bool include_end=true) Make a list of nbins + 1 values equally spaced in the CDF of x^n between start and end inclusive. Definition MathUtils.hh:410
vector< double > bwdbnspace(size_t nbins, double start, double end, double mu, double gamma, bool include_end=true) Make a list of nbins + 1 values spaced for equal area Breit-Wigner binning between start and end incl... Definition MathUtils.hh:426
double rapidity(double E, double pz) Calculate a rapidity value from the supplied energy E and longitudinal momentum pz. Definition MathUtils.hh:703
Error for e.g. use of invalid bin ranges. Definition Exceptions.hh:22
Error specialisation for where the problem is between the chair and the computer. Definition Exceptions.hh:55
|