CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoVertex/TrimmedKalmanVertexFinder/src/TrimmedVertexFinder.cc

Go to the documentation of this file.
00001 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/TrimmedVertexFinder.h"
00002 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00003 #include "RecoVertex/VertexTools/interface/PerigeeLinearizedTrackState.h"
00004 #include "FWCore/Utilities/interface/Exception.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 using namespace reco;
00008 
00009 TrimmedVertexFinder::TrimmedVertexFinder(
00010   const VertexFitter<5> * vf, const VertexUpdator<5> * vu, 
00011   const VertexTrackCompatibilityEstimator<5> * ve)
00012   : theFitter(vf->clone()), theUpdator(vu->clone()), 
00013     theEstimator(ve->clone()), theMinProb(0.05)
00014 {}
00015 
00016     
00017 TrimmedVertexFinder::TrimmedVertexFinder(
00018 const TrimmedVertexFinder & other) 
00019   : theFitter(other.theFitter->clone()), 
00020   theUpdator(other.theUpdator->clone()), 
00021   theEstimator(other.theEstimator->clone()), 
00022   theMinProb(other.theMinProb)
00023 {}
00024 
00025 
00026 TrimmedVertexFinder::~TrimmedVertexFinder() {
00027   delete theFitter;
00028   delete theUpdator;
00029   delete theEstimator;
00030 }
00031 
00032 std::vector<TransientVertex> 
00033 TrimmedVertexFinder::vertices(std::vector<TransientTrack> & tks ) const
00034 {
00035   // FIXME write this!!!
00036   return vertices ( tks, reco::BeamSpot(), false );
00037 }
00038 
00039 std::vector<TransientVertex> 
00040 TrimmedVertexFinder::vertices(std::vector<TransientTrack> & tks,
00041     const reco::BeamSpot & spot, bool use_spot )  const
00042 {
00043   std::vector<TransientVertex> all;
00044   if (tks.size() < 2) return all;
00045 
00046   // prepare vertex tracks and initial vertex
00047   CachingVertex<5> vtx;
00048   if ( use_spot ) 
00049   {
00050      vtx = theFitter->vertex(tks, spot );
00051   } else {
00052      vtx = theFitter->vertex(tks);
00053   }
00054   if (!vtx.isValid()) {
00055     edm::LogWarning ( "TrimmedVertexFinder" ) << "initial vertex invalid."
00056             << " vertex finding stops here.";
00057     return all;
00058   }
00059 
00060   std::vector<RefCountedVertexTrack> selected = vtx.tracks();
00061 
00062   // reject incompatible tracks starting from the worst
00063   std::vector<RefCountedVertexTrack> remain;
00064   bool found = false;
00065   while (!found && selected.size() >= 2) {
00066 
00067     // find track with worst compatibility
00068     std::vector<RefCountedVertexTrack>::iterator iWorst = theWorst(vtx, 
00069                                                               selected, 
00070                                                               theMinProb);
00071     
00072     if (iWorst != selected.end()) {
00073       // reject track
00074       remain.push_back(*iWorst);
00075       selected.erase(iWorst);
00076       
00077       if (selected.size() == 1) {
00078         // not enough tracks to build new vertices
00079         remain.push_back(selected.front());
00080       }
00081       else {
00082         // removing track from vertex
00083         // need to redo the vertex fit instead of removing the track;
00084         // successive removals lead to numerical problems
00085         // this is however quick since intermediate quantities are cached
00086         // in the RefCountedVertexTracks
00087   if ( use_spot ) // && all.size()==0 )
00088   {
00089           vtx = theFitter->vertex(selected,spot);
00090   } else {
00091           vtx = theFitter->vertex(selected);
00092   }
00093         if (!vtx.isValid()) {
00094     edm::LogWarning ( "TrimmedVertexFinder" ) << "current vertex invalid"
00095                << "vertex finding stops here.";
00096           return all;
00097         }
00098 
00099         // ref-counted tracks may have changed during the fit
00100         // if the linearization point has moved too much
00101         selected = vtx.tracks();
00102       }
00103       
00104     } else {
00105       // no incompatible track remains, make vertex
00106       found = true;
00107       int n_tracks_in_vertex = selected.size();
00108 
00109       // now return all tracks with weight < 0.5 to 'remain'.
00110       for ( std::vector< RefCountedVertexTrack >::const_iterator t=selected.begin();
00111           t!=selected.end() ; ++t )
00112       {
00113         if ( (**t).weight() < 0.5 )
00114         {
00115           /*
00116           cout << "[TrimmedVertexFinder] recycling track with weight "
00117                << (**t).weight() << endl;*/
00118           remain.push_back ( *t );
00119           n_tracks_in_vertex--; // one 'good' track less in the vertex
00120         };
00121       };
00122 
00123       if ( n_tracks_in_vertex > 1 ) {
00124         all.push_back(vtx);
00125       }
00126       else {
00127         edm::LogError ( "TrimmedVertexFinder" ) 
00128           << "found vertex has less than 2 tracks";
00129       }
00130     }
00131   }
00132 
00133   // modify list of incompatible tracks
00134   tks.clear();
00135   for (std::vector<RefCountedVertexTrack>::const_iterator i = remain.begin(); 
00136        i != remain.end(); i++) {
00137     const PerigeeLinearizedTrackState* plts = 
00138       dynamic_cast<const PerigeeLinearizedTrackState*>
00139       ((**i).linearizedTrack().get());
00140     if (plts == 0) {
00141       throw cms::Exception("TrimmedVertexFinder: can't take track from non-perigee track state");
00142     }
00143 
00144     tks.push_back(plts->track());
00145   }
00146 
00147   // return 0 or 1 vertex
00148   return all;
00149 }
00150 
00151 
00152 std::vector<TrimmedVertexFinder::RefCountedVertexTrack>::iterator 
00153 TrimmedVertexFinder::theWorst(const CachingVertex<5> & vtx, 
00154   std::vector<RefCountedVertexTrack> & vtxTracks, float cut) const
00155 {
00156 
00157   //  cout << "Cut is now " << cut << endl;
00158 
00159   // find track with worst compatibility
00160   std::vector<RefCountedVertexTrack>::iterator iWorst = vtxTracks.end();
00161   float worseChi2 = 0.;
00162   for (std::vector<RefCountedVertexTrack>::iterator itr = vtxTracks.begin();
00163        itr != vtxTracks.end(); itr++) {
00164 
00165     CachingVertex<5> newV = theUpdator->remove(vtx, *itr);
00166     if (!newV.isValid()) return itr;
00167     std::pair<bool, double> result = theEstimator->estimate(newV, *itr);
00168     if (!result.first) return itr;
00169     float chi2 = result.second;
00170 
00171     // compute number of degrees of freedom
00172     if ( chi2 > 0. ) {
00173       // we want to keep negative chi squares, and avoid calling
00174       // ChiSquaredProbability with a negative chi2.
00175       int ndf = 2;
00176       if (vtx.tracks().size() == 2) ndf = 1;
00177 
00178       float prob = ChiSquaredProbability(chi2, ndf);
00179       if (prob < cut && chi2 >= worseChi2) {
00180         // small inconsistency: sorting on chi2, not on chi2 probability
00181         // because large chi2 are not rounded off while small probabilities 
00182         // are rounded off to 0....
00183         // should not matter since at a given iteration 
00184         // all track chi2 increments have the same ndf
00185         iWorst = itr;
00186         worseChi2 = chi2;
00187       }
00188     }
00189   }
00190 
00191   return iWorst;
00192 }