Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTPixlMBForAlignmentFilter
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Mika Huhtinen
0007  *
0008  */
0009 
0010 #include "HLTPixlMBForAlignmentFilter.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 HLTPixlMBForAlignmentFilter::HLTPixlMBForAlignmentFilter(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       min_isol_(iConfig.getParameter<double>("MinIsol"))
0038 
0039 {
0040   pixlToken_ = consumes<reco::RecoChargedCandidateCollection>(pixlTag_);
0041   LogDebug("") << "MinPt cut " << min_Pt_ << "pixl: " << pixlTag_.encode();
0042   LogDebug("") << "Requesting : " << min_trks_ << " tracks from same vertex ";
0043   LogDebug("") << "Requesting tracks from same vertex eta-phi separation by " << min_sep_;
0044   LogDebug("") << "Requesting track to be isolated within cone of " << min_isol_;
0045 }
0046 
0047 HLTPixlMBForAlignmentFilter::~HLTPixlMBForAlignmentFilter() = default;
0048 
0049 void HLTPixlMBForAlignmentFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050   edm::ParameterSetDescription desc;
0051   makeHLTFilterDescription(desc);
0052   desc.add<edm::InputTag>("pixlTag", edm::InputTag("hltPixelCands"));
0053   desc.add<double>("MinPt", 5.0);
0054   desc.add<unsigned int>("MinTrks", 2);
0055   desc.add<double>("MinSep", 1.0);
0056   desc.add<double>("MinIsol", 0.05);
0057   descriptions.add("hltPixlMBForAlignmentFilter", desc);
0058 }
0059 
0060 //
0061 // member functions
0062 //
0063 
0064 // ------------ method called to produce the data  ------------
0065 bool HLTPixlMBForAlignmentFilter::hltFilter(edm::Event& iEvent,
0066                                             const edm::EventSetup& iSetup,
0067                                             trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0068   using namespace std;
0069   using namespace edm;
0070   using namespace reco;
0071   using namespace trigger;
0072 
0073   // All HLT filters must create and fill an HLT filter object,
0074   // recording any reconstructed physics objects satisfying (or not)
0075   // this HLT filter, and place it in the Event.
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   vector<double> etastore;
0086   vector<double> phistore;
0087   vector<double> ptstore;
0088   vector<int> itstore;
0089   bool accept = false;
0090   auto apixl(tracks->begin());
0091   auto epixl(tracks->end());
0092   int itrk = 0;
0093   if (tracks->size() >= min_trks_) {
0094     etastore.clear();
0095     phistore.clear();
0096     itstore.clear();
0097     for (auto ipixl = apixl; ipixl != epixl; ipixl++) {
0098       const double& etatrk1 = ipixl->momentum().eta();
0099       const double& phitrk1 = ipixl->momentum().phi();
0100       const double& pttrk1 = ipixl->pt();
0101       if (pttrk1 > min_Pt_) {
0102         //       the *store-vectors store the tracks above pt-cut
0103         //       itstore is the position in the original collection
0104         etastore.push_back(etatrk1);
0105         phistore.push_back(phitrk1);
0106         ptstore.push_back(pttrk1);
0107         itstore.push_back(itrk);
0108       }
0109       itrk++;
0110     }
0111     //   locisol is the position in the *store vectors
0112     vector<int> locisol;
0113     if (itstore.size() > 1) {
0114       // now check that tracks are isolated
0115       locisol.clear();
0116       for (unsigned int i = 0; i < itstore.size(); i++) {
0117         int nincone = 0;
0118         //       check isolation wrt ALL tracks, not only those above ptcut
0119         for (auto ipixl = apixl; ipixl != epixl; ipixl++) {
0120           double phidist = std::abs(phistore.at(i) - ipixl->momentum().phi());
0121           double etadist = std::abs(etastore.at(i) - ipixl->momentum().eta());
0122           double trkdist = sqrt(phidist * phidist + etadist * etadist);
0123           if (trkdist < min_isol_)
0124             nincone++;
0125         }
0126         //       the check above always find the track itself, so nincone never should be 0
0127         if (nincone < 2)
0128           locisol.push_back(i);
0129       }
0130     }
0131     //   now check that the selected tracks have enough mutual separation
0132     vector<int> itsep;
0133     for (unsigned int i = 0; i < locisol.size(); i++) {
0134       //     check for each so far selected track...
0135       itsep.clear();
0136       itsep.push_back(locisol.at(i));
0137       for (unsigned int j = i + 1; j < locisol.size(); j++) {
0138         //       ...if it is sufficiently separated from other selectad tracks...
0139         double phidist = phistore.at(locisol.at(i)) - phistore.at(locisol.at(j));
0140         double etadist = etastore.at(locisol.at(i)) - etastore.at(locisol.at(j));
0141         double dist = sqrt(phidist * phidist + etadist * etadist);
0142         if (dist > min_sep_) {
0143           if (itsep.size() == 1) {
0144             itsep.push_back(locisol.at(j));
0145           } else {
0146             bool is_separated = true;
0147             for (int k : itsep) {
0148               //             ...and the other ones, that are on the 'final acceptance' list already, if min_trks_ > 2
0149               double phisep = phistore.at(k) - phistore.at(locisol.at(j));
0150               double etasep = etastore.at(k) - etastore.at(locisol.at(j));
0151               double sep = sqrt(phisep * phisep + etasep * etasep);
0152               if (sep < min_sep_) {
0153                 //               this one was no good, too close to some other already accepted
0154                 is_separated = false;
0155                 break;
0156               }
0157             }
0158             if (is_separated)
0159               itsep.push_back(locisol.at(j));
0160           }
0161         }
0162         if (itsep.size() >= min_trks_) {
0163           accept = true;
0164           break;
0165         }
0166       }
0167       if (accept) {
0168         break;
0169       }
0170     }
0171     // At this point we have the indices of the accepted tracks stored in itstore
0172     // we now move them to the filterproduct
0173 
0174     if (accept) {
0175       for (int ipos : itsep) {
0176         int iaddr = itstore.at(ipos);
0177         filterproduct.addObject(TriggerTrack, RecoChargedCandidateRef(tracks, iaddr));
0178       }
0179       // std::cout << "Accept this event " << std::endl;
0180     }
0181   }
0182 
0183   //  LogDebug("") << "Number of pixel-track objects accepted:"
0184   //               << " " << npixl_tot;
0185 
0186   // return with final filter decision
0187   return accept;
0188 }
0189 
0190 // declare this class as a framework plugin
0191 #include "FWCore/Framework/interface/MakerMacros.h"
0192 DEFINE_FWK_MODULE(HLTPixlMBForAlignmentFilter);