Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:46

0001 #ifndef RecoEgamma_EgammaElectronAlgos_ConversionFinder_h
0002 #define RecoEgamma_EgammaElectronAlgos_ConversionFinder_h
0003 
0004 /** \class reco:: ConversionFinder.h RecoEgamma/EgammaElectronAlgos/interface/ConversionFinder.h
0005   *
0006   * Conversion finding and rejection code
0007   * Uses simple geometric methods to determine whether or not the
0008   * electron did indeed come from a conversion
0009   * \author Puneeth Kalavase, University Of California, Santa Barbara
0010   *
0011   *
0012   */
0013 
0014 #include "DataFormats/EgammaCandidates/interface/GsfElectronCore.h"
0015 #include "CommonTools/Utils/interface/KinematicTables.h"
0016 #include "CommonTools/Utils/interface/TrackSpecificColumns.h"
0017 
0018 #include <iostream>
0019 #include <optional>
0020 
0021 /*
0022    Class Looks for oppositely charged track in the
0023    track collection with the minimum delta cot(theta) between the track
0024    and the electron's CTF track (if it doesn't exist, we use the electron's
0025    GSF track). Calculate the dist, dcot, point of conversion and the
0026    radius of conversion for this pair and fill the ConversionInfo
0027 */
0028 
0029 namespace egamma::conv {
0030 
0031   struct ConversionInfo {
0032     const float dist = -9999.;
0033     const float dcot = -9999.;
0034     const float radiusOfConversion = -9999.;
0035     // if the partner track is found in the  GSF track collection,
0036     // this is a ref to the GSF partner track
0037     const std::optional<int> conversionPartnerCtfTkIdx = std::nullopt;
0038     // if the partner track is found in the  CTF track collection,
0039     // this is a ref to the CTF partner track
0040     const std::optional<int> conversionPartnerGsfTkIdx = std::nullopt;
0041     const int deltaMissingHits = -9999;
0042     const int flag = -9999;
0043 
0044     // flag 0: Partner track found in the CTF collection using the electron's CTF track
0045     // flag 1: Partner track found in the CTF collection using the electron's GSF track
0046     // flag 2: Partner track found in the GSF collection using the electron's CTF track
0047     // flag 3: Partner track found in the GSF collection using the electron's GSF track
0048   };
0049 
0050   using TrackTableSpecificColumns = std::tuple<edm::soa::col::Pz,
0051                                                edm::soa::col::PtError,
0052                                                edm::soa::col::MissingInnerHits,
0053                                                edm::soa::col::NumberOfValidHits,
0054                                                edm::soa::col::Charge,
0055                                                edm::soa::col::D0>;
0056   using TrackTable = edm::soa::AddColumns<edm::soa::PtEtaPhiTable, TrackTableSpecificColumns>::type;
0057   using TrackTableView = edm::soa::ViewFromTable_t<TrackTable>;
0058   using TrackRowView = TrackTable::const_iterator::value_type;
0059 
0060   std::vector<ConversionInfo> findConversions(const reco::GsfElectronCore& gsfElectron,
0061                                               TrackTableView ctfTable,
0062                                               TrackTableView gsfTable,
0063                                               float bFieldAtOrigin,
0064                                               float minFracSharedHits);
0065 
0066   //places different cuts on dist, dcot, delmissing hits and arbitration based on R = sqrt(dist*dist + dcot*dcot)
0067   ConversionInfo findBestConversionMatch(const std::vector<ConversionInfo>& v_convCandidates);
0068 
0069   // returns the "best" conversion,
0070   // bField has to be supplied in Tesla
0071   inline ConversionInfo findConversion(const reco::GsfElectronCore& gsfElectron,
0072                                        TrackTableView ctfTable,
0073                                        TrackTableView gsfTable,
0074                                        float bFieldAtOrigin,
0075                                        float minFracSharedHits = 0.45f) {
0076     return findBestConversionMatch(findConversions(gsfElectron, ctfTable, gsfTable, bFieldAtOrigin, minFracSharedHits));
0077   }
0078 
0079 }  // namespace egamma::conv
0080 
0081 #endif