CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTracker/FinalTrackSelectors/src/AnalyticalTrackSelector.cc

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/interface/AnalyticalTrackSelector.h"
00002 #include "RecoTracker/FinalTrackSelectors/src/MultiTrackSelector.h"
00003 
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include <Math/DistFunc.h>
00006 #include "TMath.h"
00007 
00008 using reco::modules::AnalyticalTrackSelector;
00009 
00010 AnalyticalTrackSelector::AnalyticalTrackSelector( const edm::ParameterSet & cfg ) : MultiTrackSelector( )
00011 {
00012     //Spoof the pset for each track selector!
00013     //Size is always 1!!!
00014     qualityToSet_.reserve(1);
00015     vtxNumber_.reserve(1);
00016     vertexCut_.reserve(1);
00017     res_par_.reserve(1);
00018     chi2n_par_.reserve(1);
00019     chi2n_no1Dmod_par_.reserve(1);
00020     d0_par1_.reserve(1);
00021     dz_par1_.reserve(1);
00022     d0_par2_.reserve(1);
00023     dz_par2_.reserve(1);
00024     applyAdaptedPVCuts_.reserve(1);
00025     max_d0_.reserve(1);
00026     max_z0_.reserve(1);
00027     nSigmaZ_.reserve(1);
00028     min_layers_.reserve(1);
00029     min_3Dlayers_.reserve(1);
00030     max_lostLayers_.reserve(1);
00031     min_hits_bypass_.reserve(1);
00032     applyAbsCutsIfNoPV_.reserve(1);
00033     max_d0NoPV_.reserve(1);
00034     max_z0NoPV_.reserve(1);
00035     preFilter_.reserve(1);
00036     max_relpterr_.reserve(1);
00037     min_nhits_.reserve(1);
00038     max_minMissHitOutOrIn_.reserve(1);
00039     max_lostHitFraction_.reserve(1);
00040     min_eta_.reserve(1);
00041     max_eta_.reserve(1);
00042 
00043     src_ = cfg.getParameter<edm::InputTag>( "src" );
00044     beamspot_ = cfg.getParameter<edm::InputTag>( "beamspot" );
00045     useVertices_ = cfg.getParameter<bool>( "useVertices" );
00046     useVtxError_ = cfg.getParameter<bool>( "useVtxError" );
00047     vertices_ = useVertices_ ? cfg.getParameter<edm::InputTag>( "vertices" ) : edm::InputTag("NONE");
00048     copyExtras_ = cfg.getUntrackedParameter<bool>("copyExtras", false);
00049     copyTrajectories_ = cfg.getUntrackedParameter<bool>("copyTrajectories", false);
00050     
00051     qualityToSet_.push_back( TrackBase::undefQuality );
00052     // parameters for vertex selection
00053     vtxNumber_.push_back( useVertices_ ? cfg.getParameter<int32_t>("vtxNumber") : 0 );
00054     vertexCut_.push_back( useVertices_ ? cfg.getParameter<std::string>("vertexCut") : "");
00055     //  parameters for adapted optimal cuts on chi2 and primary vertex compatibility
00056     res_par_.push_back(cfg.getParameter< std::vector<double> >("res_par") );
00057     chi2n_par_.push_back( cfg.getParameter<double>("chi2n_par") );
00058     chi2n_no1Dmod_par_.push_back( cfg.getParameter<double>("chi2n_no1Dmod_par") );
00059     d0_par1_.push_back(cfg.getParameter< std::vector<double> >("d0_par1"));
00060     dz_par1_.push_back(cfg.getParameter< std::vector<double> >("dz_par1"));
00061     d0_par2_.push_back(cfg.getParameter< std::vector<double> >("d0_par2"));
00062     dz_par2_.push_back(cfg.getParameter< std::vector<double> >("dz_par2"));
00063 
00064     // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
00065     applyAdaptedPVCuts_.push_back(cfg.getParameter<bool>("applyAdaptedPVCuts"));
00066     // Impact parameter absolute cuts.
00067     max_d0_.push_back(cfg.getParameter<double>("max_d0"));
00068     max_z0_.push_back(cfg.getParameter<double>("max_z0"));
00069     nSigmaZ_.push_back(cfg.getParameter<double>("nSigmaZ"));
00070     // Cuts on numbers of layers with hits/3D hits/lost hits.
00071     min_layers_.push_back(cfg.getParameter<uint32_t>("minNumberLayers") );
00072     min_3Dlayers_.push_back(cfg.getParameter<uint32_t>("minNumber3DLayers") );
00073     max_lostLayers_.push_back(cfg.getParameter<uint32_t>("maxNumberLostLayers"));
00074     min_hits_bypass_.push_back(cfg.getParameter<uint32_t>("minHitsToBypassChecks"));
00075     max_relpterr_.push_back(cfg.getParameter<double>("max_relpterr"));
00076     min_nhits_.push_back(cfg.getParameter<uint32_t>("min_nhits"));
00077     max_minMissHitOutOrIn_.push_back(
00078         cfg.existsAs<int32_t>("max_minMissHitOutOrIn") ? 
00079         cfg.getParameter<int32_t>("max_minMissHitOutOrIn") : 99);
00080     max_lostHitFraction_.push_back(
00081         cfg.existsAs<double>("max_lostHitFraction") ?
00082         cfg.getParameter<double>("max_lostHitFraction") : 1.0);
00083     min_eta_.push_back(cfg.getParameter<double>("min_eta"));
00084     max_eta_.push_back(cfg.getParameter<double>("max_eta"));
00085 
00086     // Flag to apply absolute cuts if no PV passes the selection
00087     applyAbsCutsIfNoPV_.push_back(cfg.getParameter<bool>("applyAbsCutsIfNoPV"));
00088     keepAllTracks_.push_back( cfg.exists("keepAllTracks") ?
00089                               cfg.getParameter<bool>("keepAllTracks") :
00090                               false ); 
00091  
00092     setQualityBit_.push_back( false );
00093     std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
00094     
00095     if(d0_par1_[0].size()!=2 || dz_par1_[0].size()!=2 || d0_par2_[0].size()!=2 || dz_par2_[0].size()!=2)
00096     {
00097       edm::LogError("MisConfiguration")<<"vector of size less then 2";
00098       throw; 
00099     }
00100 
00101     if (cfg.exists("qualityBit")) {
00102       std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
00103       if (qualityStr != "") {
00104         setQualityBit_[0] = true;
00105         qualityToSet_ [0] = TrackBase::qualityByName(cfg.getParameter<std::string>("qualityBit"));
00106       }
00107     }
00108   
00109     if (keepAllTracks_[0] && !setQualityBit_[0]) throw cms::Exception("Configuration") << 
00110       "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
00111     if (setQualityBit_[0] && (qualityToSet_[0] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
00112       "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
00113     if (applyAbsCutsIfNoPV_[0]) {
00114       max_d0NoPV_.push_back(cfg.getParameter<double>("max_d0NoPV"));
00115       max_z0NoPV_.push_back(cfg.getParameter<double>("max_z0NoPV"));
00116     }
00117     else{//dummy values
00118       max_d0NoPV_.push_back(0.);
00119       max_z0NoPV_.push_back(0.);
00120     }
00121     
00122     std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
00123     produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
00124     if (copyExtras_) {
00125       produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
00126       produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
00127     }
00128     if (copyTrajectories_) {
00129       produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
00130       produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
00131     }
00132 }
00133 
00134 AnalyticalTrackSelector::~AnalyticalTrackSelector() {
00135 }
00136 
00137 void AnalyticalTrackSelector::produce( edm::Event& evt, const edm::EventSetup& es ) 
00138 {
00139   using namespace std; 
00140   using namespace edm;
00141   using namespace reco;
00142 
00143   Handle<TrackCollection> hSrcTrack;
00144   Handle< vector<Trajectory> > hTraj;
00145   Handle< vector<Trajectory> > hTrajP;
00146   Handle< TrajTrackAssociationCollection > hTTAss;
00147 
00148   // looking for the beam spot
00149   edm::Handle<reco::BeamSpot> hBsp;
00150   evt.getByLabel(beamspot_, hBsp);
00151   reco::BeamSpot vertexBeamSpot;
00152   vertexBeamSpot = *hBsp;
00153         
00154   // Select good primary vertices for use in subsequent track selection
00155   edm::Handle<reco::VertexCollection> hVtx;
00156   std::vector<Point> points;
00157   std::vector<float> vterr, vzerr;
00158   if (useVertices_) {
00159       evt.getByLabel(vertices_, hVtx);
00160       selectVertices(0,*hVtx, points, vterr, vzerr);
00161       // Debug 
00162       LogDebug("SelectVertex") << points.size() << " good pixel vertices";
00163   }
00164 
00165   // Get tracks 
00166   evt.getByLabel( src_, hSrcTrack );
00167 
00168   selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
00169   rTracks_ = evt.getRefBeforePut<TrackCollection>();      
00170   if (copyExtras_) {
00171     selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
00172     selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
00173     rHits_ = evt.getRefBeforePut<TrackingRecHitCollection>();
00174     rTrackExtras_ = evt.getRefBeforePut<TrackExtraCollection>();
00175   }
00176 
00177   if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
00178 
00179   // Loop over tracks
00180   size_t current = 0;
00181   for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
00182     const Track & trk = * it;
00183     // Check if this track passes cuts
00184 
00185     LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
00186     
00187     bool ok = select(0,vertexBeamSpot, trk, points, vterr, vzerr);
00188     if (!ok) {
00189 
00190       LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
00191 
00192       if (copyTrajectories_) trackRefs_[current] = reco::TrackRef();
00193       if (!keepAllTracks_[0]) continue;
00194     }
00195     LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
00196     selTracks_->push_back( Track( trk ) ); // clone and store
00197     if (ok && setQualityBit_[0]) {
00198       selTracks_->back().setQuality(qualityToSet_[0]);
00199       if (!points.empty()) {
00200         if (qualityToSet_[0]==TrackBase::loose) selTracks_->back().setQuality(TrackBase::looseSetWithPV);
00201         else if (qualityToSet_[0]==TrackBase::highPurity) selTracks_->back().setQuality(TrackBase::highPuritySetWithPV);
00202       }
00203     }
00204     if (copyExtras_) {
00205       // TrackExtras
00206       selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
00207                                               trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
00208                                               trk.outerStateCovariance(), trk.outerDetId(),
00209                                               trk.innerStateCovariance(), trk.innerDetId(),
00210                                               trk.seedDirection(), trk.seedRef() ) );
00211       selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
00212       TrackExtra & tx = selTrackExtras_->back();
00213       tx.setResiduals(trk.residuals());
00214       // TrackingRecHits
00215       for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
00216         selHits_->push_back( (*hit)->clone() );
00217         tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
00218       }
00219     }
00220     if (copyTrajectories_) {
00221       trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
00222     }
00223   }
00224   if ( copyTrajectories_ ) {
00225     Handle< vector<Trajectory> > hTraj;
00226     Handle< TrajTrackAssociationCollection > hTTAss;
00227     evt.getByLabel(src_, hTTAss);
00228     evt.getByLabel(src_, hTraj);
00229     selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>()); 
00230     rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
00231     selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
00232     for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
00233       Ref< vector<Trajectory> > trajRef(hTraj, i);
00234       TrajTrackAssociationCollection::const_iterator match = hTTAss->find(trajRef);
00235       if (match != hTTAss->end()) {
00236         const Ref<TrackCollection> &trkRef = match->val; 
00237         short oldKey = static_cast<short>(trkRef.key());
00238         if (trackRefs_[oldKey].isNonnull()) {
00239           selTrajs_->push_back( Trajectory(*trajRef) );
00240           selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
00241         }
00242       }
00243     }
00244   }
00245 
00246   static const std::string emptyString;
00247   evt.put(selTracks_);
00248   if (copyExtras_ ) {
00249     evt.put(selTrackExtras_); 
00250     evt.put(selHits_);
00251   }
00252   if ( copyTrajectories_ ) {
00253     evt.put(selTrajs_);
00254     evt.put(selTTAss_);
00255   }
00256 }