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 
13 
17 
19 
22 
24 
27 
29 
30 //
31 // constructors and destructor
32 //
33 
35  pixlTag_ (iConfig.getParameter<edm::InputTag>("pixlTag")),
36  min_Pt_ (iConfig.getParameter<double>("MinPt")),
37  min_trks_ (iConfig.getParameter<unsigned int>("MinTrks")),
38  min_sep_ (iConfig.getParameter<double>("MinSep"))
39 
40 {
41  LogDebug("") << "MinPt cut " << min_Pt_ << "pixl: " << pixlTag_.encode();
42  LogDebug("") << "Requesting : " << min_trks_ << " tracks from same vertex ";
43  LogDebug("") << "Requesting tracks from same vertex eta-phi separation by " << min_sep_;
44 
45  //register your products
46  produces<trigger::TriggerFilterObjectWithRefs>();
47 }
48 
50 {
51 }
52 
53 //
54 // member functions
55 //
56 
57 // ------------ method called to produce the data ------------
59 {
60  using namespace std;
61  using namespace edm;
62  using namespace reco;
63  using namespace trigger;
64 
65  // All HLT filters must create and fill an HLT filter object,
66  // recording any reconstructed physics objects satisfying (or not)
67  // this HLT filter, and place it in the Event.
68 
69 
70 
71  // Specific filter code
72 
73  // get hold of products from Event
74 
76 
77  iEvent.getByLabel(pixlTag_,tracks);
78 
79  // pixel tracks
80  int npixl_tot = 0;
81  vector<double> etastore;
82  vector<double> phistore;
83  vector<int> itstore;
84  bool accept;
85  RecoChargedCandidateCollection::const_iterator apixl(tracks->begin());
86  RecoChargedCandidateCollection::const_iterator epixl(tracks->end());
87  RecoChargedCandidateCollection::const_iterator ipixl, jpixl;
88  unsigned int nsame_vtx=0;
89  int itrk = -1;
90  if (tracks->size() >= min_trks_) {
91  for (ipixl=apixl; ipixl!=epixl; ipixl++){
92  itrk++;
93  const double& ztrk1 = ipixl->vz();
94  const double& etatrk1 = ipixl->momentum().eta();
95  const double& phitrk1 = ipixl->momentum().phi();
96  nsame_vtx=1;
97  etastore.clear();
98  phistore.clear();
99  itstore.clear();
100  etastore.push_back(etatrk1);
101  phistore.push_back(phitrk1);
102  itstore.push_back(itrk);
103  if (fabs(ztrk1) < 15.0) {
104  // check this track against all others to see if others start from same point
105  int jtrk=-1;
106  for (jpixl=apixl; jpixl!=epixl; jpixl++) {
107  jtrk++;
108  if (jpixl==ipixl) continue;
109  const double& ztrk2 = jpixl->vz();
110  const double& etatrk2 = jpixl->momentum().eta();
111  const double& phitrk2 = jpixl->momentum().phi();
112  double eta_dist=etatrk2-etatrk1;
113  double phi_dist=phitrk2-phitrk1;
114  double etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
115  if (fabs(ztrk2-ztrk1) < 1.0 && etaphi_dist > min_sep_) {
116  if (min_trks_ <= 2 || itstore.size() <= 1) {
117  etastore.push_back(etatrk2);
118  phistore.push_back(phitrk2);
119  itstore.push_back(jtrk);
120  nsame_vtx++;
121  } else {
122  // check also separation to already found 'second' tracks
123  LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
124  bool isok = true;
125  for (unsigned int k=1; k < itstore.size(); k++) {
126  eta_dist=etatrk2-etastore.at(k);
127  phi_dist=phitrk2-phistore.at(k);
128  etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
129  if (etaphi_dist < min_sep_) {
130  isok=false;
131  break;
132  }
133  }
134  if (isok) {
135  etastore.push_back(etatrk2);
136  phistore.push_back(phitrk2);
137  itstore.push_back(jtrk);
138  nsame_vtx++;
139  }
140  }
141  }
142  if (nsame_vtx >= min_trks_) break;
143  }
144  }
145  npixl_tot++;
146 
147  if (nsame_vtx >= min_trks_) break;
148  }
149 
150  // final filter decision:
151  // request at least min_trks_ tracks compatible with vertex-region
152  accept = (nsame_vtx >= min_trks_ ) ;
153 
154  } else {
155  accept = false;
156  }
157 
158  // At this point we have the indices of the accepted tracks stored in itstore
159  // we now move them to the filterobject
160 
161  // The filter object
162  auto_ptr<TriggerFilterObjectWithRefs> filterobject (new TriggerFilterObjectWithRefs(path(),module()));
163 
164  if (accept) {
165  for (unsigned int ipos=0; ipos < itstore.size(); ipos++) {
166  int iaddr=itstore.at(ipos);
167  filterobject->addObject(TriggerTrack,RecoChargedCandidateRef(tracks,iaddr));
168  }
169  }
170  // put filter object into the Event
171  iEvent.put(filterobject);
172 
173 
174  LogDebug("") << "Number of pixel-track objects accepted:"
175  << " " << npixl_tot;
176 
177  // return with final filter decision
178  return accept;
179 
180 }
#define LogDebug(id)
int module() const
Definition: HLTadd.h:12
HLTPixlMBFilt(const edm::ParameterSet &)
edm::InputTag pixlTag_
Definition: HLTPixlMBFilt.h:30
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:21
std::string encode() const
Definition: InputTag.cc:72
int path() const
Definition: HLTadd.h:3
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
T sqrt(T t)
Definition: SSEVec.h:28
unsigned int min_trks_
Definition: HLTPixlMBFilt.h:33
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
int k[5][pyjets_maxn]
tuple tracks
Definition: testEve_cfg.py:39
virtual bool filter(edm::Event &, const edm::EventSetup &)