CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackerTrackHitFilter.cc
Go to the documentation of this file.
5 
9 
15 
22 
30 
31 //for S/N cut
38 //for angle cut
44 #include "TMath.h"
45 
49 
50 #include <boost/regex.hpp>
51 #include <map>
52 //#include <math>
53 
75 namespace reco { namespace modules {
77  public:
78  TrackerTrackHitFilter(const edm::ParameterSet &iConfig) ;
79  virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
80  int checkHit(const edm::EventSetup &iSetup,const DetId &detid, const TrackingRecHit * hit);
81  void produceFromTrajectory( const edm::EventSetup &iSetup, const Trajectory *itt, std::vector<TrackingRecHit *>&hits);
82  void produceFromTrack( const edm::EventSetup &iSetup, const Track *itt, std::vector<TrackingRecHit *>&hits);
83  private:
84  class Rule {
85  public:
86  // parse a rule from a string
87  Rule(const std::string &str) ;
88  // check this DetId, update the verdict if the rule has anything to say about it
89  void apply(DetId detid, bool & verdict) const {
90  // check detector
91  if (detid.subdetId() == subdet_) {
92  // check layer
93  if ( (layer_ == 0) || (layer_ == layer(detid)) ) {
94  // override verdict
95  verdict = keep_;
96  // std::cout<<"Verdict is "<<verdict<<" for subdet "<<subdet_<<", layer "<<layer_<<"! "<<std::endl;
97  }
98  // else std::cout<<"No, sorry, wrong layer.Retry..."<<std::endl;
99  }
100  // else{ std::cout<<"No, sorry, wrong subdet.Retry..."<<std::endl;}
101  }
102  private:
103  int subdet_;
104 
105  int layer_;
106  bool keep_;
107  int layer(DetId detid) const ;
108 
109  };
110 
112 
113  int iRun;
114  int iEvt;
115 
116  size_t minimumHits_;
117 
122 
124  std::string CMNSubtractionMode_;
125  std::vector<bool> subdetStoN_;//(6); //,std::bool(false));
126  std::vector<double> subdetStoNlowcut_;//(6,-1.0);
127  std::vector<double> subdetStoNhighcut_;//(6,-1.0);
128  bool checkStoN( const edm::EventSetup &iSetup,const DetId &id, const TrackingRecHit* therechit);
129  void parseStoN(const std::string &str);
130 
131  std::vector<uint32_t> detsToIgnore_;
132  std::vector<Rule> rules_;
133 
137  bool checkHitAngle(const TrajectoryMeasurement &meas);
141  std::vector<int32_t> pxlTPLqBin_;
142 
145 
148 
150 
151 
152 
155  int layerFromId (const DetId& id) const;
156  // bool checkOverlapHit();
157 
158  TrackCandidate makeCandidate(const reco::Track &tk, std::vector<TrackingRecHit *>::iterator hitsBegin, std::vector<TrackingRecHit *>::iterator hitsEnd) ;
159  //const TransientTrackingRecHitBuilder *RHBuilder;
160  }; // class
161 
162 
163 
164 TrackerTrackHitFilter::Rule::Rule(const std::string &str) {
165  static boost::regex rule("(keep|drop)\\s+([A-Z]+)(\\s+(\\d+))?");
166  boost::cmatch match;
167  std::string match_1;
168  std::string match_2;
169  std::string match_3;
170  // match and check it works
171  if (!regex_match(str.c_str(), match, rule)) {
172  throw cms::Exception("Configuration") << "Rule '" << str << "' not understood.\n";
173  }
174  else{
175  std::cout<<"*** Rule Command given to TrackerTrackHitFilter:\t"<<str<<std::endl;
176 
177  }
178  // Set up fields:
179  // rule type
180  keep_ = (strncmp(match[1].first,"keep",4 )== 0);
181 
182  // subdet
183  subdet_ = -1;
184  if (strncmp(match[2].first, "PXB", 3) == 0) subdet_ = PixelSubdetector::PixelBarrel;
185  else if (strncmp(match[2].first, "PXE", 3) == 0) subdet_ = PixelSubdetector::PixelEndcap;
186  else if (strncmp(match[2].first, "TIB", 3) == 0) subdet_ = StripSubdetector::TIB;
187  else if (strncmp(match[2].first, "TID", 3) == 0) subdet_ = StripSubdetector::TID;
188  else if (strncmp(match[2].first, "TOB", 3) == 0) subdet_ = StripSubdetector::TOB;
189  else if (strncmp(match[2].first, "TEC", 3) == 0) subdet_ = StripSubdetector::TEC;
190  if (subdet_ == -1) {
191  throw cms::Exception("Configuration") << "Detector '" << match[2].first << "' not understood. Should be PXB, PXE, TIB, TID, TOB, TEC.\n";
192  }
193  // layer (if present)
194  if (match[4].first != match[4].second) {
195  layer_ = atoi(match[4].first);
196  } else {
197  layer_ = 0;
198  }
199 }//end Rule::Rule
200 
202  switch (detid.subdetId()) {
203  case PixelSubdetector::PixelBarrel: return PXBDetId(detid).layer();
204  case PixelSubdetector::PixelEndcap: return PXFDetId(detid).disk();
205  case StripSubdetector::TIB: return TIBDetId(detid).layer();
206  case StripSubdetector::TID: return TIDDetId(detid).wheel();
207  case StripSubdetector::TOB: return TOBDetId(detid).layer();
208  case StripSubdetector::TEC: return TECDetId(detid).wheel();
209  }
210  return -1; // never match
211 }
212 
213 void TrackerTrackHitFilter::parseStoN(const std::string &str) {
214  // match a set of capital case chars (preceded by an arbitrary number of leading blanks),
215  //followed b an arbitrary number of blanks, one or more digits (not necessary, they cannot also be,
216  // another set of blank spaces and, again another *eventual* digit
217  // static boost::regex rule("\\s+([A-Z]+)(\\s+(\\d+)(\\.)?(\\d+))?(\\s+(\\d+)(\\.)?(\\d+))?");
218  static boost::regex rule("([A-Z]+)"
219  "\\s*(\\d+\\.*\\d*)?"
220  "\\s*(\\d+\\.*\\d*)?");
221 
222 
223  boost::cmatch match;
224  std::string match_1;
225  std::string match_2;
226  std::string match_3;
227  // match and check it works
228  if (!regex_match(str.c_str(), match, rule)) {
229  throw cms::Exception("Configuration") << "Rule for S to N cut '" << str << "' not understood.\n";
230  }
231  else{
232  std::string match_0=(match[0].first,match[0].second);
233  match_1=(match[1].first,match[1].second);
234  match_2=(match[2].first,match[2].second);
235  match_3=(match[3].first,match[3].second);
236 
237  }
238 
239  int cnt=0;
240  float subdet_ind[6];
241  for(cnt=0;cnt<6;cnt++ ){
242  subdet_ind[cnt]=-1.0;
243  }
244 
245 
246  bool doALL=false;
247  std::string match_1a(match[1].first,match[1].second);
248  if (strncmp(match[1].first, "ALL", 3) == 0) doALL=true;
249  if (doALL || strncmp(match[1].first, "PXB", 3) == 0) subdet_ind[0] = +1.0;
250  if (doALL || strncmp(match[1].first, "PXE", 3) == 0) subdet_ind[1] = +1.0;
251  if (doALL || strncmp(match[1].first, "TIB", 3) == 0) subdet_ind[2] = +1.0;
252  if (doALL || strncmp(match[1].first, "TID", 3) == 0) subdet_ind[3] = +1.0;
253  if (doALL || strncmp(match[1].first, "TOB", 3) == 0) subdet_ind[4] = +1.0;
254  if (doALL || strncmp(match[1].first, "TEC", 3) == 0) subdet_ind[5] = +1.0;
255 
256  for(cnt=0;cnt<6;cnt++ ){//loop on subdets
257  if(subdet_ind[cnt]>0.0){
258  subdetStoN_[cnt]=true;
259  if (match[2].first != match[2].second) {
260  subdetStoNlowcut_[cnt] = atof(match[2].first);
261  }
262  if (match[3].first != match[3].second ) {
263  subdetStoNhighcut_[cnt] = atof(match[3].first);
264  }
265  std::cout<<"Setting thresholds*&^ for subdet #"<<cnt+1<<" = "<<subdetStoNlowcut_[cnt]<<" - "<<subdetStoNhighcut_[cnt]<<std::endl;
266  }
267  }
268 
269  bool correct_regex=false;
270  for(cnt=0;cnt<6;cnt++ ){//check that the regex was correct
271  if(subdetStoN_[cnt])correct_regex=true;
272  }
273 
274  if (!correct_regex) {
275  throw cms::Exception("Configuration") << "Detector '" << match_1a << "' not understood in parseStoN. Should be PXB, PXE, TIB, TID, TOB, TEC.\n";
276  }
277 
278  //std::cout<<"Reached end of parseStoN"<<std::endl;
279 }//end parseStoN
280 
281 
283  src_(iConfig.getParameter<edm::InputTag>("src")),
284  minimumHits_(iConfig.getParameter<uint32_t>("minimumHits")),
285  replaceWithInactiveHits_(iConfig.getParameter<bool>("replaceWithInactiveHits")),
286  stripFrontInvalidHits_(iConfig.getParameter<bool>("stripFrontInvalidHits")),
287  stripBackInvalidHits_( iConfig.getParameter<bool>("stripBackInvalidHits") ),
288  stripAllInvalidHits_( iConfig.getParameter<bool>("stripAllInvalidHits") ),
289  rejectBadStoNHits_( iConfig.getParameter<bool>("rejectBadStoNHits") ),
290  CMNSubtractionMode_( iConfig.getParameter<std::string>("CMNSubtractionMode") ),
291  detsToIgnore_( iConfig.getParameter<std::vector<uint32_t> >("detsToIgnore") ),
292  useTrajectories_( iConfig.getParameter<bool>("useTrajectories") ),
293  rejectLowAngleHits_( iConfig.getParameter<bool>("rejectLowAngleHits") ),
294  TrackAngleCut_( iConfig.getParameter<double>("TrackAngleCut") ),
295  checkPXLQuality_(iConfig.getParameter<bool>("usePixelQualityFlag") ),
296  pxlTPLProbXY_(iConfig.getParameter<double>("PxlTemplateProbXYCut")),
297  pxlTPLProbXYQ_(iConfig.getParameter<double>("PxlTemplateProbXYChargeCut")),
298  pxlTPLqBin_(iConfig.getParameter<std::vector<int32_t> >("PxlTemplateqBinCut")),
299  PXLcorrClusChargeCut_(iConfig.getParameter<double>("PxlCorrClusterChargeCut")),
300  tagOverlaps_( iConfig.getParameter<bool>("tagOverlaps") )
301 {
302 
303 
304  // sanity check
306  throw cms::Exception("Configuration") << "Inconsistent Configuration: you can't set both 'stripAllInvalidHits' and 'replaceWithInactiveHits' to true\n";
307  }
309  throw cms::Exception("Configuration") << "Wrong configuration of TrackerTrackHitFilter. You cannot apply the cut on the track angle without using Trajectories!\n";
310  }
312  throw cms::Exception("Configuration") << "Wrong configuration of TrackerTrackHitFilter. You cannot apply the cut on the corrected pixel cluster charge without using Trajectories!\n";
313  }
314 
315  if(pxlTPLqBin_.size()>2){
316  std::cout<<"Warning from TrackerTrackHitFilter: vector with qBin cuts has size > 2. Additional items will be ignored."<<std::endl;
317  }
318 
319 
320 
321  // read and parse commands
322  std::vector<std::string> str_rules = iConfig.getParameter<std::vector<std::string> >("commands");
323  rules_.reserve(str_rules.size());
324  for (std::vector<std::string>::const_iterator it = str_rules.begin(), ed = str_rules.end(); it != ed; ++it) {
325  rules_.push_back(Rule(*it));
326  }
327 
328  if(rejectBadStoNHits_){//commands for S/N cut
329 
330  subdetStoN_.reserve(6);
331  subdetStoNlowcut_.reserve(6);
332  subdetStoNhighcut_.reserve(6);
333  int cnt=0;
334  for(cnt=0;cnt<6;cnt++ ){
335  subdetStoN_[cnt]=false;
336  subdetStoNlowcut_[cnt]=-1.0;
337  subdetStoNhighcut_[cnt]=-1.0;
338  }
339 
340  std::vector<std::string> str_StoNrules = iConfig.getParameter<std::vector<std::string> >("StoNcommands");
341  for (std::vector<std::string>::const_iterator str_StoN = str_StoNrules.begin(); str_StoN != str_StoNrules.end(); ++str_StoN) {
342  parseStoN(*str_StoN);
343  }
345  std::cout<<"Finished parsing S/N. Applying following cuts to subdets:";
346  for(cnt=0;cnt<6;cnt++ ){
348  std::cout<<"Subdet #"<<cnt+1<<" -> "<<subdetStoNlowcut_[cnt]<<" , "<<subdetStoNhighcut_[cnt];
349  }
350  }//end if rejectBadStoNHits_
351 
352 
353  if(rejectLowAngleHits_ ) std::cout<<"\nApplying cut on angle track = "<<TrackAngleCut_<<std::endl;
354 
355 
356  // sort detids to ignore
357  std::sort(detsToIgnore_.begin(), detsToIgnore_.end());
358 
359  // issue the produce<>
360  produces<TrackCandidateCollection>();
361 }
362 
363 void
365 {
366  //Dump Run and Event
367  iRun=iEvent.id().run();
368  iEvt=iEvent.id().event();
369 
370  // read with View, so we can read also a TrackRefVector
373 
374  if(useTrajectories_)iEvent.getByLabel(src_, assoMap);
375  else iEvent.getByLabel(src_, tracks);
376 
377  // read from EventSetup
380  //iSetup.get<TransientRecHitRecord>().get("WithTrackAngle",theBuilder);
381  // RHBuilder= theBuilder.product();
382 
383 
384  // prepare output collection
385  size_t candcollsize;
386  if(useTrajectories_)candcollsize= assoMap->size();
387  else candcollsize=tracks->size();
388  std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection());
389 
390  output->reserve(candcollsize);
391 
392  // working area and tools
393  std::vector<TrackingRecHit *> hits;
394 
395  if(useTrajectories_){
396  for (TrajTrackAssociationCollection::const_iterator itass = assoMap->begin(); itass != assoMap->end(); ++itass){
397 
398  const edm::Ref<std::vector<Trajectory> >traj = itass->key;//trajectory in the collection
399  const reco::TrackRef tkref = itass->val;//associated track track in the collection
400  //std::cout<<"The hit collection has size "<<hits.size()<<" (should be 0) while the track contains initially "<< tkref->recHitsEnd() - tkref->recHitsBegin()<<std::endl;
401 
402  const Track *trk = &(*tkref);
403  const Trajectory * myTrajectory= &(*traj);
404  produceFromTrajectory(iSetup,myTrajectory,hits);
405 
406  std::vector<TrackingRecHit *>::iterator begin = hits.begin(), end = hits.end();
407 
408  // strip invalid hits at the beginning
410  while ( (begin != end) && ( (*begin)->isValid() == false ) ) ++begin;
411  }
412  // strip invalid hits at the end
413  if (stripBackInvalidHits_ && (begin != end)) {
414  --end;
415  while ( (begin != end) && ( (*end)->isValid() == false ) ) --end;
416  ++end;
417  }
418 
419  // if we still have some hits build the track candidate
421  int nvalidhits=0;
422  for( std::vector<TrackingRecHit *>::iterator ithit=begin;ithit!=end;++ithit){
423  if( (*ithit)->isValid())nvalidhits++;
424  }
425  if(nvalidhits >= int(minimumHits_)){
426  output->push_back( makeCandidate ( *trk, begin, end ) );
427  }
428  nvalidhits=0;
429  }
430  else{//all invalid hits have been already kicked out
431  if ((end - begin) >= int(minimumHits_)) {
432  output->push_back( makeCandidate ( *trk, begin, end ) );
433  }
434  }
435 
436 
437  // if we still have some hits build the track candidate
438  //if ((end - begin) >= int(minimumHits_)) {
439  // output->push_back( makeCandidate ( *trk, begin, end ) );
440  //}
441 
442 
443  // now delete the hits not used by the candidate
444  for (begin = hits.begin(), end = hits.end(); begin != end; ++begin) {
445  if (*begin) delete *begin;
446  }
447  hits.clear();
448  } // loop on trajectories
449 
450 
451  }
452  else{ //use plain tracks
453 
454  // loop on tracks
455  for (std::vector<reco::Track>::const_iterator ittrk = tracks->begin(), edtrk = tracks->end(); ittrk != edtrk; ++ittrk) {
456 
457  // std::cout<<"The hit collection has size "<<hits.size()<<" (should be 0) while the track contains initially "<< ittrk->recHitsEnd() - ittrk->recHitsBegin()<<std::endl;
458 
459  const Track *trk = &(*ittrk);
460 
461  produceFromTrack(iSetup,trk,hits);
462  //-----------------------
463  /*
464  std::cout<<"Hit collection in output has size "<<hits.size()<<". Dumping hit positions..."<<std::endl;
465  for (std::vector<TrackingRecHit *>::iterator ith = hits.begin(), edh = hits.end(); ith != edh; ++ith) {
466  const TrackingRecHit *myhit = *(ith);
467  TransientTrackingRecHit::RecHitPointer ttrh;
468  float radius=0.0;
469  float xx=-999.0,yy=-999.0,zz=-999.0;
470  unsigned int myid=0;
471  if(myhit->isValid()){
472  ttrh = RHBuilder->build(myhit);
473  xx=ttrh->globalPosition().x();
474  yy=ttrh->globalPosition().y();
475  zz=ttrh->globalPosition().z();
476  radius = sqrt( pow(xx,2)+pow(yy,2) );
477  myid=myhit->geographicalId().rawId();
478  }
479  std::cout<<"-$-$ OUTPUT Hit position: ( "<<xx<<" , " <<yy<<" , " <<zz<<" ) , RADIUS = " <<radius<<" on DetID= "<< myid<<std::endl;
480  }//end loop on hits
481  */
482  //-----------------------
483 
484 
485  std::vector<TrackingRecHit *>::iterator begin = hits.begin(), end = hits.end();
486  // std::cout << "Back in the main producer (TRK), the final hit collection has size " << hits.size() << std::endl;
487  // strip invalid hits at the beginning
489  while ( (begin != end) && ( (*begin)->isValid() == false ) ) ++begin;
490  }
491  // strip invalid hits at the end
492  if (stripBackInvalidHits_ && (begin != end)) {
493  --end;
494  while ( (begin != end) && ( (*end)->isValid() == false ) ) --end;
495  ++end;
496  }
497 
498  // if we still have some hits build the track candidate
500  int nvalidhits=0;
501  for( std::vector<TrackingRecHit *>::iterator ithit=begin;ithit!=end;++ithit){
502  if( (*ithit)->isValid())nvalidhits++;
503  }
504  if(nvalidhits >= int(minimumHits_)){
505  output->push_back( makeCandidate ( *ittrk, begin, end ) );
506  }
507 
508  }
509  else{//all invalid hits have been already kicked out
510  if ((end - begin) >= int(minimumHits_)) {
511  output->push_back( makeCandidate ( *ittrk, begin, end ) );
512  }
513  }
514 
515  // now delete the hits not used by the candidate
516  for (begin = hits.begin(), end = hits.end(); begin != end; ++begin) {
517  if (*begin) delete *begin;
518  }
519  hits.clear();
520  } // loop on tracks
521  }//end else useTracks
522 
523  // std::cout<<"OUTPUT SIZE: "<<output->size()<<std::endl;
524 
525  iEvent.put(output);
526 }
527 
529 TrackerTrackHitFilter::makeCandidate(const reco::Track &tk, std::vector<TrackingRecHit *>::iterator hitsBegin, std::vector<TrackingRecHit *>::iterator hitsEnd) {
530 
531 
534  if ( pdir == anyDirection ) throw cms::Exception("UnimplementedFeature") << "Cannot work with tracks that have 'anyDirecton' \n";
535 
536  // double innerP=sqrt( pow(tk.innerMomentum().X(),2)+pow(tk.innerMomentum().Y(),2)+pow(tk.innerMomentum().Z(),2) );
537  // if ( (pdir == alongMomentum) == ( innerP >= tk.outerP() ) ) {
538 
539  if ( (pdir == alongMomentum) == ( (tk.outerPosition()-tk.innerPosition()).Dot(tk.momentum()) >= 0 ) ) {
540  // use inner state
542  state = trajectoryStateTransform::persistentState( originalTsosIn, DetId(tk.innerDetId()) );
543  } else {
544  // use outer state
546  state = trajectoryStateTransform::persistentState( originalTsosOut, DetId(tk.outerDetId()) );
547  }
548  TrajectorySeed seed(state, TrackCandidate::RecHitContainer(), pdir);
550  ownHits.reserve(hitsEnd - hitsBegin);
551  for ( ; hitsBegin != hitsEnd; ++hitsBegin) {
552  //if(! (*hitsBegin)->isValid() ) std::cout<<"Putting in the trackcandidate an INVALID HIT !"<<std::endl;
553  ownHits.push_back( *hitsBegin );
554  }
555 
556  TrackCandidate cand(ownHits, seed, state, tk.seedRef());
557 
558  return cand;
559 }
560 
561 void TrackerTrackHitFilter::produceFromTrack(const edm::EventSetup &iSetup, const Track *itt , std::vector<TrackingRecHit *>&hits){
562 
563  // loop on tracks
564  hits.clear(); // extra safety
565 
566  for (trackingRecHit_iterator ith = itt->recHitsBegin(), edh = itt->recHitsEnd(); ith != edh; ++ith) {
567  const TrackingRecHit * hit = ith->get(); // ith is an iterator on edm::Ref to rechit
568 
569  DetId detid = hit->geographicalId();
570 
571  //check that the hit is a real hit and not a constraint
572  if(hit->isValid() && hit==0 && detid.rawId()==0) continue;
573 
574  int verdict=checkHit(iSetup,detid,hit);
575  if (verdict == 0) {
576  // just copy the hit
577  hits.push_back(hit->clone());
578  }
579  else if(verdict<-2){//hit rejected because did not pass the selections
580  // still, if replaceWithInactiveHits is true we have to put a new hit
582  hits.push_back(new InvalidTrackingRecHit(detid, TrackingRecHit::inactive));
583  }
584  }
585  else if(verdict==-2) hits.push_back(hit->clone());//hit not in the tracker
586  else if(verdict==-1){ //hit not valid
587  if (!stripAllInvalidHits_) {
588  hits.push_back(hit->clone());
589  }
590  }
591  } // loop on hits
592 
593 }//end TrackerTrackHitFilter::produceFromTrack
594 
595 
596 void TrackerTrackHitFilter::produceFromTrajectory(const edm::EventSetup &iSetup, const Trajectory *itt, std::vector<TrackingRecHit *>&hits){
597  hits.clear(); // extra safety
598  nOverlaps=0;
599 
600 
601  std::vector<TrajectoryMeasurement> tmColl =itt->measurements();
602 
603  //---OverlapBegin needed eventually for overlaps, but I must create them here in any case
604  const TrajectoryMeasurement* previousTM(0);
605  DetId previousId(0);
606  //int previousLayer(-1);
608 
609  int constrhits=0;
610 
611  for(std::vector<TrajectoryMeasurement>::const_iterator itTrajMeas = tmColl.begin(); itTrajMeas!=tmColl.end(); itTrajMeas++){
612  TransientTrackingRecHit::ConstRecHitPointer hitpointer = itTrajMeas->recHit();
613 
614  //check that the hit is a real hit and not a constraint
615  if(hitpointer->isValid() && hitpointer->hit()==0){constrhits++; continue;}
616 
617  const TrackingRecHit *hit=((*hitpointer).hit());
618  DetId detid = hit->geographicalId();
619  int verdict=checkHit(iSetup,detid,hit);
620 
621  if (verdict == 0) {
622  if( rejectLowAngleHits_ && !checkHitAngle(*itTrajMeas) ){//check angle of track on module if requested
623  verdict=-6;//override previous verdicts
624  }
625  }
626 
627  /*
628  //this has been included in checkHitAngle(*itTrajMeas)
629  if (verdict == 0) {
630  if( PXLcorrClusChargeCut_>0.0 && !checkPXLCorrClustCharge(*itTrajMeas) ){//check angle of track on module if requested
631  verdict=-7;//override previous verdicts
632  }
633  }
634  */
635 
636 
637  if(verdict==0){// Hit TAKEN !!!!
638 
639  if(tagOverlaps_){
640  //std::cout<<"Looking for overlaps in Run="<<iRun<<" , Event ="<<iEvt<<std::flush;
641 
642  int layer(layerFromId(detid));//layer 1-4=TIB, layer 5-10=TOB
643  int subDet = detid.subdetId();
644  //std::cout << " Check Subdet #" <<subDet << ", layer = " <<layer<<" stereo: "<< ((subDet > 2)?(SiStripDetId(detid).stereo()):2);
645 
646  if ( ( previousTM!=0 )&& (layer!=-1 )) {
647  //std::cout<<"A previous TM exists! "<<std::endl;
648  for (std::vector<TrajectoryMeasurement>::const_iterator itmCompare =itTrajMeas-1;itmCompare >= tmColl.begin() && itmCompare > itTrajMeas - 4;--itmCompare){
649 
650  DetId compareId = itmCompare->recHit()->geographicalId();
651  if ( subDet != compareId.subdetId() ||
652  layer != layerFromId(compareId)) break;
653  if (!itmCompare->recHit()->isValid()) continue;
654  if ( (subDet<=2) ||
655  (subDet > 2 && SiStripDetId(detid).stereo()==SiStripDetId(compareId).stereo()))
656  {//if either pixel or strip stereo module
657  // overlapHits.push_back(std::make_pair(&(*itmCompare),&(*itm)));
658  //std::cout<< "Adding pair "<< ((subDet >2)?(SiStripDetId(detid).stereo()):2)
659  // << " from SubDet = "<<subDet<<" , layer = " << layer<<" Run:"<<iRun<<"\tEv: "<<iEvt<<"\tId1: "<<compareId.rawId()<<"\tId2: "<<detid.rawId()<<std::endl;
660  // if(abs(compareId.rawId()-detid.rawId())==1)std::cout<<"These two are from the same det! Id1= "<<detid.rawId()<<" has stereo type "<<SiStripDetId(detid).stereo() <<"\tId2: "<<compareId.rawId()<<" has stereo type "<<SiStripDetId(compareId).stereo()<<std::endl;
662  // if(detid.rawId()<compareId.rawId()){
663  // std::cout<< "+++ "<< "\t"<<iRun<<"\t"<<iEvt<<"\t"<<detid.rawId()<<"\t"<<compareId.rawId()<<std::endl;
664  // }
665  //else std::cout<< "+++ "<< "\t"<<iRun<<"\t"<<iEvt<<"\t"<<compareId.rawId()<<"\t"<<detid.rawId()<<std::endl;
666 
667  nOverlaps++;
668  break;
669  }
670  }//end second loop on TM for overlap tagging
671 
672  }//end if ( (layer!=-1 )&&(acceptLayer[subDet]))
673 
674  previousTM = &(* itTrajMeas);
675  previousId = detid;
676  //previousLayer = layer;
677  }//end if look for overlaps
679 
680  hits.push_back(hit->clone()); //just copy it
681  }//end if HIT TAKEN
682  else if(verdict<-2){//hit rejected because did not pass the selections
683  // still, if replaceWithInactiveHits is true we have to put a new hit
685  hits.push_back(new InvalidTrackingRecHit(detid, TrackingRecHit::inactive));
686  }
687  }
688  else if(verdict==-2) hits.push_back(hit->clone());//hit not in the tracker
689  else if(verdict==-1){ //hit not valid
690  if (!stripAllInvalidHits_) {
691  hits.push_back(hit->clone());
692  }
693  }
694  } // loop on hits
695 
696  std::reverse(hits.begin(),hits.end());
697  // std::cout<<"Finished producefromTrajecotries. Nhits in final coll"<<hits.size() <<" Nconstraints="<<constrhits<<std::endl;
698 } //end TrackerTrackHitFilter::produceFromTrajectories
699 
701 
702 
703  int hitresult=0;
704  if (hit->isValid()) {
705 
706  if (detid.det() == DetId::Tracker) { // check for tracker hits
707  bool verdict = true;
708  // first check at structure level
709  for (std::vector<Rule>::const_iterator itr = rules_.begin(), edr = rules_.end(); itr != edr; ++itr) {
710  itr->apply(detid, verdict);
711  }
712 
713  // if the hit is good, check again at module level
714  if ( verdict){
715  if(std::binary_search(detsToIgnore_.begin(), detsToIgnore_.end(), detid.rawId())) {
716  hitresult=-4;
717  }
718  }
719  else hitresult=-3;
720  //if the hit is in the desired part of the det, check other things
721  if( hitresult==0 && rejectBadStoNHits_){
722  if(!checkStoN(iSetup, detid, hit))hitresult=-5;
723  }//end if S/N is ok
724  }//end hit in tracker
725  else hitresult =-2;
726  }//end hit is valid
727  else hitresult = -1; //invalid hit
728  return hitresult;
729 }//end TrackerTrackHitFilter::checkHit()
730 
731 
732 
733 bool TrackerTrackHitFilter::checkStoN(const edm::EventSetup &iSetup, const DetId &id, const TrackingRecHit* therechit){
734 
735  bool keepthishit=true;
736  // const uint32_t& recHitDetId = id.rawId();
737 
738  //check StoN only if subdet was set by the user
739  // int subdet_cnt=0;
740  int subdet_cnt=id.subdetId();
741 
742  // for(subdet_cnt=0;subdet_cnt<6; ++subdet_cnt){
743 
744  // if( subdetStoN_[subdet_cnt-1]&& (id.subdetId()==subdet_cnt) ){//check that hit is in a det belonging to a subdet where we decided to apply a S/N cut
745 
746 
747  if(subdet_cnt>2){ //SiStrip
748  if( subdetStoN_[subdet_cnt-1]){//check that hit is in a det belonging to a subdet where we decided to apply a S/N cut
749  const std::type_info &type = typeid(*therechit);
750  const SiStripCluster* cluster;
751  if (type == typeid(SiStripRecHit2D)) {
752  const SiStripRecHit2D* hit = dynamic_cast<const SiStripRecHit2D*>(therechit);
753  if (hit!=0) cluster = &*(hit->cluster());
754  else{
755  edm::LogError("TrackerTrackHitFilter")<< "TrackerTrackHitFilter::checkStoN : Unknown valid tracker hit in subdet " << id.subdetId()<< "(detID="<<id.rawId()<<")\n ";
756  keepthishit = false;
757  }
758  }
759  else if (type == typeid(SiStripRecHit1D)) {
760  const SiStripRecHit1D* hit = dynamic_cast<const SiStripRecHit1D*>(therechit);
761  if (hit!=0) cluster = &*(hit->cluster());
762  else{
763  edm::LogError("TrackerTrackHitFilter")<< "TrackerTrackHitFilter::checkStoN : Unknown valid tracker hit in subdet " << id.subdetId()<< "(detID="<<id.rawId()<<")\n ";
764  keepthishit = false;
765  }
766  }
767  //the following two cases should not happen anymore since CMSSW > 2_0_X because of hit splitting in stereo modules
768  //const SiStripMatchedRecHit2D* matchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(therechit);
769  //const ProjectedSiStripRecHit2D* unmatchedhit = dynamic_cast<const ProjectedSiStripRecHit2D*>(therechit);
770  else{
771  throw cms::Exception("Unknown RecHit Type") << "RecHit of type " << type.name() <<
772  " not supported. (use c++filt to demangle the name)";
773  }
774 
775  if(keepthishit){
776  SiStripClusterInfo clusterInfo = SiStripClusterInfo( *cluster, iSetup, id.rawId());
777  if ( (subdetStoNlowcut_[subdet_cnt-1]>0) && (clusterInfo.signalOverNoise() < subdetStoNlowcut_[subdet_cnt-1]) ) keepthishit = false;
778  if ( (subdetStoNhighcut_[subdet_cnt-1]>0) && (clusterInfo.signalOverNoise() > subdetStoNhighcut_[subdet_cnt-1]) ) keepthishit = false;
779  //if(!keepthishit)std::cout<<"Hit rejected because of bad S/N: "<<clusterInfo.signalOverNoise()<<std::endl;
780  }
781 
782  }//end if subdetStoN_[subdet_cnt]&&...
783 
784  }//end if subdet_cnt >2
785  else if (subdet_cnt<=2){//pixel
786  //pixels have naturally a very low noise (because of their low capacitance). So the S/N cut is
787  //irrelevant in this case. Leave it dummy
788  keepthishit = true;
789 
790  /**************
791  * Cut on cluster charge corr by angle embedded in the checkHitAngle() function
792  *
793  *************/
794 
795  if(checkPXLQuality_){
796  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(therechit);
797  if(pixelhit!=0){
798  //std::cout << "ClusterCharge=" <<std::flush<<pixelhit->cluster()->charge() << std::flush;
799  float xyprob=pixelhit->clusterProbability(0);//x-y combined log_e probability of the pixel cluster
800  //singl x- and y-prob not stored sicne CMSSW 3_9_0
801  float xychargeprob =pixelhit->clusterProbability(1);//xy-prob * charge prob
802  // float chargeprob =pixelhit->clusterProbability(2);//charge prob
803  bool haspassed_tplreco= pixelhit->hasFilledProb(); //the cluster was associted to a template
804  int qbin =pixelhit->qBin(); //R==meas_charge/pred_charge: Qbin=0 ->R>1.5 , =1->1<R<1.5 ,=2->0.85<R<1 ,
805  // Qbin=3->0.95*Qminpred<R<0.85 ,=4->, =5->meas_charge<0.95*Qminpred
806 
807  // if(haspassed_tplreco) std::cout<<" CLUSTPROB=\t"<<xprob<<"\t"<<yprob<<"\t"<<combprob<<"\t"<<qbin<<std::endl;
808  // else std::cout<<"CLUSTPROBNOTDEF=\t"<<xprob<<"\t"<<yprob<<"\t"<<combprob<<"\t"<<qbin<<std::endl;
809 
810  keepthishit = false;
811  // std::cout<<"yyyyy "<<qbin<<" "<<xprob<<" "<<yprob<<std::endl;
812  if( haspassed_tplreco && xyprob>pxlTPLProbXY_ && xychargeprob>pxlTPLProbXYQ_ && qbin>pxlTPLqBin_[0] && qbin<=pxlTPLqBin_[1] )keepthishit = true;
813 
814  }
815  else {std::cout<<"HIT IN PIXEL ("<<subdet_cnt <<") but PixelRecHit is EMPTY!!!"<<std::endl;}
816  }//end if check pixel quality flag
817  }
818  // else throw cms::Exception("TrackerTrackHitFilter") <<"Loop over subdetector out of range when applying the S/N cut: "<<subdet_cnt;
819 
820  // }//end loop on subdets
821 
822 
823  return keepthishit;
824 }//end CheckStoN
825 
826 
827 
829 
830  bool angle_ok=false;
831  bool corrcharge_ok=true;
833  /*
834  edm::LogDebug("TrackerTrackHitFilter")<<"TSOS parameters: ";
835  edm::LogDebug("TrackerTrackHitFilter") <<"Global momentum: "<<tsos.globalMomentum().x()<<" "<<tsos.globalMomentum().y()<<" "<<tsos.globalMomentum().z();
836  edm::LogDebug("TrackerTrackHitFilter") <<"Local momentum: "<<tsos.localMomentum().x()<<" "<<tsos.localMomentum().y()<<" "<<tsos.localMomentum().z();
837  edm::LogDebug("TrackerTrackHitFilter") <<"Track charge: " <<tsos.charge();
838  edm::LogDebug("TrackerTrackHitFilter")<<"Local position: " <<tsos.localPosition().x()<<" "<<tsos.localPosition().y()<<" "<<tsos.localPosition().z();
839  */
840  if(tsos.isValid()){
841  //check the angle of this tsos
842  float mom_x=tsos.localDirection().x();
843  float mom_y=tsos.localDirection().y();
844  float mom_z=tsos.localDirection().z();
845  //we took LOCAL momentum, i.e. respect to surface. Thus the plane is z=0
846  float angle=TMath::ASin(TMath::Abs(mom_z) / sqrt(pow(mom_x,2)+pow(mom_y,2)+pow(mom_z,2) ) );
847  if(!rejectLowAngleHits_ || angle>=TrackAngleCut_) angle_ok=true;// keep this hit
848  // else std::cout<<"Hit rejected because angle is "<< angle<<" ( <"<<TrackAngleCut_<<" )"<<std::endl;
849 
850  if(angle_ok && PXLcorrClusChargeCut_>0.0){
851  //
852  //get the hit from the TM and check that it is in the pixel
854  if(hitpointer->isValid()){
855  const TrackingRecHit *hit=(*hitpointer).hit();
856  if( (hit->geographicalId()).subdetId()<=2 ){//do it only for pixel hits
857  corrcharge_ok=false;
858  float clust_alpha= atan2( mom_z, mom_x );
859  float clust_beta= atan2( mom_z, mom_y );
860 
861  //Now get the cluster charge
862 
863  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
864  float clust_charge=pixelhit->cluster()->charge();
865  float corr_clust_charge=clust_charge * sqrt( 1.0 / ( 1.0/pow( tan(clust_alpha), 2 ) +
866  1.0/pow( tan(clust_beta ), 2 ) +
867  1.0 )
868  );
869  //std::cout<<"xxxxx "<<clust_charge<<" "<<corr_clust_charge<<" " <<pixelhit->qBin()<<" "<<pixelhit->clusterProbability(1)<<" "<<pixelhit->clusterProbability(2)<< std::endl;
870  if(corr_clust_charge>PXLcorrClusChargeCut_){
871  corrcharge_ok=true;
872  }
873  } //end if hit is in pixel
874  }//end if hit is valid
875 
876  }//check corr cluster charge for pixel hits
877 
878  }//end if TSOS is valid
879  else{
880  edm::LogWarning("TrackerTrackHitFilter") <<"TSOS not valid ! Impossible to calculate track angle.";
881  }
882 
883  return angle_ok&&corrcharge_ok;
884 }//end TrackerTrackHitFilter::checkHitAngle
885 
886 
888  /*
889  Code taken from DPGAnalysis/SiPixelTools/plugins/PixelNtuplizer_RealData.cc
890  */
891 
892  bool corrcharge_ok=false;
893  //get the hit from the TM and check that it is in the pixel
895  if(!hitpointer->isValid()) return corrcharge_ok;
896  const TrackingRecHit *hit=(*hitpointer).hit();
897  if( (hit->geographicalId()).subdetId()>2 ){//SiStrip hit, skip
898  return corrcharge_ok;
899  }
900 
901  TrajectoryStateOnSurface tsos = meas.updatedState();
902  if(tsos.isValid()){
903  float mom_x=tsos.localDirection().x();
904  float mom_y=tsos.localDirection().y();
905  float mom_z=tsos.localDirection().z();
906  float clust_alpha= atan2( mom_z, mom_x );
907  float clust_beta= atan2( mom_z, mom_y );
908 
909  //Now get the cluster charge
910 
911  const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(hit);
912  float clust_charge=pixelhit->cluster()->charge();
913  float corr_clust_charge=clust_charge * sqrt( 1.0 / ( 1.0/pow( tan(clust_alpha), 2 ) +
914  1.0/pow( tan(clust_beta ), 2 ) +
915  1.0 )
916  );
917  if(corr_clust_charge>PXLcorrClusChargeCut_)corrcharge_ok=true;
918 
919  }//end if TSOS is valid
920  return corrcharge_ok;
921 
922 }//end TrackerTrackHitFilter::checkPXLCorrClustCharge
923 
924 
925 
927 {
928  if ( id.subdetId()== int(PixelSubdetector::PixelBarrel) ) {
929  PXBDetId tobId(id);
930  return tobId.layer();
931  }
932  else if ( id.subdetId()== int(PixelSubdetector::PixelEndcap) ) {
933  PXFDetId tobId(id);
934  return tobId.disk() + (3*(tobId.side()-1));
935  }
936  else if ( id.subdetId()==StripSubdetector::TIB ) {
937  TIBDetId tibId(id);
938  return tibId.layer();
939  }
940  else if ( id.subdetId()==StripSubdetector::TOB ) {
941  TOBDetId tobId(id);
942  return tobId.layer();
943  }
944  else if ( id.subdetId()==StripSubdetector::TEC ) {
945  TECDetId tobId(id);
946  return tobId.wheel() + (9*(tobId.side()-1));
947  }
948  else if ( id.subdetId()==StripSubdetector::TID ) {
949  TIDDetId tobId(id);
950  return tobId.wheel() + (3*(tobId.side()-1));
951  }
952  return -1;
953 
954 }
955 
956 }} //namespaces
957 
958 
959 // ========= MODULE DEF ==============
962 
RunNumber_t run() const
Definition: EventID.h:42
ClusterRef cluster() const
type
Definition: HCALResponse.h:22
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
bool hasFilledProb() const
float clusterProbability(unsigned int flags=0) const
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:150
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
uint32_t stereo() const
Definition: SiStripDetId.h:162
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
LocalVector localDirection() const
TrackerTrackHitFilter(const edm::ParameterSet &iConfig)
std::vector< TrackCandidate > TrackCandidateCollection
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
T y() const
Definition: PV3DBase.h:62
bool checkStoN(const edm::EventSetup &iSetup, const DetId &id, const TrackingRecHit *therechit)
ConstRecHitPointer recHit() const
PropagationDirection
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
void produceFromTrajectory(const edm::EventSetup &iSetup, const Trajectory *itt, std::vector< TrackingRecHit * > &hits)
unsigned int side() const
positive or negative id
Definition: TECDetId.h:47
void apply(DetId detid, bool &verdict) const
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
U second(std::pair< T, U > const &p)
bool checkPXLCorrClustCharge(const TrajectoryMeasurement &meas)
DataContainer const & measurements() const
Definition: Trajectory.h:203
float signalOverNoise() const
void push_back(D *&d)
Definition: OwnVector.h:273
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
T sqrt(T t)
Definition: SSEVec.h:46
void produceFromTrack(const edm::EventSetup &iSetup, const Track *itt, std::vector< TrackingRecHit * > &hits)
T z() const
Definition: PV3DBase.h:63
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
ClusterRef cluster() const
TrajectoryStateOnSurface updatedState() const
#define end
Definition: vmac.h:38
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
bool first
Definition: L1TdeRCT.cc:94
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
virtual TrackingRecHit * clone() const =0
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
int qBin() const
Definition: SiPixelRecHit.h:96
Definition: DetId.h:20
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:112
edm::ESHandle< TransientTrackingRecHitBuilder > theBuilder
int checkHit(const edm::EventSetup &iSetup, const DetId &detid, const TrackingRecHit *hit)
tuple tracks
Definition: testEve_cfg.py:39
edm::ESHandle< TrackerGeometry > theGeometry
unsigned int side() const
positive or negative id
Definition: TIDDetId.h:45
ClusterRef cluster() const
Definition: SiPixelRecHit.h:41
const T & get() const
Definition: EventSetup.h:55
key_type key() const
Accessor for product key.
Definition: Ref.h:266
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
bool isValid() const
void parseStoN(const std::string &str)
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
char state
Definition: procUtils.cc:75
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
edm::EventID id() const
Definition: EventBase.h:56
#define begin
Definition: vmac.h:31
bool checkHitAngle(const TrajectoryMeasurement &meas)
unsigned int side() const
positive or negative id
Definition: PXFDetId.h:38
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:105
tuple cout
Definition: gather_cfg.py:121
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
DetId geographicalId() const
Detector det() const
get the detector field from this detid
Definition: DetId.h:37
TrackCandidate makeCandidate(const reco::Track &tk, std::vector< TrackingRecHit * >::iterator hitsBegin, std::vector< TrackingRecHit * >::iterator hitsEnd)
edm::ESHandle< MagneticField > theMagField
T x() const
Definition: PV3DBase.h:61
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void reserve(size_t)
Definition: OwnVector.h:267
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50
Our base class.
Definition: SiPixelRecHit.h:22
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65