CMS 3D CMS Logo

ConversionFinder.cc
Go to the documentation of this file.
6 #include "TMath.h"
7 
8 namespace egammaTools {
9 
11 
13  const reco::Track* candPartnerTk,
14  const double bFieldAtOrigin);
15 
16  const reco::Track* getElectronTrack(const reco::GsfElectron&, const float minFracSharedHits = 0.45);
17 
18  const reco::Track* getElectronTrack(const reco::GsfElectronCore&, const float minFracSharedHits = 0.45);
19 
20  bool isFromConversion(const ConversionInfo& convInfo, double maxAbsDist, double maxAbsDcot) {
21  return (std::abs(convInfo.dist) < maxAbsDist) && (std::abs(convInfo.dcot) < maxAbsDcot);
22  }
23 
24  //-----------------------------------------------------------------------------
26  const edm::Handle<reco::TrackCollection>& ctftracks_h,
27  const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
28  const double bFieldAtOrigin,
29  const double minFracSharedHits) {
30  std::vector<ConversionInfo> temp =
31  getConversionInfos(*gsfElectron.core(), ctftracks_h, gsftracks_h, bFieldAtOrigin, minFracSharedHits);
33  }
34  //-----------------------------------------------------------------------------
36  const edm::Handle<reco::TrackCollection>& ctftracks_h,
37  const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
38  const double bFieldAtOrigin,
39  const double minFracSharedHits) {
40  std::vector<ConversionInfo> temp =
41  getConversionInfos(gsfElectron, ctftracks_h, gsftracks_h, bFieldAtOrigin, minFracSharedHits);
43  }
44 
45  //-----------------------------------------------------------------------------
46  std::vector<ConversionInfo> getConversionInfos(const reco::GsfElectronCore& gsfElectron,
47  const edm::Handle<reco::TrackCollection>& ctftracks_h,
48  const edm::Handle<reco::GsfTrackCollection>& gsftracks_h,
49  const double bFieldAtOrigin,
50  const double minFracSharedHits) {
51  using namespace reco;
52  using namespace std;
53  using namespace edm;
54 
55  //get the track collections
56  const TrackCollection* ctftracks = ctftracks_h.product();
57  const GsfTrackCollection* gsftracks = gsftracks_h.product();
58 
59  //get the references to the gsf and ctf tracks that are made
60  //by the electron
61  const reco::TrackRef el_ctftrack = gsfElectron.ctfTrack();
62  const reco::GsfTrackRef& el_gsftrack = gsfElectron.gsfTrack();
63 
64  //protect against the wrong collection being passed to the function
65  if (el_ctftrack.isNonnull() && el_ctftrack.id() != ctftracks_h.id())
66  throw cms::Exception("ConversionFinderError")
67  << "ProductID of ctf track collection does not match ProductID of electron's CTF track! \n";
68  if (el_gsftrack.isNonnull() && el_gsftrack.id() != gsftracks_h.id())
69  throw cms::Exception("ConversionFinderError")
70  << "ProductID of gsf track collection does not match ProductID of electron's GSF track! \n";
71 
72  //make p4s for the electron's tracks for use later
73  LorentzVector el_ctftrack_p4;
74  if (el_ctftrack.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits)
75  el_ctftrack_p4 = LorentzVector(el_ctftrack->px(), el_ctftrack->py(), el_ctftrack->pz(), el_ctftrack->p());
76  LorentzVector el_gsftrack_p4(el_gsftrack->px(), el_gsftrack->py(), el_gsftrack->pz(), el_gsftrack->p());
77 
78  //the electron's CTF track must share at least 45% of the inner hits
79  //with the electron's GSF track
80  int ctfidx = -999.;
81  int gsfidx = -999.;
82  if (el_ctftrack.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits)
83  ctfidx = static_cast<int>(el_ctftrack.key());
84 
85  gsfidx = static_cast<int>(el_gsftrack.key());
86 
87  //these vectors are for those candidate partner tracks that pass our cuts
88  vector<ConversionInfo> v_candidatePartners;
89  //track indices required to make references
90  int ctftk_i = 0;
91  int gsftk_i = 0;
92 
93  //loop over the CTF tracks and try to find the partner track
94  for (TrackCollection::const_iterator ctftk = ctftracks->begin(); ctftk != ctftracks->end(); ctftk++, ctftk_i++) {
95  if (ctftk_i == ctfidx)
96  continue;
97 
98  //candidate track's p4
99  LorentzVector ctftk_p4 = LorentzVector(ctftk->px(), ctftk->py(), ctftk->pz(), ctftk->p());
100 
101  //apply quality cuts to remove bad tracks
102  if (ctftk->ptError() / ctftk->pt() > 0.05)
103  continue;
104  if (ctftk->numberOfValidHits() < 5)
105  continue;
106 
107  if (el_ctftrack.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits &&
108  fabs(ctftk_p4.Pt() - el_ctftrack->pt()) / el_ctftrack->pt() < 0.2)
109  continue;
110 
111  //use the electron's CTF track, if not null, to search for the partner track
112  //look only in a cone of 0.5 to save time, and require that the track is opp. sign
113  if (el_ctftrack.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits &&
114  deltaR(el_ctftrack_p4, ctftk_p4) < 0.5 && (el_ctftrack->charge() + ctftk->charge() == 0)) {
115  ConversionInfo convInfo = getConversionInfo((const reco::Track*)(el_ctftrack.get()), &(*ctftk), bFieldAtOrigin);
116 
117  //need to add the track reference information for completeness
118  //because the overloaded fnc above does not make a trackRef
119  int deltaMissingHits = ctftk->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
120  el_ctftrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
121 
122  v_candidatePartners.push_back({convInfo.dist,
123  convInfo.dcot,
124  convInfo.radiusOfConversion,
125  convInfo.pointOfConversion,
126  TrackRef(ctftracks_h, ctftk_i),
127  GsfTrackRef(),
128  deltaMissingHits,
129  0});
130 
131  } //using the electron's CTF track
132 
133  //now we check using the electron's gsf track
134  if (deltaR(el_gsftrack_p4, ctftk_p4) < 0.5 && (el_gsftrack->charge() + ctftk->charge() == 0) &&
135  el_gsftrack->ptError() / el_gsftrack->pt() < 0.25) {
136  int deltaMissingHits = ctftk->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
137  el_gsftrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
138 
139  ConversionInfo convInfo = getConversionInfo((const reco::Track*)(el_gsftrack.get()), &(*ctftk), bFieldAtOrigin);
140 
141  v_candidatePartners.push_back({convInfo.dist,
142  convInfo.dcot,
143  convInfo.radiusOfConversion,
144  convInfo.pointOfConversion,
145  TrackRef(ctftracks_h, ctftk_i),
146  GsfTrackRef(),
147  deltaMissingHits,
148  1});
149  } //using the electron's GSF track
150 
151  } //loop over the CTF track collection
152 
153  //------------------------------------------------------ Loop over GSF collection ----------------------------------//
154  for (GsfTrackCollection::const_iterator gsftk = gsftracks->begin(); gsftk != gsftracks->end(); gsftk++, gsftk_i++) {
155  //reject the electron's own gsfTrack
156  if (gsfidx == gsftk_i)
157  continue;
158 
159  LorentzVector gsftk_p4 = LorentzVector(gsftk->px(), gsftk->py(), gsftk->pz(), gsftk->p());
160 
161  //apply quality cuts to remove bad tracks
162  if (gsftk->ptError() / gsftk->pt() > 0.5)
163  continue;
164  if (gsftk->numberOfValidHits() < 5)
165  continue;
166 
167  if (fabs(gsftk->pt() - el_gsftrack->pt()) / el_gsftrack->pt() < 0.25)
168  continue;
169 
170  //try using the electron's CTF track first if it exists
171  //look only in a cone of 0.5 around the electron's track
172  //require opposite sign
173  if (el_ctftrack.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits &&
174  deltaR(el_ctftrack_p4, gsftk_p4) < 0.5 && (el_ctftrack->charge() + gsftk->charge() == 0)) {
175  int deltaMissingHits = gsftk->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
176  el_ctftrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
177 
178  ConversionInfo convInfo =
179  getConversionInfo((const reco::Track*)(el_ctftrack.get()), (const reco::Track*)(&(*gsftk)), bFieldAtOrigin);
180  //fill the Ref info
181  v_candidatePartners.push_back({convInfo.dist,
182  convInfo.dcot,
183  convInfo.radiusOfConversion,
184  convInfo.pointOfConversion,
185  TrackRef(),
186  GsfTrackRef(gsftracks_h, gsftk_i),
187  deltaMissingHits,
188  2});
189  }
190 
191  //use the electron's gsf track
192  if (deltaR(el_gsftrack_p4, gsftk_p4) < 0.5 && (el_gsftrack->charge() + gsftk->charge() == 0) &&
193  (el_gsftrack->ptError() / el_gsftrack_p4.pt() < 0.5)) {
194  ConversionInfo convInfo =
195  getConversionInfo((const reco::Track*)(el_gsftrack.get()), (const reco::Track*)(&(*gsftk)), bFieldAtOrigin);
196  //fill the Ref info
197 
198  int deltaMissingHits = gsftk->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
199  el_gsftrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
200 
201  v_candidatePartners.push_back({convInfo.dist,
202  convInfo.dcot,
203  convInfo.radiusOfConversion,
204  convInfo.pointOfConversion,
205  TrackRef(),
206  GsfTrackRef(gsftracks_h, gsftk_i),
207  deltaMissingHits,
208  3});
209  }
210  } //loop over the gsf track collection
211 
212  return v_candidatePartners;
213  }
214 
215  //-------------------------------------------------------------------------------------
217  const reco::Track* candPartnerTk,
218  const double bFieldAtOrigin) {
219  using namespace reco;
220 
221  //now calculate the conversion related information
222  LorentzVector el_tk_p4(el_track->px(), el_track->py(), el_track->pz(), el_track->p());
223  double elCurvature = -0.3 * bFieldAtOrigin * (el_track->charge() / el_tk_p4.pt()) / 100.;
224  double rEl = fabs(1. / elCurvature);
225  double xEl = -1 * (1. / elCurvature - el_track->d0()) * sin(el_tk_p4.phi());
226  double yEl = (1. / elCurvature - el_track->d0()) * cos(el_tk_p4.phi());
227 
228  LorentzVector cand_p4 =
229  LorentzVector(candPartnerTk->px(), candPartnerTk->py(), candPartnerTk->pz(), candPartnerTk->p());
230  double candCurvature = -0.3 * bFieldAtOrigin * (candPartnerTk->charge() / cand_p4.pt()) / 100.;
231  double rCand = fabs(1. / candCurvature);
232  double xCand = -1 * (1. / candCurvature - candPartnerTk->d0()) * sin(cand_p4.phi());
233  double yCand = (1. / candCurvature - candPartnerTk->d0()) * cos(cand_p4.phi());
234 
235  double d = sqrt(pow(xEl - xCand, 2) + pow(yEl - yCand, 2));
236  double dist = d - (rEl + rCand);
237  double dcot = 1. / tan(el_tk_p4.theta()) - 1. / tan(cand_p4.theta());
238 
239  //get the point of conversion
240  double xa1 = xEl + (xCand - xEl) * rEl / d;
241  double xa2 = xCand + (xEl - xCand) * rCand / d;
242  double ya1 = yEl + (yCand - yEl) * rEl / d;
243  double ya2 = yCand + (yEl - yCand) * rCand / d;
244 
245  double x = .5 * (xa1 + xa2);
246  double y = .5 * (ya1 + ya2);
247  double rconv = sqrt(pow(x, 2) + pow(y, 2));
248  double z =
249  el_track->dz() + rEl * el_track->pz() * TMath::ACos(1 - pow(rconv, 2) / (2. * pow(rEl, 2))) / el_track->pt();
250 
251  math::XYZPoint convPoint(x, y, z);
252 
253  //now assign a sign to the radius of conversion
254  float tempsign = el_track->px() * x + el_track->py() * y;
255  tempsign = tempsign / fabs(tempsign);
256  rconv = tempsign * rconv;
257 
258  //return an instance of ConversionInfo, but with a NULL track refs
259  return ConversionInfo{dist, dcot, rconv, convPoint, TrackRef(), GsfTrackRef(), -9999, -9999};
260  }
261 
262  //-------------------------------------------------------------------------------------
263  const reco::Track* getElectronTrack(const reco::GsfElectron& electron, const float minFracSharedHits) {
264  if (electron.closestCtfTrackRef().isNonnull() && electron.shFracInnerHits() > minFracSharedHits)
265  return (const reco::Track*)electron.closestCtfTrackRef().get();
266 
267  return (const reco::Track*)(electron.gsfTrack().get());
268  }
269 
270  //------------------------------------------------------------------------------------
271 
272  //takes in a vector of candidate conversion partners
273  //and arbitrates between them returning the one with the
274  //smallest R=sqrt(dist*dist + dcot*dcot)
275  ConversionInfo arbitrateConversionPartnersbyR(const std::vector<ConversionInfo>& v_convCandidates) {
276  if (v_convCandidates.size() == 1)
277  return v_convCandidates.at(0);
278 
279  double R = sqrt(pow(v_convCandidates.at(0).dist, 2) + pow(v_convCandidates.at(0).dcot, 2));
280 
281  int iArbitrated = 0;
282  int i = 0;
283 
284  for (auto const& temp : v_convCandidates) {
285  double temp_R = sqrt(pow(temp.dist, 2) + pow(temp.dcot, 2));
286  if (temp_R < R) {
287  R = temp_R;
288  iArbitrated = i;
289  }
290  ++i;
291  }
292 
293  return v_convCandidates.at(iArbitrated);
294  }
295 
296  //------------------------------------------------------------------------------------
297  ConversionInfo findBestConversionMatch(const std::vector<ConversionInfo>& v_convCandidates) {
298  using namespace std;
299 
300  if (v_convCandidates.empty())
301  return ConversionInfo{-9999.,
302  -9999.,
303  -9999.,
304  math::XYZPoint(-9999., -9999., -9999),
305  reco::TrackRef(),
307  -9999,
308  -9999};
309 
310  if (v_convCandidates.size() == 1)
311  return v_convCandidates.at(0);
312 
313  vector<ConversionInfo> v_0;
314  vector<ConversionInfo> v_1;
315  vector<ConversionInfo> v_2;
316  vector<ConversionInfo> v_3;
317  //loop over the candidates
318  for (unsigned int i = 1; i < v_convCandidates.size(); i++) {
319  ConversionInfo temp = v_convCandidates.at(i);
320 
321  if (temp.flag == 0) {
322  bool isConv = false;
323  if (fabs(temp.dist) < 0.02 && fabs(temp.dcot) < 0.02 && temp.deltaMissingHits < 3 &&
324  temp.radiusOfConversion > -2)
325  isConv = true;
326  if (sqrt(pow(temp.dist, 2) + pow(temp.dcot, 2)) < 0.05 && temp.deltaMissingHits < 2 &&
327  temp.radiusOfConversion > -2)
328  isConv = true;
329 
330  if (isConv)
331  v_0.push_back(temp);
332  }
333 
334  if (temp.flag == 1) {
335  if (sqrt(pow(temp.dist, 2) + pow(temp.dcot, 2)) < 0.05 && temp.deltaMissingHits < 2 &&
336  temp.radiusOfConversion > -2)
337  v_1.push_back(temp);
338  }
339  if (temp.flag == 2) {
340  if (sqrt(pow(temp.dist, 2) + pow(temp.dcot * temp.dcot, 2)) < 0.05 && temp.deltaMissingHits < 2 &&
341  temp.radiusOfConversion > -2)
342  v_2.push_back(temp);
343  }
344  if (temp.flag == 3) {
345  if (sqrt(temp.dist * temp.dist + temp.dcot * temp.dcot) < 0.05 && temp.deltaMissingHits < 2 &&
346  temp.radiusOfConversion > -2)
347  v_3.push_back(temp);
348  }
349 
350  } //candidate conversion loop
351 
352  //now do some arbitration
353 
354  //give preference to conversion partners found in the CTF collection
355  //using the electron's CTF track
356  if (!v_0.empty())
357  return arbitrateConversionPartnersbyR(v_0);
358 
359  if (!v_1.empty())
360  return arbitrateConversionPartnersbyR(v_1);
361 
362  if (!v_2.empty())
363  return arbitrateConversionPartnersbyR(v_2);
364 
365  if (!v_3.empty())
366  return arbitrateConversionPartnersbyR(v_3);
367 
368  //if we get here, we didn't find a candidate conversion partner that
369  //satisfied even the loose selections
370  //return the the closest partner by R
371  return arbitrateConversionPartnersbyR(v_convCandidates);
372  }
373 
374  //------------------------------------------------------------------------------------
375 
376  //------------------------------------------------------------------------------------
377  // Exists here for backwards compatibility only. Provides only the dist and dcot
378  std::pair<double, double> getConversionInfo(LorentzVector trk1_p4,
379  int trk1_q,
380  float trk1_d0,
381  LorentzVector trk2_p4,
382  int trk2_q,
383  float trk2_d0,
384  float bFieldAtOrigin) {
385  double tk1Curvature = -0.3 * bFieldAtOrigin * (trk1_q / trk1_p4.pt()) / 100.;
386  double rTk1 = fabs(1. / tk1Curvature);
387  double xTk1 = -1. * (1. / tk1Curvature - trk1_d0) * sin(trk1_p4.phi());
388  double yTk1 = (1. / tk1Curvature - trk1_d0) * cos(trk1_p4.phi());
389 
390  double tk2Curvature = -0.3 * bFieldAtOrigin * (trk2_q / trk2_p4.pt()) / 100.;
391  double rTk2 = fabs(1. / tk2Curvature);
392  double xTk2 = -1. * (1. / tk2Curvature - trk2_d0) * sin(trk2_p4.phi());
393  double yTk2 = (1. / tk2Curvature - trk2_d0) * cos(trk2_p4.phi());
394 
395  double dist = sqrt(pow(xTk1 - xTk2, 2) + pow(yTk1 - yTk2, 2));
396  dist = dist - (rTk1 + rTk2);
397 
398  double dcot = 1. / tan(trk1_p4.theta()) - 1. / tan(trk2_p4.theta());
399 
400  return std::make_pair(dist, dcot);
401  }
402 
403  //-------------------------------------- Also for backwards compatibility reasons ---------------------------------------------------------------
405  const edm::Handle<reco::TrackCollection>& track_h,
406  const double bFieldAtOrigin,
407  const double minFracSharedHits) {
408  using namespace reco;
409  using namespace std;
410  using namespace edm;
411 
412  const TrackCollection* ctftracks = track_h.product();
413  const reco::TrackRef el_ctftrack = gsfElectron.closestCtfTrackRef();
414  int ctfidx = -999.;
415  int flag = -9999.;
416  if (el_ctftrack.isNonnull() && gsfElectron.shFracInnerHits() > minFracSharedHits) {
417  ctfidx = static_cast<int>(el_ctftrack.key());
418  flag = 0;
419  } else
420  flag = 1;
421 
422  /*
423  determine whether we're going to use the CTF track or the GSF track
424  using the electron's CTF track to find the dist, dcot has been shown
425  to reduce the inefficiency
426  */
427  const reco::Track* el_track = getElectronTrack(gsfElectron, minFracSharedHits);
428  LorentzVector el_tk_p4(el_track->px(), el_track->py(), el_track->pz(), el_track->p());
429 
430  int tk_i = 0;
431  double mindcot = 9999.;
432  //make a null Track Ref
433  TrackRef candCtfTrackRef = TrackRef();
434 
435  for (TrackCollection::const_iterator tk = ctftracks->begin(); tk != ctftracks->end(); tk++, tk_i++) {
436  //if the general Track is the same one as made by the electron, skip it
437  if (tk_i == ctfidx)
438  continue;
439 
440  LorentzVector tk_p4 = LorentzVector(tk->px(), tk->py(), tk->pz(), tk->p());
441 
442  //look only in a cone of 0.5
443  double dR = deltaR(el_tk_p4, tk_p4);
444  if (dR > 0.5)
445  continue;
446 
447  //require opp. sign -> Should we use the majority logic??
448  if (tk->charge() + el_track->charge() != 0)
449  continue;
450 
451  double dcot = fabs(1. / tan(tk_p4.theta()) - 1. / tan(el_tk_p4.theta()));
452  if (dcot < mindcot) {
453  mindcot = dcot;
454  candCtfTrackRef = reco::TrackRef(track_h, tk_i);
455  }
456  } //track loop
457 
458  if (!candCtfTrackRef.isNonnull())
459  return ConversionInfo{-9999.,
460  -9999.,
461  -9999.,
462  math::XYZPoint(-9999., -9999., -9999),
463  reco::TrackRef(),
465  -9999,
466  -9999};
467 
468  //now calculate the conversion related information
469  double elCurvature = -0.3 * bFieldAtOrigin * (el_track->charge() / el_tk_p4.pt()) / 100.;
470  double rEl = fabs(1. / elCurvature);
471  double xEl = -1 * (1. / elCurvature - el_track->d0()) * sin(el_tk_p4.phi());
472  double yEl = (1. / elCurvature - el_track->d0()) * cos(el_tk_p4.phi());
473 
474  LorentzVector cand_p4 =
475  LorentzVector(candCtfTrackRef->px(), candCtfTrackRef->py(), candCtfTrackRef->pz(), candCtfTrackRef->p());
476  double candCurvature = -0.3 * bFieldAtOrigin * (candCtfTrackRef->charge() / cand_p4.pt()) / 100.;
477  double rCand = fabs(1. / candCurvature);
478  double xCand = -1 * (1. / candCurvature - candCtfTrackRef->d0()) * sin(cand_p4.phi());
479  double yCand = (1. / candCurvature - candCtfTrackRef->d0()) * cos(cand_p4.phi());
480 
481  double d = sqrt(pow(xEl - xCand, 2) + pow(yEl - yCand, 2));
482  double dist = d - (rEl + rCand);
483  double dcot = 1. / tan(el_tk_p4.theta()) - 1. / tan(cand_p4.theta());
484 
485  //get the point of conversion
486  double xa1 = xEl + (xCand - xEl) * rEl / d;
487  double xa2 = xCand + (xEl - xCand) * rCand / d;
488  double ya1 = yEl + (yCand - yEl) * rEl / d;
489  double ya2 = yCand + (yEl - yCand) * rCand / d;
490 
491  double x = .5 * (xa1 + xa2);
492  double y = .5 * (ya1 + ya2);
493  double rconv = sqrt(pow(x, 2) + pow(y, 2));
494  double z =
495  el_track->dz() + rEl * el_track->pz() * TMath::ACos(1 - pow(rconv, 2) / (2. * pow(rEl, 2))) / el_track->pt();
496 
497  math::XYZPoint convPoint(x, y, z);
498 
499  //now assign a sign to the radius of conversion
500  float tempsign = el_track->px() * x + el_track->py() * y;
501  tempsign = tempsign / fabs(tempsign);
502  rconv = tempsign * rconv;
503 
504  int deltaMissingHits = -9999;
505 
506  deltaMissingHits = candCtfTrackRef->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
508 
509  return ConversionInfo{dist, dcot, rconv, convPoint, candCtfTrackRef, GsfTrackRef(), deltaMissingHits, flag};
510  }
511 
512 } // namespace egammaTools
513 
514 //-------------------------------------------------------------------------------------
ConversionInfo
Definition: ConversionInfo.h:8
reco::GsfTrackRef
edm::Ref< GsfTrackCollection > GsfTrackRef
persistent reference to a GsfTrack
Definition: GsfTrackFwd.h:13
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
egammaTools::arbitrateConversionPartnersbyR
ConversionInfo arbitrateConversionPartnersbyR(const std::vector< ConversionInfo > &v_convCandidates)
Definition: ConversionFinder.cc:275
reco::TrackBase::p
double p() const
momentum vector magnitude
Definition: TrackBase.h:605
ConversionInfo::radiusOfConversion
const double radiusOfConversion
Definition: ConversionInfo.h:11
edm
HLT enums.
Definition: AlignableModifier.h:19
egammaTools::getConversionInfos
std::vector< ConversionInfo > getConversionInfos(const reco::GsfElectronCore &, const edm::Handle< reco::TrackCollection > &ctftracks_h, const edm::Handle< reco::GsfTrackCollection > &gsftracks_h, const double bFieldAtOrigin, const double minFracSharedHits=0.45)
Definition: ConversionFinder.cc:46
edm::Ref::get
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
reco::TrackBase::px
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:611
reco::GsfElectronCore::ctfTrack
TrackRef ctfTrack() const
Definition: GsfElectronCore.h:52
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
reco::HitPattern::numberOfLostHits
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:860
edm::Handle< reco::TrackCollection >
groupFilesInBlocks.temp
list temp
Definition: groupFilesInBlocks.py:142
egammaTools::LorentzVector
math::XYZTLorentzVector LorentzVector
Definition: ConversionFinder.cc:10
edm::Ref< TrackCollection >
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
reco::GsfElectronCore::ctfGsfOverlap
float ctfGsfOverlap() const
Definition: GsfElectronCore.h:55
reco::GsfElectron::core
virtual GsfElectronCoreRef core() const
Definition: GsfElectron.cc:8
deltaR.h
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:608
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
reco::TrackBase::py
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:614
egammaTools::findBestConversionMatch
ConversionInfo findBestConversionMatch(const std::vector< ConversionInfo > &v_convCandidates)
Definition: ConversionFinder.cc:297
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
metsig::electron
Definition: SignAlgoResolutions.h:48
reco::Track
Definition: Track.h:27
egammaTools::getConversionInfo
ConversionInfo getConversionInfo(const reco::GsfElectronCore &, const edm::Handle< reco::TrackCollection > &ctftracks_h, const edm::Handle< reco::GsfTrackCollection > &gsftracks_h, const double bFieldAtOrigin, const double minFracSharedHits=0.45)
Definition: ConversionFinder.cc:35
reco::TrackBase::charge
int charge() const
track electric charge
Definition: TrackBase.h:581
reco::TrackBase::dz
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:602
reco::GsfElectron
Definition: GsfElectron.h:35
egammaTools::isFromConversion
bool isFromConversion(const ConversionInfo &, double maxAbsDist=0.02, double maxAbsDcot=0.02)
Definition: ConversionFinder.cc:20
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
reco::TrackRef
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
egammaTools
Definition: ConversionFinder.h:34
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
ConversionInfo::pointOfConversion
const math::XYZPoint pointOfConversion
Definition: ConversionInfo.h:12
edm::Ref::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
reco::GsfElectronCore
Definition: GsfElectronCore.h:32
funct::tan
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
ConversionInfo::dist
const double dist
Definition: ConversionInfo.h:9
GsfTrack.h
reco::TrackBase::d0
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:596
reco::GsfElectron::closestCtfTrackRef
virtual TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:188
ConversionFinder.h
reco::TrackBase::hitPattern
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:489
edm::Ref::id
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
reco::GsfTrackCollection
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Definition: GsfTrackFwd.h:9
std
Definition: JetResolutionObject.h:76
GsfTrackFwd.h
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::GsfElectronCore::gsfTrack
const GsfTrackRef & gsfTrack() const
Definition: GsfElectronCore.h:48
reco::HitPattern::MISSING_INNER_HITS
Definition: HitPattern.h:155
edm::Ref::key
key_type key() const
Accessor for product key.
Definition: Ref.h:250
egammaTools::getElectronTrack
const reco::Track * getElectronTrack(const reco::GsfElectron &, const float minFracSharedHits=0.45)
Definition: ConversionFinder.cc:263
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
ztail.d
d
Definition: ztail.py:151
reco::TrackBase::pz
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:617
cms::Exception
Definition: Exception.h:70
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
dttmaxenums::R
Definition: DTTMax.h:29
edm::HandleBase::id
ProductID id() const
Definition: HandleBase.cc:13
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::GsfElectron::shFracInnerHits
float shFracInnerHits() const
Definition: GsfElectron.h:187
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
ConversionInfo::dcot
const double dcot
Definition: ConversionInfo.h:10