CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/src/MultiTrackSelector.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/Common/interface/ValueMap.h"
00005 #include <Math/DistFunc.h>
00006 #include "TMath.h"
00007 
00008 using reco::modules::MultiTrackSelector;
00009 
00010 MultiTrackSelector::MultiTrackSelector( const edm::ParameterSet & cfg ) :
00011   src_( cfg.getParameter<edm::InputTag>( "src" ) ),
00012   beamspot_( cfg.getParameter<edm::InputTag>( "beamspot" ) ),
00013   useVertices_( cfg.getParameter<bool>( "useVertices" ) ),
00014   useVtxError_( cfg.getParameter<bool>( "useVtxError" ) ),
00015   vertices_( useVertices_ ? cfg.getParameter<edm::InputTag>( "vertices" ) : edm::InputTag("NONE"))
00016   // now get the pset for each selector
00017 {
00018 
00019   std::vector<edm::ParameterSet> trkSelectors( cfg.getParameter<std::vector< edm::ParameterSet> >("trackSelectors") );
00020   qualityToSet_.reserve(trkSelectors.size());
00021   vtxNumber_.reserve(trkSelectors.size());
00022   vertexCut_.reserve(trkSelectors.size());
00023   res_par_.reserve(trkSelectors.size());
00024   chi2n_par_.reserve(trkSelectors.size());
00025   chi2n_no1Dmod_par_.reserve(trkSelectors.size());
00026   d0_par1_.reserve(trkSelectors.size());
00027   dz_par1_.reserve(trkSelectors.size());
00028   d0_par2_.reserve(trkSelectors.size());
00029   dz_par2_.reserve(trkSelectors.size());
00030   applyAdaptedPVCuts_.reserve(trkSelectors.size());
00031   max_d0_.reserve(trkSelectors.size());
00032   max_z0_.reserve(trkSelectors.size());
00033   nSigmaZ_.reserve(trkSelectors.size());
00034   min_layers_.reserve(trkSelectors.size());
00035   min_3Dlayers_.reserve(trkSelectors.size());
00036   max_lostLayers_.reserve(trkSelectors.size());
00037   min_hits_bypass_.reserve(trkSelectors.size());
00038   applyAbsCutsIfNoPV_.reserve(trkSelectors.size());
00039   max_d0NoPV_.reserve(trkSelectors.size());
00040   max_z0NoPV_.reserve(trkSelectors.size());
00041   preFilter_.reserve(trkSelectors.size());
00042   max_relpterr_.reserve(trkSelectors.size());
00043   min_nhits_.reserve(trkSelectors.size());
00044 
00045   for ( unsigned int i=0; i<trkSelectors.size(); i++) {
00046 
00047     qualityToSet_.push_back( TrackBase::undefQuality );
00048     // parameters for vertex selection
00049     vtxNumber_.push_back( useVertices_ ? trkSelectors[i].getParameter<int32_t>("vtxNumber") : 0 );
00050     vertexCut_.push_back( useVertices_ ? trkSelectors[i].getParameter<std::string>("vertexCut") : 0);
00051     //  parameters for adapted optimal cuts on chi2 and primary vertex compatibility
00052     res_par_.push_back(trkSelectors[i].getParameter< std::vector<double> >("res_par") );
00053     chi2n_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_par") );
00054     chi2n_no1Dmod_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_no1Dmod_par") );
00055     d0_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par1"));
00056     dz_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par1"));
00057     d0_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par2"));
00058     dz_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par2"));
00059     // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
00060     applyAdaptedPVCuts_.push_back(trkSelectors[i].getParameter<bool>("applyAdaptedPVCuts"));
00061     // Impact parameter absolute cuts.
00062     max_d0_.push_back(trkSelectors[i].getParameter<double>("max_d0"));
00063     max_z0_.push_back(trkSelectors[i].getParameter<double>("max_z0"));
00064     nSigmaZ_.push_back(trkSelectors[i].getParameter<double>("nSigmaZ"));
00065     // Cuts on numbers of layers with hits/3D hits/lost hits.
00066     min_layers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumberLayers") );
00067     min_3Dlayers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumber3DLayers") );
00068     max_lostLayers_.push_back(trkSelectors[i].getParameter<uint32_t>("maxNumberLostLayers"));
00069     min_hits_bypass_.push_back(trkSelectors[i].getParameter<uint32_t>("minHitsToBypassChecks"));
00070     // Flag to apply absolute cuts if no PV passes the selection
00071     applyAbsCutsIfNoPV_.push_back(trkSelectors[i].getParameter<bool>("applyAbsCutsIfNoPV"));
00072     keepAllTracks_.push_back( trkSelectors[i].getParameter<bool>("keepAllTracks")); 
00073     max_relpterr_.push_back(trkSelectors[i].getParameter<double>("max_relpterr"));
00074     min_nhits_.push_back(trkSelectors[i].getParameter<uint32_t>("min_nhits"));
00075 
00076   
00077     setQualityBit_.push_back( false );
00078     std::string qualityStr = trkSelectors[i].getParameter<std::string>("qualityBit");
00079     if (qualityStr != "") {
00080       setQualityBit_[i] = true;
00081       qualityToSet_[i]  = TrackBase::qualityByName(trkSelectors[i].getParameter<std::string>("qualityBit"));
00082     }
00083   
00084     if (setQualityBit_[i] && (qualityToSet_[i] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
00085     "You can't set the quality bit " << trkSelectors[i].getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
00086 
00087     if (applyAbsCutsIfNoPV_[i]) {
00088       max_d0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_d0NoPV"));
00089       max_z0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_z0NoPV"));
00090     }
00091     else{//dummy values
00092       max_d0NoPV_.push_back(0.);
00093       max_z0NoPV_.push_back(0.);
00094     }
00095   
00096     name_.push_back( trkSelectors[i].getParameter<std::string>("name") );
00097 
00098     preFilter_[i]=trkSelectors.size(); // no prefilter
00099 
00100     std::string pfName=trkSelectors[i].getParameter<std::string>("preFilterName");
00101     if (pfName!="") {
00102       bool foundPF=false;
00103       for ( unsigned int j=0; j<i; j++) 
00104         if (name_[j]==pfName ) {
00105           foundPF=true;
00106           preFilter_[i]=j;
00107         }
00108       if ( !foundPF)
00109         throw cms::Exception("Configuration") << "Invalid prefilter name in MultiTrackSelector " 
00110                                               << trkSelectors[i].getParameter<std::string>("preFilterName");
00111           
00112     }
00113 
00114     //    produces<std::vector<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
00115     produces<edm::ValueMap<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
00116   }
00117 }
00118 
00119 MultiTrackSelector::~MultiTrackSelector() {
00120 }
00121 
00122 void MultiTrackSelector::produce( edm::Event& evt, const edm::EventSetup& es ) 
00123 {
00124   using namespace std; 
00125   using namespace edm;
00126   using namespace reco;
00127 
00128   // Get tracks 
00129   Handle<TrackCollection> hSrcTrack;
00130   evt.getByLabel( src_, hSrcTrack );
00131 
00132   // looking for the beam spot
00133   edm::Handle<reco::BeamSpot> hBsp;
00134   evt.getByLabel(beamspot_, hBsp);
00135   reco::BeamSpot vertexBeamSpot;
00136   vertexBeamSpot = *hBsp;
00137         
00138   // Select good primary vertices for use in subsequent track selection
00139   edm::Handle<reco::VertexCollection> hVtx;
00140   if (useVertices_) evt.getByLabel(vertices_, hVtx);
00141 
00142   unsigned int trkSize=hSrcTrack->size();
00143   std::vector<int> selTracksSave( qualityToSet_.size()*trkSize,0);
00144   
00145   //loop over all of the selectors...
00146   for (unsigned int i=0; i<qualityToSet_.size(); i++) {  
00147     std::vector<int> selTracks(trkSize,0);
00148     auto_ptr<edm::ValueMap<int> > selTracksValueMap = auto_ptr<edm::ValueMap<int> >(new edm::ValueMap<int>);
00149     edm::ValueMap<int>::Filler filler(*selTracksValueMap);
00150 
00151     std::vector<Point> points;
00152     std::vector<double> vterr;
00153     std::vector<double> vzerr;
00154     if (useVertices_) selectVertices(i,*hVtx, points, vterr, vzerr);
00155 
00156     // Loop over tracks
00157     size_t current = 0;
00158     for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
00159       const Track & trk = * it;
00160       // Check if this track passes cuts
00161 
00162       LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
00163 
00164       //already removed
00165       bool ok=true;
00166       if (preFilter_[i]<i && selTracksSave[preFilter_[i]*trkSize+current] < 0) {
00167         selTracks.at(current)=-1;
00168         ok=false;
00169         if ( !keepAllTracks_[i]) 
00170           continue;
00171       }
00172       else {
00173         ok = select(i,vertexBeamSpot, trk, points, vterr, vzerr);
00174         if (!ok) { 
00175           LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
00176           if (!keepAllTracks_[i]) { 
00177             selTracks.at(current)=-1;
00178             continue;
00179           }
00180         }
00181         else
00182           LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
00183       }
00184 
00185       if (preFilter_[i]<i ) {
00186         selTracks.at(current)=selTracksSave[preFilter_[i]*trkSize+current];
00187       }
00188       else {
00189         selTracks.at(current)=trk.qualityMask();
00190       }
00191       if ( ok && setQualityBit_[i]) {
00192         selTracks.at(current)= (selTracks.at(current) | (1<<qualityToSet_[i]));
00193         if (!points.empty()) {
00194           if (qualityToSet_[i]==TrackBase::loose) {
00195             selTracks.at(current)=(selTracks.at(current) | (1<<TrackBase::looseSetWithPV));
00196           }
00197           else if (qualityToSet_[i]==TrackBase::highPurity) {
00198             selTracks.at(current)=(selTracks.at(current) | (1<<TrackBase::highPuritySetWithPV));
00199           }
00200         }
00201       }
00202     }
00203     for ( unsigned int j=0; j< trkSize; j++ ) selTracksSave[j+i*trkSize]=selTracks.at(j);
00204     filler.insert(hSrcTrack, selTracks.begin(),selTracks.end());
00205     filler.fill();
00206 
00207     //    evt.put(selTracks,name_[i]);
00208     evt.put(selTracksValueMap,name_[i]);
00209   }
00210 }
00211 
00212 
00213  bool MultiTrackSelector::select(unsigned int tsNum, 
00214                                  const reco::BeamSpot &vertexBeamSpot, 
00215                                  const reco::Track &tk, 
00216                                  const std::vector<Point> &points,
00217                                  std::vector<double> &vterr,
00218                                  std::vector<double> &vzerr) {
00219   // Decide if the given track passes selection cuts.
00220 
00221   using namespace std; 
00222   
00223   if(tk.found()>=min_hits_bypass_[tsNum]) return true;
00224   if ( tk.ndof() < 1E-5 ) return false;
00225 
00226   // Cuts on numbers of layers with hits/3D hits/lost hits.
00227   uint32_t nlayers     = tk.hitPattern().trackerLayersWithMeasurement();
00228   uint32_t nlayers3D   = tk.hitPattern().pixelLayersWithMeasurement() +
00229     tk.hitPattern().numberOfValidStripLayersWithMonoAndStereo();
00230   uint32_t nlayersLost = tk.hitPattern().trackerLayersWithoutMeasurement();
00231   LogDebug("TrackSelection") << "cuts on nlayers: " << nlayers << " " << nlayers3D << " " << nlayersLost << " vs " 
00232                              << min_layers_[tsNum] << " " << min_3Dlayers_[tsNum] << " " << max_lostLayers_[tsNum];
00233   if (nlayers < min_layers_[tsNum]) return false;
00234   if (nlayers3D < min_3Dlayers_[tsNum]) return false;
00235   if (nlayersLost > max_lostLayers_[tsNum]) return false;
00236   LogTrace("TrackSelection") << "cuts on nlayers passed";
00237 
00238   double chi2n =  tk.normalizedChi2();
00239   double chi2n_no1Dmod = chi2n;
00240 
00241   int count1dhits = 0;
00242   for (trackingRecHit_iterator ith = tk.recHitsBegin(), edh = tk.recHitsEnd(); ith != edh; ++ith) {
00243     const TrackingRecHit * hit = ith->get();
00244     if (hit->isValid()) {
00245       if (typeid(*hit) == typeid(SiStripRecHit1D)) ++count1dhits;
00246     }
00247   }
00248   if (count1dhits > 0) {
00249     double chi2 = tk.chi2();
00250     double ndof = tk.ndof();
00251     chi2n = (chi2+count1dhits)/double(ndof+count1dhits);
00252   }
00253   // For each 1D rechit, the chi^2 and ndof is increased by one.  This is a way of retaining approximately
00254   // the same normalized chi^2 distribution as with 2D rechits.
00255   if (chi2n > chi2n_par_[tsNum]*nlayers) return false;
00256 
00257   if (chi2n_no1Dmod > chi2n_no1Dmod_par_[tsNum]*nlayers) return false;
00258 
00259   // Get track parameters
00260   double pt = tk.pt(), eta = tk.eta();
00261 
00262   //cuts on relative error on pt and number of valid hits
00263   double relpterr = tk.ptError()/max(pt,1e-9);
00264   uint32_t nhits = tk.numberOfValidHits();
00265   if(relpterr > max_relpterr_[tsNum]) return false;
00266   if(nhits < min_nhits_[tsNum]) return false;
00267 
00268 
00269   //other track parameters
00270   double d0 = -tk.dxy(vertexBeamSpot.position()), d0E =  tk.d0Error(),
00271     dz = tk.dz(vertexBeamSpot.position()), dzE =  tk.dzError();
00272 
00273   // parametrized d0 resolution for the track pt
00274   double nomd0E = sqrt(res_par_[tsNum][0]*res_par_[tsNum][0]+(res_par_[tsNum][1]/max(pt,1e-9))*(res_par_[tsNum][1]/max(pt,1e-9)));
00275   // parametrized z0 resolution for the track pt and eta
00276   double nomdzE = nomd0E*(std::cosh(eta));
00277 
00278   double dzCut = min( pow(dz_par1_[tsNum][0]*nlayers,dz_par1_[tsNum][1])*nomdzE, 
00279                       pow(dz_par2_[tsNum][0]*nlayers,dz_par2_[tsNum][1])*dzE );
00280   double d0Cut = min( pow(d0_par1_[tsNum][0]*nlayers,d0_par1_[tsNum][1])*nomd0E, 
00281                       pow(d0_par2_[tsNum][0]*nlayers,d0_par2_[tsNum][1])*d0E );
00282 
00283 
00284   // ---- PrimaryVertex compatibility cut
00285   bool primaryVertexZCompatibility(false);   
00286   bool primaryVertexD0Compatibility(false);   
00287 
00288   if (points.empty()) { //If not primaryVertices are reconstructed, check just the compatibility with the BS
00289     //z0 within (n sigma + dzCut) of the beam spot z, if no good vertex is found
00290     if ( abs(dz) < hypot(vertexBeamSpot.sigmaZ()*nSigmaZ_[tsNum],dzCut) ) primaryVertexZCompatibility = true;  
00291     // d0 compatibility with beam line
00292     if (abs(d0) < d0Cut) primaryVertexD0Compatibility = true;     
00293   }
00294 
00295   int iv=0;
00296   for (std::vector<Point>::const_iterator point = points.begin(), end = points.end(); point != end; ++point) {
00297     LogTrace("TrackSelection") << "Test track w.r.t. vertex with z position " << point->z();
00298     if(primaryVertexZCompatibility && primaryVertexD0Compatibility) break;
00299     double dzPV = tk.dz(*point); //re-evaluate the dz with respect to the vertex position
00300     double d0PV = tk.dxy(*point); //re-evaluate the dxy with respect to the vertex position
00301     if(useVtxError_){
00302        double dzErrPV = sqrt(dzE*dzE+vzerr[iv]*vzerr[iv]); // include vertex error in z
00303        double d0ErrPV = sqrt(d0E*d0E+vterr[iv]*vterr[iv]); // include vertex error in xy
00304        iv++;
00305        if (abs(dzPV) < dz_par1_[tsNum][0]*pow(nlayers,dz_par1_[tsNum][1])*nomdzE &&
00306            abs(dzPV) < dz_par2_[tsNum][0]*pow(nlayers,dz_par2_[tsNum][1])*dzErrPV &&
00307            abs(dzPV) < max_z0_[tsNum])  primaryVertexZCompatibility = true;
00308        if (abs(d0PV) < d0_par1_[tsNum][0]*pow(nlayers,d0_par1_[tsNum][1])*nomd0E &&
00309            abs(d0PV) < d0_par2_[tsNum][0]*pow(nlayers,d0_par2_[tsNum][1])*d0ErrPV &&
00310            abs(d0PV) < max_d0_[tsNum]) primaryVertexD0Compatibility = true; 
00311     }else{
00312        if (abs(dzPV) < dzCut)  primaryVertexZCompatibility = true;
00313        if (abs(d0PV) < d0Cut) primaryVertexD0Compatibility = true;     
00314     }
00315     LogTrace("TrackSelection") << "distances " << dzPV << " " << d0PV << " vs " << dzCut << " " << d0Cut;
00316   }
00317 
00318   if (points.empty() && applyAbsCutsIfNoPV_[tsNum]) {
00319     if ( abs(dz) > max_z0NoPV_[tsNum] || abs(d0) > max_d0NoPV_[tsNum]) return false;
00320   }  else {
00321     // Absolute cuts on all tracks impact parameters with respect to beam-spot.
00322     // If BS is not compatible, verify if at least the reco-vertex is compatible (useful for incorrect BS settings)
00323     if (abs(d0) > max_d0_[tsNum] && !primaryVertexD0Compatibility) return false;
00324     LogTrace("TrackSelection") << "absolute cuts on d0 passed";
00325     if (abs(dz) > max_z0_[tsNum] && !primaryVertexZCompatibility) return false;
00326     LogTrace("TrackSelection") << "absolute cuts on dz passed";
00327   }
00328 
00329   LogTrace("TrackSelection") << "cuts on PV: apply adapted PV cuts? " << applyAdaptedPVCuts_[tsNum] 
00330                              << " d0 compatibility? " << primaryVertexD0Compatibility  
00331                              << " z compatibility? " << primaryVertexZCompatibility ;
00332 
00333   if (applyAdaptedPVCuts_[tsNum]) {
00334     return (primaryVertexD0Compatibility && primaryVertexZCompatibility);
00335   } else {
00336     return true;     
00337   }
00338 
00339 }
00340 
00341  void MultiTrackSelector::selectVertices(unsigned int tsNum, 
00342                                          const reco::VertexCollection &vtxs, 
00343                                          std::vector<Point> &points,
00344                                          std::vector<double> &vterr, 
00345                                          std::vector<double> &vzerr) {
00346   // Select good primary vertices
00347   using namespace reco;
00348   int32_t toTake = vtxNumber_[tsNum]; 
00349   for (VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
00350 
00351     LogDebug("SelectVertex") << " select vertex with z position " << it->z() << " " 
00352                              << it->chi2() << " " << it->ndof() << " " << TMath::Prob(it->chi2(), static_cast<int32_t>(it->ndof()));
00353     Vertex vtx = *it;
00354     bool pass = vertexCut_[tsNum]( vtx );
00355     if( pass ) { 
00356       points.push_back(it->position()); 
00357       vterr.push_back(sqrt(it->yError()*it->xError()));
00358       vzerr.push_back(it->zError());
00359       LogTrace("SelectVertex") << " SELECTED vertex with z position " << it->z();
00360       toTake--; if (toTake == 0) break;
00361     }
00362   }
00363 }