CMS 3D CMS Logo

MultiVertexBSeeder.cc

Go to the documentation of this file.
00001 #include "RecoVertex/MultiVertexFit/interface/MultiVertexBSeeder.h"
00002 // #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
00003 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00004 // #include "RecoVertex/VertexTools/interface/BeamSpot.h"
00005 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00006 #include "RecoVertex/VertexTools/interface/FsmwModeFinder3d.h"
00007 #include "CommonTools/Clustering1D/interface/FsmwClusterizer1D.h"
00008 // #include "CommonTools/Clustering1D/interface/MultiClusterizer1D.h"
00009 #include "CommonTools/Clustering1D/interface/OutermostClusterizer1D.h"
00010 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 // #define MVBS_DEBUG
00014 #ifdef MVBS_DEBUG
00015 #include <map>
00016 #include "RecoVertex/MultiVertexFit/interface/DebuggingHarvester.h"
00017 #endif
00018 
00019 using namespace std;
00020 
00021 namespace {
00022   vector < reco::TransientTrack > convert ( 
00023       const vector < const reco::TransientTrack * > & ptrs )
00024   {
00025     vector < reco::TransientTrack > ret;
00026     for ( vector< const reco::TransientTrack * >::const_iterator i=ptrs.begin(); 
00027           i!=ptrs.end() ; ++i )
00028     {
00029       ret.push_back ( **i );
00030     }
00031     return ret;
00032   }
00033 
00034   int verbose()
00035   {
00036     return 0;
00037   }
00038 
00039   inline GlobalPoint & operator += ( GlobalPoint & a, const GlobalPoint & b )
00040   {
00041     a = GlobalPoint ( a.x() + b.x(), a.y() + b.y(), a.z() + b.z() );
00042     return a;
00043   }
00044 
00045   inline GlobalPoint & operator /= ( GlobalPoint & a, float b )
00046   {
00047     a = GlobalPoint ( a.x() / b, a.y() / b, a.z() / b );
00048     return a;
00049   }
00050 
00051   bool element ( const reco::TransientTrack & rt, const TransientVertex & rv )
00052   {
00053     const vector < reco::TransientTrack > trks = rv.originalTracks();
00054     for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); i!=trks.end() ; ++i )
00055     {
00056       if ( rt == ( *i ) ) return true;
00057     };
00058     return false;
00059   }
00060 
00061   GlobalPoint toPoint ( const GlobalVector & v )
00062   {
00063     return GlobalPoint ( v.x(), v.y(), v.z() );
00064   }
00065 
00066   GlobalVector toVector ( const GlobalPoint & v )
00067   {
00068     return GlobalVector ( v.x(), v.y(), v.z() );
00069   }
00070 
00071 
00072   GlobalPoint computeJetOrigin ( const vector < reco::TransientTrack > & trks )
00073   {
00074     FsmwModeFinder3d f;
00075     vector< ModeFinder3d::PointAndDistance> input;
00076     for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
00077           i!=trks.end() ; ++i )
00078     {
00079       input.push_back ( ModeFinder3d::PointAndDistance 
00080                         ( i->impactPointState().globalPosition(), 1. )  );
00081     }
00082     return f(input);
00083   }
00084 
00085   GlobalVector computeJetDirection ( const vector < reco::TransientTrack > & trks )
00086   {
00087     FsmwModeFinder3d f;
00088     vector< ModeFinder3d::PointAndDistance> input;
00089     for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
00090           i!=trks.end() ; ++i )
00091     {
00092       input.push_back ( ModeFinder3d::PointAndDistance ( 
00093             toPoint ( i->impactPointState().globalMomentum() ), 1. )  );
00094     }
00095     GlobalPoint pt ( f(input) );
00096     pt/=pt.mag();
00097     return toVector ( pt );
00098   }
00099 
00100   GlobalTrajectoryParameters computeJetTrajectory ( const vector < reco::TransientTrack > & trks )
00101   {
00108     if ( trks.size() == 0 ) return GlobalTrajectoryParameters();
00109 
00110     GlobalVector mom = computeJetDirection ( trks );
00111     GlobalPoint pos = computeJetOrigin ( trks );
00112 
00113     GlobalTrajectoryParameters ret ( pos, mom, 0,
00114         &(trks[0].impactPointState().globalParameters().magneticField()) );
00115     #ifdef MVBS_DEBUG
00116     DebuggingHarvester("out.txt").save ( ret , "p<sub>tot</sub>" );
00117     #endif
00118     return ret;
00119   }
00120 
00121   vector < Cluster1D < reco::TransientTrack > > computeIPs (
00122       const vector < reco::TransientTrack > & trks )
00123   {
00124     GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
00125     FreeTrajectoryState axis ( jet );
00126     TwoTrackMinimumDistance ttmd;
00127     vector < Cluster1D < reco::TransientTrack > > pts;
00128     for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); 
00129           i!=trks.end() ; ++i )
00130     {
00131       bool status = ttmd.calculate( axis,*( i->impactPointState().freeState() ) );
00132       pair < GlobalPoint, GlobalPoint > pt = ttmd.points();
00133       double d = ( pt.first - pt.second ).mag();
00134       double w = 1. / ( 0.002 + d ); // hard coded weights
00135       double s = ( pt.first - axis.position() ).mag();
00136       Measurement1D ms ( s, 1.0 );
00137       vector < const reco::TransientTrack * > trk;
00138       trk.push_back ( &(*i) );
00139       pts.push_back ( Cluster1D < reco::TransientTrack > ( ms, trk, w ) );
00140     }
00141     /*
00142     #ifdef MVBS_DEBUG
00143     map < string, harvest::MultiType > attrs;
00144     attrs["point:mag"]=0.5;
00145     attrs["point:color"]="khaki";
00146     DebuggingHarvester("out.txt").save ( pts, jet, attrs, "ip" );
00147     #endif
00148     */
00149     return pts;
00150   }
00151 
00152   GlobalPoint computePos ( const GlobalTrajectoryParameters & jet,
00153       double s )
00154   {
00155     GlobalPoint ret = jet.position();
00156     ret += s * jet.momentum();
00157     return ret;
00158   }
00159 
00160   TransientVertex pseudoVertexFit ( const Cluster1D < reco::TransientTrack > & src,
00161       bool ascending=false, bool kalmanfit=false )
00162   {
00163     // cout << "[MultiVertexBSeeder] debug: pseudoVertexFit with " << flush;
00164     vector < const reco::TransientTrack * > trkptrs=src.tracks();
00165     vector < reco::TransientTrack > trks = convert ( trkptrs );
00166     // cout << trks.size() << " tracks.";
00167     GlobalPoint gp;
00168     GlobalError ge;
00169     if ( kalmanfit )
00170     {
00171         TransientVertex v = KalmanVertexFitter().vertex ( trks );
00172         gp=v.position();
00173     }
00174     TransientVertex ret = TransientVertex ( gp, ge, trks, -1. );
00175     TransientVertex::TransientTrackToFloatMap mp;
00176     float w=1.0; float r=0.5;
00177     if ( ascending ) { w = pow ( (float) 0.5, (int) (trks.size()-1) ); r=2.0; };
00178     for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); 
00179          i!=trks.end() ; ++i )
00180     {
00181       mp[*i]=w;
00182       w*=r;
00183     }
00184     ret.weightMap ( mp );
00185     // cout << "[MultiVertexBSeeder] debug: return pseudoVertexFit with " << endl;
00186     return ret;
00187   }
00188 
00189   /* TransientVertex kalmanVertexFit ( const Cluster1D < reco::TransientTrack > & src )
00190   {
00191     KalmanVertexFitter fitter;
00192       vector < const reco::TransientTrack * > trkptrs=src.tracks();
00193       vector < reco::TransientTrack > trks = convert ( trkptrs );
00194       return fitter.vertex ( trks );
00195     return TransientVertex();
00196   }*/
00197 }
00198 
00199 MultiVertexBSeeder * MultiVertexBSeeder::clone() const
00200 {
00201   return new MultiVertexBSeeder ( * this );
00202 }
00203 
00204 MultiVertexBSeeder::MultiVertexBSeeder ( double nsigma ) :
00205   theNSigma ( nsigma ) {}
00206 
00207 vector < TransientVertex > MultiVertexBSeeder::vertices (
00208     const vector < reco::TransientTrack > & trks,
00209     const reco::BeamSpot & s) const
00210 {
00211   return vertices ( trks );
00212 }
00213 
00214 vector < TransientVertex > MultiVertexBSeeder::vertices (
00215     const vector < reco::TransientTrack > & trks ) const
00216 {
00217   vector < Cluster1D < reco::TransientTrack > > ips = computeIPs ( trks );
00218   /*
00219   FsmwClusterizer1D < reco::TransientTrack > fsmw ( .4, theNSigma );
00220   MultiClusterizer1D<reco::TransientTrack> finder ( fsmw );
00221   */
00222   OutermostClusterizer1D<reco::TransientTrack> finder;
00223 
00224   pair < vector < Cluster1D<reco::TransientTrack> >, 
00225          vector < const reco::TransientTrack * > > res;
00226   res = finder ( ips );
00227   #ifdef MVBS_DEBUG
00228   // need to compute jet trajectory again :-(
00229   GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
00230   map < string, harvest::MultiType > attrs;
00231   attrs["point:mag"]=0.75;
00232   attrs["point:color"]="blue";
00233   attrs["point:name"]="mode";
00234   DebuggingHarvester("out.txt").save ( res.first, jet, attrs, "mode" );
00235   #endif
00236 
00237   vector < TransientVertex > ret;
00238   for ( vector< Cluster1D < reco::TransientTrack > >::const_iterator i=res.first.begin(); 
00239         i!=res.first.end() ; ++i )
00240   {
00241     ret.push_back ( pseudoVertexFit ( *i, i!=res.first.begin(), 
00242                     /* kalman fit*/ true ) );
00243   }
00244   return ret;
00245 }
00246 
00247 #ifdef MVBS_DEBUG
00248 #undef MVBS_DEBUG
00249 #endif

Generated on Tue Jun 9 17:46:11 2009 for CMSSW by  doxygen 1.5.4