CMS 3D CMS Logo

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