CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackVertexArbitratration.cc
Go to the documentation of this file.
3 using namespace reco;
4 
6  primaryVertexCollection (params.getParameter<edm::InputTag>("primaryVertices")),
7  secondaryVertexCollection (params.getParameter<edm::InputTag>("secondaryVertices")),
8  trackCollection (params.getParameter<edm::InputTag>("tracks")),
9  beamSpotCollection (params.getParameter<edm::InputTag>("beamSpot")),
10  dRCut (params.getParameter<double>("dRCut")),
11  distCut (params.getParameter<double>("distCut")),
12  sigCut (params.getParameter<double>("sigCut")),
13  dLenFraction (params.getParameter<double>("dLenFraction")),
14  fitterSigmacut (params.getParameter<double>("fitterSigmacut")),
15  fitterTini (params.getParameter<double>("fitterTini")),
16  fitterRatio (params.getParameter<double>("fitterRatio")),
17  trackMinLayers (params.getParameter<int32_t>("trackMinLayers")),
18  trackMinPt (params.getParameter<double>("trackMinPt")),
19  trackMinPixels (params.getParameter<int32_t>("trackMinPixels"))
20 {
21 
22 }
23 
25 {
26  if (track->hitPattern().trackerLayersWithMeasurement() < trackMinLayers)
27  return false;
28  if (track->pt() < trackMinPt)
29  return false;
30  if (track->hitPattern().numberOfValidPixelHits() < trackMinPixels)
31  return false;
32 
33  return true;
34 }
35 
36 
37 
38 
39 
42  const reco::Vertex &pv,
45  reco::VertexCollection & secondaryVertices)
46 {
47  using namespace reco;
48 
49  //std::cout << "PV: " << pv.position() << std::endl;
50  VertexDistance3D dist;
51  AdaptiveVertexFitter theAdaptiveFitter(
57 
58 
59 
60  reco::VertexCollection recoVertices;
61  VertexDistance3D vdist;
62 
63  std::map<unsigned int, Measurement1D> cachedIP;
64  for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices.begin();
65  sv != secondaryVertices.end(); ++sv) {
66 /* recoVertices->push_back(*sv);
67 
68 
69  for(std::vector<reco::Vertex>::iterator sv = recoVertices->begin();
70  sv != recoVertices->end(); ++sv) {
71 */
72  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
73  GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z());
74  GlobalVector flightDir = ssv-ppv;
75 // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
76  Measurement1D dlen= vdist.distance(pv,*sv);
77  std::vector<reco::TransientTrack> selTracks;
78  for(unsigned int itrack = 0; itrack < selectedTracks.size(); itrack++){
79  TrackRef ref = (selectedTracks)[itrack];
80 
81  //for(TrackCollection::const_iterator track = tracks->begin();
82  // track != tracks->end(); ++track) {
83 
84  // TrackRef ref(tracks, track - tracks->begin());
85  if (!trackFilterArbitrator(ref)) continue;
86 
87  TransientTrack tt = trackBuilder->build(ref);
88  tt.setBeamSpot(*beamSpot);
89  float w = sv->trackWeight(ref);
90  Measurement1D ipv;
91  if(cachedIP.find(itrack)!=cachedIP.end()) {
92  ipv=cachedIP[itrack];
93  } else {
94  std::pair<bool,Measurement1D> ipvp = IPTools::absoluteImpactParameter3D(tt,pv);
95  cachedIP[itrack]=ipvp.second;
96  ipv=ipvp.second;
97  }
98  std::pair<bool,Measurement1D> isv = IPTools::absoluteImpactParameter3D(tt,*sv);
99  float dR = Geom::deltaR(flightDir,tt.track()); //.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
100 
101  if( w > 0 || ( isv.second.significance() < sigCut && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction ) )
102  {
103 
104  if(( isv.second.value() < ipv.value() ) && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction
105  && dR < dRCut )
106  {
107 #ifdef VTXDEBUG
108  if(w > 0.5) std::cout << " = ";
109  else std::cout << " + ";
110 #endif
111  selTracks.push_back(tt);
112  } else
113  {
114 #ifdef VTXDEBUG
115  if(w > 0.5 && isv.second.value() > ipv.value() ) std::cout << " - ";
116  else std::cout << " ";
117 #endif
118  //add also the tracks used in previous fitting that are still closer to Sv than Pv
119  if(w > 0.5 && isv.second.value() <= ipv.value() && dR < dRCut) {
120  selTracks.push_back(tt);
121 #ifdef VTXDEBUG
122  std::cout << " = ";
123 #endif
124  }
125  if(w > 0.5 && isv.second.value() <= ipv.value() && dR >= dRCut) {
126 #ifdef VTXDEBUG
127  std::cout << " - ";
128 #endif
129 
130  }
131 
132 
133  }
134 #ifdef VTXDEBUG
135 
136  std::cout << "t : " << itrack << " ref " << ref.key() << " w: " << w
137  << " svip: " << isv.second.significance() << " " << isv.second.value()
138  << " pvip: " << ipv.significance() << " " << ipv.value() << " res " << tt.track().residualX(0) << "," << tt.track().residualY(0)
139 // << " tpvip: " << itpv.second.significance() << " " << itpv.second.value() << " dr: " << dR << std::endl;
140 #endif
141 
142  }
143  else
144  {
145 #ifdef VTXDEBUG
146 
147  std::cout << " . t : " << itrack << " ref " << ref.key()<< " w: " << w
148  << " svip: " << isv.second.significance() << " " << isv.second.value()
149  << " pvip: " << ipv.significance() << " " << ipv.value() << " dr: " << dR << std::endl;
150 #endif
151 
152  }
153  }
154 
155  if(selTracks.size() >= 2)
156  {
157  TransientVertex singleFitVertex;
158  singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv);
159  if(singleFitVertex.isValid()) { recoVertices.push_back(reco::Vertex(singleFitVertex));
160 
161 #ifdef VTXDEBUG
162  const reco::Vertex & extVertex = recoVertices.back();
163  GlobalVector vtxDir = GlobalVector(extVertex.p4().X(),extVertex.p4().Y(),extVertex.p4().Z());
164 
165 
166  std::cout << " pointing : " << Geom::deltaR(extVertex.position() - pv.position(), vtxDir) << std::endl;
167 #endif
168  }
169  }
170  }
171  return recoVertices;
172 }
173 
174 
175 
176 
177 
178 
179 
180 
181 
void setBeamSpot(const reco::BeamSpot &beamSpot)
reco::VertexCollection trackVertexArbitrator(edm::Handle< reco::BeamSpot > &beamSpot, const reco::Vertex &pv, edm::ESHandle< TransientTrackBuilder > &trackBuilder, const edm::RefVector< reco::TrackCollection > &selectedTracks, reco::VertexCollection &secondaryVertices)
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const Point & position() const
position
Definition: Vertex.h:106
double residualX(int position) const
Definition: Track.cc:14
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
bool trackFilterArbitrator(const reco::TrackRef &track) const
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double significance() const
Definition: Measurement1D.h:32
const Track & track() const
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:113
key_type key() const
Accessor for product key.
Definition: Ref.h:266
double value() const
Definition: Measurement1D.h:28
double residualY(int position) const
Definition: Track.cc:19
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
tuple cout
Definition: gather_cfg.py:121
T w() const
bool isValid() const
TrackVertexArbitration(const edm::ParameterSet &params)
Global3DVector GlobalVector
Definition: GlobalVector.h:10