CMS 3D CMS Logo

HcalIsoTrkAnalyzer.cc
Go to the documentation of this file.
1 // system include files
2 #include <atomic>
3 #include <cmath>
4 #include <memory>
5 #include <string>
6 #include <vector>
7 
8 // Root objects
9 #include "TROOT.h"
10 #include "TSystem.h"
11 #include "TFile.h"
12 #include "TProfile.h"
13 #include "TDirectory.h"
14 #include "TTree.h"
15 #include "TLorentzVector.h"
16 #include "TInterpreter.h"
17 
21 
22 //Tracks
27 // Vertices
31 
32 //Triggers
39 
44 
48 
49 //Generator information
53 
56 
64 
71 
84 
85 //#define EDM_ML_DEBUG
86 
87 class HcalIsoTrkAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns,edm::one::SharedResources> {
88 
89 public:
90  explicit HcalIsoTrkAnalyzer(edm::ParameterSet const&);
91  ~HcalIsoTrkAnalyzer() override {}
92 
93  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
94 
95 private:
96  void analyze(edm::Event const&, edm::EventSetup const&) override;
97  void beginJob() override;
98  void beginRun(edm::Run const&, edm::EventSetup const&) override;
99  void endRun(edm::Run const&, edm::EventSetup const&) override;
100 
101  std::array<int,3> fillTree(std::vector< math::XYZTLorentzVector>& vecL1,
102  std::vector< math::XYZTLorentzVector>& vecL3,
103  math::XYZPoint& leadPV,
104  std::vector<spr::propagatedTrackDirection>& trkCaloDirections,
105  std::vector<spr::propagatedTrackID>& trkCaloDets,
106  const CaloGeometry* geo, const CaloTopology* topo,
107  const HcalTopology* theHBHETopology,
108  const EcalChannelStatus* theEcalChStatus,
109  const EcalSeverityLevelAlgo* theEcalSevlv,
110  edm::Handle<EcalRecHitCollection>& barrelRecHitsHandle,
111  edm::Handle<EcalRecHitCollection>& endcapRecHitsHandle,
115  const HcalRespCorrs* respCorrs);
117  double trackP(const reco::Track* ,
119  double rhoh(const edm::Handle<CaloTowerCollection>&);
120  void storeEnergy(int indx, const HcalRespCorrs* respCorrs,
121  const std::vector<DetId>& ids,
122  std::vector<double>& edet, double & eHcal,
123  std::vector<unsigned int> *detIds,
124  std::vector<double> *hitEnergies);
125 
129  const std::vector<std::string> trigNames_;
137  const double pTrackLow_, pTrackHigh_;
139  const int useRaw_, dataType_, mode_;
143  const double hitEthrEE2_, hitEthrEE3_;
147  unsigned int nRun_, nLow_, nHigh_;
150  std::vector<double> etabins_, phibins_;
164 
165  TTree *tree, *tree2;
166  unsigned int t_RunNo, t_EventNo;
178  std::vector<unsigned int> *t_DetIds, *t_DetIds1, *t_DetIds3;
179  std::vector<double> *t_HitEnergies, *t_HitEnergies1, *t_HitEnergies3;
180  std::vector<bool> *t_trgbits, *t_hltbits;
181  bool t_L1Bit;
184  std::vector<int> *t_ietaAll, *t_ietaGood, *t_trackType;
185 };
186 
187 static const bool useL1EventSetup(false);
188 static const bool useL1GtTriggerMenuLite(true);
189 
191  trigNames_(iConfig.getParameter<std::vector<std::string> >("triggers")),
192  theTrackQuality_(iConfig.getParameter<std::string>("trackQuality")),
193  processName_(iConfig.getParameter<std::string>("processName")),
194  l1Filter_(iConfig.getParameter<std::string>("l1Filter")),
195  l2Filter_(iConfig.getParameter<std::string>("l2Filter")),
196  l3Filter_(iConfig.getParameter<std::string>("l3Filter")),
197  a_coneR_(iConfig.getParameter<double>("coneRadius")),
198  a_mipR_(iConfig.getParameter<double>("coneRadiusMIP")),
199  pTrackMin_(iConfig.getParameter<double>("minimumTrackP")),
200  eEcalMax_(iConfig.getParameter<double>("maximumEcalEnergy")),
201  maxRestrictionP_(iConfig.getParameter<double>("maxTrackP")),
202  slopeRestrictionP_(iConfig.getParameter<double>("slopeTrackP")),
203  hcalScale_(iConfig.getUntrackedParameter<double>("hHcalScale",1.0)),
204  eIsolate1_(iConfig.getParameter<double>("isolationEnergyTight")),
205  eIsolate2_(iConfig.getParameter<double>("isolationEnergyLoose")),
206  pTrackLow_(iConfig.getParameter<double>("momentumLow")),
207  pTrackHigh_(iConfig.getParameter<double>("momentumHigh")),
208  prescaleLow_(iConfig.getParameter<int>("prescaleLow")),
209  prescaleHigh_(iConfig.getParameter<int>("prescaleHigh")),
210  useRaw_(iConfig.getUntrackedParameter<int>("useRaw",0)),
211  dataType_(iConfig.getUntrackedParameter<int>("dataType",0)),
212  mode_(iConfig.getUntrackedParameter<int>("outMode",11)),
213  ignoreTrigger_(iConfig.getUntrackedParameter<bool>("ignoreTriggers",false)),
214  useL1Trigger_(iConfig.getUntrackedParameter<bool>("useL1Trigger",false)),
215  unCorrect_(iConfig.getUntrackedParameter<bool>("unCorrect",false)),
216  collapseDepth_(iConfig.getUntrackedParameter<bool>("collapseDepth",false)),
217  hitEthrEB_(iConfig.getParameter<double>("EBHitEnergyThreshold") ),
218  hitEthrEE0_(iConfig.getParameter<double>("EEHitEnergyThreshold0") ),
219  hitEthrEE1_(iConfig.getParameter<double>("EEHitEnergyThreshold1") ),
220  hitEthrEE2_(iConfig.getParameter<double>("EEHitEnergyThreshold2") ),
221  hitEthrEE3_(iConfig.getParameter<double>("EEHitEnergyThreshold3") ),
222  triggerEvent_(iConfig.getParameter<edm::InputTag>("labelTriggerEvent")),
223  theTriggerResultsLabel_(iConfig.getParameter<edm::InputTag>("labelTriggerResult")),
224  labelGenTrack_(iConfig.getParameter<std::string>("labelTrack")),
225  labelRecVtx_(iConfig.getParameter<std::string>("labelVertex")),
226  labelEB_(iConfig.getParameter<std::string>("labelEBRecHit")),
227  labelEE_(iConfig.getParameter<std::string>("labelEERecHit")),
228  labelHBHE_(iConfig.getParameter<std::string>("labelHBHERecHit")),
229  labelTower_(iConfig.getParameter<std::string>("labelCaloTower")),
230  l1TrigName_(iConfig.getUntrackedParameter<std::string>("l1TrigName","L1_SingleJet60")),
231  nRun_(0), nLow_(0), nHigh_(0), hdc_(nullptr) {
232 
233  usesResource(TFileService::kSharedResource);
234 
235  //now do whatever initialization is needed
236  const double isolationRadius(28.9), innerR(10.0), outerR(30.0);
238  selectionParameter_.minPt = iConfig.getParameter<double>("minTrackPt");;
239  selectionParameter_.minQuality = trackQuality_;
240  selectionParameter_.maxDxyPV = iConfig.getParameter<double>("maxDxyPV");
241  selectionParameter_.maxDzPV = iConfig.getParameter<double>("maxDzPV");
242  selectionParameter_.maxChi2 = iConfig.getParameter<double>("maxChi2");
243  selectionParameter_.maxDpOverP = iConfig.getParameter<double>("maxDpOverP");
244  selectionParameter_.minOuterHit = iConfig.getParameter<int>("minOuterHit");
245  selectionParameter_.minLayerCrossed = iConfig.getParameter<int>("minLayerCrossed");
246  selectionParameter_.maxInMiss = iConfig.getParameter<int>("maxInMiss");
247  selectionParameter_.maxOutMiss = iConfig.getParameter<int>("maxOutMiss");
248  a_charIsoR_ = a_coneR_ + isolationRadius;
249  a_coneR1_ = a_coneR_ + innerR;
250  a_coneR2_ = a_coneR_ + outerR;
251  // Different isolation cuts are described in DN-2016/029
252  // Tight cut uses 2 GeV; Loose cut uses 10 GeV
253  // Eta dependent cut uses (maxRestrictionP_ * exp(|ieta|*log(2.5)/18))
254  // with the factor for exponential slopeRestrictionP_ = log(2.5)/18
255  // maxRestrictionP_ = 8 GeV as came from a study
256  std::string labelBS = iConfig.getParameter<std::string>("labelBeamSpot");
257  std::string modnam = iConfig.getUntrackedParameter<std::string>("moduleName","");
258  std::string prdnam = iConfig.getUntrackedParameter<std::string>("producerName","");
259  edm::InputTag algTag = iConfig.getParameter<edm::InputTag>("algInputTag");
260  edm::InputTag extTag = iConfig.getParameter<edm::InputTag>("extInputTag");
261  l1GtUtils_ = new l1t::L1TGlobalUtil(iConfig, consumesCollector(), *this, algTag, extTag);
262  // define tokens for access
263  tok_trigEvt_ = consumes<trigger::TriggerEvent>(triggerEvent_);
264  tok_trigRes_ = consumes<edm::TriggerResults>(theTriggerResultsLabel_);
265  tok_bs_ = consumes<reco::BeamSpot>(labelBS);
266  tok_genTrack_ = consumes<reco::TrackCollection>(labelGenTrack_);
267  tok_ew_ = consumes<GenEventInfoProduct>(edm::InputTag("generator"));
268  tok_parts_ = consumes<reco::GenParticleCollection>(edm::InputTag("genParticles"));
269  tok_cala_ = consumes<CaloTowerCollection>(labelTower_);
270  tok_alg_ = consumes<BXVector<GlobalAlgBlk>>(algTag);
271 
272  if (modnam.empty()) {
273  tok_recVtx_ = consumes<reco::VertexCollection>(labelRecVtx_);
274  tok_EB_ = consumes<EcalRecHitCollection>(edm::InputTag("ecalRecHit",labelEB_));
275  tok_EE_ = consumes<EcalRecHitCollection>(edm::InputTag("ecalRecHit",labelEE_));
276  tok_hbhe_ = consumes<HBHERecHitCollection>(labelHBHE_);
277  edm::LogVerbatim("HcalIsoTrack") << "Labels used " << triggerEvent_ << " "
278  << theTriggerResultsLabel_ << " "
279  << labelBS << " " << labelRecVtx_ << " "
280  << labelGenTrack_ << " "
281  << edm::InputTag("ecalRecHit",labelEB_)
282  << " "
283  << edm::InputTag("ecalRecHit",labelEE_)
284  << " " << labelHBHE_ << " " << labelTower_;
285  } else {
286  tok_recVtx_ = consumes<reco::VertexCollection>(edm::InputTag(modnam,labelRecVtx_,prdnam));
287  tok_EB_ = consumes<EcalRecHitCollection>(edm::InputTag(modnam,labelEB_,prdnam));
288  tok_EE_ = consumes<EcalRecHitCollection>(edm::InputTag(modnam,labelEE_,prdnam));
289  tok_hbhe_ = consumes<HBHERecHitCollection>(edm::InputTag(modnam,labelHBHE_,prdnam));
290  edm::LogVerbatim("HcalIsoTrack") << "Labels used " << triggerEvent_ << " "
291  << theTriggerResultsLabel_ << " "
292  << labelBS << " "
293  << edm::InputTag(modnam,labelRecVtx_,prdnam)
294  << " " << labelGenTrack_ << " "
295  << edm::InputTag(modnam,labelEB_,prdnam)
296  << " "
297  << edm::InputTag(modnam,labelEE_,prdnam)
298  << " "
299  << edm::InputTag(modnam,labelHBHE_,prdnam)
300  << " " << labelTower_;
301  }
302 
303  edm::LogVerbatim("HcalIsoTrack") <<"Parameters read from config file \n"
304  <<"\t minPt " << selectionParameter_.minPt
305  <<"\t theTrackQuality " << theTrackQuality_
306  <<"\t minQuality " << selectionParameter_.minQuality
307  <<"\t maxDxyPV " << selectionParameter_.maxDxyPV
308  <<"\t maxDzPV " << selectionParameter_.maxDzPV
309  <<"\t maxChi2 " << selectionParameter_.maxChi2
310  <<"\t maxDpOverP " << selectionParameter_.maxDpOverP
311  <<"\t minOuterHit " << selectionParameter_.minOuterHit
312  <<"\t minLayerCrossed " << selectionParameter_.minLayerCrossed
313  <<"\t maxInMiss " << selectionParameter_.maxInMiss
314  <<"\t maxOutMiss " << selectionParameter_.maxOutMiss
315  <<"\t a_coneR " << a_coneR_
316  << ":" << a_coneR1_ << ":" << a_coneR2_
317  <<"\t a_charIsoR " << a_charIsoR_
318  <<"\t a_mipR " << a_mipR_
319  <<"\n pTrackMin_ " << pTrackMin_
320  <<"\t eEcalMax_ " << eEcalMax_
321  <<"\t maxRestrictionP_ "<< maxRestrictionP_
322  <<"\t slopeRestrictionP_ " << slopeRestrictionP_
323  <<"\t eIsolateStrong_ " << eIsolate1_
324  <<"\t eIsolateSoft_ " << eIsolate2_
325  <<"\t hcalScale_ " << hcalScale_
326  <<"\n\t momentumLow_ " << pTrackLow_
327  <<"\t prescaleLow_ " << prescaleLow_
328  <<"\t momentumHigh_ " << pTrackHigh_
329  <<"\t prescaleHigh_ " << prescaleHigh_
330  <<"\n\t useRaw_ " << useRaw_
331  <<"\t ignoreTrigger_ " << ignoreTrigger_
332  <<"\n\t useL1Trigegr_ " << useL1Trigger_
333  <<"\t dataType_ " << dataType_
334  <<"\t mode_ " << mode_
335  <<"\t unCorrect_ " << unCorrect_
336  <<"\t collapseDepth_ " << collapseDepth_
337  <<"\t L1TrigName_ " << l1TrigName_;
338  edm::LogVerbatim("HcalIsoTrack") << "Process " << processName_
339  << " L1Filter:" << l1Filter_ << " L2Filter:"
340  << l2Filter_ << " L3Filter:" << l3Filter_;
341  for (unsigned int k=0; k<trigNames_.size(); ++k) {
342  edm::LogVerbatim("HcalIsoTrack") << "Trigger[" << k << "] "
343  << trigNames_[k];
344  }
345 
346  for (int i=0; i<10;i++) phibins_.push_back(-M_PI+0.1*(2*i+1)*M_PI);
347  for (int i=0; i<8; ++i) etabins_.push_back(-2.1+0.6*i);
348  etadist_ = etabins_[1]-etabins_[0];
349  phidist_ = phibins_[1]-phibins_[0];
350  etahalfdist_ = 0.5*etadist_;
351  phihalfdist_ = 0.5*phidist_;
352  edm::LogVerbatim("HcalIsoTrack") << "EtaDist " << etadist_ << " "
353  << etahalfdist_ << " PhiDist " << phidist_
354  << " " <<phihalfdist_;
355  unsigned int k1(0), k2(0);
356  for (auto phi : phibins_) {
357  edm::LogVerbatim("HcalIsoTrack") << "phibin_[" << k1 << "] " << phi; ++k1;
358  }
359  for (auto eta : etabins_) {
360  edm::LogVerbatim("HcalIsoTrack") << "etabin_[" << k2 << "] " << eta; ++k2;
361  }
362 }
363 
365 
366  t_Run = iEvent.id().run();
367  t_Event = iEvent.id().event();
369 #ifdef EDM_ML_DEBUG
370  edm::LogVerbatim("HcalIsoTrack") << "Run " << t_Run << " Event " << t_Event
371  << " type " << t_DataType << " Luminosity "
372  << iEvent.luminosityBlock() << " Bunch "
373  << iEvent.bunchCrossing();
374 #endif
375  //Get magnetic field and ECAL channel status
377  iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
378  const MagneticField *bField = bFieldH.product();
379 
381  iSetup.get<EcalChannelStatusRcd>().get(ecalChStatus);
382  const EcalChannelStatus* theEcalChStatus = ecalChStatus.product();
383 
385  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
386  const EcalSeverityLevelAlgo* theEcalSevlv = sevlv.product();
387 
388  // get handles to calogeometry and calotopology
390  iSetup.get<CaloGeometryRecord>().get(pG);
391  const CaloGeometry* geo = pG.product();
392 
393  edm::ESHandle<CaloTopology> theCaloTopology;
394  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
395  const CaloTopology *caloTopology = theCaloTopology.product();
396 
398  iSetup.get<HcalRecNumberingRecord>().get(htopo);
399  const HcalTopology* theHBHETopology = htopo.product();
400 
402  iSetup.get<HcalRespCorrsRcd>().get(resp);
403  HcalRespCorrs* respCorrs = new HcalRespCorrs(*resp.product());
404  respCorrs->setTopo(theHBHETopology);
405 
406  //=== genParticle information
408  iEvent.getByToken(tok_parts_, genParticles);
409 
410  bool okC(true);
411  //Get track collection
413  iEvent.getByToken(tok_genTrack_, trkCollection);
414  reco::TrackCollection::const_iterator trkItr;
415  if (!trkCollection.isValid()) {
416  edm::LogWarning("HcalIsoTrack") << "Cannot access the collection "
417  << labelGenTrack_;
418  okC = false;
419  }
420 
421  //event weight for FLAT sample
422  t_EventWeight = 1.0;
424  iEvent.getByToken(tok_ew_, genEventInfo);
425  if (genEventInfo.isValid()) t_EventWeight = genEventInfo->weight();
426 
427  //Define the best vertex and the beamspot
429  iEvent.getByToken(tok_recVtx_, recVtxs);
430  edm::Handle<reco::BeamSpot> beamSpotH;
431  iEvent.getByToken(tok_bs_, beamSpotH);
432  math::XYZPoint leadPV(0,0,0);
433  t_goodPV = t_nVtx = 0;
434  if (recVtxs.isValid() && !(recVtxs->empty())) {
435  t_nVtx = recVtxs->size();
436  for (unsigned int k=0; k<recVtxs->size(); ++k) {
437  if (!((*recVtxs)[k].isFake()) && ((*recVtxs)[k].ndof() > 4)) {
438  if (t_goodPV == 0) leadPV = math::XYZPoint((*recVtxs)[k].x(),(*recVtxs)[k].y(),(*recVtxs)[k].z());
439  t_goodPV++;
440  }
441  }
442  }
443  if (t_goodPV == 0 && beamSpotH.isValid()) {
444  leadPV = beamSpotH->position();
445  }
447 #ifdef EDM_ML_DEBUG
448  edm::LogVerbatim("HcalIsoTrack") << "Primary Vertex " << leadPV << " out of "
449  << t_goodPV << " vertex";
450  if (beamSpotH.isValid()) {
451  edm::LogVerbatim("HcalIsoTrack") << " Beam Spot " << beamSpotH->position();
452  }
453 #endif
454  // RecHits
455  edm::Handle<EcalRecHitCollection> barrelRecHitsHandle;
456  iEvent.getByToken(tok_EB_, barrelRecHitsHandle);
457  if (!barrelRecHitsHandle.isValid()) {
458  edm::LogWarning("HcalIsoTrack") << "Cannot access the collection "
459  << labelEB_;
460  okC = false;
461  }
462  edm::Handle<EcalRecHitCollection> endcapRecHitsHandle;
463  iEvent.getByToken(tok_EE_, endcapRecHitsHandle);
464  if (!endcapRecHitsHandle.isValid()) {
465  edm::LogWarning("HcalIsoTrack") << "Cannot access the collection "
466  << labelEE_;
467  okC = false;
468  }
470  iEvent.getByToken(tok_hbhe_, hbhe);
471  if (!hbhe.isValid()) {
472  edm::LogWarning("HcalIsoTrack") << "Cannot access the collection "
473  << labelHBHE_;
474  okC = false;
475  }
477  iEvent.getByToken(tok_cala_, caloTower);
478 
479  //Propagate tracks to calorimeter surface)
480  std::vector<spr::propagatedTrackDirection> trkCaloDirections;
481  spr::propagateCALO(trkCollection, geo, bField, theTrackQuality_,
482  trkCaloDirections, false);
483  std::vector<spr::propagatedTrackID> trkCaloDets;
484  spr::propagateCALO(trkCollection, geo, bField, theTrackQuality_,
485  trkCaloDets, false);
486  std::vector<math::XYZTLorentzVector> vecL1, vecL3;
487  t_RunNo = iEvent.id().run();
488  t_EventNo = iEvent.id().event();
489  t_Tracks = trkCollection->size();
490  t_TracksProp = trkCaloDirections.size();
491  t_ietaAll->clear(); t_ietaGood->clear(); t_trackType->clear();
492  t_trgbits->clear(); t_hltbits->clear();
493 #ifdef EDM_ML_DEBUG
494  edm::LogVerbatim("HcalIsoTrack") << "# of propagated tracks " << t_TracksProp
495  << " out of " << t_Tracks << " with Trigger "
496  << ignoreTrigger_;
497 #endif
498 
499  //Trigger
500  t_trgbits->assign(trigNames_.size(),false);
501  t_hltbits->assign(trigNames_.size(),false);
503  t_L1Bit = true;
504  t_TrigPass = false;
505 
506  //L1
507  l1GtUtils_->retrieveL1(iEvent,iSetup,tok_alg_);
508  const std::vector<std::pair<std::string, bool> > & finalDecisions = l1GtUtils_->decisionsFinal();
509  for (const auto& decision : finalDecisions) {
510  if (decision.first.find(l1TrigName_) != std::string::npos) {
511  t_L1Bit = decision.second;
512  break;
513  }
514  }
515 #ifdef EDM_ML_DEBUG
516  edm::LogVerbatim("HcalIsoTrack") << "Trigger Information for " << l1TrigName_
517  << " is " << t_L1Bit << " from a list of "
518  << finalDecisions.size() << " decisions";
519 #endif
520 
521  //HLT
523  iEvent.getByToken(tok_trigRes_, triggerResults);
524  if (triggerResults.isValid()) {
525  const edm::TriggerNames & triggerNames = iEvent.triggerNames(*triggerResults);
526  const std::vector<std::string> & names = triggerNames.triggerNames();
527  if (!trigNames_.empty()) {
528  for (unsigned int iHLT=0; iHLT<triggerResults->size(); iHLT++) {
529  int hlt = triggerResults->accept(iHLT);
530  for (unsigned int i=0; i<trigNames_.size(); ++i) {
531  if (names[iHLT].find(trigNames_[i]) != std::string::npos) {
532  t_trgbits->at(i) = (hlt>0);
533  t_hltbits->at(i) = (hlt>0);
534  if (hlt>0) t_TrigPass = true;
535 #ifdef EDM_ML_DEBUG
536  edm::LogVerbatim("HcalIsoTrack") << "This trigger "
537  << names[iHLT] << " Flag "
538  << hlt << ":" << t_trgbits->at(i);
539 #endif
540  }
541  }
542  }
543  }
544  }
545 #ifdef EDM_ML_DEBUG
546  edm::LogVerbatim("HcalIsoTrack") << "HLT Information shows " << t_TrigPass
547  << ":" << trigNames_.empty() << ":" << okC;
548 #endif
549 
550  std::array<int,3> ntksave{ {0,0,0} };
551  if (ignoreTrigger_ || useL1Trigger_) {
552  t_l1pt = t_l1eta = t_l1phi = 0;
553  t_l3pt = t_l3eta = t_l3phi = 0;
554  if (ignoreTrigger_ || t_L1Bit)
555  ntksave = fillTree(vecL1, vecL3, leadPV, trkCaloDirections, trkCaloDets,
556  geo, caloTopology, theHBHETopology, theEcalChStatus,
557  theEcalSevlv, barrelRecHitsHandle, endcapRecHitsHandle,
558  hbhe, caloTower, genParticles, respCorrs);
559  t_TracksSaved = ntksave[0];
560  t_TracksLoose = ntksave[1];
561  t_TracksTight = ntksave[2];
562  } else {
563  trigger::TriggerEvent triggerEvent;
564  edm::Handle<trigger::TriggerEvent> triggerEventHandle;
565  iEvent.getByToken(tok_trigEvt_, triggerEventHandle);
566  if (!triggerEventHandle.isValid()) {
567  edm::LogWarning("HcalIsoTrack") << "Error! Can't get the product "
568  << triggerEvent_.label() ;
569  } else if (okC) {
570  triggerEvent = *(triggerEventHandle.product());
571  const trigger::TriggerObjectCollection& TOC(triggerEvent.getObjects());
572  bool done(false);
573  if (triggerResults.isValid()) {
574  std::vector<std::string> modules;
575  const edm::TriggerNames & triggerNames = iEvent.triggerNames(*triggerResults);
576  const std::vector<std::string> & names = triggerNames.triggerNames();
577  for (unsigned int iHLT=0; iHLT<triggerResults->size(); iHLT++) {
578  bool ok = (t_TrigPass) || (trigNames_.empty());
579  if (ok) {
580  unsigned int triggerindx = hltConfig_.triggerIndex(names[iHLT]);
581  const std::vector<std::string>& moduleLabels(hltConfig_.moduleLabels(triggerindx));
582  std::vector<math::XYZTLorentzVector> vecL2;
583  vecL1.clear(); vecL3.clear();
584  //loop over all trigger filters in event (i.e. filters passed)
585  for (unsigned int ifilter=0; ifilter<triggerEvent.sizeFilters();
586  ++ifilter) {
587  std::vector<int> Keys;
588  std::string label = triggerEvent.filterTag(ifilter).label();
589  //loop over keys to objects passing this filter
590  for (unsigned int imodule=0; imodule<moduleLabels.size(); imodule++) {
591  if (label.find(moduleLabels[imodule]) != std::string::npos) {
592 #ifdef EDM_ML_DEBUG
593  edm::LogVerbatim("HcalIsoTrack") << "FilterName " << label;
594 #endif
595  for (unsigned int ifiltrKey=0; ifiltrKey<triggerEvent.filterKeys(ifilter).size(); ++ifiltrKey) {
596  Keys.push_back(triggerEvent.filterKeys(ifilter)[ifiltrKey]);
597  const trigger::TriggerObject& TO(TOC[Keys[ifiltrKey]]);
598  math::XYZTLorentzVector v4(TO.px(), TO.py(), TO.pz(), TO.energy());
599  if (label.find(l2Filter_) != std::string::npos) {
600  vecL2.push_back(v4);
601  } else if (label.find(l3Filter_) != std::string::npos) {
602  vecL3.push_back(v4);
603  } else if ((label.find(l1Filter_) != std::string::npos) ||
604  (l1Filter_.empty())) {
605  vecL1.push_back(v4);
606  }
607 #ifdef EDM_ML_DEBUG
608  edm::LogVerbatim("HcalIsoTrack") << "key " << ifiltrKey
609  << " : pt " << TO.pt()
610  << " eta " << TO.eta()
611  << " phi " << TO.phi()
612  << " mass " << TO.mass()
613  <<" Id " << TO.id();
614 #endif
615  }
616 #ifdef EDM_ML_DEBUG
617  edm::LogVerbatim("HcalIsoTrack") << "sizes " << vecL1.size()
618  <<":" << vecL2.size() << ":"
619  << vecL3.size();
620 #endif
621  }
622  }
623  }
625  math::XYZTLorentzVector mindRvec1;
626  double mindR1(999);
627  for (unsigned int i=0; i<vecL2.size(); i++) {
628  double dr = dR(vecL1[0],vecL2[i]);
629 #ifdef EDM_ML_DEBUG
630  edm::LogVerbatim("HcalIsoTrack") << "lvl2[" << i << "] dR " << dr;
631 #endif
632  if (dr<mindR1) {
633  mindR1 = dr;
634  mindRvec1 = vecL2[i];
635  }
636  }
637 #ifdef EDM_ML_DEBUG
638  edm::LogVerbatim("HcalIsoTrack") << "L2 object closest to L1 "
639  << mindRvec1 << " at Dr "
640  << mindR1;
641 #endif
642 
643  if (!vecL1.empty()) {
644  t_l1pt = vecL1[0].pt();
645  t_l1eta = vecL1[0].eta();
646  t_l1phi = vecL1[0].phi();
647  } else {
648  t_l1pt = t_l1eta = t_l1phi = 0;
649  }
650  if (!vecL3.empty()) {
651  t_l3pt = vecL3[0].pt();
652  t_l3eta = vecL3[0].eta();
653  t_l3phi = vecL3[0].phi();
654  } else {
655  t_l3pt = t_l3eta = t_l3phi = 0;
656  }
657  // Now fill in the tree for each selected track
658  if (!done) {
659  ntksave = fillTree(vecL1, vecL3, leadPV, trkCaloDirections,
660  trkCaloDets, geo, caloTopology,theHBHETopology,
661  theEcalChStatus, theEcalSevlv,
662  barrelRecHitsHandle, endcapRecHitsHandle,
663  hbhe, caloTower, genParticles, respCorrs);
664  t_TracksSaved += ntksave[0];
665  t_TracksLoose += ntksave[1];
666  t_TracksTight += ntksave[2];
667  done = true;
668  }
669  }
670  }
671  }
672  }
673  }
674 #ifdef EDM_ML_DEBUG
675  edm::LogVerbatim("HcalIsoTrack") << "Final results on selected tracks "
676  << t_TracksSaved << ":" << t_TracksLoose
677  << ":" << t_TracksTight;
678 #endif
680  tree2->Fill();
681 }
682 
684 
685  tree = fs->make<TTree>("CalibTree", "CalibTree");
686 
687  tree->Branch("t_Run", &t_Run, "t_Run/I");
688  tree->Branch("t_Event", &t_Event, "t_Event/I");
689  tree->Branch("t_DataType", &t_DataType, "t_DataType/I");
690  tree->Branch("t_ieta", &t_ieta, "t_ieta/I");
691  tree->Branch("t_iphi", &t_iphi, "t_iphi/I");
692  tree->Branch("t_EventWeight", &t_EventWeight, "t_EventWeight/D");
693  tree->Branch("t_nVtx", &t_nVtx, "t_nVtx/I");
694  tree->Branch("t_nTrk", &t_nTrk, "t_nTrk/I");
695  tree->Branch("t_goodPV", &t_goodPV, "t_goodPV/I");
696  if (((mode_/10)%10) == 1) {
697  tree->Branch("t_l1pt", &t_l1pt, "t_l1pt/D");
698  tree->Branch("t_l1eta", &t_l1eta, "t_l1eta/D");
699  tree->Branch("t_l1phi", &t_l1phi, "t_l1phi/D");
700  tree->Branch("t_l3pt", &t_l3pt, "t_l3pt/D");
701  tree->Branch("t_l3eta", &t_l3eta, "t_l3eta/D");
702  tree->Branch("t_l3phi", &t_l3phi, "t_l3phi/D");
703  }
704  tree->Branch("t_p", &t_p, "t_p/D");
705  tree->Branch("t_pt", &t_pt, "t_pt/D");
706  tree->Branch("t_phi", &t_phi, "t_phi/D");
707  tree->Branch("t_mindR1", &t_mindR1, "t_mindR1/D");
708  tree->Branch("t_mindR2", &t_mindR2, "t_mindR2/D");
709  tree->Branch("t_eMipDR", &t_eMipDR, "t_eMipDR/D");
710  tree->Branch("t_eHcal", &t_eHcal, "t_eHcal/D");
711  tree->Branch("t_eHcal10", &t_eHcal10, "t_eHcal10/D");
712  tree->Branch("t_eHcal30", &t_eHcal30, "t_eHcal30/D");
713  tree->Branch("t_hmaxNearP", &t_hmaxNearP, "t_hmaxNearP/D");
714  tree->Branch("t_emaxNearP", &t_emaxNearP, "t_emaxNearP/D");
715  tree->Branch("t_eAnnular", &t_eAnnular, "t_eAnnular/D");
716  tree->Branch("t_hAnnular", &t_hAnnular, "t_hAnnular/D");
717  tree->Branch("t_rhoh", &t_rhoh, "t_rhoh/D");
718  tree->Branch("t_selectTk", &t_selectTk, "t_selectTk/O");
719  tree->Branch("t_qltyFlag", &t_qltyFlag, "t_qltyFlag/O");
720  tree->Branch("t_qltyMissFlag",&t_qltyMissFlag,"t_qltyMissFlag/O");
721  tree->Branch("t_qltyPVFlag", &t_qltyPVFlag, "t_qltyPVFlag/O");
722  tree->Branch("t_gentrackP", &t_gentrackP, "t_gentrackP/D");
723 
724  t_DetIds = new std::vector<unsigned int>();
725  t_DetIds1 = new std::vector<unsigned int>();
726  t_DetIds3 = new std::vector<unsigned int>();
727  t_HitEnergies = new std::vector<double>();
728  t_HitEnergies1 = new std::vector<double>();
729  t_HitEnergies3 = new std::vector<double>();
730  t_trgbits = new std::vector<bool>();
731  tree->Branch("t_DetIds", "std::vector<unsigned int>", &t_DetIds);
732  tree->Branch("t_HitEnergies", "std::vector<double>", &t_HitEnergies);
733  if (((mode_/10)%10) == 1) {
734  tree->Branch("t_trgbits", "std::vector<bool>", &t_trgbits);
735  }
736  if ((mode_%10) == 1) {
737  tree->Branch("t_DetIds1", "std::vector<unsigned int>", &t_DetIds1);
738  tree->Branch("t_DetIds3", "std::vector<unsigned int>", &t_DetIds3);
739  tree->Branch("t_HitEnergies1", "std::vector<double>", &t_HitEnergies1);
740  tree->Branch("t_HitEnergies3", "std::vector<double>", &t_HitEnergies3);
741  }
742  tree2 = fs->make<TTree>("EventInfo", "Event Information");
743 
744  tree2->Branch("t_RunNo", &t_RunNo, "t_RunNo/i");
745  tree2->Branch("t_EventNo", &t_EventNo, "t_EventNo/i");
746  tree2->Branch("t_Tracks", &t_Tracks, "t_Tracks/I");
747  tree2->Branch("t_TracksProp", &t_TracksProp, "t_TracksProp/I");
748  tree2->Branch("t_TracksSaved", &t_TracksSaved, "t_TracksSaved/I");
749  tree2->Branch("t_TracksLoose", &t_TracksLoose, "t_TracksLoose/I");
750  tree2->Branch("t_TracksTight", &t_TracksTight, "t_TracksTight/I");
751  tree2->Branch("t_TrigPass", &t_TrigPass, "t_TrigPass/O");
752  tree2->Branch("t_TrigPassSel", &t_TrigPassSel, "t_TrigPassSel/O");
753  tree2->Branch("t_L1Bit", &t_L1Bit, "t_L1Bit/O");
754  tree2->Branch("t_allvertex", &t_allvertex, "t_allvertex/I");
755  t_hltbits = new std::vector<bool>();
756  t_ietaAll = new std::vector<int>();
757  t_ietaGood = new std::vector<int>();
758  t_trackType = new std::vector<int>();
759  tree2->Branch("t_ietaAll", "std::vector<int>", &t_ietaAll);
760  tree2->Branch("t_ietaGood", "std::vector<int>", &t_ietaGood);
761  tree2->Branch("t_trackType", "std::vector<int>", &t_trackType);
762  tree2->Branch("t_hltbits", "std::vector<bool>", &t_hltbits);
763 }
764 
765 
766 // ------------ method called when starting to processes a run ------------
767 void HcalIsoTrkAnalyzer::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
769  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
770  hdc_ = pHRNDC.product();
771 
772  bool changed_(true);
773  bool flag = hltConfig_.init(iRun,iSetup,processName_,changed_);
774  edm::LogVerbatim("HcalIsoTrack") << "Run[" << nRun_ << "] " << iRun.run()
775  << " process " << processName_
776  << " init flag " << flag << " change flag "
777  << changed_;
778  // check if trigger names in (new) config
779  if (changed_) {
780  changed_ = false;
781 #ifdef EDM_ML_DEBUG
782  edm::LogVerbatim("HcalIsoTrack") << "New trigger menu found !!!";
783 #endif
784  const unsigned int n(hltConfig_.size());
785  for (unsigned itrig=0; itrig<trigNames_.size(); itrig++) {
786  unsigned int triggerindx = hltConfig_.triggerIndex(trigNames_[itrig]);
787  if (triggerindx >= n) {
788  edm::LogWarning("HcalIsoTrack") << trigNames_[itrig] << " "
789  << triggerindx << " does not exist in "
790  << "the current menu";
791 #ifdef EDM_ML_DEBUG
792  } else {
793  edm::LogVerbatim("HcalIsoTrack") << trigNames_[itrig] << " "
794  << triggerindx << " exists";
795 #endif
796  }
797  }
798  }
799 }
800 
801 // ------------ method called when ending the processing of a run ------------
803  nRun_++;
804  edm::LogVerbatim("HcalIsoTrack") << "endRun[" << nRun_ << "] " << iRun.run();
805 }
806 
809  std::vector<std::string> trig = {"HLT_PFJet40","HLT_PFJet60","HLT_PFJet80",
810  "HLT_PFJet140","HLT_PFJet200","HLT_PFJet260",
811  "HLT_PFJet320","HLT_PFJet400","HLT_PFJet450",
812  "HLT_PFJet500"};
813  desc.add<std::vector<std::string>>("triggers",trig);
814  desc.add<std::string>("processName","HLT");
815  desc.add<std::string>("l1Filter","");
816  desc.add<std::string>("l2Filter","L2Filter");
817  desc.add<std::string>("l3Filter","Filter");
818  // following 10 parameters are parameters to select good tracks
819  desc.add<std::string>("trackQuality","highPurity");
820  desc.add<double>("minTrackPt",1.0);
821  desc.add<double>("maxDxyPV",0.02);
822  desc.add<double>("maxDzPV",0.02);
823  desc.add<double>("maxChi2",5.0);
824  desc.add<double>("maxDpOverP",0.1);
825  desc.add<int>("minOuterHit",4);
826  desc.add<int>("minLayerCrossed",8);
827  desc.add<int>("maxInMiss",0);
828  desc.add<int>("maxOutMiss",0);
829  // Minimum momentum of selected isolated track and signal zone
830  desc.add<double>("minimumTrackP",20.0);
831  desc.add<double>("coneRadius",34.98);
832  // signal zone in ECAL and MIP energy cutoff
833  desc.add<double>("coneRadiusMIP",14.0);
834  desc.add<double>("maximumEcalEnergy",2.0);
835  // following 4 parameters are for isolation cuts and described in the code
836  desc.add<double>("maxTrackP",8.0);
837  desc.add<double>("slopeTrackP",0.05090504066);
838  desc.add<double>("isolationEnergyTight",2.0);
839  desc.add<double>("isolationEnergyLoose",10.0);
840  // energy thershold for ECAL (from Egamma group)
841  desc.add<double>("EBHitEnergyThreshold",0.10);
842  desc.add<double>("EEHitEnergyThreshold0",-41.0664);
843  desc.add<double>("EEHitEnergyThreshold1",68.7950);
844  desc.add<double>("EEHitEnergyThreshold2",-38.1483);
845  desc.add<double>("EEHitEnergyThreshold3",7.04303);
846  // prescale factors
847  desc.add<double>("momentumLow", 40.0);
848  desc.add<double>("momentumHigh",60.0);
849  desc.add<int>("prescaleLow", 1);
850  desc.add<int>("prescaleHigh",1);
851  // various labels for collections used in the code
852  desc.add<edm::InputTag>("labelTriggerEvent",edm::InputTag("hltTriggerSummaryAOD","","HLT"));
853  desc.add<edm::InputTag>("labelTriggerResult",edm::InputTag("TriggerResults","","HLT"));
854  desc.add<std::string>("labelTrack","generalTracks");
855  desc.add<std::string>("labelVertex","offlinePrimaryVertices");
856  desc.add<std::string>("labelEBRecHit","EcalRecHitsEB");
857  desc.add<std::string>("labelEERecHit","EcalRecHitsEE");
858  desc.add<std::string>("labelHBHERecHit","hbhereco");
859  desc.add<std::string>("labelBeamSpot","offlineBeamSpot");
860  desc.add<std::string>("labelCaloTower","towerMaker");
861  desc.add<edm::InputTag>("algInputTag", edm::InputTag("gtStage2Digis"));
862  desc.add<edm::InputTag>("extInputTag", edm::InputTag("gtStage2Digis"));
863  desc.addUntracked<std::string>("moduleName","");
864  desc.addUntracked<std::string>("producerName","");
865  // Various flags used for selecting tracks, choice of energy Method2/0
866  // Data type 0/1 for single jet trigger or others
867  desc.addUntracked<int>("useRaw",0);
868  desc.addUntracked<bool>("ignoreTriggers",false);
869  desc.addUntracked<bool>("useL1Trigger",false);
870  desc.addUntracked<double>("hcalScale",1.0);
871  desc.addUntracked<int>("dataType",0);
872  desc.addUntracked<int>("outMode",11);
873  desc.addUntracked<bool>("unCorrect",false);
874  desc.addUntracked<bool>("collapseDepth",false);
875  desc.addUntracked<std::string>("l1TrigName","L1_SingleJet60");
876  descriptions.add("HcalIsoTrkAnalyzer",desc);
877 }
878 
879 std::array<int,3> HcalIsoTrkAnalyzer::fillTree(std::vector< math::XYZTLorentzVector>& vecL1,
880  std::vector< math::XYZTLorentzVector>& vecL3,
881  math::XYZPoint& leadPV,
882  std::vector<spr::propagatedTrackDirection>& trkCaloDirections,
883  std::vector<spr::propagatedTrackID>& trkCaloDets,
884  const CaloGeometry* geo,
885  const CaloTopology* caloTopology,
886  const HcalTopology* theHBHETopology,
887  const EcalChannelStatus* theEcalChStatus,
888  const EcalSeverityLevelAlgo* theEcalSevlv,
889  edm::Handle<EcalRecHitCollection>& barrelRecHitsHandle,
890  edm::Handle<EcalRecHitCollection>& endcapRecHitsHandle,
894  const HcalRespCorrs* respCorrs) {
895 
896  int nSave(0), nLoose(0), nTight(0);
897  //Loop over tracks
898  std::vector<spr::propagatedTrackDirection>::const_iterator trkDetItr;
899  unsigned int nTracks(0), nselTracks(0);
900  t_nTrk = trkCaloDirections.size();
901  t_rhoh = (tower.isValid()) ? rhoh(tower) : 0;
902  for (trkDetItr = trkCaloDirections.begin(),nTracks=0;
903  trkDetItr != trkCaloDirections.end(); trkDetItr++,nTracks++) {
904  const reco::Track* pTrack = &(*(trkDetItr->trkItr));
905  math::XYZTLorentzVector v4(pTrack->px(), pTrack->py(),
906  pTrack->pz(), pTrack->p());
907 #ifdef EDM_ML_DEBUG
908  edm::LogVerbatim("HcalIsoTrack") << "This track : " << nTracks
909  << " (pt|eta|phi|p) :" << pTrack->pt()
910  << "|" << pTrack->eta() << "|"
911  << pTrack->phi() << "|" <<pTrack->p();
912 #endif
913  t_mindR2 = 999;
914  for (unsigned int k=0; k<vecL3.size(); ++k) {
915  double dr = dR(vecL3[k],v4);
916  if (dr<t_mindR2) {
917  t_mindR2 = dr;
918  }
919  }
920  t_mindR1 = (!vecL1.empty()) ? dR(vecL1[0],v4) : 999;
921 #ifdef EDM_ML_DEBUG
922  edm::LogVerbatim("HcalIsoTrack") << "Closest L3 object at dr :"
923  << t_mindR2 << " and from L1 " << t_mindR1;
924 #endif
925  t_ieta = t_iphi = 0;
926  if (trkDetItr->okHCAL) {
927  HcalDetId detId = (HcalDetId)(trkDetItr->detIdHCAL);
928  t_ieta = detId.ieta();
929  t_iphi = detId.iphi();
930  if (t_p > 40.0 && t_p <= 60.0) t_ietaAll->emplace_back(t_ieta);
931  }
932  //Selection of good track
933  t_selectTk = spr::goodTrack(pTrack,leadPV,selectionParameter_,false);
935  oneCutParameters.maxDxyPV = 10;
936  oneCutParameters.maxDzPV = 100;
937  oneCutParameters.maxInMiss = 2;
938  oneCutParameters.maxOutMiss= 2;
939  bool qltyFlag = spr::goodTrack(pTrack,leadPV,oneCutParameters,false);
940  oneCutParameters = selectionParameter_;
941  oneCutParameters.maxDxyPV = 10;
942  oneCutParameters.maxDzPV = 100;
943  t_qltyMissFlag = spr::goodTrack(pTrack,leadPV,oneCutParameters,false);
944  oneCutParameters = selectionParameter_;
945  oneCutParameters.maxInMiss = 2;
946  oneCutParameters.maxOutMiss= 2;
947  t_qltyPVFlag = spr::goodTrack(pTrack,leadPV,oneCutParameters,false);
948  double eIsolation = maxRestrictionP_*exp(slopeRestrictionP_*std::abs((double)t_ieta));
949  if (eIsolation < eIsolate1_) eIsolation = eIsolate1_;
950  if (eIsolation < eIsolate2_) eIsolation = eIsolate2_;
951 #ifdef EDM_ML_DEBUG
952  edm::LogVerbatim("HcalIsoTrack") << "qltyFlag|okECAL|okHCAL : "
953  << qltyFlag << "|" << trkDetItr->okECAL
954  << "|" << trkDetItr->okHCAL
955  << " eIsolation " << eIsolation;
956 #endif
957  t_qltyFlag = (qltyFlag && trkDetItr->okECAL && trkDetItr->okHCAL);
958  if (t_qltyFlag) {
959  nselTracks++;
960  int nNearTRKs(0);
961  std::vector<DetId> eIds;
962  std::vector<double> eHit;
963  t_eMipDR = spr::eCone_ecal(geo, barrelRecHitsHandle,
964  endcapRecHitsHandle, trkDetItr->pointHCAL,
965  trkDetItr->pointECAL, a_mipR_,
966  trkDetItr->directionECAL, eIds, eHit);
967  double eEcal(0);
968  for (unsigned int k=0; k<eIds.size(); ++k) {
969  const GlobalPoint& pos = geo->getPosition(eIds[k]);
970  double eta = std::abs(pos.eta());
971  double eThr = (eIds[k].subdetId() == EcalBarrel) ? hitEthrEB_ :
973  if (eHit[k] > eThr) eEcal += eHit[k];
974  }
975 #ifdef EDM_ML_DEBUG
976  edm::LogVerbatim("HcalIsoTrack") << "eMIP before and after: " << t_eMipDR
977  << ":" << eEcal;
978 #endif
979  t_eMipDR = eEcal;
980  t_emaxNearP = spr::chargeIsolationEcal(nTracks, trkCaloDets, geo,
981  caloTopology, 15,15);
982  const DetId cellE(trkDetItr->detIdECAL);
983  std::pair<double, bool> e11x11P =
984  spr::eECALmatrix(cellE,barrelRecHitsHandle,endcapRecHitsHandle,
985  *theEcalChStatus,geo,caloTopology,theEcalSevlv,
986  5,5,-100.0,-100.0,-100.0,100.0);
987  std::pair<double, bool> e15x15P =
988  spr::eECALmatrix(cellE,barrelRecHitsHandle,endcapRecHitsHandle,
989  *theEcalChStatus,geo,caloTopology,theEcalSevlv,7,7,
990  -100.0,-100.0,-100.0,100.0);
991  if (e11x11P.second && e15x15P.second) {
992  t_eAnnular = (e15x15P.first - e11x11P.first);
993  } else {
994  t_eAnnular =-(e15x15P.first - e11x11P.first);
995  }
996  t_hmaxNearP = spr::chargeIsolationCone(nTracks, trkCaloDirections,
997  a_charIsoR_, nNearTRKs, false);
998  const DetId cellH(trkDetItr->detIdHCAL);
999  double h5x5 = spr::eHCALmatrix(theHBHETopology,cellH,hbhe,2,2,false,true,
1000  -100.0,-100.0,-100.0,-100.0,-100.0,100.0);
1001  double h7x7 = spr::eHCALmatrix(theHBHETopology,cellH,hbhe,3,3,false,true,
1002  -100.0,-100.0,-100.0,-100.0,-100.0,100.0);
1003  t_hAnnular = h7x7 - h5x5;
1004 #ifdef EDM_ML_DEBUG
1005  edm::LogVerbatim("HcalIsoTrack") << "max p Near (Ecal) " << t_emaxNearP
1006  << " (Hcal) " << t_hmaxNearP
1007  << " Annular E (Ecal) " << e11x11P.first
1008  << ":" << e15x15P.first << ":"
1009  << t_eAnnular << " (Hcal) " << h5x5
1010  << ":" << h7x7 << ":" << t_hAnnular;
1011 #endif
1012  t_gentrackP = trackP(pTrack, genParticles);
1013  if (t_eMipDR < eEcalMax_ && t_hmaxNearP < eIsolation) {
1014  t_DetIds->clear(); t_HitEnergies->clear();
1015  t_DetIds1->clear(); t_HitEnergies1->clear();
1016  t_DetIds3->clear(); t_HitEnergies3->clear();
1017  int nRecHits(-999), nRecHits1(-999), nRecHits3(-999);
1018  std::vector<DetId> ids, ids1, ids3;
1019  std::vector<double> edet0, edet1, edet3;
1020  t_eHcal = spr::eCone_hcal(geo, hbhe, trkDetItr->pointHCAL,
1021  trkDetItr->pointECAL, a_coneR_,
1022  trkDetItr->directionHCAL,nRecHits,
1023  ids, edet0, useRaw_);
1024  storeEnergy(0,respCorrs,ids, edet0,t_eHcal, t_DetIds, t_HitEnergies);
1025 
1026  //----- hcal energy in the extended cone 1 (a_coneR+10) --------------
1027  t_eHcal10 = spr::eCone_hcal(geo, hbhe, trkDetItr->pointHCAL,
1028  trkDetItr->pointECAL, a_coneR1_,
1029  trkDetItr->directionHCAL,nRecHits1,
1030  ids1, edet1, useRaw_);
1031  storeEnergy(1,respCorrs,ids1,edet1,t_eHcal10,t_DetIds1,t_HitEnergies1);
1032 
1033  //----- hcal energy in the extended cone 3 (a_coneR+30) --------------
1034  t_eHcal30 = spr::eCone_hcal(geo, hbhe, trkDetItr->pointHCAL,
1035  trkDetItr->pointECAL, a_coneR2_,
1036  trkDetItr->directionHCAL,nRecHits3,
1037  ids3, edet3, useRaw_);
1038  storeEnergy(3,respCorrs,ids3,edet3,t_eHcal30,t_DetIds3,t_HitEnergies3);
1039 
1040  t_p = pTrack->p();
1041  t_pt = pTrack->pt();
1042  t_phi = pTrack->phi();
1043 
1044 #ifdef EDM_ML_DEBUG
1045  edm::LogVerbatim("HcalIsoTrack") << "This track : " << nTracks
1046  << " (pt|eta|phi|p) :" << t_pt
1047  << "|" << pTrack->eta() << "|"
1048  << t_phi << "|" << t_p
1049  << " Generator Level p "
1050  << t_gentrackP;
1051  edm::LogVerbatim("HcalIsoTrack") << "e_MIP " << t_eMipDR
1052  << " Chg Isolation " << t_hmaxNearP
1053  << " eHcal" << t_eHcal << " ieta "
1054  << t_ieta << " Quality "
1055  << t_qltyMissFlag << ":"
1056  << t_qltyPVFlag << ":" << t_selectTk;
1057  for (unsigned int ll=0; ll<t_DetIds->size(); ll++) {
1058  edm::LogVerbatim("HcalIsoTrack") << "det id is = " << t_DetIds->at(ll)
1059  << " hit enery is = "
1060  << t_HitEnergies->at(ll);
1061  }
1062  for (unsigned int ll=0; ll<t_DetIds1->size(); ll++) {
1063  edm::LogVerbatim("HcalIsoTrack") << "det id is = " << t_DetIds1->at(ll)
1064  << " hit enery is = "
1065  << t_HitEnergies1->at(ll);
1066  }
1067  for (unsigned int ll=0; ll<t_DetIds3->size(); ll++) {
1068  edm::LogVerbatim("HcalIsoTrack") << "det id is = " << t_DetIds3->at(ll)
1069  << " hit enery is = "
1070  << t_HitEnergies3->at(ll);
1071  }
1072 #endif
1073  bool accept(false);
1074  if (t_p>pTrackMin_) {
1075  if (t_p<pTrackLow_) {
1076  ++nLow_;
1077  if (prescaleLow_ <= 1) accept = true;
1078  else if (nLow_%prescaleLow_ == 1) accept = true;
1079  } else if (t_p>pTrackHigh_) {
1080  ++nHigh_;
1081  if (prescaleHigh_ <= 1) accept = true;
1082  else if (nHigh_%prescaleHigh_ == 1) accept = true;
1083  } else {
1084  accept = true;
1085  }
1086  }
1087  if (accept) {
1088  tree->Fill();
1089  nSave++;
1090  int type(0);
1091  if (t_eMipDR < 1.0) {
1092  if (t_hmaxNearP < eIsolate2_) { ++nLoose; type = 1;}
1093  if (t_hmaxNearP < eIsolate1_) { ++nTight; type = 2;}
1094  }
1095  if (t_p > 40.0 && t_p <= 60.0 && t_selectTk) {
1096  t_ietaGood->emplace_back(t_ieta);
1097  t_trackType->emplace_back(type);
1098  }
1099 #ifdef EDM_ML_DEBUG
1100  for (unsigned int k=0; k<t_trgbits->size(); k++) {
1101  edm::LogVerbatim("HcalIsoTrack") << "trigger bit is = "
1102  << t_trgbits->at(k);
1103  }
1104 #endif
1105  }
1106  }
1107  }
1108  }
1109  std::array<int,3> i3{ {nSave,nLoose,nTight} };
1110  return i3;
1111 }
1112 
1114  return reco::deltaR(vec1.eta(),vec1.phi(),vec2.eta(),vec2.phi());
1115 }
1116 
1119 
1120  double pmom = -1.0;
1121  if (genParticles.isValid()) {
1122  double mindR(999.9);
1123  for (const auto &p : (*genParticles)) {
1124  double dR = reco::deltaR(pTrack->eta(), pTrack->phi(),
1125  p.momentum().Eta(), p.momentum().Phi());
1126  if (dR < mindR) {
1127  mindR = dR; pmom = p.momentum().R();
1128  }
1129  }
1130  }
1131  return pmom;
1132 }
1133 
1135 
1136  std::vector<double> sumPFNallSMDQH2;
1137  sumPFNallSMDQH2.reserve(phibins_.size()*etabins_.size());
1138 
1139  for (auto eta : etabins_) {
1140  for (auto phi : phibins_) {
1141  double hadder = 0;
1142  for (const auto & pf_it : (*tower)) {
1143  if (fabs(eta-pf_it.eta())>etahalfdist_) continue;
1144  if (fabs(reco::deltaPhi(phi,pf_it.phi()))>phihalfdist_) continue;
1145  hadder += pf_it.hadEt();
1146  }
1147  sumPFNallSMDQH2.emplace_back(hadder);
1148  }
1149  }
1150 
1151  double evt_smdq(0);
1152  std::sort(sumPFNallSMDQH2.begin(),sumPFNallSMDQH2.end());
1153  if (sumPFNallSMDQH2.size()%2) evt_smdq = sumPFNallSMDQH2[(sumPFNallSMDQH2.size()-1)/2];
1154  else evt_smdq = (sumPFNallSMDQH2[sumPFNallSMDQH2.size()/2]+sumPFNallSMDQH2[(sumPFNallSMDQH2.size()-2)/2])/2.;
1155  double rhoh = evt_smdq/(etadist_*phidist_);
1156 #ifdef EDM_ML_DEBUG
1157  edm::LogVerbatim("HcalIsoTrack") << "Rho " << evt_smdq << ":" << rhoh;
1158 #endif
1159  return rhoh;
1160 }
1161 
1162 void HcalIsoTrkAnalyzer::storeEnergy(int indx, const HcalRespCorrs* respCorrs,
1163  const std::vector<DetId>& ids,
1164  std::vector<double>& edet,
1165  double & eHcal,
1166  std::vector<unsigned int> *detIds,
1167  std::vector<double> *hitEnergies) {
1168  double ehcal(0);
1169  if (unCorrect_) {
1170  for (unsigned int k=0; k<ids.size(); ++k) {
1171  double corr = (respCorrs->getValues(ids[k]))->getValue();
1172  if (corr != 0) edet[k] /= corr;
1173  ehcal += edet[k];
1174  }
1175  } else {
1176  for (const auto& en : edet) ehcal += en;
1177  }
1178  if (std::abs(ehcal-eHcal) > 0.001)
1179  edm::LogWarning("HcalIsoTrack") << "Check inconsistent energies: " << indx
1180  << " " << eHcal << ":" << ehcal
1181  << " from " << ids.size() << " cells";
1182  eHcal = hcalScale_*ehcal;
1183 
1184  if (collapseDepth_) {
1185  std::map<HcalDetId,double> hitMap;
1186  for (unsigned int k=0; k<ids.size(); ++k) {
1187  HcalDetId id = hdc_->mergedDepthDetId(HcalDetId(ids[k]));
1188  auto itr = hitMap.find(id);
1189  if (itr == hitMap.end()) {
1190  hitMap[id] = edet[k];
1191  } else {
1192  (itr->second) += edet[k];
1193  }
1194  }
1195  detIds->reserve(hitMap.size()); hitEnergies->reserve(hitMap.size());
1196  for (const auto& hit : hitMap) {
1197  detIds->emplace_back(hit.first.rawId());
1198  hitEnergies->emplace_back(hit.second);
1199  }
1200  } else {
1201  detIds->reserve(ids.size()); hitEnergies->reserve(ids.size());
1202  for (unsigned int k=0; k<ids.size(); ++k) {
1203  detIds->emplace_back(ids[k].rawId());
1204  hitEnergies->emplace_back(edet[k]);
1205  }
1206  }
1207 #ifdef EDM_ML_DEBUG
1208  edm::LogVerbatim("HcalIsoTrack") << "Input to storeEnergy with "
1209  << ids.size() << " cells";
1210  for (unsigned int k=0; k<ids.size(); ++k)
1211  edm::LogVerbatim("HcalIsoTrack") << "Hit [" << k << "] "
1212  << HcalDetId(ids[k]) << " E " << edet[k];
1213  edm::LogVerbatim("HcalIsoTrack") << "Output of storeEnergy with "
1214  << detIds->size() << " cells and Etot "
1215  << eHcal;
1216  for (unsigned int k=0; k<detIds->size(); ++k)
1217  edm::LogVerbatim("HcalIsoTrack") << "Hit [" << k << "] "
1218  << HcalDetId((*detIds)[k]) << " E "
1219  << (*hitEnergies)[k];
1220 #endif
1221 }
1222 
1223 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
unsigned int size() const
number of trigger paths in trigger table
static const std::string kSharedResource
Definition: TFileService.h:76
const std::string labelHBHE_
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
const std::vector< std::string > trigNames_
double p() const
momentum vector magnitude
Definition: TrackBase.h:615
double eCone_hcal(const CaloGeometry *geo, edm::Handle< T > &hits, const GlobalPoint &hpoint1, const GlobalPoint &point1, double dR, const GlobalVector &trackMom, int &nRecHits, double hbThr=-100, double heThr=-100, double hfThr=-100, double hoThr=-100, double tMin=-500, double tMax=500, int detOnly=-1, int useRaw=0, bool debug=false)
type
Definition: HCALResponse.h:21
spr::trackSelectionParameters selectionParameter_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
std::vector< spr::propagatedTrackID > propagateCALO(edm::Handle< reco::TrackCollection > &trkCollection, const CaloGeometry *geo, const MagneticField *bField, const std::string &theTrackQuality, bool debug=false)
~HcalIsoTrkAnalyzer() override
const unsigned int nTracks(const reco::Vertex &sv)
const std::string l1TrigName_
const double maxRestrictionP_
edm::EDGetTokenT< EcalRecHitCollection > tok_EE_
int id() const
getters
Definition: TriggerObject.h:55
RunNumber_t run() const
Definition: RunBase.h:40
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:135
const std::string labelRecVtx_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
float phi() const
Definition: TriggerObject.h:58
TrackQuality
track quality
Definition: TrackBase.h:151
HcalDetId mergedDepthDetId(const HcalDetId &id) const
const std::string labelTower_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< int > * t_ietaAll
bool accept() const
Has at least one path accepted the event?
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
int bunchCrossing() const
Definition: EventBase.h:66
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:63
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< BXVector< GlobalAlgBlk > > tok_alg_
const Item * getValues(DetId fId, bool throwOnFail=true) const
edm::EDGetTokenT< EcalRecHitCollection > tok_EB_
float energy() const
Definition: TriggerObject.h:65
std::vector< double > * t_HitEnergies
double chargeIsolationCone(unsigned int trkIndex, std::vector< spr::propagatedTrackDirection > &trkDirs, double dR, int &nNearTRKs, bool debug=false)
void beginRun(edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< reco::TrackCollection > tok_genTrack_
float eta() const
Definition: TriggerObject.h:57
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:627
double weight() const
void storeEnergy(int indx, const HcalRespCorrs *respCorrs, const std::vector< DetId > &ids, std::vector< double > &edet, double &eHcal, std::vector< unsigned int > *detIds, std::vector< double > *hitEnergies)
std::array< int, 3 > fillTree(std::vector< math::XYZTLorentzVector > &vecL1, std::vector< math::XYZTLorentzVector > &vecL3, math::XYZPoint &leadPV, std::vector< spr::propagatedTrackDirection > &trkCaloDirections, std::vector< spr::propagatedTrackID > &trkCaloDets, const CaloGeometry *geo, const CaloTopology *topo, const HcalTopology *theHBHETopology, const EcalChannelStatus *theEcalChStatus, const EcalSeverityLevelAlgo *theEcalSevlv, edm::Handle< EcalRecHitCollection > &barrelRecHitsHandle, edm::Handle< EcalRecHitCollection > &endcapRecHitsHandle, edm::Handle< HBHERecHitCollection > &hbhe, edm::Handle< CaloTowerCollection > &towerHandle, edm::Handle< reco::GenParticleCollection > &genParticles, const HcalRespCorrs *respCorrs)
#define nullptr
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::vector< double > etabins_
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
std::vector< unsigned int > * t_DetIds
const std::string names[nVars_]
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
std::vector< double > * t_HitEnergies1
void analyze(edm::Event const &, edm::EventSetup const &) override
static const bool useL1GtTriggerMenuLite(true)
edm::EDGetTokenT< reco::GenParticleCollection > tok_parts_
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:22
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
static const bool useL1EventSetup(false)
double chargeIsolationEcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, bool debug=false)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< trigger::TriggerEvent > tok_trigEvt_
std::vector< double > * t_HitEnergies3
const std::string labelGenTrack_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
edm::EDGetTokenT< reco::BeamSpot > tok_bs_
double trackP(const reco::Track *, const edm::Handle< reco::GenParticleCollection > &)
edm::EDGetTokenT< GenEventInfoProduct > tok_ew_
edm::Service< TFileService > fs
bool goodTrack(const reco::Track *pTrack, math::XYZPoint leadPV, trackSelectionParameters parameters, bool debug=false)
const TriggerObjectCollection & getObjects() const
Definition: TriggerEvent.h:98
std::vector< double > vec1
Definition: HCALResponse.h:15
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
l1t::L1TGlobalUtil * l1GtUtils_
unsigned int size() const
Get number of paths stored.
double pt() const
track transverse momentum
Definition: TrackBase.h:621
std::vector< bool > * t_hltbits
int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:74
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double eCone_ecal(const CaloGeometry *geo, edm::Handle< T > &barrelhits, edm::Handle< T > &endcaphits, const GlobalPoint &hpoint1, const GlobalPoint &point1, double dR, const GlobalVector &trackMom, int &nRecHits, double ebThr=-100, double eeThr=-100, double tMin=-500, double tMax=500, bool debug=false)
static std::string const triggerResults
Definition: EdmProvDump.cc:42
ParameterDescriptionBase * add(U const &iLabel, T const &value)
HLTConfigProvider hltConfig_
bool isValid() const
Definition: HandleBase.h:74
const std::string l3Filter_
JetCorrectorParameters corr
Definition: classes.h:5
#define M_PI
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:639
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
int k[5][pyjets_maxn]
const std::vector< std::string > & moduleLabels(unsigned int trigger) const
label(s) of module(s) on a trigger path
std::vector< bool > * t_trgbits
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
const double slopeRestrictionP_
Definition: DetId.h:18
const edm::InputTag filterTag(trigger::size_type index) const
Definition: TriggerEvent.h:103
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
const std::string l2Filter_
T const * product() const
Definition: Handle.h:81
const std::string processName_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< size_type > Keys
const std::string labelEB_
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
std::vector< double > phibins_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double rhoh(const edm::Handle< CaloTowerCollection > &)
double dR(math::XYZTLorentzVector &, math::XYZTLorentzVector &)
const edm::InputTag theTriggerResultsLabel_
std::string const & label() const
Definition: InputTag.h:36
T eta() const
Definition: PV3DBase.h:76
edm::EventID id() const
Definition: EventBase.h:60
void endRun(edm::Run const &, edm::EventSetup const &) override
HLT enums.
HcalIsoTrkAnalyzer(edm::ParameterSet const &)
const HcalDDDRecConstants * hdc_
T get() const
Definition: EventSetup.h:62
reco::TrackBase::TrackQuality minQuality
const std::vector< std::pair< std::string, bool > > & decisionsFinal()
const JetExtendedData & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
const Point & position() const
position
Definition: BeamSpot.h:62
std::vector< unsigned int > * t_DetIds1
edm::EDGetTokenT< edm::TriggerResults > tok_trigRes_
Definition: tree.py:1
const edm::InputTag triggerEvent_
float mass() const
Definition: TriggerObject.h:59
std::vector< int > * t_ietaGood
void beginJob() override
edm::EDGetTokenT< reco::VertexCollection > tok_recVtx_
edm::EDGetTokenT< CaloTowerCollection > tok_cala_
std::vector< vec1 > vec2
Definition: HCALResponse.h:16
std::vector< unsigned int > * t_DetIds3
T const * product() const
Definition: ESHandle.h:86
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
const std::string labelEE_
void setTopo(const HcalTopology *topo)
void retrieveL1(const edm::Event &iEvent, const edm::EventSetup &evSetup)
initialize the class (mainly reserve)
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:633
Definition: Run.h:44
const std::string l1Filter_
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)
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)
const std::string theTrackQuality_
std::vector< int > * t_trackType