CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 }
47 
48 void
52  desc.add<edm::InputTag>("pixlTag",edm::InputTag("hltPixelCands"));
53  desc.add<double>("MinPt",0.);
54  desc.add<unsigned int>("MinTrks",2);
55  desc.add<double>("MinSep",1.);
56  descriptions.add("hltPixlMBFilt",desc);
57 }
58 
59 //
60 // member functions
61 //
62 
63 // ------------ method called to produce the data ------------
65 {
66  using namespace std;
67  using namespace edm;
68  using namespace reco;
69  using namespace trigger;
70 
71  // All HLT filters must create and fill an HLT filter object,
72  // recording any reconstructed physics objects satisfying (or not)
73  // this HLT filter, and place it in the Event.
74 
75 
76 
77  // Specific filter code
78 
79  // get hold of products from Event
80 
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  RecoChargedCandidateCollection::const_iterator apixl(tracks->begin());
91  RecoChargedCandidateCollection::const_iterator 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  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  jtrk++;
113  if (jpixl==ipixl) continue;
114  const double& ztrk2 = jpixl->vz();
115  const double& etatrk2 = jpixl->momentum().eta();
116  const double& phitrk2 = jpixl->momentum().phi();
117  double eta_dist=etatrk2-etatrk1;
118  double phi_dist=phitrk2-phitrk1;
119  double etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
120  if (fabs(ztrk2-ztrk1) < 1.0 && etaphi_dist > min_sep_) {
121  if (min_trks_ <= 2 || itstore.size() <= 1) {
122  etastore.push_back(etatrk2);
123  phistore.push_back(phitrk2);
124  itstore.push_back(jtrk);
125  nsame_vtx++;
126  } else {
127  // check also separation to already found 'second' tracks
128  LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
129  bool isok = true;
130  for (unsigned int k=1; k < itstore.size(); k++) {
131  eta_dist=etatrk2-etastore.at(k);
132  phi_dist=phitrk2-phistore.at(k);
133  etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
134  if (etaphi_dist < min_sep_) {
135  isok=false;
136  break;
137  }
138  }
139  if (isok) {
140  etastore.push_back(etatrk2);
141  phistore.push_back(phitrk2);
142  itstore.push_back(jtrk);
143  nsame_vtx++;
144  }
145  }
146  }
147  if (nsame_vtx >= min_trks_) break;
148  }
149  }
150  npixl_tot++;
151 
152  if (nsame_vtx >= min_trks_) break;
153  }
154 
155  // final filter decision:
156  // request at least min_trks_ tracks compatible with vertex-region
157  accept = (nsame_vtx >= min_trks_ ) ;
158 
159  } else {
160  accept = false;
161  }
162 
163  // At this point we have the indices of the accepted tracks stored in itstore
164  // we now move them to the filterproduct
165 
166  if (accept) {
167  for (unsigned int ipos=0; ipos < itstore.size(); ipos++) {
168  int iaddr=itstore.at(ipos);
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:434
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:26
std::string encode() const
Definition: InputTag.cc:164
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:243
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
T sqrt(T t)
Definition: SSEVec.h:48
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
tuple tracks
Definition: testEve_cfg.py:39
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)