Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:15

0001 // -*- C++ -*-
0002 //
0003 // Package:    DQMOffline/CalibCalo
0004 // Class:      DQMHcalIsoTrackAlCaReco
0005 //
0006 /**\class DQMHcalIsoTrackAlCaReco DQMHcalIsoTrackAlCaReco.cc
0007  DQMOffline/CalibCalo/src/DQMHcalIsoTrackAlCaReco.cc
0008 
0009  Description: <one line class summary>
0010 
0011  Implementation:
0012      <Notes on implementation>
0013 */
0014 //
0015 // Original Author:  Grigory SAFRONOV
0016 //         Created:  Tue Oct  14 16:10:31 CEST 2008
0017 //         Modified: Tue Mar   3 16:10:31 CEST 2015
0018 //
0019 //
0020 
0021 // system include files
0022 #include <cmath>
0023 #include <fstream>
0024 #include <vector>
0025 
0026 // user include files
0027 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0028 #include "DataFormats/HcalIsolatedTrack/interface/HcalIsolatedTrackCandidate.h"
0029 #include "DataFormats/HcalIsolatedTrack/interface/HcalIsolatedTrackCandidateFwd.h"
0030 #include "DataFormats/Math/interface/deltaR.h"
0031 
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/Frameworkfwd.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 #include "FWCore/ServiceRegistry/interface/Service.h"
0038 
0039 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0040 #include "DQMServices/Core/interface/DQMStore.h"
0041 
0042 class DQMHcalIsoTrackAlCaReco : public DQMEDAnalyzer {
0043 public:
0044   explicit DQMHcalIsoTrackAlCaReco(const edm::ParameterSet &);
0045   ~DQMHcalIsoTrackAlCaReco() override;
0046 
0047   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0048   void analyze(const edm::Event &, const edm::EventSetup &) override;
0049 
0050 private:
0051   std::string folderName_;
0052   std::vector<std::string> l1FilterTag_, hltFilterTag_;
0053   std::vector<int> type_;
0054   edm::InputTag labelTrigger_, labelTrack_;
0055   edm::EDGetTokenT<trigger::TriggerEvent> tokTrigger_;
0056   edm::EDGetTokenT<reco::HcalIsolatedTrackCandidateCollection> tokTrack_;
0057 
0058   double pThr_;
0059 
0060   std::vector<MonitorElement *> hL1Pt_, hL1Eta_, hL1phi_;
0061   std::vector<MonitorElement *> hHltP_, hHltEta_, hHltPhi_;
0062   MonitorElement *hL3Dr_, *hL3Rat_;
0063   std::vector<MonitorElement *> hOffP_;
0064   MonitorElement *hMaxP_, *hEnEcal_, *hIeta_, *hIphi_;
0065 
0066   int nTotal_, nHLTaccepts_;
0067   std::vector<double> etaRange_;
0068   std::vector<unsigned int> indexH_;
0069   std::vector<bool> ifL3_;
0070 };
0071 
0072 DQMHcalIsoTrackAlCaReco::DQMHcalIsoTrackAlCaReco(const edm::ParameterSet &iConfig) {
0073   folderName_ = iConfig.getParameter<std::string>("FolderName");
0074   l1FilterTag_ = iConfig.getParameter<std::vector<std::string>>("L1FilterLabel");
0075   hltFilterTag_ = iConfig.getParameter<std::vector<std::string>>("HltFilterLabels");
0076   type_ = iConfig.getParameter<std::vector<int>>("TypeFilter");
0077   labelTrigger_ = iConfig.getParameter<edm::InputTag>("TriggerLabel");
0078   labelTrack_ = iConfig.getParameter<edm::InputTag>("TracksLabel");
0079   pThr_ = iConfig.getUntrackedParameter<double>("pThrL3", 0);
0080 
0081   nTotal_ = nHLTaccepts_ = 0;
0082   tokTrigger_ = consumes<trigger::TriggerEvent>(labelTrigger_);
0083   tokTrack_ = consumes<reco::HcalIsolatedTrackCandidateCollection>(labelTrack_);
0084   LogDebug("HcalIsoTrack") << "Folder " << folderName_ << " Input Tag for Trigger " << labelTrigger_ << " track "
0085                            << labelTrack_ << " threshold " << pThr_ << " with " << l1FilterTag_.size()
0086                            << " level 1 and " << hltFilterTag_.size() << " hlt filter tags"
0087                            << "\n";
0088   for (unsigned int k = 0; k < l1FilterTag_.size(); ++k)
0089     LogDebug("HcalIsoTrack") << "L1FilterTag[" << k << "] " << l1FilterTag_[k] << "\n";
0090   for (unsigned int k = 0; k < hltFilterTag_.size(); ++k)
0091     LogDebug("HcalIsoTrack") << "HLTFilterTag[" << k << "] " << hltFilterTag_[k] << "\n";
0092 }
0093 
0094 DQMHcalIsoTrackAlCaReco::~DQMHcalIsoTrackAlCaReco() {}
0095 
0096 void DQMHcalIsoTrackAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0097   nTotal_++;
0098   bool accept(false);
0099   edm::Handle<trigger::TriggerEvent> trEv;
0100   iEvent.getByToken(tokTrigger_, trEv);
0101 
0102   edm::Handle<reco::HcalIsolatedTrackCandidateCollection> recoIsoTracks;
0103   iEvent.getByToken(tokTrack_, recoIsoTracks);
0104   LogDebug("HcalIsoTrack") << "Gets Trigger information with " << trEv.isValid() << " and offline tracks with "
0105                            << recoIsoTracks.isValid() << "\n";
0106 
0107   if (trEv.isValid()) {
0108     const trigger::TriggerObjectCollection &TOCol(trEv->getObjects());
0109     const trigger::size_type nFilt(trEv->sizeFilters());
0110     // plots for L1 trigger
0111     for (unsigned int k = 0; k < l1FilterTag_.size(); ++k) {
0112       trigger::Keys KEYSl1;
0113       double etaTrigl1(-10000), phiTrigl1(-10000), ptMaxl1(0);
0114       for (trigger::size_type iFilt = 0; iFilt != nFilt; iFilt++) {
0115         LogDebug("HcalIsoTrack") << trEv->filterTag(iFilt).label() << " find for " << l1FilterTag_[k] << " gives "
0116                                  << (trEv->filterTag(iFilt).label()).find(l1FilterTag_[k]) << "\n";
0117         if ((trEv->filterTag(iFilt).label()).find(l1FilterTag_[k]) != std::string::npos) {
0118           KEYSl1 = trEv->filterKeys(iFilt);
0119           trigger::size_type nRegl1 = KEYSl1.size();
0120           LogDebug("HcalIsoTrack") << "# of objects " << nRegl1 << "\n";
0121           for (trigger::size_type iReg = 0; iReg < nRegl1; iReg++) {
0122             const trigger::TriggerObject &TObj(TOCol[KEYSl1[iReg]]);
0123             LogDebug("HcalIsoTrack") << "Object[" << iReg << "] with pt " << TObj.pt() << " " << TObj.eta() << " "
0124                                      << TObj.phi() << "\n";
0125             if (TObj.pt() > ptMaxl1) {
0126               etaTrigl1 = TObj.eta();
0127               phiTrigl1 = TObj.phi();
0128               ptMaxl1 = TObj.pt();
0129             }
0130           }
0131         }
0132       }
0133       LogDebug("HcalIsoTrack") << "For L1 trigger type " << k << " pt " << ptMaxl1 << " eta " << etaTrigl1 << " phi "
0134                                << phiTrigl1 << "\n";
0135       if (ptMaxl1 > 0) {
0136         hL1Pt_[k]->Fill(ptMaxl1);
0137         hL1Eta_[k]->Fill(etaTrigl1);
0138         hL1phi_[k]->Fill(phiTrigl1);
0139       }
0140     }
0141     // Now make plots for hlt objects
0142     trigger::Keys KEYS;
0143     for (unsigned l = 0; l < hltFilterTag_.size(); l++) {
0144       for (trigger::size_type iFilt = 0; iFilt != nFilt; iFilt++) {
0145         LogDebug("HcalIsoTrack") << trEv->filterTag(iFilt).label() << " find for " << hltFilterTag_[l] << " gives "
0146                                  << (trEv->filterTag(iFilt).label()).find(hltFilterTag_[l]) << "\n";
0147         if ((trEv->filterTag(iFilt).label()).find(hltFilterTag_[l]) != std::string::npos) {
0148           KEYS = trEv->filterKeys(iFilt);
0149           trigger::size_type nReg = KEYS.size();
0150           LogDebug("HcalIsoTrack") << "# of objects for HLT " << nReg << "\n";
0151           // checks with IsoTrack trigger results
0152           for (trigger::size_type iReg = 0; iReg < nReg; iReg++) {
0153             const trigger::TriggerObject &TObj(TOCol[KEYS[iReg]]);
0154             LogDebug("HcalIsoTrack") << "HLT Filter Tag " << l << " trigger " << iFilt << " object " << iReg << " p "
0155                                      << TObj.p() << " pointer " << indexH_[l] << ":" << hHltP_[indexH_[l]] << ":"
0156                                      << hHltEta_[indexH_[l]] << ":" << hHltPhi_[indexH_[l]] << "\n";
0157             if (TObj.p() > pThr_) {
0158               hHltP_[indexH_[l]]->Fill(TObj.p());
0159               hHltEta_[indexH_[l]]->Fill(TObj.eta());
0160               hHltPhi_[indexH_[l]]->Fill(TObj.phi());
0161               if (ifL3_[l])
0162                 accept = true;
0163               if (recoIsoTracks.isValid() && ifL3_[l]) {
0164                 double minRecoL3dist(1000), pt(1000);
0165                 reco::HcalIsolatedTrackCandidateCollection::const_iterator mrtr;
0166                 for (mrtr = recoIsoTracks->begin(); mrtr != recoIsoTracks->end(); mrtr++) {
0167                   double R = deltaR(mrtr->eta(), mrtr->phi(), TObj.eta(), TObj.phi());
0168                   if (R < minRecoL3dist) {
0169                     minRecoL3dist = R;
0170                     pt = mrtr->pt();
0171                   }
0172                 }
0173                 LogDebug("HcalIsoTrack") << "Minimum R " << minRecoL3dist << " pt " << pt << ":" << TObj.pt() << "\n";
0174                 hL3Dr_->Fill(minRecoL3dist);
0175                 if (minRecoL3dist < 0.02)
0176                   hL3Rat_->Fill(TObj.pt() / pt);
0177               }
0178             }
0179           }
0180         }
0181       }
0182     }
0183   }
0184 
0185   // general distributions
0186   if (recoIsoTracks.isValid()) {
0187     for (reco::HcalIsolatedTrackCandidateCollection::const_iterator itr = recoIsoTracks->begin();
0188          itr != recoIsoTracks->end();
0189          itr++) {
0190       hMaxP_->Fill(itr->maxP());
0191       hEnEcal_->Fill(itr->energyEcal());
0192       std::pair<int, int> etaphi = itr->towerIndex();
0193       hIeta_->Fill(etaphi.first);
0194       hIphi_->Fill(etaphi.second);
0195       LogDebug("HcalIsoTrack") << "Reco track p " << itr->p() << " eta|phi " << etaphi.first << "|" << etaphi.second
0196                                << " maxP " << itr->maxP() << " EcalE " << itr->energyEcal() << " pointers " << hHltP_[3]
0197                                << ":" << hHltEta_[3] << ":" << hHltPhi_[3] << "\n";
0198       if (itr->p() >= pThr_) {
0199         hHltP_[3]->Fill(itr->p());
0200         hHltEta_[3]->Fill(itr->eta());
0201         hHltPhi_[3]->Fill(itr->phi());
0202       }
0203       double etaAbs = std::abs(itr->eta());
0204       hOffP_[0]->Fill(itr->p());
0205       for (unsigned int l = 1; l < etaRange_.size(); l++) {
0206         if (etaAbs >= etaRange_[l - 1] && etaAbs < etaRange_[l]) {
0207           LogDebug("HcalIsoTrack") << "Range " << l << " p " << itr->p() << " pointer " << hOffP_[l];
0208           hOffP_[l]->Fill(itr->p());
0209           break;
0210         }
0211       }
0212     }
0213   }
0214 
0215   if (accept)
0216     nHLTaccepts_++;
0217   LogDebug("HcalIsoTrack") << "Accept " << accept << "\n";
0218 }
0219 
0220 void DQMHcalIsoTrackAlCaReco::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0221   iBooker.setCurrentFolder(folderName_);
0222   LogDebug("HcalIsoTrack") << "Set the folder to " << folderName_ << "\n";
0223   char name[100], title[200];
0224   for (unsigned int k = 0; k < l1FilterTag_.size(); ++k) {
0225     sprintf(name, "hp%s", l1FilterTag_[k].c_str());
0226     sprintf(title, "p_T of L1 object for %s", l1FilterTag_[k].c_str());
0227     hL1Pt_.push_back(iBooker.book1D(name, title, 1000, 0, 1000));
0228     hL1Pt_[k]->setAxisTitle("p_T (GeV)", 1);
0229     sprintf(name, "heta%s", l1FilterTag_[k].c_str());
0230     sprintf(title, "#eta of L1 object for %s", l1FilterTag_[k].c_str());
0231     hL1Eta_.push_back(iBooker.book1D(name, title, 100, -2.5, 2.5));
0232     hL1Eta_[k]->setAxisTitle("#eta", 1);
0233     sprintf(name, "hphi%s", l1FilterTag_[k].c_str());
0234     sprintf(title, "#phi of L1 object for %s", l1FilterTag_[k].c_str());
0235     hL1phi_.push_back(iBooker.book1D(name, title, 100, -3.2, 3.2));
0236     hL1phi_[k]->setAxisTitle("#phi", 1);
0237   }
0238 
0239   std::string types[4] = {"L2", "L2x", "L3", "Off"};
0240   for (unsigned int l = 0; l < 4; l++) {
0241     sprintf(name, "hp%s", types[l].c_str());
0242     sprintf(title, "Momentum of %s object", types[l].c_str());
0243     hHltP_.push_back(iBooker.book1D(name, title, 200, 0, 1000));
0244     hHltP_[l]->setAxisTitle("p (GeV)", 1);
0245     sprintf(name, "heta%s", types[l].c_str());
0246     sprintf(title, "#eta of %s object", types[l].c_str());
0247     hHltEta_.push_back(iBooker.book1D(name, title, 16, -2, 2));
0248     hHltEta_[l]->setAxisTitle("#eta", 1);
0249     sprintf(name, "hphi%s", types[l].c_str());
0250     sprintf(title, "#phi of %s object", types[l].c_str());
0251     hHltPhi_.push_back(iBooker.book1D(name, title, 16, -3.2, 3.2));
0252     hHltPhi_[l]->setAxisTitle("#phi", 1);
0253   }
0254   sprintf(title, "Distance of offline track from L3 object");
0255   hL3Dr_ = (iBooker.book1D("hDRL3", title, 40, 0, 0.2));
0256   hL3Dr_->setAxisTitle("R(#eta,#phi)", 1);
0257   sprintf(title, "Ratio of p L3/Offline");
0258   hL3Rat_ = (iBooker.book1D("hRatL3", title, 500, 0, 3));
0259   indexH_.clear();
0260   ifL3_.clear();
0261   for (unsigned int l = 0; l < hltFilterTag_.size(); l++) {
0262     unsigned int indx = (type_[l] >= 0 && type_[l] < 3) ? type_[l] : 0;
0263     indexH_.push_back(indx);
0264     ifL3_.push_back(indx == 2);
0265     LogDebug("HcalIsoTrack") << "Filter[" << l << "] " << hltFilterTag_[l] << " type " << type_[l] << " index "
0266                              << indexH_[l] << " L3? " << ifL3_[l] << "\n";
0267   }
0268 
0269   double etaV[6] = {0.0, 0.5, 1.0, 1.5, 2.0, 2.5};
0270   for (unsigned int k = 0; k < 6; ++k) {
0271     sprintf(name, "hOffP%d", k);
0272     if (k == 0) {
0273       sprintf(title, "p of AlCaReco object (All)");
0274     } else {
0275       sprintf(title, "p of AlCaReco object (%3.1f < |#eta| < %3.1f)", etaV[k - 1], etaV[k]);
0276     }
0277     etaRange_.push_back(etaV[k]);
0278     hOffP_.push_back(iBooker.book1D(name, title, 1000, 0, 1000));
0279     hOffP_[k]->setAxisTitle("E (GeV)", 1);
0280   }
0281   hMaxP_ = iBooker.book1D("hChgIsol", "Energy for charge isolation", 110, -10, 100);
0282   hMaxP_->setAxisTitle("p (GeV)", 1);
0283   hEnEcal_ = iBooker.book1D("hEnEcal", "Energy in ECAL", 100, 0, 20);
0284   hEnEcal_->setAxisTitle("E (GeV)", 1);
0285   hIeta_ = iBooker.book1D("hIEta", "i#eta for HCAL tower", 90, -45, 45);
0286   hIeta_->setAxisTitle("i#eta", 1);
0287   hIphi_ = iBooker.book1D("hIPhi", "i#phi for HCAL tower", 72, 0, 72);
0288   hIphi_->setAxisTitle("i#phi", 1);
0289 }
0290 
0291 DEFINE_FWK_MODULE(DQMHcalIsoTrackAlCaReco);