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  max3DIPSignificance(params.getParameter<double>("seedMax3DIPSignificance")),
8  max3DIPValue(params.getParameter<double>("seedMax3DIPValue")),
9  min3DIPSignificance(params.getParameter<double>("seedMin3DIPSignificance")),
10  min3DIPValue(params.getParameter<double>("seedMin3DIPValue")),
11  clusterMaxDistance(params.getParameter<double>("clusterMaxDistance")),
12  clusterMaxSignificance(params.getParameter<double>("clusterMaxSignificance")), //3
13  distanceRatio(params.getParameter<double>("distanceRatio")),//was clusterScale/densityFactor
14  clusterMinAngleCosine(params.getParameter<double>("clusterMinAngleCosine")) //0.0
15 
16 {
17 
18 }
19 
20 std::pair<std::vector<reco::TransientTrack>,GlobalPoint> TracksClusteringFromDisplacedSeed::nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const
21 {
22  VertexDistance3D distanceComputer;
23  GlobalPoint pv(primaryVertex.position().x(),primaryVertex.position().y(),primaryVertex.position().z());
24  std::vector<reco::TransientTrack> result;
26  GlobalPoint seedingPoint;
27  float sumWeights=0;
28  std::pair<bool,Measurement1D> ipSeed = IPTools::absoluteImpactParameter3D(seed,primaryVertex);
29  float pvDistance = ipSeed.second.value();
30  for(std::vector<reco::TransientTrack>::const_iterator tt = tracks.begin();tt!=tracks.end(); ++tt ) {
31 
32  if(*tt==seed) continue;
33 
34  if(dist.calculate(tt->impactPointState(),seed.impactPointState()))
35  {
36  GlobalPoint ttPoint = dist.points().first;
37  GlobalError ttPointErr = tt->impactPointState().cartesianError().position();
38  GlobalPoint seedPosition = dist.points().second;
39  GlobalError seedPositionErr = seed.impactPointState().cartesianError().position();
40  Measurement1D m = distanceComputer.distance(VertexState(seedPosition,seedPositionErr), VertexState(ttPoint, ttPointErr));
41  GlobalPoint cp(dist.crossingPoint());
42 
43 
44  float distanceFromPV = (dist.points().second-pv).mag();
45  float distance = dist.distance();
46  GlobalVector trackDir2D(tt->impactPointState().globalDirection().x(),tt->impactPointState().globalDirection().y(),0.);
48  //SK:UNUSED// float dotprodTrackSeed2D = trackDir2D.unit().dot(seedDir2D.unit());
49 
50  float dotprodTrack = (dist.points().first-pv).unit().dot(tt->impactPointState().globalDirection().unit());
51  float dotprodSeed = (dist.points().second-pv).unit().dot(seed.impactPointState().globalDirection().unit());
52 
53  float w = distanceFromPV*distanceFromPV/(pvDistance*distance);
54  bool selected = (m.significance() < clusterMaxSignificance &&
55  dotprodSeed > clusterMinAngleCosine && //Angles between PV-PCAonSeed vectors and seed directions
56  dotprodTrack > clusterMinAngleCosine && //Angles between PV-PCAonTrack vectors and track directions
57 // dotprodTrackSeed2D > clusterMinAngleCosine && //Angle between track and seed
58  // distance*clusterScale*tracks.size() < (distanceFromPV+pvDistance)*(distanceFromPV+pvDistance)/pvDistance && // cut scaling with track density
59  distance*distanceRatio < distanceFromPV && // cut scaling with track density
60  distance < clusterMaxDistance); // absolute distance cut
61 
62 #ifdef VTXDEBUG
63  std::cout << tt->trackBaseRef().key() << " : " << (selected?"+":" ")<< " " << m.significance() << " < " << clusterMaxSignificance << " && " <<
64  dotprodSeed << " > " << clusterMinAngleCosine << " && " <<
65  dotprodTrack << " > " << clusterMinAngleCosine << " && " <<
66  dotprodTrackSeed2D << " > " << clusterMinAngleCosine << " && " <<
67  distance*distanceRatio << " < " << distanceFromPV << " crossingtoPV: " << distanceFromPV << " dis*scal " << distance*distanceRatio << " < " << distanceFromPV << " dist: " << distance << " < " << clusterMaxDistance << std::endl; // cut scaling with track density
68 #endif
69  if(selected)
70  {
71  result.push_back(*tt);
72  seedingPoint = GlobalPoint(cp.x()*w+seedingPoint.x(),cp.y()*w+seedingPoint.y(),cp.z()*w+seedingPoint.z());
73  sumWeights+=w;
74  }
75  }
76  }
77 
78  seedingPoint = GlobalPoint(seedingPoint.x()/sumWeights,seedingPoint.y()/sumWeights,seedingPoint.z()/sumWeights);
79  return std::pair<std::vector<reco::TransientTrack>,GlobalPoint>(result,seedingPoint);
80 
81 }
82 
83 
84 
85 
86 
87 std::vector<TracksClusteringFromDisplacedSeed::Cluster> TracksClusteringFromDisplacedSeed::clusters(
88  const reco::Vertex &pv,
89  const std::vector<reco::TransientTrack> & selectedTracks
90  )
91 {
92  using namespace reco;
93  std::vector<TransientTrack> seeds;
94  for(std::vector<TransientTrack>::const_iterator it = selectedTracks.begin(); it != selectedTracks.end(); it++){
95  std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*it,pv);
96  if(ip.first && ip.second.value() >= min3DIPValue && ip.second.significance() >= min3DIPSignificance && ip.second.value() <= max3DIPValue && ip.second.significance() <= max3DIPSignificance)
97  {
98 #ifdef VTXDEBUG
99  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;
100 #endif
101  seeds.push_back(*it);
102  }
103 
104  }
105 
106  std::vector< Cluster > clusters;
107  int i = 0;
108  for(std::vector<TransientTrack>::const_iterator s = seeds.begin();
109  s != seeds.end(); ++s, ++i)
110  {
111 #ifdef VTXDEBUG
112  std::cout << "Seed N. "<<i << std::endl;
113 #endif // VTXDEBUG
114  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> ntracks = nearTracks(*s,selectedTracks,pv);
115 // std::cout << ntracks.first.size() << " " << ntracks.first.size() << std::endl;
116 // if(ntracks.first.size() == 0 || ntracks.first.size() > maxNTracks ) continue;
117  ntracks.first.push_back(*s);
118  Cluster aCl;
119  aCl.seedingTrack = *s;
120  aCl.seedPoint = ntracks.second;
121  aCl.tracks = ntracks.first;
122  clusters.push_back(aCl);
123  }
124 
125 return clusters;
126 }
127 
virtual float distance() const
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
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
const double w
Definition: UKUtility.cc:23
std::vector< Cluster > clusters(const reco::Vertex &pv, const std::vector< reco::TransientTrack > &selectedTracks)
virtual GlobalPoint crossingPoint() const
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:99
tuple result
Definition: mps_fire.py:84
string unit
Definition: csvLumiCalc.py:46
T z() const
Definition: PV3DBase.h:64
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:145
TrajectoryStateOnSurface impactPointState() const
T x() const
Definition: PV3DBase.h:62
virtual std::pair< GlobalPoint, GlobalPoint > points() const
GlobalVector globalDirection() const