ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/Validation/src/JetValidationMod.cc
Revision: 1.1
Committed: Wed Oct 15 06:05:02 2008 UTC (16 years, 7 months ago) by loizides
Content type: text/plain
Branch: MAIN
Log Message:
Added MitPhysics.

File Contents

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