CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/CommonTools/RecoAlgos/interface/RecoTrackSelector.h

Go to the documentation of this file.
00001 #ifndef RecoSelectors_RecoTrackSelector_h
00002 #define RecoSelectors_RecoTrackSelector_h
00003 /* \class RecoTrackSelector
00004  *
00005  * \author Giuseppe Cerati, INFN
00006  *
00007  *  $Date: 2010/02/11 00:10:50 $
00008  *  $Revision: 1.3 $
00009  *
00010  */
00011 #include "DataFormats/TrackReco/interface/Track.h"
00012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00013 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00014 #include "FWCore/Utilities/interface/InputTag.h"
00015 
00016 class RecoTrackSelector {
00017  public:
00018   typedef reco::TrackCollection collection;
00019   typedef std::vector<const reco::Track *> container;
00020   typedef container::const_iterator const_iterator;
00021 
00023   RecoTrackSelector() {}
00024   RecoTrackSelector ( const edm::ParameterSet & cfg ) :
00025     ptMin_(cfg.getParameter<double>("ptMin")),
00026     minRapidity_(cfg.getParameter<double>("minRapidity")),
00027     maxRapidity_(cfg.getParameter<double>("maxRapidity")),
00028     tip_(cfg.getParameter<double>("tip")),
00029     lip_(cfg.getParameter<double>("lip")),
00030     minHit_(cfg.getParameter<int>("minHit")),
00031     min3DHit_(cfg.getParameter<int>("min3DHit")),
00032     maxChi2_(cfg.getParameter<double>("maxChi2")),
00033     bsSrc_(cfg.getParameter<edm::InputTag>("beamSpot")) 
00034     {
00035       std::vector<std::string> quality = cfg.getParameter<std::vector<std::string> >("quality");
00036       for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
00037       std::vector<std::string> algorithm = cfg.getParameter<std::vector<std::string> >("algorithm");
00038       for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
00039     }
00040   
00041   RecoTrackSelector ( double ptMin, double minRapidity, double maxRapidity,
00042                       double tip, double lip, int minHit, int min3DHit, double maxChi2, 
00043                       std::vector<std::string> quality , std::vector<std::string> algorithm ) :
00044     ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
00045     tip_( tip ), lip_( lip ), minHit_( minHit ), min3DHit_( min3DHit), maxChi2_( maxChi2 ) 
00046     { 
00047       for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
00048       for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
00049     }
00050 
00051   const_iterator begin() const { return selected_.begin(); }
00052   const_iterator end() const { return selected_.end(); }
00053   
00054   void select( const edm::Handle<collection>& c, const edm::Event & event, const edm::EventSetup&) {
00055     selected_.clear();
00056     edm::Handle<reco::BeamSpot> beamSpot;
00057     event.getByLabel(bsSrc_,beamSpot); 
00058     for( reco::TrackCollection::const_iterator trk = c->begin(); 
00059          trk != c->end(); ++ trk )
00060       if ( operator()(*trk,beamSpot.product()) ) {
00061         selected_.push_back( & * trk );
00062       }
00063   }
00064 
00066   bool operator()( const reco::Track & t, const reco::BeamSpot* bs) {
00067     bool quality_ok = true;
00068     if (quality_.size()!=0) {
00069       quality_ok = false;
00070       for (unsigned int i = 0; i<quality_.size();++i) {
00071         if (t.quality(quality_[i])){
00072           quality_ok = true;
00073           break;          
00074         }
00075       }
00076     }
00077     bool algo_ok = true;
00078     if (algorithm_.size()!=0) {
00079       if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
00080     }
00081     return
00082       (t.hitPattern().trackerLayersWithMeasurement() >= minHit_ &&
00083        t.hitPattern().pixelLayersWithMeasurement() +
00084        t.hitPattern().numberOfValidStripLayersWithMonoAndStereo() >= min3DHit_ &&
00085        fabs(t.pt()) >= ptMin_ &&
00086        t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
00087        fabs(t.dxy(bs->position())) <= tip_ &&
00088        fabs(t.dsz(bs->position())) <= lip_  &&
00089        t.normalizedChi2()<=maxChi2_ &&
00090        quality_ok &&
00091        algo_ok);
00092   }
00093 
00094   size_t size() const { return selected_.size(); }
00095   
00096  protected:
00097   double ptMin_;
00098   double minRapidity_;
00099   double maxRapidity_;
00100   double tip_;
00101   double lip_;
00102   int    minHit_;
00103   int    min3DHit_;
00104   double maxChi2_;
00105   std::vector<reco::TrackBase::TrackQuality> quality_;
00106   std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
00107   edm::InputTag bsSrc_;
00108   container selected_;
00109 };
00110 
00111 #endif