Go to the documentation of this file.00001 #include "RecoVertex/MultiVertexFit/interface/MultiVertexBSeeder.h"
00002
00003 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00004
00005 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00006 #include "RecoVertex/VertexTools/interface/FsmwModeFinder3d.h"
00007 #include "CommonTools/Clustering1D/interface/FsmwClusterizer1D.h"
00008
00009 #include "CommonTools/Clustering1D/interface/OutermostClusterizer1D.h"
00010 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012
00013
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 if (status) {
00133 pair < GlobalPoint, GlobalPoint > pt = ttmd.points();
00134 double d = ( pt.first - pt.second ).mag();
00135 double w = 1. / ( 0.002 + d );
00136 double s = ( pt.first - axis.position() ).mag();
00137 Measurement1D ms ( s, 1.0 );
00138 vector < const reco::TransientTrack * > trk;
00139 trk.push_back ( &(*i) );
00140 pts.push_back ( Cluster1D < reco::TransientTrack > ( ms, trk, w ) );
00141 }
00142 }
00143
00144
00145
00146
00147
00148
00149
00150
00151 return pts;
00152 }
00153
00154 GlobalPoint computePos ( const GlobalTrajectoryParameters & jet,
00155 double s )
00156 {
00157 GlobalPoint ret = jet.position();
00158 ret += s * jet.momentum();
00159 return ret;
00160 }
00161
00162 TransientVertex pseudoVertexFit ( const Cluster1D < reco::TransientTrack > & src,
00163 bool ascending=false, bool kalmanfit=false )
00164 {
00165
00166 vector < const reco::TransientTrack * > trkptrs=src.tracks();
00167 vector < reco::TransientTrack > trks = convert ( trkptrs );
00168
00169 GlobalPoint gp;
00170 GlobalError ge;
00171 if ( kalmanfit )
00172 {
00173 TransientVertex v = KalmanVertexFitter().vertex ( trks );
00174 gp=v.position();
00175 }
00176 TransientVertex ret = TransientVertex ( gp, ge, trks, -1. );
00177 TransientVertex::TransientTrackToFloatMap mp;
00178 float w=1.0; float r=0.5;
00179 if ( ascending ) { w = pow ( (float) 0.5, (int) (trks.size()-1) ); r=2.0; };
00180 for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
00181 i!=trks.end() ; ++i )
00182 {
00183 mp[*i]=w;
00184 w*=r;
00185 }
00186 ret.weightMap ( mp );
00187
00188 return ret;
00189 }
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199 }
00200
00201 MultiVertexBSeeder * MultiVertexBSeeder::clone() const
00202 {
00203 return new MultiVertexBSeeder ( * this );
00204 }
00205
00206 MultiVertexBSeeder::MultiVertexBSeeder ( double nsigma ) :
00207 theNSigma ( nsigma ) {}
00208
00209 vector < TransientVertex > MultiVertexBSeeder::vertices (
00210 const vector < reco::TransientTrack > & trks,
00211 const reco::BeamSpot & s) const
00212 {
00213 return vertices ( trks );
00214 }
00215
00216 vector < TransientVertex > MultiVertexBSeeder::vertices (
00217 const vector < reco::TransientTrack > & trks ) const
00218 {
00219 vector < Cluster1D < reco::TransientTrack > > ips = computeIPs ( trks );
00220
00221
00222
00223
00224 OutermostClusterizer1D<reco::TransientTrack> finder;
00225
00226 pair < vector < Cluster1D<reco::TransientTrack> >,
00227 vector < const reco::TransientTrack * > > res;
00228 res = finder ( ips );
00229 #ifdef MVBS_DEBUG
00230
00231 GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
00232 map < string, harvest::MultiType > attrs;
00233 attrs["point:mag"]=0.75;
00234 attrs["point:color"]="blue";
00235 attrs["point:name"]="mode";
00236 DebuggingHarvester("out.txt").save ( res.first, jet, attrs, "mode" );
00237 #endif
00238
00239 vector < TransientVertex > ret;
00240 for ( vector< Cluster1D < reco::TransientTrack > >::const_iterator i=res.first.begin();
00241 i!=res.first.end() ; ++i )
00242 {
00243 ret.push_back ( pseudoVertexFit ( *i, i!=res.first.begin(),
00244 true ) );
00245 }
00246 return ret;
00247 }
00248
00249 #ifdef MVBS_DEBUG
00250 #undef MVBS_DEBUG
00251 #endif