CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoVertex/TrimmedKalmanVertexFinder/src/ConfigurableTrimmedVertexFinder.cc

Go to the documentation of this file.
00001 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/ConfigurableTrimmedVertexFinder.h"
00002 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00003 
00004 using namespace reco;
00005 
00006 ConfigurableTrimmedVertexFinder::ConfigurableTrimmedVertexFinder(
00007   const VertexFitter<5> * vf, 
00008   const VertexUpdator<5> * vu, 
00009   const VertexTrackCompatibilityEstimator<5> * ve) 
00010   : theClusterFinder(vf, vu, ve), theVtxFitProbCut(0.01), 
00011     theTrackCompatibilityToPV(0.05), theTrackCompatibilityToSV(0.01), 
00012     theMaxNbOfVertices(0)
00013 {
00014   // default pt cut is 1.5 GeV
00015   theFilter.setPtCut(1.5);
00016 }
00017 
00018 void ConfigurableTrimmedVertexFinder::setParameters ( const edm::ParameterSet & s )
00019 {
00020   theFilter.setPtCut(s.getParameter<double>("ptCut"));
00021   theTrackCompatibilityToPV = s.getParameter<double>("trackCompatibilityToPVcut");
00022   theTrackCompatibilityToSV = s.getParameter<double>("trackCompatibilityToSVcut");
00023   theVtxFitProbCut = s.getParameter<double>("vtxFitProbCut"); 
00024   theMaxNbOfVertices =  s.getParameter<int>("maxNbOfVertices");
00025 }
00026 
00027 
00028 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
00029   const std::vector<TransientTrack> & tracks) const
00030 {
00031   std::vector<TransientTrack> remaining;
00032 
00033   return vertices(tracks, remaining, reco::BeamSpot(), false );
00034 
00035 }
00036 
00037 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
00038   const std::vector<TransientTrack> & tracks, const reco::BeamSpot & spot ) const
00039 {
00040   std::vector<TransientTrack> remaining;
00041   return vertices ( tracks, remaining, spot, true );
00042 }
00043 
00044 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertices(
00045   const std::vector<TransientTrack> & tracks, std::vector<TransientTrack> & unused,
00046   const reco::BeamSpot & spot, bool use_spot ) const
00047 {
00048   resetEvent(tracks);
00049   analyseInputTracks(tracks);
00050 
00051   std::vector<TransientTrack> filtered;
00052   for (std::vector<TransientTrack>::const_iterator it = tracks.begin();
00053        it != tracks.end(); it++) {
00054     if (theFilter(*it)) { 
00055       filtered.push_back(*it);
00056     }
00057     else {
00058       unused.push_back(*it);
00059     }
00060   }
00061 
00062   std::vector<TransientVertex> all = vertexCandidates(filtered, unused,
00063       spot, use_spot );
00064 
00065   analyseVertexCandidates(all);
00066 
00067   std::vector<TransientVertex> sel = clean(all);
00068 
00069   analyseFoundVertices(sel);
00070 
00071   return sel;
00072 
00073 }
00074 
00075 
00076 std::vector<TransientVertex> ConfigurableTrimmedVertexFinder::vertexCandidates(
00077   const std::vector<TransientTrack> & tracks, std::vector<TransientTrack> & unused,
00078   const reco::BeamSpot & spot, bool use_spot ) const 
00079 {
00080 
00081   std::vector<TransientVertex> cand;
00082 
00083   std::vector<TransientTrack> remain = tracks;
00084 
00085   while (true) {
00086 
00087     float tkCompCut = (cand.size() == 0 ? 
00088                        theTrackCompatibilityToPV 
00089                        : theTrackCompatibilityToSV);
00090 
00091     //    std::cout << "PVR:compat cut " << tkCompCut << std::endl;
00092     theClusterFinder.setTrackCompatibilityCut(tkCompCut);
00093     //    std::cout << "PVCF:compat cut after setting " 
00094     //   << theClusterFinder.trackCompatibilityCut() << std::endl;
00095 
00096     std::vector<TransientVertex> newVertices;
00097     if ( cand.size() == 0 && use_spot )
00098     {
00099       newVertices = theClusterFinder.vertices(remain, spot );
00100     } else {
00101       newVertices = theClusterFinder.vertices(remain);
00102     }
00103     if (newVertices.empty()) break;
00104 
00105     analyseClusterFinder(newVertices, remain);
00106     
00107     for (std::vector<TransientVertex>::const_iterator iv = newVertices.begin();
00108          iv != newVertices.end(); iv++) {
00109       if ( iv->originalTracks().size() > 1 ) {
00110         cand.push_back(*iv);
00111       } 
00112       else {
00113         // candidate has too few tracks - get them back into the vector
00114         for ( std::vector< TransientTrack >::const_iterator trk
00115                 = iv->originalTracks().begin();
00116               trk != iv->originalTracks().end(); ++trk ) {
00117           unused.push_back ( *trk );
00118         }
00119       }
00120     }
00121 
00122     // when max number of vertices reached, stop
00123     if (theMaxNbOfVertices != 0) {
00124       if (cand.size() >= (unsigned int) theMaxNbOfVertices) break;
00125     }
00126   }
00127 
00128   for (std::vector<TransientTrack>::const_iterator it = remain.begin();
00129        it != remain.end(); it++) {
00130     unused.push_back(*it);
00131   }
00132 
00133   return cand;
00134 }
00135 
00136 
00137 std::vector<TransientVertex> 
00138 ConfigurableTrimmedVertexFinder::clean(const std::vector<TransientVertex> & candidates) const
00139 {
00140   std::vector<TransientVertex> sel;
00141   for (std::vector<TransientVertex>::const_iterator i = candidates.begin(); 
00142        i != candidates.end(); i++) {
00143 
00144     if (ChiSquaredProbability((*i).totalChiSquared(), (*i).degreesOfFreedom())
00145         > theVtxFitProbCut) { sel.push_back(*i); }
00146   }
00147 
00148   return sel;
00149 }