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 
47 {
48 }
49 
50 //
51 // member functions
52 //
53 
54 // ------------ method called to produce the data ------------
56 {
57  using namespace std;
58  using namespace edm;
59  using namespace reco;
60  using namespace trigger;
61 
62  // All HLT filters must create and fill an HLT filter object,
63  // recording any reconstructed physics objects satisfying (or not)
64  // this HLT filter, and place it in the Event.
65 
66 
67 
68  // Specific filter code
69 
70  // get hold of products from Event
71 
73 
74  iEvent.getByLabel(pixlTag_,tracks);
75 
76  // pixel tracks
77  int npixl_tot = 0;
78  vector<double> etastore;
79  vector<double> phistore;
80  vector<int> itstore;
81  bool accept;
82  RecoChargedCandidateCollection::const_iterator apixl(tracks->begin());
83  RecoChargedCandidateCollection::const_iterator epixl(tracks->end());
84  RecoChargedCandidateCollection::const_iterator ipixl, jpixl;
85  unsigned int nsame_vtx=0;
86  int itrk = -1;
87  if (tracks->size() >= min_trks_) {
88  for (ipixl=apixl; ipixl!=epixl; ipixl++){
89  itrk++;
90  const double& ztrk1 = ipixl->vz();
91  const double& etatrk1 = ipixl->momentum().eta();
92  const double& phitrk1 = ipixl->momentum().phi();
93  nsame_vtx=1;
94  etastore.clear();
95  phistore.clear();
96  itstore.clear();
97  etastore.push_back(etatrk1);
98  phistore.push_back(phitrk1);
99  itstore.push_back(itrk);
100  if (fabs(ztrk1) < 15.0) {
101  // check this track against all others to see if others start from same point
102  int jtrk=-1;
103  for (jpixl=apixl; jpixl!=epixl; jpixl++) {
104  jtrk++;
105  if (jpixl==ipixl) continue;
106  const double& ztrk2 = jpixl->vz();
107  const double& etatrk2 = jpixl->momentum().eta();
108  const double& phitrk2 = jpixl->momentum().phi();
109  double eta_dist=etatrk2-etatrk1;
110  double phi_dist=phitrk2-phitrk1;
111  double etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
112  if (fabs(ztrk2-ztrk1) < 1.0 && etaphi_dist > min_sep_) {
113  if (min_trks_ <= 2 || itstore.size() <= 1) {
114  etastore.push_back(etatrk2);
115  phistore.push_back(phitrk2);
116  itstore.push_back(jtrk);
117  nsame_vtx++;
118  } else {
119  // check also separation to already found 'second' tracks
120  LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
121  bool isok = true;
122  for (unsigned int k=1; k < itstore.size(); k++) {
123  eta_dist=etatrk2-etastore.at(k);
124  phi_dist=phitrk2-phistore.at(k);
125  etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
126  if (etaphi_dist < min_sep_) {
127  isok=false;
128  break;
129  }
130  }
131  if (isok) {
132  etastore.push_back(etatrk2);
133  phistore.push_back(phitrk2);
134  itstore.push_back(jtrk);
135  nsame_vtx++;
136  }
137  }
138  }
139  if (nsame_vtx >= min_trks_) break;
140  }
141  }
142  npixl_tot++;
143 
144  if (nsame_vtx >= min_trks_) break;
145  }
146 
147  // final filter decision:
148  // request at least min_trks_ tracks compatible with vertex-region
149  accept = (nsame_vtx >= min_trks_ ) ;
150 
151  } else {
152  accept = false;
153  }
154 
155  // At this point we have the indices of the accepted tracks stored in itstore
156  // we now move them to the filterproduct
157 
158  if (accept) {
159  for (unsigned int ipos=0; ipos < itstore.size(); ipos++) {
160  int iaddr=itstore.at(ipos);
162  }
163  }
164 
165  LogDebug("") << "Number of pixel-track objects accepted:"
166  << " " << npixl_tot;
167 
168  // return with final filter decision
169  return accept;
170 
171 }
#define LogDebug(id)
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:22
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)
T sqrt(T t)
Definition: SSEVec.h:48
unsigned int min_trks_
Definition: HLTPixlMBFilt.h:33
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
int k[5][pyjets_maxn]
tuple tracks
Definition: testEve_cfg.py:39