CMS 3D CMS Logo

HcalHBHEMuonAnalyzer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <iostream>
3 #include <fstream>
4 #include <vector>
5 #include <TFile.h>
6 #include <TTree.h>
7 #include "TPRegexp.h"
8 
9 // user include files
19 
27 
29 
39 
42 
45 
49 
53 
57 
69 
70 //#define EDM_ML_DEBUG
71 
72 class HcalHBHEMuonAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns,edm::one::SharedResources> {
73 
74 public:
75  explicit HcalHBHEMuonAnalyzer(const edm::ParameterSet&);
76 
77  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
78 
79 private:
80  void beginJob() override;
81  void analyze(edm::Event const&, edm::EventSetup const&) override;
82  void beginRun(edm::Run const&, edm::EventSetup const&) override;
83  void endRun(edm::Run const&, edm::EventSetup const&) override {}
84  void clearVectors();
85  int matchId(const HcalDetId&, const HcalDetId&);
86  double activeLength(const DetId&);
87  bool isGoodVertex(const reco::Vertex& vtx);
88  double respCorr(const DetId& id);
89  double gainFactor(const edm::ESHandle<HcalDbService>&, const HcalDetId& id);
90  int depth16HE(int ieta, int iphi);
91  bool goodCell(const HcalDetId& hcid, const reco::Track* pTrack,
92  const CaloGeometry* geo, const MagneticField* bField);
93 
94  // ----------member data ---------------------------
100  const std::vector<std::string> triggers_;
101  const int verbosity_, useRaw_;
107 
111 
118 
120  static const int depthMax_ = 7;
121  TTree *tree_;
123  unsigned int goodVertex_;
125  std::vector<bool> muon_is_tight_, muon_is_medium_;
126  std::vector<double> ptGlob_, etaGlob_, phiGlob_, energyMuon_, pMuon_;
129  std::vector<bool> innerTrack_, outerTrack_, globalTrack_;
130  std::vector<double> chiTracker_, dxyTracker_, dzTracker_;
132  std::vector<double> tight_validFraction_, outerTrackChi_;
136  std::vector<int> globalMuonHits_, matchedStat_;
138  std::vector<double> isolationR04_, isolationR03_;
139  std::vector<double> ecalEnergy_, hcalEnergy_, hoEnergy_;
140  std::vector<bool> matchedId_, hcalHot_;
141  std::vector<double> ecal3x3Energy_, hcal1x1Energy_;
142  std::vector<unsigned int> ecalDetId_, hcalDetId_, ehcalDetId_;
143  std::vector<int> hcal_ieta_, hcal_iphi_;
144  std::vector<double> hcalDepthEnergy_[depthMax_];
145  std::vector<double> hcalDepthActiveLength_[depthMax_];
146  std::vector<double> hcalDepthEnergyHot_[depthMax_];
147  std::vector<double> hcalDepthActiveLengthHot_[depthMax_];
148  std::vector<double> hcalDepthChargeHot_[depthMax_];
149  std::vector<double> hcalDepthChargeHotBG_[depthMax_];
150  std::vector<double> hcalDepthEnergyCorr_[depthMax_];
151  std::vector<double> hcalDepthEnergyHotCorr_[depthMax_];
152  std::vector<bool> hcalDepthMatch_[depthMax_];
153  std::vector<bool> hcalDepthMatchHot_[depthMax_];
155  std::vector<std::string> all_triggers_;
156  std::vector<int> hltresults_;
157 
158  std::vector<HcalDDDRecConstants::HcalActiveLength> actHB, actHE;
159  std::map<DetId,double> corrValue_;
161 
162 };
163 
165  hlTriggerResults_(iConfig.getParameter<edm::InputTag>("hlTriggerResults")),
166  labelEBRecHit_(iConfig.getParameter<edm::InputTag>("labelEBRecHit")),
167  labelEERecHit_(iConfig.getParameter<edm::InputTag>("labelEERecHit")),
168  labelHBHERecHit_(iConfig.getParameter<edm::InputTag>("labelHBHERecHit")),
169  labelVtx_(iConfig.getParameter<std::string>("labelVertex")),
170  labelMuon_(iConfig.getParameter<std::string>("labelMuon")),
171  fileInCorr_(iConfig.getUntrackedParameter<std::string>("fileInCorr","")),
172  triggers_(iConfig.getParameter<std::vector<std::string>>("triggers")),
173  verbosity_(iConfig.getUntrackedParameter<int>("verbosity",0)),
174  useRaw_(iConfig.getParameter<int>("useRaw")),
175  unCorrect_(iConfig.getParameter<bool>("unCorrect")),
176  collapseDepth_(iConfig.getParameter<bool>("collapseDepth")),
177  isItPlan1_(iConfig.getParameter<bool>("isItPlan1")),
178  ignoreHECorr_(iConfig.getUntrackedParameter<bool>("ignoreHECorr",false)),
179  isItPreRecHit_(iConfig.getUntrackedParameter<bool>("isItPreRecHit",false)),
180  getCharge_(iConfig.getParameter<bool>("getCharge")),
181  writeRespCorr_(iConfig.getUntrackedParameter<bool>("writeRespCorr",false)),
183 
184  usesResource(TFileService::kSharedResource);
185  //now do what ever initialization is needed
186  kount_ = 0;
187  maxDepth_ = iConfig.getUntrackedParameter<int>("maxDepth",4);
189  else if (maxDepth_ < 1) maxDepth_ = 4;
190  std::string modnam = iConfig.getUntrackedParameter<std::string>("moduleName","");
191  std::string procnm = iConfig.getUntrackedParameter<std::string>("processName","");
192 
194  tok_trigRes_ = consumes<edm::TriggerResults>(hlTriggerResults_);
195  tok_EB_ = consumes<EcalRecHitCollection>(labelEBRecHit_);
196  tok_EE_ = consumes<EcalRecHitCollection>(labelEERecHit_);
197  tok_HBHE_ = consumes<HBHERecHitCollection>(labelHBHERecHit_);
198  if (modnam.empty()) {
199  tok_Vtx_ = consumes<reco::VertexCollection>(labelVtx_);
200  tok_Muon_ = consumes<reco::MuonCollection>(labelMuon_);
201  edm::LogVerbatim("HBHEMuon") << "Labels used: Trig " << hlTriggerResults_
202  << " Vtx " << labelVtx_ << " EB "
203  << labelEBRecHit_ << " EE "
204  << labelEERecHit_ << " HBHE "
205  << labelHBHERecHit_ << " MU " << labelMuon_;
206  } else {
207  tok_Vtx_ = consumes<reco::VertexCollection>(edm::InputTag(modnam,labelVtx_,procnm));
208  tok_Muon_ = consumes<reco::MuonCollection>(edm::InputTag(modnam,labelMuon_,procnm));
209  edm::LogVerbatim("HBHEMuon") << "Labels used Trig " << hlTriggerResults_
210  << "\n Vtx " << edm::InputTag(modnam,labelVtx_,procnm)
211  << "\n EB " << labelEBRecHit_
212  << "\n EE " << labelEERecHit_
213  << "\n HBHE " << labelHBHERecHit_
214  << "\n MU " << edm::InputTag(modnam,labelMuon_,procnm);
215  }
216 
217  if (!fileInCorr_.empty()) {
218  std::ifstream infile(fileInCorr_.c_str());
219  if (infile.is_open()) {
220  while (true) {
221  unsigned int id;
222  double cfac;
223  infile >> id >> cfac;
224  if (!infile.good()) break;
225  corrValue_[DetId(id)] = cfac;
226  }
227  infile.close();
228  }
229  }
230  useMyCorr_ = (!corrValue_.empty());
231  edm::LogVerbatim("HBHEMuon") << "Flags used: UseRaw " << useRaw_
232  << " GetCharge " << getCharge_ << " UnCorrect "
233  << unCorrect_ << " IgnoreHECorr "
234  << ignoreHECorr_ << " CollapseDepth "
235  << collapseDepth_ << ":" << mergedDepth_
236  << " IsItPlan1 " << isItPlan1_
237  << " IsItPreRecHit " << isItPreRecHit_
238  << " UseMyCorr " << useMyCorr_;
239 }
240 
241 //
242 // member functions
243 //
244 
245 // ------------ method called for each event ------------
247  ++kount_;
248  clearVectors();
249  runNumber_ = iEvent.id().run();
250  eventNumber_ = iEvent.id().event();
251  lumiNumber_ = iEvent.id().luminosityBlock();
252  bxNumber_ = iEvent.bunchCrossing();
253 #ifdef EDM_ML_DEBUG
254  edm::LogVerbatim("HBHEMuon") << "Run " << runNumber_ << " Event "
255  << eventNumber_ << " Lumi " << lumiNumber_
256  << " BX " << bxNumber_ << std::endl;
257 #endif
259  iEvent.getByToken(tok_trigRes_, _Triggers);
260 #ifdef EDM_ML_DEBUG
261  if ((verbosity_/10000)%10>0)
262  edm::LogVerbatim("HBHEMuon") << "Size of all triggers "
263  << all_triggers_.size() << std::endl;
264 #endif
265  int Ntriggers = all_triggers_.size();
266 #ifdef EDM_ML_DEBUG
267  if ((verbosity_/10000)%10>0)
268  edm::LogVerbatim("HBHEMuon") << "Size of HLT MENU: " << _Triggers->size()
269  << std::endl;
270 #endif
271  if (_Triggers.isValid()) {
272  const edm::TriggerNames &triggerNames_ = iEvent.triggerNames(*_Triggers);
273  std::vector<int> index;
274  for (int i=0; i<Ntriggers; i++) {
275  index.push_back(triggerNames_.triggerIndex(all_triggers_[i]));
276  int triggerSize = int( _Triggers->size());
277 #ifdef EDM_ML_DEBUG
278  if ((verbosity_/10000)%10>0)
279  edm::LogVerbatim("HBHEMuon") << "outside loop " << index[i]
280  << "\ntriggerSize " << triggerSize
281  << std::endl;
282 #endif
283  if (index[i] < triggerSize) {
284  hltresults_.push_back(_Triggers->accept(index[i]));
285 #ifdef EDM_ML_DEBUG
286  if ((verbosity_/10000)%10>0)
287  edm::LogVerbatim("HBHEMuon") << "Trigger_info " << triggerSize
288  << " triggerSize " << index[i]
289  << " trigger_index " << hltresults_.at(i)
290  << " hltresult" << std::endl;
291 #endif
292  } else {
293  if ((verbosity_/10000)%10>0)
294  edm::LogVerbatim("HBHEMuon") << "Requested HLT path \""
295  << "\" does not exist\n";
296  }
297  }
298  }
299 
300  // get handles to calogeometry and calotopology
302  iSetup.get<CaloGeometryRecord>().get(pG);
303  const CaloGeometry* geo = pG.product();
304 
306  iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
307  const MagneticField* bField = bFieldH.product();
308 
310  iSetup.get<EcalChannelStatusRcd>().get(ecalChStatus);
311  const EcalChannelStatus* theEcalChStatus = ecalChStatus.product();
312 
314  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
315 
316  edm::ESHandle<CaloTopology> theCaloTopology;
317  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
318  const CaloTopology *caloTopology = theCaloTopology.product();
319 
320  edm::ESHandle<HcalDbService> conditions;
321  iSetup.get<HcalDbRecord>().get(conditions);
322 
323  // Relevant blocks from iEvent
325  iEvent.getByToken(tok_Vtx_, vtx);
326 
327  edm::Handle<EcalRecHitCollection> barrelRecHitsHandle;
328  iEvent.getByToken(tok_EB_, barrelRecHitsHandle);
329  edm::Handle<EcalRecHitCollection> endcapRecHitsHandle;
330  iEvent.getByToken(tok_EE_, endcapRecHitsHandle);
331 
333  iEvent.getByToken(tok_HBHE_, hbhe);
334 
336  iEvent.getByToken(tok_Muon_, _Muon);
337 
338  // require a good vertex
339  math::XYZPoint pvx;
340  goodVertex_ = 0;
341  if (!vtx.isValid()) {
342 #ifdef EDM_ML_DEBUG
343  edm::LogVerbatim("HBHEMuon") << "No Good Vertex found == Reject\n";
344 #endif
345  return;
346  }
347  reco::VertexCollection::const_iterator firstGoodVertex = vtx->end();
348  for (reco::VertexCollection::const_iterator it = vtx->begin();
349  it != vtx->end(); it++) {
350  if (isGoodVertex(*it)) {
351  if (firstGoodVertex == vtx->end()) firstGoodVertex = it;
352  ++goodVertex_;
353  }
354  }
355  if (firstGoodVertex != vtx->end()) pvx = firstGoodVertex->position();
356 
357  bool accept(false);
358  if (_Muon.isValid() && barrelRecHitsHandle.isValid() &&
359  endcapRecHitsHandle.isValid() && hbhe.isValid()) {
360  for (reco::MuonCollection::const_iterator RecMuon = _Muon->begin(); RecMuon!= _Muon->end(); ++RecMuon) {
361 
362  muon_is_good_.push_back(RecMuon->isPFMuon());
363  muon_global_.push_back(RecMuon->isGlobalMuon());
364  muon_tracker_.push_back(RecMuon->isTrackerMuon());
365  ptGlob_.push_back((RecMuon)->pt());
366  etaGlob_.push_back(RecMuon->eta());
367  phiGlob_.push_back(RecMuon->phi());
368  energyMuon_.push_back(RecMuon->energy());
369  pMuon_.push_back(RecMuon->p());
370 #ifdef EDM_ML_DEBUG
371  edm::LogVerbatim("HBHEMuon") << "Energy:" << RecMuon->energy() << " P:"
372  << RecMuon->p() << std::endl;
373 #endif
374  muon_is_tight_.push_back(muon::isTightMuon(*RecMuon,*firstGoodVertex));
375  muon_is_medium_.push_back(muon::isMediumMuon(*RecMuon));
376  muon_trkKink.push_back(RecMuon->combinedQuality().trkKink);
377  muon_chi2LocalPosition.push_back(RecMuon->combinedQuality().chi2LocalPosition);
378  muon_segComp.push_back(muon::segmentCompatibility(*RecMuon));
379  // acessing tracker hits info
380  if (RecMuon->track().isNonnull()) {
381  trackerLayer_.push_back(RecMuon->track()->hitPattern().trackerLayersWithMeasurement());
382  } else {
383  trackerLayer_.push_back(-1);
384  }
385  if (RecMuon->innerTrack().isNonnull()) {
386  innerTrack_.push_back(true);
387  numPixelLayers_.push_back(RecMuon->innerTrack()->hitPattern().pixelLayersWithMeasurement());
388  chiTracker_.push_back(RecMuon->innerTrack()->normalizedChi2());
389  dxyTracker_.push_back(fabs(RecMuon->innerTrack()->dxy(pvx)));
390  dzTracker_.push_back(fabs(RecMuon->innerTrack()->dz(pvx)));
391  innerTrackpt_.push_back(RecMuon->innerTrack()->pt());
392  innerTracketa_.push_back(RecMuon->innerTrack()->eta());
393  innerTrackphi_.push_back(RecMuon->innerTrack()->phi());
394  tight_PixelHits_.push_back(RecMuon->innerTrack()->hitPattern().numberOfValidPixelHits());
395  tight_validFraction_.push_back(RecMuon->innerTrack()->validFraction());
396  } else {
397  innerTrack_.push_back(false);
398  numPixelLayers_.push_back(0);
399  chiTracker_.push_back(0);
400  dxyTracker_.push_back(0);
401  dzTracker_.push_back(0);
402  innerTrackpt_.push_back(0);
403  innerTracketa_.push_back(0);
404  innerTrackphi_.push_back(0);
405  tight_PixelHits_.push_back(0);
406  tight_validFraction_.push_back(-99);
407  }
408  // outer track info
409  if (RecMuon->outerTrack().isNonnull()) {
410  outerTrack_.push_back(true);
411  outerTrackPt_.push_back(RecMuon->outerTrack()->pt());
412  outerTrackEta_.push_back(RecMuon->outerTrack()->eta());
413  outerTrackPhi_.push_back(RecMuon->outerTrack()->phi());
414  outerTrackChi_.push_back(RecMuon->outerTrack()->normalizedChi2());
415  outerTrackHits_.push_back(RecMuon->outerTrack()->numberOfValidHits());
416  outerTrackRHits_.push_back(RecMuon->outerTrack()->recHitsSize());
417  } else {
418  outerTrack_.push_back(false);
419  outerTrackPt_.push_back(0);
420  outerTrackEta_.push_back(0);
421  outerTrackPhi_.push_back(0);
422  outerTrackChi_.push_back(0);
423  outerTrackHits_.push_back(0);
424  outerTrackRHits_.push_back(0);
425  }
426  // Tight Muon cuts
427  if (RecMuon->globalTrack().isNonnull()) {
428  globalTrack_.push_back(true);
429  chiGlobal_.push_back(RecMuon->globalTrack()->normalizedChi2());
430  globalMuonHits_.push_back(RecMuon->globalTrack()->hitPattern().numberOfValidMuonHits());
431  matchedStat_.push_back(RecMuon->numberOfMatchedStations());
432  globalTrckPt_.push_back(RecMuon->globalTrack()->pt());
433  globalTrckEta_.push_back(RecMuon->globalTrack()->eta());
434  globalTrckPhi_.push_back(RecMuon->globalTrack()->phi());
435  tight_TransImpara_.push_back(fabs(RecMuon->muonBestTrack()->dxy(pvx)));
436  tight_LongPara_.push_back(fabs(RecMuon->muonBestTrack()->dz(pvx)));
437  } else {
438  globalTrack_.push_back(false);
439  chiGlobal_.push_back(0);
440  globalMuonHits_.push_back(0);
441  matchedStat_.push_back(0);
442  globalTrckPt_.push_back(0);
443  globalTrckEta_.push_back(0);
444  globalTrckPhi_.push_back(0);
445  tight_TransImpara_.push_back(0);
446  tight_LongPara_.push_back(0);
447  }
448 
449  isolationR04_.push_back(((RecMuon->pfIsolationR04().sumChargedHadronPt + std::max(0.,RecMuon->pfIsolationR04().sumNeutralHadronEt + RecMuon->pfIsolationR04().sumPhotonEt - (0.5 *RecMuon->pfIsolationR04().sumPUPt))) / RecMuon->pt()) );
450 
451  isolationR03_.push_back(((RecMuon->pfIsolationR03().sumChargedHadronPt + std::max(0.,RecMuon->pfIsolationR03().sumNeutralHadronEt + RecMuon->pfIsolationR03().sumPhotonEt - (0.5 * RecMuon->pfIsolationR03().sumPUPt))) / RecMuon->pt()));
452 
453  ecalEnergy_.push_back(RecMuon->calEnergy().emS9);
454  hcalEnergy_.push_back(RecMuon->calEnergy().hadS9);
455  hoEnergy_.push_back(RecMuon->calEnergy().hoS9);
456 
457  double eEcal(0), eHcal(0), activeLengthTot(0), activeLengthHotTot(0);
458  double eHcalDepth[depthMax_], eHcalDepthHot[depthMax_];
459  double eHcalDepthC[depthMax_], eHcalDepthHotC[depthMax_];
460  double cHcalDepthHot[depthMax_], cHcalDepthHotBG[depthMax_];
461  double activeL[depthMax_], activeHotL[depthMax_];
462  bool matchDepth[depthMax_], matchDepthHot[depthMax_];
463  HcalDetId eHcalDetId[depthMax_];
464  unsigned int isHot(0);
465  bool tmpmatch(false);
466  int ieta(-1000), iphi(-1000);
467  for (int i=0; i<depthMax_; ++i) {
468  eHcalDepth[i] = eHcalDepthHot[i] = 0;
469  eHcalDepthC[i] = eHcalDepthHotC[i] = 0;
470  cHcalDepthHot[i] = cHcalDepthHotBG[i]= 0;
471  activeL[i] = activeHotL[i] = 0;
472  matchDepth[i] = matchDepthHot[i] = true;
473  }
474  if (RecMuon->innerTrack().isNonnull()) {
475  const reco::Track* pTrack = (RecMuon->innerTrack()).get();
476  spr::propagatedTrackID trackID = spr::propagateCALO(pTrack, geo, bField, (((verbosity_/100)%10>0)));
477  if ((RecMuon->p()>10.0) && (trackID.okHCAL)) accept = true;
478 
479  ecalDetId_.push_back((trackID.detIdECAL)());
480  hcalDetId_.push_back((trackID.detIdHCAL)());
481  ehcalDetId_.push_back((trackID.detIdEHCAL)());
482 
484  std::pair<bool,HcalDetId> info = spr::propagateHCALBack(pTrack, geo, bField, (((verbosity_/100)%10>0)));
485  if (info.first) {
486  check = info.second;
487  }
488 
489  bool okE = trackID.okECAL;
490  if (okE) {
491  const DetId isoCell(trackID.detIdECAL);
492  std::pair<double,bool> e3x3 = spr::eECALmatrix(isoCell,barrelRecHitsHandle,endcapRecHitsHandle,*theEcalChStatus,geo,caloTopology,sevlv.product(),1,1,-100.0,-100.0,-500.0,500.0,false);
493  eEcal = e3x3.first;
494  okE = e3x3.second;
495  }
496 #ifdef EDM_ML_DEBUG
497  edm::LogVerbatim("HBHEMuon") << "Propagate Track to ECAL: " << okE
498  << ":" << trackID.okECAL << " E " <<eEcal;
499 #endif
500 
501  if (trackID.okHCAL) {
502  DetId closestCell(trackID.detIdHCAL);
503  HcalDetId hcidt(closestCell.rawId());
504  if ((hcidt.ieta() == check.ieta()) && (hcidt.iphi() == check.iphi()))
505  tmpmatch = true;
506 #ifdef EDM_ML_DEBUG
507  edm::LogVerbatim("HBHEMuon") << "Front " << hcidt << " Back "
508  << info.first << ":" << check
509  << " Match " << tmpmatch;
510 #endif
511 
512  HcalSubdetector subdet = hcidt.subdet();
513  ieta = hcidt.ieta();
514  iphi = hcidt.iphi();
515  bool hborhe = (std::abs(ieta) == 16);
516 
517  eHcal = spr::eHCALmatrix(theHBHETopology_, closestCell, hbhe,0,0, false, true, -100.0, -100.0, -100.0, -100.0, -500.,500.,useRaw_);
518  std::vector<std::pair<double,int> > ehdepth;
519  spr::energyHCALCell((HcalDetId)closestCell, hbhe, ehdepth, maxDepth_, -100.0, -100.0, -100.0, -100.0, -500.0, 500.0, useRaw_, depth16HE(ieta,iphi), (((verbosity_/1000)%10)>0));
520  for (int i=0; i<depthMax_; ++i) eHcalDetId[i] = HcalDetId();
521  for (unsigned int i=0; i<ehdepth.size(); ++i) {
522  HcalSubdetector subdet0 = (hborhe) ? ((ehdepth[i].second >= depth16HE(ieta,iphi)) ? HcalEndcap : HcalBarrel) : subdet;
523  HcalDetId hcid0(subdet0,ieta,iphi,ehdepth[i].second);
524  double actL = activeLength(DetId(hcid0));
525  double ene = ehdepth[i].first;
526  bool tmpC(false);
527  if (ene > 0.0) {
528  if (!(theHBHETopology_->validHcal(hcid0))) {
529  edm::LogWarning("HBHEMuon") << "(1) Invalid ID " << hcid0
530  << " with E = " << ene;
531  edm::LogWarning("HBHEMuon") << HcalDetId(closestCell)
532  << " with " << ehdepth.size()
533  << " depths:";
534  for (const auto& ehd : ehdepth)
535  edm::LogWarning("HBHEMuon") << " " << ehd.second << ":"
536  << ehd.first;
537  } else {
538  tmpC = goodCell(hcid0, pTrack, geo, bField);
539  double enec(ene);
540  if (unCorrect_) {
541  double corr = (ignoreHECorr_ && (subdet0==HcalEndcap)) ? 1.0 : respCorr(DetId(hcid0));
542  if (corr != 0) ene /= corr;
543 #ifdef EDM_ML_DEBUG
544  HcalDetId id = (isItPlan1_ && isItPreRecHit_) ? hdc_->mergedDepthDetId(hcid0) : hcid0;
545  edm::LogVerbatim("HBHEMuon") << hcid0 << ":" << id << " Corr "
546  << corr;
547 #endif
548  }
549  int depth = ehdepth[i].second - 1;
550  if (collapseDepth_) {
551  HcalDetId id = hdc_->mergedDepthDetId(hcid0);
552  depth = id.depth() - 1;
553  }
554  eHcalDepth[depth] += ene;
555  eHcalDepthC[depth]+= enec;
556  activeL[depth] += actL;
557  activeLengthTot += actL;
558  matchDepth[depth] = (matchDepth[depth] && tmpC);
559 #ifdef EDM_ML_DEBUG
560  if ((verbosity_%10) > 0)
561  edm::LogVerbatim("HBHEMuon") << hcid0 << " E " << ene << ":"
562  << enec << " L " << actL
563  << " Match " << tmpC;
564 #endif
565  }
566  }
567  }
568 #ifdef EDM_ML_DEBUG
569  if ((verbosity_%10) > 0) {
570  edm::LogVerbatim("HBHEMuon") << hcidt << " Match " << tmpmatch
571  << " Depths " << ehdepth.size();
572  for (unsigned int k=0; k<ehdepth.size(); ++k)
573  edm::LogVerbatim("HBHEMuon") << " [" << k << ":"
574  << ehdepth[k].second << "] "
575  << matchDepth[k];
576  }
577 #endif
578  HcalDetId hotCell;
579  spr::eHCALmatrix(geo, theHBHETopology_, closestCell, hbhe, 1,1, hotCell, false, useRaw_, false);
580  isHot = matchId(closestCell,hotCell);
581  if (hotCell != HcalDetId()) {
582  subdet = HcalDetId(hotCell).subdet();
583  ieta = HcalDetId(hotCell).ieta();
584  iphi = HcalDetId(hotCell).iphi();
585  hborhe = (std::abs(ieta) == 16);
586  std::vector<std::pair<double,int> > ehdepth;
587  spr::energyHCALCell(hotCell, hbhe, ehdepth, maxDepth_, -100.0, -100.0, -100.0, -100.0, -500.0, 500.0, useRaw_, depth16HE(ieta,iphi), false);//(((verbosity_/1000)%10)>0 ));
588  for (int i=0; i<depthMax_; ++i) eHcalDetId[i] = HcalDetId();
589  for (unsigned int i=0; i<ehdepth.size(); ++i) {
590  HcalSubdetector subdet0 = (hborhe) ? ((ehdepth[i].second >= depth16HE(ieta,iphi)) ? HcalEndcap : HcalBarrel) : subdet;
591  HcalDetId hcid0(subdet0,ieta,iphi,ehdepth[i].second);
592  double actL = activeLength(DetId(hcid0));
593  double ene = ehdepth[i].first;
594  bool tmpC(false);
595  if (ene > 0.0) {
596  if (!(theHBHETopology_->validHcal(hcid0))) {
597  edm::LogWarning("HBHEMuon") << "(2) Invalid ID " << hcid0
598  << " with E = " << ene;
599  edm::LogWarning("HBHEMuon") << HcalDetId(hotCell)
600  << " with " << ehdepth.size()
601  << " depths:";
602  for (const auto& ehd : ehdepth)
603  edm::LogWarning("HBHEMuon") << " " << ehd.second << ":"
604  << ehd.first;
605  } else {
606  tmpC = goodCell(hcid0, pTrack, geo, bField);
607  double chg(ene), enec(ene);
608  if (unCorrect_) {
609  double corr = (ignoreHECorr_ && (subdet0==HcalEndcap)) ? 1.0 : respCorr(DetId(hcid0));
610  if (corr != 0) ene /= corr;
611 #ifdef EDM_ML_DEBUG
612  HcalDetId id = (isItPlan1_ && isItPreRecHit_) ? hdc_->mergedDepthDetId(hcid0) : hcid0;
613  edm::LogVerbatim("HBHEMuon") << hcid0 << ":" << id
614  << " Corr " << corr << " E "
615  << ene << ":" << enec;
616 #endif
617  }
618  if (getCharge_) {
619  double gain = gainFactor(conditions,hcid0);
620  if (gain != 0) chg /= gain;
621 #ifdef EDM_ML_DEBUG
622  edm::LogVerbatim("HBHEMuon") << hcid0 << " Gain " << gain
623  << " C " << chg;
624 #endif
625  }
626  int depth = ehdepth[i].second - 1;
627  if (collapseDepth_) {
628  HcalDetId id = hdc_->mergedDepthDetId(hcid0);
629  depth = id.depth() - 1;
630  }
631  eHcalDepthHot[depth] += ene;
632  eHcalDepthHotC[depth] += enec;
633  cHcalDepthHot[depth] += chg;
634  activeHotL[depth] += actL;
635  activeLengthHotTot += actL;
636  matchDepthHot[depth] = (matchDepthHot[depth] && tmpC);
637 #ifdef EDM_ML_DEBUG
638  if ((verbosity_%10) > 0)
639  edm::LogVerbatim("HBHEMuon") << hcid0 << " depth " << depth
640  << " E " << ene << ":" << enec
641  << " C " << chg << " L "
642  << actL << " Match " << tmpC;
643 #endif
644  }
645  }
646  }
647 
648  HcalDetId oppCell(subdet,-ieta,iphi,HcalDetId(hotCell).depth());
649  std::vector<std::pair<double,int> > ehdeptho;
650  spr::energyHCALCell(oppCell, hbhe, ehdeptho, maxDepth_, -100.0, -100.0, -100.0, -100.0, -500.0, 500.0, useRaw_, depth16HE(-ieta,iphi), false); //(((verbosity_/1000)%10)>0));
651  for (unsigned int i=0; i<ehdeptho.size(); ++i) {
652  HcalSubdetector subdet0 = (hborhe) ? ((ehdeptho[i].second >= depth16HE(-ieta,iphi)) ? HcalEndcap : HcalBarrel) : subdet;
653  HcalDetId hcid0(subdet0,-ieta,iphi,ehdeptho[i].second);
654  double ene = ehdeptho[i].first;
655  if (ene > 0.0) {
656  if (!(theHBHETopology_->validHcal(hcid0))) {
657  edm::LogWarning("HBHEMuon") << "(3) Invalid ID " << hcid0
658  << " with E = " << ene;
659  edm::LogWarning("HBHEMuon") << oppCell << " with "
660  << ehdeptho.size() << " depths:";
661  for (const auto& ehd : ehdeptho)
662  edm::LogWarning("HBHEMuon") << " " << ehd.second << ":"
663  << ehd.first;
664  } else {
665  double chg(ene);
666  if (unCorrect_) {
667  double corr = (ignoreHECorr_ && (subdet0==HcalEndcap)) ? 1.0 : respCorr(DetId(hcid0));
668  if (corr != 0) ene /= corr;
669 #ifdef EDM_ML_DEBUG
670  HcalDetId id = (isItPlan1_ && isItPreRecHit_) ? hdc_->mergedDepthDetId(hcid0) : hcid0;
671  edm::LogVerbatim("HBHEMuon") << hcid0 << ":" << id
672  << " Corr " << corr << " E "
673  << ene << ":"
674  << ehdeptho[i].first;
675 #endif
676  }
677  if (getCharge_) {
678  double gain = gainFactor(conditions,hcid0);
679  if (gain != 0) chg /= gain;
680 #ifdef EDM_ML_DEBUG
681  edm::LogVerbatim("HBHEMuon") << hcid0 << " Gain " << gain
682  << " C " << chg;
683 #endif
684  }
685  int depth = ehdeptho[i].second - 1;
686  if (collapseDepth_) {
687  HcalDetId id = hdc_->mergedDepthDetId(hcid0);
688  depth = id.depth() - 1;
689  }
690  cHcalDepthHotBG[depth] += chg;
691 #ifdef EDM_ML_DEBUG
692  if ((verbosity_%10) > 0)
693  edm::LogVerbatim("HBHEMuon") << hcid0 << " Depth " << depth
694  << " E " << ene << " C "
695  << chg;
696 #endif
697  }
698  }
699  }
700  }
701  }
702 #ifdef EDM_ML_DEBUG
703  edm::LogVerbatim("HBHEMuon") << "Propagate Track to HCAL: "
704  << trackID.okHCAL << " Match " << tmpmatch
705  << " Hot " << isHot << " Energy "
706  << eHcal << std::endl;
707 #endif
708 
709  } else {
710  ecalDetId_.push_back(0);
711  hcalDetId_.push_back(0);
712  ehcalDetId_.push_back(0);
713  }
714 
715  matchedId_.push_back(tmpmatch);
716  ecal3x3Energy_.push_back(eEcal);
717  hcal1x1Energy_.push_back(eHcal);
718  hcal_ieta_.push_back(ieta);
719  hcal_iphi_.push_back(iphi);
720  for (int i=0; i<depthMax_; ++i) {
721  hcalDepthEnergy_[i].push_back(eHcalDepth[i]);
722  hcalDepthActiveLength_[i].push_back(activeL[i]);
723  hcalDepthEnergyHot_[i].push_back(eHcalDepthHot[i]);
724  hcalDepthActiveLengthHot_[i].push_back(activeHotL[i]);
725  hcalDepthEnergyCorr_[i].push_back(eHcalDepthC[i]);
726  hcalDepthEnergyHotCorr_[i].push_back(eHcalDepthHotC[i]);
727  hcalDepthChargeHot_[i].push_back(cHcalDepthHot[i]);
728  hcalDepthChargeHotBG_[i].push_back(cHcalDepthHotBG[i]);
729  hcalDepthMatch_[i].push_back(matchDepth[i]);
730  hcalDepthMatchHot_[i].push_back(matchDepthHot[i]);
731  }
732  hcalActiveLength_.push_back(activeLengthTot);
733  hcalHot_.push_back(isHot);
734  hcalActiveLengthHot_.push_back(activeLengthHotTot);
735  }
736  }
737  if (accept) {
738 #ifdef EDM_ML_DEBUG
739  for (unsigned int i=0; i<hcal_ieta_.size(); ++i)
740  edm::LogVerbatim("HBHEMuon") << "[" << i << "] ieta/iphi for entry to "
741  << "HCAL has value of " << hcal_ieta_[i]
742  << ":" << hcal_iphi_[i];
743 #endif
744  tree_->Fill();
745  }
746 }
747 
748 // ------------ method called once each job just before starting event loop ------------
750 
751  tree_ = fs->make<TTree>("TREE", "TREE");
752  tree_->Branch("Event_No", &eventNumber_);
753  tree_->Branch("Run_No", &runNumber_);
754  tree_->Branch("LumiNumber", &lumiNumber_);
755  tree_->Branch("BXNumber", &bxNumber_);
756  tree_->Branch("GoodVertex", &goodVertex_);
757  tree_->Branch("PF_Muon", &muon_is_good_);
758  tree_->Branch("Global_Muon", &muon_global_);
759  tree_->Branch("Tracker_muon", &muon_tracker_);
760  tree_->Branch("MuonIsTight", &muon_is_tight_);
761  tree_->Branch("MuonIsMedium", &muon_is_medium_);
762  tree_->Branch("pt_of_muon", &ptGlob_);
763  tree_->Branch("eta_of_muon", &etaGlob_);
764  tree_->Branch("phi_of_muon", &phiGlob_);
765  tree_->Branch("energy_of_muon", &energyMuon_);
766  tree_->Branch("p_of_muon", &pMuon_);
767  tree_->Branch("muon_trkKink", &muon_trkKink);
768  tree_->Branch("muon_chi2LocalPosition", &muon_chi2LocalPosition);
769  tree_->Branch("muon_segComp", &muon_segComp);
770 
771  tree_->Branch("TrackerLayer", &trackerLayer_);
772  tree_->Branch("NumPixelLayers", &numPixelLayers_);
773  tree_->Branch("InnerTrackPixelHits", &tight_PixelHits_);
774  tree_->Branch("innerTrack", &innerTrack_);
775  tree_->Branch("chiTracker", &chiTracker_);
776  tree_->Branch("DxyTracker", &dxyTracker_);
777  tree_->Branch("DzTracker", &dzTracker_);
778  tree_->Branch("innerTrackpt", &innerTrackpt_);
779  tree_->Branch("innerTracketa", &innerTracketa_);
780  tree_->Branch("innerTrackphi", &innerTrackphi_);
781  tree_->Branch("tight_validFraction", &tight_validFraction_);
782 
783  tree_->Branch("OuterTrack", &outerTrack_);
784  tree_->Branch("OuterTrackChi", &outerTrackChi_);
785  tree_->Branch("OuterTrackPt", &outerTrackPt_);
786  tree_->Branch("OuterTrackEta", &outerTrackEta_);
787  tree_->Branch("OuterTrackPhi", &outerTrackPhi_);
788  tree_->Branch("OuterTrackHits", &outerTrackHits_);
789  tree_->Branch("OuterTrackRHits", &outerTrackRHits_);
790 
791  tree_->Branch("GlobalTrack", &globalTrack_);
792  tree_->Branch("GlobalTrckPt", &globalTrckPt_);
793  tree_->Branch("GlobalTrckEta", &globalTrckEta_);
794  tree_->Branch("GlobalTrckPhi", &globalTrckPhi_);
795  tree_->Branch("Global_Muon_Hits", &globalMuonHits_);
796  tree_->Branch("MatchedStations", &matchedStat_);
797  tree_->Branch("GlobTrack_Chi", &chiGlobal_);
798  tree_->Branch("Tight_LongitudinalImpactparameter",&tight_LongPara_);
799  tree_->Branch("Tight_TransImpactparameter", &tight_TransImpara_);
800 
801  tree_->Branch("IsolationR04", &isolationR04_);
802  tree_->Branch("IsolationR03", &isolationR03_);
803  tree_->Branch("ecal_3into3", &ecalEnergy_);
804  tree_->Branch("hcal_3into3", &hcalEnergy_);
805  tree_->Branch("tracker_3into3", &hoEnergy_);
806 
807  tree_->Branch("matchedId", &matchedId_);
808  tree_->Branch("hcal_cellHot", &hcalHot_);
809 
810  tree_->Branch("ecal_3x3", &ecal3x3Energy_);
811  tree_->Branch("hcal_1x1", &hcal1x1Energy_);
812  tree_->Branch("ecal_detID", &ecalDetId_);
813  tree_->Branch("hcal_detID", &hcalDetId_);
814  tree_->Branch("ehcal_detID", &ehcalDetId_);
815  tree_->Branch("hcal_ieta", &hcal_ieta_);
816  tree_->Branch("hcal_iphi", &hcal_iphi_);
817 
818  char name[100];
819  for (int k=0; k<maxDepth_; ++k) {
820  sprintf (name, "hcal_edepth%d", (k+1));
821  tree_->Branch(name, &hcalDepthEnergy_[k]);
822  sprintf (name, "hcal_activeL%d", (k+1));
823  tree_->Branch(name, &hcalDepthActiveLength_[k]);
824  sprintf (name, "hcal_edepthHot%d", (k+1));
825  tree_->Branch(name, &hcalDepthEnergyHot_[k]);
826  sprintf (name, "hcal_activeHotL%d", (k+1));
827  tree_->Branch(name, &hcalDepthActiveLengthHot_[k]);
828  sprintf (name, "hcal_cdepthHot%d", (k+1));
829  tree_->Branch(name, &hcalDepthChargeHot_[k]);
830  sprintf (name, "hcal_cdepthHotBG%d", (k+1));
831  tree_->Branch(name, &hcalDepthChargeHotBG_[k]);
832  sprintf (name, "hcal_edepthCorrect%d", (k+1));
833  tree_->Branch(name, &hcalDepthEnergyCorr_[k]);
834  sprintf (name, "hcal_edepthHotCorrect%d", (k+1));
835  tree_->Branch(name, &hcalDepthEnergyHotCorr_[k]);
836  sprintf (name, "hcal_depthMatch%d", (k+1));
837  tree_->Branch(name, &hcalDepthMatch_[k]);
838  sprintf (name, "hcal_depthMatchHot%d", (k+1));
839  tree_->Branch(name, &hcalDepthMatchHot_[k]);
840  }
841 
842  tree_->Branch("activeLength", &hcalActiveLength_);
843  tree_->Branch("activeLengthHot", &hcalActiveLengthHot_);
844 
845  tree_->Branch("hltresults", &hltresults_);
846  tree_->Branch("all_triggers", &all_triggers_);
847 }
848 
849 // ------------ method called when starting to processes a run ------------
850 void HcalHBHEMuonAnalyzer::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
851 
853  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
854  hdc_ = pHRNDC.product();
855  actHB.clear();
856  actHE.clear();
857  actHB = hdc_->getThickActive(0);
858  actHE = hdc_->getThickActive(1);
859 #ifdef EDM_ML_DEBUG
860  unsigned int k1(0), k2(0);
861  edm::LogVerbatim("HBHEMuon") << actHB.size() << " Active Length for HB";
862  for (const auto& act : actHB) {
863  edm::LogVerbatim("HBHEMuon") << "[" << k1 << "] ieta " << act.ieta
864  << " depth " << act.depth << " zside "
865  << act.zside << " type " << act.stype
866  << " phi " << act.iphis.size() << ":"
867  << act.iphis[0] << " L " << act.thick;
868  HcalDetId hcid1(HcalBarrel,(act.ieta)*(act.zside),act.iphis[0],act.depth);
869  HcalDetId hcid2 = mergedDepth_ ? hdc_->mergedDepthDetId(hcid1) : hcid1;
870  edm::LogVerbatim("HBHEMuon") << hcid1 << " | " << hcid2 << " L "
871  << activeLength(DetId(hcid2));
872  ++k1;
873  }
874  edm::LogVerbatim("HBHEMuon") << actHE.size() << " Active Length for HE";
875  for (const auto& act : actHE) {
876  edm::LogVerbatim("HBHEMuon") << "[" << k2 << "] ieta " << act.ieta
877  << " depth " << act.depth << " zside "
878  << act.zside << " type " << act.stype
879  << " phi " << act.iphis.size() << ":"
880  << act.iphis[0] << " L " << act.thick;
881  HcalDetId hcid1(HcalEndcap,(act.ieta)*(act.zside),act.iphis[0],act.depth);
882  HcalDetId hcid2 = mergedDepth_ ? hdc_->mergedDepthDetId(hcid1) : hcid1;
883  edm::LogVerbatim("HBHEMuon") << hcid1 << " | " << hcid2 << " L "
884  << activeLength(DetId(hcid2));
885  ++k2;
886  }
887 #endif
888 
889  bool changed = true;
890  all_triggers_.clear();
891  if (hltConfig_.init(iRun, iSetup, "HLT" , changed)) {
892  // if init returns TRUE, initialisation has succeeded!
893 #ifdef EDM_ML_DEBUG
894  edm::LogVerbatim("HBHEMuon") << "HLT config with process name "
895  << "HLT" << " successfully extracted"
896  << std::endl;
897 #endif
898  unsigned int ntriggers = hltConfig_.size();
899  for (unsigned int t=0;t<ntriggers;++t) {
901  for (unsigned int ik=0; ik<6; ++ik) {
902  if (hltname.find(triggers_[ik])!=std::string::npos ){
903  all_triggers_.push_back(hltname);
904  break;
905  }
906  }
907  }//loop over ntriggers
908  edm::LogVerbatim("HBHEMuon") << "All triggers size in begin run "
909  << all_triggers_.size() << std::endl;
910  } else {
911  edm::LogError("HBHEMuon") << "Error! HLT config extraction with process "
912  << "name HLT failed";
913  }
914 
916  iSetup.get<HcalRecNumberingRecord>().get(htopo);
917  theHBHETopology_ = htopo.product();
918 
920  iSetup.get<HcalRespCorrsRcd>().get(resp);
921  respCorrs_ = new HcalRespCorrs(*resp.product());
923 
924  // Write correction factors for all HB/HE events
925  if (writeRespCorr_) {
927  iSetup.get<CaloGeometryRecord>().get(pG);
928  const CaloGeometry* geo = pG.product();
930  const std::vector<DetId>& ids = gHcal->getValidDetIds(DetId::Hcal,0);
931  edm::LogVerbatim("HBHEMuon") << "\nTable of Correction Factors for Run "
932  << iRun.run() << "\n";
933  for (auto const& id: ids) {
934  if ((id.det() == DetId::Hcal) &&
935  ((id.subdetId() == HcalBarrel) || (id.subdetId() == HcalEndcap))) {
936  edm::LogVerbatim("HBHEMuon") << HcalDetId(id) << " " << id.rawId() <<" "
937  << (respCorrs_->getValues(id))->getValue();
938  }
939  }
940  }
941 }
942 
943 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
946  desc.add<edm::InputTag>("hlTriggerResults",edm::InputTag("TriggerResults","","HLT"));
947  desc.add<edm::InputTag>("labelEBRecHit",edm::InputTag("ecalRecHit","EcalRecHitsEB"));
948  desc.add<edm::InputTag>("labelEERecHit",edm::InputTag("ecalRecHit","EcalRecHitsEE"));
949  desc.add<edm::InputTag>("labelHBHERecHit",edm::InputTag("hbhereco"));
950  desc.add<std::string>("labelVertex","offlinePrimaryVertices");
951  desc.add<std::string>("labelMuon","muons");
952  std::vector<std::string> trig = {"HLT_IsoMu17","HLT_IsoMu20",
953  "HLT_IsoMu24","HLT_IsoMu27",
954  "HLT_Mu45","HLT_Mu50"};
955  desc.add<std::vector<std::string>>("triggers",trig);
956  desc.addUntracked<int>("verbosity",0);
957  desc.add<int>("useRaw",0);
958  desc.add<bool>("unCorrect",false);
959  desc.add<bool>("getCharge",false);
960  desc.add<bool>("collapseDepth",false);
961  desc.add<bool>("isItPlan1",false);
962  desc.addUntracked<bool>("ignoreHECorr",false);
963  desc.addUntracked<bool>("isItPreRecHit",false);
964  desc.addUntracked<std::string>("moduleName","");
965  desc.addUntracked<std::string>("processName","");
966  desc.addUntracked<int>("maxDepth",4);
967  desc.addUntracked<std::string>("fileInCorr","");
968  desc.addUntracked<bool>("writeRespCorr",false);
969  descriptions.add("hcalHBHEMuon",desc);
970 }
971 
974  eventNumber_ = -99999;
975  runNumber_ = -99999;
976  lumiNumber_ = -99999;
977  bxNumber_ = -99999;
978  goodVertex_ = -99999;
979 
980  muon_is_good_.clear();
981  muon_global_.clear();
982  muon_tracker_.clear();
983  ptGlob_.clear();
984  etaGlob_.clear();
985  phiGlob_.clear();
986  energyMuon_.clear();
987  pMuon_.clear();
988  muon_trkKink.clear();
989  muon_chi2LocalPosition.clear();
990  muon_segComp.clear();
991  muon_is_tight_.clear();
992  muon_is_medium_.clear();
993 
994  trackerLayer_.clear();
995  numPixelLayers_.clear();
996  tight_PixelHits_.clear();
997  innerTrack_.clear();
998  chiTracker_.clear();
999  dxyTracker_.clear();
1000  dzTracker_.clear();
1001  innerTrackpt_.clear();
1002  innerTracketa_.clear();
1003  innerTrackphi_.clear();
1004  tight_validFraction_.clear();
1005 
1006  outerTrack_.clear();
1007  outerTrackPt_.clear();
1008  outerTrackEta_.clear();
1009  outerTrackPhi_.clear();
1010  outerTrackHits_.clear();
1011  outerTrackRHits_.clear();
1012  outerTrackChi_.clear();
1013 
1014  globalTrack_.clear();
1015  globalTrckPt_.clear();
1016  globalTrckEta_.clear();
1017  globalTrckPhi_.clear();
1018  globalMuonHits_.clear();
1019  matchedStat_.clear();
1020  chiGlobal_.clear();
1021  tight_LongPara_.clear();
1022  tight_TransImpara_.clear();
1023 
1024  isolationR04_.clear();
1025  isolationR03_.clear();
1026  ecalEnergy_.clear();
1027  hcalEnergy_.clear();
1028  hoEnergy_.clear();
1029  matchedId_.clear();
1030  hcalHot_.clear();
1031  ecal3x3Energy_.clear();
1032  hcal1x1Energy_.clear();
1033  ecalDetId_.clear();
1034  hcalDetId_.clear();
1035  ehcalDetId_.clear();
1036  hcal_ieta_.clear();
1037  hcal_iphi_.clear();
1038  for (int i=0; i<maxDepth_; ++i) {
1039  hcalDepthEnergy_[i].clear();
1040  hcalDepthActiveLength_[i].clear();
1041  hcalDepthEnergyHot_[i].clear();
1042  hcalDepthActiveLengthHot_[i].clear();
1043  hcalDepthChargeHot_[i].clear();
1044  hcalDepthChargeHotBG_[i].clear();
1045  hcalDepthEnergyCorr_[i].clear();
1046  hcalDepthEnergyHotCorr_[i].clear();
1047  hcalDepthMatch_[i].clear();
1048  hcalDepthMatchHot_[i].clear();
1049  }
1050  hcalActiveLength_.clear();
1051  hcalActiveLengthHot_.clear();
1052  hltresults_.clear();
1053 }
1054 
1056 
1057  HcalDetId kd1(id1.subdet(),id1.ieta(),id1.iphi(),1);
1058  HcalDetId kd2(id1.subdet(),id2.ieta(),id2.iphi(),1);
1059  int match = ((kd1 == kd2) ? 1 : 0);
1060  return match;
1061 }
1062 
1064  HcalDetId id(id_);
1065  int ieta = id.ietaAbs();
1066  int zside= id.zside();
1067  int iphi = id.iphi();
1068  std::vector<int> dpths;
1069  if (mergedDepth_) {
1070  std::vector<HcalDetId> ids;
1071  hdc_->unmergeDepthDetId(id,ids);
1072  for (auto idh : ids)
1073  dpths.emplace_back(idh.depth());
1074  } else {
1075  dpths.emplace_back(id.depth());
1076  }
1077  double lx(0);
1078  if (id.subdet() == HcalBarrel) {
1079  for (unsigned int i=0; i<actHB.size(); ++i) {
1080  if ((ieta == actHB[i].ieta) && (zside == actHB[i].zside) &&
1081  (std::find(dpths.begin(),dpths.end(),actHB[i].depth) != dpths.end())&&
1082  (std::find(actHB[i].iphis.begin(),actHB[i].iphis.end(),iphi) !=
1083  actHB[i].iphis.end())) {
1084  lx += actHB[i].thick;
1085  }
1086  }
1087  } else {
1088  for (unsigned int i=0; i<actHE.size(); ++i) {
1089  if ((ieta == actHE[i].ieta) && (zside == actHE[i].zside) &&
1090  (std::find(dpths.begin(),dpths.end(),actHE[i].depth) != dpths.end())&&
1091  (std::find(actHE[i].iphis.begin(),actHE[i].iphis.end(),iphi) !=
1092  actHE[i].iphis.end())) {
1093  lx += actHE[i].thick;
1094  }
1095  }
1096  }
1097  return lx;
1098 }
1099 
1101  if (vtx.isFake()) return false;
1102  if (vtx.ndof() < 4) return false;
1103  if (vtx.position().Rho() > 2.) return false;
1104  if (fabs(vtx.position().Z()) > 24.) return false;
1105  return true;
1106 }
1107 
1109  double cfac(1.0);
1110  if (useMyCorr_) {
1111  auto itr = corrValue_.find(id);
1112  if (itr != corrValue_.end()) cfac = itr->second;
1113  } else if (respCorrs_ != nullptr) {
1114  cfac = (respCorrs_->getValues(id))->getValue();
1115  }
1116  return cfac;
1117 }
1118 
1120  double gain(0.0);
1121  const HcalCalibrations& calibs=conditions->getHcalCalibrations(id);
1122  for (int capid=0; capid<4; ++capid)
1123  gain += (0.25*calibs.respcorrgain(capid));
1124  return gain;
1125 }
1126 
1127 int HcalHBHEMuonAnalyzer::depth16HE(int ieta, int iphi) {
1128  // Transition between HB/HE is special
1129  // For Run 1 or for Plan1 standard reconstruction it is 3
1130  // For runs beyond 2018 or in Plan1 for HEP17 it is 4
1131  int zside = (ieta > 0) ? 1 : -1;
1132  int depth = theHBHETopology_->dddConstants()->getMinDepth(1,16,iphi,zside);
1133  if (isItPlan1_ && (!isItPreRecHit_)) depth = 3;
1134 #ifdef EDM_ML_DEBUG
1135  edm::LogVerbatim("HBHEMuon") << "Plan1 " << isItPlan1_ << " PreRecHit "
1136  << isItPreRecHit_ << " phi " << iphi
1137  << " depth " << depth;
1138 #endif
1139  return depth;
1140 }
1141 
1143  const reco::Track* pTrack,
1144  const CaloGeometry* geo,
1145  const MagneticField* bField) {
1146 
1147  std::pair<double,double> rz = hdc_->getRZ(hcid);
1148  bool typeRZ = (hcid.subdet() == HcalEndcap) ? false : true;
1149  bool match = spr::propagateHCAL(pTrack, geo, bField, typeRZ, rz, (((verbosity_/10000)%10)>0));
1150  return match;
1151 }
1152 
1153 //define this as a plug-in
1155 
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
std::vector< double > tight_validFraction_
std::vector< int > hltresults_
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
EventNumber_t event() const
Definition: EventID.h:41
std::vector< double > innerTrackphi_
HcalHBHEMuonAnalyzer(const edm::ParameterSet &)
const std::string labelMuon_
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)
const HcalTopology * theHBHETopology_
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:167
static const TGPicture * info(bool iBackgroundIsBlack)
std::vector< std::string > all_triggers_
std::vector< double > dzTracker_
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:146
RunNumber_t run() const
Definition: RunBase.h:40
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void analyze(edm::Event const &, edm::EventSetup const &) override
const std::string & triggerName(unsigned int triggerIndex) const
std::vector< double > outerTrackChi_
bool isMediumMuon(const reco::Muon &, bool run2016_hip_mitigation=false)
const float chg[109]
Definition: CoreSimTrack.cc:5
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
HcalDetId mergedDepthDetId(const HcalDetId &id) const
std::vector< double > hcalDepthChargeHotBG_[depthMax_]
std::vector< double > tight_LongPara_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::pair< math::XYZPoint, bool > propagateHCAL(const reco::Track *, const MagneticField *, bool debug=false)
edm::EDGetTokenT< EcalRecHitCollection > tok_EE_
std::vector< unsigned int > ecalDetId_
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
std::vector< double > pMuon_
bool accept() const
Has at least one path accepted the event?
std::vector< double > globalTrckPhi_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
std::vector< bool > muon_is_tight_
int bunchCrossing() const
Definition: EventBase.h:66
double respCorr(const DetId &id)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const Item * getValues(DetId fId, bool throwOnFail=true) const
int zside(DetId const &)
std::vector< HcalDDDRecConstants::HcalActiveLength > actHE
bool validHcal(const HcalDetId &id) const
const edm::InputTag hlTriggerResults_
std::vector< double > isolationR03_
std::vector< int > outerTrackHits_
#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
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
const Point & position() const
position
Definition: Vertex.h:109
double gainFactor(const edm::ESHandle< HcalDbService > &, const HcalDetId &id)
std::vector< bool > muon_global_
edm::EDGetTokenT< edm::TriggerResults > tok_trigRes_
std::vector< float > muon_trkKink
std::vector< double > hcalDepthEnergy_[depthMax_]
const edm::InputTag labelEERecHit_
std::vector< double > outerTrackPt_
std::map< DetId, double > corrValue_
U second(std::pair< T, U > const &p)
const edm::InputTag labelEBRecHit_
edm::EDGetTokenT< reco::VertexCollection > tok_Vtx_
void beginRun(edm::Run const &, edm::EventSetup const &) override
std::vector< double > globalTrckPt_
int depth() const
get the tower depth
Definition: HcalDetId.h:166
std::vector< int > tight_PixelHits_
int iEvent
Definition: GenABIO.cc:230
std::vector< bool > muon_tracker_
std::vector< double > ecalEnergy_
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
std::vector< int > numPixelLayers_
spr::propagatedTrackDirection propagateHCALBack(unsigned int thisTrk, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, const CaloGeometry *geo, const MagneticField *bField, bool debug=false)
std::vector< double > hcalActiveLengthHot_
std::vector< double > innerTracketa_
std::vector< double > globalTrckEta_
std::vector< double > isolationR04_
unsigned int size() const
Get number of paths stored.
std::vector< double > hcalDepthActiveLengthHot_[depthMax_]
std::vector< double > hcalDepthEnergyHotCorr_[depthMax_]
int getMinDepth(const int &itype, const int &ieta, const int &iphi, const int &zside) const
double getRZ(const int &subdet, const int &ieta, const int &depth) const
std::vector< double > ptGlob_
std::vector< unsigned int > ehcalDetId_
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float segmentCompatibility(const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
std::vector< double > hcalEnergy_
std::vector< int > hcal_iphi_
const std::string labelVtx_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< double > outerTrackEta_
bool isValid() const
Definition: HandleBase.h:74
double activeLength(const DetId &)
std::vector< bool > globalTrack_
std::vector< double > hcal1x1Energy_
edm::EDGetTokenT< reco::MuonCollection > tok_Muon_
bool goodCell(const HcalDetId &hcid, const reco::Track *pTrack, const CaloGeometry *geo, const MagneticField *bField)
std::vector< int > outerTrackRHits_
double ndof() const
Definition: Vertex.h:105
std::vector< bool > hcalHot_
JetCorrectorParameters corr
Definition: classes.h:5
const edm::InputTag labelHBHERecHit_
int k[5][pyjets_maxn]
std::vector< bool > muon_is_good_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
Definition: DetId.h:18
std::vector< double > chiGlobal_
std::vector< HcalDDDRecConstants::HcalActiveLength > actHB
std::vector< float > muon_segComp
int depth16HE(int ieta, int iphi)
std::vector< unsigned int > hcalDetId_
bool isFake() const
Definition: Vertex.h:72
std::vector< double > hcalDepthEnergyHot_[depthMax_]
std::vector< int > matchedStat_
std::vector< bool > matchedId_
void endRun(edm::Run const &, edm::EventSetup const &) override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< double > hcalDepthEnergyCorr_[depthMax_]
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< int > hcal_ieta_
std::vector< double > hcalDepthActiveLength_[depthMax_]
std::vector< double > hoEnergy_
int matchId(const HcalDetId &, const HcalDetId &)
const HcalDDDRecConstants * hdc_
edm::EventID id() const
Definition: EventBase.h:60
std::vector< double > phiGlob_
HLT enums.
HLTConfigProvider hltConfig_
T get() const
Definition: EventSetup.h:68
edm::EDGetTokenT< EcalRecHitCollection > tok_EB_
std::vector< bool > muon_is_medium_
std::vector< int > globalMuonHits_
std::vector< double > hcalActiveLength_
std::vector< double > chiTracker_
std::vector< double > hcalDepthChargeHot_[depthMax_]
const JetExtendedData & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
std::vector< bool > hcalDepthMatchHot_[depthMax_]
std::vector< HcalActiveLength > getThickActive(const int &type) const
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< int > trackerLayer_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< bool > outerTrack_
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
def check(config)
Definition: trackerTree.py:14
edm::EDGetTokenT< HBHERecHitCollection > tok_HBHE_
const std::string fileInCorr_
edm::Service< TFileService > fs
T const * product() const
Definition: ESHandle.h:84
std::vector< double > dxyTracker_
void energyHCALCell(HcalDetId detId, edm::Handle< T > &hits, std::vector< std::pair< double, int > > &energyCell, int maxDepth=1, double hbThr=-100, double heThr=-100, double hfThr=-100, double hoThr=-100, double tMin=-500, double tMax=500, int useRaw=0, int depthHE=3, bool debug=false)
constexpr double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
std::vector< double > ecal3x3Energy_
void setTopo(const HcalTopology *topo)
std::vector< double > tight_TransImpara_
const std::vector< std::string > triggers_
std::vector< double > energyMuon_
std::vector< double > etaGlob_
Definition: Run.h:44
bool isGoodVertex(const reco::Vertex &vtx)
std::vector< double > outerTrackPhi_
std::vector< double > innerTrackpt_
std::vector< bool > innerTrack_
std::vector< bool > hcalDepthMatch_[depthMax_]
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)
std::vector< float > muon_chi2LocalPosition