CMS 3D CMS Logo

TrackVertexArbitration.h
Go to the documentation of this file.
1 #ifndef TrackVertexArbitration_H
2 #define TrackVertexArbitration_H
3 #include <memory>
4 #include <set>
5 #include <unordered_map>
6 
12 
17 
24 
30 
39 
41 
43 //#include "DataFormats/Math/interface/deltaR.h"
44 
45 //#define VTXDEBUG
46 
47 namespace svhelper {
48  inline double cov33(const reco::Vertex &sv) { return sv.covariance(3, 3); }
49  inline double cov33(const reco::VertexCompositePtrCandidate &sv) { return sv.vertexCovariance(3, 3); }
50 } // namespace svhelper
51 
52 template <class VTX>
54 public:
56 
58  const reco::Vertex &pv,
59  std::vector<reco::TransientTrack> &selectedTracks,
60  std::vector<VTX> &secondaryVertices);
61 
62 private:
64 
69  double dRCut;
70  double distCut;
71  double sigCut;
72  double dLenFraction;
73  double fitterSigmacut; // = 3.0;
74  double fitterTini; // = 256.;
75  double fitterRatio; // = 0.25;
77  double trackMinPt;
80  double sign;
81 };
82 
84 template <class VTX>
86  : primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
87  secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
88  trackCollection(params.getParameter<edm::InputTag>("tracks")),
89  beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")),
90  dRCut(params.getParameter<double>("dRCut")),
91  distCut(params.getParameter<double>("distCut")),
92  sigCut(params.getParameter<double>("sigCut")),
93  dLenFraction(params.getParameter<double>("dLenFraction")),
94  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
95  fitterTini(params.getParameter<double>("fitterTini")),
96  fitterRatio(params.getParameter<double>("fitterRatio")),
97  trackMinLayers(params.getParameter<int32_t>("trackMinLayers")),
98  trackMinPt(params.getParameter<double>("trackMinPt")),
99  trackMinPixels(params.getParameter<int32_t>("trackMinPixels")),
100  maxTimeSignificance(params.getParameter<double>("maxTimeSignificance")) {
101  sign = 1.0;
102  if (dRCut < 0) {
103  sign = -1.0;
104  }
105  dRCut *= dRCut;
106 }
107 template <class VTX>
109  if (!track.isValid())
110  return false;
111  if (track.track().hitPattern().trackerLayersWithMeasurement() < trackMinLayers)
112  return false;
113  if (track.track().pt() < trackMinPt)
114  return false;
115  if (track.track().hitPattern().numberOfValidPixelHits() < trackMinPixels)
116  return false;
117 
118  return true;
119 }
120 
121 inline float trackWeight(const reco::Vertex &sv, const reco::TransientTrack &track) {
122  return sv.trackWeight(track.trackBaseRef());
123 }
126  dynamic_cast<const reco::CandidatePtrTransientTrack *>(tt.basicTransientTrack());
127  if (cptt == nullptr)
128  edm::LogError("DynamicCastingFailed") << "Casting of TransientTrack to CandidatePtrTransientTrack failed!";
129  else {
130  const reco::CandidatePtr &c = cptt->candidate();
131  if (std::find(sv.daughterPtrVector().begin(), sv.daughterPtrVector().end(), c) != sv.daughterPtrVector().end())
132  return 1.0;
133  else
134  return 0;
135  }
136  return 0;
137 }
138 
139 template <class VTX>
141  const reco::Vertex &pv,
142  std::vector<reco::TransientTrack> &selectedTracks,
143  std::vector<VTX> &secondaryVertices) {
144  using namespace reco;
145 
146  //std::cout << "PV: " << pv.position() << std::endl;
147  VertexDistance3D dist;
153 
154  std::vector<VTX> recoVertices;
155  VertexDistance3D vdist;
156  GlobalPoint ppv(pv.position().x(), pv.position().y(), pv.position().z());
157 
158  std::unordered_map<unsigned int, Measurement1D> cachedIP;
159  for (typename std::vector<VTX>::const_iterator sv = secondaryVertices.begin(); sv != secondaryVertices.end(); ++sv) {
160  const double svTime(sv->t()), svTimeCov(svhelper::cov33(*sv));
161  const bool svHasTime = svTimeCov > 0.;
162 
163  GlobalPoint ssv(sv->position().x(), sv->position().y(), sv->position().z());
164  GlobalVector flightDir = ssv - ppv;
165  // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
166  Measurement1D dlen =
168  std::vector<reco::TransientTrack> selTracks;
169  for (unsigned int itrack = 0; itrack < selectedTracks.size(); itrack++) {
170  TransientTrack &tt = (selectedTracks)[itrack];
171  if (!trackFilterArbitrator(tt))
172  continue;
173  tt.setBeamSpot(*beamSpot);
174  float w = trackWeight(*sv, tt);
175  Measurement1D ipv;
176  if (cachedIP.count(itrack)) {
177  ipv = cachedIP[itrack];
178  } else {
179  std::pair<bool, Measurement1D> ipvp = IPTools::absoluteImpactParameter3D(tt, pv);
180  cachedIP[itrack] = ipvp.second;
181  ipv = ipvp.second;
182  }
183 
184  AnalyticalImpactPointExtrapolator extrapolator(tt.field());
186  extrapolator.extrapolate(tt.impactPointState(), RecoVertex::convertPos(sv->position()));
187  if (!tsos.isValid())
188  continue;
189  GlobalPoint refPoint = tsos.globalPosition();
190  GlobalError refPointErr = tsos.cartesianError().position();
191  Measurement1D isv =
193  VertexState(refPoint, refPointErr));
194 
195  float dR = Geom::deltaR2(((sign > 0) ? flightDir : -flightDir),
196  tt.track()); //.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
197 
198  double timeSig = 0.;
199  if (svHasTime && edm::isFinite(tt.timeExt())) {
200  double tError = std::sqrt(std::pow(tt.dtErrorExt(), 2) + svTimeCov);
201  timeSig = std::abs(tt.timeExt() - svTime) / tError;
202  }
203 
204  if (w > 0 || (isv.significance() < sigCut && isv.value() < distCut && isv.value() < dlen.value() * dLenFraction &&
205  timeSig < maxTimeSignificance)) {
206  if ((isv.value() < ipv.value()) && isv.value() < distCut && isv.value() < dlen.value() * dLenFraction &&
207  dR < dRCut && timeSig < maxTimeSignificance) {
208 #ifdef VTXDEBUG
209  if (w > 0.5)
210  std::cout << " = ";
211  else
212  std::cout << " + ";
213 #endif
214  selTracks.push_back(tt);
215  } else {
216 #ifdef VTXDEBUG
217  if (w > 0.5 && isv.value() > ipv.value())
218  std::cout << " - ";
219  else
220  std::cout << " ";
221 #endif
222  //add also the tracks used in previous fitting that are still closer to Sv than Pv
223  if (w > 0.5 && isv.value() <= ipv.value() && dR < dRCut && timeSig < maxTimeSignificance) {
224  selTracks.push_back(tt);
225 #ifdef VTXDEBUG
226  std::cout << " = ";
227 #endif
228  }
229  if (w > 0.5 && isv.value() <= ipv.value() && dR >= dRCut) {
230 #ifdef VTXDEBUG
231  std::cout << " - ";
232 #endif
233  }
234  }
235 #ifdef VTXDEBUG
236 
237  std::cout << "t : " << itrack << " ref " << ref.key() << " w: " << w << " svip: " << isv.significance() << " "
238  << isv.value() << " pvip: " << ipv.significance() << " " << ipv.value() << " res "
239  << tt.track().residualX(0) << "," << tt.track().residualY(0)
240 // << " tpvip: " << itpv.second.significance() << " " << itpv.second.value() << " dr: " << dR << std::endl;
241 #endif
242 
243  } else {
244 #ifdef VTXDEBUG
245 
246  std::cout << " . t : " << itrack << " ref " << ref.key() << " w: " << w
247  << " svip: " << isv.second.significance() << " " << isv.second.value()
248  << " pvip: " << ipv.significance() << " " << ipv.value() << " dr: " << dR << std::endl;
249 #endif
250  }
251  }
252 
253  if (selTracks.size() >= 2) {
254  TransientVertex singleFitVertex;
255  singleFitVertex = theAdaptiveFitter.vertex(selTracks, ssv);
256 
257  if (singleFitVertex.isValid()) {
258  if (pv.covariance(3, 3) > 0.)
259  svhelper::updateVertexTime(singleFitVertex);
260  recoVertices.push_back(VTX(singleFitVertex));
261 
262 #ifdef VTXDEBUG
263  const VTX &extVertex = recoVertices.back();
264  GlobalVector vtxDir = GlobalVector(extVertex.p4().X(), extVertex.p4().Y(), extVertex.p4().Z());
265 
266  std::cout << " pointing : " << Geom::deltaR(extVertex.position() - pv.position(), vtxDir) << std::endl;
267 #endif
268  }
269  }
270  }
271  return recoVertices;
272 }
273 
274 #endif
reco::Vertex::Point convertPos(const GlobalPoint &p)
T w() const
CandidatePtr candidate() const override
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:38
Log< level::Error, false > LogError
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
constexpr bool isFinite(T x)
Definition: TTTypes.h:54
bool isValid() const
T sqrt(T t)
Definition: SSEVec.h:23
edm::InputTag secondaryVertexCollection
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackCollection
Definition: JetHT_cfg.py:51
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
std::vector< VTX > trackVertexArbitrator(edm::Handle< reco::BeamSpot > &beamSpot, const reco::Vertex &pv, std::vector< reco::TransientTrack > &selectedTracks, std::vector< VTX > &secondaryVertices)
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
fixed size matrix
bool trackFilterArbitrator(const reco::TransientTrack &track) const
HLT enums.
double cov33(const reco::Vertex &sv)
TrackVertexArbitration(const edm::ParameterSet &params)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
float trackWeight(const reco::Vertex &sv, const reco::TransientTrack &track)
Global3DVector GlobalVector
Definition: GlobalVector.h:10