CMS 3D CMS Logo

L1ValidatorHists.cc
Go to the documentation of this file.
2 
3 //#include <DataFormats/HepMCCandidate/interface/GenParticle.h>
4 
6 
7 /*#define BOOKHISTS(TYPE) \
8 TYPE ## _N_Pt = new TH2F(#TYPE "_N_Pt", #TYPE " Number", 20, 0, 200); \
9 TYPE ## _N_Eta = new TH2F(#TYPE "_N_Eta", #TYPE " Number", 20, -4, 4); \
10 TYPE ## _Eff_Pt = new TH2F(#TYPE "_Eff_Pt", #TYPE " Number", 20, 0, 200); \
11 TYPE ## _Eff_Eta = new TH2F(#TYPE "_Eff_Eta", #TYPE " Number", 20, -4, 4); \
12 TYPE ## _dR = new TH2F(#TYPE "_dR", #TYPE " Number", 20, 0, 1); \
13 TYPE ## _dPt = new TH2F(#TYPE "_dPt", #TYPE " Number", 20, -1, 1);
14 */
16 // Name[0]="IsoEG"; // Run I legacy
17 // Name[1]="NonIsoEG";
18 // Name[2]="CenJet";
19 // Name[3]="ForJet";
20 // Name[4]="TauJet";
21 // Name[5]="Muon";
22  Name[0]="Egamma";
23  Name[1]="Jet";
24  Name[2]="Tau";
25  Name[3]="Muon";
26 
27 }
29 
31  NEvents=0;
32 
33  float ptbins[14] = { 0,5,10,15,20,25,30,35, 40, 50, 60, 80, 120, 160};
34  int Nptbin = 13;
35 
36  for(int i=0; i<Type::Number; i++){
37  N[i] = iBooker.book1D( (Name[i]+"_N").c_str(), ("L1 " + Name[i]+" Number with BX=0").c_str(), 16, -0.5, 15.5);
38 
39  Eff_Pt[i] = iBooker.book1D( (Name[i]+"_Eff_Pt").c_str(), (Name[i]+" Efficiency vs Pt; Gen p_{T} [GeV]; L1T Efficiency").c_str(), Nptbin, ptbins);
40  Eff_Pt_Denom[i] = iBooker.book1D( (Name[i]+"_Eff_Pt_Denom").c_str(), (Name[i]+" Efficiency vs Pt Denom; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
41  Eff_Pt_Nomin[i] = iBooker.book1D( (Name[i]+"_Eff_Pt_Nomin").c_str(), (Name[i]+" Efficiency vs Pt Nomin; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
42  Eff_Eta[i] = iBooker.book1D( (Name[i]+"_Eff_Eta").c_str(), (Name[i]+" Efficiency vs #eta (Gen p_{T} > 10GeV); Gen #eta; L1T Efficiency").c_str(), 80, -4, 4);
43  Eff_Eta_Denom[i] = iBooker.book1D( (Name[i]+"_Eff_Eta_Denom").c_str(), (Name[i]+" Efficiency vs #eta Denom; Gen #eta; Entries").c_str(), 80, -4, 4);
44  Eff_Eta_Nomin[i] = iBooker.book1D( (Name[i]+"_Eff_Eta_Nomin").c_str(), (Name[i]+" Efficiency vs #eta Nomin; Gen #eta; Entries").c_str(), 80, -4, 4);
45  TurnOn_15[i] = iBooker.book1D( (Name[i]+"_TurnOn_15").c_str(), (Name[i]+" Turn On (15 GeV); Gen p_{T} [GeV]; L1T Efficiency").c_str(), Nptbin, ptbins);
46  TurnOn_15_Denom[i] = iBooker.book1D( (Name[i]+"_TurnOn_15_Denom").c_str(), (Name[i]+" Turn On (15 GeV) Denom; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
47  TurnOn_15_Nomin[i] = iBooker.book1D( (Name[i]+"_TurnOn_15_Nomin").c_str(), (Name[i]+" Turn On (15 GeV) Nomin; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
48  TurnOn_30[i] = iBooker.book1D( (Name[i]+"_TurnOn_30").c_str(), (Name[i]+" Turn On (30 GeV); Gen p_{T} [GeV]; L1T Efficiency").c_str(), Nptbin, ptbins);
49  TurnOn_30_Denom[i] = iBooker.book1D( (Name[i]+"_TurnOn_30_Denom").c_str(), (Name[i]+" Turn On (30 GeV) Denom; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
50  TurnOn_30_Nomin[i] = iBooker.book1D( (Name[i]+"_TurnOn_30_Nomin").c_str(), (Name[i]+" Turn On (30 GeV) Nomin; Gen p_{T} [GeV]; Entries").c_str(), Nptbin, ptbins);
51  dR[i] = iBooker.book1D( (Name[i]+"_dR").c_str(), (Name[i]+" #DeltaR; #DeltaR(L1 object, Gen object); Entries").c_str(), 40, 0, 1);
52  dR_vs_Pt[i] = iBooker.book2D( (Name[i]+"_dR_vs_Pt").c_str(), (Name[i]+" #DeltaR vs p_{T}; Gen p_{T} [GeV]; #DeltaR(L1 object, Gen object); Entries").c_str(), 12, 0, 120, 40, 0, 1);
53  dPt[i] = iBooker.book1D( (Name[i]+"_dPt").c_str(), (Name[i]+" #Deltap_{T}; (p_{T}^{L1}-p_{T}^{Gen})/p_{T}^{Gen}; Entries").c_str(), 100, -2, 2);
54  dPt_vs_Pt[i] = iBooker.book2D( (Name[i]+"_dPt_vs_Pt").c_str(), (Name[i]+" #Deltap_{T} vs p_{T}; Gen p_{T} [GeV]; (p_{T}^{L1}-p_{T}^{Gen})/p_{T}^{Gen}; Entries").c_str(), 12, 0, 120, 40, -2, 2);
55  }
56 
57 }
58 
60  double GenPartPt = GenPart->pt();
61  // fill the overflow in the last bin
62  if(GenPart->pt()>=160.0) GenPartPt = 159.0;
63  if(L1Part==nullptr) {
64  Eff_Pt_Denom[i]->Fill(GenPartPt);
65  if(GenPart->pt()>10)Eff_Eta_Denom[i]->Fill(GenPart->eta());
66  TurnOn_15_Denom[i]->Fill(GenPartPt);
67  TurnOn_30_Denom[i]->Fill(GenPartPt);
68  } else {
69  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), L1Part->eta(), L1Part->phi());
70  bool matched = idR < 0.15;
71  Eff_Pt_Denom[i]->Fill(GenPartPt);
72  if(GenPart->pt()>10)Eff_Eta_Denom[i]->Fill(GenPart->eta());
73  if(matched)Eff_Pt_Nomin[i]->Fill(GenPartPt);
74  if(matched && GenPart->pt()>10)Eff_Eta_Nomin[i]->Fill(GenPart->eta());
75  TurnOn_15_Denom[i]->Fill(GenPartPt);
76  TurnOn_30_Denom[i]->Fill(GenPartPt);
77  if(L1Part->pt()>15 && matched) TurnOn_15_Nomin[i]->Fill(GenPartPt);
78  if(L1Part->pt()>30 && matched) TurnOn_30_Nomin[i]->Fill(GenPartPt);
79  dR[i]->Fill(idR);
80  dPt[i]->Fill( (L1Part->pt()-GenPart->pt()) / GenPart->pt() );
81  dR_vs_Pt[i]->Fill(GenPart->pt(), idR);
82  dPt_vs_Pt[i]->Fill( GenPart->pt(), (L1Part->pt()-GenPart->pt()) / GenPart->pt() );
83  }
84 }
85 
87  N[i]->Fill(Number);
88 }
89 
91  for(int i=0; i<Type::Number; i++){
92  N[i]->getTH1()->Write();
93  Eff_Pt[i]->getTH1()->Write();
94  Eff_Pt_Denom[i]->getTH1()->Write();
95  Eff_Pt_Nomin[i]->getTH1()->Write();
96  Eff_Eta[i]->getTH1()->Write();
97  Eff_Eta_Denom[i]->getTH1()->Write();
98  Eff_Eta_Nomin[i]->getTH1()->Write();
99  TurnOn_15[i]->getTH1()->Write();
100  TurnOn_15_Denom[i]->getTH1()->Write();
101  TurnOn_15_Nomin[i]->getTH1()->Write();
102  TurnOn_30[i]->getTH1()->Write();
103  TurnOn_30_Denom[i]->getTH1()->Write();
104  TurnOn_30_Nomin[i]->getTH1()->Write();
105  dR[i]->getTH1()->Write();
106  dPt[i]->getTH1()->Write();
107  dR_vs_Pt[i]->getTH2F()->Write();
108  dPt_vs_Pt[i]->getTH2F()->Write();
109  }
110 }
111 
112 /*void L1ValidatorHists::NormalizeSlices(TH2F *Hist){
113  int NBinsX = Hist->GetNbinsX();
114  int NBinsY = Hist->GetNbinsY();
115  for(int i=0; i<NBinsX+2; i++){
116  float Total = Hist->Integral(i, i, 0, -1);
117  if(Total == 0) continue;
118  for(int j=0; j<NBinsY+2; j++){
119  Hist->SetBinContent(i,j, Hist->GetBinContent(i,j)/Total);
120  }
121  }
122 }
123 */
void Fill(int, const reco::LeafCandidate *, const reco::LeafCandidate *)
MonitorElement * dR_vs_Pt[Type::Number]
double eta() const final
momentum pseudorapidity
MonitorElement * Eff_Eta[Type::Number]
MonitorElement * N[Type::Number]
MonitorElement * TurnOn_15[Type::Number]
MonitorElement * Eff_Pt_Denom[Type::Number]
MonitorElement * TurnOn_15_Denom[Type::Number]
double pt() const final
transverse momentum
MonitorElement * Eff_Eta_Nomin[Type::Number]
void Fill(long long x)
void FillNumber(int, int)
MonitorElement * TurnOn_15_Nomin[Type::Number]
MonitorElement * TurnOn_30_Denom[Type::Number]
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * dPt[Type::Number]
TH1 * getTH1(void) const
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
MonitorElement * TurnOn_30[Type::Number]
MonitorElement * dR[Type::Number]
MonitorElement * TurnOn_30_Nomin[Type::Number]
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * dPt_vs_Pt[Type::Number]
MonitorElement * Eff_Pt_Nomin[Type::Number]
void Book(DQMStore::IBooker &)
TH2F * getTH2F(void) const
MonitorElement * Eff_Eta_Denom[Type::Number]
double phi() const final
momentum azimuthal angle
std::string Name[4]
MonitorElement * Eff_Pt[Type::Number]