CMS 3D CMS Logo

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