CMS 3D CMS Logo

TracksClusteringFromDisplacedSeed.cc
Go to the documentation of this file.
2 //#define VTXDEBUG 1
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  maxTimeSignificance(params.getParameter<double>("maxTimeSignificance"))
16 
17 {
18 
19 }
20 
21 std::pair<std::vector<reco::TransientTrack>,GlobalPoint> TracksClusteringFromDisplacedSeed::nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const
22 {
23  VertexDistance3D distanceComputer;
24  GlobalPoint pv(primaryVertex.position().x(),primaryVertex.position().y(),primaryVertex.position().z());
25  std::vector<reco::TransientTrack> result;
27  GlobalPoint seedingPoint;
28  float sumWeights=0;
29  std::pair<bool,Measurement1D> ipSeed = IPTools::absoluteImpactParameter3D(seed,primaryVertex);
30  float pvDistance = ipSeed.second.value();
31  for(std::vector<reco::TransientTrack>::const_iterator tt = tracks.begin();tt!=tracks.end(); ++tt ) {
32 
33  if(*tt==seed) continue;
34 
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  double timeSig = 0.;
45  if( primaryVertex.covariance(3,3) > 0. &&
46  edm::isFinite(seed.timeExt()) && edm::isFinite(tt->timeExt()) ) {
47  // apply only if time available and being used in vertexing
48  const double tError = std::sqrt( std::pow(seed.dtErrorExt(),2) + std::pow(tt->dtErrorExt(),2) );
49  timeSig = std::abs( seed.timeExt() - tt->timeExt() ) / tError;
50  }
51 
52  float distanceFromPV = (dist.points().second-pv).mag();
53  float distance = dist.distance();
54  GlobalVector trackDir2D(tt->impactPointState().globalDirection().x(),tt->impactPointState().globalDirection().y(),0.);
56  //SK:UNUSED// float dotprodTrackSeed2D = trackDir2D.unit().dot(seedDir2D.unit());
57 
58  float dotprodTrack = (dist.points().first-pv).unit().dot(tt->impactPointState().globalDirection().unit());
59  float dotprodSeed = (dist.points().second-pv).unit().dot(seed.impactPointState().globalDirection().unit());
60 
61  float w = distanceFromPV*distanceFromPV/(pvDistance*distance);
62  bool selected = (m.significance() < clusterMaxSignificance &&
63  ((clusterMinAngleCosine > 0) ? (dotprodSeed > clusterMinAngleCosine) : (dotprodSeed < clusterMinAngleCosine)) && //Angles between PV-PCAonSeed vectors and seed directions
64  ((clusterMinAngleCosine > 0) ? (dotprodTrack > clusterMinAngleCosine) : (dotprodTrack < clusterMinAngleCosine)) && //Angles between PV-PCAonTrack vectors and track directions
65  //dotprodTrackSeed2D > clusterMinAngleCosine && //Angle between track and seed
66  //distance*clusterScale*tracks.size() < (distanceFromPV+pvDistance)*(distanceFromPV+pvDistance)/pvDistance && // cut scaling with track density
67  distance*distanceRatio < distanceFromPV && // cut scaling with track density
68  distance < clusterMaxDistance &&
69  timeSig < maxTimeSignificance); // absolute distance cut
70 
71 #ifdef VTXDEBUG
72  std::cout << tt->trackBaseRef().key() << " : " << (selected?"+":" ")<< " " << m.significance() << " < " << clusterMaxSignificance << " && " <<
73  dotprodSeed << " > " << clusterMinAngleCosine << " && " <<
74  dotprodTrack << " > " << clusterMinAngleCosine << " && " <<
75  dotprodTrackSeed2D << " > " << clusterMinAngleCosine << " && " <<
76  distance*distanceRatio << " < " << distanceFromPV << " crossingtoPV: " << distanceFromPV << " dis*scal " << distance*distanceRatio << " < " << distanceFromPV << " dist: " << distance << " < " << clusterMaxDistance <<
77  << "timeSig: " << timeSig << std::endl; // cut scaling with track density
78 #endif
79  if(selected)
80  {
81  result.push_back(*tt);
82  seedingPoint = GlobalPoint(cp.x()*w+seedingPoint.x(),cp.y()*w+seedingPoint.y(),cp.z()*w+seedingPoint.z());
83  sumWeights+=w;
84  }
85  }
86  }
87 
88  seedingPoint = GlobalPoint(seedingPoint.x()/sumWeights,seedingPoint.y()/sumWeights,seedingPoint.z()/sumWeights);
89  return std::pair<std::vector<reco::TransientTrack>,GlobalPoint>(result,seedingPoint);
90 
91 }
92 
93 
94 
95 
96 
97 std::vector<TracksClusteringFromDisplacedSeed::Cluster> TracksClusteringFromDisplacedSeed::clusters(
98  const reco::Vertex &pv,
99  const std::vector<reco::TransientTrack> & selectedTracks
100  )
101 {
102  using namespace reco;
103  std::vector<TransientTrack> seeds;
104  for(std::vector<TransientTrack>::const_iterator it = selectedTracks.begin(); it != selectedTracks.end(); it++){
105  std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*it,pv);
106  if(ip.first && ip.second.value() >= min3DIPValue && ip.second.significance() >= min3DIPSignificance && ip.second.value() <= max3DIPValue && ip.second.significance() <= max3DIPSignificance)
107  {
108 #ifdef VTXDEBUG
109  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;
110 #endif
111  seeds.push_back(*it);
112  }
113 
114  }
115 
116  std::vector< Cluster > clusters;
117  int i = 0;
118  for(std::vector<TransientTrack>::const_iterator s = seeds.begin();
119  s != seeds.end(); ++s, ++i)
120  {
121 #ifdef VTXDEBUG
122  std::cout << "Seed N. "<<i << std::endl;
123 #endif // VTXDEBUG
124  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> ntracks = nearTracks(*s,selectedTracks,pv);
125 // std::cout << ntracks.first.size() << " " << ntracks.first.size() << std::endl;
126 // if(ntracks.first.size() == 0 || ntracks.first.size() > maxNTracks ) continue;
127  ntracks.first.push_back(*s);
128  Cluster aCl;
129  aCl.seedingTrack = *s;
130  aCl.seedPoint = ntracks.second;
131  aCl.tracks = ntracks.first;
132  clusters.push_back(aCl);
133  }
134 
135 return clusters;
136 }
137 
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > nearTracks(const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
std::pair< GlobalPoint, GlobalPoint > points() const override
const double w
Definition: UKUtility.cc:23
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
double timeExt() const
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
const Point & position() const
position
Definition: Vertex.h:109
bool isFinite(T x)
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GlobalPoint crossingPoint() const override
float y() const
Same as rapidity().
Definition: CaloParticle.h:176
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double significance() const
Definition: Measurement1D.h:32
float distance() const override
const GlobalError position() const
Position error submatrix.
double dtErrorExt() const
TracksClusteringFromDisplacedSeed(const edm::ParameterSet &params)
fixed size matrix
TrajectoryStateOnSurface impactPointState() const
T x() const
Definition: PV3DBase.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
GlobalVector globalDirection() const