File indexing completed on 2024-04-06 12:24:48
0001 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0002 #include "DataFormats/TrackReco/interface/Track.h"
0003 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "FWCore/Utilities/interface/EDGetToken.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "RecoEgamma/EgammaElectronAlgos/interface/ConversionFinder.h"
0015
0016 namespace egamma::conv {
0017
0018 std::ostream& operator<<(std::ostream& os, ConversionInfo const& conv) {
0019 auto partnerIdx = conv.conversionPartnerCtfTkIdx.value_or(conv.conversionPartnerGsfTkIdx.value_or(-9999));
0020
0021 return os << " - flag: " << conv.flag << "\n - partnerIdx:" << partnerIdx << "\n - dist: " << conv.dist
0022 << "\n - dcot: " << conv.dcot << "\n - radius: " << conv.radiusOfConversion
0023 << "\n - deltaMissingHits: " << conv.deltaMissingHits;
0024 }
0025
0026 }
0027
0028 class TestGsfElectronConversionFinder : public edm::one::EDAnalyzer<> {
0029 public:
0030 explicit TestGsfElectronConversionFinder(const edm::ParameterSet&);
0031
0032 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0033
0034 private:
0035 void analyze(const edm::Event&, const edm::EventSetup&) override;
0036
0037 const edm::EDGetTokenT<edm::View<reco::GsfElectron>> gsfElectronsToken_;
0038 };
0039
0040 TestGsfElectronConversionFinder::TestGsfElectronConversionFinder(const edm::ParameterSet& pset)
0041 : gsfElectronsToken_{consumes<edm::View<reco::GsfElectron>>(pset.getParameter<edm::InputTag>("gsfElectrons"))} {}
0042
0043 void TestGsfElectronConversionFinder::analyze(const edm::Event& event, const edm::EventSetup&) {
0044
0045 constexpr float magneticFieldInTesla = 4.0f;
0046
0047
0048 auto const& gsfElectrons = event.get(gsfElectronsToken_);
0049
0050
0051 std::optional<egamma::conv::TrackTable> ctfTrackTable = std::nullopt;
0052 std::optional<egamma::conv::TrackTable> gsfTrackTable = std::nullopt;
0053
0054 int gsfElectronIdx{0};
0055
0056
0057 for (auto const& gsfElectron : gsfElectrons) {
0058 auto const& gsfTrack = gsfElectron.core()->gsfTrack();
0059 auto const& ctfTrack = gsfElectron.core()->ctfTrack();
0060 edm::LogVerbatim("TestConversionFinder") << "===================================" << std::endl;
0061 edm::LogVerbatim("TestConversionFinder")
0062 << "event " << event.id().event() << ", gsfElectron " << gsfElectronIdx << std::endl;
0063 edm::LogVerbatim("TestConversionFinder") << "-------------------------" << std::endl;
0064 if (gsfTrackTable == std::nullopt) {
0065 edm::Handle<reco::GsfTrackCollection> originalGsfTracksHandle;
0066 event.get(gsfTrack.id(), originalGsfTracksHandle);
0067 gsfTrackTable = egamma::conv::TrackTable(*originalGsfTracksHandle);
0068 }
0069 if (ctfTrackTable == std::nullopt) {
0070 if (!ctfTrack.isNull()) {
0071 edm::Handle<reco::TrackCollection> originalCtfTracksHandle;
0072 event.get(ctfTrack.id(), originalCtfTracksHandle);
0073 ctfTrackTable = egamma::conv::TrackTable(*originalCtfTracksHandle);
0074 }
0075 }
0076
0077
0078
0079
0080 auto conversions = egamma::conv::findConversions(*gsfElectron.core(),
0081 ctfTrackTable.value_or(egamma::conv::TrackTable{}),
0082 gsfTrackTable.value(),
0083 magneticFieldInTesla,
0084 0.45f);
0085 auto bestConversion = egamma::conv::findBestConversionMatch(conversions);
0086
0087
0088 int conversionIdx{0};
0089 for (auto const& conversion : conversions) {
0090 edm::LogVerbatim("TestConversionFinder") << "conversion " << conversionIdx << std::endl;
0091 edm::LogVerbatim("TestConversionFinder") << conversion << std::endl;
0092 edm::LogVerbatim("TestConversionFinder") << "-------------------------" << std::endl;
0093 ++conversionIdx;
0094 }
0095 edm::LogVerbatim("TestConversionFinder") << "best conversion " << std::endl;
0096 edm::LogVerbatim("TestConversionFinder") << bestConversion << std::endl;
0097
0098 ++gsfElectronIdx;
0099 }
0100 }
0101
0102 void TestGsfElectronConversionFinder::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0103 edm::ParameterSetDescription desc;
0104 desc.add<edm::InputTag>("gsfElectrons", {"gedGsfElectrons"});
0105 descriptions.addDefault(desc);
0106 }
0107
0108 #include "FWCore/Framework/interface/MakerMacros.h"
0109 DEFINE_FWK_MODULE(TestGsfElectronConversionFinder);