1 #ifndef PhysicsTools_Utilities_interface_Lumi3DReWeighting_cc
2 #define PhysicsTools_Utilities_interface_Lumi3DReWeighting_cc
25 #include "TStopwatch.h"
45 : generatedFileName_(generatedFile),
46 dataFileName_(dataFile),
47 GenHistName_(GenHistName),
48 DataHistName_(DataHistName),
49 weightFileName_(WeightOutputFile) {
53 std::shared_ptr<TH1> Data_temp =
56 std::shared_ptr<TH1> MC_temp =
71 const std::vector<float> &Lumi_distr,
79 Int_t NMCBins = MC_distr.size();
81 MC_distr_ = std::shared_ptr<TH1>(
new TH1F(
"MC_distr",
"MC dist", NMCBins, 0.,
float(NMCBins)));
83 Int_t NDBins = Lumi_distr.size();
85 Data_distr_ = std::shared_ptr<TH1>(
new TH1F(
"Data_distr",
"Data dist", NDBins, 0.,
float(NDBins)));
87 for (
int ibin = 1; ibin < NMCBins + 1; ++ibin) {
88 MC_distr_->SetBinContent(ibin, MC_distr[ibin - 1]);
91 for (
int ibin = 1; ibin < NDBins + 1; ++ibin) {
92 Data_distr_->SetBinContent(ibin, Lumi_distr[ibin - 1]);
98 if (fabs(1.0 - deltaH) > 0.001) {
102 if (fabs(1.0 - deltaMC) > 0.001) {
110 int npm1 =
min(pv1, 49);
111 int np0 =
min(pv2, 49);
112 int npp1 =
min(pv3, 49);
127 std::vector<PileupSummaryInfo>::const_iterator PVI;
133 for (PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
134 int BX = PVI->getBunchCrossing();
137 npm1 = PVI->getPU_NumInteractions();
140 np0 = PVI->getPU_NumInteractions();
143 npp1 = PVI->getPU_NumInteractions();
147 npm1 =
min(npm1, 49);
149 npp1 =
min(npp1, 49);
173 std::shared_ptr<TH1> Data_temp =
176 std::shared_ptr<TH1> MC_temp =
195 TH3D *WHist =
new TH3D(
"WHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
196 TH3D *DHist =
new TH3D(
"DHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
197 TH3D *MHist =
new TH3D(
"MHist",
"3D weights", 50, 0., 50., 50, 0., 50., 50, 0., 50.);
202 std::cout <<
" MC and Data distributions are not initialized! You must call the Lumi3DReWeighting constructor. "
208 double MC_ints[50][50][50];
209 double Data_ints[50][50][50];
211 for (
int i = 0;
i < 50;
i++) {
212 for (
int j = 0;
j < 50;
j++) {
213 for (
int k = 0;
k < 50;
k++) {
214 MC_ints[
i][
j][
k] = 0.;
215 Data_ints[
i][
j][
k] = 0.;
227 for (
int i = 1;
i < 50; ++
i) {
234 double probi, probj, probk;
242 for (
int jbin = 1; jbin < NMCbin + 1; jbin++) {
244 xweight =
MC_distr_->GetBinContent(jbin);
254 throw cms::Exception(
"BadInputValue") <<
" Your histogram generates MC luminosity values less than zero!"
255 <<
" Please Check. Terminating." << std::endl;
266 for (
int i = 1;
i < 50; ++
i) {
273 for (
int i = 0;
i < 50;
i++) {
275 for (
int j = 0;
j < 50;
j++) {
277 for (
int k = 0;
k < 50;
k++) {
280 MC_ints[
i][
j][
k] = MC_ints[
i][
j][
k] + probi * probj * probk * xweight;
288 for (
int jbin = 1; jbin < NDatabin + 1; jbin++) {
295 throw cms::Exception(
"BadInputValue") <<
" Your histogram generates Data luminosity values less than zero!"
296 <<
" Please Check. Terminating." << std::endl;
307 for (
int i = 1;
i < 50; ++
i) {
314 for (
int i = 0;
i < 50;
i++) {
316 for (
int j = 0;
j < 50;
j++) {
318 for (
int k = 0;
k < 50;
k++) {
321 Data_ints[
i][
j][
k] = Data_ints[
i][
j][
k] + probi * probj * probk * xweight;
327 for (
int i = 0;
i < 50;
i++) {
329 for (
int j = 0;
j < 50;
j++) {
330 for (
int k = 0;
k < 50;
k++) {
331 if ((MC_ints[
i][
j][
k]) > 0.) {
337 DHist->SetBinContent(
i + 1,
j + 1,
k + 1, Data_ints[
i][
j][
k]);
338 MHist->SetBinContent(
i + 1,
j + 1,
k + 1, MC_ints[
i][
j][
k]);
346 std::cout <<
" 3D Weight Matrix initialized! " << std::endl;
362 TFile *
infile =
new TFile(WeightFileName.c_str());
363 TH1F *WHist = (TH1F *)
infile->Get(
"WHist");
367 throw cms::Exception(
"HistogramNotFound") <<
" Could not find the histogram WHist in the file "
368 <<
"in the file " << WeightFileName <<
"." << std::endl;
371 for (
int i = 0;
i < 50;
i++) {
373 for (
int j = 0;
j < 50;
j++) {
374 for (
int k = 0;
k < 50;
k++) {
382 std::cout <<
" 3D Weight Matrix initialized! " << std::endl;
388 TFile *infileMC =
new TFile(MCWeightFileName.c_str());
389 TH3D *MHist = (TH3D *)infileMC->Get(
"MHist");
393 throw cms::Exception(
"HistogramNotFound") <<
" Could not find the histogram MHist in the file "
394 <<
"in the file " << MCWeightFileName <<
"." << std::endl;
397 TFile *infileD =
new TFile(DataWeightFileName.c_str());
398 TH3D *DHist = (TH3D *)infileD->Get(
"DHist");
402 throw cms::Exception(
"HistogramNotFound") <<
" Could not find the histogram DHist in the file "
403 <<
"in the file " << DataWeightFileName <<
"." << std::endl;
406 for (
int i = 0;
i < 50;
i++) {
407 for (
int j = 0;
j < 50;
j++) {
408 for (
int k = 0;
k < 50;
k++) {
409 Weight3D_[
i][
j][
k] = DHist->GetBinContent(
i + 1,
j + 1,
k + 1) / MHist->GetBinContent(
i + 1,
j + 1,
k + 1);
414 std::cout <<
" 3D Weight Matrix initialized! " << std::endl;