CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoTracker/FinalTrackSelectors/src/CosmicTrackSelector.cc

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/interface/CosmicTrackSelector.h"
00002 
00003 #include <Math/DistFunc.h>
00004 #include "TMath.h"
00005 
00006 using reco::modules::CosmicTrackSelector;
00007 
00008 CosmicTrackSelector::CosmicTrackSelector( const edm::ParameterSet & cfg ) :
00009   src_( cfg.getParameter<edm::InputTag>( "src" ) ),
00010   beamspot_( cfg.getParameter<edm::InputTag>( "beamspot" ) ),
00011   copyExtras_(cfg.getUntrackedParameter<bool>("copyExtras", false)),
00012   copyTrajectories_(cfg.getUntrackedParameter<bool>("copyTrajectories", false)),
00013   keepAllTracks_( cfg.exists("keepAllTracks") ?
00014                   cfg.getParameter<bool>("keepAllTracks") :
00015                   false ),  // as this is what you expect from a well behaved selector
00016   setQualityBit_( false ),
00017   qualityToSet_( TrackBase::undefQuality ),
00018   chi2n_par_( cfg.getParameter<double>("chi2n_par") ),
00019   // Impact parameter absolute cuts.
00020   max_d0_(cfg.getParameter<double>("max_d0")),
00021   max_z0_(cfg.getParameter<double>("max_z0")),
00022   // Track parameter cuts.
00023   min_pt_(cfg.getParameter<double>("min_pt")),
00024   max_eta_(cfg.getParameter<double>("max_eta")),
00025   // Cut on number of valid hits
00026   min_nHit_(cfg.getParameter<uint32_t>("min_nHit")),
00027   // Cut on number of valid hits
00028   min_nPixelHit_(cfg.getParameter<uint32_t>("min_nPixelHit")),
00029   // Cuts on numbers of layers with hits/3D hits/lost hits.
00030   min_layers_(cfg.getParameter<uint32_t>("minNumberLayers") ),
00031   min_3Dlayers_(cfg.getParameter<uint32_t>("minNumber3DLayers") ),
00032   max_lostLayers_(cfg.getParameter<uint32_t>("maxNumberLostLayers") )
00033 {
00034   if (cfg.exists("qualityBit")) {
00035     std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
00036     if (qualityStr != "") {
00037       setQualityBit_ = true;
00038       qualityToSet_  = TrackBase::qualityByName(cfg.getParameter<std::string>("qualityBit"));
00039     }
00040   }
00041   if (keepAllTracks_ && !setQualityBit_) throw cms::Exception("Configuration") << 
00042     "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
00043   if (setQualityBit_ && (qualityToSet_ == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
00044     "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
00045   
00046   std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
00047   produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
00048   if (copyExtras_) {
00049     produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
00050     produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
00051   }
00052   if (copyTrajectories_) {
00053     produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
00054     produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
00055   }
00056   
00057 }
00058 
00059 CosmicTrackSelector::~CosmicTrackSelector() {
00060 }
00061 
00062 void CosmicTrackSelector::produce( edm::Event& evt, const edm::EventSetup& es ) 
00063 {
00064   using namespace std; 
00065   using namespace edm;
00066   using namespace reco;
00067   
00068   Handle<TrackCollection> hSrcTrack;
00069   Handle< vector<Trajectory> > hTraj;
00070   Handle< vector<Trajectory> > hTrajP;
00071   Handle< TrajTrackAssociationCollection > hTTAss;
00072   
00073   // looking for the beam spot
00074   edm::Handle<reco::BeamSpot> hBsp;
00075   evt.getByLabel(beamspot_, hBsp);
00076   reco::BeamSpot vertexBeamSpot;
00077   vertexBeamSpot = *hBsp;
00078   
00079   // Get tracks 
00080   evt.getByLabel( src_, hSrcTrack );
00081   
00082   selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
00083   rTracks_ = evt.getRefBeforePut<TrackCollection>();      
00084   if (copyExtras_) {
00085     selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
00086     selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
00087     rHits_ = evt.getRefBeforePut<TrackingRecHitCollection>();
00088     rTrackExtras_ = evt.getRefBeforePut<TrackExtraCollection>();
00089   }
00090   
00091   if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
00092   
00093   // Loop over tracks
00094   size_t current = 0;
00095   for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
00096     const Track & trk = * it;
00097     // Check if this track passes cuts
00098     bool ok = select(vertexBeamSpot, trk);
00099     if (!ok) {
00100       if (copyTrajectories_) trackRefs_[current] = reco::TrackRef();
00101       if (!keepAllTracks_) continue;
00102     }
00103     selTracks_->push_back( Track( trk ) ); // clone and store
00104     if (ok && setQualityBit_) selTracks_->back().setQuality(qualityToSet_);
00105     if (copyExtras_) {
00106       // TrackExtras
00107       selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
00108                                               trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
00109                                               trk.outerStateCovariance(), trk.outerDetId(),
00110                                               trk.innerStateCovariance(), trk.innerDetId(),
00111                                               trk.seedDirection(), trk.seedRef() ) );
00112       selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
00113       TrackExtra & tx = selTrackExtras_->back();
00114       tx.setResiduals(trk.residuals());
00115       // TrackingRecHits
00116       for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
00117         selHits_->push_back( (*hit)->clone() );
00118         tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
00119       }
00120     }
00121     if (copyTrajectories_) {
00122       trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
00123     }
00124     }
00125   if ( copyTrajectories_ ) {
00126     Handle< vector<Trajectory> > hTraj;
00127     Handle< TrajTrackAssociationCollection > hTTAss;
00128     evt.getByLabel(src_, hTTAss);
00129     evt.getByLabel(src_, hTraj);
00130     selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>()); 
00131     rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
00132     selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
00133     for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
00134       Ref< vector<Trajectory> > trajRef(hTraj, i);
00135       TrajTrackAssociationCollection::const_iterator match = hTTAss->find(trajRef);
00136       if (match != hTTAss->end()) {
00137         const Ref<TrackCollection> &trkRef = match->val; 
00138         short oldKey = static_cast<short>(trkRef.key());
00139         if (trackRefs_[oldKey].isNonnull()) {
00140           selTrajs_->push_back( Trajectory(*trajRef) );
00141           selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
00142         }
00143       }
00144     }
00145   }
00146   
00147   static const std::string emptyString;
00148   evt.put(selTracks_);
00149   if (copyExtras_ ) {
00150     evt.put(selTrackExtras_); 
00151     evt.put(selHits_);
00152   }
00153   if ( copyTrajectories_ ) {
00154     evt.put(selTrajs_);
00155     evt.put(selTTAss_);
00156   }
00157 }
00158 
00159 
00160 bool CosmicTrackSelector::select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk) {
00161   // Decide if the given track passes selection cuts.
00162   
00163   using namespace std; 
00164   
00165   // Cuts on numbers of layers with hits/3D hits/lost hits.
00166   uint32_t nlayers     = tk.hitPattern().trackerLayersWithMeasurement();
00167   uint32_t nlayers3D   = tk.hitPattern().pixelLayersWithMeasurement() +
00168     tk.hitPattern().numberOfValidStripLayersWithMonoAndStereo();
00169   uint32_t nlayersLost = tk.hitPattern().trackerLayersWithoutMeasurement();
00170   
00171   // Get the number of valid hits and PixelHits
00172   uint32_t nHit = 0;
00173   uint32_t nPixelHit = 0;
00174   for ( trackingRecHit_iterator recHit = tk.recHitsBegin(); recHit != tk.recHitsEnd(); ++recHit ) {
00175     if ( !((*recHit)->isValid()) ) continue; 
00176     ++nHit;
00177     DetId id((*recHit)->geographicalId());
00178     if ( (unsigned int)id.subdetId() == PixelSubdetector::PixelBarrel 
00179          || (unsigned int)id.subdetId() == PixelSubdetector::PixelEndcap )
00180       ++nPixelHit;
00181   }  
00182   
00183   // Cut on the number of valid hits
00184   if (nHit < min_nHit_) return false;
00185   // Cut on the number of valid Pixel hits
00186   if (nPixelHit < min_nPixelHit_) return false;
00187   if (nlayers < min_layers_) return false;
00188   if (nlayers3D < min_3Dlayers_) return false;
00189   if (nlayersLost > max_lostLayers_) return false;
00190   
00191   // Get track parameters
00192   double pt = tk.pt(),eta = tk.eta(), chi2n =  tk.normalizedChi2();
00193   double d0 = -tk.dxy(vertexBeamSpot.position()), dz = tk.dz();
00194   
00195   // Absolute cuts on all tracks impact parameters with respect to beam-spot.
00196   if (abs(d0) > max_d0_) return false;
00197   if (abs(dz) > max_z0_) return false;
00198   
00199   // optimized cuts adapted to the track eta, pt and  chiquare/ndof 
00200   if (abs(eta) > max_eta_) return false;
00201   if (pt < min_pt_) return false;
00202   if (chi2n > chi2n_par_*nlayers) return false;
00203 
00204   
00205   else    
00206     return true;
00207   
00208 }
00209 
00210