CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiVertexBSeeder.cc
Go to the documentation of this file.
2 // #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
4 // #include "RecoVertex/VertexTools/interface/BeamSpot.h"
8 // #include "CommonTools/Clustering1D/interface/MultiClusterizer1D.h"
12 
13 // #define MVBS_DEBUG
14 #ifdef MVBS_DEBUG
15 #include <map>
16 #include "RecoVertex/MultiVertexFit/interface/DebuggingHarvester.h"
17 #endif
18 
19 using namespace std;
20 
21 namespace {
22  vector < reco::TransientTrack > convert (
23  const vector < const reco::TransientTrack * > & ptrs )
24  {
25  vector < reco::TransientTrack > ret;
26  for ( vector< const reco::TransientTrack * >::const_iterator i=ptrs.begin();
27  i!=ptrs.end() ; ++i )
28  {
29  ret.push_back ( **i );
30  }
31  return ret;
32  }
33 
34  int verbose()
35  {
36  return 0;
37  }
38 
39  inline GlobalPoint & operator += ( GlobalPoint & a, const GlobalPoint & b )
40  {
41  a = GlobalPoint ( a.x() + b.x(), a.y() + b.y(), a.z() + b.z() );
42  return a;
43  }
44 
45  inline GlobalPoint & operator /= ( GlobalPoint & a, float b )
46  {
47  a = GlobalPoint ( a.x() / b, a.y() / b, a.z() / b );
48  return a;
49  }
50 
51  bool element ( const reco::TransientTrack & rt, const TransientVertex & rv )
52  {
53  const vector < reco::TransientTrack > trks = rv.originalTracks();
54  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); i!=trks.end() ; ++i )
55  {
56  if ( rt == ( *i ) ) return true;
57  };
58  return false;
59  }
60 
61  GlobalPoint toPoint ( const GlobalVector & v )
62  {
63  return GlobalPoint ( v.x(), v.y(), v.z() );
64  }
65 
66  GlobalVector toVector ( const GlobalPoint & v )
67  {
68  return GlobalVector ( v.x(), v.y(), v.z() );
69  }
70 
71 
72  GlobalPoint computeJetOrigin ( const vector < reco::TransientTrack > & trks )
73  {
75  vector< ModeFinder3d::PointAndDistance> input;
76  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
77  i!=trks.end() ; ++i )
78  {
79  input.push_back ( ModeFinder3d::PointAndDistance
80  ( i->impactPointState().globalPosition(), 1. ) );
81  }
82  return f(input);
83  }
84 
85  GlobalVector computeJetDirection ( const vector < reco::TransientTrack > & trks )
86  {
88  vector< ModeFinder3d::PointAndDistance> input;
89  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
90  i!=trks.end() ; ++i )
91  {
92  input.push_back ( ModeFinder3d::PointAndDistance (
93  toPoint ( i->impactPointState().globalMomentum() ), 1. ) );
94  }
95  GlobalPoint pt ( f(input) );
96  pt/=pt.mag();
97  return toVector ( pt );
98  }
99 
100  GlobalTrajectoryParameters computeJetTrajectory ( const vector < reco::TransientTrack > & trks )
101  {
108  if ( trks.size() == 0 ) return GlobalTrajectoryParameters();
109 
110  GlobalVector mom = computeJetDirection ( trks );
111  GlobalPoint pos = computeJetOrigin ( trks );
112 
113  GlobalTrajectoryParameters ret ( pos, mom, 0,
114  &(trks[0].impactPointState().globalParameters().magneticField()) );
115  #ifdef MVBS_DEBUG
116  DebuggingHarvester("out.txt").save ( ret , "p<sub>tot</sub>" );
117  #endif
118  return ret;
119  }
120 
121  vector < Cluster1D < reco::TransientTrack > > computeIPs (
122  const vector < reco::TransientTrack > & trks )
123  {
124  GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
125  FreeTrajectoryState axis ( jet );
126  TwoTrackMinimumDistance ttmd;
127  vector < Cluster1D < reco::TransientTrack > > pts;
128  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
129  i!=trks.end() ; ++i )
130  {
131  bool status = ttmd.calculate( axis,*( i->impactPointState().freeState() ) );
132  if (status) {
133  pair < GlobalPoint, GlobalPoint > pt = ttmd.points();
134  double d = ( pt.first - pt.second ).mag();
135  double w = 1. / ( 0.002 + d ); // hard coded weights
136  double s = ( pt.first - axis.position() ).mag();
137  Measurement1D ms ( s, 1.0 );
138  vector < const reco::TransientTrack * > trk;
139  trk.push_back ( &(*i) );
140  pts.push_back ( Cluster1D < reco::TransientTrack > ( ms, trk, w ) );
141  }
142  }
143  /*
144  #ifdef MVBS_DEBUG
145  map < string, harvest::MultiType > attrs;
146  attrs["point:mag"]=0.5;
147  attrs["point:color"]="khaki";
148  DebuggingHarvester("out.txt").save ( pts, jet, attrs, "ip" );
149  #endif
150  */
151  return pts;
152  }
153 
154  GlobalPoint computePos ( const GlobalTrajectoryParameters & jet,
155  double s )
156  {
157  GlobalPoint ret = jet.position();
158  ret += s * jet.momentum();
159  return ret;
160  }
161 
162  TransientVertex pseudoVertexFit ( const Cluster1D < reco::TransientTrack > & src,
163  bool ascending=false, bool kalmanfit=false )
164  {
165  // cout << "[MultiVertexBSeeder] debug: pseudoVertexFit with " << flush;
166  vector < const reco::TransientTrack * > trkptrs=src.tracks();
167  vector < reco::TransientTrack > trks = convert ( trkptrs );
168  // cout << trks.size() << " tracks.";
169  GlobalPoint gp;
170  GlobalError ge;
171  if ( kalmanfit )
172  {
174  gp=v.position();
175  }
176  TransientVertex ret = TransientVertex ( gp, ge, trks, -1. );
178  float w=1.0; float r=0.5;
179  if ( ascending ) { w = pow ( (float) 0.5, (int) (trks.size()-1) ); r=2.0; };
180  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
181  i!=trks.end() ; ++i )
182  {
183  mp[*i]=w;
184  w*=r;
185  }
186  ret.weightMap ( mp );
187  // cout << "[MultiVertexBSeeder] debug: return pseudoVertexFit with " << endl;
188  return ret;
189  }
190 
191  /* TransientVertex kalmanVertexFit ( const Cluster1D < reco::TransientTrack > & src )
192  {
193  KalmanVertexFitter fitter;
194  vector < const reco::TransientTrack * > trkptrs=src.tracks();
195  vector < reco::TransientTrack > trks = convert ( trkptrs );
196  return fitter.vertex ( trks );
197  return TransientVertex();
198  }*/
199 }
200 
202 {
203  return new MultiVertexBSeeder ( * this );
204 }
205 
207  theNSigma ( nsigma ) {}
208 
209 vector < TransientVertex > MultiVertexBSeeder::vertices (
210  const vector < reco::TransientTrack > & trks,
211  const reco::BeamSpot & s) const
212 {
213  return vertices ( trks );
214 }
215 
216 vector < TransientVertex > MultiVertexBSeeder::vertices (
217  const vector < reco::TransientTrack > & trks ) const
218 {
219  vector < Cluster1D < reco::TransientTrack > > ips = computeIPs ( trks );
220  /*
221  FsmwClusterizer1D < reco::TransientTrack > fsmw ( .4, theNSigma );
222  MultiClusterizer1D<reco::TransientTrack> finder ( fsmw );
223  */
225 
226  pair < vector < Cluster1D<reco::TransientTrack> >,
227  vector < const reco::TransientTrack * > > res;
228  res = finder ( ips );
229  #ifdef MVBS_DEBUG
230  // need to compute jet trajectory again :-(
231  GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
232  map < string, harvest::MultiType > attrs;
233  attrs["point:mag"]=0.75;
234  attrs["point:color"]="blue";
235  attrs["point:name"]="mode";
236  DebuggingHarvester("out.txt").save ( res.first, jet, attrs, "mode" );
237  #endif
238 
239  vector < TransientVertex > ret;
240  for ( vector< Cluster1D < reco::TransientTrack > >::const_iterator i=res.first.begin();
241  i!=res.first.end() ; ++i )
242  {
243  ret.push_back ( pseudoVertexFit ( *i, i!=res.first.begin(),
244  /* kalman fit*/ true ) );
245  }
246  return ret;
247 }
248 
249 #ifdef MVBS_DEBUG
250 #undef MVBS_DEBUG
251 #endif
int i
Definition: DBlmapReader.cc:9
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
std::pair< GlobalPoint, float > PointAndDistance
Definition: ModeFinder3d.h:17
void convert(uint32 i, char_uint32 v)
Definition: MsgTools.h:46
static std::string const input
Definition: EdmProvDump.cc:44
static const double pts[33]
Definition: Constants.h:30
std::vector< reco::TransientTrack > const & originalTracks() const
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
MultiVertexBSeeder * clone() const
double f[11][100]
std::vector< const T * > tracks() const
Definition: Cluster1D.h:49
Basic3DVector & operator/=(T t)
Scaling by a scalar value (division)
double b
Definition: hdecay.h:120
susybsm::MuonSegment ms
Definition: classes.h:31
std::vector< T > toVector(boost::python::list &l)
Definition: PythonWrapper.h:28
double a
Definition: hdecay.h:121
T w() const
tuple status
Definition: ntuplemaker.py:245
MultiVertexBSeeder(double nsigma=50.)
T x() const
Definition: PV3DBase.h:62
TransientTrackToFloatMap weightMap() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Basic3DVector & operator+=(const Basic3DVector< U > &p)
Global3DVector GlobalVector
Definition: GlobalVector.h:10
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &) const