1 #ifndef PhysicsTools_Utilities_interface_LumiReWeighting_h 2 #define PhysicsTools_Utilities_interface_LumiReWeighting_h 25 #include "TStopwatch.h" 38 class PoissonMeanShifter {
48 static const double p0_minus[20] = {1., 1., 1., 1., 1., 1., 1., 1., 1., 1.,
49 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
50 static const double p1_minus[20] = {-0.677786, -0.619614, -0.49465, -0.357963, -0.238359, -0.110002, 0.0348629,
51 0.191263, 0.347648, 0.516615, 0.679646, 0.836673, 0.97764, 1.135,
52 1.29922, 1.42467, 1.55901, 1.61762, 1.67275, 1.96008};
53 static const double p2_minus[20] = {
54 0.526164, 0.251816, 0.11049, 0.026917, -0.0464692, -0.087022, -0.0931581, -0.0714295, -0.0331772, 0.0347473,
55 0.108658, 0.193048, 0.272314, 0.376357, 0.4964, 0.58854, 0.684959, 0.731063, 0.760044, 1.02386};
57 static const double p1_expoM[5] = {1.63363e-03, 6.79290e-04, 3.69900e-04, 2.24349e-04, 9.87156e-06};
59 static const double p2_expoM[5] = {2.64692, 3.26585, 3.53229, 4.18035, 5.64027};
61 static const double p0_plus[20] = {1., 1., 1., 1., 1., 1., 1., 1., 1., 1.,
62 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
63 static const double p1_plus[20] = {-0.739059, -0.594445, -0.477276, -0.359707, -0.233573, -0.103458, 0.0373401,
64 0.176571, 0.337617, 0.499074, 0.675126, 0.840522, 1.00917, 1.15847,
65 1.23816, 1.44271, 1.52982, 1.46385, 1.5802, 0.988689};
66 static const double p2_plus[20] = {0.208068, 0.130033, 0.0850356, 0.0448344, 0.000749832,
67 -0.0331347, -0.0653281, -0.0746009, -0.0800667, -0.0527636,
68 -0.00402649, 0.103338, 0.261261, 0.491084, 0.857966,
69 1.19495, 1.75071, 2.65559, 3.35433, 5.48835};
71 static const double p1_expoP[5] = {1.42463e-01, 4.18966e-02, 1.12697e-01, 1.66197e-01, 1.50768e-01};
73 static const double p2_expoP[5] = {1.98758, 2.27217, 2.26799, 2.38455, 2.52428};
77 for (
int ibin = 0; ibin < 20; ibin++) {
79 Pweight_[ibin] = p0_minus[ibin] + p1_minus[ibin] * Shift + p2_minus[ibin] * Shift * Shift;
81 Pweight_[ibin] = p0_plus[ibin] + p1_plus[ibin] * Shift + p2_plus[ibin] * Shift * Shift;
87 for (
int ibin = 20; ibin < 25; ibin++) {
89 Pweight_[ibin] = p1_expoM[ibin - 20] *
exp(p2_expoM[ibin - 20] * Shift);
91 Pweight_[ibin] = p1_expoP[ibin - 20] *
exp(p2_expoP[ibin - 20] * Shift);
97 if (ibin < 25 && ibin >= 0) {
105 int ibin =
int(pvnum);
107 if (ibin < 25 && ibin >= 0) {
148 std::cout <<
" Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
152 for (
int ibin = 1; ibin <
NBins + 1; ++ibin) {
153 std::cout <<
" " << ibin - 1 <<
" " <<
weights_->GetBinContent(ibin) << std::endl;
161 LumiReWeighting(
const std::vector<float>& MC_distr,
const std::vector<float>& Lumi_distr) {
168 if (MC_distr.size() != Lumi_distr.size()) {
169 std::cerr <<
"ERROR: LumiReWeighting: input vectors have different sizes. Quitting... \n";
173 Int_t
NBins = MC_distr.size();
179 TH1F* den =
new TH1F(
"lumidenom",
"lumidenom",
NBins, -0.5,
float(
NBins) - 0.5);
181 for (
int ibin = 1; ibin <
NBins + 1; ++ibin) {
182 weights_->SetBinContent(ibin, Lumi_distr[ibin - 1]);
183 Data_distr_->SetBinContent(ibin, Lumi_distr[ibin - 1]);
184 den->SetBinContent(ibin, MC_distr[ibin - 1]);
185 MC_distr_->SetBinContent(ibin, MC_distr[ibin - 1]);
190 float deltaH =
weights_->Integral();
191 if (fabs(1.0 - deltaH) > 0.02) {
195 float deltaMC = den->Integral();
196 if (fabs(1.0 - deltaMC) > 0.02) {
197 den->Scale(1.0 / deltaMC);
203 std::cout <<
" Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
205 for (
int ibin = 1; ibin <
NBins + 1; ++ibin) {
206 std::cout <<
" " << ibin - 1 <<
" " <<
weights_->GetBinContent(ibin) << std::endl;
217 TH3D* WHist =
new TH3D(
"WHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
218 TH3D* DHist =
new TH3D(
"DHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
219 TH3D* MHist =
new TH3D(
"MHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
224 std::cout <<
" MC and Data distributions are not initialized! You must call the LumiReWeighting constructor. " 230 double MC_ints[50][50][50];
231 double Data_ints[50][50][50];
233 for (
int i = 0;
i < 50;
i++) {
234 for (
int j = 0;
j < 50;
j++) {
235 for (
int k = 0;
k < 50;
k++) {
236 MC_ints[
i][
j][
k] = 0.;
237 Data_ints[
i][
j][
k] = 0.;
249 for (
int i = 1;
i < 50; ++
i) {
256 double probi, probj, probk;
263 for (
int jbin = 1; jbin < NMCbin + 1; jbin++) {
265 xweight =
MC_distr_->GetBinContent(jbin);
274 std::cout <<
"LumiReweighting:BadInputValue" 275 <<
" Your histogram generates MC luminosity values less than zero!" 276 <<
" Please Check. Terminating." << std::endl;
287 for (
int i = 1;
i < 50; ++
i) {
293 for (
int i = 0;
i < 50;
i++) {
295 for (
int j = 0;
j < 50;
j++) {
297 for (
int k = 0;
k < 50;
k++) {
300 MC_ints[
i][
j][
k] = MC_ints[
i][
j][
k] + probi * probj * probk * xweight;
308 for (
int jbin = 1; jbin < NDatabin + 1; jbin++) {
314 std::cout <<
"LumiReweighting:BadInputValue" 315 <<
" Your histogram generates MC luminosity values less than zero!" 316 <<
" Please Check. Terminating." << std::endl;
327 for (
int i = 1;
i < 50; ++
i) {
334 for (
int i = 0;
i < 50;
i++) {
336 for (
int j = 0;
j < 50;
j++) {
338 for (
int k = 0;
k < 50;
k++) {
341 Data_ints[
i][
j][
k] = Data_ints[
i][
j][
k] + probi * probj * probk * xweight;
347 for (
int i = 0;
i < 50;
i++) {
349 for (
int j = 0;
j < 50;
j++) {
350 for (
int k = 0;
k < 50;
k++) {
351 if ((MC_ints[
i][
j][
k]) > 0.) {
357 DHist->SetBinContent(
i + 1,
j + 1,
k + 1, Data_ints[
i][
j][
k]);
358 MHist->SetBinContent(
i + 1,
j + 1,
k + 1, MC_ints[
i][
j][
k]);
365 if (!WeightOutputFile.empty()) {
366 std::cout <<
" 3D Weight Matrix initialized! " << std::endl;
367 std::cout <<
" Writing weights to file " << WeightOutputFile <<
" for re-use... " << std::endl;
369 TFile*
outfile =
new TFile(WeightOutputFile.c_str(),
"RECREATE");
382 TFile*
infile =
new TFile(WeightFileName.c_str());
383 TH1F* WHist = (TH1F*)
infile->Get(
"WHist");
387 std::cout <<
" Could not find the histogram WHist in the file " 388 <<
"in the file " << WeightFileName <<
"." << std::endl;
392 for (
int i = 0;
i < 50;
i++) {
393 for (
int j = 0;
j < 50;
j++) {
394 for (
int k = 0;
k < 50;
k++) {
400 std::cout <<
" 3D Weight Matrix initialized! " << std::endl;
413 static const double weight_24[25] = {0, 0, 0, 0, 2.46277e-06,
414 2.95532e-05, 0.000104668, 0.000401431, 0.00130034, 0.00342202,
415 0.00818132, 0.0175534, 0.035784, 0.0650836, 0.112232,
416 0.178699, 0.268934, 0.380868, 0.507505, 0.640922,
417 0.768551, 0.877829, 0.958624, 0.99939, 1};
419 static const double weight_23[25] = {0, 1.20628e-06, 1.20628e-06, 2.41255e-06, 1.20628e-05,
420 6.39326e-05, 0.000252112, 0.000862487, 0.00244995, 0.00616527,
421 0.0140821, 0.0293342, 0.0564501, 0.100602, 0.164479,
422 0.252659, 0.36268, 0.491427, 0.627979, 0.75918,
423 0.873185, 0.957934, 0.999381, 1, 0.957738};
425 static const double weight_22[25] = {
426 0, 0, 0, 5.88636e-06, 3.0609e-05, 0.000143627, 0.000561558, 0.00173059, 0.00460078,
427 0.0110616, 0.0238974, 0.0475406, 0.0875077, 0.148682, 0.235752, 0.343591, 0.473146, 0.611897,
428 0.748345, 0.865978, 0.953199, 0.997848, 1, 0.954245, 0.873688};
430 static const double weight_21[25] = {
431 0, 0, 1.15381e-06, 8.07665e-06, 7.1536e-05, 0.000280375, 0.00107189, 0.00327104, 0.00809396,
432 0.0190978, 0.0401894, 0.0761028, 0.13472, 0.216315, 0.324649, 0.455125, 0.598241, 0.739215,
433 0.861866, 0.953911, 0.998918, 1, 0.956683, 0.872272, 0.76399};
435 static const double weight_20[25] = {
436 0, 0, 1.12532e-06, 2.58822e-05, 0.000145166, 0.000633552, 0.00215048, 0.00592816, 0.0145605,
437 0.0328367, 0.0652649, 0.11893, 0.19803, 0.305525, 0.436588, 0.581566, 0.727048, 0.8534,
438 0.949419, 0.999785, 1, 0.953008, 0.865689, 0.753288, 0.62765};
439 static const double weight_19[25] = {
440 0, 0, 1.20714e-05, 5.92596e-05, 0.000364337, 0.00124994, 0.00403953, 0.0108149, 0.025824,
441 0.0544969, 0.103567, 0.17936, 0.283532, 0.416091, 0.562078, 0.714714, 0.846523, 0.947875,
442 1, 0.999448, 0.951404, 0.859717, 0.742319, 0.613601, 0.48552};
444 static const double weight_18[25] = {
445 0, 3.20101e-06, 2.88091e-05, 0.000164319, 0.000719161, 0.00250106, 0.00773685, 0.0197513, 0.0443693,
446 0.0885998, 0.159891, 0.262607, 0.392327, 0.543125, 0.69924, 0.837474, 0.943486, 0.998029,
447 1, 0.945937, 0.851807, 0.729309, 0.596332, 0.467818, 0.350434};
449 static const double weight_17[25] = {
450 1.03634e-06, 7.25437e-06, 4.97443e-05, 0.000340956, 0.00148715, 0.00501485, 0.0143067, 0.034679, 0.0742009,
451 0.140287, 0.238288, 0.369416, 0.521637, 0.682368, 0.828634, 0.939655, 1, 0.996829,
452 0.94062, 0.841575, 0.716664, 0.582053, 0.449595, 0.331336, 0.234332};
454 static const double weight_16[25] = {
455 4.03159e-06, 2.41895e-05, 0.000141106, 0.00081942, 0.00314565, 0.00990662, 0.026293, 0.0603881, 0.120973,
456 0.214532, 0.343708, 0.501141, 0.665978, 0.820107, 0.938149, 1, 0.99941, 0.940768,
457 0.837813, 0.703086, 0.564023, 0.42928, 0.312515, 0.216251, 0.14561};
459 static const double weight_15[25] = {
460 9.76084e-07, 5.07564e-05, 0.000303562, 0.00174036, 0.00617959, 0.0188579, 0.047465, 0.101656, 0.189492,
461 0.315673, 0.474383, 0.646828, 0.809462, 0.934107, 0.998874, 1, 0.936163, 0.827473,
462 0.689675, 0.544384, 0.40907, 0.290648, 0.198861, 0.12951, 0.0808051};
464 static const double weight_14[25] = {
465 1.13288e-05, 0.000124617, 0.000753365, 0.00345056, 0.0123909, 0.0352712, 0.0825463, 0.16413, 0.287213,
466 0.44615, 0.625826, 0.796365, 0.930624, 0.999958, 1, 0.934414, 0.816456, 0.672939,
467 0.523033, 0.386068, 0.269824, 0.180342, 0.114669, 0.0698288, 0.0406496};
469 static const double weight_13[25] = {
470 2.54296e-05, 0.000261561, 0.00167018, 0.00748083, 0.0241308, 0.0636801, 0.138222, 0.255814, 0.414275,
471 0.600244, 0.779958, 0.92256, 0.999155, 1, 0.927126, 0.804504, 0.651803, 0.497534,
472 0.35976, 0.245834, 0.160904, 0.0991589, 0.0585434, 0.0332437, 0.0180159};
474 static const double weight_12[25] = {
475 5.85742e-05, 0.000627706, 0.00386677, 0.0154068, 0.0465892, 0.111683, 0.222487, 0.381677, 0.5719,
476 0.765001, 0.915916, 1, 0.999717, 0.921443, 0.791958, 0.632344, 0.475195, 0.334982,
477 0.223666, 0.141781, 0.0851538, 0.048433, 0.0263287, 0.0133969, 0.00696683};
479 static const double weight_11[25] = {
480 0.00015238, 0.00156064, 0.00846044, 0.0310939, 0.0856225, 0.187589, 0.343579, 0.541892, 0.74224,
481 0.909269, 0.998711, 1, 0.916889, 0.77485, 0.608819, 0.447016, 0.307375, 0.198444,
482 0.121208, 0.070222, 0.0386492, 0.0201108, 0.0100922, 0.00484937, 0.00222458};
484 static const double weight_10[25] = {
485 0.000393044, 0.00367001, 0.0179474, 0.060389, 0.151477, 0.302077, 0.503113, 0.720373, 0.899568,
486 1, 0.997739, 0.909409, 0.75728, 0.582031, 0.415322, 0.277663, 0.174147, 0.102154,
487 0.0566719, 0.0298642, 0.0147751, 0.00710995, 0.00319628, 0.00140601, 0.000568796};
489 static const double weight_9[25] = {
490 0.00093396, 0.00854448, 0.0380306, 0.113181, 0.256614, 0.460894, 0.690242, 0.888781, 1,
491 0.998756, 0.899872, 0.735642, 0.552532, 0.382726, 0.246114, 0.147497, 0.0825541, 0.0441199,
492 0.0218157, 0.0103578, 0.00462959, 0.0019142, 0.000771598, 0.000295893, 0.000111529};
494 static const double weight_8[25] = {
495 0.00240233, 0.0192688, 0.0768653, 0.205008, 0.410958, 0.65758, 0.875657, 0.999886, 1,
496 0.889476, 0.711446, 0.517781, 0.345774, 0.212028, 0.121208, 0.0644629, 0.0324928, 0.0152492,
497 0.00673527, 0.0028547, 0.00117213, 0.000440177, 0.000168471, 5.80689e-05, 1.93563e-05};
499 static const double weight_7[25] = {0.00617233, 0.0428714, 0.150018, 0.350317, 0.612535,
500 0.856525, 0.999923, 1, 0.87544, 0.679383,
501 0.478345, 0.303378, 0.176923, 0.0950103, 0.0476253,
502 0.0222211, 0.00972738, 0.00392962, 0.0015258, 0.000559168,
503 0.000183928, 6.77983e-05, 1.67818e-05, 7.38398e-06, 6.71271e-07};
505 static const double weight_6[25] = {0.0154465, 0.0923472, 0.277322, 0.55552, 0.833099,
506 0.999035, 1, 0.855183, 0.641976, 0.428277,
507 0.256804, 0.139798, 0.0700072, 0.0321586, 0.0137971,
508 0.00544756, 0.00202316, 0.000766228, 0.000259348, 8.45836e-05,
509 1.80362e-05, 8.70713e-06, 3.73163e-06, 6.21938e-07, 0};
511 static const double weight_5[25] = {0.0382845, 0.191122, 0.478782, 0.797314, 1,
512 0.997148, 0.831144, 0.59461, 0.371293, 0.205903,
513 0.103102, 0.0471424, 0.0194997, 0.00749415, 0.00273709,
514 0.000879189, 0.000286049, 0.000102364, 1.70606e-05, 3.98081e-06,
515 2.27475e-06, 0, 0, 0, 0};
517 static const double weight_4[25] = {0.0941305, 0.373824, 0.750094, 1, 0.997698,
518 0.800956, 0.532306, 0.304597, 0.152207, 0.0676275,
519 0.0270646, 0.00975365, 0.00326077, 0.00101071, 0.000301781,
520 7.41664e-05, 1.58563e-05, 3.58045e-06, 1.02299e-06, 0,
521 5.11493e-07, 0, 0, 0, 0};
523 static const double weight_3[25] = {0.222714, 0.667015, 1, 0.999208, 0.750609,
524 0.449854, 0.224968, 0.0965185, 0.0361225, 0.012084,
525 0.00359618, 0.000977166, 0.000239269, 6.29422e-05, 1.16064e-05,
526 1.78559e-06, 0, 4.46398e-07, 0, 0,
529 static const double weight_2[25] = {
530 0.499541, 0.999607, 1, 0.666607, 0.333301, 0.13279, 0.0441871, 0.0127455, 0.00318434,
531 0.00071752, 0.000132204, 2.69578e-05, 5.16999e-06, 2.21571e-06, 0, 0, 0, 0,
532 0, 0, 0, 0, 0, 0, 0};
534 static const double weight_1[25] = {
535 0.999165, 1, 0.499996, 0.166868, 0.0414266, 0.00831053, 0.00137472, 0.000198911, 2.66302e-05,
536 2.44563e-06, 2.71737e-07, 2.71737e-07, 0, 0, 0, 0, 0, 0,
537 0, 0, 0, 0, 0, 0, 0};
539 static const double weight_0[25] = {1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
543 const double* WeightPtr =
nullptr;
545 for (
int iint = 0; iint < 25; ++iint) {
547 WeightPtr = weight_0;
549 WeightPtr = weight_1;
551 WeightPtr = weight_2;
553 WeightPtr = weight_3;
555 WeightPtr = weight_4;
557 WeightPtr = weight_5;
559 WeightPtr = weight_6;
561 WeightPtr = weight_7;
563 WeightPtr = weight_8;
565 WeightPtr = weight_9;
567 WeightPtr = weight_10;
569 WeightPtr = weight_11;
571 WeightPtr = weight_12;
573 WeightPtr = weight_13;
575 WeightPtr = weight_14;
577 WeightPtr = weight_15;
579 WeightPtr = weight_16;
581 WeightPtr = weight_17;
583 WeightPtr = weight_18;
585 WeightPtr = weight_19;
587 WeightPtr = weight_20;
589 WeightPtr = weight_21;
591 WeightPtr = weight_22;
593 WeightPtr = weight_23;
595 WeightPtr = weight_24;
597 for (
int ibin = 0; ibin < 25; ++ibin) {
621 int npm1 =
min(pv1, 34);
622 int np0 =
min(pv2, 34);
623 int npp1 =
min(pv3, 34);
629 static const double Correct_Weights2011[25] = {
631 5.30031, 2.07903, 1.40729, 1.27687, 1.0702, 0.902094, 0.902345, 0.931449, 0.78202,
632 0.824686, 0.837735, 0.910261, 1.01394, 1.1599, 1.12778, 1.58423, 1.78868, 1.58296,
633 2.3291, 3.86641, 0, 0, 0, 0, 0};
636 std::cout <<
" **** Warning: Out-of-time pileup reweighting appropriate only for PU_S3 **** " << std::endl;
637 std::cout <<
" **** will be applied **** " << std::endl;
646 if (npv_in_time < 0) {
647 std::cerr <<
" no in-time beam crossing found\n! ";
648 std::cerr <<
" Returning event weight=0\n! ";
651 if (npv_m50nsBX < 0) {
652 std::cerr <<
" no out-of-time beam crossing found\n! ";
653 std::cerr <<
" Returning event weight=0\n! ";
657 int bin =
weights_->GetXaxis()->FindBin(npv_in_time);
659 double inTimeWeight =
weights_->GetBinContent(
bin);
661 double TotalWeight = 1.0;
663 TotalWeight = inTimeWeight *
WeightOOTPU_[
bin - 1][npv_m50nsBX] * Correct_Weights2011[
bin - 1];
double WeightOOTPU_[25][25]
LumiReWeighting(std::string generatedFile, std::string dataFile, std::string GenHistName, std::string DataHistName)
void weight3D_init(float ScaleFactor, std::string WeightOutputFile="")
double weight(float n_int)
double ShiftWeight(int ibin)
void weight3D_set(std::string WeightFileName)
LumiReWeighting(const std::vector< float > &MC_distr, const std::vector< float > &Lumi_distr)
std::string DataHistName_
double ShiftWeight(float pvnum)
double weightOOT(int npv_in_time, int npv_m50nsBX)
std::string dataFileName_
double Weight3D_[50][50][50]
double ITweight3BX(float ave_int)
double weight3D(int pv1, int pv2, int pv3)
PoissonMeanShifter(float Shift)
int factorial(int n)
factorial function
std::string generatedFileName_