ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/Validation/src/JetValidationMod.cc
Revision: 1.4
Committed: Mon Jun 15 15:00:23 2009 UTC (15 years, 10 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_009c, Mit_009b
Changes since 1.3: +7 -5 lines
Log Message:
Added proper fwd defs plus split up complilation of MitAna/DataTree LinkDefs.

File Contents

# User Rev Content
1 loizides 1.4 // $Id: JetValidationMod.cc,v 1.3 2009/04/07 15:37:10 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     fIC5GenJetRecoJetDeltaPhi->Fill(MathUtils::DeltaPhi(GoodIC5GenJets->At(i)->Phi(),iter->second->Phi()));
123    
124     if (fabs(iter->second->Eta()) <= 1.4)
125     fIC5JetResponseVsGenJetPtInBarrel->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
126     else if (fabs(iter->second->Eta()) <= 3.0)
127     fIC5JetResponseVsGenJetPtInEndcap->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
128     else if (fabs(iter->second->Eta()) <= 5.0)
129     fIC5JetResponseVsGenJetPtForward->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),GoodIC5GenJets->At(i)->Pt());
130    
131     //avoid the low pt region
132     if (iter->second->Pt() > 30)
133     fIC5JetResponseVsCaloJetEta->Fill(iter->second->Et()/GoodIC5GenJets->At(i)->Et(),iter->second->Eta());
134    
135     if (fabs(iter->second->Eta()) <= 3.0) {
136     fIC5CentralGenJetRecoJetDeltaR->Fill(dR);
137     } else if (fabs(iter->second->Eta()) <= 5.0) {
138     fIC5ForwardGenJetRecoJetDeltaR->Fill(dR);
139     }
140    
141     //Number of Matched CaloJet vs Genjet pt and eta
142     fIC5NMatchedCaloJetsVsGenJetPt->Fill(GoodIC5GenJets->At(i)->Pt());
143     fIC5NMatchedCaloJetsVsGenJetEta->Fill(GoodIC5GenJets->At(i)->Eta());
144     if (GoodIC5GenJets->At(i)->Pt() > 20.0 && GoodIC5GenJets->At(i)->Pt() < 30.0)
145     fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30->Fill(GoodIC5GenJets->At(i)->Eta());
146     if (GoodIC5GenJets->At(i)->Pt() > 30.0 && GoodIC5GenJets->At(i)->Pt() < 40.0)
147     fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40->Fill(GoodIC5GenJets->At(i)->Eta());
148     if (GoodIC5GenJets->At(i)->Pt() > 60.0 && GoodIC5GenJets->At(i)->Pt() < 80.0)
149     fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80->Fill(GoodIC5GenJets->At(i)->Eta());
150    
151     fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt->Fill(iter->second->Pt() * iter->second->L2RelativeCorrectionScale() * iter->second->L3AbsoluteCorrectionScale() / GoodIC5GenJets->At(i)->Pt(), GoodIC5GenJets->At(i)->Pt());
152     }
153     }
154    
155     for (UInt_t i=0; i<GoodIC5Jets->GetEntries(); i++) {
156     map<Jet*, GenJet*>::iterator iter =
157     IC5JetToIC5GenJetMap.find(GoodIC5Jets->At(i));
158     //unmatched ones
159     if (iter == IC5JetToIC5GenJetMap.end()) {
160     fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt->Fill(GoodIC5Jets->At(i)->Pt() * GoodIC5Jets->At(i)->L2RelativeCorrectionScale() * GoodIC5Jets->At(i)->L3AbsoluteCorrectionScale());
161     fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta->Fill(GoodIC5Jets->At(i)->Eta());
162     }
163     }
164    
165    
166     //************************************************************************************************
167     //
168     // For SisCone 0.5 Jets
169     //
170     //************************************************************************************************
171    
172     ObjArray<Jet> *GoodSC5Jets = new ObjArray<Jet>;
173     LoadBranch(fSC5JetName);
174     for (UInt_t i=0; i<fSC5Jets->GetEntries(); ++i) {
175     if(fabs(fSC5Jets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
176     if(fSC5Jets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
177     GoodSC5Jets->Add(fSC5Jets->At(i));
178     fSC5CaloJetsPt->Fill(fSC5Jets->At(i)->Pt());
179     fSC5CaloJetsEta->Fill(fSC5Jets->At(i)->Eta());
180     }
181     fSC5NCaloJets->Fill(GoodSC5Jets->GetEntries());
182    
183     ObjArray<GenJet> *GoodSC5GenJets = new ObjArray<GenJet>;
184     LoadBranch(fSC5GenJetName);
185     for (UInt_t i=0; i<fSC5GenJets->GetEntries(); ++i) {
186     if(fabs(fSC5GenJets->At(i)->Eta()) >= 5.0) continue; // minimum eta requirement
187     if(fSC5GenJets->At(i)->Pt() <= 5) continue; // minimum Pt requirement
188     GoodSC5GenJets->Add(fSC5GenJets->At(i));
189     }
190    
191     std::map <GenJet*, Jet*> SC5GenJetToSC5JetMap;
192     for (UInt_t i=0; i<GoodSC5GenJets->GetEntries(); ++i) {
193 loizides 1.3 Double_t minDR = 5000;
194 loizides 1.1 int closestMatchIndex = -1;
195     for (UInt_t j=0; j<GoodSC5Jets->GetEntries(); j++) {
196 loizides 1.3 Double_t tempDR = MathUtils::DeltaR(GoodSC5GenJets->At(i)->Mom(), GoodSC5Jets->At(j)->Mom());
197 loizides 1.1 if (tempDR < minDR) {
198     closestMatchIndex = j;
199     minDR = tempDR;
200     }
201     }
202     if (closestMatchIndex > -1)
203 loizides 1.3 if (fPrintDebug) cerr << "GoodSC5GenJets->At(i) " << GoodSC5GenJets->At(i)->Pt() << " "
204     << GoodSC5GenJets->At(i)->Eta() << " " << GoodSC5GenJets->At(i)->Phi()
205     << " GoodSC5Jets->At(" << closestMatchIndex << ") "
206     << GoodSC5Jets->At(closestMatchIndex)->Pt() << " "
207     << GoodSC5Jets->At(closestMatchIndex)->Eta() << " "
208     << GoodSC5Jets->At(closestMatchIndex)->Phi() << " DR = "
209     << minDR << endl;
210 loizides 1.1 if (minDR < 0.5)
211     {
212     SC5GenJetToSC5JetMap[GoodSC5GenJets->At(i)] = GoodSC5Jets->At(closestMatchIndex);
213     }
214     }
215     std::map <Jet*, GenJet*> SC5JetToSC5GenJetMap;
216     for (UInt_t i=0; i<GoodSC5Jets->GetEntries(); ++i) {
217 loizides 1.3 Double_t minDR = 5000;
218 loizides 1.1 int closestMatchIndex = -1;
219     for (UInt_t j=0; j<GoodSC5GenJets->GetEntries(); j++) {
220 loizides 1.3 Double_t tempDR = MathUtils::DeltaR(GoodSC5Jets->At(i)->Mom(), GoodSC5GenJets->At(j)->Mom());
221 loizides 1.1 if (tempDR < minDR) {
222     closestMatchIndex = j;
223     minDR = tempDR;
224     }
225     }
226     if (minDR < 0.5)
227     SC5JetToSC5GenJetMap[GoodSC5Jets->At(i)] = GoodSC5GenJets->At(closestMatchIndex);
228     }
229     for (UInt_t i=0; i<GoodSC5GenJets->GetEntries(); i++) {
230    
231     fSC5NGenJetsVsGenJetPt->Fill(GoodSC5GenJets->At(i)->Pt());
232     fSC5NGenJetsVsGenJetEta->Fill(GoodSC5GenJets->At(i)->Eta());
233     if (GoodSC5GenJets->At(i)->Pt() > 20.0 && GoodSC5GenJets->At(i)->Pt() < 30.0)
234     fSC5NGenJetsVsGenJetEta_Pt20To30->Fill(GoodSC5GenJets->At(i)->Eta());
235     if (GoodSC5GenJets->At(i)->Pt() > 30.0 && GoodSC5GenJets->At(i)->Pt() < 40.0)
236     fSC5NGenJetsVsGenJetEta_Pt30To40->Fill(GoodSC5GenJets->At(i)->Eta());
237     if (GoodSC5GenJets->At(i)->Pt() > 60.0 && GoodSC5GenJets->At(i)->Pt() < 80.0)
238     fSC5NGenJetsVsGenJetEta_Pt60To80->Fill(GoodSC5GenJets->At(i)->Eta());
239    
240     map<GenJet*, Jet*>::iterator iter =
241     SC5GenJetToSC5JetMap.find(GoodSC5GenJets->At(i));
242     if (iter != SC5GenJetToSC5JetMap.end()) {
243 loizides 1.3 Double_t dR = MathUtils::DeltaR(GoodSC5GenJets->At(i)->Mom(), iter->second->Mom());
244 loizides 1.1
245     fSC5GenJetRecoJetDeltaR->Fill(dR);
246     fSC5GenJetRecoJetDeltaEta->Fill(fabs(GoodSC5GenJets->At(i)->Eta() - iter->second->Eta()));
247 loizides 1.3 fSC5GenJetRecoJetDeltaPhi->Fill(MathUtils::DeltaPhi(GoodSC5GenJets->At(i)->Phi(),
248     iter->second->Phi()));
249 loizides 1.1
250     if (fabs(iter->second->Eta()) <= 1.4)
251 loizides 1.3 fSC5JetResponseVsGenJetPtInBarrel->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),
252     GoodSC5GenJets->At(i)->Pt());
253 loizides 1.1 else if (fabs(iter->second->Eta()) <= 3.0)
254 loizides 1.3 fSC5JetResponseVsGenJetPtInEndcap->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),
255     GoodSC5GenJets->At(i)->Pt());
256 loizides 1.1 else if (fabs(iter->second->Eta()) <= 5.0)
257 loizides 1.3 fSC5JetResponseVsGenJetPtForward->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),
258     GoodSC5GenJets->At(i)->Pt());
259 loizides 1.1 //avoid the low pt region
260     if (iter->second->Pt() > 30)
261 loizides 1.3 fSC5JetResponseVsCaloJetEta->Fill(iter->second->Et()/GoodSC5GenJets->At(i)->Et(),
262     iter->second->Eta());
263 loizides 1.1
264     if (fabs(iter->second->Eta()) <= 3.0) {
265     fSC5CentralGenJetRecoJetDeltaR->Fill(dR);
266     } else if (fabs(iter->second->Eta()) <= 5.0) {
267     fSC5ForwardGenJetRecoJetDeltaR->Fill(dR);
268     }
269    
270     //Number of Matched CaloJet vs Genjet pt and eta
271     fSC5NMatchedCaloJetsVsGenJetPt->Fill(GoodSC5GenJets->At(i)->Pt());
272     fSC5NMatchedCaloJetsVsGenJetEta->Fill(GoodSC5GenJets->At(i)->Eta());
273     if (GoodSC5GenJets->At(i)->Pt() > 20.0 && GoodSC5GenJets->At(i)->Pt() < 30.0)
274     fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30->Fill(GoodSC5GenJets->At(i)->Eta());
275     if (GoodSC5GenJets->At(i)->Pt() > 30.0 && GoodSC5GenJets->At(i)->Pt() < 40.0)
276     fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40->Fill(GoodSC5GenJets->At(i)->Eta());
277     if (GoodSC5GenJets->At(i)->Pt() > 60.0 && GoodSC5GenJets->At(i)->Pt() < 80.0)
278     fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80->Fill(GoodSC5GenJets->At(i)->Eta());
279 loizides 1.3 fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt->Fill(
280     iter->second->Pt() * iter->second->L2RelativeCorrectionScale() *
281     iter->second->L3AbsoluteCorrectionScale() / GoodSC5GenJets->At(i)->Pt(),
282     GoodSC5GenJets->At(i)->Pt());
283 loizides 1.1 }
284     }
285    
286     for (UInt_t i=0; i<GoodSC5Jets->GetEntries(); i++) {
287     map<Jet*, GenJet*>::iterator iter =
288     SC5JetToSC5GenJetMap.find(GoodSC5Jets->At(i));
289     //unmatched ones
290     if (iter == SC5JetToSC5GenJetMap.end()) {
291 loizides 1.3 fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt->Fill(
292     GoodSC5Jets->At(i)->Pt() * GoodSC5Jets->At(i)->L2RelativeCorrectionScale() *
293     GoodSC5Jets->At(i)->L3AbsoluteCorrectionScale());
294 loizides 1.1 fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta->Fill(GoodSC5Jets->At(i)->Eta());
295     }
296     }
297    
298     }
299    
300     //--------------------------------------------------------------------------------------------------
301     void JetValidationMod::SlaveBegin()
302     {
303     // Run startup code on the computer (slave) doing the actual analysis. Here,
304     // we typically initialize histograms and other analysis objects and request
305     // branches. For this module, we request a branch of the MitTree.
306    
307     ReqBranch(fIC5GenJetName, fIC5GenJets);
308     ReqBranch(fSC5GenJetName, fSC5GenJets);
309     ReqBranch(fIC5JetName, fIC5Jets);
310     ReqBranch(fSC5JetName, fSC5Jets);
311    
312 loizides 1.2 fIC5GenJetRecoJetDeltaR =
313     new TH1D("hIC5GenJetRecoJetDeltaR",";IC5GenJetRecoJetDeltaR;#",100,0.,1.0);
314     fIC5GenJetRecoJetDeltaEta =
315     new TH1D("hIC5GenJetRecoJetDeltaEta",";IC5GenJetRecoJetDeltaEta;#",100,0.,0.5);
316     fIC5GenJetRecoJetDeltaPhi =
317     new TH1D("hIC5GenJetRecoJetDeltaPhi",";IC5GenJetRecoJetDeltaPhi;#",100,0.,0.5);
318 loizides 1.1
319 loizides 1.3 fIC5JetResponseVsGenJetPtInBarrel = new TH2D("hIC5JetResponseVsGenJetPtInBarrel",
320     ";IC5JetResponseVsGenJetPtInBarrel;#",
321     500,0.,5.0,5000,0,5000);
322     fIC5JetResponseVsGenJetPtInEndcap = new TH2D("hIC5JetResponseVsGenJetPtInEndcap",
323     ";IC5JetResponseVsGenJetPtInEndcap;#",
324     500,0.,5.0,5000,0,5000);
325     fIC5JetResponseVsGenJetPtForward = new TH2D("hIC5JetResponseVsGenJetPtForward",
326     ";IC5JetResponseVsGenJetPtForward;#",
327     500,0.,5.0,5000,0,5000);
328     fIC5JetResponseVsCaloJetEta = new TH2D("hIC5JetResponseVsCaloJetEta",
329     ";IC5JetResponseVsCaloJetEta;#",500,0.,5.0,100,-5,5);
330 loizides 1.2 fIC5CentralGenJetRecoJetDeltaR =
331     new TH1D("hIC5CentralGenJetRecoJetDeltaR",";IC5CentralGenJetRecoJetDeltaR;#",100,0.,1.0);
332 loizides 1.3 fIC5ForwardGenJetRecoJetDeltaR = new TH1D("hIC5ForwardGenJetRecoJetDeltaR",
333     ";IC5ForwardGenJetRecoJetDeltaR;#",100,0.,1.0);
334     fIC5NMatchedCaloJetsVsGenJetPt = new TH1D("hIC5NMatchedCaloJetsVsGenJetPt",
335     ";IC5GenJetPt;#",5000,0.,5000);
336     fIC5NMatchedCaloJetsVsGenJetEta = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta",
337     ";IC5GenJetEta;#",100,-5.0,5.0);
338     fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt20To30",
339     ";IC5GenJetEta;#",100,-5.0,5.0);
340     fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt30To40",
341     ";IC5GenJetEta;#",100,-5.0,5.0);
342     fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80 = new TH1D("hIC5NMatchedCaloJetsVsGenJetEta_Pt60To80",
343     ";IC5GenJetEta;#",100,-5.0,5.0);
344 loizides 1.1 fIC5NGenJetsVsGenJetPt = new TH1D("hIC5NGenJetsVsGenJetPt",";IC5GenJetPt;#",5000,0.,5000);
345     fIC5NGenJetsVsGenJetEta = new TH1D("hIC5NGenJetsVsGenJetEta",";IC5GenJetEta;#",100,-5.,5.0);
346 loizides 1.3 fIC5NGenJetsVsGenJetEta_Pt20To30 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt20To30",
347     ";IC5GenJetEta;#",100,-5.,5.0);
348     fIC5NGenJetsVsGenJetEta_Pt30To40 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt30To40",
349     ";IC5GenJetEta;#",100,-5.,5.0);
350     fIC5NGenJetsVsGenJetEta_Pt60To80 = new TH1D("hIC5NGenJetsVsGenJetEta_Pt60To80",
351     ";IC5GenJetEta;#",100,-5.,5.0);
352 loizides 1.1 fIC5CaloJetsPt = new TH1D("hIC5CaloJetsPt",";IC5CaloJetPt;#",5000,0.,5000);
353     fIC5CaloJetsEta = new TH1D("hIC5CaloJetsEta",";IC5CaloJetEta;#",100,-5.0,5.0);
354 loizides 1.3 fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt =
355     new TH1D("hIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt",
356     ";IC5CorrectedCaloJetPt;#",5000,0.0,5000.0);
357     fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta =
358     new TH1D("hIC5NUnmatchedCalojetsVsCorrectedCaloJetEta",
359     ";IC5CorrectedCaloJetEta;#",100,-5.0,5.0);
360     fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt =
361     new TH2D("hIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt",
362     ";IC5GenJetPt;#",1000,0,1000.0, 100,0.0,5.0);
363 loizides 1.1 fIC5NCaloJets = new TH1D("hIC5NCaloJets",";IC5NCaloJets;#",100,0.,15);
364 loizides 1.3 fSC5GenJetRecoJetDeltaR = new TH1D("hSC5GenJetRecoJetDeltaR",
365     ";SC5GenJetRecoJetDeltaR;#",100,0.,1.0);
366     fSC5GenJetRecoJetDeltaEta = new TH1D("hSC5GenJetRecoJetDeltaEta",
367     ";SC5GenJetRecoJetDeltaEta;#",100,0.,0.5);
368     fSC5GenJetRecoJetDeltaPhi = new TH1D("hSC5GenJetRecoJetDeltaPhi",
369     ";SC5GenJetRecoJetDeltaPhi;#",100,0.,0.5);
370     fSC5JetResponseVsGenJetPtInBarrel =
371     new TH2D("hSC5JetResponseVsGenJetPtInBarrel",
372     ";SC5JetResponseVsGenJetPtInBarrel;#",500,0.,5.0,5000,0,5000);
373     fSC5JetResponseVsGenJetPtInEndcap =
374     new TH2D("hSC5JetResponseVsGenJetPtInEndcap",
375     ";SC5JetResponseVsGenJetPtInEndcap;#",500,0.,5.0,5000,0,5000);
376     fSC5JetResponseVsGenJetPtForward =
377     new TH2D("hSC5JetResponseVsGenJetPtForward",
378     ";SC5JetResponseVsGenJetPtForward;#",500,0.,5.0,5000,0,5000);
379    
380     fSC5JetResponseVsCaloJetEta = new TH2D("hSC5JetResponseVsCaloJetEta",
381     ";SC5JetResponseVsCaloJetEta;#",500,0.,5.0,100,-5,5);
382     fSC5CentralGenJetRecoJetDeltaR =
383     new TH1D("hSC5CentralGenJetRecoJetDeltaR",
384     ";SC5CentralGenJetRecoJetDeltaR;#",100,0.,1.0);
385     fSC5ForwardGenJetRecoJetDeltaR =
386     new TH1D("hSC5ForwardGenJetRecoJetDeltaR",
387     ";SC5ForwardGenJetRecoJetDeltaR;#",100,0.,1.0);
388     fSC5NMatchedCaloJetsVsGenJetPt =
389     new TH1D("hSC5NMatchedCaloJetsVsGenJetPt",";SC5GenJetPt;#",5000,0.,5000);
390     fSC5NMatchedCaloJetsVsGenJetEta =
391     new TH1D("hSC5NMatchedCaloJetsVsGenJetEta",";SC5GenJetEta;#",100,-5.0,5.0);
392     fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30 =
393     new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt20To30",";SC5GenJetEta;#",100,-5.0,5.0);
394     fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40 =
395     new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt30To40",";SC5GenJetEta;#",100,-5.0,5.0);
396     fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80 =
397     new TH1D("hSC5NMatchedCaloJetsVsGenJetEta_Pt60To80",";SC5GenJetEta;#",100,-5.0,5.0);
398 loizides 1.1 fSC5NGenJetsVsGenJetPt = new TH1D("hSC5NGenJetsVsGenJetPt",";SC5GenJetPt;#",5000,0.,5000);
399     fSC5NGenJetsVsGenJetEta = new TH1D("hSC5NGenJetsVsGenJetEta",";SC5GenJetEta;#",100,-5.,5.0);
400 loizides 1.3 fSC5NGenJetsVsGenJetEta_Pt20To30 =
401     new TH1D("hSC5NGenJetsVsGenJetEta_Pt20To30",";SC5GenJetEta;#",100,-5.,5.0);
402     fSC5NGenJetsVsGenJetEta_Pt30To40 =
403     new TH1D("hSC5NGenJetsVsGenJetEta_Pt30To40",";SC5GenJetEta;#",100,-5.,5.0);
404     fSC5NGenJetsVsGenJetEta_Pt60To80 =
405     new TH1D("hSC5NGenJetsVsGenJetEta_Pt60To80",";SC5GenJetEta;#",100,-5.,5.0);
406 loizides 1.1 fSC5CaloJetsPt = new TH1D("hSC5CaloJetsPt",";SC5CaloJetPt;#",5000,0.,5000);
407     fSC5CaloJetsEta = new TH1D("hSC5CaloJetsEta",";SC5CaloJetEta;#",100,-5.0,5.0);
408 loizides 1.3 fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt =
409     new TH1D("hSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt",
410     ";SC5CorrectedCaloJetPt;#",5000,0.0,5000.0);
411     fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta =
412     new TH1D("hSC5NUnmatchedCalojetsVsCorrectedCaloJetEta",
413     ";SC5CorrectedCaloJetEta;#",100,-5.0,5.0);
414     fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt =
415     new TH2D("hSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt",
416     ";SC5GenJetPt;#",1000,0,1000.0, 100,0.0,5.0);
417 loizides 1.1 fSC5NCaloJets = new TH1D("hSC5NCaloJets",";SC5NCaloJets;#",100,0.,15);
418    
419     AddOutput(fIC5GenJetRecoJetDeltaR);
420     AddOutput(fIC5GenJetRecoJetDeltaEta);
421     AddOutput(fIC5GenJetRecoJetDeltaPhi);
422     AddOutput(fIC5JetResponseVsGenJetPtInBarrel);
423     AddOutput(fIC5JetResponseVsGenJetPtInEndcap);
424     AddOutput(fIC5JetResponseVsGenJetPtForward);
425     AddOutput(fIC5JetResponseVsCaloJetEta);
426     AddOutput(fIC5CentralGenJetRecoJetDeltaR);
427     AddOutput(fIC5ForwardGenJetRecoJetDeltaR);
428     AddOutput(fIC5NMatchedCaloJetsVsGenJetPt);
429     AddOutput(fIC5NMatchedCaloJetsVsGenJetEta);
430     AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt20To30);
431     AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt30To40);
432     AddOutput(fIC5NMatchedCaloJetsVsGenJetEta_Pt60To80);
433     AddOutput(fIC5NGenJetsVsGenJetPt);
434     AddOutput(fIC5NGenJetsVsGenJetEta);
435     AddOutput(fIC5NGenJetsVsGenJetEta_Pt20To30);
436     AddOutput(fIC5NGenJetsVsGenJetEta_Pt30To40);
437     AddOutput(fIC5NGenJetsVsGenJetEta_Pt60To80);
438     AddOutput(fIC5CaloJetsPt);
439     AddOutput(fIC5CaloJetsEta);
440     AddOutput(fIC5NUnmatchedCaloJetsVsCorrectedCaloJetPt);
441     AddOutput(fIC5NUnmatchedCalojetsVsCorrectedCaloJetEta);
442     AddOutput(fIC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt);
443     AddOutput(fIC5NCaloJets);
444     AddOutput(fSC5GenJetRecoJetDeltaR);
445     AddOutput(fSC5GenJetRecoJetDeltaEta);
446     AddOutput(fSC5GenJetRecoJetDeltaPhi);
447     AddOutput(fSC5JetResponseVsGenJetPtInBarrel);
448     AddOutput(fSC5JetResponseVsGenJetPtInEndcap);
449     AddOutput(fSC5JetResponseVsGenJetPtForward);
450     AddOutput(fSC5JetResponseVsCaloJetEta);
451     AddOutput(fSC5CentralGenJetRecoJetDeltaR);
452     AddOutput(fSC5ForwardGenJetRecoJetDeltaR);
453     AddOutput(fSC5NMatchedCaloJetsVsGenJetPt);
454     AddOutput(fSC5NMatchedCaloJetsVsGenJetEta);
455     AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt20To30);
456     AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt30To40);
457     AddOutput(fSC5NMatchedCaloJetsVsGenJetEta_Pt60To80);
458     AddOutput(fSC5NGenJetsVsGenJetPt);
459     AddOutput(fSC5NGenJetsVsGenJetEta);
460     AddOutput(fSC5NGenJetsVsGenJetEta_Pt20To30);
461     AddOutput(fSC5NGenJetsVsGenJetEta_Pt30To40);
462     AddOutput(fSC5NGenJetsVsGenJetEta_Pt60To80);
463     AddOutput(fSC5CaloJetsPt);
464     AddOutput(fSC5CaloJetsEta);
465     AddOutput(fSC5NUnmatchedCaloJetsVsCorrectedCaloJetPt);
466     AddOutput(fSC5NUnmatchedCalojetsVsCorrectedCaloJetEta);
467     AddOutput(fSC5CorrPtCaloJetsOverGenJetsPtVsGenJetPt);
468     AddOutput(fSC5NCaloJets);
469     }