Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:42

0001 /** \class HLTPixlMBFilt
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Mika Huhtinen
0007  *
0008  */
0009 
0010 #include "HLTPixlMBFilt.h"
0011 
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 
0018 #include "DataFormats/Common/interface/Handle.h"
0019 
0020 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0021 
0022 #include "DataFormats/TrackReco/interface/Track.h"
0023 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0024 
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 
0027 //
0028 // constructors and destructor
0029 //
0030 
0031 HLTPixlMBFilt::HLTPixlMBFilt(const edm::ParameterSet& iConfig)
0032     : HLTFilter(iConfig),
0033       pixlTag_(iConfig.getParameter<edm::InputTag>("pixlTag")),
0034       min_Pt_(iConfig.getParameter<double>("MinPt")),
0035       min_trks_(iConfig.getParameter<unsigned int>("MinTrks")),
0036       min_sep_(iConfig.getParameter<double>("MinSep"))
0037 
0038 {
0039   pixlToken_ = consumes<reco::RecoChargedCandidateCollection>(pixlTag_);
0040   LogDebug("") << "MinPt cut " << min_Pt_ << "pixl: " << pixlTag_.encode();
0041   LogDebug("") << "Requesting : " << min_trks_ << " tracks from same vertex ";
0042   LogDebug("") << "Requesting tracks from same vertex eta-phi separation by " << min_sep_;
0043 }
0044 
0045 HLTPixlMBFilt::~HLTPixlMBFilt() = default;
0046 
0047 void HLTPixlMBFilt::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0048   edm::ParameterSetDescription desc;
0049   makeHLTFilterDescription(desc);
0050   desc.add<edm::InputTag>("pixlTag", edm::InputTag("hltPixelCands"));
0051   desc.add<double>("MinPt", 0.);
0052   desc.add<unsigned int>("MinTrks", 2);
0053   desc.add<double>("MinSep", 1.);
0054   descriptions.add("hltPixlMBFilt", desc);
0055 }
0056 
0057 //
0058 // member functions
0059 //
0060 
0061 // ------------ method called to produce the data  ------------
0062 bool HLTPixlMBFilt::hltFilter(edm::Event& iEvent,
0063                               const edm::EventSetup& iSetup,
0064                               trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0065   using namespace std;
0066   using namespace edm;
0067   using namespace reco;
0068   using namespace trigger;
0069 
0070   // All HLT filters must create and fill an HLT filter object,
0071   // recording any reconstructed physics objects satisfying (or not)
0072   // this HLT filter, and place it in the Event.
0073   if (saveTags()) {
0074     filterproduct.addCollectionTag(pixlTag_);
0075   }
0076 
0077   // Specific filter code
0078 
0079   // get hold of products from Event
0080 
0081   Handle<RecoChargedCandidateCollection> tracks;
0082   iEvent.getByToken(pixlToken_, tracks);
0083 
0084   // pixel tracks
0085   int npixl_tot = 0;
0086   vector<double> etastore;
0087   vector<double> phistore;
0088   vector<int> itstore;
0089   bool accept;
0090   auto apixl(tracks->begin());
0091   auto epixl(tracks->end());
0092   RecoChargedCandidateCollection::const_iterator ipixl, jpixl;
0093   unsigned int nsame_vtx = 0;
0094   int itrk = -1;
0095   if (tracks->size() >= min_trks_) {
0096     for (ipixl = apixl; ipixl != epixl; ipixl++) {
0097       if (ipixl->pt() < min_Pt_)
0098         continue;
0099       itrk++;
0100       const double& ztrk1 = ipixl->vz();
0101       const double& etatrk1 = ipixl->momentum().eta();
0102       const double& phitrk1 = ipixl->momentum().phi();
0103       nsame_vtx = 1;
0104       etastore.clear();
0105       phistore.clear();
0106       itstore.clear();
0107       etastore.push_back(etatrk1);
0108       phistore.push_back(phitrk1);
0109       itstore.push_back(itrk);
0110       if (fabs(ztrk1) < 15.0) {
0111         //  check this track against all others to see if others start from same point
0112         int jtrk = -1;
0113         for (jpixl = apixl; jpixl != epixl; jpixl++) {
0114           if (jpixl->pt() < min_Pt_)
0115             continue;
0116           jtrk++;
0117           if (jpixl == ipixl)
0118             continue;
0119           const double& ztrk2 = jpixl->vz();
0120           const double& etatrk2 = jpixl->momentum().eta();
0121           const double& phitrk2 = jpixl->momentum().phi();
0122           double eta_dist = etatrk2 - etatrk1;
0123           double phi_dist = phitrk2 - phitrk1;
0124           double etaphi_dist = sqrt(eta_dist * eta_dist + phi_dist * phi_dist);
0125           if (fabs(ztrk2 - ztrk1) < 1.0 && etaphi_dist > min_sep_) {
0126             if (min_trks_ <= 2 || itstore.size() <= 1) {
0127               etastore.push_back(etatrk2);
0128               phistore.push_back(phitrk2);
0129               itstore.push_back(jtrk);
0130               nsame_vtx++;
0131             } else {
0132               // check also separation to already found 'second' tracks
0133               LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
0134               bool isok = true;
0135               for (unsigned int k = 1; k < itstore.size(); k++) {
0136                 eta_dist = etatrk2 - etastore.at(k);
0137                 phi_dist = phitrk2 - phistore.at(k);
0138                 etaphi_dist = sqrt(eta_dist * eta_dist + phi_dist * phi_dist);
0139                 if (etaphi_dist < min_sep_) {
0140                   isok = false;
0141                   break;
0142                 }
0143               }
0144               if (isok) {
0145                 etastore.push_back(etatrk2);
0146                 phistore.push_back(phitrk2);
0147                 itstore.push_back(jtrk);
0148                 nsame_vtx++;
0149               }
0150             }
0151           }
0152           if (nsame_vtx >= min_trks_)
0153             break;
0154         }
0155       }
0156       npixl_tot++;
0157 
0158       if (nsame_vtx >= min_trks_)
0159         break;
0160     }
0161 
0162     //   final filter decision:
0163     //   request at least min_trks_ tracks compatible with vertex-region
0164     accept = (nsame_vtx >= min_trks_);
0165 
0166   } else {
0167     accept = false;
0168   }
0169 
0170   // At this point we have the indices of the accepted tracks stored in itstore
0171   // we now move them to the filterproduct
0172 
0173   if (accept) {
0174     for (int iaddr : itstore) {
0175       filterproduct.addObject(TriggerTrack, RecoChargedCandidateRef(tracks, iaddr));
0176     }
0177   }
0178 
0179   LogDebug("") << "Number of pixel-track objects accepted:"
0180                << " " << npixl_tot;
0181 
0182   // return with final filter decision
0183   return accept;
0184 }
0185 
0186 // declare this class as a framework plugin
0187 #include "FWCore/Framework/interface/MakerMacros.h"
0188 DEFINE_FWK_MODULE(HLTPixlMBFilt);