CMS 3D CMS Logo

HLTPixlMBForAlignmentFilter.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  min_isol_ (iConfig.getParameter<double>("MinIsol"))
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  LogDebug("") << "Requesting track to be isolated within cone of " << min_isol_;
44 }
45 
47 
48 void
52  desc.add<edm::InputTag>("pixlTag",edm::InputTag("hltPixelCands"));
53  desc.add<double>("MinPt",5.0);
54  desc.add<unsigned int>("MinTrks",2);
55  desc.add<double>("MinSep",1.0);
56  desc.add<double>("MinIsol",0.05);
57  descriptions.add("hltPixlMBForAlignmentFilter",desc);
58 }
59 
60 //
61 // member functions
62 //
63 
64 // ------------ method called to produce the data ------------
66 {
67  using namespace std;
68  using namespace edm;
69  using namespace reco;
70  using namespace trigger;
71 
72  // All HLT filters must create and fill an HLT filter object,
73  // recording any reconstructed physics objects satisfying (or not)
74  // this HLT filter, and place it in the Event.
75 
76 
77 
78  // Specific filter code
79 
80  // get hold of products from Event
81 
83  iEvent.getByToken(pixlToken_,tracks);
84 
85  // pixel tracks
86  vector<double> etastore;
87  vector<double> phistore;
88  vector<double> ptstore;
89  vector<int> itstore;
90  bool accept = false;
91  auto apixl(tracks->begin());
92  auto epixl(tracks->end());
93  int itrk = 0;
94  double zvtxfit = 0.0;
95  double zvtxfit2 = 0.0;
96  if (tracks->size() >= min_trks_) {
97  etastore.clear();
98  phistore.clear();
99  itstore.clear();
100  for (auto ipixl=apixl; ipixl!=epixl; ipixl++){
101  const double& ztrk1 = ipixl->vz();
102  const double& etatrk1 = ipixl->momentum().eta();
103  const double& phitrk1 = ipixl->momentum().phi();
104  const double& pttrk1 = ipixl->pt();
105  zvtxfit = zvtxfit + ztrk1;
106  zvtxfit2 = zvtxfit2 + ztrk1 * ztrk1;
107  if (pttrk1 > min_Pt_) {
108 // the *store-vectors store the tracks above pt-cut
109 // itstore is the position in the original collection
110  etastore.push_back(etatrk1);
111  phistore.push_back(phitrk1);
112  ptstore.push_back(pttrk1);
113  itstore.push_back(itrk);
114  }
115  itrk++;
116  }
117  if (itrk > 0) {
118 // implement proper vertex fit here ?
119  zvtxfit = zvtxfit / itrk;
120  zvtxfit2 = zvtxfit2 / itrk;
121  zvtxfit2 = sqrt(zvtxfit2 - zvtxfit*zvtxfit);
122  }
123 // locisol is the position in the *store vectors
124  vector<int> locisol;
125  if (itstore.size() > 1) {
126  // now check that tracks are isolated
127  locisol.clear();
128  for (unsigned int i=0; i<itstore.size(); i++) {
129  int nincone=0;
130 // check isolation wrt ALL tracks, not only those above ptcut
131  for (auto ipixl=apixl; ipixl!=epixl; ipixl++){
132  double phidist=std::abs( phistore.at(i) - ipixl->momentum().phi() );
133  double etadist=std::abs( etastore.at(i) - ipixl->momentum().eta() );
134  double trkdist = sqrt(phidist*phidist + etadist*etadist);
135  if (trkdist < min_isol_) nincone++;
136  }
137 // the check above always find the track itself, so nincone never should be 0
138  if (nincone < 2) locisol.push_back(i);
139  }
140  }
141 // now check that the selected tracks have enough mutual separation
142  vector<int> itsep;
143  for (unsigned int i=0; i<locisol.size(); i++) {
144 // check for each so far selected track...
145  itsep.clear();
146  itsep.push_back(locisol.at(i));
147  for (unsigned int j=i+1; j<locisol.size(); j++) {
148 // ...if it is sufficiently separated from other selectad tracks...
149  double phidist = phistore.at(locisol.at(i))-phistore.at(locisol.at(j));
150  double etadist = etastore.at(locisol.at(i))-etastore.at(locisol.at(j));
151  double dist = sqrt(phidist*phidist + etadist*etadist);
152  if (dist > min_sep_) {
153  if (itsep.size() == 1) {
154  itsep.push_back(locisol.at(j));
155  } else {
156  bool is_separated = true;
157  for (int k : itsep){
158 // ...and the other ones, that are on the 'final acceptance' list already, if min_trks_ > 2
159  double phisep = phistore.at(k)-phistore.at(locisol.at(j));
160  double etasep = etastore.at(k)-etastore.at(locisol.at(j));
161  double sep = sqrt(phisep*phisep + etasep*etasep);
162  if (sep < min_sep_) {
163 // this one was no good, too close to some other already accepted
164  is_separated = false;
165  break;
166  }
167  }
168  if (is_separated) itsep.push_back(locisol.at(j));
169  }
170  }
171  if (itsep.size() >= min_trks_) {
172  accept = true;
173  break;
174  }
175  }
176  if (accept) {
177  break;
178  }
179  }
180  // At this point we have the indices of the accepted tracks stored in itstore
181  // we now move them to the filterproduct
182 
183  if (accept) {
184  for (int ipos : itsep) {
185  int iaddr=itstore.at(ipos);
186  filterproduct.addObject(TriggerTrack,RecoChargedCandidateRef(tracks,iaddr));
187  }
188  // std::cout << "Accept this event " << std::endl;
189  }
190  }
191 
192 
193 // LogDebug("") << "Number of pixel-track objects accepted:"
194 // << " " << npixl_tot;
195 
196  // return with final filter decision
197  return accept;
198 
199 }
200 
201 // declare this class as a framework plugin
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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:159
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HLTPixlMBForAlignmentFilter(const edm::ParameterSet &)
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > pixlToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
int k[5][pyjets_maxn]
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
~HLTPixlMBForAlignmentFilter() override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override