ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MPIAnalyzer/src/MPIntuple.cc
Revision: 1.37
Committed: Wed Apr 27 18:46:34 2011 UTC (14 years ago) by naodell
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Changes since 1.36: +5 -5 lines
Error occurred while calculating annotation data.
Log Message:
Fixed bug in jet-vertex association

File Contents

# Content
1 // Original Author: "Nathaniel Odell"
2 // Secondary Author: Steven Won
3 // With contributions from: Andrey Pozdnyakov
4 // Created: Thurs April 22 2010
5 // $Id: MPIntuple.cc,v 1.36 2011/04/27 15:54:46 naodell Exp $
6
7 // system include files
8 #include <memory>
9 #include <string>
10
11 // user include files
12 #include "FWCore/Framework/interface/Frameworkfwd.h"
13 #include "FWCore/Framework/interface/EDAnalyzer.h"
14 #include "FWCore/Framework/interface/ESHandle.h"
15 #include "FWCore/Framework/interface/EventSetup.h"
16 #include "FWCore/Framework/interface/Event.h"
17 #include "FWCore/Framework/interface/MakerMacros.h"
18 #include "FWCore/Framework/interface/LuminosityBlock.h"
19 #include "FWCore/ParameterSet/interface/ParameterSet.h"
20
21 #include "Geometry/Records/interface/CaloGeometryRecord.h"
22 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
23 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
24 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
25 #include "DataFormats/GeometryVector/interface/LocalVector.h"
26 #include "DataFormats/Math/interface/Point3D.h"
27 #include "DataFormats/Math/interface/Vector3D.h"
28 #include "DataFormats/Math/interface/LorentzVector.h"
29 #include "DataFormats/Math/interface/deltaPhi.h"
30
31 // Libraries for objects
32 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
33 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
34 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
35 #include "DataFormats/JetReco/interface/CaloJet.h"
36 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
37 #include "DataFormats/JetReco/interface/GenJet.h"
38 #include "DataFormats/JetReco/interface/GenJetCollection.h"
39 #include "DataFormats/JetReco/interface/PFJet.h"
40 #include "DataFormats/JetReco/interface/PFJetCollection.h"
41 #include "DataFormats/METReco/interface/PFMET.h"
42 #include "DataFormats/METReco/interface/PFMETCollection.h"
43 #include "DataFormats/METReco/interface/PFMETFwd.h"
44 #include "DataFormats/METReco/interface/MET.h"
45 #include "DataFormats/METReco/interface/METCollection.h"
46 #include "DataFormats/METReco/interface/METFwd.h"
47 #include "DataFormats/MuonReco/interface/Muon.h"
48 #include "DataFormats/MuonReco/interface/MuonFwd.h"
49 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
50 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
51 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
52 #include "DataFormats/EgammaCandidates/interface/Photon.h"
53 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
54 #include "DataFormats/TrackReco/interface/Track.h"
55 #include "DataFormats/TrackReco/interface/TrackFwd.h"
56 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
57 #include "DataFormats/VertexReco/interface/Vertex.h"
58 #include "DataFormats/VertexReco/interface/VertexFwd.h"
59 #include "DataFormats/BTauReco/interface/JetTag.h"
60 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
61 #include "DataFormats/Luminosity/interface/LumiSummary.h"
62 #include "DataFormats/Common/interface/ValueMap.h"
63 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
64 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
65 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
66 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
67
68 //#include "RecoVertex/PrimaryVertexProducer/interface/VertexHigherPtSquared.h"
69
70 // JEC
71 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
72 #include "JetMETCorrections/Objects/interface/JetCorrectionsRecord.h"
73 #include "CondFormats/JetMETObjects/interface/JetCorrectionUncertainty.h"
74 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
75
76 // ntuple storage classes
77 #include "TCJet.h"
78 #include "TCMET.h"
79 #include "TCPrimaryVtx.h"
80 #include "TCGenJet.h"
81 #include "TCMuon.h"
82 #include "TCElectron.h"
83
84 // Need for HLT trigger info:
85 #include "FWCore/Common/interface/TriggerNames.h"
86 #include "DataFormats/Common/interface/TriggerResults.h"
87 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
88
89 //Root stuff
90 #include "TROOT.h"
91 #include "TH1.h"
92 #include "TH2.h"
93 #include "TProfile.h"
94 #include "TFile.h"
95 #include "TTree.h"
96 #include "TString.h"
97 #include "TObject.h"
98 #include "TObjArray.h"
99 #include "TClonesArray.h"
100 #include "TRefArray.h"
101 #include "TLorentzVector.h"
102 #include "TVector3.h"
103
104 //MC stuff
105 //#include "HepMC/GenEvent"
106
107 using namespace edm;
108 using namespace std;
109 using namespace reco;
110
111 //
112 // class declaration
113 //
114
115 class MPIntuple : public edm::EDAnalyzer {
116 public:
117 explicit MPIntuple(const edm::ParameterSet&);
118 ~MPIntuple();
119
120
121 private:
122 virtual void beginJob() ;
123 virtual void beginRun(const edm::Run&, const edm::EventSetup&) ;
124 virtual void analyze(const edm::Event&, const edm::EventSetup&);
125 virtual void endLuminosityBlock(const edm::LuminosityBlock&,const edm::EventSetup&);
126 virtual void endRun(const edm::Run&, const edm::Event&, const edm::EventSetup&);
127 virtual void endJob() ;
128
129 bool genParticleMatch(const HepMC::GenParticle& genParticle, const reco::Candidate& motherParticle);
130 bool triggerDecision(edm::Handle<edm::TriggerResults>& hltR, int iTrigger);
131 float sumPtSquared(const Vertex& v);
132
133 // ----------member data ---------------------------
134
135
136 int eventNumber, runNumber, lumiSection, bunchCross;
137 float ptHat, qScale, crossSection, evtWeight;
138 float lumiDeadCount, lumiLiveFrac, intDelLumi;
139
140 TTree* sTree;
141 TFile* ntupleFile;
142 edm::InputTag recoJetTag_;
143 edm::InputTag recoMETTag_;
144 edm::InputTag genJetTag_;
145 edm::InputTag muonTag_;
146 edm::InputTag electronTag_;
147 edm::InputTag primaryVtxTag_;
148 edm::InputTag triggerResultsTag_;
149 edm::InputTag electronIDMap_;
150 bool doGenJets_;
151 bool doPFJets_;
152 bool triggerHLT_;
153 bool isRealData;
154 int nJets_;
155 string rootfilename;
156
157 //Physics object containers
158 TClonesArray* recoJets;
159 TClonesArray* recoMET;
160 TClonesArray* genJets;
161 TClonesArray* primaryVtx;
162 TClonesArray* recoMuons;
163 TClonesArray* recoElectrons;
164
165 //GenParticles
166 TClonesArray* hardPartonP4;
167 int partonPdgId[4];
168
169 //Triggers
170 string hlTriggerResults_, hltProcess_, triggerName_;
171 TriggerNames triggerNames;
172 HLTConfigProvider hltConfig_;
173 vector<string> hlNames;
174 vector<string> mpiTriggers_;
175 unsigned int triggerStatus;
176 unsigned int hltPrescale[32];
177
178
179 //Histograms
180 TH1D * h1_ptHat;
181
182 TH1D * h1_fracAssociatedTracks;
183 TH1D * h1_meanJetTrackZ;
184 TH1D * h1_trackDxy;
185 TH1D * h1_jetVertexZ;
186 TH1D * h1_associatedSumPt;
187 TH1D * h1_associatedVertexIndex;
188 TH2F * h2_nAssTracksVsJetPt;
189 TProfile * p1_nVtcs;
190
191 TH1D * h1_partonJetdR;
192 TH1D * h1_partonJet2dR;
193 TH1D * h1_partonJetdpT;
194 TH1D * h1_partonJet2dpT;
195 TH1D * h1_partonJetMatch;
196 TH1D * h1_partonJetMatch_dR5;
197 TH1D * h1_partonJetMatch_dR9;
198
199 };
200
201 MPIntuple::MPIntuple(const edm::ParameterSet& iConfig)
202 {
203 recoJetTag_ = iConfig.getUntrackedParameter<edm::InputTag>("RecoJetTag");
204 recoMETTag_ = iConfig.getUntrackedParameter<edm::InputTag>("RecoMETTag");
205 muonTag_ = iConfig.getUntrackedParameter<edm::InputTag>("MuonTag");
206 electronTag_ = iConfig.getUntrackedParameter<edm::InputTag>("ElectronTag");
207 genJetTag_ = iConfig.getUntrackedParameter<edm::InputTag>("GenJetTag");
208 primaryVtxTag_ = iConfig.getUntrackedParameter<edm::InputTag>("PrimaryVtxTag");
209 electronIDMap_ = iConfig.getParameter<edm::InputTag>("electronIDMap");
210 nJets_ = iConfig.getUntrackedParameter<int>("nJets");
211 doGenJets_ = iConfig.getUntrackedParameter<bool>("doGenJets");
212 doPFJets_ = iConfig.getUntrackedParameter<bool>("doPFJets");
213 triggerHLT_ = iConfig.getUntrackedParameter<bool>("triggerHLT");
214 hlTriggerResults_ = iConfig.getUntrackedParameter<string>("HLTriggerResults","TriggerResults");
215 hltProcess_ = iConfig.getUntrackedParameter<string>("hltName");
216 mpiTriggers_ = iConfig.getUntrackedParameter<vector<string> >("triggers");
217 rootfilename = iConfig.getUntrackedParameter<string>("rootfilename");
218 }
219
220 MPIntuple::~MPIntuple()
221 {
222
223 }
224
225 //
226 // member functions
227 //
228
229 // ------------ method called to for each event ------------
230 void MPIntuple::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
231 {
232
233 eventNumber = iEvent.id().event();
234 runNumber = iEvent.id().run();
235 lumiSection = (unsigned int)iEvent.getLuminosityBlock().luminosityBlock();
236 bunchCross = (unsigned int)iEvent.bunchCrossing();
237 isRealData = iEvent.isRealData();
238
239 edm::Handle<reco::BeamSpot> beamSpotHandle;
240 iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
241 reco::BeamSpot vertexBeamSpot = *beamSpotHandle;
242
243 int pfCount = 0;
244 int genCount = 0;
245 int vtxCount = 0;
246 float primaryVertexZ = -999;
247
248 //////////////////////////
249 //Get vertex information//
250 //////////////////////////
251
252 Handle<reco::VertexCollection> primaryVtcs;
253 iEvent.getByLabel(primaryVtxTag_, primaryVtcs);
254
255 for(VertexCollection::const_iterator vtx_iter = primaryVtcs->begin(); vtx_iter!= primaryVtcs->end(); ++vtx_iter){
256 reco::Vertex myVtx = reco::Vertex(*vtx_iter);
257 if(!myVtx.isValid() || myVtx.isFake()) continue;
258 TCPrimaryVtx* vtxCon = new ((*primaryVtx)[vtxCount]) TCPrimaryVtx;
259 vtxCon->SetPosition(myVtx.x(), myVtx.y(), myVtx.z());
260 vtxCon->SetNDof(myVtx.ndof());
261 vtxCon->SetChi2(myVtx.chi2());
262 vtxCon->SetNTrks(myVtx.tracksSize());
263 vtxCon->SetSumPt2Trks(sumPtSquared(myVtx));
264 if(vtxCount == 0) primaryVertexZ = myVtx.z();
265 ++vtxCount;
266 }
267
268 p1_nVtcs->Fill(runNumber, vtxCount);
269
270 ///////////////////////
271 //get jet information//
272 ///////////////////////
273
274 if(doPFJets_){
275
276 edm::Handle<reco::JetTagCollection> bTagHandle1;
277 iEvent.getByLabel("trackCountingHighEffBJetTags", bTagHandle1);
278 const reco::JetTagCollection & bTags1 = *(bTagHandle1.product());
279
280 edm::Handle<reco::JetTagCollection> bTagHandle2;
281 iEvent.getByLabel("trackCountingHighPurBJetTags", bTagHandle2);
282 const reco::JetTagCollection & bTags2 = *(bTagHandle2.product());
283 reco::JetTagCollection::const_iterator jet_it_1, jet_it_2;
284
285 edm::ESHandle<JetCorrectorParametersCollection> JetCorParColl;
286 iSetup.get<JetCorrectionsRecord>().get("AK5PF",JetCorParColl);
287 JetCorrectorParameters const & JetCorPar = (*JetCorParColl)["Uncertainty"];
288 JetCorrectionUncertainty *jecUnc = new JetCorrectionUncertainty(JetCorPar);
289
290 //const JetCorrector* correctorL1 = JetCorrector::getJetCorrector("ak5PFL1Offset",iSetup);
291 const JetCorrector* correctorL2 = JetCorrector::getJetCorrector("ak5PFL2Relative",iSetup);
292 const JetCorrector* correctorL3 = JetCorrector::getJetCorrector("ak5PFL3Absolute",iSetup);
293 const JetCorrector* correctorRes = JetCorrector::getJetCorrector("ak5PFResidual", iSetup);
294
295 Handle<reco::PFJetCollection> PFJets;
296 iEvent.getByLabel(recoJetTag_, PFJets);
297
298
299 for (PFJetCollection::const_iterator jet_iter = PFJets->begin(); jet_iter!= PFJets->end(); ++jet_iter) {
300
301 reco::PFJet myJet = reco::PFJet(*jet_iter);
302 reco::PFJet corJet = reco::PFJet(*jet_iter);
303
304 float scale2 = correctorL2->correction(corJet);
305 corJet.scaleEnergy(scale2);
306 float scale3 = correctorL3->correction(corJet);
307 corJet.scaleEnergy(scale3);
308
309 if (corJet.pt() > 5) {
310
311 TCJet* jetCon = new ((*recoJets)[pfCount]) TCJet;
312
313 jetCon->SetP4(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
314 jetCon->SetVtx(-999.0, -999.0, -999.0);
315 jetCon->SetChHadFrac(myJet.chargedHadronEnergyFraction());
316 jetCon->SetNeuHadFrac(myJet.neutralHadronEnergyFraction());
317 jetCon->SetChEmFrac(myJet.chargedEmEnergyFraction());
318 jetCon->SetNeuEmFrac(myJet.neutralEmEnergyFraction());
319 jetCon->SetNumConstit(myJet.chargedMultiplicity() + myJet.neutralMultiplicity());
320 jetCon->SetNumChPart(myJet.chargedMultiplicity());
321
322 // Get b-tag information
323 for (jet_it_1 = bTags1.begin(); jet_it_1 != bTags1.end(); jet_it_1++) {
324 if (sqrt(pow(jet_it_1->first->eta() - corJet.eta(), 2) + pow(deltaPhi(jet_it_1->first->phi(),corJet.phi()), 2)) == 0.) {
325 jetCon->SetBDiscrTrkCountHiEff(jet_it_1->second);
326 break;
327 }
328 }
329
330 for (jet_it_2 = bTags2.begin(); jet_it_2 != bTags2.end(); jet_it_2++) {
331 if (sqrt(pow(jet_it_2->first->eta() - corJet.eta(), 2) + pow(deltaPhi(jet_it_2->first->phi(),corJet.phi()), 2)) == 0.) {
332 jetCon->SetBDiscrTrkCountHiPure(jet_it_2->second);
333 break;
334 }
335 }
336
337
338 //add more corrections
339
340 jetCon->SetJetCorr(2, scale2);
341 jetCon->SetJetCorr(3, scale3);
342
343 if (isRealData) {
344 float scaleRes = correctorRes->correction(corJet);
345 jetCon->SetJetCorr(4, scaleRes);
346
347 jecUnc->setJetEta(corJet.eta());
348 jecUnc->setJetPt(corJet.pt());
349 jetCon->SetUncertaintyJES(jecUnc->getUncertainty(true));
350 }
351
352 /////////////////////////
353 //get associated tracks//
354 /////////////////////////
355
356 const reco::TrackRefVector &tracks = myJet.getTrackRefs();
357
358 vector<TVector3> vtxPositionCollection;
359 vector<float> associatedTrackSumPt;
360 vector<const reco::Track*> jetTrackAddresses;
361 float sumTrackX, sumTrackY, sumTrackZ, sumTrackIP, sumTrackPt;
362 int nJetTracks, nVertexTracks, nAssociatedTracks, vertexIndex;
363 int vCount = 0;
364
365 nJetTracks = nVertexTracks = nAssociatedTracks = 0;
366 sumTrackX = sumTrackY = sumTrackZ = sumTrackIP = sumTrackPt = 0;
367
368 if(fabs(myJet.eta()) < 2.5){
369
370 for (TrackRefVector::const_iterator iTrack = tracks.begin(); iTrack != tracks.end(); ++iTrack) {
371 const reco::Track &myJetTrack = **iTrack;
372
373 sumTrackPt += myJetTrack.pt();
374 sumTrackX += myJetTrack.vx();
375 sumTrackY += myJetTrack.vy();
376 sumTrackZ += myJetTrack.vz();
377 sumTrackIP += myJetTrack.dxy(vertexBeamSpot.position());
378 jetTrackAddresses.push_back(&myJetTrack);
379 ++nJetTracks;
380 }
381
382 if (nJetTracks > 0) {
383 jetCon->SetVtx(sumTrackX/nJetTracks, sumTrackY/nJetTracks, sumTrackZ/nJetTracks);
384 h1_meanJetTrackZ->Fill(sumTrackZ/nJetTracks);
385 h1_trackDxy->Fill(sumTrackIP/nJetTracks);
386 }
387 if(jetTrackAddresses.size() > 0){
388
389 for (VertexCollection::const_iterator vtx_iter = primaryVtcs->begin(); vtx_iter!= primaryVtcs->end(); ++vtx_iter) {
390 reco::Vertex myVtx = reco::Vertex(*vtx_iter);
391 if(!myVtx.isValid() || myVtx.isFake()) continue;
392 TVector3 *iVtxPosition = new TVector3(myVtx.x(), myVtx.y(), myVtx.z());
393 vtxPositionCollection.push_back(*iVtxPosition);
394 associatedTrackSumPt.push_back(0);
395
396 for(Vertex::trackRef_iterator iTrackRef = myVtx.tracks_begin(); iTrackRef != myVtx.tracks_end(); ++iTrackRef){
397 const edm::RefToBase<reco::Track> &myTrackRef = *iTrackRef;
398 if(myTrackRef.isAvailable()){
399 const reco::Track &myVertexTrack = *myTrackRef.get();
400
401 for(vector<const reco::Track*>::const_iterator iTrackAddress = jetTrackAddresses.begin(); iTrackAddress != jetTrackAddresses.end(); ++iTrackAddress){
402 if (*iTrackAddress == &myVertexTrack) {
403 associatedTrackSumPt.at(vCount) += myVertexTrack.pt()/sumTrackPt;
404 ++nAssociatedTracks;
405 }
406 }
407 }
408 }
409 ++vCount;
410 }
411
412 float maxSumPtFraction = 0;
413 vCount = vertexIndex = 0;
414
415 for (vector<float>::const_iterator iTrackSumPt = associatedTrackSumPt.begin(); iTrackSumPt != associatedTrackSumPt.end(); ++iTrackSumPt) {
416 if (*iTrackSumPt > maxSumPtFraction) {
417 maxSumPtFraction = *iTrackSumPt;
418 vertexIndex = vCount + 1;
419 }
420 ++vCount;
421 }
422 if (vertexIndex > 0) {
423 h1_jetVertexZ->Fill(vtxPositionCollection[vertexIndex-1].z());
424 }
425
426 jetCon->SetVtxSumPtFrac(maxSumPtFraction);
427 jetCon->SetVtxSumPt(sumTrackPt);
428 jetCon->SetVtxTrackFrac((float)nAssociatedTracks/(float)nJetTracks);
429 jetCon->SetVtxNTracks(nJetTracks);
430 jetCon->SetVtxIndex(vertexIndex);
431 }
432 }
433 ++pfCount;
434 }
435 }
436 }
437
438
439
440 //////////////////
441 // Get MET info //
442 //////////////////
443
444
445 Handle<PFMETCollection> MET;
446 iEvent.getByLabel(recoMETTag_, MET);
447
448 int metCount = 0;
449 for (PFMETCollection::const_iterator iMET = MET->begin(); iMET != MET->end(); ++iMET) {
450 TCMET* metCon = new ((*recoMET)[metCount]) TCMET;
451 metCon->SetSumEt(iMET->sumEt());
452 metCon->SetMet(iMET->et());
453 metCon->SetPhi(iMET->phi());
454 metCon->SetPhotonEtFraction(iMET->photonEtFraction());
455 metCon->SetElectronEtFraction(iMET->electronEtFraction());
456 metCon->SetMuonEtFraction(iMET->muonEtFraction());
457 metCon->SetNeutralHadronEtFraction(iMET->neutralHadronEtFraction());
458 metCon->SetChargedHadronEtFraction(iMET->chargedHadronEtFraction());
459 metCon->SetHFHadronEtFraction(iMET->HFHadronEtFraction());
460 metCon->SetHFEMEtFraction(iMET->HFEMEtFraction());
461 }
462
463 ///////////////////
464 // Get muon info //
465 ///////////////////
466
467
468 Handle<MuonCollection> muons;
469 iEvent.getByLabel(muonTag_, muons);
470
471 int muCount = 0;
472 for (MuonCollection::const_iterator mu = muons->begin(); mu != muons->end(); ++mu) {
473 if (mu->isGlobalMuon() && mu->pt() > 8){
474 TCMuon* lepCon = new ((*recoMuons)[muCount]) TCMuon;
475 lepCon->Setp4(mu->px(), mu->py(), mu->pz(), mu->p());
476 lepCon->SetVtx(mu->globalTrack()->vx(),mu->globalTrack()->vy(),mu->globalTrack()->vz());
477 lepCon->SetEta(mu->eta());
478 lepCon->SetPhi(mu->phi());
479 lepCon->SetCharge(mu->charge());
480 lepCon->SetisGLB(mu->isGlobalMuon());
481 lepCon->SetisTRK(mu->isTrackerMuon());
482 lepCon->Setdxy(mu->globalTrack()->dxy(vertexBeamSpot.position()));
483 lepCon->SetnPXLHits(mu->globalTrack()->hitPattern().numberOfValidPixelHits());
484 lepCon->SetnTRKHits(mu->globalTrack()->hitPattern().numberOfValidTrackerHits());
485 lepCon->SetnValidMuHits(mu->globalTrack()->hitPattern().numberOfValidMuonHits());
486 lepCon->SetnMatchSeg(mu->numberOfMatches());
487 lepCon->SetNormChi2(mu->globalTrack()->normalizedChi2());
488 lepCon->SetCaloComp(mu->caloCompatibility());
489 lepCon->SetSegComp(muon::segmentCompatibility(*mu));
490 lepCon->SetEMIso(mu->isolationR03().emEt);
491 lepCon->SetHADIso(mu->isolationR03().hadEt);
492 lepCon->SetTRKIso(mu->isolationR03().sumPt);
493 muCount++;
494 }
495 }
496
497
498 ///////////////////////
499 // Get electron info //
500 ///////////////////////
501
502
503 Handle<edm::ValueMap<float> > eIDValueMap;
504 iEvent.getByLabel( electronIDMap_ , eIDValueMap );
505 const edm::ValueMap<float> & eIDmap = * eIDValueMap ;
506
507 Handle<GsfElectronCollection> electrons;
508 iEvent.getByLabel(electronTag_, electrons);
509
510 int elCount = 0;
511 for (unsigned int i = 0; i < electrons->size(); i++){
512 edm::Ref<reco::GsfElectronCollection> electronRef(electrons,i);
513 if ( eIDmap[electronRef] == 7 && electronRef->pt() > 15){
514 TCElectron* lepCon = new ((*recoElectrons)[elCount]) TCElectron;
515 lepCon->Setp4(electronRef->px(),electronRef->py(),electronRef->pz(),electronRef->p());
516 lepCon->SetVtx(electronRef->gsfTrack()->vx(),electronRef->gsfTrack()->vy(),electronRef->gsfTrack()->vz());
517 lepCon->SetCharge(electronRef->charge());
518 lepCon->SetEta(electronRef->eta());
519 lepCon->SetPhi(electronRef->phi());
520 lepCon->Setdxy(electronRef->gsfTrack()->dxy(vertexBeamSpot.position()));
521 lepCon->SetNormChi2(electronRef->gsfTrack()->normalizedChi2());
522 lepCon->SetEMIso(electronRef->dr03EcalRecHitSumEt());
523 lepCon->SetHADIso(electronRef->dr03HcalTowerSumEt());
524 lepCon->SetTRKIso(electronRef->dr03TkSumPt());
525 elCount++;
526 }
527 }
528
529
530 ////////////////////////
531 // Get gen-level info //
532 ////////////////////////
533
534
535 if (!isRealData) {
536
537 Handle<HepMCProduct > genEvtHandle;
538 iEvent.getByLabel( "generator", genEvtHandle) ;
539 const HepMC::GenEvent* Evt = genEvtHandle->GetEvent();
540
541 Handle<GenEventInfoProduct> GenEventInfoHandle;
542 iEvent.getByLabel("generator", GenEventInfoHandle);
543
544 Handle<reco::GenJetCollection> GenJets;
545 iEvent.getByLabel(genJetTag_, GenJets);
546
547 ptHat = qScale = -1; crossSection = 0;
548
549 if (GenEventInfoHandle.isValid()) {
550 qScale = GenEventInfoHandle->qScale();
551 ptHat = (GenEventInfoHandle->hasBinningValues() ? GenEventInfoHandle->binningValues()[0] : 0.0);
552 evtWeight = GenEventInfoHandle->weight();
553
554 if (evtWeight != 1) cout<<"Event weight: "<<evtWeight<<endl;
555 //if (qScale != ptHat) cout<<"q scale: "<<qScale<<endl;
556 h1_ptHat->Fill(ptHat);
557 }
558 vector<HepMC::GenParticle*> genPartons;
559 int hardCount = 0;
560
561 //cout<<"Event: "<<Evt->signal_process_id()<<endl;
562
563 for (HepMC::GenEvent::particle_const_iterator iGenParticle = Evt->particles_begin(); iGenParticle != Evt->particles_end(); ++iGenParticle) {
564 HepMC::GenParticle *myGenPart = *iGenParticle;
565 if (myGenPart->status() == 23) {
566 new ((*hardPartonP4)[hardCount]) TLorentzVector(myGenPart->momentum().px(), myGenPart->momentum().py(), myGenPart->momentum().pz(), myGenPart->momentum().e());
567 partonPdgId[hardCount] = myGenPart->pdg_id();
568
569 //genPartons.push_back(myGenPart);
570 //float genPt = sqrt(pow(myGenPart->momentum().px(), 2) + pow(myGenPart->momentum().py(), 2));
571 ////cout<<"Hard scattering parton "<<hardCount<<" : ("<<genPt<<", "<<myGenPart->momentum().phi()<<", "<<myGenPart->momentum().eta()<<"), vtx : "<<myGenPart->production_vertex()->barcode()<<", pdgId : "<<myGenPart->pdg_id()<<endl;
572 //
573 //for(HepMC::GenVertex::particle_iterator iGenDescend = myGenPart->end_vertex()->particles_begin(HepMC::descendants);
574 // iGenDescend != myGenPart->end_vertex()->particles_end(HepMC::descendants); ++iGenDescend) {
575 // HepMC::GenParticle *myGenDescend = *iGenDescend;
576 // if (myGenDescend->status() == 71){
577 // float descPt = sqrt(pow(myGenDescend->momentum().px(), 2) + pow(myGenDescend->momentum().py(), 2));
578 // //cout<<"\t End partons"<<" : ("<<descPt<<", "<<myGenDescend->momentum().phi()<<", "<<myGenDescend->momentum().eta()<<"), vtx : "<<", pdgId : "<<myGenDescend->pdg_id()<<endl;
579 // }
580 //}
581 //cout<<hardCount<<endl;
582 ++hardCount;
583 }
584 }
585
586 for (GenJetCollection::const_iterator jet_iter = GenJets->begin(); jet_iter!= GenJets->end(); ++jet_iter) {
587 reco::GenJet myJet = reco::GenJet(*jet_iter);
588 if (myJet.pt() > 15 && fabs(myJet.eta()) < 2.5) {
589
590 TCGenJet* jetCon = new ((*genJets)[genCount]) TCGenJet;
591 jetCon->SetP4(myJet.px(), myJet.py(), myJet.pz(), myJet.energy());
592 jetCon->SetHadEnergy(myJet.hadEnergy());
593 jetCon->SetEmEnergy(myJet.emEnergy());
594 jetCon->SetInvEnergy(myJet.invisibleEnergy());
595 jetCon->SetAuxEnergy(myJet.auxiliaryEnergy());
596 jetCon->SetNumConstit(myJet.getGenConstituents().size());
597
598 //cout<<"~~JET "<<genCount+1<<"~~"<<endl;
599 //cout<<"(pt, phi, eta) = ("<<myJet.pt()<<", "<<myJet.phi()<<", "<<myJet.eta()<<")"<<endl;
600
601 int pCount, nearestIndex, dR5Index, dR9Index;
602 float nearestDR, nearest2DR, nearestDR5, nearestDR9;
603 float nearestDPT, nearest2DPT, nearestDPT5, nearestDPT9;
604 pCount = nearestIndex = dR5Index = dR9Index = 0;
605 nearestDR = nearest2DR = nearestDR5 = nearestDR9 = 10.;
606 nearestDPT = nearest2DPT = nearestDPT5 = nearestDPT9 = 15.;
607
608 for (vector<HepMC::GenParticle*>::const_iterator parton_iter = genPartons.begin(); parton_iter != genPartons.end(); ++parton_iter) {
609 HepMC::GenParticle *myParton = *parton_iter;
610 float partonPt = sqrt(pow(myParton->momentum().px(), 2) + pow(myParton->momentum().py(), 2));
611 float dPhi = myJet.phi() - myParton->momentum().phi();
612 if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
613 float dR = sqrt(pow(myJet.eta() - myParton->momentum().eta(), 2) + pow(dPhi, 2));
614
615 if (dR < nearestDR) {
616 nearest2DR = nearestDR;
617 nearestDR = dR;
618 nearest2DPT = nearestDPT;
619 nearestDPT = fabs(partonPt -myJet.pt())/myJet.pt();
620 nearestIndex = pCount + 1;
621 }
622 if (dR < 0.5) {
623 nearestDR5 = dR;
624 nearestDPT5 = fabs(partonPt -myJet.pt())/myJet.pt();
625 if (dR5Index == 0) {
626 dR5Index = pCount + 1;
627 } else {
628 dR5Index = 5;
629 }
630 }
631 if (dR < 0.9) {
632 nearestDR9 = dR;
633 nearestDPT9 = fabs(partonPt -myJet.pt())/myJet.pt();
634 if (dR9Index == 0) {
635 dR9Index = pCount + 1;
636 } else {
637 dR9Index = 5;
638 }
639 }
640 //cout<<"dR jet-parton: "<<dR<<" | dpT: "<<partonPt - myJet.pt()<<endl;
641 ++pCount;
642 }
643
644 h1_partonJetMatch->Fill(nearestIndex);
645 h1_partonJetMatch_dR5->Fill(dR5Index);
646 h1_partonJetMatch_dR9->Fill(dR9Index);
647 h1_partonJetdR->Fill(nearestDR);
648 h1_partonJet2dR->Fill(nearest2DR);
649 h1_partonJetdpT->Fill(nearestDPT);
650 h1_partonJet2dpT->Fill(nearest2DPT);
651 //h1_partonJetMatch_dR5->Fill();
652 //h1_partonJetMatch_dR9->Fill();
653
654 ++genCount;
655 }
656 }
657 //cout<<"\nEND EVENT\n"<<endl;
658 }
659
660 ///////////////////////////
661 //get trigger information//
662 ///////////////////////////
663
664 if (triggerHLT_) {
665
666 edm::Handle<TriggerResults> hltR;
667 triggerResultsTag_ = InputTag(hlTriggerResults_,"",hltProcess_);
668 iEvent.getByLabel(triggerResultsTag_,hltR);
669
670 const TriggerNames & triggerNames = iEvent.triggerNames(*hltR);
671 hlNames=triggerNames.triggerNames();
672
673 triggerStatus = 0x0;
674
675 //for(int i = 0; i < (int)mpiTriggers_.size(); ++i) hltPrescale[i] = hltConfig_.prescaleValue(iEvent, iSetup, mpiTriggers_[i]);
676
677 for (int i=0; i < (int)hlNames.size(); ++i) {
678 if (!triggerDecision(hltR, i)) continue;
679 for (int j = 0; j < (int)mpiTriggers_.size(); ++j){
680 if (hlNames[i].compare(0, mpiTriggers_[j].length(),mpiTriggers_[j]) == 0) {
681 triggerStatus |= 0x01 << j;
682 }
683 }
684 }
685 }
686
687 if((pfCount >= nJets_ || genCount >= nJets_) && vtxCount > 0) sTree -> Fill();
688
689 recoJets->Clear();
690 recoMET->Clear();
691 recoElectrons->Clear();
692 recoMuons->Clear();
693 primaryVtx->Clear();
694 genJets->Clear();
695 hardPartonP4->Clear();
696 }
697
698 // ------------ method called once each job just before starting event loop ------------
699 void MPIntuple::beginJob()
700 {
701 ntupleFile = new TFile(rootfilename.c_str(), "RECREATE");
702 sTree = new TTree("mpiTree", "Tree for Jets");
703
704 h1_ptHat = new TH1D("h1_ptHat", "ptHat", 15, 10.0, 160.0);
705 h1_fracAssociatedTracks = new TH1D("h1_fracAssociatedTracks", "Fraction of associated jet tracks", 20, 0.0, 1.0);
706 h1_trackDxy = new TH1D("h1_trackDxy", "dxy for all associated tracks", 50, -0.4, 0.4);
707 h1_meanJetTrackZ = new TH1D("h1_meanJetTrackZ", "<z> from jet tracks", 50, -25.0, 25.0);
708 h1_jetVertexZ = new TH1D("h1_jetVertexZ", "z for vertex associated to jet", 50, -25.0, 25.0);
709 h1_associatedSumPt = new TH1D("h1_associatedSumPt", "ratio of sum p_{T} between associatedTracks and jetTracks", 20, 0.0, 1.0);
710 h1_associatedVertexIndex = new TH1D("h1_associatedVertex", "number of jets associated to primary vertex", 5, -0.5, 4.5);
711 h1_partonJetMatch = new TH1D("h1_partonJetMatch", "jet association to hardest scattering", 4, 0.5, 4.5);
712 h1_partonJetMatch_dR5 = new TH1D("h1_partonJetMatch_dR5", "jet association to hardest scattering (dR < 0.5)", 6, -0.5, 5.5);
713 h1_partonJetMatch_dR9 = new TH1D("h1_partonJetMatch_dR9", "jet association to hardest scattering (dR < 0.9)", 6, -0.5, 5.5);
714 h1_partonJetdR = new TH1D("h1_partonJetdR", "jet-nearest parton dR", 50, 0., 5.);
715 h1_partonJet2dR = new TH1D("h1_partonJet2dR", "jet-2nd nearest parton dR", 50, 0., 5.);
716 h1_partonJetdpT = new TH1D("h1_partonJetDpT", "jet-nearest parton dpT", 20, 0., 1.);
717 h1_partonJet2dpT = new TH1D("h1_partonJet2DpT", "jet-2nd nearest parton dpT", 20, 0., 1.);
718 h2_nAssTracksVsJetPt = new TH2F("h2_nAssTracksVsJetPt", "nTracks vs Jet p_{T}", 50, 0, 100, 20, -0.5, 19.5);
719 p1_nVtcs = new TProfile("p1_nVtcs", "<N> per run", 5200.0, 132300.0, 137500.0, 0.0, 6.0);
720
721 recoJets = new TClonesArray("TCJet");
722 recoMET = new TClonesArray("TCMET");
723 recoElectrons = new TClonesArray("TCElectron");
724 recoMuons = new TClonesArray("TCMuon");
725 genJets = new TClonesArray("TCGenJet");
726 primaryVtx = new TClonesArray("TCPrimaryVtx");
727 hardPartonP4 = new TClonesArray("TLorentzVector");
728
729 sTree->Branch("recoJets",&recoJets, 6400, 0);
730 sTree->Branch("recoElectrons",&recoElectrons, 6400, 0);
731 sTree->Branch("recoMuons",&recoMuons, 6400, 0);
732 sTree->Branch("genJets",&genJets, 6400, 0);
733 sTree->Branch("primaryVtx",&primaryVtx, 6400, 0);
734 sTree->Branch("recoMET",&recoMET, 6400, 0);
735 sTree->Branch("hardPartonP4",&hardPartonP4, 6400, 0);
736 sTree->Branch("partonPdgId",partonPdgId, "partonPdgId[32]/i");
737 sTree->Branch("eventNumber",&eventNumber, "eventNumber/I");
738 sTree->Branch("runNumber",&runNumber, "runNumber/I");
739 sTree->Branch("lumiSection",&lumiSection, "lumiSection/I");
740 sTree->Branch("triggerStatus",&triggerStatus, "triggerStatus/i");
741 sTree->Branch("hltPrescale",hltPrescale, "hltPrescale[32]/i");
742 sTree->Branch("isRealData",&isRealData, "isRealData/i");
743 sTree->Branch("bunchCross",&bunchCross, "bunchCross/i");
744 sTree->Branch("lumiDeadCount",&lumiDeadCount, "lumiDeadCount/f");
745 sTree->Branch("lumiLiveFrac",&lumiLiveFrac, "lumiLiveFrac/f");
746 sTree->Branch("intDelLumi",&intDelLumi, "intDelLumi/f");
747 sTree->Branch("ptHat",&ptHat, "ptHat/f");
748 sTree->Branch("qScale", &qScale, "qScale/f");
749 sTree->Branch("evtWeight", &evtWeight, "evtWeight/f");
750 sTree->Branch("crossSection", &crossSection, "crossSection/f");
751 }
752
753 void MPIntuple::beginRun(const edm::Run& iRun, const edm::EventSetup& iEvent)
754 {
755 bool changed = true;
756 hltConfig_.init(iRun, iEvent, hltProcess_, changed);
757 }
758
759 void MPIntuple::endLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup& iEvent)
760 {
761 edm::Handle<LumiSummary> lumiSummary;
762 iLumi.getByLabel("lumiProducer", lumiSummary);
763
764 lumiDeadCount = lumiSummary->deadcount();
765 lumiLiveFrac = lumiSummary->liveFrac();
766 intDelLumi = lumiSummary->avgInsDelLumi()*93.244;
767
768 //cout<<iLumi.id().luminosityBlock()<<endl;
769 //cout<<"\t Dead Count = "<<lumiSummary->deadcount()<<endl;
770 //cout<<"\t Fraction of dead time = "<<1 - lumiSummary->liveFrac()<<endl;
771 //cout<<"\t Integrated luminosity = "<<lumiSummary->avgInsDelLumi()*93.244<<endl;
772 //cout<<"\t Dead time corrected luminosity = "<<lumiSummary->avgInsDelLumi()*lumiSummary->liveFrac()*93.244<<endl;
773 }
774
775 void MPIntuple::endRun(const edm::Run& iRun, const edm::Event& iEvent, const edm::EventSetup& iSetup)
776 {
777 Handle<GenRunInfoProduct> GenRunInfoHandle;
778 iEvent.getByLabel("generator", GenRunInfoHandle);
779
780 if (GenRunInfoHandle.isValid()) {
781 crossSection = GenRunInfoHandle->crossSection();
782 }
783
784 for(int i = 0; i < (int)mpiTriggers_.size(); ++i) hltPrescale[i] = hltConfig_.prescaleValue(iEvent, iSetup, mpiTriggers_[i]);
785 //for (int i = 0; i < (int)mpiTriggers_.size(); ++i) cout << mpiTriggers_[i] << " prescale = " << hltPrescale[i] <<endl;
786 }
787 // ------------ method called once each job just after ending the event loop ------------
788 void MPIntuple::endJob()
789 {
790 ntupleFile->cd();
791
792 h1_ptHat->Write();
793 h1_fracAssociatedTracks->Write();
794 h1_meanJetTrackZ->Write();
795 h1_trackDxy->Write();
796 h1_jetVertexZ->Write();
797 h1_associatedSumPt->Write();
798 h1_associatedVertexIndex->Write();
799 h1_partonJetMatch->Write();
800 h1_partonJetdR->Write();
801 h1_partonJet2dR->Write();
802 h1_partonJetdpT->Write();
803 h1_partonJet2dpT->Write();
804 h1_partonJetMatch_dR5->Write();
805 h1_partonJetMatch_dR9->Write();
806 h2_nAssTracksVsJetPt->Write();
807 p1_nVtcs->Write();
808
809 ntupleFile->Write();
810 ntupleFile->Close();
811 }
812
813 bool MPIntuple::genParticleMatch(const HepMC::GenParticle& genParticle, const reco::Candidate& motherParticle)
814 {
815 float genPhi, momPhi, genEta, momEta;
816 genPhi = genParticle.momentum().phi();
817 genEta = genParticle.momentum().eta();
818 momPhi = motherParticle.phi();
819 momEta = motherParticle.eta();
820
821 float dPhi = momPhi - genPhi;
822 if (dPhi > TMath::Pi()) dPhi = 2*TMath::Pi() - dPhi;
823 float dR = sqrt(pow((genPhi - momPhi), 2) + pow((genEta - momEta), 2));
824 if (dR < 0.00001) return true;
825 else return false;
826 }
827
828 bool MPIntuple::triggerDecision(edm::Handle<edm::TriggerResults> &hltR, int iTrigger)
829 {
830 bool triggerPassed = false;
831 if(hltR->wasrun(iTrigger) &&
832 hltR->accept(iTrigger) &&
833 !hltR->error(iTrigger) ){
834 triggerPassed = true;
835 }
836 return triggerPassed;
837 }
838
839 float MPIntuple::sumPtSquared(const Vertex & v)
840 {
841 float sum = 0.;
842 float pT;
843 for (Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); it++) {
844 pT = (**it).pt();
845 float epT=(**it).ptError(); pT=pT>epT ? pT-epT : 0;
846
847 sum += pT*pT;
848 }
849 return sum;
850 }
851
852 //define this as a plug-in
853 DEFINE_FWK_MODULE(MPIntuple);