CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ConversionFinder.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaElectronAlgos_ConversionFinder_h
2 #define RecoEgamma_EgammaElectronAlgos_ConversionFinder_h
3 
17 
18 #include <iostream>
19 #include <optional>
20 
21 /*
22  Class Looks for oppositely charged track in the
23  track collection with the minimum delta cot(theta) between the track
24  and the electron's CTF track (if it doesn't exist, we use the electron's
25  GSF track). Calculate the dist, dcot, point of conversion and the
26  radius of conversion for this pair and fill the ConversionInfo
27 */
28 
29 namespace egamma::conv {
30 
31  struct ConversionInfo {
32  const float dist = -9999.;
33  const float dcot = -9999.;
34  const float radiusOfConversion = -9999.;
35  // if the partner track is found in the GSF track collection,
36  // this is a ref to the GSF partner track
37  const std::optional<int> conversionPartnerCtfTkIdx = std::nullopt;
38  // if the partner track is found in the CTF track collection,
39  // this is a ref to the CTF partner track
40  const std::optional<int> conversionPartnerGsfTkIdx = std::nullopt;
41  const int deltaMissingHits = -9999;
42  const int flag = -9999;
43 
44  // flag 0: Partner track found in the CTF collection using the electron's CTF track
45  // flag 1: Partner track found in the CTF collection using the electron's GSF track
46  // flag 2: Partner track found in the GSF collection using the electron's CTF track
47  // flag 3: Partner track found in the GSF collection using the electron's GSF track
48  };
49 
59 
60  std::vector<ConversionInfo> findConversions(const reco::GsfElectronCore& gsfElectron,
61  TrackTableView ctfTable,
62  TrackTableView gsfTable,
63  float bFieldAtOrigin,
64  float minFracSharedHits);
65 
66  //places different cuts on dist, dcot, delmissing hits and arbitration based on R = sqrt(dist*dist + dcot*dcot)
67  ConversionInfo findBestConversionMatch(const std::vector<ConversionInfo>& v_convCandidates);
68 
69  // returns the "best" conversion,
70  // bField has to be supplied in Tesla
72  TrackTableView ctfTable,
73  TrackTableView gsfTable,
74  float bFieldAtOrigin,
75  float minFracSharedHits = 0.45f) {
76  return findBestConversionMatch(findConversions(gsfElectron, ctfTable, gsfTable, bFieldAtOrigin, minFracSharedHits));
77  }
78 
79 } // namespace egamma::conv
80 
81 #endif
const std::optional< int > conversionPartnerGsfTkIdx
TrackTable::const_iterator::value_type TrackRowView
ConversionInfo findBestConversionMatch(const std::vector< ConversionInfo > &v_convCandidates)
ConversionInfo findConversion(const reco::GsfElectronCore &gsfElectron, TrackTableView ctfTable, TrackTableView gsfTable, float bFieldAtOrigin, float minFracSharedHits=0.45f)
typename ViewFromTable< T >::type ViewFromTable_t
Definition: TableView.h:104
edm::soa::ViewFromTable_t< TrackTable > TrackTableView
std::vector< ConversionInfo > findConversions(const reco::GsfElectronCore &gsfElectron, TrackTableView ctfTable, TrackTableView gsfTable, float bFieldAtOrigin, float minFracSharedHits)
const std::optional< int > conversionPartnerCtfTkIdx
EPOS::IO_EPOS conv
edm::soa::AddColumns< edm::soa::PtEtaPhiTable, TrackTableSpecificColumns >::type TrackTable
std::tuple< edm::soa::col::Pz, edm::soa::col::PtError, edm::soa::col::MissingInnerHits, edm::soa::col::NumberOfValidHits, edm::soa::col::Charge, edm::soa::col::D0 > TrackTableSpecificColumns