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 #ifndef __clang__
35  inline GlobalPoint & operator += ( GlobalPoint & a, const GlobalPoint & b )
36  {
37  a = GlobalPoint ( a.x() + b.x(), a.y() + b.y(), a.z() + b.z() );
38  return a;
39  }
40 #endif
41 
42  inline GlobalPoint & operator /= ( GlobalPoint & a, float b )
43  {
44  a = GlobalPoint ( a.x() / b, a.y() / b, a.z() / b );
45  return a;
46  }
47 
48 #ifndef __clang__
49  inline bool element ( const reco::TransientTrack & rt, const TransientVertex & rv )
50  {
51  const vector < reco::TransientTrack > trks = rv.originalTracks();
52  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); i!=trks.end() ; ++i )
53  {
54  if ( rt == ( *i ) ) return true;
55  };
56  return false;
57  }
58 #endif
59 
60  GlobalPoint toPoint ( const GlobalVector & v )
61  {
62  return GlobalPoint ( v.x(), v.y(), v.z() );
63  }
64 
65  GlobalVector toVector ( const GlobalPoint & v )
66  {
67  return GlobalVector ( v.x(), v.y(), v.z() );
68  }
69 
70 
71  GlobalPoint computeJetOrigin ( const vector < reco::TransientTrack > & trks )
72  {
74  vector< ModeFinder3d::PointAndDistance> input;
75  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
76  i!=trks.end() ; ++i )
77  {
78  input.push_back ( ModeFinder3d::PointAndDistance
79  ( i->impactPointState().globalPosition(), 1. ) );
80  }
81  return f(input);
82  }
83 
84  GlobalVector computeJetDirection ( const vector < reco::TransientTrack > & trks )
85  {
87  vector< ModeFinder3d::PointAndDistance> input;
88  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
89  i!=trks.end() ; ++i )
90  {
91  input.push_back ( ModeFinder3d::PointAndDistance (
92  toPoint ( i->impactPointState().globalMomentum() ), 1. ) );
93  }
94  GlobalPoint pt ( f(input) );
95  pt/=pt.mag();
96  return toVector ( pt );
97  }
98 
99  GlobalTrajectoryParameters computeJetTrajectory ( const vector < reco::TransientTrack > & trks )
100  {
107  if ( trks.size() == 0 ) return GlobalTrajectoryParameters();
108 
109  GlobalVector mom = computeJetDirection ( trks );
110  GlobalPoint pos = computeJetOrigin ( trks );
111 
112  GlobalTrajectoryParameters ret ( pos, mom, 0,
113  &(trks[0].impactPointState().globalParameters().magneticField()) );
114  #ifdef MVBS_DEBUG
115  DebuggingHarvester("out.txt").save ( ret , "p<sub>tot</sub>" );
116  #endif
117  return ret;
118  }
119 
120  vector < Cluster1D < reco::TransientTrack > > computeIPs (
121  const vector < reco::TransientTrack > & trks )
122  {
123  GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
124  FreeTrajectoryState axis ( jet );
126  vector < Cluster1D < reco::TransientTrack > > pts;
127  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
128  i!=trks.end() ; ++i )
129  {
130  bool status = ttmd.calculate( axis,*( i->impactPointState().freeState() ) );
131  if (status) {
132  pair < GlobalPoint, GlobalPoint > pt = ttmd.points();
133  double d = ( pt.first - pt.second ).mag();
134  double w = 1. / ( 0.002 + d ); // hard coded weights
135  double s = ( pt.first - axis.position() ).mag();
136  Measurement1D ms ( s, 1.0 );
137  vector < const reco::TransientTrack * > trk;
138  trk.push_back ( &(*i) );
139  pts.push_back ( Cluster1D < reco::TransientTrack > ( ms, trk, w ) );
140  }
141  }
142  /*
143  #ifdef MVBS_DEBUG
144  map < string, harvest::MultiType > attrs;
145  attrs["point:mag"]=0.5;
146  attrs["point:color"]="khaki";
147  DebuggingHarvester("out.txt").save ( pts, jet, attrs, "ip" );
148  #endif
149  */
150  return pts;
151  }
152 
153 #ifndef __clang__
154  inline GlobalPoint computePos ( const GlobalTrajectoryParameters & jet,
155  double s )
156  {
157  GlobalPoint ret = jet.position();
158  ret += s * jet.momentum();
159  return ret;
160  }
161 #endif
162 
163  TransientVertex pseudoVertexFit ( const Cluster1D < reco::TransientTrack > & src,
164  bool ascending=false, bool kalmanfit=false )
165  {
166  // cout << "[MultiVertexBSeeder] debug: pseudoVertexFit with " << flush;
167  vector < const reco::TransientTrack * > trkptrs=src.tracks();
168  vector < reco::TransientTrack > trks = convert ( trkptrs );
169  // cout << trks.size() << " tracks.";
170  GlobalPoint gp;
171  GlobalError ge;
172  if ( kalmanfit )
173  {
175  gp=v.position();
176  }
177  TransientVertex ret = TransientVertex ( gp, ge, trks, -1. );
179  float w=1.0; float r=0.5;
180  if ( ascending ) { w = pow ( (float) 0.5, (int) (trks.size()-1) ); r=2.0; };
181  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
182  i!=trks.end() ; ++i )
183  {
184  mp[*i]=w;
185  w*=r;
186  }
187  ret.weightMap ( mp );
188  // cout << "[MultiVertexBSeeder] debug: return pseudoVertexFit with " << endl;
189  return ret;
190  }
191 
192  /* TransientVertex kalmanVertexFit ( const Cluster1D < reco::TransientTrack > & src )
193  {
194  KalmanVertexFitter fitter;
195  vector < const reco::TransientTrack * > trkptrs=src.tracks();
196  vector < reco::TransientTrack > trks = convert ( trkptrs );
197  return fitter.vertex ( trks );
198  return TransientVertex();
199  }*/
200 }
201 
203 {
204  return new MultiVertexBSeeder ( * this );
205 }
206 
208  theNSigma ( nsigma ) {}
209 
210 vector < TransientVertex > MultiVertexBSeeder::vertices (
211  const vector < reco::TransientTrack > & trks,
212  const reco::BeamSpot & s) const
213 {
214  return vertices ( trks );
215 }
216 
217 vector < TransientVertex > MultiVertexBSeeder::vertices (
218  const vector < reco::TransientTrack > & trks ) const
219 {
220  vector < Cluster1D < reco::TransientTrack > > ips = computeIPs ( trks );
221  /*
222  FsmwClusterizer1D < reco::TransientTrack > fsmw ( .4, theNSigma );
223  MultiClusterizer1D<reco::TransientTrack> finder ( fsmw );
224  */
226 
227  pair < vector < Cluster1D<reco::TransientTrack> >,
228  vector < const reco::TransientTrack * > > res;
229  res = finder ( ips );
230  #ifdef MVBS_DEBUG
231  // need to compute jet trajectory again :-(
232  GlobalTrajectoryParameters jet = computeJetTrajectory ( trks );
233  map < string, harvest::MultiType > attrs;
234  attrs["point:mag"]=0.75;
235  attrs["point:color"]="blue";
236  attrs["point:name"]="mode";
237  DebuggingHarvester("out.txt").save ( res.first, jet, attrs, "mode" );
238  #endif
239 
240  vector < TransientVertex > ret;
241  for ( vector< Cluster1D < reco::TransientTrack > >::const_iterator i=res.first.begin();
242  i!=res.first.end() ; ++i )
243  {
244  ret.push_back ( pseudoVertexFit ( *i, i!=res.first.begin(),
245  /* kalman fit*/ true ) );
246  }
247  return ret;
248 }
249 
250 #ifdef MVBS_DEBUG
251 #undef MVBS_DEBUG
252 #endif
int i
Definition: DBlmapReader.cc:9
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
const double w
Definition: UKUtility.cc:23
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
static std::string const input
Definition: EdmProvDump.cc:43
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
tuple status
Definition: ntuplemaker.py:245
MultiVertexBSeeder(double nsigma=50.)
T x() const
Definition: PV3DBase.h:62
virtual std::pair< GlobalPoint, GlobalPoint > points() const
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