ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/Validation/src/JetValidationMod.cc
Revision: 1.5
Committed: Fri Jul 10 14:04:44 2009 UTC (15 years, 9 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_029c, Mit_029b, Mit_029a, Mit_028a, Mit_028, Mit_027, Mit_027a, Mit_025e, Mit_025d, Mit_025c, Mit_025b, Mit_025a, Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024, Mit_023, Mit_022a, Mit_022, Mit_020d, TMit_020d, Mit_020c, Mit_021, Mit_021pre2, Mit_021pre1, Mit_020b, Mit_020a, Mit_020, Mit_020pre1, Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, Mit_013pre1, Mit_012i, Mit_012g, Mit_012f, Mit_012e, Mit_012d, Mit_012c, Mit_012b, Mit_012a, Mit_012, Mit_011a, Mit_011, Mit_010a, Mit_010, HEAD
Changes since 1.4: +19 -8 lines
Log Message:
Cleanup

File Contents

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