CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
IsolatedTracksNxN.cc
Go to the documentation of this file.
1 // -*- C++ -*
2 //
3 // Package: IsolatedParticles
4 // Class: IsolatedTracksNxN
5 //
13 //
14 // Original Author: Seema Sharma
15 // Created: Mon Aug 10 15:30:40 CST 2009
16 //
17 //
18 
19 // system include files
20 #include <cmath>
21 #include <map>
22 #include <memory>
23 #include <string>
24 #include <vector>
25 
26 // user include files
27 #include <Math/GenVector/VectorUtil.h>
28 
29 // root objects
30 #include "TROOT.h"
31 #include "TSystem.h"
32 #include "TFile.h"
33 #include "TH1F.h"
34 #include "TH2F.h"
35 #include "TProfile.h"
36 #include "TDirectory.h"
37 #include "TTree.h"
38 
47 
49 //L1 trigger Menus etc
58 
63 // muons and tracks
68 // Vertices
72 // Calorimeters
79 // Trigger
87 //L1 objects
94 // Jets in the event
98 
101 
109 // TFile Service
112 
113 // ecal / hcal
121 
124 
128 
129 // SimHit + SimTrack
135 
136 // track associator
140 
141 class IsolatedTracksNxN : public edm::one::EDAnalyzer<edm::one::SharedResources> {
142 public:
143  explicit IsolatedTracksNxN(const edm::ParameterSet &);
144  ~IsolatedTracksNxN() override;
145 
146  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
147 
148 private:
149  void beginJob() override;
150  void analyze(const edm::Event &, const edm::EventSetup &) override;
151  void endJob() override;
152 
153  void printTrack(const reco::Track *pTrack);
154 
155  void bookHistograms();
156 
157  void clearTreeVectors();
158 
159 private:
160  std::unique_ptr<L1GtUtils> m_l1GtUtils;
161  static constexpr size_t nL1BitsMax = 128;
163 
164  // map of trig bit, algo name and num events passed
165  std::map<std::pair<unsigned int, std::string>, int> l1AlgoMap_;
166  std::vector<unsigned int> m_triggerMaskAlgoTrig;
167 
168  const bool doMC_, writeAllTracks_;
169  const int myverbose_;
170  const double pvTracksPtMin_;
171  const int debugTrks_;
173  const double minTrackP_, maxTrackEta_;
175  const double tMinE_, tMaxE_, tMinH_, tMaxH_;
176  bool initL1_;
178 
182 
186 
189 
193 
201 
209  static constexpr size_t NPBins = 15;
210  static constexpr size_t NEtaBins = 3;
212 
217  TH1F *h_nTracks;
224 
225  TTree *tree_;
226 
229  std::vector<std::string> *t_L1AlgoNames;
230  std::vector<int> *t_L1PreScale;
231  int t_L1Decision[128];
232 
233  std::vector<double> *t_L1CenJetPt, *t_L1CenJetEta, *t_L1CenJetPhi;
234  std::vector<double> *t_L1FwdJetPt, *t_L1FwdJetEta, *t_L1FwdJetPhi;
235  std::vector<double> *t_L1TauJetPt, *t_L1TauJetEta, *t_L1TauJetPhi;
236  std::vector<double> *t_L1MuonPt, *t_L1MuonEta, *t_L1MuonPhi;
237  std::vector<double> *t_L1IsoEMPt, *t_L1IsoEMEta, *t_L1IsoEMPhi;
239  std::vector<double> *t_L1METPt, *t_L1METEta, *t_L1METPhi;
240 
241  std::vector<double> *t_PVx, *t_PVy, *t_PVz, *t_PVTracksSumPt;
244  std::vector<int> *t_PVNTracksHP, *t_PVNTracksHPWt;
245 
246  std::vector<double> *t_jetPt, *t_jetEta, *t_jetPhi;
247  std::vector<double> *t_nTrksJetCalo, *t_nTrksJetVtx;
248 
250  std::vector<double> *t_trackPtAll;
252 
253  std::vector<double> *t_trackP, *t_trackPt, *t_trackEta, *t_trackPhi;
255  std::vector<double> *t_trackDxy, *t_trackDxyBS, *t_trackDz, *t_trackDzBS;
256  std::vector<double> *t_trackDxyPV, *t_trackDzPV;
257  std::vector<double> *t_trackChiSq;
258  std::vector<int> *t_trackPVIdx;
259 
260  std::vector<int> *t_NLayersCrossed, *t_trackNOuterHits;
261  std::vector<int> *t_trackHitsTOB, *t_trackHitsTEC;
268  std::vector<double> *t_trackOutPosOutHitDr, *t_trackL;
269 
270  std::vector<double> *t_maxNearP31x31;
271  std::vector<double> *t_maxNearP21x21;
272  std::vector<int> *t_ecalSpike11x11;
273 
274  std::vector<double> *t_e7x7, *t_e9x9, *t_e11x11, *t_e15x15;
280 
281  std::vector<double> *t_esimPdgId, *t_simTrackP, *t_trkEcalEne;
282  std::vector<double> *t_esim7x7, *t_esim9x9, *t_esim11x11, *t_esim15x15;
288 
290  std::vector<double> *t_h3x3, *t_h5x5, *t_h7x7;
291  std::vector<double> *t_h3x3Sig, *t_h5x5Sig, *t_h7x7Sig;
292  std::vector<int> *t_infoHcal;
293 
294  std::vector<double> *t_trkHcalEne;
295  std::vector<double> *t_hsim3x3, *t_hsim5x5, *t_hsim7x7;
297  std::vector<double> *t_hsim3x3Rest, *t_hsim5x5Rest, *t_hsim7x7Rest;
301 };
302 
303 static const bool useL1EventSetup(true);
304 static const bool useL1GtTriggerMenuLite(true);
305 
307  : trackerHitAssociatorConfig_(consumesCollector()),
308  doMC_(iConfig.getUntrackedParameter<bool>("doMC", false)),
309  writeAllTracks_(iConfig.getUntrackedParameter<bool>("writeAllTracks", false)),
310  myverbose_(iConfig.getUntrackedParameter<int>("verbosity", 5)),
311  pvTracksPtMin_(iConfig.getUntrackedParameter<double>("pvTracksPtMin", 1.0)),
312  debugTrks_(iConfig.getUntrackedParameter<int>("debugTracks", 0)),
313  printTrkHitPattern_(iConfig.getUntrackedParameter<bool>("printTrkHitPattern", false)),
314  minTrackP_(iConfig.getUntrackedParameter<double>("minTrackP", 1.0)),
315  maxTrackEta_(iConfig.getUntrackedParameter<double>("maxTrackEta", 5.0)),
316  debugL1Info_(iConfig.getUntrackedParameter<bool>("debugL1Info", false)),
317  L1TriggerAlgoInfo_(iConfig.getUntrackedParameter<bool>("l1TriggerAlgoInfo", false)),
318  tMinE_(iConfig.getUntrackedParameter<double>("timeMinCutECAL", -500.)),
319  tMaxE_(iConfig.getUntrackedParameter<double>("timeMaxCutECAL", 500.)),
320  tMinH_(iConfig.getUntrackedParameter<double>("timeMinCutHCAL", -500.)),
321  tMaxH_(iConfig.getUntrackedParameter<double>("timeMaxCutHCAL", 500.)),
322  t_L1AlgoNames(nullptr),
323  t_L1PreScale(nullptr),
324  t_L1CenJetPt(nullptr),
325  t_L1CenJetEta(nullptr),
326  t_L1CenJetPhi(nullptr),
327  t_L1FwdJetPt(nullptr),
328  t_L1FwdJetEta(nullptr),
329  t_L1FwdJetPhi(nullptr),
330  t_L1TauJetPt(nullptr),
331  t_L1TauJetEta(nullptr),
332  t_L1TauJetPhi(nullptr),
333  t_L1MuonPt(nullptr),
334  t_L1MuonEta(nullptr),
335  t_L1MuonPhi(nullptr),
336  t_L1IsoEMPt(nullptr),
337  t_L1IsoEMEta(nullptr),
338  t_L1IsoEMPhi(nullptr),
339  t_L1NonIsoEMPt(nullptr),
340  t_L1NonIsoEMEta(nullptr),
341  t_L1NonIsoEMPhi(nullptr),
342  t_L1METPt(nullptr),
343  t_L1METEta(nullptr),
344  t_L1METPhi(nullptr),
345  t_PVx(nullptr),
346  t_PVy(nullptr),
347  t_PVz(nullptr),
348  t_PVTracksSumPt(nullptr),
349  t_PVTracksSumPtWt(nullptr),
350  t_PVTracksSumPtHP(nullptr),
351  t_PVTracksSumPtHPWt(nullptr),
352  t_PVisValid(nullptr),
353  t_PVNTracks(nullptr),
354  t_PVNTracksWt(nullptr),
355  t_PVndof(nullptr),
356  t_PVNTracksHP(nullptr),
357  t_PVNTracksHPWt(nullptr),
358  t_jetPt(nullptr),
359  t_jetEta(nullptr),
360  t_jetPhi(nullptr),
361  t_nTrksJetCalo(nullptr),
362  t_nTrksJetVtx(nullptr),
363  t_trackPAll(nullptr),
364  t_trackEtaAll(nullptr),
365  t_trackPhiAll(nullptr),
366  t_trackPdgIdAll(nullptr),
367  t_trackPtAll(nullptr),
368  t_trackDxyAll(nullptr),
369  t_trackDzAll(nullptr),
370  t_trackDxyPVAll(nullptr),
371  t_trackDzPVAll(nullptr),
372  t_trackChiSqAll(nullptr),
373  t_trackP(nullptr),
374  t_trackPt(nullptr),
375  t_trackEta(nullptr),
376  t_trackPhi(nullptr),
377  t_trackEcalEta(nullptr),
378  t_trackEcalPhi(nullptr),
379  t_trackHcalEta(nullptr),
380  t_trackHcalPhi(nullptr),
381  t_trackDxy(nullptr),
382  t_trackDxyBS(nullptr),
383  t_trackDz(nullptr),
384  t_trackDzBS(nullptr),
385  t_trackDxyPV(nullptr),
386  t_trackDzPV(nullptr),
387  t_trackChiSq(nullptr),
388  t_trackPVIdx(nullptr),
389  t_NLayersCrossed(nullptr),
390  t_trackNOuterHits(nullptr),
391  t_trackHitsTOB(nullptr),
392  t_trackHitsTEC(nullptr),
393  t_trackHitInMissTOB(nullptr),
394  t_trackHitInMissTEC(nullptr),
395  t_trackHitInMissTIB(nullptr),
396  t_trackHitInMissTID(nullptr),
397  t_trackHitInMissTIBTID(nullptr),
398  t_trackHitOutMissTOB(nullptr),
399  t_trackHitOutMissTEC(nullptr),
400  t_trackHitOutMissTIB(nullptr),
401  t_trackHitOutMissTID(nullptr),
402  t_trackHitOutMissTOBTEC(nullptr),
403  t_trackHitInMeasTOB(nullptr),
404  t_trackHitInMeasTEC(nullptr),
405  t_trackHitInMeasTIB(nullptr),
406  t_trackHitInMeasTID(nullptr),
407  t_trackHitOutMeasTOB(nullptr),
408  t_trackHitOutMeasTEC(nullptr),
409  t_trackHitOutMeasTIB(nullptr),
410  t_trackHitOutMeasTID(nullptr),
411  t_trackOutPosOutHitDr(nullptr),
412  t_trackL(nullptr),
413  t_maxNearP31x31(nullptr),
414  t_maxNearP21x21(nullptr),
415  t_ecalSpike11x11(nullptr),
416  t_e7x7(nullptr),
417  t_e9x9(nullptr),
418  t_e11x11(nullptr),
419  t_e15x15(nullptr),
420  t_e7x7_10Sig(nullptr),
421  t_e9x9_10Sig(nullptr),
422  t_e11x11_10Sig(nullptr),
423  t_e15x15_10Sig(nullptr),
424  t_e7x7_15Sig(nullptr),
425  t_e9x9_15Sig(nullptr),
426  t_e11x11_15Sig(nullptr),
427  t_e15x15_15Sig(nullptr),
428  t_e7x7_20Sig(nullptr),
429  t_e9x9_20Sig(nullptr),
430  t_e11x11_20Sig(nullptr),
431  t_e15x15_20Sig(nullptr),
432  t_e7x7_25Sig(nullptr),
433  t_e9x9_25Sig(nullptr),
434  t_e11x11_25Sig(nullptr),
435  t_e15x15_25Sig(nullptr),
436  t_e7x7_30Sig(nullptr),
437  t_e9x9_30Sig(nullptr),
438  t_e11x11_30Sig(nullptr),
439  t_e15x15_30Sig(nullptr),
440  t_esimPdgId(nullptr),
441  t_simTrackP(nullptr),
442  t_trkEcalEne(nullptr),
443  t_esim7x7(nullptr),
444  t_esim9x9(nullptr),
445  t_esim11x11(nullptr),
446  t_esim15x15(nullptr),
447  t_esim7x7Matched(nullptr),
448  t_esim9x9Matched(nullptr),
449  t_esim11x11Matched(nullptr),
450  t_esim15x15Matched(nullptr),
451  t_esim7x7Rest(nullptr),
452  t_esim9x9Rest(nullptr),
453  t_esim11x11Rest(nullptr),
454  t_esim15x15Rest(nullptr),
455  t_esim7x7Photon(nullptr),
456  t_esim9x9Photon(nullptr),
457  t_esim11x11Photon(nullptr),
458  t_esim15x15Photon(nullptr),
459  t_esim7x7NeutHad(nullptr),
460  t_esim9x9NeutHad(nullptr),
461  t_esim11x11NeutHad(nullptr),
462  t_esim15x15NeutHad(nullptr),
463  t_esim7x7CharHad(nullptr),
464  t_esim9x9CharHad(nullptr),
465  t_esim11x11CharHad(nullptr),
466  t_esim15x15CharHad(nullptr),
467  t_maxNearHcalP3x3(nullptr),
468  t_maxNearHcalP5x5(nullptr),
469  t_maxNearHcalP7x7(nullptr),
470  t_h3x3(nullptr),
471  t_h5x5(nullptr),
472  t_h7x7(nullptr),
473  t_h3x3Sig(nullptr),
474  t_h5x5Sig(nullptr),
475  t_h7x7Sig(nullptr),
476  t_infoHcal(nullptr),
477  t_trkHcalEne(nullptr),
478  t_hsim3x3(nullptr),
479  t_hsim5x5(nullptr),
480  t_hsim7x7(nullptr),
481  t_hsim3x3Matched(nullptr),
482  t_hsim5x5Matched(nullptr),
483  t_hsim7x7Matched(nullptr),
484  t_hsim3x3Rest(nullptr),
485  t_hsim5x5Rest(nullptr),
486  t_hsim7x7Rest(nullptr),
487  t_hsim3x3Photon(nullptr),
488  t_hsim5x5Photon(nullptr),
489  t_hsim7x7Photon(nullptr),
490  t_hsim3x3NeutHad(nullptr),
491  t_hsim5x5NeutHad(nullptr),
492  t_hsim7x7NeutHad(nullptr),
493  t_hsim3x3CharHad(nullptr),
494  t_hsim5x5CharHad(nullptr),
495  t_hsim7x7CharHad(nullptr) {
496  if (L1TriggerAlgoInfo_) {
497  m_l1GtUtils = std::make_unique<L1GtUtils>(
499  }
500 
501  usesResource(TFileService::kSharedResource);
502 
503  //now do what ever initialization is needed
504 
505  edm::InputTag L1extraTauJetSource_ = iConfig.getParameter<edm::InputTag>("l1extraTauJetSource");
506  edm::InputTag L1extraCenJetSource_ = iConfig.getParameter<edm::InputTag>("l1extraCenJetSource");
507  edm::InputTag L1extraFwdJetSource_ = iConfig.getParameter<edm::InputTag>("l1extraFwdJetSource");
508  edm::InputTag L1extraMuonSource_ = iConfig.getParameter<edm::InputTag>("l1extraMuonSource");
509  edm::InputTag L1extraIsoEmSource_ = iConfig.getParameter<edm::InputTag>("l1extraIsoEmSource");
510  edm::InputTag L1extraNonIsoEmSource_ = iConfig.getParameter<edm::InputTag>("l1extraNonIsoEmSource");
511  edm::InputTag L1GTReadoutRcdSource_ = iConfig.getParameter<edm::InputTag>("l1GTReadoutRcdSource");
512  edm::InputTag L1GTObjectMapRcdSource_ = iConfig.getParameter<edm::InputTag>("l1GTObjectMapRcdSource");
513  edm::InputTag JetSrc_ = iConfig.getParameter<edm::InputTag>("jetSource");
514  edm::InputTag JetExtender_ = iConfig.getParameter<edm::InputTag>("jetExtender");
515  edm::InputTag HBHERecHitSource_ = iConfig.getParameter<edm::InputTag>("hbheRecHitSource");
516 
517  // define tokens for access
518  tok_L1extTauJet_ = consumes<l1extra::L1JetParticleCollection>(L1extraTauJetSource_);
519  tok_L1extCenJet_ = consumes<l1extra::L1JetParticleCollection>(L1extraCenJetSource_);
520  tok_L1extFwdJet_ = consumes<l1extra::L1JetParticleCollection>(L1extraFwdJetSource_);
521  tok_L1extMu_ = consumes<l1extra::L1MuonParticleCollection>(L1extraMuonSource_);
522  tok_L1extIsoEm_ = consumes<l1extra::L1EmParticleCollection>(L1extraIsoEmSource_);
523  tok_L1extNoIsoEm_ = consumes<l1extra::L1EmParticleCollection>(L1extraNonIsoEmSource_);
524  tok_jets_ = consumes<reco::CaloJetCollection>(JetSrc_);
525  tok_hbhe_ = consumes<HBHERecHitCollection>(HBHERecHitSource_);
526 
527  tok_genTrack_ = consumes<reco::TrackCollection>(edm::InputTag("generalTracks"));
528  tok_recVtx_ = consumes<reco::VertexCollection>(edm::InputTag("offlinePrimaryVertices"));
529  tok_bs_ = consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"));
530  tok_EB_ = consumes<EcalRecHitCollection>(edm::InputTag("ecalRecHit", "EcalRecHitsEB"));
531  tok_EE_ = consumes<EcalRecHitCollection>(edm::InputTag("ecalRecHit", "EcalRecHitsEE"));
532 
533  tok_simTk_ = consumes<edm::SimTrackContainer>(edm::InputTag("g4SimHits"));
534  tok_simVtx_ = consumes<edm::SimVertexContainer>(edm::InputTag("g4SimHits"));
535  tok_caloEB_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", "EcalHitsEB"));
536  tok_caloEE_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", "EcalHitsEE"));
537  tok_caloHH_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", "HcalHits"));
538 
539  if (myverbose_ >= 0) {
540  edm::LogVerbatim("IsoTrack") << "Parameters read from config file \n"
541  << " doMC " << doMC_ << "\t myverbose " << myverbose_ << "\t minTrackP "
542  << minTrackP_ << "\t maxTrackEta " << maxTrackEta_ << "\t tMinE " << tMinE_
543  << "\t tMaxE " << tMaxE_ << "\t tMinH " << tMinH_ << "\t tMaxH " << tMaxH_
544  << "\n debugL1Info " << debugL1Info_ << "\t L1TriggerAlgoInfo " << L1TriggerAlgoInfo_
545  << "\n";
546  }
547 
548  tok_geom_ = esConsumes<CaloGeometry, CaloGeometryRecord>();
549  tok_caloTopology_ = esConsumes<CaloTopology, CaloTopologyRecord>();
550  tok_topo_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
551  tok_magField_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
552  tok_ecalChStatus_ = esConsumes<EcalChannelStatus, EcalChannelStatusRcd>();
553  tok_sevlv_ = esConsumes<EcalSeverityLevelAlgo, EcalSeverityLevelAlgoRcd>();
554  tok_htmap_ = esConsumes<EcalTrigTowerConstituentsMap, IdealGeometryRecord>();
555 }
556 
558  delete t_PVx;
559  delete t_PVy;
560  delete t_PVz;
561  delete t_PVisValid;
562  delete t_PVndof;
563  delete t_PVNTracks;
564  delete t_PVNTracksWt;
565  delete t_PVTracksSumPt;
566  delete t_PVTracksSumPtWt;
567  delete t_PVNTracksHP;
568  delete t_PVNTracksHPWt;
569  delete t_PVTracksSumPtHP;
570  delete t_PVTracksSumPtHPWt;
571  delete t_L1AlgoNames;
572  delete t_L1PreScale;
573  delete t_L1CenJetPt;
574  delete t_L1CenJetEta;
575  delete t_L1CenJetPhi;
576  delete t_L1FwdJetPt;
577  delete t_L1FwdJetEta;
578  delete t_L1FwdJetPhi;
579  delete t_L1TauJetPt;
580  delete t_L1TauJetEta;
581  delete t_L1TauJetPhi;
582  delete t_L1MuonPt;
583  delete t_L1MuonEta;
584  delete t_L1MuonPhi;
585  delete t_L1IsoEMPt;
586  delete t_L1IsoEMEta;
587  delete t_L1IsoEMPhi;
588  delete t_L1NonIsoEMPt;
589  delete t_L1NonIsoEMEta;
590  delete t_L1NonIsoEMPhi;
591  delete t_L1METPt;
592  delete t_L1METEta;
593  delete t_L1METPhi;
594  delete t_jetPt;
595  delete t_jetEta;
596  delete t_jetPhi;
597  delete t_nTrksJetCalo;
598  delete t_nTrksJetVtx;
599  delete t_trackPAll;
600  delete t_trackEtaAll;
601  delete t_trackPhiAll;
602  delete t_trackPdgIdAll;
603  delete t_trackPtAll;
604  delete t_trackDxyAll;
605  delete t_trackDzAll;
606  delete t_trackDxyPVAll;
607  delete t_trackDzPVAll;
608  delete t_trackChiSqAll;
609  delete t_trackP;
610  delete t_trackPt;
611  delete t_trackEta;
612  delete t_trackPhi;
613  delete t_trackEcalEta;
614  delete t_trackEcalPhi;
615  delete t_trackHcalEta;
616  delete t_trackHcalPhi;
617  delete t_trackNOuterHits;
618  delete t_NLayersCrossed;
619  delete t_trackDxy;
620  delete t_trackDxyBS;
621  delete t_trackDz;
622  delete t_trackDzBS;
623  delete t_trackDxyPV;
624  delete t_trackDzPV;
625  delete t_trackPVIdx;
626  delete t_trackChiSq;
627  delete t_trackHitsTOB;
628  delete t_trackHitsTEC;
629  delete t_trackHitInMissTOB;
630  delete t_trackHitInMissTEC;
631  delete t_trackHitInMissTIB;
632  delete t_trackHitInMissTID;
633  delete t_trackHitInMissTIBTID;
634  delete t_trackHitOutMissTOB;
635  delete t_trackHitOutMissTEC;
636  delete t_trackHitOutMissTIB;
637  delete t_trackHitOutMissTID;
639  delete t_trackHitInMeasTOB;
640  delete t_trackHitInMeasTEC;
641  delete t_trackHitInMeasTIB;
642  delete t_trackHitInMeasTID;
643  delete t_trackHitOutMeasTOB;
644  delete t_trackHitOutMeasTEC;
645  delete t_trackHitOutMeasTIB;
646  delete t_trackHitOutMeasTID;
647  delete t_trackOutPosOutHitDr;
648  delete t_trackL;
649  delete t_maxNearP31x31;
650  delete t_maxNearP21x21;
651  delete t_ecalSpike11x11;
652  delete t_e7x7;
653  delete t_e9x9;
654  delete t_e11x11;
655  delete t_e15x15;
656  delete t_e7x7_10Sig;
657  delete t_e9x9_10Sig;
658  delete t_e11x11_10Sig;
659  delete t_e15x15_10Sig;
660  delete t_e7x7_15Sig;
661  delete t_e9x9_15Sig;
662  delete t_e11x11_15Sig;
663  delete t_e15x15_15Sig;
664  delete t_e7x7_20Sig;
665  delete t_e9x9_20Sig;
666  delete t_e11x11_20Sig;
667  delete t_e15x15_20Sig;
668  delete t_e7x7_25Sig;
669  delete t_e9x9_25Sig;
670  delete t_e11x11_25Sig;
671  delete t_e15x15_25Sig;
672  delete t_e7x7_30Sig;
673  delete t_e9x9_30Sig;
674  delete t_e11x11_30Sig;
675  delete t_e15x15_30Sig;
676  delete t_esim7x7;
677  delete t_esim9x9;
678  delete t_esim11x11;
679  delete t_esim15x15;
680  delete t_esim7x7Matched;
681  delete t_esim9x9Matched;
682  delete t_esim11x11Matched;
683  delete t_esim15x15Matched;
684  delete t_esim7x7Rest;
685  delete t_esim9x9Rest;
686  delete t_esim11x11Rest;
687  delete t_esim15x15Rest;
688  delete t_esim7x7Photon;
689  delete t_esim9x9Photon;
690  delete t_esim11x11Photon;
691  delete t_esim15x15Photon;
692  delete t_esim7x7NeutHad;
693  delete t_esim9x9NeutHad;
694  delete t_esim11x11NeutHad;
695  delete t_esim15x15NeutHad;
696  delete t_esim7x7CharHad;
697  delete t_esim9x9CharHad;
698  delete t_esim11x11CharHad;
699  delete t_esim15x15CharHad;
700  delete t_trkEcalEne;
701  delete t_simTrackP;
702  delete t_esimPdgId;
703  delete t_maxNearHcalP3x3;
704  delete t_maxNearHcalP5x5;
705  delete t_maxNearHcalP7x7;
706  delete t_h3x3;
707  delete t_h5x5;
708  delete t_h7x7;
709  delete t_h3x3Sig;
710  delete t_h5x5Sig;
711  delete t_h7x7Sig;
712  delete t_infoHcal;
713  delete t_trkHcalEne;
714  delete t_hsim3x3;
715  delete t_hsim5x5;
716  delete t_hsim7x7;
717  delete t_hsim3x3Matched;
718  delete t_hsim5x5Matched;
719  delete t_hsim7x7Matched;
720  delete t_hsim3x3Rest;
721  delete t_hsim5x5Rest;
722  delete t_hsim7x7Rest;
723  delete t_hsim3x3Photon;
724  delete t_hsim5x5Photon;
725  delete t_hsim7x7Photon;
726  delete t_hsim3x3NeutHad;
727  delete t_hsim5x5NeutHad;
728  delete t_hsim7x7NeutHad;
729  delete t_hsim3x3CharHad;
730  delete t_hsim5x5CharHad;
731  delete t_hsim7x7CharHad;
732 }
733 
736  desc.addUntracked<bool>("doMC", false);
737  desc.addUntracked<bool>("writeAllTracks", false);
738  desc.addUntracked<int>("verbosity", 1);
739  desc.addUntracked<double>("pvTracksPtMin", 0.200);
740  desc.addUntracked<int>("debugTracks", 0);
741  desc.addUntracked<bool>("printTrkHitPattern", true);
742  desc.addUntracked<double>("minTrackP", 1.0);
743  desc.addUntracked<double>("maxTrackEta", 2.6);
744  desc.addUntracked<bool>("debugL1Info", false);
745  desc.addUntracked<bool>("l1TriggerAlgoInfo", false);
746  desc.add<edm::InputTag>("l1extraTauJetSource", edm::InputTag("l1extraParticles", "Tau"));
747  desc.add<edm::InputTag>("l1extraCenJetSource", edm::InputTag("l1extraParticles", "Central"));
748  desc.add<edm::InputTag>("l1extraFwdJetSource", edm::InputTag("l1extraParticles", "Forward"));
749  desc.add<edm::InputTag>("l1extraMuonSource", edm::InputTag("l1extraParticles"));
750  desc.add<edm::InputTag>("l1extraIsoEmSource", edm::InputTag("l1extraParticles", "Isolated"));
751  desc.add<edm::InputTag>("l1extraNonIsoEmSource", edm::InputTag("l1extraParticles", "NonIsolated"));
752  desc.add<edm::InputTag>("l1GTReadoutRcdSource", edm::InputTag("gtDigis"));
753  desc.add<edm::InputTag>("l1GTObjectMapRcdSource", edm::InputTag("hltL1GtObjectMap"));
754  desc.add<edm::InputTag>("jetSource", edm::InputTag("iterativeCone5CaloJets"));
755  desc.add<edm::InputTag>("jetExtender", edm::InputTag("iterativeCone5JetExtender"));
756  desc.add<edm::InputTag>("hbheRecHitSource", edm::InputTag("hbhereco"));
757  desc.addUntracked<double>("maxNearTrackPT", 1.0);
758  desc.addUntracked<double>("timeMinCutECAL", -500.0);
759  desc.addUntracked<double>("timeMaxCutECAL", 500.0);
760  desc.addUntracked<double>("timeMinCutHCAL", -500.0);
761  desc.addUntracked<double>("timeMaxCutHCAL", 500.0);
762  descriptions.add("isolatedTracksNxN", desc);
763 }
764 
766  bool haveIsoTrack = false;
767 
768  const MagneticField *bField = &iSetup.getData(tok_magField_);
769 
771 
772  t_RunNo = iEvent.id().run();
773  t_EvtNo = iEvent.id().event();
774  t_Lumi = iEvent.luminosityBlock();
775  t_Bunch = iEvent.bunchCrossing();
776 
777  ++nEventProc_;
778 
780  iEvent.getByToken(tok_genTrack_, trkCollection);
781  if (debugTrks_ > 1) {
782  edm::LogVerbatim("IsoTrack") << "Track Collection: ";
783  edm::LogVerbatim("IsoTrack") << "Number of Tracks " << trkCollection->size();
784  }
785  std::string theTrackQuality = "highPurity";
786  reco::TrackBase::TrackQuality trackQuality_ = reco::TrackBase::qualityByName(theTrackQuality);
787 
788  //===================== save L1 Trigger information =======================
789  if (L1TriggerAlgoInfo_) {
790  m_l1GtUtils->getL1GtRunCache(iEvent, iSetup, useL1EventSetup, useL1GtTriggerMenuLite);
791 
792  int iErrorCode = -1;
793  int l1ConfCode = -1;
794  const bool l1Conf = m_l1GtUtils->availableL1Configuration(iErrorCode, l1ConfCode);
795  if (!l1Conf) {
796  edm::LogVerbatim("IsoTrack")
797  << "\nL1 configuration code:" << l1ConfCode << "\nNo valid L1 trigger configuration available."
798  << "\nSee text above for error code interpretation"
799  << "\nNo return here, in order to test each method, protected against configuration error.";
800  }
801 
802  const L1GtTriggerMenu *m_l1GtMenu = m_l1GtUtils->ptrL1TriggerMenuEventSetup(iErrorCode);
803  const AlgorithmMap &algorithmMap = m_l1GtMenu->gtAlgorithmMap();
804  const std::string &menuName = m_l1GtMenu->gtTriggerMenuName();
805 
806  if (!initL1_) {
807  initL1_ = true;
808  edm::LogVerbatim("IsoTrack") << "menuName " << menuName;
809  for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
810  std::string algName = itAlgo->first;
811  int algBitNumber = (itAlgo->second).algoBitNumber();
812  l1AlgoMap_.insert(std::pair<std::pair<unsigned int, std::string>, int>(
813  std::pair<unsigned int, std::string>(algBitNumber, algName), 0));
814  }
815  std::map<std::pair<unsigned int, std::string>, int>::iterator itr;
816  for (itr = l1AlgoMap_.begin(); itr != l1AlgoMap_.end(); itr++) {
817  edm::LogVerbatim("IsoTrack") << " ********** " << (itr->first).first << " " << (itr->first).second << " "
818  << itr->second;
819  }
820  }
821 
822  std::vector<int> algbits;
823  for (CItAlgo itAlgo = algorithmMap.begin(); itAlgo != algorithmMap.end(); itAlgo++) {
824  std::string algName = itAlgo->first;
825  int algBitNumber = (itAlgo->second).algoBitNumber();
826  bool decision = m_l1GtUtils->decision(iEvent, itAlgo->first, iErrorCode);
827  int preScale = m_l1GtUtils->prescaleFactor(iEvent, itAlgo->first, iErrorCode);
828 
829  // save the algo names which fired
830  if (decision) {
831  l1AlgoMap_[std::pair<unsigned int, std::string>(algBitNumber, algName)] += 1;
832  h_L1AlgoNames->Fill(algBitNumber);
833  t_L1AlgoNames->push_back(itAlgo->first);
834  t_L1PreScale->push_back(preScale);
835  t_L1Decision[algBitNumber] = 1;
836  algbits.push_back(algBitNumber);
837  }
838  }
839 
840  if (debugL1Info_) {
841  for (unsigned int ii = 0; ii < t_L1AlgoNames->size(); ii++) {
842  edm::LogVerbatim("IsoTrack") << ii << " " << (*t_L1AlgoNames)[ii] << " " << (*t_L1PreScale)[ii] << " "
843  << algbits[ii];
844  }
845  for (int i = 0; i < 128; ++i) {
846  edm::LogVerbatim("IsoTrack") << "L1Decision: " << i << ":" << t_L1Decision[i];
847  }
848  }
849 
850  // L1Taus
852  iEvent.getByToken(tok_L1extTauJet_, l1TauHandle);
853  l1extra::L1JetParticleCollection::const_iterator itr;
854  int iL1Obj = 0;
855  for (itr = l1TauHandle->begin(), iL1Obj = 0; itr != l1TauHandle->end(); ++itr, iL1Obj++) {
856  if (iL1Obj < 1) {
857  t_L1TauJetPt->push_back(itr->pt());
858  t_L1TauJetEta->push_back(itr->eta());
859  t_L1TauJetPhi->push_back(itr->phi());
860  }
861  if (debugL1Info_) {
862  edm::LogVerbatim("IsoTrack") << "tauJ p/pt " << itr->momentum() << " " << itr->pt() << " eta/phi "
863  << itr->eta() << " " << itr->phi();
864  }
865  }
866 
867  // L1 Central Jets
869  iEvent.getByToken(tok_L1extCenJet_, l1CenJetHandle);
870  for (itr = l1CenJetHandle->begin(), iL1Obj = 0; itr != l1CenJetHandle->end(); ++itr, iL1Obj++) {
871  if (iL1Obj < 1) {
872  t_L1CenJetPt->push_back(itr->pt());
873  t_L1CenJetEta->push_back(itr->eta());
874  t_L1CenJetPhi->push_back(itr->phi());
875  }
876  if (debugL1Info_) {
877  edm::LogVerbatim("IsoTrack") << "cenJ p/pt " << itr->momentum() << " " << itr->pt() << " eta/phi "
878  << itr->eta() << " " << itr->phi();
879  }
880  }
881 
882  // L1 Forward Jets
884  iEvent.getByToken(tok_L1extFwdJet_, l1FwdJetHandle);
885  for (itr = l1FwdJetHandle->begin(), iL1Obj = 0; itr != l1FwdJetHandle->end(); ++itr, iL1Obj++) {
886  if (iL1Obj < 1) {
887  t_L1FwdJetPt->push_back(itr->pt());
888  t_L1FwdJetEta->push_back(itr->eta());
889  t_L1FwdJetPhi->push_back(itr->phi());
890  }
891  if (debugL1Info_) {
892  edm::LogVerbatim("IsoTrack") << "fwdJ p/pt " << itr->momentum() << " " << itr->pt() << " eta/phi "
893  << itr->eta() << " " << itr->phi();
894  }
895  }
896 
897  // L1 Isolated EM onjects
898  l1extra::L1EmParticleCollection::const_iterator itrEm;
900  iEvent.getByToken(tok_L1extIsoEm_, l1IsoEmHandle);
901  for (itrEm = l1IsoEmHandle->begin(), iL1Obj = 0; itrEm != l1IsoEmHandle->end(); ++itrEm, iL1Obj++) {
902  if (iL1Obj < 1) {
903  t_L1IsoEMPt->push_back(itrEm->pt());
904  t_L1IsoEMEta->push_back(itrEm->eta());
905  t_L1IsoEMPhi->push_back(itrEm->phi());
906  }
907  if (debugL1Info_) {
908  edm::LogVerbatim("IsoTrack") << "isoEm p/pt " << itrEm->momentum() << " " << itrEm->pt() << " eta/phi "
909  << itrEm->eta() << " " << itrEm->phi();
910  }
911  }
912 
913  // L1 Non-Isolated EM onjects
915  iEvent.getByToken(tok_L1extNoIsoEm_, l1NonIsoEmHandle);
916  for (itrEm = l1NonIsoEmHandle->begin(), iL1Obj = 0; itrEm != l1NonIsoEmHandle->end(); ++itrEm, iL1Obj++) {
917  if (iL1Obj < 1) {
918  t_L1NonIsoEMPt->push_back(itrEm->pt());
919  t_L1NonIsoEMEta->push_back(itrEm->eta());
920  t_L1NonIsoEMPhi->push_back(itrEm->phi());
921  }
922  if (debugL1Info_) {
923  edm::LogVerbatim("IsoTrack") << "nonIsoEm p/pt " << itrEm->momentum() << " " << itrEm->pt() << " eta/phi "
924  << itrEm->eta() << " " << itrEm->phi();
925  }
926  }
927 
928  // L1 Muons
929  l1extra::L1MuonParticleCollection::const_iterator itrMu;
931  iEvent.getByToken(tok_L1extMu_, l1MuHandle);
932  for (itrMu = l1MuHandle->begin(), iL1Obj = 0; itrMu != l1MuHandle->end(); ++itrMu, iL1Obj++) {
933  if (iL1Obj < 1) {
934  t_L1MuonPt->push_back(itrMu->pt());
935  t_L1MuonEta->push_back(itrMu->eta());
936  t_L1MuonPhi->push_back(itrMu->phi());
937  }
938  if (debugL1Info_) {
939  edm::LogVerbatim("IsoTrack") << "l1muon p/pt " << itrMu->momentum() << " " << itrMu->pt() << " eta/phi "
940  << itrMu->eta() << " " << itrMu->phi();
941  }
942  }
943  }
944 
945  //============== store the information about all the Non-Fake vertices ===============
946 
948  iEvent.getByToken(tok_recVtx_, recVtxs);
949 
950  std::vector<reco::Track> svTracks;
951  math::XYZPoint leadPV(0, 0, 0);
952  double sumPtMax = -1.0;
953  for (unsigned int ind = 0; ind < recVtxs->size(); ind++) {
954  if (!((*recVtxs)[ind].isFake())) {
955  double vtxTrkSumPt = 0.0, vtxTrkSumPtWt = 0.0;
956  int vtxTrkNWt = 0;
957  double vtxTrkSumPtHP = 0.0, vtxTrkSumPtHPWt = 0.0;
958  int vtxTrkNHP = 0, vtxTrkNHPWt = 0;
959 
960  reco::Vertex::trackRef_iterator vtxTrack = (*recVtxs)[ind].tracks_begin();
961 
962  for (vtxTrack = (*recVtxs)[ind].tracks_begin(); vtxTrack != (*recVtxs)[ind].tracks_end(); vtxTrack++) {
963  if ((*vtxTrack)->pt() < pvTracksPtMin_)
964  continue;
965 
966  bool trkQuality = (*vtxTrack)->quality(trackQuality_);
967 
968  vtxTrkSumPt += (*vtxTrack)->pt();
969 
970  vtxTrkSumPt += (*vtxTrack)->pt();
971  if (trkQuality) {
972  vtxTrkSumPtHP += (*vtxTrack)->pt();
973  vtxTrkNHP++;
974  }
975 
976  double weight = (*recVtxs)[ind].trackWeight(*vtxTrack);
977  h_PVTracksWt->Fill(weight);
978  if (weight > 0.5) {
979  vtxTrkSumPtWt += (*vtxTrack)->pt();
980  vtxTrkNWt++;
981  if (trkQuality) {
982  vtxTrkSumPtHPWt += (*vtxTrack)->pt();
983  vtxTrkNHPWt++;
984  }
985  }
986  }
987 
988  if (vtxTrkSumPt > sumPtMax) {
989  sumPtMax = vtxTrkSumPt;
990  leadPV = math::XYZPoint((*recVtxs)[ind].x(), (*recVtxs)[ind].y(), (*recVtxs)[ind].z());
991  }
992 
993  t_PVx->push_back((*recVtxs)[ind].x());
994  t_PVy->push_back((*recVtxs)[ind].y());
995  t_PVz->push_back((*recVtxs)[ind].z());
996  t_PVisValid->push_back((*recVtxs)[ind].isValid());
997  t_PVNTracks->push_back((*recVtxs)[ind].tracksSize());
998  t_PVndof->push_back((*recVtxs)[ind].ndof());
999  t_PVNTracksWt->push_back(vtxTrkNWt);
1000  t_PVTracksSumPt->push_back(vtxTrkSumPt);
1001  t_PVTracksSumPtWt->push_back(vtxTrkSumPtWt);
1002 
1003  t_PVNTracksHP->push_back(vtxTrkNHP);
1004  t_PVNTracksHPWt->push_back(vtxTrkNHPWt);
1005  t_PVTracksSumPtHP->push_back(vtxTrkSumPtHP);
1006  t_PVTracksSumPtHPWt->push_back(vtxTrkSumPtHPWt);
1007 
1008  if (myverbose_ == 4) {
1009  edm::LogVerbatim("IsoTrack") << "PV " << ind << " isValid " << (*recVtxs)[ind].isValid() << " isFake "
1010  << (*recVtxs)[ind].isFake() << " hasRefittedTracks() " << ind << " "
1011  << (*recVtxs)[ind].hasRefittedTracks() << " refittedTrksSize "
1012  << (*recVtxs)[ind].refittedTracks().size() << " tracksSize() "
1013  << (*recVtxs)[ind].tracksSize() << " sumPt " << vtxTrkSumPt;
1014  }
1015  } // if vtx is not Fake
1016  } // loop over PVs
1017  //===================================================================================
1018 
1019  // Get the beamspot
1020  edm::Handle<reco::BeamSpot> beamSpotH;
1021  iEvent.getByToken(tok_bs_, beamSpotH);
1022  math::XYZPoint bspot;
1023  bspot = (beamSpotH.isValid()) ? beamSpotH->position() : math::XYZPoint(0, 0, 0);
1024 
1025  //=====================================================================
1026 
1028  iEvent.getByToken(tok_jets_, jets);
1029  // edm::Handle<reco::JetExtendedAssociation::Container> jetExtender;
1030  // iEvent.getByLabel(JetExtender_,jetExtender);
1031 
1032  for (unsigned int ijet = 0; ijet < (*jets).size(); ijet++) {
1033  t_jetPt->push_back((*jets)[ijet].pt());
1034  t_jetEta->push_back((*jets)[ijet].eta());
1035  t_jetPhi->push_back((*jets)[ijet].phi());
1036  t_nTrksJetVtx->push_back(-1.0);
1037  t_nTrksJetCalo->push_back(-1.0);
1038  }
1039 
1040  //===================================================================================
1041 
1042  // get handles to calogeometry and calotopology
1043  const CaloGeometry *geo = &iSetup.getData(tok_geom_);
1044  const CaloTopology *caloTopology = &iSetup.getData(tok_caloTopology_);
1045  const HcalTopology *theHBHETopology = &iSetup.getData(tok_topo_);
1046 
1047  edm::Handle<EcalRecHitCollection> barrelRecHitsHandle;
1048  edm::Handle<EcalRecHitCollection> endcapRecHitsHandle;
1049  iEvent.getByToken(tok_EB_, barrelRecHitsHandle);
1050  iEvent.getByToken(tok_EE_, endcapRecHitsHandle);
1051 
1052  // Retrieve the good/bad ECAL channels from the DB
1053  const EcalChannelStatus *theEcalChStatus = &iSetup.getData(tok_ecalChStatus_);
1054  const EcalSeverityLevelAlgo *sevlv = &iSetup.getData(tok_sevlv_);
1055 
1056  // Retrieve trigger tower map
1057  const EcalTrigTowerConstituentsMap &ttMap = iSetup.getData(tok_htmap_);
1058 
1060  iEvent.getByToken(tok_hbhe_, hbhe);
1061  if (!hbhe.isValid()) {
1062  ++nbad_;
1063  if (nbad_ < 10)
1064  edm::LogVerbatim("IsoTrack") << "No HBHE rechit collection";
1065  return;
1066  }
1067  const HBHERecHitCollection Hithbhe = *(hbhe.product());
1068 
1069  //get Handles to SimTracks and SimHits
1071  if (doMC_)
1072  iEvent.getByToken(tok_simTk_, SimTk);
1073 
1075  if (doMC_)
1076  iEvent.getByToken(tok_simVtx_, SimVtx);
1077 
1078  //get Handles to PCaloHitContainers of eb/ee/hbhe
1080  if (doMC_)
1081  iEvent.getByToken(tok_caloEB_, pcaloeb);
1082 
1084  if (doMC_)
1085  iEvent.getByToken(tok_caloEE_, pcaloee);
1086 
1088  if (doMC_)
1089  iEvent.getByToken(tok_caloHH_, pcalohh);
1090 
1091  //associates tracker rechits/simhits to a track
1092  std::unique_ptr<TrackerHitAssociator> associate;
1093  if (doMC_)
1094  associate = std::make_unique<TrackerHitAssociator>(iEvent, trackerHitAssociatorConfig_);
1095 
1096  //===================================================================================
1097 
1098  h_nTracks->Fill(trkCollection->size());
1099 
1100  int nTracks = 0;
1101 
1102  t_nTracks = trkCollection->size();
1103 
1104  // get the list of DetIds closest to the impact point of track on surface calorimeters
1105  std::vector<spr::propagatedTrackID> trkCaloDets;
1106  spr::propagateCALO(trkCollection, geo, bField, theTrackQuality, trkCaloDets, false);
1107  std::vector<spr::propagatedTrackID>::const_iterator trkDetItr;
1108 
1109  if (myverbose_ > 2) {
1110  for (trkDetItr = trkCaloDets.begin(); trkDetItr != trkCaloDets.end(); trkDetItr++) {
1111  edm::LogVerbatim("IsoTrack") << trkDetItr->trkItr->p() << " " << trkDetItr->trkItr->eta() << " "
1112  << trkDetItr->okECAL << " " << trkDetItr->okHCAL;
1113  if (trkDetItr->okECAL) {
1114  if (trkDetItr->detIdECAL.subdetId() == EcalBarrel)
1115  edm::LogVerbatim("IsoTrack") << (EBDetId)trkDetItr->detIdECAL;
1116  else
1117  edm::LogVerbatim("IsoTrack") << (EEDetId)trkDetItr->detIdECAL;
1118  }
1119  if (trkDetItr->okHCAL)
1120  edm::LogVerbatim("IsoTrack") << (HcalDetId)trkDetItr->detIdHCAL;
1121  }
1122  }
1123 
1124  int nvtxTracks = 0;
1125  for (trkDetItr = trkCaloDets.begin(), nTracks = 0; trkDetItr != trkCaloDets.end(); trkDetItr++, nTracks++) {
1126  const reco::Track *pTrack = &(*(trkDetItr->trkItr));
1127 
1128  // find vertex index the track is associated with
1129  int pVtxTkId = -1;
1130  for (unsigned int ind = 0; ind < recVtxs->size(); ind++) {
1131  if (!((*recVtxs)[ind].isFake())) {
1132  reco::Vertex::trackRef_iterator vtxTrack = (*recVtxs)[ind].tracks_begin();
1133  for (vtxTrack = (*recVtxs)[ind].tracks_begin(); vtxTrack != (*recVtxs)[ind].tracks_end(); vtxTrack++) {
1134  const edm::RefToBase<reco::Track> pvtxTrack = (*vtxTrack);
1135  if (pTrack == pvtxTrack.get()) {
1136  pVtxTkId = ind;
1137  break;
1138  if (myverbose_ > 2) {
1139  if (pTrack->pt() > 1.0) {
1140  edm::LogVerbatim("IsoTrack") << "Debug the track association with vertex ";
1141  edm::LogVerbatim("IsoTrack") << pTrack << " " << pvtxTrack.get();
1142  edm::LogVerbatim("IsoTrack")
1143  << " trkVtxIndex " << nvtxTracks << " vtx " << ind << " pt " << pTrack->pt() << " eta "
1144  << pTrack->eta() << " " << pTrack->pt() - pvtxTrack->pt() << " "
1145  << pTrack->eta() - pvtxTrack->eta();
1146  nvtxTracks++;
1147  }
1148  }
1149  }
1150  }
1151  }
1152  }
1153 
1154  const reco::HitPattern &hitp = pTrack->hitPattern();
1155  int nLayersCrossed = hitp.trackerLayersWithMeasurement();
1156  int nOuterHits = hitp.stripTOBLayersWithMeasurement() + hitp.stripTECLayersWithMeasurement();
1157 
1158  bool ifGood = pTrack->quality(trackQuality_);
1159  double pt1 = pTrack->pt();
1160  double p1 = pTrack->p();
1161  double eta1 = pTrack->momentum().eta();
1162  double phi1 = pTrack->momentum().phi();
1163  double etaEcal1 = trkDetItr->etaECAL;
1164  double phiEcal1 = trkDetItr->phiECAL;
1165  double etaHcal1 = trkDetItr->etaHCAL;
1166  double phiHcal1 = trkDetItr->phiHCAL;
1167  double dxy1 = pTrack->dxy();
1168  double dz1 = pTrack->dz();
1169  double chisq1 = pTrack->normalizedChi2();
1170  double dxybs1 = beamSpotH.isValid() ? pTrack->dxy(bspot) : pTrack->dxy();
1171  double dzbs1 = beamSpotH.isValid() ? pTrack->dz(bspot) : pTrack->dz();
1172  double dxypv1 = pTrack->dxy();
1173  double dzpv1 = pTrack->dz();
1174  if (pVtxTkId >= 0) {
1175  math::XYZPoint thisTkPV =
1176  math::XYZPoint((*recVtxs)[pVtxTkId].x(), (*recVtxs)[pVtxTkId].y(), (*recVtxs)[pVtxTkId].z());
1177  dxypv1 = pTrack->dxy(thisTkPV);
1178  dzpv1 = pTrack->dz(thisTkPV);
1179  }
1180 
1181  h_recEtaPt_0->Fill(eta1, pt1);
1182  h_recEtaP_0->Fill(eta1, p1);
1183  h_recPt_0->Fill(pt1);
1184  h_recP_0->Fill(p1);
1185  h_recEta_0->Fill(eta1);
1186  h_recPhi_0->Fill(phi1);
1187 
1188  if (ifGood && nLayersCrossed > 7) {
1189  h_recEtaPt_1->Fill(eta1, pt1);
1190  h_recEtaP_1->Fill(eta1, p1);
1191  h_recPt_1->Fill(pt1);
1192  h_recP_1->Fill(p1);
1193  h_recEta_1->Fill(eta1);
1194  h_recPhi_1->Fill(phi1);
1195  }
1196 
1197  if (!ifGood)
1198  continue;
1199 
1200  if (writeAllTracks_ && p1 > 2.0 && nLayersCrossed > 7) {
1201  t_trackPAll->push_back(p1);
1202  t_trackEtaAll->push_back(eta1);
1203  t_trackPhiAll->push_back(phi1);
1204  t_trackPtAll->push_back(pt1);
1205  t_trackDxyAll->push_back(dxy1);
1206  t_trackDzAll->push_back(dz1);
1207  t_trackDxyPVAll->push_back(dxypv1);
1208  t_trackDzPVAll->push_back(dzpv1);
1209  t_trackChiSqAll->push_back(chisq1);
1210  }
1211  if (doMC_) {
1212  edm::SimTrackContainer::const_iterator matchedSimTrkAll =
1213  spr::matchedSimTrack(iEvent, SimTk, SimVtx, pTrack, *associate, false);
1214  if (writeAllTracks_ && matchedSimTrkAll != SimTk->end())
1215  t_trackPdgIdAll->push_back(matchedSimTrkAll->type());
1216  }
1217 
1218  if (pt1 > minTrackP_ && std::abs(eta1) < maxTrackEta_ && trkDetItr->okECAL) {
1219  double maxNearP31x31 = 999.0, maxNearP25x25 = 999.0, maxNearP21x21 = 999.0, maxNearP15x15 = 999.0;
1220  maxNearP31x31 = spr::chargeIsolationEcal(nTracks, trkCaloDets, geo, caloTopology, 15, 15);
1221  maxNearP25x25 = spr::chargeIsolationEcal(nTracks, trkCaloDets, geo, caloTopology, 12, 12);
1222  maxNearP21x21 = spr::chargeIsolationEcal(nTracks, trkCaloDets, geo, caloTopology, 10, 10);
1223  maxNearP15x15 = spr::chargeIsolationEcal(nTracks, trkCaloDets, geo, caloTopology, 7, 7);
1224 
1225  int iTrkEtaBin = -1, iTrkMomBin = -1;
1226  for (unsigned int ieta = 0; ieta < NEtaBins; ieta++) {
1227  if (std::abs(eta1) > genPartEtaBins[ieta] && std::abs(eta1) < genPartEtaBins[ieta + 1])
1228  iTrkEtaBin = ieta;
1229  }
1230  for (unsigned int ipt = 0; ipt < NPBins; ipt++) {
1231  if (p1 > genPartPBins[ipt] && p1 < genPartPBins[ipt + 1])
1232  iTrkMomBin = ipt;
1233  }
1234  if (iTrkMomBin >= 0 && iTrkEtaBin >= 0) {
1235  h_maxNearP31x31[iTrkMomBin][iTrkEtaBin]->Fill(maxNearP31x31);
1236  h_maxNearP25x25[iTrkMomBin][iTrkEtaBin]->Fill(maxNearP25x25);
1237  h_maxNearP21x21[iTrkMomBin][iTrkEtaBin]->Fill(maxNearP21x21);
1238  h_maxNearP15x15[iTrkMomBin][iTrkEtaBin]->Fill(maxNearP15x15);
1239  }
1240  if (maxNearP31x31 < 0.0 && nLayersCrossed > 7 && nOuterHits > 4) {
1241  h_recEtaPt_2->Fill(eta1, pt1);
1242  h_recEtaP_2->Fill(eta1, p1);
1243  h_recPt_2->Fill(pt1);
1244  h_recP_2->Fill(p1);
1245  h_recEta_2->Fill(eta1);
1246  h_recPhi_2->Fill(phi1);
1247  }
1248 
1249  // if charge isolated in ECAL, store the further quantities
1250  if (maxNearP31x31 < 1.0) {
1251  haveIsoTrack = true;
1252 
1253  // get the matching simTrack
1254  double simTrackP = -1;
1255  if (doMC_) {
1256  edm::SimTrackContainer::const_iterator matchedSimTrk =
1257  spr::matchedSimTrack(iEvent, SimTk, SimVtx, pTrack, *associate, false);
1258  if (matchedSimTrk != SimTk->end())
1259  simTrackP = matchedSimTrk->momentum().P();
1260  }
1261 
1262  // get ECal Tranverse Profile
1263  std::pair<double, bool> e7x7P, e9x9P, e11x11P, e15x15P;
1264  std::pair<double, bool> e7x7_10SigP, e9x9_10SigP, e11x11_10SigP, e15x15_10SigP;
1265  std::pair<double, bool> e7x7_15SigP, e9x9_15SigP, e11x11_15SigP, e15x15_15SigP;
1266  std::pair<double, bool> e7x7_20SigP, e9x9_20SigP, e11x11_20SigP, e15x15_20SigP;
1267  std::pair<double, bool> e7x7_25SigP, e9x9_25SigP, e11x11_25SigP, e15x15_25SigP;
1268  std::pair<double, bool> e7x7_30SigP, e9x9_30SigP, e11x11_30SigP, e15x15_30SigP;
1269 
1270  spr::caloSimInfo simInfo3x3, simInfo5x5, simInfo7x7, simInfo9x9;
1271  spr::caloSimInfo simInfo11x11, simInfo13x13, simInfo15x15, simInfo21x21, simInfo25x25, simInfo31x31;
1272  double trkEcalEne = 0;
1273 
1274  const DetId isoCell = trkDetItr->detIdECAL;
1275  e7x7P = spr::eECALmatrix(isoCell,
1276  barrelRecHitsHandle,
1277  endcapRecHitsHandle,
1278  *theEcalChStatus,
1279  geo,
1280  caloTopology,
1281  sevlv,
1282  3,
1283  3,
1284  -100.0,
1285  -100.0,
1286  tMinE_,
1287  tMaxE_);
1288  e9x9P = spr::eECALmatrix(isoCell,
1289  barrelRecHitsHandle,
1290  endcapRecHitsHandle,
1291  *theEcalChStatus,
1292  geo,
1293  caloTopology,
1294  sevlv,
1295  4,
1296  4,
1297  -100.0,
1298  -100.0,
1299  tMinE_,
1300  tMaxE_);
1301  e11x11P = spr::eECALmatrix(isoCell,
1302  barrelRecHitsHandle,
1303  endcapRecHitsHandle,
1304  *theEcalChStatus,
1305  geo,
1306  caloTopology,
1307  sevlv,
1308  5,
1309  5,
1310  -100.0,
1311  -100.0,
1312  tMinE_,
1313  tMaxE_);
1314  e15x15P = spr::eECALmatrix(isoCell,
1315  barrelRecHitsHandle,
1316  endcapRecHitsHandle,
1317  *theEcalChStatus,
1318  geo,
1319  caloTopology,
1320  sevlv,
1321  7,
1322  7,
1323  -100.0,
1324  -100.0,
1325  tMinE_,
1326  tMaxE_);
1327 
1328  e7x7_10SigP = spr::eECALmatrix(isoCell,
1329  barrelRecHitsHandle,
1330  endcapRecHitsHandle,
1331  *theEcalChStatus,
1332  geo,
1333  caloTopology,
1334  sevlv,
1335  3,
1336  3,
1337  0.030,
1338  0.150,
1339  tMinE_,
1340  tMaxE_);
1341  e9x9_10SigP = spr::eECALmatrix(isoCell,
1342  barrelRecHitsHandle,
1343  endcapRecHitsHandle,
1344  *theEcalChStatus,
1345  geo,
1346  caloTopology,
1347  sevlv,
1348  4,
1349  4,
1350  0.030,
1351  0.150,
1352  tMinE_,
1353  tMaxE_);
1354  e11x11_10SigP = spr::eECALmatrix(isoCell,
1355  barrelRecHitsHandle,
1356  endcapRecHitsHandle,
1357  *theEcalChStatus,
1358  geo,
1359  caloTopology,
1360  sevlv,
1361  5,
1362  5,
1363  0.030,
1364  0.150,
1365  tMinE_,
1366  tMaxE_);
1367  e15x15_10SigP = spr::eECALmatrix(isoCell,
1368  barrelRecHitsHandle,
1369  endcapRecHitsHandle,
1370  *theEcalChStatus,
1371  geo,
1372  caloTopology,
1373  sevlv,
1374  7,
1375  7,
1376  0.030,
1377  0.150,
1378  tMinE_,
1379  tMaxE_);
1380 
1381  e7x7_15SigP = spr::eECALmatrix(isoCell,
1382  barrelRecHitsHandle,
1383  endcapRecHitsHandle,
1384  *theEcalChStatus,
1385  geo,
1386  caloTopology,
1387  sevlv,
1388  ttMap,
1389  3,
1390  3,
1391  0.20,
1392  0.45,
1393  tMinE_,
1394  tMaxE_);
1395  e9x9_15SigP = spr::eECALmatrix(isoCell,
1396  barrelRecHitsHandle,
1397  endcapRecHitsHandle,
1398  *theEcalChStatus,
1399  geo,
1400  caloTopology,
1401  sevlv,
1402  ttMap,
1403  4,
1404  4,
1405  0.20,
1406  0.45,
1407  tMinE_,
1408  tMaxE_);
1409  e11x11_15SigP = spr::eECALmatrix(isoCell,
1410  barrelRecHitsHandle,
1411  endcapRecHitsHandle,
1412  *theEcalChStatus,
1413  geo,
1414  caloTopology,
1415  sevlv,
1416  ttMap,
1417  5,
1418  5,
1419  0.20,
1420  0.45,
1421  tMinE_,
1422  tMaxE_);
1423  e15x15_15SigP = spr::eECALmatrix(isoCell,
1424  barrelRecHitsHandle,
1425  endcapRecHitsHandle,
1426  *theEcalChStatus,
1427  geo,
1428  caloTopology,
1429  sevlv,
1430  ttMap,
1431  7,
1432  7,
1433  0.20,
1434  0.45,
1435  tMinE_,
1436  tMaxE_,
1437  false);
1438 
1439  e7x7_20SigP = spr::eECALmatrix(isoCell,
1440  barrelRecHitsHandle,
1441  endcapRecHitsHandle,
1442  *theEcalChStatus,
1443  geo,
1444  caloTopology,
1445  sevlv,
1446  3,
1447  3,
1448  0.060,
1449  0.300,
1450  tMinE_,
1451  tMaxE_);
1452  e9x9_20SigP = spr::eECALmatrix(isoCell,
1453  barrelRecHitsHandle,
1454  endcapRecHitsHandle,
1455  *theEcalChStatus,
1456  geo,
1457  caloTopology,
1458  sevlv,
1459  4,
1460  4,
1461  0.060,
1462  0.300,
1463  tMinE_,
1464  tMaxE_);
1465  e11x11_20SigP = spr::eECALmatrix(isoCell,
1466  barrelRecHitsHandle,
1467  endcapRecHitsHandle,
1468  *theEcalChStatus,
1469  geo,
1470  caloTopology,
1471  sevlv,
1472  5,
1473  5,
1474  0.060,
1475  0.300,
1476  tMinE_,
1477  tMaxE_);
1478  e15x15_20SigP = spr::eECALmatrix(isoCell,
1479  barrelRecHitsHandle,
1480  endcapRecHitsHandle,
1481  *theEcalChStatus,
1482  geo,
1483  caloTopology,
1484  sevlv,
1485  7,
1486  7,
1487  0.060,
1488  0.300,
1489  tMinE_,
1490  tMaxE_);
1491 
1492  e7x7_25SigP = spr::eECALmatrix(isoCell,
1493  barrelRecHitsHandle,
1494  endcapRecHitsHandle,
1495  *theEcalChStatus,
1496  geo,
1497  caloTopology,
1498  sevlv,
1499  3,
1500  3,
1501  0.075,
1502  0.375,
1503  tMinE_,
1504  tMaxE_);
1505  e9x9_25SigP = spr::eECALmatrix(isoCell,
1506  barrelRecHitsHandle,
1507  endcapRecHitsHandle,
1508  *theEcalChStatus,
1509  geo,
1510  caloTopology,
1511  sevlv,
1512  4,
1513  4,
1514  0.075,
1515  0.375,
1516  tMinE_,
1517  tMaxE_);
1518  e11x11_25SigP = spr::eECALmatrix(isoCell,
1519  barrelRecHitsHandle,
1520  endcapRecHitsHandle,
1521  *theEcalChStatus,
1522  geo,
1523  caloTopology,
1524  sevlv,
1525  5,
1526  5,
1527  0.075,
1528  0.375,
1529  tMinE_,
1530  tMaxE_);
1531  e15x15_25SigP = spr::eECALmatrix(isoCell,
1532  barrelRecHitsHandle,
1533  endcapRecHitsHandle,
1534  *theEcalChStatus,
1535  geo,
1536  caloTopology,
1537  sevlv,
1538  7,
1539  7,
1540  0.075,
1541  0.375,
1542  tMinE_,
1543  tMaxE_);
1544 
1545  e7x7_30SigP = spr::eECALmatrix(isoCell,
1546  barrelRecHitsHandle,
1547  endcapRecHitsHandle,
1548  *theEcalChStatus,
1549  geo,
1550  caloTopology,
1551  sevlv,
1552  3,
1553  3,
1554  0.090,
1555  0.450,
1556  tMinE_,
1557  tMaxE_);
1558  e9x9_30SigP = spr::eECALmatrix(isoCell,
1559  barrelRecHitsHandle,
1560  endcapRecHitsHandle,
1561  *theEcalChStatus,
1562  geo,
1563  caloTopology,
1564  sevlv,
1565  4,
1566  4,
1567  0.090,
1568  0.450,
1569  tMinE_,
1570  tMaxE_);
1571  e11x11_30SigP = spr::eECALmatrix(isoCell,
1572  barrelRecHitsHandle,
1573  endcapRecHitsHandle,
1574  *theEcalChStatus,
1575  geo,
1576  caloTopology,
1577  sevlv,
1578  5,
1579  5,
1580  0.090,
1581  0.450,
1582  tMinE_,
1583  tMaxE_);
1584  e15x15_30SigP = spr::eECALmatrix(isoCell,
1585  barrelRecHitsHandle,
1586  endcapRecHitsHandle,
1587  *theEcalChStatus,
1588  geo,
1589  caloTopology,
1590  sevlv,
1591  7,
1592  7,
1593  0.090,
1594  0.450,
1595  tMinE_,
1596  tMaxE_);
1597  if (myverbose_ == 2) {
1598  edm::LogVerbatim("IsoTrack") << "clean ecal rechit ";
1599  edm::LogVerbatim("IsoTrack") << "e7x7 " << e7x7P.first << " e9x9 " << e9x9P.first << " e11x11 "
1600  << e11x11P.first << " e15x15 " << e15x15P.first;
1601  edm::LogVerbatim("IsoTrack") << "e7x7_10Sig " << e7x7_10SigP.first << " e11x11_10Sig " << e11x11_10SigP.first
1602  << " e15x15_10Sig " << e15x15_10SigP.first;
1603  }
1604 
1605  if (doMC_) {
1606  // check the energy from SimHits
1608  iEvent, isoCell, geo, caloTopology, pcaloeb, pcaloee, SimTk, SimVtx, pTrack, *associate, 1, 1, simInfo3x3);
1610  iEvent, isoCell, geo, caloTopology, pcaloeb, pcaloee, SimTk, SimVtx, pTrack, *associate, 2, 2, simInfo5x5);
1612  iEvent, isoCell, geo, caloTopology, pcaloeb, pcaloee, SimTk, SimVtx, pTrack, *associate, 3, 3, simInfo7x7);
1614  iEvent, isoCell, geo, caloTopology, pcaloeb, pcaloee, SimTk, SimVtx, pTrack, *associate, 4, 4, simInfo9x9);
1615  spr::eECALSimInfo(iEvent,
1616  isoCell,
1617  geo,
1618  caloTopology,
1619  pcaloeb,
1620  pcaloee,
1621  SimTk,
1622  SimVtx,
1623  pTrack,
1624  *associate,
1625  5,
1626  5,
1627  simInfo11x11);
1628  spr::eECALSimInfo(iEvent,
1629  isoCell,
1630  geo,
1631  caloTopology,
1632  pcaloeb,
1633  pcaloee,
1634  SimTk,
1635  SimVtx,
1636  pTrack,
1637  *associate,
1638  6,
1639  6,
1640  simInfo13x13);
1641  spr::eECALSimInfo(iEvent,
1642  isoCell,
1643  geo,
1644  caloTopology,
1645  pcaloeb,
1646  pcaloee,
1647  SimTk,
1648  SimVtx,
1649  pTrack,
1650  *associate,
1651  7,
1652  7,
1653  simInfo15x15,
1654  150.0,
1655  false);
1656  spr::eECALSimInfo(iEvent,
1657  isoCell,
1658  geo,
1659  caloTopology,
1660  pcaloeb,
1661  pcaloee,
1662  SimTk,
1663  SimVtx,
1664  pTrack,
1665  *associate,
1666  10,
1667  10,
1668  simInfo21x21);
1669  spr::eECALSimInfo(iEvent,
1670  isoCell,
1671  geo,
1672  caloTopology,
1673  pcaloeb,
1674  pcaloee,
1675  SimTk,
1676  SimVtx,
1677  pTrack,
1678  *associate,
1679  12,
1680  12,
1681  simInfo25x25);
1682  spr::eECALSimInfo(iEvent,
1683  isoCell,
1684  geo,
1685  caloTopology,
1686  pcaloeb,
1687  pcaloee,
1688  SimTk,
1689  SimVtx,
1690  pTrack,
1691  *associate,
1692  15,
1693  15,
1694  simInfo31x31);
1695 
1696  trkEcalEne = spr::eCaloSimInfo(iEvent, geo, pcaloeb, pcaloee, SimTk, SimVtx, pTrack, *associate);
1697  if (myverbose_ == 1) {
1698  edm::LogVerbatim("IsoTrack") << "Track momentum " << pt1;
1699  edm::LogVerbatim("IsoTrack") << "ecal siminfo ";
1700  edm::LogVerbatim("IsoTrack") << "simInfo3x3: eTotal " << simInfo3x3.eTotal << " eMatched "
1701  << simInfo3x3.eMatched << " eRest " << simInfo3x3.eRest << " eGamma "
1702  << simInfo3x3.eGamma << " eNeutralHad " << simInfo3x3.eNeutralHad
1703  << " eChargedHad " << simInfo3x3.eChargedHad;
1704  edm::LogVerbatim("IsoTrack") << "simInfo5x5: eTotal " << simInfo5x5.eTotal << " eMatched "
1705  << simInfo5x5.eMatched << " eRest " << simInfo5x5.eRest << " eGamma "
1706  << simInfo5x5.eGamma << " eNeutralHad " << simInfo5x5.eNeutralHad
1707  << " eChargedHad " << simInfo5x5.eChargedHad;
1708  edm::LogVerbatim("IsoTrack") << "simInfo7x7: eTotal " << simInfo7x7.eTotal << " eMatched "
1709  << simInfo7x7.eMatched << " eRest " << simInfo7x7.eRest << " eGamma "
1710  << simInfo7x7.eGamma << " eNeutralHad " << simInfo7x7.eNeutralHad
1711  << " eChargedHad " << simInfo7x7.eChargedHad;
1712  edm::LogVerbatim("IsoTrack") << "simInfo9x9: eTotal " << simInfo9x9.eTotal << " eMatched "
1713  << simInfo9x9.eMatched << " eRest " << simInfo9x9.eRest << " eGamma "
1714  << simInfo9x9.eGamma << " eNeutralHad " << simInfo9x9.eNeutralHad
1715  << " eChargedHad " << simInfo9x9.eChargedHad;
1716  edm::LogVerbatim("IsoTrack") << "simInfo11x11: eTotal " << simInfo11x11.eTotal << " eMatched "
1717  << simInfo11x11.eMatched << " eRest " << simInfo11x11.eRest << " eGamma "
1718  << simInfo11x11.eGamma << " eNeutralHad " << simInfo11x11.eNeutralHad
1719  << " eChargedHad " << simInfo11x11.eChargedHad;
1720  edm::LogVerbatim("IsoTrack") << "simInfo15x15: eTotal " << simInfo15x15.eTotal << " eMatched "
1721  << simInfo15x15.eMatched << " eRest " << simInfo15x15.eRest << " eGamma "
1722  << simInfo15x15.eGamma << " eNeutralHad " << simInfo15x15.eNeutralHad
1723  << " eChargedHad " << simInfo15x15.eChargedHad;
1724  edm::LogVerbatim("IsoTrack") << "simInfo31x31: eTotal " << simInfo31x31.eTotal << " eMatched "
1725  << simInfo31x31.eMatched << " eRest " << simInfo31x31.eRest << " eGamma "
1726  << simInfo31x31.eGamma << " eNeutralHad " << simInfo31x31.eNeutralHad
1727  << " eChargedHad " << simInfo31x31.eChargedHad;
1728  edm::LogVerbatim("IsoTrack") << "trkEcalEne" << trkEcalEne;
1729  }
1730  }
1731 
1732  // ======= Get HCAL information
1733  double hcalScale = 1.0;
1734  if (std::abs(pTrack->eta()) < 1.4) {
1735  hcalScale = 120.0;
1736  } else {
1737  hcalScale = 135.0;
1738  }
1739 
1740  double maxNearHcalP3x3 = -1, maxNearHcalP5x5 = -1, maxNearHcalP7x7 = -1;
1741  maxNearHcalP3x3 = spr::chargeIsolationHcal(nTracks, trkCaloDets, theHBHETopology, 1, 1);
1742  maxNearHcalP5x5 = spr::chargeIsolationHcal(nTracks, trkCaloDets, theHBHETopology, 2, 2);
1743  maxNearHcalP7x7 = spr::chargeIsolationHcal(nTracks, trkCaloDets, theHBHETopology, 3, 3);
1744 
1745  double h3x3 = 0, h5x5 = 0, h7x7 = 0;
1746  double h3x3Sig = 0, h5x5Sig = 0, h7x7Sig = 0;
1747  double trkHcalEne = 0;
1748  spr::caloSimInfo hsimInfo3x3, hsimInfo5x5, hsimInfo7x7;
1749 
1750  if (trkDetItr->okHCAL) {
1751  const DetId ClosestCell(trkDetItr->detIdHCAL);
1752  // bool includeHO=false, bool algoNew=true, bool debug=false
1753  h3x3 = spr::eHCALmatrix(
1754  theHBHETopology, ClosestCell, hbhe, 1, 1, false, true, -100.0, -100.0, -100.0, -100.0, tMinH_, tMaxH_);
1755  h5x5 = spr::eHCALmatrix(
1756  theHBHETopology, ClosestCell, hbhe, 2, 2, false, true, -100.0, -100.0, -100.0, -100.0, tMinH_, tMaxH_);
1757  h7x7 = spr::eHCALmatrix(
1758  theHBHETopology, ClosestCell, hbhe, 3, 3, false, true, -100.0, -100.0, -100.0, -100.0, tMinH_, tMaxH_);
1759  h3x3Sig = spr::eHCALmatrix(
1760  theHBHETopology, ClosestCell, hbhe, 1, 1, false, true, 0.7, 0.8, -100.0, -100.0, tMinH_, tMaxH_);
1761  h5x5Sig = spr::eHCALmatrix(
1762  theHBHETopology, ClosestCell, hbhe, 2, 2, false, true, 0.7, 0.8, -100.0, -100.0, tMinH_, tMaxH_);
1763  h7x7Sig = spr::eHCALmatrix(
1764  theHBHETopology, ClosestCell, hbhe, 3, 3, false, true, 0.7, 0.8, -100.0, -100.0, tMinH_, tMaxH_);
1765 
1766  if (myverbose_ == 2) {
1767  edm::LogVerbatim("IsoTrack") << "HCAL 3x3 " << h3x3 << " " << h3x3Sig << " 5x5 " << h5x5 << " " << h5x5Sig
1768  << " 7x7 " << h7x7 << " " << h7x7Sig;
1769  }
1770 
1771  if (doMC_) {
1773  iEvent, theHBHETopology, ClosestCell, geo, pcalohh, SimTk, SimVtx, pTrack, *associate, 1, 1, hsimInfo3x3);
1775  iEvent, theHBHETopology, ClosestCell, geo, pcalohh, SimTk, SimVtx, pTrack, *associate, 2, 2, hsimInfo5x5);
1776  spr::eHCALSimInfo(iEvent,
1777  theHBHETopology,
1778  ClosestCell,
1779  geo,
1780  pcalohh,
1781  SimTk,
1782  SimVtx,
1783  pTrack,
1784  *associate,
1785  3,
1786  3,
1787  hsimInfo7x7,
1788  150.0,
1789  false,
1790  false);
1791  trkHcalEne = spr::eCaloSimInfo(iEvent, geo, pcalohh, SimTk, SimVtx, pTrack, *associate);
1792  if (myverbose_ == 1) {
1793  edm::LogVerbatim("IsoTrack") << "Hcal siminfo ";
1794  edm::LogVerbatim("IsoTrack")
1795  << "hsimInfo3x3: eTotal " << hsimInfo3x3.eTotal << " eMatched " << hsimInfo3x3.eMatched << " eRest "
1796  << hsimInfo3x3.eRest << " eGamma " << hsimInfo3x3.eGamma << " eNeutralHad " << hsimInfo3x3.eNeutralHad
1797  << " eChargedHad " << hsimInfo3x3.eChargedHad;
1798  edm::LogVerbatim("IsoTrack")
1799  << "hsimInfo5x5: eTotal " << hsimInfo5x5.eTotal << " eMatched " << hsimInfo5x5.eMatched << " eRest "
1800  << hsimInfo5x5.eRest << " eGamma " << hsimInfo5x5.eGamma << " eNeutralHad " << hsimInfo5x5.eNeutralHad
1801  << " eChargedHad " << hsimInfo5x5.eChargedHad;
1802  edm::LogVerbatim("IsoTrack")
1803  << "hsimInfo7x7: eTotal " << hsimInfo7x7.eTotal << " eMatched " << hsimInfo7x7.eMatched << " eRest "
1804  << hsimInfo7x7.eRest << " eGamma " << hsimInfo7x7.eGamma << " eNeutralHad " << hsimInfo7x7.eNeutralHad
1805  << " eChargedHad " << hsimInfo7x7.eChargedHad;
1806  edm::LogVerbatim("IsoTrack") << "trkHcalEne " << trkHcalEne;
1807  }
1808  }
1809 
1810  // debug the ecal and hcal matrix
1811  if (myverbose_ == 4) {
1812  edm::LogVerbatim("IsoTrack") << "Run " << iEvent.id().run() << " Event " << iEvent.id().event();
1813  std::vector<std::pair<DetId, double> > v7x7 =
1814  spr::eHCALmatrixCell(theHBHETopology, ClosestCell, hbhe, 3, 3, false, false);
1815  double sumv = 0.0;
1816 
1817  for (unsigned int iv = 0; iv < v7x7.size(); iv++) {
1818  sumv += v7x7[iv].second;
1819  }
1820  edm::LogVerbatim("IsoTrack") << "h7x7 " << h7x7 << " v7x7 " << sumv << " in " << v7x7.size();
1821  for (unsigned int iv = 0; iv < v7x7.size(); iv++) {
1822  HcalDetId id = v7x7[iv].first;
1823  edm::LogVerbatim("IsoTrack") << " Cell " << iv << " 0x" << std::hex << v7x7[iv].first() << std::dec << " "
1824  << id << " Energy " << v7x7[iv].second;
1825  }
1826  }
1827  }
1828  if (doMC_) {
1829  trkHcalEne = spr::eCaloSimInfo(iEvent, geo, pcalohh, SimTk, SimVtx, pTrack, *associate);
1830  }
1831 
1832  // ====================================================================================================
1833  // get diff between track outermost hit position and the propagation point at outermost surface of tracker
1834  std::pair<math::XYZPoint, double> point2_TK0 = spr::propagateTrackerEnd(pTrack, bField, false);
1835  math::XYZPoint diff(pTrack->outerPosition().X() - point2_TK0.first.X(),
1836  pTrack->outerPosition().Y() - point2_TK0.first.Y(),
1837  pTrack->outerPosition().Z() - point2_TK0.first.Z());
1838  double trackOutPosOutHitDr = diff.R();
1839  double trackL = point2_TK0.second;
1840  if (myverbose_ == 5) {
1841  edm::LogVerbatim("IsoTrack") << " propagted " << point2_TK0.first << " " << point2_TK0.first.eta() << " "
1842  << point2_TK0.first.phi();
1843  edm::LogVerbatim("IsoTrack") << " outerPosition() " << pTrack->outerPosition() << " "
1844  << pTrack->outerPosition().eta() << " " << pTrack->outerPosition().phi();
1845  edm::LogVerbatim("IsoTrack") << "diff " << diff << " diffR " << diff.R() << " diffR/L "
1846  << diff.R() / point2_TK0.second;
1847  }
1848 
1849  for (unsigned int ind = 0; ind < recVtxs->size(); ind++) {
1850  if (!((*recVtxs)[ind].isFake())) {
1851  reco::Vertex::trackRef_iterator vtxTrack = (*recVtxs)[ind].tracks_begin();
1852  if (reco::deltaR(eta1, phi1, (*vtxTrack)->eta(), (*vtxTrack)->phi()) < 0.01)
1853  t_trackPVIdx->push_back(ind);
1854  else
1855  t_trackPVIdx->push_back(-1);
1856  }
1857  }
1858 
1859  // Fill the tree Branches here
1860  t_trackPVIdx->push_back(pVtxTkId);
1861  t_trackP->push_back(p1);
1862  t_trackPt->push_back(pt1);
1863  t_trackEta->push_back(eta1);
1864  t_trackPhi->push_back(phi1);
1865  t_trackEcalEta->push_back(etaEcal1);
1866  t_trackEcalPhi->push_back(phiEcal1);
1867  t_trackHcalEta->push_back(etaHcal1);
1868  t_trackHcalPhi->push_back(phiHcal1);
1869  t_trackDxy->push_back(dxy1);
1870  t_trackDz->push_back(dz1);
1871  t_trackDxyBS->push_back(dxybs1);
1872  t_trackDzBS->push_back(dzbs1);
1873  t_trackDxyPV->push_back(dxypv1);
1874  t_trackDzPV->push_back(dzpv1);
1875  t_trackChiSq->push_back(chisq1);
1876  t_trackNOuterHits->push_back(nOuterHits);
1877  t_NLayersCrossed->push_back(nLayersCrossed);
1878 
1879  t_trackHitsTOB->push_back(hitp.stripTOBLayersWithMeasurement());
1880  t_trackHitsTEC->push_back(hitp.stripTECLayersWithMeasurement());
1887 
1894 
1903  t_trackOutPosOutHitDr->push_back(trackOutPosOutHitDr);
1904  t_trackL->push_back(trackL);
1905 
1906  t_maxNearP31x31->push_back(maxNearP31x31);
1907  t_maxNearP21x21->push_back(maxNearP21x21);
1908 
1909  t_ecalSpike11x11->push_back(e11x11P.second);
1910  t_e7x7->push_back(e7x7P.first);
1911  t_e9x9->push_back(e9x9P.first);
1912  t_e11x11->push_back(e11x11P.first);
1913  t_e15x15->push_back(e15x15P.first);
1914 
1915  t_e7x7_10Sig->push_back(e7x7_10SigP.first);
1916  t_e9x9_10Sig->push_back(e9x9_10SigP.first);
1917  t_e11x11_10Sig->push_back(e11x11_10SigP.first);
1918  t_e15x15_10Sig->push_back(e15x15_10SigP.first);
1919  t_e7x7_15Sig->push_back(e7x7_15SigP.first);
1920  t_e9x9_15Sig->push_back(e9x9_15SigP.first);
1921  t_e11x11_15Sig->push_back(e11x11_15SigP.first);
1922  t_e15x15_15Sig->push_back(e15x15_15SigP.first);
1923  t_e7x7_20Sig->push_back(e7x7_20SigP.first);
1924  t_e9x9_20Sig->push_back(e9x9_20SigP.first);
1925  t_e11x11_20Sig->push_back(e11x11_20SigP.first);
1926  t_e15x15_20Sig->push_back(e15x15_20SigP.first);
1927  t_e7x7_25Sig->push_back(e7x7_25SigP.first);
1928  t_e9x9_25Sig->push_back(e9x9_25SigP.first);
1929  t_e11x11_25Sig->push_back(e11x11_25SigP.first);
1930  t_e15x15_25Sig->push_back(e15x15_25SigP.first);
1931  t_e7x7_30Sig->push_back(e7x7_30SigP.first);
1932  t_e9x9_30Sig->push_back(e9x9_30SigP.first);
1933  t_e11x11_30Sig->push_back(e11x11_30SigP.first);
1934  t_e15x15_30Sig->push_back(e15x15_30SigP.first);
1935 
1936  if (doMC_) {
1937  t_esim7x7->push_back(simInfo7x7.eTotal);
1938  t_esim9x9->push_back(simInfo9x9.eTotal);
1939  t_esim11x11->push_back(simInfo11x11.eTotal);
1940  t_esim15x15->push_back(simInfo15x15.eTotal);
1941 
1942  t_esim7x7Matched->push_back(simInfo7x7.eMatched);
1943  t_esim9x9Matched->push_back(simInfo9x9.eMatched);
1944  t_esim11x11Matched->push_back(simInfo11x11.eMatched);
1945  t_esim15x15Matched->push_back(simInfo15x15.eMatched);
1946 
1947  t_esim7x7Rest->push_back(simInfo7x7.eRest);
1948  t_esim9x9Rest->push_back(simInfo9x9.eRest);
1949  t_esim11x11Rest->push_back(simInfo11x11.eRest);
1950  t_esim15x15Rest->push_back(simInfo15x15.eRest);
1951 
1952  t_esim7x7Photon->push_back(simInfo7x7.eGamma);
1953  t_esim9x9Photon->push_back(simInfo9x9.eGamma);
1954  t_esim11x11Photon->push_back(simInfo11x11.eGamma);
1955  t_esim15x15Photon->push_back(simInfo15x15.eGamma);
1956 
1957  t_esim7x7NeutHad->push_back(simInfo7x7.eNeutralHad);
1958  t_esim9x9NeutHad->push_back(simInfo9x9.eNeutralHad);
1959  t_esim11x11NeutHad->push_back(simInfo11x11.eNeutralHad);
1960  t_esim15x15NeutHad->push_back(simInfo15x15.eNeutralHad);
1961 
1962  t_esim7x7CharHad->push_back(simInfo7x7.eChargedHad);
1963  t_esim9x9CharHad->push_back(simInfo9x9.eChargedHad);
1964  t_esim11x11CharHad->push_back(simInfo11x11.eChargedHad);
1965  t_esim15x15CharHad->push_back(simInfo15x15.eChargedHad);
1966 
1967  t_trkEcalEne->push_back(trkEcalEne);
1968  t_simTrackP->push_back(simTrackP);
1969  t_esimPdgId->push_back(simInfo11x11.pdgMatched);
1970  }
1971 
1972  t_maxNearHcalP3x3->push_back(maxNearHcalP3x3);
1973  t_maxNearHcalP5x5->push_back(maxNearHcalP5x5);
1974  t_maxNearHcalP7x7->push_back(maxNearHcalP7x7);
1975 
1976  t_h3x3->push_back(h3x3);
1977  t_h5x5->push_back(h5x5);
1978  t_h7x7->push_back(h7x7);
1979  t_h3x3Sig->push_back(h3x3Sig);
1980  t_h5x5Sig->push_back(h5x5Sig);
1981  t_h7x7Sig->push_back(h7x7Sig);
1982 
1983  t_infoHcal->push_back(trkDetItr->okHCAL);
1984  if (doMC_) {
1985  t_trkHcalEne->push_back(hcalScale * trkHcalEne);
1986 
1987  t_hsim3x3->push_back(hcalScale * hsimInfo3x3.eTotal);
1988  t_hsim5x5->push_back(hcalScale * hsimInfo5x5.eTotal);
1989  t_hsim7x7->push_back(hcalScale * hsimInfo7x7.eTotal);
1990 
1991  t_hsim3x3Matched->push_back(hcalScale * hsimInfo3x3.eMatched);
1992  t_hsim5x5Matched->push_back(hcalScale * hsimInfo5x5.eMatched);
1993  t_hsim7x7Matched->push_back(hcalScale * hsimInfo7x7.eMatched);
1994 
1995  t_hsim3x3Rest->push_back(hcalScale * hsimInfo3x3.eRest);
1996  t_hsim5x5Rest->push_back(hcalScale * hsimInfo5x5.eRest);
1997  t_hsim7x7Rest->push_back(hcalScale * hsimInfo7x7.eRest);
1998 
1999  t_hsim3x3Photon->push_back(hcalScale * hsimInfo3x3.eGamma);
2000  t_hsim5x5Photon->push_back(hcalScale * hsimInfo5x5.eGamma);
2001  t_hsim7x7Photon->push_back(hcalScale * hsimInfo7x7.eGamma);
2002 
2003  t_hsim3x3NeutHad->push_back(hcalScale * hsimInfo3x3.eNeutralHad);
2004  t_hsim5x5NeutHad->push_back(hcalScale * hsimInfo5x5.eNeutralHad);
2005  t_hsim7x7NeutHad->push_back(hcalScale * hsimInfo7x7.eNeutralHad);
2006 
2007  t_hsim3x3CharHad->push_back(hcalScale * hsimInfo3x3.eChargedHad);
2008  t_hsim5x5CharHad->push_back(hcalScale * hsimInfo5x5.eChargedHad);
2009  t_hsim7x7CharHad->push_back(hcalScale * hsimInfo7x7.eChargedHad);
2010  }
2011 
2012  } // if loosely isolated track
2013  } // check p1/eta
2014  } // loop over track collection
2015 
2016  if (haveIsoTrack)
2017  tree_->Fill();
2018 }
2019 
2020 // ----- method called once each job just before starting event loop ----
2022  nEventProc_ = 0;
2023  nbad_ = 0;
2024  initL1_ = false;
2025  double tempgen_TH[NPBins + 1] = {
2026  0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 9.0, 11.0, 15.0, 20.0, 30.0, 50.0, 75.0, 100.0};
2027 
2028  for (unsigned int i = 0; i < NPBins + 1; i++)
2029  genPartPBins[i] = tempgen_TH[i];
2030 
2031  double tempgen_Eta[NEtaBins + 1] = {0.0, 1.131, 1.653, 2.172};
2032 
2033  for (unsigned int i = 0; i < NEtaBins + 1; i++)
2034  genPartEtaBins[i] = tempgen_Eta[i];
2035 
2036  bookHistograms();
2037 }
2038 
2039 // ----- method called once each job just after ending the event loop ----
2041  if (L1TriggerAlgoInfo_) {
2042  std::map<std::pair<unsigned int, std::string>, int>::iterator itr;
2043  for (itr = l1AlgoMap_.begin(); itr != l1AlgoMap_.end(); itr++) {
2044  edm::LogVerbatim("IsoTrack") << " ****endjob**** " << (itr->first).first << " " << (itr->first).second << " "
2045  << itr->second;
2046  int ibin = (itr->first).first;
2047  TString name((itr->first).second);
2048  h_L1AlgoNames->GetXaxis()->SetBinLabel(ibin + 1, name);
2049  }
2050  edm::LogVerbatim("IsoTrack") << "Number of Events Processed " << nEventProc_;
2051  }
2052 }
2053 
2054 //========================================================================================================
2055 
2057  t_PVx->clear();
2058  t_PVy->clear();
2059  t_PVz->clear();
2060  t_PVisValid->clear();
2061  t_PVndof->clear();
2062  t_PVNTracks->clear();
2063  t_PVNTracksWt->clear();
2064  t_PVTracksSumPt->clear();
2065  t_PVTracksSumPtWt->clear();
2066  t_PVNTracksHP->clear();
2067  t_PVNTracksHPWt->clear();
2068  t_PVTracksSumPtHP->clear();
2069  t_PVTracksSumPtHPWt->clear();
2070 
2071  for (int i = 0; i < 128; i++)
2072  t_L1Decision[i] = 0;
2073  t_L1AlgoNames->clear();
2074  t_L1PreScale->clear();
2075 
2076  t_L1CenJetPt->clear();
2077  t_L1CenJetEta->clear();
2078  t_L1CenJetPhi->clear();
2079  t_L1FwdJetPt->clear();
2080  t_L1FwdJetEta->clear();
2081  t_L1FwdJetPhi->clear();
2082  t_L1TauJetPt->clear();
2083  t_L1TauJetEta->clear();
2084  t_L1TauJetPhi->clear();
2085  t_L1MuonPt->clear();
2086  t_L1MuonEta->clear();
2087  t_L1MuonPhi->clear();
2088  t_L1IsoEMPt->clear();
2089  t_L1IsoEMEta->clear();
2090  t_L1IsoEMPhi->clear();
2091  t_L1NonIsoEMPt->clear();
2092  t_L1NonIsoEMEta->clear();
2093  t_L1NonIsoEMPhi->clear();
2094  t_L1METPt->clear();
2095  t_L1METEta->clear();
2096  t_L1METPhi->clear();
2097 
2098  t_jetPt->clear();
2099  t_jetEta->clear();
2100  t_jetPhi->clear();
2101  t_nTrksJetCalo->clear();
2102  t_nTrksJetVtx->clear();
2103 
2104  t_trackPAll->clear();
2105  t_trackEtaAll->clear();
2106  t_trackPhiAll->clear();
2107  t_trackPdgIdAll->clear();
2108  t_trackPtAll->clear();
2109  t_trackDxyAll->clear();
2110  t_trackDzAll->clear();
2111  t_trackDxyPVAll->clear();
2112  t_trackDzPVAll->clear();
2113  t_trackChiSqAll->clear();
2114 
2115  t_trackP->clear();
2116  t_trackPt->clear();
2117  t_trackEta->clear();
2118  t_trackPhi->clear();
2119  t_trackEcalEta->clear();
2120  t_trackEcalPhi->clear();
2121  t_trackHcalEta->clear();
2122  t_trackHcalPhi->clear();
2123  t_NLayersCrossed->clear();
2124  t_trackNOuterHits->clear();
2125  t_trackDxy->clear();
2126  t_trackDxyBS->clear();
2127  t_trackDz->clear();
2128  t_trackDzBS->clear();
2129  t_trackDxyPV->clear();
2130  t_trackDzPV->clear();
2131  t_trackChiSq->clear();
2132  t_trackPVIdx->clear();
2133  t_trackHitsTOB->clear();
2134  t_trackHitsTEC->clear();
2135  t_trackHitInMissTOB->clear();
2136  t_trackHitInMissTEC->clear();
2137  t_trackHitInMissTIB->clear();
2138  t_trackHitInMissTID->clear();
2139  t_trackHitInMissTIBTID->clear();
2140  t_trackHitOutMissTOB->clear();
2141  t_trackHitOutMissTEC->clear();
2142  t_trackHitOutMissTIB->clear();
2143  t_trackHitOutMissTID->clear();
2144  t_trackHitOutMissTOBTEC->clear();
2145  t_trackHitInMeasTOB->clear();
2146  t_trackHitInMeasTEC->clear();
2147  t_trackHitInMeasTIB->clear();
2148  t_trackHitInMeasTID->clear();
2149  t_trackHitOutMeasTOB->clear();
2150  t_trackHitOutMeasTEC->clear();
2151  t_trackHitOutMeasTIB->clear();
2152  t_trackHitOutMeasTID->clear();
2153  t_trackOutPosOutHitDr->clear();
2154  t_trackL->clear();
2155 
2156  t_maxNearP31x31->clear();
2157  t_maxNearP21x21->clear();
2158 
2159  t_ecalSpike11x11->clear();
2160  t_e7x7->clear();
2161  t_e9x9->clear();
2162  t_e11x11->clear();
2163  t_e15x15->clear();
2164 
2165  t_e7x7_10Sig->clear();
2166  t_e9x9_10Sig->clear();
2167  t_e11x11_10Sig->clear();
2168  t_e15x15_10Sig->clear();
2169  t_e7x7_15Sig->clear();
2170  t_e9x9_15Sig->clear();
2171  t_e11x11_15Sig->clear();
2172  t_e15x15_15Sig->clear();
2173  t_e7x7_20Sig->clear();
2174  t_e9x9_20Sig->clear();
2175  t_e11x11_20Sig->clear();
2176  t_e15x15_20Sig->clear();
2177  t_e7x7_25Sig->clear();
2178  t_e9x9_25Sig->clear();
2179  t_e11x11_25Sig->clear();
2180  t_e15x15_25Sig->clear();
2181  t_e7x7_30Sig->clear();
2182  t_e9x9_30Sig->clear();
2183  t_e11x11_30Sig->clear();
2184  t_e15x15_30Sig->clear();
2185 
2186  if (doMC_) {
2187  t_simTrackP->clear();
2188  t_esimPdgId->clear();
2189  t_trkEcalEne->clear();
2190 
2191  t_esim7x7->clear();
2192  t_esim9x9->clear();
2193  t_esim11x11->clear();
2194  t_esim15x15->clear();
2195 
2196  t_esim7x7Matched->clear();
2197  t_esim9x9Matched->clear();
2198  t_esim11x11Matched->clear();
2199  t_esim15x15Matched->clear();
2200 
2201  t_esim7x7Rest->clear();
2202  t_esim9x9Rest->clear();
2203  t_esim11x11Rest->clear();
2204  t_esim15x15Rest->clear();
2205 
2206  t_esim7x7Photon->clear();
2207  t_esim9x9Photon->clear();
2208  t_esim11x11Photon->clear();
2209  t_esim15x15Photon->clear();
2210 
2211  t_esim7x7NeutHad->clear();
2212  t_esim9x9NeutHad->clear();
2213  t_esim11x11NeutHad->clear();
2214  t_esim15x15NeutHad->clear();
2215 
2216  t_esim7x7CharHad->clear();
2217  t_esim9x9CharHad->clear();
2218  t_esim11x11CharHad->clear();
2219  t_esim15x15CharHad->clear();
2220  }
2221 
2222  t_maxNearHcalP3x3->clear();
2223  t_maxNearHcalP5x5->clear();
2224  t_maxNearHcalP7x7->clear();
2225 
2226  t_h3x3->clear();
2227  t_h5x5->clear();
2228  t_h7x7->clear();
2229  t_h3x3Sig->clear();
2230  t_h5x5Sig->clear();
2231  t_h7x7Sig->clear();
2232 
2233  t_infoHcal->clear();
2234 
2235  if (doMC_) {
2236  t_trkHcalEne->clear();
2237 
2238  t_hsim3x3->clear();
2239  t_hsim5x5->clear();
2240  t_hsim7x7->clear();
2241  t_hsim3x3Matched->clear();
2242  t_hsim5x5Matched->clear();
2243  t_hsim7x7Matched->clear();
2244  t_hsim3x3Rest->clear();
2245  t_hsim5x5Rest->clear();
2246  t_hsim7x7Rest->clear();
2247  t_hsim3x3Photon->clear();
2248  t_hsim5x5Photon->clear();
2249  t_hsim7x7Photon->clear();
2250  t_hsim3x3NeutHad->clear();
2251  t_hsim5x5NeutHad->clear();
2252  t_hsim7x7NeutHad->clear();
2253  t_hsim3x3CharHad->clear();
2254  t_hsim5x5CharHad->clear();
2255  t_hsim7x7CharHad->clear();
2256  }
2257 }
2258 
2260  char hname[100], htit[100];
2261 
2263  TFileDirectory dir = fs->mkdir("nearMaxTrackP");
2264 
2265  for (unsigned int ieta = 0; ieta < NEtaBins; ieta++) {
2266  double lowEta = -5.0, highEta = 5.0;
2267  lowEta = genPartEtaBins[ieta];
2268  highEta = genPartEtaBins[ieta + 1];
2269 
2270  for (unsigned int ipt = 0; ipt < NPBins; ipt++) {
2271  double lowP = 0.0, highP = 300.0;
2272  lowP = genPartPBins[ipt];
2273  highP = genPartPBins[ipt + 1];
2274  sprintf(hname, "h_maxNearP31x31_ptBin%i_etaBin%i", ipt, ieta);
2275  sprintf(htit, "maxNearP in 31x31 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP);
2276  h_maxNearP31x31[ipt][ieta] = dir.make<TH1F>(hname, htit, 220, -2.0, 20.0);
2277  h_maxNearP31x31[ipt][ieta]->Sumw2();
2278  sprintf(hname, "h_maxNearP25x25_ptBin%i_etaBin%i", ipt, ieta);
2279  sprintf(htit, "maxNearP in 25x25 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP);
2280  h_maxNearP25x25[ipt][ieta] = dir.make<TH1F>(hname, htit, 220, -2.0, 20.0);
2281  h_maxNearP25x25[ipt][ieta]->Sumw2();
2282  sprintf(hname, "h_maxNearP21x21_ptBin%i_etaBin%i", ipt, ieta);
2283  sprintf(htit, "maxNearP in 21x21 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP);
2284  h_maxNearP21x21[ipt][ieta] = dir.make<TH1F>(hname, htit, 220, -2.0, 20.0);
2285  h_maxNearP21x21[ipt][ieta]->Sumw2();
2286  sprintf(hname, "h_maxNearP15x15_ptBin%i_etaBin%i", ipt, ieta);
2287  sprintf(htit, "maxNearP in 15x15 (%3.2f<|#eta|<%3.2f), (%2.0f<trkP<%3.0f)", lowEta, highEta, lowP, highP);
2288  h_maxNearP15x15[ipt][ieta] = dir.make<TH1F>(hname, htit, 220, -2.0, 20.0);
2289  h_maxNearP15x15[ipt][ieta]->Sumw2();
2290  }
2291  }
2292 
2293  h_L1AlgoNames = fs->make<TH1I>("h_L1AlgoNames", "h_L1AlgoNames:Bin Labels", 128, -0.5, 127.5);
2294 
2295  // Reconstructed Tracks
2296 
2297  h_PVTracksWt = fs->make<TH1F>("h_PVTracksWt", "h_PVTracksWt", 600, -0.1, 1.1);
2298 
2299  h_nTracks = fs->make<TH1F>("h_nTracks", "h_nTracks", 1000, -0.5, 999.5);
2300 
2301  sprintf(hname, "h_recEtaPt_0");
2302  sprintf(htit, "h_recEtaPt (all tracks Eta vs pT)");
2303  h_recEtaPt_0 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2304 
2305  sprintf(hname, "h_recEtaP_0");
2306  sprintf(htit, "h_recEtaP (all tracks Eta vs pT)");
2307  h_recEtaP_0 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2308 
2309  h_recPt_0 = fs->make<TH1F>("h_recPt_0", "Pt (all tracks)", 15, genPartPBins);
2310  h_recP_0 = fs->make<TH1F>("h_recP_0", "P (all tracks)", 15, genPartPBins);
2311  h_recEta_0 = fs->make<TH1F>("h_recEta_0", "Eta (all tracks)", 60, -3.0, 3.0);
2312  h_recPhi_0 = fs->make<TH1F>("h_recPhi_0", "Phi (all tracks)", 100, -3.2, 3.2);
2313  //-------------------------
2314  sprintf(hname, "h_recEtaPt_1");
2315  sprintf(htit, "h_recEtaPt (all good tracks Eta vs pT)");
2316  h_recEtaPt_1 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2317 
2318  sprintf(hname, "h_recEtaP_1");
2319  sprintf(htit, "h_recEtaP (all good tracks Eta vs pT)");
2320  h_recEtaP_1 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2321 
2322  h_recPt_1 = fs->make<TH1F>("h_recPt_1", "Pt (all good tracks)", 15, genPartPBins);
2323  h_recP_1 = fs->make<TH1F>("h_recP_1", "P (all good tracks)", 15, genPartPBins);
2324  h_recEta_1 = fs->make<TH1F>("h_recEta_1", "Eta (all good tracks)", 60, -3.0, 3.0);
2325  h_recPhi_1 = fs->make<TH1F>("h_recPhi_1", "Phi (all good tracks)", 100, -3.2, 3.2);
2326  //-------------------------
2327  sprintf(hname, "h_recEtaPt_2");
2328  sprintf(htit, "h_recEtaPt (charge isolation Eta vs pT)");
2329  h_recEtaPt_2 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2330 
2331  sprintf(hname, "h_recEtaP_2");
2332  sprintf(htit, "h_recEtaP (charge isolation Eta vs pT)");
2333  h_recEtaP_2 = fs->make<TH2F>(hname, htit, 30, -3.0, 3.0, 15, genPartPBins);
2334 
2335  h_recPt_2 = fs->make<TH1F>("h_recPt_2", "Pt (charge isolation)", 15, genPartPBins);
2336  h_recP_2 = fs->make<TH1F>("h_recP_2", "P (charge isolation)", 15, genPartPBins);
2337  h_recEta_2 = fs->make<TH1F>("h_recEta_2", "Eta (charge isolation)", 60, -3.0, 3.0);
2338  h_recPhi_2 = fs->make<TH1F>("h_recPhi_2", "Phi (charge isolation)", 100, -3.2, 3.2);
2339 
2340  tree_ = fs->make<TTree>("tree", "tree");
2341  tree_->SetAutoSave(10000);
2342 
2343  tree_->Branch("t_EvtNo", &t_EvtNo, "t_EvtNo/I");
2344  tree_->Branch("t_RunNo", &t_RunNo, "t_RunNo/I");
2345  tree_->Branch("t_Lumi", &t_Lumi, "t_Lumi/I");
2346  tree_->Branch("t_Bunch", &t_Bunch, "t_Bunch/I");
2347 
2348  t_PVx = new std::vector<double>();
2349  t_PVy = new std::vector<double>();
2350  t_PVz = new std::vector<double>();
2351  t_PVisValid = new std::vector<int>();
2352  t_PVndof = new std::vector<int>();
2353  t_PVNTracks = new std::vector<int>();
2354  t_PVNTracksWt = new std::vector<int>();
2355  t_PVTracksSumPt = new std::vector<double>();
2356  t_PVTracksSumPtWt = new std::vector<double>();
2357  t_PVNTracksHP = new std::vector<int>();
2358  t_PVNTracksHPWt = new std::vector<int>();
2359  t_PVTracksSumPtHP = new std::vector<double>();
2360  t_PVTracksSumPtHPWt = new std::vector<double>();
2361 
2362  tree_->Branch("PVx", "std::vector<double>", &t_PVx);
2363  tree_->Branch("PVy", "std::vector<double>", &t_PVy);
2364  tree_->Branch("PVz", "std::vector<double>", &t_PVz);
2365  tree_->Branch("PVisValid", "std::vector<int>", &t_PVisValid);
2366  tree_->Branch("PVndof", "std::vector<int>", &t_PVndof);
2367  tree_->Branch("PVNTracks", "std::vector<int>", &t_PVNTracks);
2368  tree_->Branch("PVNTracksWt", "std::vector<int>", &t_PVNTracksWt);
2369  tree_->Branch("t_PVTracksSumPt", "std::vector<double>", &t_PVTracksSumPt);
2370  tree_->Branch("t_PVTracksSumPtWt", "std::vector<double>", &t_PVTracksSumPtWt);
2371  tree_->Branch("PVNTracksHP", "std::vector<int>", &t_PVNTracksHP);
2372  tree_->Branch("PVNTracksHPWt", "std::vector<int>", &t_PVNTracksHPWt);
2373  tree_->Branch("t_PVTracksSumPtHP", "std::vector<double>", &t_PVTracksSumPtHP);
2374  tree_->Branch("t_PVTracksSumPtHPWt", "std::vector<double>", &t_PVTracksSumPtHPWt);
2375 
2376  //----- L1Trigger information
2377  for (int i = 0; i < 128; i++)
2378  t_L1Decision[i] = 0;
2379  t_L1AlgoNames = new std::vector<std::string>();
2380  t_L1PreScale = new std::vector<int>();
2381  t_L1CenJetPt = new std::vector<double>();
2382  t_L1CenJetEta = new std::vector<double>();
2383  t_L1CenJetPhi = new std::vector<double>();
2384  t_L1FwdJetPt = new std::vector<double>();
2385  t_L1FwdJetEta = new std::vector<double>();
2386  t_L1FwdJetPhi = new std::vector<double>();
2387  t_L1TauJetPt = new std::vector<double>();
2388  t_L1TauJetEta = new std::vector<double>();
2389  t_L1TauJetPhi = new std::vector<double>();
2390  t_L1MuonPt = new std::vector<double>();
2391  t_L1MuonEta = new std::vector<double>();
2392  t_L1MuonPhi = new std::vector<double>();
2393  t_L1IsoEMPt = new std::vector<double>();
2394  t_L1IsoEMEta = new std::vector<double>();
2395  t_L1IsoEMPhi = new std::vector<double>();
2396  t_L1NonIsoEMPt = new std::vector<double>();
2397  t_L1NonIsoEMEta = new std::vector<double>();
2398  t_L1NonIsoEMPhi = new std::vector<double>();
2399  t_L1METPt = new std::vector<double>();
2400  t_L1METEta = new std::vector<double>();
2401  t_L1METPhi = new std::vector<double>();
2402 
2403  tree_->Branch("t_L1Decision", t_L1Decision, "t_L1Decision[128]/I");
2404  tree_->Branch("t_L1AlgoNames", "std::vector<string>", &t_L1AlgoNames);
2405  tree_->Branch("t_L1PreScale", "std::vector<int>", &t_L1PreScale);
2406  tree_->Branch("t_L1CenJetPt", "std::vector<double>", &t_L1CenJetPt);
2407  tree_->Branch("t_L1CenJetEta", "std::vector<double>", &t_L1CenJetEta);
2408  tree_->Branch("t_L1CenJetPhi", "std::vector<double>", &t_L1CenJetPhi);
2409  tree_->Branch("t_L1FwdJetPt", "std::vector<double>", &t_L1FwdJetPt);
2410  tree_->Branch("t_L1FwdJetEta", "std::vector<double>", &t_L1FwdJetEta);
2411  tree_->Branch("t_L1FwdJetPhi", "std::vector<double>", &t_L1FwdJetPhi);
2412  tree_->Branch("t_L1TauJetPt", "std::vector<double>", &t_L1TauJetPt);
2413  tree_->Branch("t_L1TauJetEta", "std::vector<double>", &t_L1TauJetEta);
2414  tree_->Branch("t_L1TauJetPhi", "std::vector<double>", &t_L1TauJetPhi);
2415  tree_->Branch("t_L1MuonPt", "std::vector<double>", &t_L1MuonPt);
2416  tree_->Branch("t_L1MuonEta", "std::vector<double>", &t_L1MuonEta);
2417  tree_->Branch("t_L1MuonPhi", "std::vector<double>", &t_L1MuonPhi);
2418  tree_->Branch("t_L1IsoEMPt", "std::vector<double>", &t_L1IsoEMPt);
2419  tree_->Branch("t_L1IsoEMEta", "std::vector<double>", &t_L1IsoEMEta);
2420  tree_->Branch("t_L1IsoEMPhi", "std::vector<double>", &t_L1IsoEMPhi);
2421  tree_->Branch("t_L1NonIsoEMPt", "std::vector<double>", &t_L1NonIsoEMPt);
2422  tree_->Branch("t_L1NonIsoEMEta", "std::vector<double>", &t_L1NonIsoEMEta);
2423  tree_->Branch("t_L1NonIsoEMPhi", "std::vector<double>", &t_L1NonIsoEMPhi);
2424  tree_->Branch("t_L1METPt", "std::vector<double>", &t_L1METPt);
2425  tree_->Branch("t_L1METEta", "std::vector<double>", &t_L1METEta);
2426  tree_->Branch("t_L1METPhi", "std::vector<double>", &t_L1METPhi);
2427 
2428  t_jetPt = new std::vector<double>();
2429  t_jetEta = new std::vector<double>();
2430  t_jetPhi = new std::vector<double>();
2431  t_nTrksJetCalo = new std::vector<double>();
2432  t_nTrksJetVtx = new std::vector<double>();
2433  tree_->Branch("t_jetPt", "std::vector<double>", &t_jetPt);
2434  tree_->Branch("t_jetEta", "std::vector<double>", &t_jetEta);
2435  tree_->Branch("t_jetPhi", "std::vector<double>", &t_jetPhi);
2436  tree_->Branch("t_nTrksJetCalo", "std::vector<double>", &t_nTrksJetCalo);
2437  tree_->Branch("t_nTrksJetVtx", "std::vector<double>", &t_nTrksJetVtx);
2438 
2439  t_trackPAll = new std::vector<double>();
2440  t_trackEtaAll = new std::vector<double>();
2441  t_trackPhiAll = new std::vector<double>();
2442  t_trackPdgIdAll = new std::vector<double>();
2443  t_trackPtAll = new std::vector<double>();
2444  t_trackDxyAll = new std::vector<double>();
2445  t_trackDzAll = new std::vector<double>();
2446  t_trackDxyPVAll = new std::vector<double>();
2447  t_trackDzPVAll = new std::vector<double>();
2448  t_trackChiSqAll = new std::vector<double>();
2449  tree_->Branch("t_trackPAll", "std::vector<double>", &t_trackPAll);
2450  tree_->Branch("t_trackPhiAll", "std::vector<double>", &t_trackPhiAll);
2451  tree_->Branch("t_trackEtaAll", "std::vector<double>", &t_trackEtaAll);
2452  tree_->Branch("t_trackPtAll", "std::vector<double>", &t_trackPtAll);
2453  tree_->Branch("t_trackDxyAll", "std::vector<double>", &t_trackDxyAll);
2454  tree_->Branch("t_trackDzAll", "std::vector<double>", &t_trackDzAll);
2455  tree_->Branch("t_trackDxyPVAll", "std::vector<double>", &t_trackDxyPVAll);
2456  tree_->Branch("t_trackDzPVAll", "std::vector<double>", &t_trackDzPVAll);
2457  tree_->Branch("t_trackChiSqAll", "std::vector<double>", &t_trackChiSqAll);
2458  //tree_->Branch("t_trackPdgIdAll", "std::vector<double>", &t_trackPdgIdAll);
2459 
2460  t_trackP = new std::vector<double>();
2461  t_trackPt = new std::vector<double>();
2462  t_trackEta = new std::vector<double>();
2463  t_trackPhi = new std::vector<double>();
2464  t_trackEcalEta = new std::vector<double>();
2465  t_trackEcalPhi = new std::vector<double>();
2466  t_trackHcalEta = new std::vector<double>();
2467  t_trackHcalPhi = new std::vector<double>();
2468  t_trackNOuterHits = new std::vector<int>();
2469  t_NLayersCrossed = new std::vector<int>();
2470  t_trackDxy = new std::vector<double>();
2471  t_trackDxyBS = new std::vector<double>();
2472  t_trackDz = new std::vector<double>();
2473  t_trackDzBS = new std::vector<double>();
2474  t_trackDxyPV = new std::vector<double>();
2475  t_trackDzPV = new std::vector<double>();
2476  t_trackPVIdx = new std::vector<int>();
2477  t_trackChiSq = new std::vector<double>();
2478  t_trackHitsTOB = new std::vector<int>();
2479  t_trackHitsTEC = new std::vector<int>();
2480  t_trackHitInMissTOB = new std::vector<int>();
2481  t_trackHitInMissTEC = new std::vector<int>();
2482  t_trackHitInMissTIB = new std::vector<int>();
2483  t_trackHitInMissTID = new std::vector<int>();
2484  t_trackHitOutMissTOB = new std::vector<int>();
2485  t_trackHitOutMissTEC = new std::vector<int>();
2486  t_trackHitOutMissTIB = new std::vector<int>();
2487  t_trackHitOutMissTID = new std::vector<int>();
2488  t_trackHitInMissTIBTID = new std::vector<int>();
2489  t_trackHitOutMissTOB = new std::vector<int>();
2490  t_trackHitOutMissTEC = new std::vector<int>();
2491  t_trackHitOutMissTIB = new std::vector<int>();
2492  t_trackHitOutMissTID = new std::vector<int>();
2493  t_trackHitOutMissTOBTEC = new std::vector<int>();
2494  t_trackHitInMeasTOB = new std::vector<int>();
2495  t_trackHitInMeasTEC = new std::vector<int>();
2496  t_trackHitInMeasTIB = new std::vector<int>();
2497  t_trackHitInMeasTID = new std::vector<int>();
2498  t_trackHitOutMeasTOB = new std::vector<int>();
2499  t_trackHitOutMeasTEC = new std::vector<int>();
2500  t_trackHitOutMeasTIB = new std::vector<int>();
2501  t_trackHitOutMeasTID = new std::vector<int>();
2502  t_trackOutPosOutHitDr = new std::vector<double>();
2503  t_trackL = new std::vector<double>();
2504 
2505  tree_->Branch("t_trackP", "std::vector<double>", &t_trackP);
2506  tree_->Branch("t_trackPt", "std::vector<double>", &t_trackPt);
2507  tree_->Branch("t_trackEta", "std::vector<double>", &t_trackEta);
2508  tree_->Branch("t_trackPhi", "std::vector<double>", &t_trackPhi);
2509  tree_->Branch("t_trackEcalEta", "std::vector<double>", &t_trackEcalEta);
2510  tree_->Branch("t_trackEcalPhi", "std::vector<double>", &t_trackEcalPhi);
2511  tree_->Branch("t_trackHcalEta", "std::vector<double>", &t_trackHcalEta);
2512  tree_->Branch("t_trackHcalPhi", "std::vector<double>", &t_trackHcalPhi);
2513 
2514  tree_->Branch("t_trackNOuterHits", "std::vector<int>", &t_trackNOuterHits);
2515  tree_->Branch("t_NLayersCrossed", "std::vector<int>", &t_NLayersCrossed);
2516  tree_->Branch("t_trackHitsTOB", "std::vector<int>", &t_trackHitsTOB);
2517  tree_->Branch("t_trackHitsTEC", "std::vector<int>", &t_trackHitsTEC);
2518  tree_->Branch("t_trackHitInMissTOB", "std::vector<int>", &t_trackHitInMissTOB);
2519  tree_->Branch("t_trackHitInMissTEC", "std::vector<int>", &t_trackHitInMissTEC);
2520  tree_->Branch("t_trackHitInMissTIB", "std::vector<int>", &t_trackHitInMissTIB);
2521  tree_->Branch("t_trackHitInMissTID", "std::vector<int>", &t_trackHitInMissTID);
2522  tree_->Branch("t_trackHitInMissTIBTID", "std::vector<int>", &t_trackHitInMissTIBTID);
2523  tree_->Branch("t_trackHitOutMissTOB", "std::vector<int>", &t_trackHitOutMissTOB);
2524  tree_->Branch("t_trackHitOutMissTEC", "std::vector<int>", &t_trackHitOutMissTEC);
2525  tree_->Branch("t_trackHitOutMissTIB", "std::vector<int>", &t_trackHitOutMissTIB);
2526  tree_->Branch("t_trackHitOutMissTID", "std::vector<int>", &t_trackHitOutMissTID);
2527  tree_->Branch("t_trackHitOutMissTOBTEC", "std::vector<int>", &t_trackHitOutMissTOBTEC);
2528  tree_->Branch("t_trackHitInMeasTOB", "std::vector<int>", &t_trackHitInMeasTOB);
2529  tree_->Branch("t_trackHitInMeasTEC", "std::vector<int>", &t_trackHitInMeasTEC);
2530  tree_->Branch("t_trackHitInMeasTIB", "std::vector<int>", &t_trackHitInMeasTIB);
2531  tree_->Branch("t_trackHitInMeasTID", "std::vector<int>", &t_trackHitInMeasTID);
2532  tree_->Branch("t_trackHitOutMeasTOB", "std::vector<int>", &t_trackHitOutMeasTOB);
2533  tree_->Branch("t_trackHitOutMeasTEC", "std::vector<int>", &t_trackHitOutMeasTEC);
2534  tree_->Branch("t_trackHitOutMeasTIB", "std::vector<int>", &t_trackHitOutMeasTIB);
2535  tree_->Branch("t_trackHitOutMeasTID", "std::vector<int>", &t_trackHitOutMeasTID);
2536  tree_->Branch("t_trackOutPosOutHitDr", "std::vector<double>", &t_trackOutPosOutHitDr);
2537  tree_->Branch("t_trackL", "std::vector<double>", &t_trackL);
2538 
2539  tree_->Branch("t_trackDxy", "std::vector<double>", &t_trackDxy);
2540  tree_->Branch("t_trackDxyBS", "std::vector<double>", &t_trackDxyBS);
2541  tree_->Branch("t_trackDz", "std::vector<double>", &t_trackDz);
2542  tree_->Branch("t_trackDzBS", "std::vector<double>", &t_trackDzBS);
2543  tree_->Branch("t_trackDxyPV", "std::vector<double>", &t_trackDxyPV);
2544  tree_->Branch("t_trackDzPV", "std::vector<double>", &t_trackDzPV);
2545  tree_->Branch("t_trackChiSq", "std::vector<double>", &t_trackChiSq);
2546  tree_->Branch("t_trackPVIdx", "std::vector<int>", &t_trackPVIdx);
2547 
2548  t_maxNearP31x31 = new std::vector<double>();
2549  t_maxNearP21x21 = new std::vector<double>();
2550 
2551  tree_->Branch("t_maxNearP31x31", "std::vector<double>", &t_maxNearP31x31);
2552  tree_->Branch("t_maxNearP21x21", "std::vector<double>", &t_maxNearP21x21);
2553 
2554  t_ecalSpike11x11 = new std::vector<int>();
2555  t_e7x7 = new std::vector<double>();
2556  t_e9x9 = new std::vector<double>();
2557  t_e11x11 = new std::vector<double>();
2558  t_e15x15 = new std::vector<double>();
2559 
2560  tree_->Branch("t_ecalSpike11x11", "std::vector<int>", &t_ecalSpike11x11);
2561  tree_->Branch("t_e7x7", "std::vector<double>", &t_e7x7);
2562  tree_->Branch("t_e9x9", "std::vector<double>", &t_e9x9);
2563  tree_->Branch("t_e11x11", "std::vector<double>", &t_e11x11);
2564  tree_->Branch("t_e15x15", "std::vector<double>", &t_e15x15);
2565 
2566  t_e7x7_10Sig = new std::vector<double>();
2567  t_e9x9_10Sig = new std::vector<double>();
2568  t_e11x11_10Sig = new std::vector<double>();
2569  t_e15x15_10Sig = new std::vector<double>();
2570  t_e7x7_15Sig = new std::vector<double>();
2571  t_e9x9_15Sig = new std::vector<double>();
2572  t_e11x11_15Sig = new std::vector<double>();
2573  t_e15x15_15Sig = new std::vector<double>();
2574  t_e7x7_20Sig = new std::vector<double>();
2575  t_e9x9_20Sig = new std::vector<double>();
2576  t_e11x11_20Sig = new std::vector<double>();
2577  t_e15x15_20Sig = new std::vector<double>();
2578  t_e7x7_25Sig = new std::vector<double>();
2579  t_e9x9_25Sig = new std::vector<double>();
2580  t_e11x11_25Sig = new std::vector<double>();
2581  t_e15x15_25Sig = new std::vector<double>();
2582  t_e7x7_30Sig = new std::vector<double>();
2583  t_e9x9_30Sig = new std::vector<double>();
2584  t_e11x11_30Sig = new std::vector<double>();
2585  t_e15x15_30Sig = new std::vector<double>();
2586 
2587  tree_->Branch("t_e7x7_10Sig", "std::vector<double>", &t_e7x7_10Sig);
2588  tree_->Branch("t_e9x9_10Sig", "std::vector<double>", &t_e9x9_10Sig);
2589  tree_->Branch("t_e11x11_10Sig", "std::vector<double>", &t_e11x11_10Sig);
2590  tree_->Branch("t_e15x15_10Sig", "std::vector<double>", &t_e15x15_10Sig);
2591  tree_->Branch("t_e7x7_15Sig", "std::vector<double>", &t_e7x7_15Sig);
2592  tree_->Branch("t_e9x9_15Sig", "std::vector<double>", &t_e9x9_15Sig);
2593  tree_->Branch("t_e11x11_15Sig", "std::vector<double>", &t_e11x11_15Sig);
2594  tree_->Branch("t_e15x15_15Sig", "std::vector<double>", &t_e15x15_15Sig);
2595  tree_->Branch("t_e7x7_20Sig", "std::vector<double>", &t_e7x7_20Sig);
2596  tree_->Branch("t_e9x9_20Sig", "std::vector<double>", &t_e9x9_20Sig);
2597  tree_->Branch("t_e11x11_20Sig", "std::vector<double>", &t_e11x11_20Sig);
2598  tree_->Branch("t_e15x15_20Sig", "std::vector<double>", &t_e15x15_20Sig);
2599  tree_->Branch("t_e7x7_25Sig", "std::vector<double>", &t_e7x7_25Sig);
2600  tree_->Branch("t_e9x9_25Sig", "std::vector<double>", &t_e9x9_25Sig);
2601  tree_->Branch("t_e11x11_25Sig", "std::vector<double>", &t_e11x11_25Sig);
2602  tree_->Branch("t_e15x15_25Sig", "std::vector<double>", &t_e15x15_25Sig);
2603  tree_->Branch("t_e7x7_30Sig", "std::vector<double>", &t_e7x7_30Sig);
2604  tree_->Branch("t_e9x9_30Sig", "std::vector<double>", &t_e9x9_30Sig);
2605  tree_->Branch("t_e11x11_30Sig", "std::vector<double>", &t_e11x11_30Sig);
2606  tree_->Branch("t_e15x15_30Sig", "std::vector<double>", &t_e15x15_30Sig);
2607 
2608  if (doMC_) {
2609  t_esim7x7 = new std::vector<double>();
2610  t_esim9x9 = new std::vector<double>();
2611  t_esim11x11 = new std::vector<double>();
2612  t_esim15x15 = new std::vector<double>();
2613 
2614  t_esim7x7Matched = new std::vector<double>();
2615  t_esim9x9Matched = new std::vector<double>();
2616  t_esim11x11Matched = new std::vector<double>();
2617  t_esim15x15Matched = new std::vector<double>();
2618 
2619  t_esim7x7Rest = new std::vector<double>();
2620  t_esim9x9Rest = new std::vector<double>();
2621  t_esim11x11Rest = new std::vector<double>();
2622  t_esim15x15Rest = new std::vector<double>();
2623 
2624  t_esim7x7Photon = new std::vector<double>();
2625  t_esim9x9Photon = new std::vector<double>();
2626  t_esim11x11Photon = new std::vector<double>();
2627  t_esim15x15Photon = new std::vector<double>();
2628 
2629  t_esim7x7NeutHad = new std::vector<double>();
2630  t_esim9x9NeutHad = new std::vector<double>();
2631  t_esim11x11NeutHad = new std::vector<double>();
2632  t_esim15x15NeutHad = new std::vector<double>();
2633 
2634  t_esim7x7CharHad = new std::vector<double>();
2635  t_esim9x9CharHad = new std::vector<double>();
2636  t_esim11x11CharHad = new std::vector<double>();
2637  t_esim15x15CharHad = new std::vector<double>();
2638 
2639  t_trkEcalEne = new std::vector<double>();
2640  t_simTrackP = new std::vector<double>();
2641  t_esimPdgId = new std::vector<double>();
2642 
2643  tree_->Branch("t_esim7x7", "std::vector<double>", &t_esim7x7);
2644  tree_->Branch("t_esim9x9", "std::vector<double>", &t_esim9x9);
2645  tree_->Branch("t_esim11x11", "std::vector<double>", &t_esim11x11);
2646  tree_->Branch("t_esim15x15", "std::vector<double>", &t_esim15x15);
2647 
2648  tree_->Branch("t_esim7x7Matched", "std::vector<double>", &t_esim7x7Matched);
2649  tree_->Branch("t_esim9x9Matched", "std::vector<double>", &t_esim9x9Matched);
2650  tree_->Branch("t_esim11x11Matched", "std::vector<double>", &t_esim11x11Matched);
2651  tree_->Branch("t_esim15x15Matched", "std::vector<double>", &t_esim15x15Matched);
2652 
2653  tree_->Branch("t_esim7x7Rest", "std::vector<double>", &t_esim7x7Rest);
2654  tree_->Branch("t_esim9x9Rest", "std::vector<double>", &t_esim9x9Rest);
2655  tree_->Branch("t_esim11x11Rest", "std::vector<double>", &t_esim11x11Rest);
2656  tree_->Branch("t_esim15x15Rest", "std::vector<double>", &t_esim15x15Rest);
2657 
2658  tree_->Branch("t_esim7x7Photon", "std::vector<double>", &t_esim7x7Photon);
2659  tree_->Branch("t_esim9x9Photon", "std::vector<double>", &t_esim9x9Photon);
2660  tree_->Branch("t_esim11x11Photon", "std::vector<double>", &t_esim11x11Photon);
2661  tree_->Branch("t_esim15x15Photon", "std::vector<double>", &t_esim15x15Photon);
2662 
2663  tree_->Branch("t_esim7x7NeutHad", "std::vector<double>", &t_esim7x7NeutHad);
2664  tree_->Branch("t_esim9x9NeutHad", "std::vector<double>", &t_esim9x9NeutHad);
2665  tree_->Branch("t_esim11x11NeutHad", "std::vector<double>", &t_esim11x11NeutHad);
2666  tree_->Branch("t_esim15x15NeutHad", "std::vector<double>", &t_esim15x15NeutHad);
2667 
2668  tree_->Branch("t_esim7x7CharHad", "std::vector<double>", &t_esim7x7CharHad);
2669  tree_->Branch("t_esim9x9CharHad", "std::vector<double>", &t_esim9x9CharHad);
2670  tree_->Branch("t_esim11x11CharHad", "std::vector<double>", &t_esim11x11CharHad);
2671  tree_->Branch("t_esim15x15CharHad", "std::vector<double>", &t_esim15x15CharHad);
2672 
2673  tree_->Branch("t_trkEcalEne", "std::vector<double>", &t_trkEcalEne);
2674  tree_->Branch("t_simTrackP", "std::vector<double>", &t_simTrackP);
2675  tree_->Branch("t_esimPdgId", "std::vector<double>", &t_esimPdgId);
2676  }
2677 
2678  t_maxNearHcalP3x3 = new std::vector<double>();
2679  t_maxNearHcalP5x5 = new std::vector<double>();
2680  t_maxNearHcalP7x7 = new std::vector<double>();
2681  t_h3x3 = new std::vector<double>();
2682  t_h5x5 = new std::vector<double>();
2683  t_h7x7 = new std::vector<double>();
2684  t_h3x3Sig = new std::vector<double>();
2685  t_h5x5Sig = new std::vector<double>();
2686  t_h7x7Sig = new std::vector<double>();
2687  t_infoHcal = new std::vector<int>();
2688 
2689  if (doMC_) {
2690  t_trkHcalEne = new std::vector<double>();
2691  t_hsim3x3 = new std::vector<double>();
2692  t_hsim5x5 = new std::vector<double>();
2693  t_hsim7x7 = new std::vector<double>();
2694  t_hsim3x3Matched = new std::vector<double>();
2695  t_hsim5x5Matched = new std::vector<double>();
2696  t_hsim7x7Matched = new std::vector<double>();
2697  t_hsim3x3Rest = new std::vector<double>();
2698  t_hsim5x5Rest = new std::vector<double>();
2699  t_hsim7x7Rest = new std::vector<double>();
2700  t_hsim3x3Photon = new std::vector<double>();
2701  t_hsim5x5Photon = new std::vector<double>();
2702  t_hsim7x7Photon = new std::vector<double>();
2703  t_hsim3x3NeutHad = new std::vector<double>();
2704  t_hsim5x5NeutHad = new std::vector<double>();
2705  t_hsim7x7NeutHad = new std::vector<double>();
2706  t_hsim3x3CharHad = new std::vector<double>();
2707  t_hsim5x5CharHad = new std::vector<double>();
2708  t_hsim7x7CharHad = new std::vector<double>();
2709  }
2710 
2711  tree_->Branch("t_maxNearHcalP3x3", "std::vector<double>", &t_maxNearHcalP3x3);
2712  tree_->Branch("t_maxNearHcalP5x5", "std::vector<double>", &t_maxNearHcalP5x5);
2713  tree_->Branch("t_maxNearHcalP7x7", "std::vector<double>", &t_maxNearHcalP7x7);
2714  tree_->Branch("t_h3x3", "std::vector<double>", &t_h3x3);
2715  tree_->Branch("t_h5x5", "std::vector<double>", &t_h5x5);
2716  tree_->Branch("t_h7x7", "std::vector<double>", &t_h7x7);
2717  tree_->Branch("t_h3x3Sig", "std::vector<double>", &t_h3x3Sig);
2718  tree_->Branch("t_h5x5Sig", "std::vector<double>", &t_h5x5Sig);
2719  tree_->Branch("t_h7x7Sig", "std::vector<double>", &t_h7x7Sig);
2720  tree_->Branch("t_infoHcal", "std::vector<int>", &t_infoHcal);
2721 
2722  if (doMC_) {
2723  tree_->Branch("t_trkHcalEne", "std::vector<double>", &t_trkHcalEne);
2724  tree_->Branch("t_hsim3x3", "std::vector<double>", &t_hsim3x3);
2725  tree_->Branch("t_hsim5x5", "std::vector<double>", &t_hsim5x5);
2726  tree_->Branch("t_hsim7x7", "std::vector<double>", &t_hsim7x7);
2727  tree_->Branch("t_hsim3x3Matched", "std::vector<double>", &t_hsim3x3Matched);
2728  tree_->Branch("t_hsim5x5Matched", "std::vector<double>", &t_hsim5x5Matched);
2729  tree_->Branch("t_hsim7x7Matched", "std::vector<double>", &t_hsim7x7Matched);
2730  tree_->Branch("t_hsim3x3Rest", "std::vector<double>", &t_hsim3x3Rest);
2731  tree_->Branch("t_hsim5x5Rest", "std::vector<double>", &t_hsim5x5Rest);
2732  tree_->Branch("t_hsim7x7Rest", "std::vector<double>", &t_hsim7x7Rest);
2733  tree_->Branch("t_hsim3x3Photon", "std::vector<double>", &t_hsim3x3Photon);
2734  tree_->Branch("t_hsim5x5Photon", "std::vector<double>", &t_hsim5x5Photon);
2735  tree_->Branch("t_hsim7x7Photon", "std::vector<double>", &t_hsim7x7Photon);
2736  tree_->Branch("t_hsim3x3NeutHad", "std::vector<double>", &t_hsim3x3NeutHad);
2737  tree_->Branch("t_hsim5x5NeutHad", "std::vector<double>", &t_hsim5x5NeutHad);
2738  tree_->Branch("t_hsim7x7NeutHad", "std::vector<double>", &t_hsim7x7NeutHad);
2739  tree_->Branch("t_hsim3x3CharHad", "std::vector<double>", &t_hsim3x3CharHad);
2740  tree_->Branch("t_hsim5x5CharHad", "std::vector<double>", &t_hsim5x5CharHad);
2741  tree_->Branch("t_hsim7x7CharHad", "std::vector<double>", &t_hsim7x7CharHad);
2742  }
2743  tree_->Branch("t_nTracks", &t_nTracks, "t_nTracks/I");
2744 }
2745 
2747  std::string theTrackQuality = "highPurity";
2748  reco::TrackBase::TrackQuality trackQuality_ = reco::TrackBase::qualityByName(theTrackQuality);
2749 
2750  edm::LogVerbatim("IsoTrack") << " Reference Point " << pTrack->referencePoint() << "\n TrackMmentum "
2751  << pTrack->momentum() << " (pt,eta,phi)(" << pTrack->pt() << "," << pTrack->eta() << ","
2752  << pTrack->phi() << ")"
2753  << " p " << pTrack->p() << "\n Normalized chi2 " << pTrack->normalizedChi2()
2754  << " charge " << pTrack->charge() << " qoverp() " << pTrack->qoverp() << "+-"
2755  << pTrack->qoverpError() << " d0 " << pTrack->d0() << "\n NValidHits "
2756  << pTrack->numberOfValidHits() << " NLostHits " << pTrack->numberOfLostHits()
2757  << " TrackQuality " << pTrack->qualityName(trackQuality_) << " "
2758  << pTrack->quality(trackQuality_);
2759 
2760  if (printTrkHitPattern_) {
2761  const reco::HitPattern &p = pTrack->hitPattern();
2762 
2763  edm::LogVerbatim("IsoTrack") << "default ";
2764  for (int i = 0; i < p.numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
2766  }
2767  edm::LogVerbatim("IsoTrack") << "trackerMissingInnerHits() ";
2768  for (int i = 0; i < p.numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS); i++) {
2770  }
2771  edm::LogVerbatim("IsoTrack") << "trackerMissingOuterHits() ";
2772  for (int i = 0; i < p.numberOfAllHits(reco::HitPattern::MISSING_OUTER_HITS); i++) {
2774  }
2775 
2776  edm::LogVerbatim("IsoTrack") << "\n \t trackerLayersWithMeasurement() " << p.trackerLayersWithMeasurement()
2777  << "\n \t pixelLayersWithMeasurement() " << p.pixelLayersWithMeasurement()
2778  << "\n \t stripLayersWithMeasurement() " << p.stripLayersWithMeasurement()
2779  << "\n \t pixelBarrelLayersWithMeasurement() " << p.pixelBarrelLayersWithMeasurement()
2780  << "\n \t pixelEndcapLayersWithMeasurement() " << p.pixelEndcapLayersWithMeasurement()
2781  << "\n \t stripTIBLayersWithMeasurement() " << p.stripTIBLayersWithMeasurement()
2782  << "\n \t stripTIDLayersWithMeasurement() " << p.stripTIDLayersWithMeasurement()
2783  << "\n \t stripTOBLayersWithMeasurement() " << p.stripTOBLayersWithMeasurement()
2784  << "\n \t stripTECLayersWithMeasurement() " << p.stripTECLayersWithMeasurement();
2785  }
2786 }
2787 
2788 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
std::vector< double > * t_hsim7x7CharHad
double qoverp() const
q / p
Definition: TrackBase.h:599
static const std::string kSharedResource
Definition: TFileService.h:76
double p() const
momentum vector magnitude
Definition: TrackBase.h:631
std::vector< double > * t_hsim3x3Matched
Log< level::Info, true > LogVerbatim
int stripTOBLayersWithMeasurement() const
Definition: HitPattern.cc:618
std::vector< double > * t_esim7x7CharHad
edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_L1extCenJet_
const std::string & gtTriggerMenuName() const
std::vector< double > * t_e15x15
value_type const * get() const
Definition: RefToBase.h:209
std::vector< double > * t_nTrksJetCalo
const bool L1TriggerAlgoInfo_
EventNumber_t event() const
Definition: EventID.h:40
std::vector< double > * t_maxNearHcalP7x7
const Point & referencePoint() const
Reference point on the track.
Definition: TrackBase.h:667
std::vector< double > * t_trackPt
std::vector< double > * t_trackPAll
const bool printTrkHitPattern_
std::vector< double > * t_e7x7
static constexpr size_t NEtaBins
std::vector< spr::propagatedTrackID > propagateCALO(edm::Handle< reco::TrackCollection > &trkCollection, const CaloGeometry *geo, const MagneticField *bField, const std::string &theTrackQuality, bool debug=false)
std::vector< double > * t_trackHcalPhi
std::vector< double > * t_trackDxyPVAll
~IsolatedTracksNxN() override
static std::string qualityName(TrackQuality)
Definition: TrackBase.h:572
std::vector< int > * t_trackHitOutMeasTEC
std::vector< double > * t_trackPdgIdAll
void eCaloSimInfo(std::vector< DetId > vdets, const CaloGeometry *geo, edm::Handle< T > &hitsEB, edm::Handle< T > &hitsEE, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, edm::SimTrackContainer::const_iterator trkInfo, caloSimInfo &info, double timeCut=150, bool debug=false)
std::vector< double > * t_trackOutPosOutHitDr
std::vector< double > * t_maxNearHcalP5x5
std::vector< double > * t_L1NonIsoEMPt
std::vector< int > * t_trackHitOutMeasTIB
int32_t *__restrict__ iv
std::vector< double > * t_esim7x7
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
std::vector< int > * t_trackHitInMissTOB
std::vector< double > * t_trackDz
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:611
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
std::vector< double > * t_trackEcalEta
int stripTIBLayersWithMeasurement() const
Definition: HitPattern.cc:598
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
std::vector< double > * t_trackEtaAll
std::vector< double > * t_esim11x11Matched
std::vector< double > * t_L1IsoEMEta
std::vector< double > * t_PVTracksSumPt
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::vector< double > * t_e9x9_10Sig
TH1F * h_maxNearP25x25[NPBins][NEtaBins]
TrackQuality
track quality
Definition: TrackBase.h:150
std::vector< double > * t_hsim3x3CharHad
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int stripTOBLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:680
std::vector< int > * t_PVndof
std::vector< std::pair< DetId, double > > eHCALmatrixCell(const HcalTopology *topology, const DetId &det, edm::Handle< T > &hits, int ieta, int iphi, bool includeHO=false, double hbThr=-100, double heThr=-100, double hfThr=-100, double hoThr=-100, bool debug=false)
std::map< std::string, L1GtAlgorithm > AlgorithmMap
map containing the algorithms
std::vector< double > * t_hsim7x7Rest
std::vector< double > * t_L1TauJetPhi
std::vector< double > * t_PVTracksSumPtWt
std::vector< int > * t_infoHcal
static constexpr size_t nL1BitsMax
std::vector< double > * t_e11x11_15Sig
std::vector< double > * t_L1TauJetPt
std::vector< int > * t_PVNTracksWt
std::vector< double > * t_trackEcalPhi
std::vector< double > * t_esim11x11NeutHad
std::vector< double > * t_L1TauJetEta
std::vector< double > * t_e15x15_15Sig
std::vector< double > * t_L1MuonEta
std::vector< double > * t_L1MuonPt
int bunchCrossing() const
Definition: EventBase.h:64
std::vector< double > * t_h3x3
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
std::vector< double > * t_L1MuonPhi
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:801
std::vector< double > * t_trackDxyAll
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > tok_ecalChStatus_
std::vector< double > * t_maxNearHcalP3x3
IsolatedTracksNxN(const edm::ParameterSet &)
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
std::vector< double > * t_hsim3x3Photon
edm::EDGetTokenT< reco::BeamSpot > tok_bs_
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
static constexpr size_t NPBins
int ii
Definition: cuy.py:589
std::vector< double > * t_PVTracksSumPtHP
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:62
std::vector< double > * t_esim11x11Rest
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
int pixelEndcapLayersWithMeasurement() const
Definition: HitPattern.cc:587
const double maxTrackEta_
std::vector< int > * t_trackNOuterHits
std::vector< double > * t_h7x7Sig
std::vector< double > * t_e15x15_30Sig
std::vector< double > * t_trackPhiAll
std::vector< double > * t_e7x7_25Sig
double genPartPBins[NPBins+1]
std::vector< double > * t_trkHcalEne
std::vector< double > * t_e7x7_20Sig
std::vector< double > * t_e11x11_10Sig
AlgorithmMap::const_iterator CItAlgo
iterators through map containing the algorithms
std::vector< double > * t_trackPtAll
std::vector< double > * t_h3x3Sig
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::vector< double > * t_esimPdgId
double eNeutralHad
Definition: CaloSimInfo.h:46
void printHitPattern(HitCategory category, int position, std::ostream &stream) const
Definition: HitPattern.cc:824
double eChargedHad
Definition: CaloSimInfo.h:47
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
std::vector< double > * t_e11x11_20Sig
std::vector< double > * t_jetEta
std::vector< double > * t_maxNearP31x31
std::vector< int > * t_trackHitInMissTEC
double chargeIsolationEcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, bool debug=false)
std::vector< double > * t_trkEcalEne
int iEvent
Definition: GenABIO.cc:224
std::vector< double > * t_trackEta
std::vector< double > * t_esim9x9Photon
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > tok_magField_
std::vector< int > * t_trackHitInMeasTID
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
std::vector< double > * t_hsim5x5Rest
std::vector< double > * t_e11x11
std::vector< double > * t_esim15x15
std::map< std::string, double > eHCALSimInfo(const edm::Event &, const HcalTopology *topology, const DetId &det, const CaloGeometry *geo, edm::Handle< T > &hits, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const reco::Track *pTrack, TrackerHitAssociator &associate, int ieta, int iphi, double timeCut=150, bool includeHO=false, bool debug=false)
edm::EDGetTokenT< edm::SimTrackContainer > tok_simTk_
std::vector< double > * t_L1CenJetPt
std::vector< double > * t_trackChiSqAll
std::vector< double > * t_e15x15_25Sig
std::vector< int > * t_trackHitOutMissTIB
edm::EDGetTokenT< edm::PCaloHitContainer > tok_caloHH_
std::vector< double > * t_esim15x15Rest
std::vector< double > * t_esim11x11
std::vector< double > * t_hsim7x7
std::vector< int > * t_PVNTracksHPWt
int stripTIDLayersWithMeasurement() const
Definition: HitPattern.cc:608
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::vector< int > * t_trackHitOutMissTID
std::vector< int > * t_trackHitOutMissTOBTEC
vector< PseudoJet > jets
double pt() const
track transverse momentum
Definition: TrackBase.h:637
std::vector< int > * t_trackHitsTEC
std::vector< double > * t_e9x9_30Sig
edm::EDGetTokenT< edm::PCaloHitContainer > tok_caloEB_
std::vector< double > * t_L1CenJetPhi
int numberOfAllHits(HitCategory category) const
Definition: HitPattern.h:804
std::vector< double > * t_L1FwdJetEta
const TString p1
Definition: fwPaths.cc:12
std::vector< double > * t_h7x7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< int > * t_ecalSpike11x11
std::vector< double > * t_esim15x15NeutHad
std::vector< double > * t_e9x9_20Sig
std::vector< int > * t_trackHitOutMeasTOB
std::vector< int > * t_trackHitInMeasTIB
edm::SimTrackContainer::const_iterator matchedSimTrack(const edm::Event &iEvent, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const reco::Track *pTrack, TrackerHitAssociator &associate, bool debug=false)
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
edm::EDGetTokenT< reco::TrackCollection > tok_genTrack_
edm::ESGetToken< CaloTopology, CaloTopologyRecord > tok_caloTopology_
T * make(const Args &...args) const
make new ROOT object
std::vector< double > * t_trackDzPV
std::vector< int > * t_PVNTracksHP
std::vector< double > * t_L1FwdJetPhi
edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_L1extFwdJet_
edm::EDGetTokenT< l1extra::L1EmParticleCollection > tok_L1extIsoEm_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< double > * t_nTrksJetVtx
bool isValid() const
Definition: HandleBase.h:70
std::vector< double > * t_esim7x7NeutHad
std::vector< double > * t_hsim3x3NeutHad
std::vector< double > * t_L1NonIsoEMPhi
std::vector< double > * t_trackDzBS
std::vector< int > * t_trackHitOutMissTOB
std::vector< double > * t_L1CenJetEta
std::vector< double > * t_trackDzPVAll
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > tok_L1extMu_
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > tok_topo_
std::vector< double > * t_esim7x7Matched
std::vector< double > * t_esim9x9
std::vector< double > * t_hsim5x5
std::vector< double > * t_h5x5
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:732
std::vector< double > * t_trackDxyPV
std::vector< double > * t_hsim5x5NeutHad
std::vector< double > * t_trackChiSq
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
std::vector< double > * t_hsim5x5Photon
std::vector< double > * t_PVy
std::vector< unsigned int > m_triggerMaskAlgoTrig
std::vector< int > * t_L1PreScale
std::vector< double > * t_hsim3x3Rest
int stripTIDLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:670
edm::ESGetToken< EcalTrigTowerConstituentsMap, IdealGeometryRecord > tok_htmap_
Definition: DetId.h:17
std::vector< int > * t_PVisValid
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
std::vector< double > * t_esim9x9Rest
std::vector< double > * t_L1METPhi
std::vector< int > * t_trackHitInMissTIB
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
std::vector< double > * t_e9x9_15Sig
std::vector< double > * t_e7x7_15Sig
std::vector< double > * t_trackPhi
T const * product() const
Definition: Handle.h:70
std::vector< double > * t_e7x7_10Sig
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
void endJob() override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const double pvTracksPtMin_
edm::EDGetTokenT< l1extra::L1EmParticleCollection > tok_L1extNoIsoEm_
std::vector< double > * t_trackDzAll
edm::EDGetTokenT< reco::VertexCollection > tok_recVtx_
TH1F * h_maxNearP31x31[NPBins][NEtaBins]
std::vector< double > * t_esim11x11CharHad
int stripLayersWithMeasurement() const
Definition: HitPattern.h:1005
int pixelBarrelLayersWithMeasurement() const
Definition: HitPattern.cc:576
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::PCaloHitContainer > tok_caloEE_
edm::EDGetTokenT< EcalRecHitCollection > tok_EB_
TH1F * h_maxNearP15x15[NPBins][NEtaBins]
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< double > * t_e9x9_25Sig
void eECALSimInfo(const edm::Event &, const DetId &det, const CaloGeometry *geo, const CaloTopology *caloTopology, edm::Handle< T > &hitsEB, edm::Handle< T > &hitsEE, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const reco::Track *pTrack, TrackerHitAssociator &associate, int ieta, int iphi, caloSimInfo &info, double timeCut=150, bool debug=false)
std::vector< double > * t_esim7x7Photon
static const bool useL1GtTriggerMenuLite(true)
std::vector< double > * t_esim9x9NeutHad
std::vector< double > * t_L1FwdJetPt
std::vector< double > * t_simTrackP
std::vector< double > * t_hsim7x7Matched
std::vector< int > * t_trackPVIdx
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:552
std::vector< int > * t_trackHitsTOB
std::vector< double > * t_esim9x9Matched
std::vector< double > * t_hsim7x7NeutHad
std::vector< double > * t_L1IsoEMPhi
std::vector< int > * t_PVNTracks
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
std::vector< double > * t_maxNearP21x21
TrackerHitAssociator::Config trackerHitAssociatorConfig_
edm::EventID id() const
Definition: EventBase.h:59
std::vector< double > * t_L1NonIsoEMEta
std::vector< double > * t_PVz
std::pair< math::XYZPoint, double > propagateTrackerEnd(const reco::Track *, const MagneticField *, bool debug=false)
void printTrack(const reco::Track *pTrack)
std::unique_ptr< L1GtUtils > m_l1GtUtils
std::vector< double > * t_e15x15_20Sig
std::vector< double > * t_PVTracksSumPtHPWt
std::vector< double > * t_trackDxy
std::vector< std::string > * t_L1AlgoNames
std::vector< double > * t_hsim5x5CharHad
std::vector< double > * t_e15x15_10Sig
const AlgorithmMap & gtAlgorithmMap() const
get / set the algorithm map (by name)
int stripTECLayersWithMeasurement() const
Definition: HitPattern.cc:628
TH1F * h_maxNearP21x21[NPBins][NEtaBins]
int stripTECLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:690
std::vector< double > * t_jetPhi
std::vector< double > * t_trackL
std::vector< double > * t_esim7x7Rest
std::vector< double > * t_trackDxyBS
tuple cout
Definition: gather_cfg.py:144
int charge() const
track electric charge
Definition: TrackBase.h:596
std::vector< double > * t_trackHcalEta
std::vector< double > * t_L1METPt
std::map< std::pair< unsigned int, std::string >, int > l1AlgoMap_
std::vector< int > * t_NLayersCrossed
edm::EDGetTokenT< reco::CaloJetCollection > tok_jets_
int weight
Definition: histoStyle.py:51
edm::ESGetToken< EcalSeverityLevelAlgo, EcalSeverityLevelAlgoRcd > tok_sevlv_
std::vector< int > * t_trackHitInMeasTOB
std::vector< double > * t_PVx
std::vector< int > * t_trackHitInMeasTEC
void beginJob() override
std::vector< double > * t_esim15x15Photon
std::vector< double > * t_e9x9
std::vector< double > * t_hsim5x5Matched
double genPartEtaBins[NEtaBins+1]
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
std::vector< double > * t_esim15x15Matched
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)
std::vector< int > * t_trackHitOutMeasTID
std::vector< double > * t_L1METEta
std::vector< double > * t_e7x7_30Sig
std::vector< int > * t_trackHitOutMissTEC
std::vector< double > * t_e11x11_25Sig
std::vector< double > * t_hsim3x3
std::vector< double > * t_trackP
std::vector< double > * t_hsim7x7Photon
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
std::vector< double > * t_jetPt
static const bool useL1EventSetup(true)
std::vector< double > * t_L1IsoEMPt
std::vector< double > * t_esim11x11Photon
int stripTIBLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:660
std::vector< double > * t_h5x5Sig
std::vector< int > * t_trackHitInMissTID
std::vector< double > * t_e11x11_30Sig
std::vector< int > * t_trackHitInMissTIBTID
std::vector< double > * t_esim9x9CharHad
double eECALmatrix(const DetId &detId, edm::Handle< T > &hitsEB, edm::Handle< T > &hitsEE, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, double ebThr=-100, double eeThr=-100, double tMin=-500, double tMax=500, bool debug=false)
edm::EDGetTokenT< EcalRecHitCollection > tok_EE_
double eHCALmatrix(const HcalTopology *topology, const DetId &det, edm::Handle< T > &hits, int ieta, int iphi, bool includeHO=false, bool algoNew=true, double hbThr=-100, double heThr=-100, double hfThr=-100, double hoThr=-100, double tMin=-500, double tMax=500, int useRaw=0, bool debug=false)
edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_L1extTauJet_
std::vector< double > * t_esim15x15CharHad
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)