ATLAS Offline Software
Loading...
Searching...
No Matches
APWeightEntry.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5#define APWeightEntry_cxx
8
11#include <cmath>
12#include <limits>
13#include <vector>
14#include "TH1F.h"
15#include "TRandom3.h"
16#include "TMath.h"
17
18using namespace std;
19
24 m_variance(0),
27 m_sys_uncert(0),
29 m_is_trig(0),
30 m_is_nan(0),
31 m_integral(0)
32{
33 m_hist = 0;
34 m_cumul = 0;
35 m_pdf = 0;
36 m_bins = 0;
37 m_ID = 0;
38 //m_coords = 0;
39}
40
41APWeightEntry::APWeightEntry(unsigned int val_denominator, unsigned int val_numerator, double scale, bool isTrig)
44{
45 m_val_denominator = val_denominator;
46 double val_denominator_f = (double) val_denominator;
47 m_val_numerator = val_numerator;
48 double val_numerator_f = (double) val_numerator;
49 m_sys_uncert = 0.;
50 m_sys_uncert2 = 0.;
51 m_integral = 0.;
52 m_hist = 0;
53 m_cumul = 0;
54 m_bins = 0;
55 m_pdf = 0;
56 m_is_trig = isTrig;
57 m_ID = 0;
58 //m_coords = 0;
59
60 if (!isTrig) {
61 if (m_val_numerator == 0 || m_val_denominator == 0) {
63 m_variance = 0.;
66 m_is_nan = true;
67 } else {
68 m_is_nan = false;
69 m_pdf = new double[1000];
70 m_bins = new double[1001];
71
72 m_expectancy_val = scale * (val_numerator_f / val_denominator_f);
73 m_variance = (val_numerator_f / (val_denominator_f * val_denominator_f)) + ((val_numerator_f * val_numerator_f) / (val_denominator_f * val_denominator_f * val_denominator_f));
76
77 double sig_low = (std::max((double) 0.0, (m_expectancy_val - 5. * scale * sqrt(m_variance))));
78 double sig_high = (m_expectancy_val + 5. * scale * sqrt(m_variance));
79 double step = (sig_high - sig_low) * 1e-3;
80 const double inv_scale = 1. / scale;
81 double shift = val_numerator_f * log(m_expectancy_val * inv_scale)+(-val_numerator_f - val_denominator_f - 2.) * log((m_expectancy_val * inv_scale) + 1.);
82 //double shift = ln_factorialApp(m_val_numerator);
83
84 double temp_value = sig_low;
85 double arr_position = sig_low;
86
87 for (int i = 0; i < 1000; ++i) {
88 temp_value = exp(m_val_numerator * log(arr_position * inv_scale)+(-val_numerator_f - val_denominator_f - 2.) * log((arr_position * inv_scale) + 1.) - shift);
89 if (temp_value != temp_value) temp_value = 0.;
90 m_pdf[i] = temp_value;
91 m_integral += temp_value;
92 m_bins[i] = arr_position;
93 arr_position += step;
94 }
95 m_bins[1000] = arr_position;
96 }
97 } else {
98 if (m_val_denominator == 0) {
100 m_variance = 0.;
103 m_is_nan = true;
104 } else {
105 m_is_nan = false;
106 m_pdf = new double[1000];
107 m_bins = new double[1001];
108
109 const double inv_val_denominator_f = 1. / val_denominator_f;
110 m_expectancy_val = val_numerator_f * inv_val_denominator_f;
111 m_variance = (m_expectancy_val * (1. - m_expectancy_val)) * inv_val_denominator_f;
112 double var_sqrt = sqrt(m_variance);
113 double sig_low = std::max(0., m_expectancy_val - 5. * var_sqrt);
114 double sig_high = std::min(1., m_expectancy_val + 5. * var_sqrt);
115 if (val_numerator_f == 0) {
116 sig_low = 0.;
117 sig_high = std::min(1., (8. * inv_val_denominator_f));
118 }
119 if (val_numerator_f == val_denominator_f) {
120 sig_low = std::max(0., (1. - 8. * inv_val_denominator_f));
121 sig_high = 1.;
122 }
123
124 double step = fabs(sig_high - sig_low) * 1e-3;
125 double shift = 0.;
126 if (val_numerator_f == 0 || val_numerator_f == val_denominator_f) {
127 shift = 0.;
128 } else {
129 shift = val_numerator_f * log(m_expectancy_val) + (val_denominator_f - val_numerator_f) * log(1. - m_expectancy_val);
130 }
131
132 double temp_value = sig_low;
133 double arr_position = sig_low;
134
135 if (val_numerator_f == 0) {
136 for (int i = 0; i < 1000; ++i) {
137 temp_value = exp(val_denominator_f * log(1. - arr_position) - shift);
138 m_pdf[i] = temp_value;
139 m_integral += temp_value;
140 m_bins[i] = arr_position;
141 arr_position += step;
142 }
143 } else {
144 if (val_numerator_f == val_denominator_f) {
145 for (int i = 0; i < 1000; ++i) {
146 temp_value = exp(val_numerator_f * log(arr_position) - shift);
147 m_pdf[i] = temp_value;
148 m_integral += temp_value;
149 m_bins[i] = arr_position;
150 arr_position += step;
151 }
152 } else {
153 for (int i = 0; i < 1000; ++i) {
154 temp_value = exp(val_numerator_f * log(arr_position) + (val_denominator_f - val_numerator_f) * log(1. - arr_position) - shift);
155 m_pdf[i] = temp_value;
156 m_integral += temp_value;
157 m_bins[i] = arr_position;
158 arr_position += step;
159 }
160 }
161 }
162 m_bins[1000] = arr_position;
163 }
164 }
165
166}
167
168void APWeightEntry::ReadEfficiency(double efficiency, double err_low, double err_high) {
170 m_val_numerator = 0;
171 m_sys_uncert = 0.;
172 m_sys_uncert2 = 0.;
174 m_stat_uncert_low = fabs(err_low);
175 m_stat_uncert_high = fabs(err_high);
178}
179
180void APWeightEntry::SetCoordinates(const std::vector<int>& coords,
181 const std::vector<int>& n_dim_origin) {
182 m_coords.clear();
183 m_n_dim_origin.clear();
184 m_coords = coords;
185 m_n_dim_origin = n_dim_origin;
186}
187
189 delete [] m_pdf;
190 m_pdf = 0;
191 delete [] m_bins;
192 m_bins = 0;
193 delete [] m_cumul;
194 m_cumul = 0;
195 delete m_hist;
196}
197
199 return m_val_denominator;
200}
201
202unsigned int APWeightEntry::GetValNumerator() const {
203 return m_val_numerator;
204}
205
207 return m_expectancy_val;
208}
209
211 return m_variance;
212}
213
215 return m_stat_uncert_low;
216}
217
219 return m_stat_uncert_high;
220}
221
223 return m_sys_uncert;
224}
225
227 return m_sys_uncert2;
228}
229
231
233 if (m_is_nan || m_integral < std::numeric_limits<double>::epsilon()) {
234 return 1e-100;
235 }
236 if (!m_cumul) _ComputeCum();
237 double rn = random->Rndm();
238 int ibin = TMath::BinarySearch(1000, m_cumul, rn);
239 return m_bins[ibin];
240 //if (rn > m_cumul[ibin]) x += fabs(m_bins[1]-m_bins[0]) * (rn-m_cumul[ibin])/(m_cumul[ibin+1] - m_cumul[ibin]);
241}
242
244 if (!m_hist) _CreateHist();
245 return m_hist;
246}
247
248unsigned int APWeightEntry::GetID() const {
249 return m_ID;
250}
251
252const vector<int> & APWeightEntry::GetCoords() const {
253 return m_coords;
254}
255
256const vector< int > & APWeightEntry::GetOriginalDimensions() const {
257 return m_n_dim_origin;
258}
259
261 return m_is_nan;
262}
263
265 return m_is_trig;
266}
267
272
273void APWeightEntry::SetID(unsigned int id) {
274 m_ID = id;
275}
276
278 m_hist = new TH1F("", "", 1000, m_bins);
279 if (m_is_nan == false) {
280 for (int i = 0; i < 1000; ++i) {
281 m_hist->SetBinContent(i + 1, m_pdf[i]);
282 }
283 }
284}
285
287 if (m_cumul) delete [] m_cumul;
288 m_cumul = new double[1001];
289 m_cumul[0] = 0.;
290 for (int i = 0; i < 1000; ++i) {
291 m_cumul[i + 1] = m_cumul[i] + m_pdf[i];
292 }
293 const double inv_normalize = 1. / m_cumul[1000];
294 for (int i = 1; i < 1000; ++i) {
295 m_cumul[i] *= inv_normalize;
296 }
297}
Thread-local TRandom generator.
Define macros for attributes used to control the static checker.
#define ATLAS_THREAD_SAFE
double GetStatUncertHigh() const
Get upper bound of asymmetric statistical uncertainty.
void SetID(unsigned int id)
Set the internal ID (used by APReweight/APReweight2D/APReweight3D/APReweightND).
const std::vector< int > & GetOriginalDimensions() const
Returns the dimensions and amounts of bins for each dimension of the original histogram.
double GetStatUncertLow() const
Get lower bound of asymmetric statistical uncertainty.
double GetSysUncert2() const
Get absolute systematic uncertainty squared value of efficiency/weight.
virtual ~APWeightEntry()
Default destructor.
bool m_is_nan
Flag, set to true if denominator is zero.
double * m_cumul
Histograms to hold the probability distribution and the cumulative distribution.
void ReadEfficiency(double efficiency, double err_low, double err_high)
Read efficiencies and upper/lower uncertainty (if numerator/denominator not applicable (e....
double GetSysUncert() const
Get absolute systematic uncertainty value of efficiency/weight.
const std::vector< int > & GetCoords() const
Returns the coordinates of the current entry in the original histogram.
unsigned int GetID() const
Returns the internal ID (used by APReweight/APReweight2D/APReweight3D/APReweightND).
APWeightEntry()
Default constructor.
double GetRandom()
Get random number from PDF.
double m_stat_uncert_high
Holds upper bound of asymmetric statistical uncertainty.
std::vector< int > m_coords
Holds the coordinates of the current entry in the original histogram.
double GetVariance() const
Get Variance of efficiency/weight (classical binomial/poisson model).
void _CreateHist()
Creates a TH1F instance from the arrays if necessary.
unsigned int GetValNumerator() const
Get value of original numerator.
double m_integral
Holds the integral of the probability distribution.
void _ComputeCum()
Calculates the cumulative function of the pdf if necessary.
double m_sys_uncert2
Holds absolute systematic uncertainty squared value of efficiency/weight.
TH1F * GetPDF()
Returns the calculated PDF.
unsigned int m_val_numerator
Holds the value of original numerator.
bool m_is_trig
Flag, set to true if weight entry is trigger based.
double GetExpectancy() const
Get Expectancy value of efficiency/weight.
double m_sys_uncert
Holds absolute systematic uncertainty value of efficiency/weight.
double m_stat_uncert_low
Holds lower bound of asymmetric statistical uncertainty.
unsigned int GetValDenominator() const
Get value of original denominator.
unsigned int m_val_denominator
Holds the value of original denominator.
bool IsNaN() const
Returns true if instance is NaN.
void SetSystUncert(double rel_uncert)
Set the relative (!) systematic uncertainty for the efficiency/weight.
TH1F * m_hist
Holds the TH1F instance from the arrays if computed.
double m_variance
Holds Variance of efficiency/weight (classical binomial/poisson model).
void SetCoordinates(const std::vector< int > &coords, const std::vector< int > &n_dim_origin)
std::vector< int > m_n_dim_origin
Holds the amount of dimensions and bins per axis in the original histogram.
double m_expectancy_val
Holds the Expectancy value of efficiency/weight.
double * m_bins
unsigned int m_ID
Holds internal ID (used by APReweight/APReweight2D/APReweight3D/APReweightND).
bool IsTrig() const
Returns true if instance is trigger based.
Thread-local TRandom generator.
Definition TRandomTLS.h:29
void efficiency(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
STL namespace.