1 |
yilmaz |
1.1 |
// -*- C++ -*-
|
2 |
|
|
//
|
3 |
|
|
// Package: RecHitComparison
|
4 |
|
|
// Class: RecHitComparison
|
5 |
|
|
//
|
6 |
|
|
/**\class RecHitComparison RecHitComparison.cc CmsHi/RecHitComparison/src/RecHitComparison.cc
|
7 |
|
|
|
8 |
|
|
Description: [one line class summary]
|
9 |
|
|
|
10 |
|
|
Implementation:
|
11 |
|
|
[Notes on implementation]
|
12 |
|
|
*/
|
13 |
|
|
//
|
14 |
|
|
// Original Author: Yetkin Yilmaz
|
15 |
|
|
// Created: Tue Sep 7 11:38:19 EDT 2010
|
16 |
|
|
// $Id$
|
17 |
|
|
//
|
18 |
|
|
//
|
19 |
|
|
|
20 |
|
|
|
21 |
|
|
// system include files
|
22 |
|
|
#include <memory>
|
23 |
|
|
#include <vector>
|
24 |
|
|
#include <iostream>
|
25 |
|
|
|
26 |
|
|
// user include files
|
27 |
|
|
#include "FWCore/Framework/interface/Frameworkfwd.h"
|
28 |
|
|
#include "FWCore/Framework/interface/EDAnalyzer.h"
|
29 |
|
|
|
30 |
|
|
#include "FWCore/Framework/interface/Event.h"
|
31 |
|
|
#include "FWCore/Framework/interface/MakerMacros.h"
|
32 |
|
|
#include "FWCore/Framework/interface/ESHandle.h"
|
33 |
|
|
#include "FWCore/ParameterSet/interface/ParameterSet.h"
|
34 |
|
|
|
35 |
|
|
#include "DataFormats/DetId/interface/DetId.h"
|
36 |
|
|
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
|
37 |
|
|
#include "Geometry/Records/interface/IdealGeometryRecord.h"
|
38 |
|
|
#include "Geometry/Records/interface/CaloGeometryRecord.h"
|
39 |
|
|
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"
|
40 |
|
|
|
41 |
|
|
#include "CommonTools/UtilAlgos/interface/TFileService.h"
|
42 |
|
|
#include "FWCore/ServiceRegistry/interface/Service.h"
|
43 |
|
|
|
44 |
|
|
#include "DataFormats/Math/interface/deltaR.h"
|
45 |
|
|
|
46 |
|
|
#include "DataFormats/Common/interface/Handle.h"
|
47 |
|
|
#include "DataFormats/FWLite/interface/ChainEvent.h"
|
48 |
|
|
#include "FWCore/Utilities/interface/InputTag.h"
|
49 |
|
|
#include "DataFormats/HeavyIonEvent/interface/CentralityBins.h"
|
50 |
|
|
#include "DataFormats/CaloTowers/interface/CaloTower.h"
|
51 |
|
|
#include "DataFormats/HeavyIonEvent/interface/Centrality.h"
|
52 |
|
|
#include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
|
53 |
|
|
#include "DataFormats/PatCandidates/interface/Jet.h"
|
54 |
|
|
#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
|
55 |
|
|
|
56 |
|
|
#include "TNtuple.h"
|
57 |
|
|
|
58 |
|
|
using namespace std;
|
59 |
|
|
|
60 |
|
|
//
|
61 |
|
|
// class declaration
|
62 |
|
|
//
|
63 |
|
|
|
64 |
|
|
class RecHitComparison : public edm::EDAnalyzer {
|
65 |
|
|
public:
|
66 |
|
|
explicit RecHitComparison(const edm::ParameterSet&);
|
67 |
|
|
~RecHitComparison();
|
68 |
|
|
|
69 |
|
|
|
70 |
|
|
private:
|
71 |
|
|
virtual void beginJob() ;
|
72 |
|
|
virtual void analyze(const edm::Event&, const edm::EventSetup&);
|
73 |
|
|
virtual void endJob() ;
|
74 |
|
|
|
75 |
|
|
// ----------member data ---------------------------
|
76 |
|
|
|
77 |
|
|
|
78 |
|
|
edm::Handle<reco::Centrality> cent;
|
79 |
|
|
edm::Handle<vector<double> > ktRhos;
|
80 |
|
|
edm::Handle<vector<double> > akRhos;
|
81 |
|
|
|
82 |
|
|
edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > jets1;
|
83 |
|
|
edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > jets2;
|
84 |
|
|
|
85 |
|
|
// typedef edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > >::const_iterator EcalIterator;
|
86 |
|
|
typedef vector<EcalRecHit>::const_iterator EcalIterator;
|
87 |
|
|
|
88 |
|
|
edm::Handle<reco::CaloJetCollection> signalJets;
|
89 |
|
|
|
90 |
|
|
TNtuple* nthit;
|
91 |
|
|
TNtuple* ntjet;
|
92 |
|
|
|
93 |
|
|
double cone;
|
94 |
|
|
|
95 |
|
|
edm::Service<TFileService> fs;
|
96 |
|
|
const CentralityBins * cbins_;
|
97 |
|
|
const CaloGeometry *geo;
|
98 |
|
|
};
|
99 |
|
|
|
100 |
|
|
//
|
101 |
|
|
// constants, enums and typedefs
|
102 |
|
|
//
|
103 |
|
|
|
104 |
|
|
//
|
105 |
|
|
// static data member definitions
|
106 |
|
|
//
|
107 |
|
|
|
108 |
|
|
//
|
109 |
|
|
// constructors and destructor
|
110 |
|
|
//
|
111 |
|
|
RecHitComparison::RecHitComparison(const edm::ParameterSet& iConfig) :
|
112 |
|
|
cbins_(0),
|
113 |
|
|
geo(0),
|
114 |
|
|
cone(0.5)
|
115 |
|
|
{
|
116 |
|
|
//now do what ever initialization is needed
|
117 |
|
|
|
118 |
|
|
}
|
119 |
|
|
|
120 |
|
|
|
121 |
|
|
RecHitComparison::~RecHitComparison()
|
122 |
|
|
{
|
123 |
|
|
|
124 |
|
|
// do anything here that needs to be done at desctruction time
|
125 |
|
|
// (e.g. close files, deallocate resources etc.)
|
126 |
|
|
|
127 |
|
|
}
|
128 |
|
|
|
129 |
|
|
|
130 |
|
|
//
|
131 |
|
|
// member functions
|
132 |
|
|
//
|
133 |
|
|
|
134 |
|
|
// ------------ method called to for each event ------------
|
135 |
|
|
void
|
136 |
|
|
RecHitComparison::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
|
137 |
|
|
{
|
138 |
|
|
|
139 |
|
|
if(!cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
|
140 |
|
|
if(!geo){
|
141 |
|
|
edm::ESHandle<CaloGeometry> pGeo;
|
142 |
|
|
iSetup.get<CaloGeometryRecord>().get(pGeo);
|
143 |
|
|
geo = pGeo.product();
|
144 |
|
|
}
|
145 |
|
|
|
146 |
|
|
edm::InputTag jetTag1("ecalRecHit","EcalRecHitsEB","SIGNALONLY");
|
147 |
|
|
edm::InputTag jetTag2("ecalRecHit","EcalRecHitsEB","RECO");
|
148 |
|
|
|
149 |
|
|
edm::InputTag signalTag("iterativeConePu5CaloJets","","SIGNALONLY");
|
150 |
|
|
edm::InputTag centTag("hiCentrality","","RECO");
|
151 |
|
|
|
152 |
|
|
using namespace edm;
|
153 |
|
|
|
154 |
|
|
ev.getByLabel(centTag,cent);
|
155 |
|
|
ev.getByLabel(jetTag1,jets1);
|
156 |
|
|
ev.getByLabel(jetTag2,jets2);
|
157 |
|
|
ev.getByLabel(signalTag,signalJets);
|
158 |
|
|
|
159 |
|
|
double hf = cent->EtHFhitSum();
|
160 |
|
|
double sumet = cent->EtMidRapiditySum();
|
161 |
|
|
int run = ev.id().run();
|
162 |
|
|
run = 1;
|
163 |
|
|
int bin = cbins_->getBin(hf);
|
164 |
|
|
int margin = 0;
|
165 |
|
|
|
166 |
|
|
vector<double> fFull;
|
167 |
|
|
vector<double> f05;
|
168 |
|
|
vector<double> f1;
|
169 |
|
|
vector<double> f15;
|
170 |
|
|
vector<double> f2;
|
171 |
|
|
vector<double> f25;
|
172 |
|
|
vector<double> f3;
|
173 |
|
|
|
174 |
|
|
int njets = signalJets->size();
|
175 |
|
|
fFull.reserve(njets);
|
176 |
|
|
f05.reserve(njets);
|
177 |
|
|
f1.reserve(njets);
|
178 |
|
|
f15.reserve(njets);
|
179 |
|
|
f2.reserve(njets);
|
180 |
|
|
f25.reserve(njets);
|
181 |
|
|
f3.reserve(njets);
|
182 |
|
|
|
183 |
|
|
for(unsigned int j1 = 0 ; j1 < signalJets->size(); ++j1){
|
184 |
|
|
fFull.push_back(0);
|
185 |
|
|
f05.push_back(0);
|
186 |
|
|
f1.push_back(0);
|
187 |
|
|
f15.push_back(0);
|
188 |
|
|
f2.push_back(0);
|
189 |
|
|
f25.push_back(0);
|
190 |
|
|
f3.push_back(0);
|
191 |
|
|
}
|
192 |
|
|
|
193 |
|
|
for(unsigned int j1 = 0 ; j1 < jets1->size(); ++j1){
|
194 |
|
|
|
195 |
|
|
const EcalRecHit& jet1 = (*jets1)[j1];
|
196 |
|
|
double e1 = jet1.energy();
|
197 |
|
|
|
198 |
|
|
const GlobalPoint& pos1=geo->getPosition(jet1.id());
|
199 |
|
|
double eta1 = pos1.eta();
|
200 |
|
|
double phi1 = pos1.phi();
|
201 |
|
|
double et1 = e1*sin(pos1.theta());
|
202 |
|
|
|
203 |
|
|
double drjet = -1;
|
204 |
|
|
double jetpt = -1;
|
205 |
|
|
bool isjet = false;
|
206 |
|
|
int matchedJet = -1;
|
207 |
|
|
|
208 |
|
|
for(unsigned int j = 0 ; j < signalJets->size(); ++j){
|
209 |
|
|
const reco::CaloJet & jet = (*signalJets)[j];
|
210 |
|
|
double dr = reco::deltaR(eta1,phi1,jet.eta(),jet.phi());
|
211 |
|
|
if(dr < cone){
|
212 |
|
|
jetpt = jet.pt();
|
213 |
|
|
drjet = dr;
|
214 |
|
|
isjet = true;
|
215 |
|
|
matchedJet = j;
|
216 |
|
|
fFull[j] += et1;
|
217 |
|
|
|
218 |
|
|
if(et1 > 0.5){
|
219 |
|
|
f05[j] += et1;
|
220 |
|
|
}
|
221 |
|
|
if(et1 > 1.){
|
222 |
|
|
f1[j] += et1;
|
223 |
|
|
}
|
224 |
|
|
if(et1 > 1.5){
|
225 |
|
|
f15[j] += et1;
|
226 |
|
|
}
|
227 |
|
|
if(et1 > 2){
|
228 |
|
|
f2[j] += et1;
|
229 |
|
|
}
|
230 |
|
|
if(et1 > 2.5){
|
231 |
|
|
f25[j] += et1;
|
232 |
|
|
}
|
233 |
|
|
if(et1 > 3.){
|
234 |
|
|
f3[j] += et1;
|
235 |
|
|
}
|
236 |
|
|
}
|
237 |
|
|
}
|
238 |
|
|
|
239 |
|
|
GlobalPoint pos2;
|
240 |
|
|
double e2 = -1;
|
241 |
|
|
int siz = jets2->size();
|
242 |
|
|
|
243 |
|
|
EcalIterator hitit = jets2->find(jet1.id());
|
244 |
|
|
if(hitit != jets2->end()){
|
245 |
|
|
e2 = hitit->energy();
|
246 |
|
|
pos2=geo->getPosition(hitit->id());
|
247 |
|
|
}else{
|
248 |
|
|
e2 = 0;
|
249 |
|
|
pos2 = pos1;
|
250 |
|
|
}
|
251 |
|
|
|
252 |
|
|
double eta2 = pos2.eta();
|
253 |
|
|
double phi2 = pos2.eta();
|
254 |
|
|
double et2 = e2*sin(pos2.theta());
|
255 |
|
|
|
256 |
|
|
if(isjet) nthit->Fill(e1,et1,e2,et2,eta2,phi2,sumet,hf,bin,jetpt,drjet);
|
257 |
|
|
|
258 |
|
|
}
|
259 |
|
|
|
260 |
|
|
for(unsigned int j1 = 0 ; j1 < signalJets->size(); ++j1){
|
261 |
|
|
const reco::CaloJet & jet = (*signalJets)[j1];
|
262 |
|
|
double em = (jet.emEnergyInEB() + jet.emEnergyInEE()) * sin(jet.theta());
|
263 |
|
|
double emf = jet.emEnergyFraction();
|
264 |
|
|
double pt = jet.pt();
|
265 |
|
|
double eta = jet.eta();
|
266 |
|
|
ntjet->Fill(bin,pt,eta,fFull[j1],f05[j1],f1[j1],f15[j1],f2[j1],f25[j1],f3[j1],em,emf);
|
267 |
|
|
}
|
268 |
|
|
|
269 |
|
|
}
|
270 |
|
|
|
271 |
|
|
|
272 |
|
|
// ------------ method called once each job just before starting event loop ------------
|
273 |
|
|
void
|
274 |
|
|
RecHitComparison::beginJob()
|
275 |
|
|
{
|
276 |
|
|
nthit = fs->make<TNtuple>("nthit","","e1:et1:e2:et2:eta:phi:sumet:hf:bin:ptjet:drjet");
|
277 |
|
|
ntjet = fs->make<TNtuple>("ntjet","","bin:pt:eta:ethit:f05:f1:f15:f2:f25:f3:em:emf");
|
278 |
|
|
|
279 |
|
|
}
|
280 |
|
|
|
281 |
|
|
// ------------ method called once each job just after ending the event loop ------------
|
282 |
|
|
void
|
283 |
|
|
RecHitComparison::endJob() {
|
284 |
|
|
}
|
285 |
|
|
|
286 |
|
|
//define this as a plug-in
|
287 |
|
|
DEFINE_FWK_MODULE(RecHitComparison);
|