CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTPixlMBFilt.cc
Go to the documentation of this file.
1 
10 #include "HLTPixlMBFilt.h"
11 
17 
19 
21 
24 
26 
27 //
28 // constructors and destructor
29 //
30 
32  : HLTFilter(iConfig),
33  pixlTag_(iConfig.getParameter<edm::InputTag>("pixlTag")),
34  min_Pt_(iConfig.getParameter<double>("MinPt")),
35  min_trks_(iConfig.getParameter<unsigned int>("MinTrks")),
36  min_sep_(iConfig.getParameter<double>("MinSep"))
37 
38 {
39  pixlToken_ = consumes<reco::RecoChargedCandidateCollection>(pixlTag_);
40  LogDebug("") << "MinPt cut " << min_Pt_ << "pixl: " << pixlTag_.encode();
41  LogDebug("") << "Requesting : " << min_trks_ << " tracks from same vertex ";
42  LogDebug("") << "Requesting tracks from same vertex eta-phi separation by " << min_sep_;
43 }
44 
46 
50  desc.add<edm::InputTag>("pixlTag", edm::InputTag("hltPixelCands"));
51  desc.add<double>("MinPt", 0.);
52  desc.add<unsigned int>("MinTrks", 2);
53  desc.add<double>("MinSep", 1.);
54  descriptions.add("hltPixlMBFilt", desc);
55 }
56 
57 //
58 // member functions
59 //
60 
61 // ------------ method called to produce the data ------------
63  const edm::EventSetup& iSetup,
64  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
65  using namespace std;
66  using namespace edm;
67  using namespace reco;
68  using namespace trigger;
69 
70  // All HLT filters must create and fill an HLT filter object,
71  // recording any reconstructed physics objects satisfying (or not)
72  // this HLT filter, and place it in the Event.
73  if (saveTags()) {
74  filterproduct.addCollectionTag(pixlTag_);
75  }
76 
77  // Specific filter code
78 
79  // get hold of products from Event
80 
82  iEvent.getByToken(pixlToken_, tracks);
83 
84  // pixel tracks
85  int npixl_tot = 0;
86  vector<double> etastore;
87  vector<double> phistore;
88  vector<int> itstore;
89  bool accept;
90  auto apixl(tracks->begin());
91  auto epixl(tracks->end());
92  RecoChargedCandidateCollection::const_iterator ipixl, jpixl;
93  unsigned int nsame_vtx = 0;
94  int itrk = -1;
95  if (tracks->size() >= min_trks_) {
96  for (ipixl = apixl; ipixl != epixl; ipixl++) {
97  if (ipixl->pt() < min_Pt_)
98  continue;
99  itrk++;
100  const double& ztrk1 = ipixl->vz();
101  const double& etatrk1 = ipixl->momentum().eta();
102  const double& phitrk1 = ipixl->momentum().phi();
103  nsame_vtx = 1;
104  etastore.clear();
105  phistore.clear();
106  itstore.clear();
107  etastore.push_back(etatrk1);
108  phistore.push_back(phitrk1);
109  itstore.push_back(itrk);
110  if (fabs(ztrk1) < 15.0) {
111  // check this track against all others to see if others start from same point
112  int jtrk = -1;
113  for (jpixl = apixl; jpixl != epixl; jpixl++) {
114  if (jpixl->pt() < min_Pt_)
115  continue;
116  jtrk++;
117  if (jpixl == ipixl)
118  continue;
119  const double& ztrk2 = jpixl->vz();
120  const double& etatrk2 = jpixl->momentum().eta();
121  const double& phitrk2 = jpixl->momentum().phi();
122  double eta_dist = etatrk2 - etatrk1;
123  double phi_dist = phitrk2 - phitrk1;
124  double etaphi_dist = sqrt(eta_dist * eta_dist + phi_dist * phi_dist);
125  if (fabs(ztrk2 - ztrk1) < 1.0 && etaphi_dist > min_sep_) {
126  if (min_trks_ <= 2 || itstore.size() <= 1) {
127  etastore.push_back(etatrk2);
128  phistore.push_back(phitrk2);
129  itstore.push_back(jtrk);
130  nsame_vtx++;
131  } else {
132  // check also separation to already found 'second' tracks
133  LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
134  bool isok = true;
135  for (unsigned int k = 1; k < itstore.size(); k++) {
136  eta_dist = etatrk2 - etastore.at(k);
137  phi_dist = phitrk2 - phistore.at(k);
138  etaphi_dist = sqrt(eta_dist * eta_dist + phi_dist * phi_dist);
139  if (etaphi_dist < min_sep_) {
140  isok = false;
141  break;
142  }
143  }
144  if (isok) {
145  etastore.push_back(etatrk2);
146  phistore.push_back(phitrk2);
147  itstore.push_back(jtrk);
148  nsame_vtx++;
149  }
150  }
151  }
152  if (nsame_vtx >= min_trks_)
153  break;
154  }
155  }
156  npixl_tot++;
157 
158  if (nsame_vtx >= min_trks_)
159  break;
160  }
161 
162  // final filter decision:
163  // request at least min_trks_ tracks compatible with vertex-region
164  accept = (nsame_vtx >= min_trks_);
165 
166  } else {
167  accept = false;
168  }
169 
170  // At this point we have the indices of the accepted tracks stored in itstore
171  // we now move them to the filterproduct
172 
173  if (accept) {
174  for (int iaddr : itstore) {
175  filterproduct.addObject(TriggerTrack, RecoChargedCandidateRef(tracks, iaddr));
176  }
177  }
178 
179  LogDebug("") << "Number of pixel-track objects accepted:"
180  << " " << npixl_tot;
181 
182  // return with final filter decision
183  return accept;
184 }
185 
186 // declare this class as a framework plugin
HLTPixlMBFilt(const edm::ParameterSet &)
edm::InputTag pixlTag_
Definition: HLTPixlMBFilt.h:38
~HLTPixlMBFilt() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Ref< RecoChargedCandidateCollection > RecoChargedCandidateRef
reference to an object in a collection of RecoChargedCandidate objects
auto const & tracks
cannot be loose
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
std::string encode() const
Definition: InputTag.cc:159
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
int iEvent
Definition: GenABIO.cc:224
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > pixlToken_
Definition: HLTPixlMBFilt.h:39
unsigned int min_trks_
Definition: HLTPixlMBFilt.h:42
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:25
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:46
#define LogDebug(id)