CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TracksClusteringFromDisplacedSeed.cc
Go to the documentation of this file.
2 //#define VTXDEBUG 1
3 
4 
6 // maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
7  min3DIPSignificance(params.getParameter<double>("seedMin3DIPSignificance")),
8  min3DIPValue(params.getParameter<double>("seedMin3DIPValue")),
9  clusterMaxDistance(params.getParameter<double>("clusterMaxDistance")),
10  clusterMaxSignificance(params.getParameter<double>("clusterMaxSignificance")), //3
11  clusterScale(params.getParameter<double>("clusterScale")),//10.
12  clusterMinAngleCosine(params.getParameter<double>("clusterMinAngleCosine")) //0.0
13 
14 {
15 
16 }
17 
18 std::pair<std::vector<reco::TransientTrack>,GlobalPoint> TracksClusteringFromDisplacedSeed::nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const
19 {
20  VertexDistance3D distanceComputer;
21  GlobalPoint pv(primaryVertex.position().x(),primaryVertex.position().y(),primaryVertex.position().z());
22  std::vector<reco::TransientTrack> result;
23  TwoTrackMinimumDistance dist;
24  GlobalPoint seedingPoint;
25  float sumWeights=0;
26  std::pair<bool,Measurement1D> ipSeed = IPTools::absoluteImpactParameter3D(seed,primaryVertex);
27  float pvDistance = ipSeed.second.value();
28 // float densityFactor = 2./sqrt(20.*tracks.size()); // assuming all tracks being in 2 narrow jets of cone 0.3
29  float densityFactor = 2./sqrt(20.*80); // assuming 80 tracks being in 2 narrow jets of cone 0.3
30  for(std::vector<reco::TransientTrack>::const_iterator tt = tracks.begin();tt!=tracks.end(); ++tt ) {
31 
32  if(*tt==seed) continue;
33 
34  std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*tt,primaryVertex);
35  if(dist.calculate(tt->impactPointState(),seed.impactPointState()))
36  {
37  GlobalPoint ttPoint = dist.points().first;
38  GlobalError ttPointErr = tt->impactPointState().cartesianError().position();
39  GlobalPoint seedPosition = dist.points().second;
40  GlobalError seedPositionErr = seed.impactPointState().cartesianError().position();
41  Measurement1D m = distanceComputer.distance(VertexState(seedPosition,seedPositionErr), VertexState(ttPoint, ttPointErr));
42  GlobalPoint cp(dist.crossingPoint());
43 
44 
45  float distanceFromPV = (dist.points().second-pv).mag();
46  float distance = dist.distance();
47  GlobalVector trackDir2D(tt->impactPointState().globalDirection().x(),tt->impactPointState().globalDirection().y(),0.);
49  //SK:UNUSED// float dotprodTrackSeed2D = trackDir2D.unit().dot(seedDir2D.unit());
50 
51  float dotprodTrack = (dist.points().first-pv).unit().dot(tt->impactPointState().globalDirection().unit());
52  float dotprodSeed = (dist.points().second-pv).unit().dot(seed.impactPointState().globalDirection().unit());
53 
54  float w = distanceFromPV*distanceFromPV/(pvDistance*distance);
55  bool selected = (m.significance() < clusterMaxSignificance &&
56  dotprodSeed > clusterMinAngleCosine && //Angles between PV-PCAonSeed vectors and seed directions
57  dotprodTrack > clusterMinAngleCosine && //Angles between PV-PCAonTrack vectors and track directions
58 // dotprodTrackSeed2D > clusterMinAngleCosine && //Angle between track and seed
59  // distance*clusterScale*tracks.size() < (distanceFromPV+pvDistance)*(distanceFromPV+pvDistance)/pvDistance && // cut scaling with track density
60  distance*clusterScale < densityFactor*distanceFromPV && // cut scaling with track density
61  distance < clusterMaxDistance); // absolute distance cut
62 
63 #ifdef VTXDEBUG
64  std::cout << tt->trackBaseRef().key() << " : " << (selected?"+":" ")<< " " << m.significance() << " < " << clusterMaxSignificance << " && " <<
65  dotprodSeed << " > " << clusterMinAngleCosine << " && " <<
66  dotprodTrack << " > " << clusterMinAngleCosine << " && " <<
67  dotprodTrackSeed2D << " > " << clusterMinAngleCosine << " && " <<
68  distance*clusterScale << " < " << densityFactor*distanceFromPV << " crossingtoPV: " << distanceFromPV << " dis*scal " << distance*clusterScale << " < " << densityFactor*distanceFromPV << " dist: " << distance << " < " << clusterMaxDistance << std::endl; // cut scaling with track density
69 #endif
70  if(selected)
71  {
72  result.push_back(*tt);
73  seedingPoint = GlobalPoint(cp.x()*w+seedingPoint.x(),cp.y()*w+seedingPoint.y(),cp.z()*w+seedingPoint.z());
74  sumWeights+=w;
75  }
76  }
77  }
78 
79  seedingPoint = GlobalPoint(seedingPoint.x()/sumWeights,seedingPoint.y()/sumWeights,seedingPoint.z()/sumWeights);
80  return std::pair<std::vector<reco::TransientTrack>,GlobalPoint>(result,seedingPoint);
81 
82 }
83 
84 
85 
86 
87 
88 std::vector<TracksClusteringFromDisplacedSeed::Cluster> TracksClusteringFromDisplacedSeed::clusters(
89  const reco::Vertex &pv,
90  const std::vector<reco::TransientTrack> & selectedTracks
91  )
92 {
93  using namespace reco;
94  std::vector<TransientTrack> seeds;
95  for(std::vector<TransientTrack>::const_iterator it = selectedTracks.begin(); it != selectedTracks.end(); it++){
96  std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*it,pv);
97  if(ip.first && ip.second.value() >= min3DIPValue && ip.second.significance() >= min3DIPSignificance)
98  {
99 #ifdef VTXDEBUG
100  std::cout << "new seed " << it-selectedTracks.begin() << " ref " << it->trackBaseRef().key() << " " << ip.second.value() << " " << ip.second.significance() << " " << it->track().hitPattern().trackerLayersWithMeasurement() << " " << it->track().pt() << " " << it->track().eta() << std::endl;
101 #endif
102  seeds.push_back(*it);
103  }
104 
105  }
106 
107  std::vector< Cluster > clusters;
108  int i = 0;
109  for(std::vector<TransientTrack>::const_iterator s = seeds.begin();
110  s != seeds.end(); ++s, ++i)
111  {
112 #ifdef VTXDEBUG
113  std::cout << "Seed N. "<<i << std::endl;
114 #endif // VTXDEBUG
115  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> ntracks = nearTracks(*s,selectedTracks,pv);
116 // std::cout << ntracks.first.size() << " " << ntracks.first.size() << std::endl;
117 // if(ntracks.first.size() == 0 || ntracks.first.size() > maxNTracks ) continue;
118  ntracks.first.push_back(*s);
119  Cluster aCl;
120  aCl.seedingTrack = *s;
121  aCl.seedPoint = ntracks.second;
122  aCl.tracks = ntracks.first;
123  clusters.push_back(aCl);
124  }
125 
126 return clusters;
127 }
128 
virtual Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > nearTracks(const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
int i
Definition: DBlmapReader.cc:9
std::vector< Cluster > clusters(const reco::Vertex &pv, const std::vector< reco::TransientTrack > &selectedTracks)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const CartesianTrajectoryError cartesianError() const
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
T y() const
Definition: PV3DBase.h:63
const Point & position() const
position
Definition: Vertex.h:93
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double significance() const
Definition: Measurement1D.h:32
const GlobalError position() const
Position error submatrix.
tuple tracks
Definition: testEve_cfg.py:39
TracksClusteringFromDisplacedSeed(const edm::ParameterSet &params)
tuple cout
Definition: gather_cfg.py:121
TrajectoryStateOnSurface impactPointState() const
T w() const
T x() const
Definition: PV3DBase.h:62
GlobalVector globalDirection() const