CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  _dbe=dbe;
17  //_dbe->setCurrentFolder("L1T/L1Extra");
18 
19  Name[0]="IsoEG";
20  Name[1]="NonIsoEG";
21  Name[2]="CenJet";
22  Name[3]="ForJet";
23  Name[4]="TauJet";
24  Name[5]="Muon";
25 }
26 
28  NEvents=0;
29 
30  for(int i=0; i<Type::Number; i++){
31  N[i] = _dbe->book1D( (Name[i]+"_N").c_str(), (Name[i]+" Number").c_str(), 5, -0.5, 4.5);
32  N_Pt[i] = new TH1F( (Name[i]+"_N_Pt").c_str(), (Name[i]+" Number").c_str(), 20, 0, 100);
33  N_Eta[i] = new TH1F( (Name[i]+"_N_Eta").c_str(), (Name[i]+" Number").c_str(), 20, -4, 4);
34  Eff_Pt[i] = _dbe->book1D( (Name[i]+"_Eff_Pt").c_str(), (Name[i]+" Efficiency").c_str(), 20, 0, 100);
35  Eff_Eta[i] = _dbe->book1D( (Name[i]+"_Eff_Eta").c_str(), (Name[i]+" Efficiency").c_str(), 20, -4, 4);
36  TurnOn_15[i] = _dbe->book1D( (Name[i]+"_TurnOn_15").c_str(), (Name[i]+" Turn On (15 GeV)").c_str(), 20, 0, 100);
37  TurnOn_30[i] = _dbe->book1D( (Name[i]+"_TurnOn_30").c_str(), (Name[i]+" Turn On (30 GeV)").c_str(), 20, 0, 100);
38  dR[i] = _dbe->book1D( (Name[i]+"_dR").c_str(), (Name[i]+" dR").c_str(), 20, 0, 1);
39  dPt[i] = _dbe->book1D( (Name[i]+"_dPt").c_str(), (Name[i]+" dPt").c_str(), 20, -1, 1);
40 
41 // N[i]->getTH1()->Sumw2();
42 // N_Pt[i]->Sumw2();
43 // N_Eta[i]->Sumw2();
44 // Eff_Pt[i]->getTH1()->Sumw2();
45 // Eff_Eta[i]->getTH1()->Sumw2();
46 // TurnOn_15[i]->getTH1()->Sumw2();
47 // TurnOn_30[i]->getTH1()->Sumw2();
48 // dR[i]->getTH1()->Sumw2();
49 // dPt[i]->getTH1()->Sumw2();
50  }
51 
52 }
53 
54 void L1ValidatorHists::Fill(int i, const reco::LeafCandidate *GenPart, const reco::LeafCandidate *RecoPart){
55  N_Pt[i]->Fill(GenPart->pt());
56  N_Eta[i]->Fill(GenPart->eta());
57 
58  if(RecoPart==NULL) return;
59 
60  Eff_Pt[i]->Fill(GenPart->pt());
61  if(RecoPart->pt()>15) TurnOn_15[i]->Fill(GenPart->pt());
62  if(RecoPart->pt()>30) TurnOn_30[i]->Fill(GenPart->pt());
63  Eff_Eta[i]->Fill(GenPart->eta());
64  dR[i]->Fill(reco::deltaR(GenPart->eta(), GenPart->phi(), RecoPart->eta(), RecoPart->phi()));
65  dPt[i]->Fill( (RecoPart->pt()-GenPart->pt()) / GenPart->pt() );
66 }
67 
68 void L1ValidatorHists::FillNumber(int i, int Number){
69  N[i]->Fill(Number);
70 }
71 
73  for(int i=0; i<Type::Number; i++){
74  Eff_Pt[i]->getTH1F()->Divide(N_Pt[i]);
75  TurnOn_15[i]->getTH1F()->Divide(N_Pt[i]);
76  TurnOn_30[i]->getTH1F()->Divide(N_Pt[i]);
77  Eff_Eta[i]->getTH1F()->Divide(N_Eta[i]);
78  N[i]->getTH1()->Scale(1./N[i]->getEntries());
79  dR[i]->getTH1()->Scale(1./dR[i]->getEntries());
80  dPt[i]->getTH1()->Scale(1./dPt[i]->getEntries());
81  }
82 }
83 
85  for(int i=0; i<Type::Number; i++){
86  N[i]->getTH1()->Write();
87  Eff_Pt[i]->getTH1()->Write();
88  Eff_Eta[i]->getTH1()->Write();
89  TurnOn_15[i]->getTH1()->Write();
90  TurnOn_30[i]->getTH1()->Write();
91  dR[i]->getTH1()->Write();
92  dPt[i]->getTH1()->Write();
93  }
94 }
95 
96 /*void L1ValidatorHists::NormalizeSlices(TH2F *Hist){
97  int NBinsX = Hist->GetNbinsX();
98  int NBinsY = Hist->GetNbinsY();
99  for(int i=0; i<NBinsX+2; i++){
100  float Total = Hist->Integral(i, i, 0, -1);
101  if(Total == 0) continue;
102  for(int j=0; j<NBinsY+2; j++){
103  Hist->SetBinContent(i,j, Hist->GetBinContent(i,j)/Total);
104  }
105  }
106 }
107 */
void Fill(int, const reco::LeafCandidate *, const reco::LeafCandidate *)
int i
Definition: DBlmapReader.cc:9
TH1F * N_Pt[Type::Number]
MonitorElement * Eff_Eta[Type::Number]
MonitorElement * N[Type::Number]
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:964
std::string Name[6]
TH1F * N_Eta[Type::Number]
MonitorElement * TurnOn_15[Type::Number]
#define NULL
Definition: scimark2.h:8
void Fill(long long x)
void FillNumber(int, int)
L1ValidatorHists(DQMStore *dbe)
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
MonitorElement * dPt[Type::Number]
TH1 * getTH1(void) const
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
MonitorElement * TurnOn_30[Type::Number]
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
MonitorElement * dR[Type::Number]
TH1F * getTH1F(void) const
virtual float pt() const GCC11_FINAL
transverse momentum
MonitorElement * Eff_Pt[Type::Number]