CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoVertex/AdaptiveVertexFinder/src/AdaptiveVertexReconstructor.cc

Go to the documentation of this file.
00001 #include "RecoVertex/AdaptiveVertexFinder/interface/AdaptiveVertexReconstructor.h"
00002 #include "RecoVertex/VertexTools/interface/GeometricAnnealing.h"
00003 #include "FWCore/Utilities/interface/EDMException.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "RecoVertex/VertexTools/interface/DummyVertexSmoother.h"
00006 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexSmoother.h"
00007 #include <algorithm>
00008 
00009 using namespace std;
00010 
00011 TransientVertex AdaptiveVertexReconstructor::cleanUp ( const TransientVertex & old ) const
00012 {
00013   vector < reco::TransientTrack > trks = old.originalTracks();
00014   vector < reco::TransientTrack > newtrks;
00015   TransientVertex::TransientTrackToFloatMap mp;
00016   static const float minweight = 1.e-8; // discard all tracks with lower weight
00017   for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
00018         i!=trks.end() ; ++i )
00019   {
00020     if ( old.trackWeight ( *i ) > minweight )
00021     {
00022       newtrks.push_back ( *i );
00023       mp[*i]=old.trackWeight ( *i );
00024     }
00025   }
00026 
00027   TransientVertex ret;
00028 
00029   if ( old.hasPrior() )
00030   {
00031     VertexState priorstate ( old.priorPosition(), old.priorError() );
00032     ret=TransientVertex ( priorstate, old.vertexState(), newtrks,
00033         old.totalChiSquared(), old.degreesOfFreedom() );
00034   } else {
00035     ret=TransientVertex ( old.vertexState(), newtrks,
00036                           old.totalChiSquared(), old.degreesOfFreedom() );
00037   }
00038   ret.weightMap ( mp ); // set weight map
00039 
00040   if ( old.hasRefittedTracks() )
00041   {
00042     // we have refitted tracks -- copy them!
00043     vector < reco::TransientTrack > newrfs;
00044     vector < reco::TransientTrack > oldrfs=old.refittedTracks();
00045     vector < reco::TransientTrack >::const_iterator origtrkiter=trks.begin();
00046     for ( vector< reco::TransientTrack >::const_iterator i=oldrfs.begin(); i!=oldrfs.end() ; ++i )
00047     {
00048       if ( old.trackWeight ( *origtrkiter ) > minweight )
00049       {
00050         newrfs.push_back ( *i );
00051       }
00052       origtrkiter++;
00053     }
00054     if ( newrfs.size() ) ret.refittedTracks ( newrfs ); // copy refitted tracks
00055   }
00056 
00057   if ( ret.refittedTracks().size() > ret.originalTracks().size() )
00058   {
00059     edm::LogError("AdaptiveVertexReconstructor" )
00060       << "More refitted tracks than original tracks!";
00061   }
00062 
00063   return ret;
00064 }
00065 
00066 void AdaptiveVertexReconstructor::erase (
00067     const TransientVertex & newvtx,
00068     set < reco::TransientTrack > & remainingtrks,
00069     float w ) const
00070 {
00071   /*
00072    * Erase tracks that are in newvtx from remainingtrks 
00073    * But erase only if trackweight > w
00074    */
00075   const vector < reco::TransientTrack > & origtrks = newvtx.originalTracks();
00076   bool erased=false;
00077   for ( vector< reco::TransientTrack >::const_iterator i=origtrks.begin();
00078         i!=origtrks.end(); ++i )
00079   {
00080     double weight = newvtx.trackWeight ( *i );
00081     if ( weight > w )
00082     {
00083       remainingtrks.erase ( *i );
00084       erased=true;
00085     };
00086   };
00087 }
00088 
00089 AdaptiveVertexReconstructor::AdaptiveVertexReconstructor(
00090     float primcut, float seccut, float min_weight, bool smoothing ) :
00091   thePrimaryFitter ( 0 ), theSecondaryFitter ( 0 ),
00092        theMinWeight( min_weight ), theWeightThreshold ( 0.001 )
00093 {
00094   setupFitters ( primcut, 256., 0.25, seccut, 256., 0.25, smoothing );
00095 }
00096 
00097 AdaptiveVertexReconstructor::~AdaptiveVertexReconstructor()
00098 {
00099   if ( thePrimaryFitter ) delete thePrimaryFitter;
00100   if ( theSecondaryFitter ) delete theSecondaryFitter;
00101 }
00102 
00103 void AdaptiveVertexReconstructor::setupFitters ( float primcut, 
00104     float primT, float primr, float seccut, float secT,
00105     float secr, bool smoothing )
00106 {
00107   VertexSmoother<5> * smoother ;
00108   if ( smoothing )
00109   {
00110     smoother = new KalmanVertexSmoother();
00111   } else {
00112     smoother = new DummyVertexSmoother<5>();
00113   }
00114 
00115   if ( thePrimaryFitter ) delete thePrimaryFitter;
00116   if ( theSecondaryFitter ) delete theSecondaryFitter;
00117 
00118   /*
00119   edm::LogError ("AdaptiveVertexReconstructor" )
00120     << "Tini and r are hardcoded now!";
00121     */
00122   thePrimaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( primcut, primT, primr ), DefaultLinearizationPointFinder(),
00123       KalmanVertexUpdator<5>(), KalmanVertexTrackCompatibilityEstimator<5>(), *smoother );
00124   thePrimaryFitter->setWeightThreshold ( theWeightThreshold );
00125   // if the primary fails, sth is wrong, so here we set a threshold on the weight.
00126   theSecondaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( seccut, secT, secr ), DefaultLinearizationPointFinder(),
00127       KalmanVertexUpdator<5>(), 
00128       KalmanVertexTrackCompatibilityEstimator<5>(), *smoother );
00129   theSecondaryFitter->setWeightThreshold ( 0. );
00130   // need to set it or else we have 
00131   // unwanted exceptions to deal with.
00132   // cleanup can come later!
00133   delete smoother;
00134 }
00135 
00136 AdaptiveVertexReconstructor::AdaptiveVertexReconstructor( const edm::ParameterSet & m )
00137   : thePrimaryFitter(0), theSecondaryFitter(0), theMinWeight(0.5), theWeightThreshold ( 0.001 )
00138 {
00139   float primcut = 2.0;
00140   float seccut = 6.0;
00141   bool smoothing=false;
00142   // float primT = 4096.;
00143   // float primr = 0.125;
00144   float primT = 256.;
00145   float primr = 0.25;
00146   float secT = 256.;
00147   float secr = 0.25;
00148   
00149   try {
00150     primcut =  m.getParameter<double>("primcut");
00151     primT =  m.getParameter<double>("primT");
00152     primr =  m.getParameter<double>("primr");
00153     seccut =  m.getParameter<double>("seccut");
00154     secT =  m.getParameter<double>("secT");
00155     secr =  m.getParameter<double>("secr");
00156     theMinWeight = m.getParameter<double>("minweight");
00157     theWeightThreshold = m.getParameter<double>("weightthreshold");
00158     smoothing =  m.getParameter<bool>("smoothing");
00159   } catch ( edm::Exception & e ) {
00160     edm::LogError ("AdaptiveVertexReconstructor") << e.what();
00161   }
00162 
00163   setupFitters ( primcut, primT, primr, seccut, secT, secr, smoothing );
00164 }
00165 
00166 vector<TransientVertex> 
00167     AdaptiveVertexReconstructor::vertices(const vector<reco::TransientTrack> & t, 
00168         const reco::BeamSpot & s ) const
00169 {
00170   return vertices ( vector<reco::TransientTrack>(), t, s, false, true );
00171 }
00172 
00173 vector<TransientVertex> 
00174     AdaptiveVertexReconstructor::vertices(const vector<reco::TransientTrack> & primaries,
00175        const vector<reco::TransientTrack > & tracks, const reco::BeamSpot & s ) const
00176 {
00177   return vertices ( primaries, tracks, s, true, true );
00178 }
00179 
00180 
00181 vector<TransientVertex> AdaptiveVertexReconstructor::vertices (
00182     const vector<reco::TransientTrack> & tracks ) const
00183 {
00184   return vertices ( vector<reco::TransientTrack>(), tracks, reco::BeamSpot(), 
00185                     false, false );
00186 }
00187 
00188 vector<TransientVertex> AdaptiveVertexReconstructor::vertices (
00189     const vector<reco::TransientTrack> & primaries,
00190     const vector<reco::TransientTrack> & tracks,
00191     const reco::BeamSpot & s, bool has_primaries, bool usespot ) const
00192 {
00193   vector < TransientVertex > ret;
00194   set < reco::TransientTrack > remainingtrks;
00195 
00196   copy(tracks.begin(), tracks.end(), 
00197             inserter(remainingtrks, remainingtrks.begin()));
00198 
00199   int ctr=0;
00200   unsigned int n_tracks = remainingtrks.size();
00201 
00202   // cout << "[AdaptiveVertexReconstructor] DEBUG ::vertices!!" << endl;
00203   try {
00204     while ( remainingtrks.size() > 1 )
00205     {
00206       /*
00207       cout << "[AdaptiveVertexReconstructor] next round: "
00208            << remainingtrks.size() << endl;
00209            */
00210       ctr++;
00211       const AdaptiveVertexFitter * fitter = theSecondaryFitter;
00212       if ( ret.size() == 0 )
00213       {
00214         fitter = thePrimaryFitter;
00215       };
00216       vector < reco::TransientTrack > fittrks;
00217       fittrks.reserve ( remainingtrks.size() );
00218 
00219       copy(remainingtrks.begin(), remainingtrks.end(), back_inserter(fittrks));
00220 
00221       TransientVertex tmpvtx;
00222       if ( (ret.size() == 0) && has_primaries )
00223       {
00224         // add the primaries to the fitted tracks.
00225         copy ( primaries.begin(), primaries.end(), back_inserter(fittrks) );
00226       }
00227       if ( (ret.size() == 0) && usespot )
00228       {
00229         tmpvtx=fitter->vertex ( fittrks, s );
00230       } else {
00231         tmpvtx=fitter->vertex ( fittrks );
00232       }
00233       TransientVertex newvtx = cleanUp ( tmpvtx );
00234       ret.push_back ( newvtx );
00235       erase ( newvtx, remainingtrks, theMinWeight );
00236       if ( n_tracks == remainingtrks.size() )
00237       {
00238         if ( usespot )
00239         {
00240           // try once more without beamspot constraint!
00241           usespot=false;
00242           LogDebug("AdaptiveVertexReconstructor") 
00243             << "no tracks in vertex. trying again without beamspot constraint!";
00244           continue;
00245         }
00246         LogDebug("AdaptiveVertexReconstructor") << "all tracks (" << n_tracks
00247              << ") would be recycled for next fit. Trying with low threshold!";
00248         erase ( newvtx, remainingtrks, 1.e-5 );
00249         if ( n_tracks == remainingtrks.size() )
00250         {
00251           LogDebug("AdaptiveVertexReconstructor") << "low threshold didnt help! "
00252                                                   << "Discontinue procedure!";
00253           break;
00254         }
00255       };
00256 
00257       // cout << "[AdaptiveVertexReconstructor] erased" << endl;
00258       n_tracks = remainingtrks.size();
00259     };
00260   } catch ( VertexException & v ) {
00261     // Will catch all (not enough significant tracks exceptions.
00262     // in this case, the iteration can safely terminate.
00263   };
00264 
00265   return cleanUpVertices ( ret );
00266 }
00267 
00268 vector<TransientVertex> AdaptiveVertexReconstructor::cleanUpVertices (
00269     const vector < TransientVertex > & old ) const
00270 {
00271   vector < TransientVertex > ret;
00272   for ( vector< TransientVertex >::const_iterator i=old.begin(); i!=old.end() ; ++i )
00273   {
00274     if (!(i->hasTrackWeight()))
00275     { // if we dont have track weights, we take the vtx
00276       ret.push_back ( *i );
00277       continue;
00278     }
00279 
00280     // maybe this should be replaced with asking for the ndf ...
00281     // e.g. if ( ndf > - 1. )
00282     int nsig=0; // number of significant tracks. 
00283     TransientVertex::TransientTrackToFloatMap wm = i->weightMap();
00284     for ( TransientVertex::TransientTrackToFloatMap::const_iterator w=wm.begin(); w!=wm.end() ; ++w )
00285     {
00286       if (w->second > theWeightThreshold) nsig++;
00287     }
00288     if ( nsig > 1 ) ret.push_back ( *i );
00289   }
00290 
00291   return ret;
00292 }