CMS 3D CMS Logo

HLTSingleVertexPixelTrackFilter.cc
Go to the documentation of this file.
1 
3 
10 
17 
19 
20 // constructors and destructor
21 //
22 
24  pixelVerticesTag_ (iConfig.getParameter<edm::InputTag>("vertexCollection")),
25  pixelTracksTag_ (iConfig.getParameter<edm::InputTag>("trackCollection")),
26  min_Pt_ (iConfig.getParameter<double>("MinPt")),
27  max_Pt_ (iConfig.getParameter<double>("MaxPt")),
28  max_Eta_ (iConfig.getParameter<double>("MaxEta")),
29  max_Vz_ (iConfig.getParameter<double>("MaxVz")),
30  min_trks_ (iConfig.getParameter<int>("MinTrks")),
31  min_sep_ (iConfig.getParameter<double>("MinSep"))
32 {
33  pixelVerticesToken_ = consumes<reco::VertexCollection>(pixelVerticesTag_);
34  pixelTracksToken_ = consumes<reco::RecoChargedCandidateCollection>(pixelTracksTag_);
35 }
36 
38 
39 void
43  desc.add<edm::InputTag>("vertexCollection",edm::InputTag("hltPixelVerticesForMinBias"));
44  desc.add<edm::InputTag>("trackCollection",edm::InputTag("hltPixelCands"));
45  desc.add<double>("MinPt",0.2);
46  desc.add<double>("MaxPt",10000.0);
47  desc.add<double>("MaxEta",1.0);
48  desc.add<double>("MaxVz",10.0);
49  desc.add<int>("MinTrks",30);
50  desc.add<double>("MinSep",0.12);
51  descriptions.add("hltSingleVertexPixelTrackFilter",desc);
52 }
53 
54 //
55 // member functions
56 //
57 
58 // ------------ method called to produce the data ------------
60 {
61  // All HLT filters must create and fill an HLT filter object,
62  // recording any reconstructed physics objects satisfying (or not)
63  // this HLT filter, and place it in the Event.
64  if (saveTags())
65  filterproduct.addCollectionTag(pixelTracksTag_);
66 
67  // Ref to Candidate object to be recorded in filter object
69 
70  // Specific filter code
71  bool accept = false;
72 
73  int nTrackCandidate = 0;
74 
75  // get hold of products from Event
77  iEvent.getByToken( pixelVerticesToken_, vertexCollection );
78  if(vertexCollection.isValid())
79  {
80  const reco::VertexCollection * vertices = vertexCollection.product();
81  int npixelvertices = vertices->size();
82  if (npixelvertices!=0)
83  {
84  double vzmax = 0;
85  int nmax = 0;
86  reco::VertexCollection::const_iterator verticesItr;
87  for (verticesItr=vertices->begin(); verticesItr!=vertices->end(); ++verticesItr)
88  {
89  int ntracksize = verticesItr->tracksSize();
90  double vz = verticesItr->z();
91  if(fabs(vz) > max_Vz_) continue;
92  if( ntracksize > nmax)
93  {
94  vzmax = vz;
95  nmax = ntracksize;
96  }
97  }
98 
100  iEvent.getByToken(pixelTracksToken_,trackCollection);
101  if(trackCollection.isValid())
102  {
103  const reco::RecoChargedCandidateCollection * tracks = trackCollection.product();
104  reco::RecoChargedCandidateCollection::const_iterator tracksItr;
105  int icount=-1;
106  for (tracksItr=tracks->begin(); tracksItr!=tracks->end(); ++tracksItr)
107  {
108  icount++;
109  double eta = tracksItr->eta();
110  if(fabs(eta) > max_Eta_) continue;
111  double pt = tracksItr->pt();
112  if(pt < min_Pt_ || pt > max_Pt_) continue;
113  double vz = tracksItr->vz();
114  if(fabs(vz-vzmax) > min_sep_) continue;
115 
117  filterproduct.addObject(trigger::TriggerTrack, candref);
118  nTrackCandidate++;
119  }
120  }
121  }
122  }
123 
124  accept = ( nTrackCandidate >= min_trks_ );
125 
126  return accept;
127 }
128 
HLTSingleVertexPixelTrackFilter(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > pixelTracksToken_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::VertexCollection > pixelVerticesToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
T const * product() const
Definition: Handle.h:74
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.