Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 #include "HLTSingleVertexPixelTrackFilter.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 
0011 #include "DataFormats/Common/interface/Handle.h"
0012 #include "DataFormats/Common/interface/RefToBase.h"
0013 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0014 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0015 #include "DataFormats/TrackReco/interface/Track.h"
0016 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 // constructors and destructor
0021 //
0022 
0023 HLTSingleVertexPixelTrackFilter::HLTSingleVertexPixelTrackFilter(const edm::ParameterSet& iConfig)
0024     : HLTFilter(iConfig),
0025       pixelVerticesTag_(iConfig.getParameter<edm::InputTag>("vertexCollection")),
0026       pixelTracksTag_(iConfig.getParameter<edm::InputTag>("trackCollection")),
0027       min_Pt_(iConfig.getParameter<double>("MinPt")),
0028       max_Pt_(iConfig.getParameter<double>("MaxPt")),
0029       max_Eta_(iConfig.getParameter<double>("MaxEta")),
0030       max_Vz_(iConfig.getParameter<double>("MaxVz")),
0031       min_trks_(iConfig.getParameter<int>("MinTrks")),
0032       min_sep_(iConfig.getParameter<double>("MinSep")) {
0033   pixelVerticesToken_ = consumes<reco::VertexCollection>(pixelVerticesTag_);
0034   pixelTracksToken_ = consumes<reco::RecoChargedCandidateCollection>(pixelTracksTag_);
0035 }
0036 
0037 HLTSingleVertexPixelTrackFilter::~HLTSingleVertexPixelTrackFilter() = default;
0038 
0039 void HLTSingleVertexPixelTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040   edm::ParameterSetDescription desc;
0041   makeHLTFilterDescription(desc);
0042   desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVerticesForMinBias"));
0043   desc.add<edm::InputTag>("trackCollection", edm::InputTag("hltPixelCands"));
0044   desc.add<double>("MinPt", 0.2);
0045   desc.add<double>("MaxPt", 10000.0);
0046   desc.add<double>("MaxEta", 1.0);
0047   desc.add<double>("MaxVz", 10.0);
0048   desc.add<int>("MinTrks", 30);
0049   desc.add<double>("MinSep", 0.12);
0050   descriptions.add("hltSingleVertexPixelTrackFilter", desc);
0051 }
0052 
0053 //
0054 // member functions
0055 //
0056 
0057 // ------------ method called to produce the data  ------------
0058 bool HLTSingleVertexPixelTrackFilter::hltFilter(edm::Event& iEvent,
0059                                                 const edm::EventSetup& iSetup,
0060                                                 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0061   // All HLT filters must create and fill an HLT filter object,
0062   // recording any reconstructed physics objects satisfying (or not)
0063   // this HLT filter, and place it in the Event.
0064   if (saveTags())
0065     filterproduct.addCollectionTag(pixelTracksTag_);
0066 
0067   // Ref to Candidate object to be recorded in filter object
0068   edm::Ref<reco::RecoChargedCandidateCollection> candref;
0069 
0070   // Specific filter code
0071   bool accept = false;
0072 
0073   int nTrackCandidate = 0;
0074 
0075   // get hold of products from Event
0076   edm::Handle<reco::VertexCollection> vertexCollection;
0077   iEvent.getByToken(pixelVerticesToken_, vertexCollection);
0078   if (vertexCollection.isValid()) {
0079     const reco::VertexCollection* vertices = vertexCollection.product();
0080     int npixelvertices = vertices->size();
0081     if (npixelvertices != 0) {
0082       double vzmax = 0;
0083       int nmax = 0;
0084       reco::VertexCollection::const_iterator verticesItr;
0085       for (verticesItr = vertices->begin(); verticesItr != vertices->end(); ++verticesItr) {
0086         int ntracksize = verticesItr->tracksSize();
0087         double vz = verticesItr->z();
0088         if (fabs(vz) > max_Vz_)
0089           continue;
0090         if (ntracksize > nmax) {
0091           vzmax = vz;
0092           nmax = ntracksize;
0093         }
0094       }
0095 
0096       edm::Handle<reco::RecoChargedCandidateCollection> trackCollection;
0097       iEvent.getByToken(pixelTracksToken_, trackCollection);
0098       if (trackCollection.isValid()) {
0099         const reco::RecoChargedCandidateCollection* tracks = trackCollection.product();
0100         reco::RecoChargedCandidateCollection::const_iterator tracksItr;
0101         int icount = -1;
0102         for (tracksItr = tracks->begin(); tracksItr != tracks->end(); ++tracksItr) {
0103           icount++;
0104           double eta = tracksItr->eta();
0105           if (fabs(eta) > max_Eta_)
0106             continue;
0107           double pt = tracksItr->pt();
0108           if (pt < min_Pt_ || pt > max_Pt_)
0109             continue;
0110           double vz = tracksItr->vz();
0111           if (fabs(vz - vzmax) > min_sep_)
0112             continue;
0113 
0114           candref = edm::Ref<reco::RecoChargedCandidateCollection>(trackCollection, icount);
0115           filterproduct.addObject(trigger::TriggerTrack, candref);
0116           nTrackCandidate++;
0117         }
0118       }
0119     }
0120   }
0121 
0122   accept = (nTrackCandidate >= min_trks_);
0123 
0124   return accept;
0125 }
0126 
0127 DEFINE_FWK_MODULE(HLTSingleVertexPixelTrackFilter);