1 |
loizides |
1.2 |
// $Id: JetValidationMod.cc,v 1.1 2008/10/15 06:05:02 loizides Exp $
|
2 |
loizides |
1.1 |
|
3 |
|
|
#include "MitPhysics/Validation/interface/JetValidationMod.h"
|
4 |
|
|
#include <TH1D.h>
|
5 |
|
|
#include <TH2D.h>
|
6 |
|
|
#include <TH3D.h>
|
7 |
|
|
#include "MitAna/DataTree/interface/Names.h"
|
8 |
|
|
#include "MitAna/DataCont/interface/ObjArray.h"
|
9 |
|
|
#include "MitCommon/MathTools/interface/MathUtils.h"
|
10 |
|
|
#include <map>
|
11 |
|
|
|
12 |
|
|
using namespace mithep;
|
13 |
|
|
|
14 |
|
|
ClassImp(mithep::JetValidationMod)
|
15 |
|
|
|
16 |
|
|
//--------------------------------------------------------------------------------------------------
|
17 |
|
|
JetValidationMod::JetValidationMod(const char *name, const char *title) :
|
18 |
|
|
BaseMod(name,title),
|
19 |
|
|
fPrintDebug(false),
|
20 |
|
|
fIC5GenJetName("IC5GenJets"),
|
21 |
|
|
fSC5GenJetName("SC5GenJets"),
|
22 |
|
|
fIC5JetName(Names::gkCaloJetBrn),
|
23 |
|
|
fSC5JetName("SisCone5Jets"),
|
24 |
|
|
fIC5Jets(0),
|
25 |
|
|
fSC5Jets(0),
|
26 |
|
|
fIC5GenJets(0),
|
27 |
loizides |
1.2 |
fSC5GenJets(0)
|
28 |
loizides |
1.1 |
{
|
29 |
|
|
// Constructor.
|
30 |
|
|
}
|
31 |
|
|
|
32 |
|
|
//--------------------------------------------------------------------------------------------------
|
33 |
|
|
void JetValidationMod::Process()
|
34 |
|
|
{
|
35 |
|
|
// Process entries of the tree.
|
36 |
|
|
|
37 |
|
|
//************************************************************************************************
|
38 |
|
|
//
|
39 |
|
|
// For IterativeCone 0.5 Jets
|
40 |
|
|
//
|
41 |
|
|
//************************************************************************************************
|
42 |
|
|
|
43 |
|
|
ObjArray<Jet> *GoodIC5Jets = new ObjArray<Jet>;
|
44 |
|
|
LoadBranch(fIC5JetName);
|
45 |
|
|
for (UInt_t i=0; i<fIC5Jets->GetEntries(); ++i) {
|
46 |
|
|
if(fabs(fIC5Jets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
|
47 |
|
|
if(fIC5Jets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
|
48 |
|
|
GoodIC5Jets->Add(fIC5Jets->At(i));
|
49 |
|
|
fIC5CaloJetsPt->Fill(fIC5Jets->At(i)->Pt());
|
50 |
|
|
fIC5CaloJetsEta->Fill(fIC5Jets->At(i)->Eta());
|
51 |
|
|
}
|
52 |
|
|
fIC5NCaloJets->Fill(GoodIC5Jets->GetEntries());
|
53 |
|
|
|
54 |
|
|
ObjArray<GenJet> *GoodIC5GenJets = new ObjArray<GenJet>;
|
55 |
|
|
LoadBranch(fIC5GenJetName);
|
56 |
|
|
for (UInt_t i=0; i<fIC5GenJets->GetEntries(); ++i) {
|
57 |
|
|
if(fabs(fIC5GenJets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
|
58 |
|
|
if(fIC5GenJets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
|
59 |
|
|
GoodIC5GenJets->Add(fIC5GenJets->At(i));
|
60 |
|
|
}
|
61 |
|
|
|
62 |
|
|
std::map <GenJet*, Jet*> IC5GenJetToIC5JetMap;
|
63 |
|
|
for (UInt_t i=0; i<GoodIC5GenJets->GetEntries(); ++i) {
|
64 |
|
|
double minDR = 5000;
|
65 |
|
|
int closestMatchIndex = -1;
|
66 |
|
|
for (UInt_t j=0; j<GoodIC5Jets->GetEntries(); j++) {
|
67 |
|
|
double tempDR = MathUtils::DeltaR(GoodIC5GenJets->At(i)->Mom(), GoodIC5Jets->At(j)->Mom());
|
68 |
|
|
if (tempDR < minDR) {
|
69 |
|
|
closestMatchIndex = j;
|
70 |
|
|
minDR = tempDR;
|
71 |
|
|
}
|
72 |
|
|
}
|
73 |
|
|
if (closestMatchIndex > -1)
|
74 |
|
|
if (fPrintDebug) cerr << "GoodIC5GenJets->At(i) " << GoodIC5GenJets->At(i)->Pt() << " " << GoodIC5GenJets->At(i)->Eta() << " " << GoodIC5GenJets->At(i)->Phi() << " GoodIC5Jets->At(" << closestMatchIndex << ") " << GoodIC5Jets->At(closestMatchIndex)->Pt() << " " << GoodIC5Jets->At(closestMatchIndex)->Eta() << " " << GoodIC5Jets->At(closestMatchIndex)->Phi() << " DR = " << minDR << endl;
|
75 |
|
|
if (minDR < 0.5)
|
76 |
|
|
{
|
77 |
|
|
IC5GenJetToIC5JetMap[GoodIC5GenJets->At(i)] = GoodIC5Jets->At(closestMatchIndex);
|
78 |
|
|
}
|
79 |
|
|
}
|
80 |
|
|
|
81 |
|
|
std::map <Jet*, GenJet*> IC5JetToIC5GenJetMap;
|
82 |
|
|
for (UInt_t i=0; i<GoodIC5Jets->GetEntries(); ++i) {
|
83 |
|
|
double minDR = 5000;
|
84 |
|
|
int closestMatchIndex = -1;
|
85 |
|
|
for (UInt_t j=0; j<GoodIC5GenJets->GetEntries(); j++) {
|
86 |
|
|
double tempDR = MathUtils::DeltaR(GoodIC5Jets->At(i)->Mom(), GoodIC5GenJets->At(j)->Mom());
|
87 |
|
|
if (tempDR < minDR) {
|
88 |
|
|
closestMatchIndex = j;
|
89 |
|
|
minDR = tempDR;
|
90 |
|
|
}
|
91 |
|
|
}
|
92 |
|
|
if (minDR < 0.5)
|
93 |
|
|
IC5JetToIC5GenJetMap[GoodIC5Jets->At(i)] = GoodIC5GenJets->At(closestMatchIndex);
|
94 |
|
|
}
|
95 |
|
|
|
96 |
|
|
for (UInt_t i=0; i<GoodIC5GenJets->GetEntries(); i++) {
|
97 |
|
|
|
98 |
|
|
fIC5NGenJetsVsGenJetPt->Fill(GoodIC5GenJets->At(i)->Pt());
|
99 |
|
|
fIC5NGenJetsVsGenJetEta->Fill(GoodIC5GenJets->At(i)->Eta());
|
100 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 20.0 && GoodIC5GenJets->At(i)->Pt() < 30.0)
|
101 |
|
|
fIC5NGenJetsVsGenJetEta_Pt20To30->Fill(GoodIC5GenJets->At(i)->Eta());
|
102 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 30.0 && GoodIC5GenJets->At(i)->Pt() < 40.0)
|
103 |
|
|
fIC5NGenJetsVsGenJetEta_Pt30To40->Fill(GoodIC5GenJets->At(i)->Eta());
|
104 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 60.0 && GoodIC5GenJets->At(i)->Pt() < 80.0)
|
105 |
|
|
fIC5NGenJetsVsGenJetEta_Pt60To80->Fill(GoodIC5GenJets->At(i)->Eta());
|
106 |
|
|
|
107 |
|
|
map<GenJet*, Jet*>::iterator iter =
|
108 |
|
|
IC5GenJetToIC5JetMap.find(GoodIC5GenJets->At(i));
|
109 |
|
|
if (iter != IC5GenJetToIC5JetMap.end()) {
|
110 |
|
|
double dR = MathUtils::DeltaR(GoodIC5GenJets->At(i)->Mom(), iter->second->Mom());
|
111 |
|
|
|
112 |
|
|
fIC5GenJetRecoJetDeltaR->Fill(dR);
|
113 |
|
|
fIC5GenJetRecoJetDeltaEta->Fill(fabs(GoodIC5GenJets->At(i)->Eta() - iter->second->Eta()));
|
114 |
|
|
fIC5GenJetRecoJetDeltaPhi->Fill(MathUtils::DeltaPhi(GoodIC5GenJets->At(i)->Phi(),iter->second->Phi()));
|
115 |
|
|
|
116 |
|
|
if (fabs(iter->second->Eta()) <= 1.4)
|
117 |
|
|
fIC5JetResponseVsGenJetPtInBarrel->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
|
118 |
|
|
else if (fabs(iter->second->Eta()) <= 3.0)
|
119 |
|
|
fIC5JetResponseVsGenJetPtInEndcap->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
|
120 |
|
|
else if (fabs(iter->second->Eta()) <= 5.0)
|
121 |
|
|
fIC5JetResponseVsGenJetPtForward->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
|
122 |
|
|
|
123 |
|
|
//avoid the low pt region
|
124 |
|
|
if (iter->second->Pt() > 30)
|
125 |
|
|
fIC5JetResponseVsCaloJetEta->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),iter->second->Eta());
|
126 |
|
|
|
127 |
|
|
if (fabs(iter->second->Eta()) <= 3.0) {
|
128 |
|
|
fIC5CentralGenJetRecoJetDeltaR->Fill(dR);
|
129 |
|
|
} else if (fabs(iter->second->Eta()) <= 5.0) {
|
130 |
|
|
fIC5ForwardGenJetRecoJetDeltaR->Fill(dR);
|
131 |
|
|
}
|
132 |
|
|
|
133 |
|
|
//Number of Matched CaloJet vs Genjet pt and eta
|
134 |
|
|
fIC5NMatchedCaloJetsVsGenJetPt->Fill(GoodIC5GenJets->At(i)->Pt());
|
135 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta->Fill(GoodIC5GenJets->At(i)->Eta());
|
136 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 20.0 && GoodIC5GenJets->At(i)->Pt() < 30.0)
|
137 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30->Fill(GoodIC5GenJets->At(i)->Eta());
|
138 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 30.0 && GoodIC5GenJets->At(i)->Pt() < 40.0)
|
139 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40->Fill(GoodIC5GenJets->At(i)->Eta());
|
140 |
|
|
if (GoodIC5GenJets->At(i)->Pt() > 60.0 && GoodIC5GenJets->At(i)->Pt() < 80.0)
|
141 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80->Fill(GoodIC5GenJets->At(i)->Eta());
|
142 |
|
|
|
143 |
|
|
fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt->Fill(iter->second->Pt() * iter->second->L2RelativeCorrectionScale() * iter->second->L3AbsoluteCorrectionScale() / GoodIC5GenJets->At(i)->Pt(), GoodIC5GenJets->At(i)->Pt());
|
144 |
|
|
}
|
145 |
|
|
}
|
146 |
|
|
|
147 |
|
|
for (UInt_t i=0; i<GoodIC5Jets->GetEntries(); i++) {
|
148 |
|
|
map<Jet*, GenJet*>::iterator iter =
|
149 |
|
|
IC5JetToIC5GenJetMap.find(GoodIC5Jets->At(i));
|
150 |
|
|
//unmatched ones
|
151 |
|
|
if (iter == IC5JetToIC5GenJetMap.end()) {
|
152 |
|
|
fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt->Fill(GoodIC5Jets->At(i)->Pt() * GoodIC5Jets->At(i)->L2RelativeCorrectionScale() * GoodIC5Jets->At(i)->L3AbsoluteCorrectionScale());
|
153 |
|
|
fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta->Fill(GoodIC5Jets->At(i)->Eta());
|
154 |
|
|
}
|
155 |
|
|
}
|
156 |
|
|
|
157 |
|
|
|
158 |
|
|
//************************************************************************************************
|
159 |
|
|
//
|
160 |
|
|
// For SisCone 0.5 Jets
|
161 |
|
|
//
|
162 |
|
|
//************************************************************************************************
|
163 |
|
|
|
164 |
|
|
ObjArray<Jet> *GoodSC5Jets = new ObjArray<Jet>;
|
165 |
|
|
LoadBranch(fSC5JetName);
|
166 |
|
|
for (UInt_t i=0; i<fSC5Jets->GetEntries(); ++i) {
|
167 |
|
|
if(fabs(fSC5Jets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
|
168 |
|
|
if(fSC5Jets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
|
169 |
|
|
GoodSC5Jets->Add(fSC5Jets->At(i));
|
170 |
|
|
fSC5CaloJetsPt->Fill(fSC5Jets->At(i)->Pt());
|
171 |
|
|
fSC5CaloJetsEta->Fill(fSC5Jets->At(i)->Eta());
|
172 |
|
|
}
|
173 |
|
|
fSC5NCaloJets->Fill(GoodSC5Jets->GetEntries());
|
174 |
|
|
|
175 |
|
|
ObjArray<GenJet> *GoodSC5GenJets = new ObjArray<GenJet>;
|
176 |
|
|
LoadBranch(fSC5GenJetName);
|
177 |
|
|
for (UInt_t i=0; i<fSC5GenJets->GetEntries(); ++i) {
|
178 |
|
|
if(fabs(fSC5GenJets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
|
179 |
|
|
if(fSC5GenJets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
|
180 |
|
|
GoodSC5GenJets->Add(fSC5GenJets->At(i));
|
181 |
|
|
}
|
182 |
|
|
|
183 |
|
|
std::map <GenJet*, Jet*> SC5GenJetToSC5JetMap;
|
184 |
|
|
for (UInt_t i=0; i<GoodSC5GenJets->GetEntries(); ++i) {
|
185 |
|
|
double minDR = 5000;
|
186 |
|
|
int closestMatchIndex = -1;
|
187 |
|
|
for (UInt_t j=0; j<GoodSC5Jets->GetEntries(); j++) {
|
188 |
|
|
double tempDR = MathUtils::DeltaR(GoodSC5GenJets->At(i)->Mom(), GoodSC5Jets->At(j)->Mom());
|
189 |
|
|
if (tempDR < minDR) {
|
190 |
|
|
closestMatchIndex = j;
|
191 |
|
|
minDR = tempDR;
|
192 |
|
|
}
|
193 |
|
|
}
|
194 |
|
|
if (closestMatchIndex > -1)
|
195 |
|
|
if (fPrintDebug) cerr << "GoodSC5GenJets->At(i) " << GoodSC5GenJets->At(i)->Pt() << " " << GoodSC5GenJets->At(i)->Eta() << " " << GoodSC5GenJets->At(i)->Phi() << " GoodSC5Jets->At(" << closestMatchIndex << ") " << GoodSC5Jets->At(closestMatchIndex)->Pt() << " " << GoodSC5Jets->At(closestMatchIndex)->Eta() << " " << GoodSC5Jets->At(closestMatchIndex)->Phi() << " DR = " << minDR << endl;
|
196 |
|
|
if (minDR < 0.5)
|
197 |
|
|
{
|
198 |
|
|
SC5GenJetToSC5JetMap[GoodSC5GenJets->At(i)] = GoodSC5Jets->At(closestMatchIndex);
|
199 |
|
|
}
|
200 |
|
|
}
|
201 |
|
|
std::map <Jet*, GenJet*> SC5JetToSC5GenJetMap;
|
202 |
|
|
for (UInt_t i=0; i<GoodSC5Jets->GetEntries(); ++i) {
|
203 |
|
|
double minDR = 5000;
|
204 |
|
|
int closestMatchIndex = -1;
|
205 |
|
|
for (UInt_t j=0; j<GoodSC5GenJets->GetEntries(); j++) {
|
206 |
|
|
double tempDR = MathUtils::DeltaR(GoodSC5Jets->At(i)->Mom(), GoodSC5GenJets->At(j)->Mom());
|
207 |
|
|
if (tempDR < minDR) {
|
208 |
|
|
closestMatchIndex = j;
|
209 |
|
|
minDR = tempDR;
|
210 |
|
|
}
|
211 |
|
|
}
|
212 |
|
|
if (minDR < 0.5)
|
213 |
|
|
SC5JetToSC5GenJetMap[GoodSC5Jets->At(i)] = GoodSC5GenJets->At(closestMatchIndex);
|
214 |
|
|
}
|
215 |
|
|
for (UInt_t i=0; i<GoodSC5GenJets->GetEntries(); i++) {
|
216 |
|
|
|
217 |
|
|
fSC5NGenJetsVsGenJetPt->Fill(GoodSC5GenJets->At(i)->Pt());
|
218 |
|
|
fSC5NGenJetsVsGenJetEta->Fill(GoodSC5GenJets->At(i)->Eta());
|
219 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 20.0 && GoodSC5GenJets->At(i)->Pt() < 30.0)
|
220 |
|
|
fSC5NGenJetsVsGenJetEta_Pt20To30->Fill(GoodSC5GenJets->At(i)->Eta());
|
221 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 30.0 && GoodSC5GenJets->At(i)->Pt() < 40.0)
|
222 |
|
|
fSC5NGenJetsVsGenJetEta_Pt30To40->Fill(GoodSC5GenJets->At(i)->Eta());
|
223 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 60.0 && GoodSC5GenJets->At(i)->Pt() < 80.0)
|
224 |
|
|
fSC5NGenJetsVsGenJetEta_Pt60To80->Fill(GoodSC5GenJets->At(i)->Eta());
|
225 |
|
|
|
226 |
|
|
map<GenJet*, Jet*>::iterator iter =
|
227 |
|
|
SC5GenJetToSC5JetMap.find(GoodSC5GenJets->At(i));
|
228 |
|
|
if (iter != SC5GenJetToSC5JetMap.end()) {
|
229 |
|
|
double dR = MathUtils::DeltaR(GoodSC5GenJets->At(i)->Mom(), iter->second->Mom());
|
230 |
|
|
|
231 |
|
|
fSC5GenJetRecoJetDeltaR->Fill(dR);
|
232 |
|
|
fSC5GenJetRecoJetDeltaEta->Fill(fabs(GoodSC5GenJets->At(i)->Eta() - iter->second->Eta()));
|
233 |
|
|
fSC5GenJetRecoJetDeltaPhi->Fill(MathUtils::DeltaPhi(GoodSC5GenJets->At(i)->Phi(),iter->second->Phi()));
|
234 |
|
|
|
235 |
|
|
if (fabs(iter->second->Eta()) <= 1.4)
|
236 |
|
|
fSC5JetResponseVsGenJetPtInBarrel->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),GoodSC5GenJets->At(i)->Pt());
|
237 |
|
|
else if (fabs(iter->second->Eta()) <= 3.0)
|
238 |
|
|
fSC5JetResponseVsGenJetPtInEndcap->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),GoodSC5GenJets->At(i)->Pt());
|
239 |
|
|
else if (fabs(iter->second->Eta()) <= 5.0)
|
240 |
|
|
fSC5JetResponseVsGenJetPtForward->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),GoodSC5GenJets->At(i)->Pt());
|
241 |
|
|
//avoid the low pt region
|
242 |
|
|
if (iter->second->Pt() > 30)
|
243 |
|
|
fSC5JetResponseVsCaloJetEta->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),iter->second->Eta());
|
244 |
|
|
|
245 |
|
|
if (fabs(iter->second->Eta()) <= 3.0) {
|
246 |
|
|
fSC5CentralGenJetRecoJetDeltaR->Fill(dR);
|
247 |
|
|
} else if (fabs(iter->second->Eta()) <= 5.0) {
|
248 |
|
|
fSC5ForwardGenJetRecoJetDeltaR->Fill(dR);
|
249 |
|
|
}
|
250 |
|
|
|
251 |
|
|
//Number of Matched CaloJet vs Genjet pt and eta
|
252 |
|
|
fSC5NMatchedCaloJetsVsGenJetPt->Fill(GoodSC5GenJets->At(i)->Pt());
|
253 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta->Fill(GoodSC5GenJets->At(i)->Eta());
|
254 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 20.0 && GoodSC5GenJets->At(i)->Pt() < 30.0)
|
255 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30->Fill(GoodSC5GenJets->At(i)->Eta());
|
256 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 30.0 && GoodSC5GenJets->At(i)->Pt() < 40.0)
|
257 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40->Fill(GoodSC5GenJets->At(i)->Eta());
|
258 |
|
|
if (GoodSC5GenJets->At(i)->Pt() > 60.0 && GoodSC5GenJets->At(i)->Pt() < 80.0)
|
259 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80->Fill(GoodSC5GenJets->At(i)->Eta());
|
260 |
|
|
fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt->Fill(iter->second->Pt() * iter->second->L2RelativeCorrectionScale() * iter->second->L3AbsoluteCorrectionScale() / GoodSC5GenJets->At(i)->Pt(), GoodSC5GenJets->At(i)->Pt());
|
261 |
|
|
}
|
262 |
|
|
}
|
263 |
|
|
|
264 |
|
|
for (UInt_t i=0; i<GoodSC5Jets->GetEntries(); i++) {
|
265 |
|
|
map<Jet*, GenJet*>::iterator iter =
|
266 |
|
|
SC5JetToSC5GenJetMap.find(GoodSC5Jets->At(i));
|
267 |
|
|
//unmatched ones
|
268 |
|
|
if (iter == SC5JetToSC5GenJetMap.end()) {
|
269 |
|
|
fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt->Fill(GoodSC5Jets->At(i)->Pt() * GoodSC5Jets->At(i)->L2RelativeCorrectionScale() * GoodSC5Jets->At(i)->L3AbsoluteCorrectionScale());
|
270 |
|
|
fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta->Fill(GoodSC5Jets->At(i)->Eta());
|
271 |
|
|
}
|
272 |
|
|
}
|
273 |
|
|
|
274 |
|
|
}
|
275 |
|
|
|
276 |
|
|
//--------------------------------------------------------------------------------------------------
|
277 |
|
|
void JetValidationMod::SlaveBegin()
|
278 |
|
|
{
|
279 |
|
|
// Run startup code on the computer (slave) doing the actual analysis. Here,
|
280 |
|
|
// we typically initialize histograms and other analysis objects and request
|
281 |
|
|
// branches. For this module, we request a branch of the MitTree.
|
282 |
|
|
|
283 |
|
|
ReqBranch(fIC5GenJetName, fIC5GenJets);
|
284 |
|
|
ReqBranch(fSC5GenJetName, fSC5GenJets);
|
285 |
|
|
ReqBranch(fIC5JetName, fIC5Jets);
|
286 |
|
|
ReqBranch(fSC5JetName, fSC5Jets);
|
287 |
|
|
|
288 |
loizides |
1.2 |
fIC5GenJetRecoJetDeltaR =
|
289 |
|
|
new TH1D("hIC5GenJetRecoJetDeltaR",";IC5GenJetRecoJetDeltaR;#",100,0.,1.0);
|
290 |
|
|
fIC5GenJetRecoJetDeltaEta =
|
291 |
|
|
new TH1D("hIC5GenJetRecoJetDeltaEta",";IC5GenJetRecoJetDeltaEta;#",100,0.,0.5);
|
292 |
|
|
fIC5GenJetRecoJetDeltaPhi =
|
293 |
|
|
new TH1D("hIC5GenJetRecoJetDeltaPhi",";IC5GenJetRecoJetDeltaPhi;#",100,0.,0.5);
|
294 |
loizides |
1.1 |
|
295 |
|
|
fIC5JetResponseVsGenJetPtInBarrel = new TH2D("hIC5JetResponseVsGenJetPtInBarrel",";IC5JetResponseVsGenJetPtInBarrel;#",500,0.,5.0,5000,0,5000);
|
296 |
|
|
fIC5JetResponseVsGenJetPtInEndcap = new TH2D("hIC5JetResponseVsGenJetPtInEndcap",";IC5JetResponseVsGenJetPtInEndcap;#",500,0.,5.0,5000,0,5000);
|
297 |
|
|
fIC5JetResponseVsGenJetPtForward = new TH2D("hIC5JetResponseVsGenJetPtForward",";IC5JetResponseVsGenJetPtForward;#",500,0.,5.0,5000,0,5000);
|
298 |
|
|
fIC5JetResponseVsCaloJetEta = new TH2D("hIC5JetResponseVsCaloJetEta",";IC5JetResponseVsCaloJetEta;#",500,0.,5.0,100,-5,5);
|
299 |
|
|
|
300 |
loizides |
1.2 |
fIC5CentralGenJetRecoJetDeltaR =
|
301 |
|
|
new TH1D("hIC5CentralGenJetRecoJetDeltaR",";IC5CentralGenJetRecoJetDeltaR;#",100,0.,1.0);
|
302 |
loizides |
1.1 |
fIC5ForwardGenJetRecoJetDeltaR = new TH1D("hIC5ForwardGenJetRecoJetDeltaR",";IC5ForwardGenJetRecoJetDeltaR;#",100,0.,1.0);
|
303 |
|
|
|
304 |
|
|
fIC5NMatchedCaloJetsVsGenJetPt = new TH1D("hIC5NMatchedCaloJetsVsGenJetPt",";IC5GenJetPt;#",5000,0.,5000);
|
305 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta",";IC5GenJetEta;#",100,-5.0,5.0);
|
306 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt20To30",";IC5GenJetEta;#",100,-5.0,5.0);
|
307 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt30To40",";IC5GenJetEta;#",100,-5.0,5.0);
|
308 |
|
|
fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt60To80",";IC5GenJetEta;#",100,-5.0,5.0);
|
309 |
|
|
fIC5NGenJetsVsGenJetPt = new TH1D("hIC5NGenJetsVsGenJetPt",";IC5GenJetPt;#",5000,0.,5000);
|
310 |
|
|
fIC5NGenJetsVsGenJetEta = new TH1D("hIC5NGenJetsVsGenJetEta",";IC5GenJetEta;#",100,-5.,5.0);
|
311 |
|
|
fIC5NGenJetsVsGenJetEta_Pt20To30 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt20To30",";IC5GenJetEta;#",100,-5.,5.0);
|
312 |
|
|
fIC5NGenJetsVsGenJetEta_Pt30To40 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt30To40",";IC5GenJetEta;#",100,-5.,5.0);
|
313 |
|
|
fIC5NGenJetsVsGenJetEta_Pt60To80 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt60To80",";IC5GenJetEta;#",100,-5.,5.0);
|
314 |
|
|
|
315 |
|
|
fIC5CaloJetsPt = new TH1D("hIC5CaloJetsPt",";IC5CaloJetPt;#",5000,0.,5000);
|
316 |
|
|
fIC5CaloJetsEta = new TH1D("hIC5CaloJetsEta",";IC5CaloJetEta;#",100,-5.0,5.0);
|
317 |
|
|
fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt = new TH1D("hIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt",";IC5CorrectedCaloJetPt;#",5000,0.0,5000.0);
|
318 |
|
|
fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta = new TH1D("hIC5NUnmatchedCalojetsVsCorrectedCaloJetEta",";IC5CorrectedCaloJetEta;#",100,-5.0,5.0);
|
319 |
|
|
|
320 |
|
|
fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt = new TH2D("hIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt",";IC5GenJetPt;#",1000,0,1000.0, 100,0.0,5.0);
|
321 |
|
|
|
322 |
|
|
fIC5NCaloJets = new TH1D("hIC5NCaloJets",";IC5NCaloJets;#",100,0.,15);
|
323 |
|
|
|
324 |
|
|
|
325 |
|
|
fSC5GenJetRecoJetDeltaR = new TH1D("hSC5GenJetRecoJetDeltaR",";SC5GenJetRecoJetDeltaR;#",100,0.,1.0);
|
326 |
|
|
fSC5GenJetRecoJetDeltaEta = new TH1D("hSC5GenJetRecoJetDeltaEta",";SC5GenJetRecoJetDeltaEta;#",100,0.,0.5);
|
327 |
|
|
fSC5GenJetRecoJetDeltaPhi = new TH1D("hSC5GenJetRecoJetDeltaPhi",";SC5GenJetRecoJetDeltaPhi;#",100,0.,0.5);
|
328 |
|
|
|
329 |
|
|
fSC5JetResponseVsGenJetPtInBarrel = new TH2D("hSC5JetResponseVsGenJetPtInBarrel",";SC5JetResponseVsGenJetPtInBarrel;#",500,0.,5.0,5000,0,5000);
|
330 |
|
|
fSC5JetResponseVsGenJetPtInEndcap = new TH2D("hSC5JetResponseVsGenJetPtInEndcap",";SC5JetResponseVsGenJetPtInEndcap;#",500,0.,5.0,5000,0,5000);
|
331 |
|
|
fSC5JetResponseVsGenJetPtForward = new TH2D("hSC5JetResponseVsGenJetPtForward",";SC5JetResponseVsGenJetPtForward;#",500,0.,5.0,5000,0,5000);
|
332 |
|
|
|
333 |
|
|
fSC5JetResponseVsCaloJetEta = new TH2D("hSC5JetResponseVsCaloJetEta",";SC5JetResponseVsCaloJetEta;#",500,0.,5.0,100,-5,5);
|
334 |
|
|
|
335 |
|
|
|
336 |
|
|
|
337 |
|
|
fSC5CentralGenJetRecoJetDeltaR = new TH1D("hSC5CentralGenJetRecoJetDeltaR",";SC5CentralGenJetRecoJetDeltaR;#",100,0.,1.0);
|
338 |
|
|
fSC5ForwardGenJetRecoJetDeltaR = new TH1D("hSC5ForwardGenJetRecoJetDeltaR",";SC5ForwardGenJetRecoJetDeltaR;#",100,0.,1.0);
|
339 |
|
|
|
340 |
|
|
fSC5NMatchedCaloJetsVsGenJetPt = new TH1D("hSC5NMatchedCaloJetsVsGenJetPt",";SC5GenJetPt;#",5000,0.,5000);
|
341 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta = new TH1D("hSC5NMatchedCaloJetsVsGenJetEta",";SC5GenJetEta;#",100,-5.0,5.0);
|
342 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30 = new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt20To30",";SC5GenJetEta;#",100,-5.0,5.0);
|
343 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40 = new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt30To40",";SC5GenJetEta;#",100,-5.0,5.0);
|
344 |
|
|
fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80 = new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt60To80",";SC5GenJetEta;#",100,-5.0,5.0);
|
345 |
|
|
fSC5NGenJetsVsGenJetPt = new TH1D("hSC5NGenJetsVsGenJetPt",";SC5GenJetPt;#",5000,0.,5000);
|
346 |
|
|
fSC5NGenJetsVsGenJetEta = new TH1D("hSC5NGenJetsVsGenJetEta",";SC5GenJetEta;#",100,-5.,5.0);
|
347 |
|
|
fSC5NGenJetsVsGenJetEta_Pt20To30 = new TH1D("hSC5NGenJetsVsGenJetEta_Pt20To30",";SC5GenJetEta;#",100,-5.,5.0);
|
348 |
|
|
fSC5NGenJetsVsGenJetEta_Pt30To40 = new TH1D("hSC5NGenJetsVsGenJetEta_Pt30To40",";SC5GenJetEta;#",100,-5.,5.0);
|
349 |
|
|
fSC5NGenJetsVsGenJetEta_Pt60To80 = new TH1D("hSC5NGenJetsVsGenJetEta_Pt60To80",";SC5GenJetEta;#",100,-5.,5.0);
|
350 |
|
|
|
351 |
|
|
fSC5CaloJetsPt = new TH1D("hSC5CaloJetsPt",";SC5CaloJetPt;#",5000,0.,5000);
|
352 |
|
|
fSC5CaloJetsEta = new TH1D("hSC5CaloJetsEta",";SC5CaloJetEta;#",100,-5.0,5.0);
|
353 |
|
|
fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt = new TH1D("hSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt",";SC5CorrectedCaloJetPt;#",5000,0.0,5000.0);
|
354 |
|
|
fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta = new TH1D("hSC5NUnmatchedCalojetsVsCorrectedCaloJetEta",";SC5CorrectedCaloJetEta;#",100,-5.0,5.0);
|
355 |
|
|
|
356 |
|
|
fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt = new TH2D("hSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt",";SC5GenJetPt;#",1000,0,1000.0, 100,0.0,5.0);
|
357 |
|
|
|
358 |
|
|
fSC5NCaloJets = new TH1D("hSC5NCaloJets",";SC5NCaloJets;#",100,0.,15);
|
359 |
|
|
|
360 |
|
|
AddOutput(fIC5GenJetRecoJetDeltaR);
|
361 |
|
|
AddOutput(fIC5GenJetRecoJetDeltaEta);
|
362 |
|
|
AddOutput(fIC5GenJetRecoJetDeltaPhi);
|
363 |
|
|
AddOutput(fIC5JetResponseVsGenJetPtInBarrel);
|
364 |
|
|
AddOutput(fIC5JetResponseVsGenJetPtInEndcap);
|
365 |
|
|
AddOutput(fIC5JetResponseVsGenJetPtForward);
|
366 |
|
|
AddOutput(fIC5JetResponseVsCaloJetEta);
|
367 |
|
|
AddOutput(fIC5CentralGenJetRecoJetDeltaR);
|
368 |
|
|
AddOutput(fIC5ForwardGenJetRecoJetDeltaR);
|
369 |
|
|
AddOutput(fIC5NMatchedCaloJetsVsGenJetPt);
|
370 |
|
|
AddOutput(fIC5NMatchedCaloJetsVsGenJetEta);
|
371 |
|
|
AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30);
|
372 |
|
|
AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40);
|
373 |
|
|
AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80);
|
374 |
|
|
AddOutput(fIC5NGenJetsVsGenJetPt);
|
375 |
|
|
AddOutput(fIC5NGenJetsVsGenJetEta);
|
376 |
|
|
AddOutput(fIC5NGenJetsVsGenJetEta_Pt20To30);
|
377 |
|
|
AddOutput(fIC5NGenJetsVsGenJetEta_Pt30To40);
|
378 |
|
|
AddOutput(fIC5NGenJetsVsGenJetEta_Pt60To80);
|
379 |
|
|
AddOutput(fIC5CaloJetsPt);
|
380 |
|
|
AddOutput(fIC5CaloJetsEta);
|
381 |
|
|
AddOutput(fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt);
|
382 |
|
|
AddOutput(fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta);
|
383 |
|
|
AddOutput(fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt);
|
384 |
|
|
AddOutput(fIC5NCaloJets);
|
385 |
|
|
AddOutput(fSC5GenJetRecoJetDeltaR);
|
386 |
|
|
AddOutput(fSC5GenJetRecoJetDeltaEta);
|
387 |
|
|
AddOutput(fSC5GenJetRecoJetDeltaPhi);
|
388 |
|
|
AddOutput(fSC5JetResponseVsGenJetPtInBarrel);
|
389 |
|
|
AddOutput(fSC5JetResponseVsGenJetPtInEndcap);
|
390 |
|
|
AddOutput(fSC5JetResponseVsGenJetPtForward);
|
391 |
|
|
AddOutput(fSC5JetResponseVsCaloJetEta);
|
392 |
|
|
AddOutput(fSC5CentralGenJetRecoJetDeltaR);
|
393 |
|
|
AddOutput(fSC5ForwardGenJetRecoJetDeltaR);
|
394 |
|
|
AddOutput(fSC5NMatchedCaloJetsVsGenJetPt);
|
395 |
|
|
AddOutput(fSC5NMatchedCaloJetsVsGenJetEta);
|
396 |
|
|
AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30);
|
397 |
|
|
AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40);
|
398 |
|
|
AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80);
|
399 |
|
|
AddOutput(fSC5NGenJetsVsGenJetPt);
|
400 |
|
|
AddOutput(fSC5NGenJetsVsGenJetEta);
|
401 |
|
|
AddOutput(fSC5NGenJetsVsGenJetEta_Pt20To30);
|
402 |
|
|
AddOutput(fSC5NGenJetsVsGenJetEta_Pt30To40);
|
403 |
|
|
AddOutput(fSC5NGenJetsVsGenJetEta_Pt60To80);
|
404 |
|
|
AddOutput(fSC5CaloJetsPt);
|
405 |
|
|
AddOutput(fSC5CaloJetsEta);
|
406 |
|
|
AddOutput(fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt);
|
407 |
|
|
AddOutput(fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta);
|
408 |
|
|
AddOutput(fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt);
|
409 |
|
|
AddOutput(fSC5NCaloJets);
|
410 |
|
|
}
|