File indexing completed on 2024-06-22 02:24:08
0001
0002
0003
0004
0005
0006
0007
0008 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
0009 #include "DataFormats/TrackReco/interface/Track.h"
0010 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
0011 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "FWCore/Framework/interface/global/EDProducer.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0017
0018 class ElectronSeedMerger : public edm::global::EDProducer<> {
0019 public:
0020 explicit ElectronSeedMerger(const edm::ParameterSet&);
0021 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0022
0023 private:
0024 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0025
0026 const edm::EDGetTokenT<reco::ElectronSeedCollection> ecalSeedToken_;
0027 edm::EDGetTokenT<reco::ElectronSeedCollection> tkSeedToken_;
0028 };
0029
0030 using namespace edm;
0031 using namespace std;
0032 using namespace reco;
0033
0034 void ElectronSeedMerger::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0035 edm::ParameterSetDescription desc;
0036 desc.add<edm::InputTag>("EcalBasedSeeds", edm::InputTag("ecalDrivenElectronSeeds"));
0037 desc.add<edm::InputTag>("TkBasedSeeds", edm::InputTag("trackerDrivenElectronSeeds:SeedsForGsf"));
0038 descriptions.addWithDefaultLabel(desc);
0039 }
0040
0041 ElectronSeedMerger::ElectronSeedMerger(const ParameterSet& iConfig)
0042 : ecalSeedToken_{consumes<ElectronSeedCollection>(iConfig.getParameter<InputTag>("EcalBasedSeeds"))} {
0043 edm::InputTag tkSeedLabel_ = iConfig.getParameter<InputTag>("TkBasedSeeds");
0044 if (!tkSeedLabel_.label().empty())
0045 tkSeedToken_ = consumes<ElectronSeedCollection>(tkSeedLabel_);
0046
0047 produces<ElectronSeedCollection>();
0048 }
0049
0050
0051 void ElectronSeedMerger::produce(edm::StreamID, Event& iEvent, const EventSetup& iSetup) const {
0052
0053 auto const& eSeeds = iEvent.get(ecalSeedToken_);
0054
0055 ElectronSeedCollection tSeedsEmpty;
0056 auto const& tSeeds = tkSeedToken_.isUninitialized() ? tSeedsEmpty : iEvent.get(tkSeedToken_);
0057
0058
0059 auto output = std::make_unique<ElectronSeedCollection>();
0060 output->reserve(eSeeds.size() + tSeeds.size());
0061
0062
0063 vector<bool> tSeedsMatched(tSeeds.size(), false);
0064
0065
0066 for (auto newSeed : eSeeds) {
0067
0068
0069 int it = -1;
0070 for (auto const& tSeed : tSeeds) {
0071 it++;
0072
0073
0074 unsigned int hitShared = 0;
0075 unsigned int hitSeed = 0;
0076 for (auto const& eh : newSeed.recHits()) {
0077 if (!eh.isValid())
0078 continue;
0079 hitSeed++;
0080
0081 for (auto const& th : tSeed.recHits()) {
0082 if (!th.isValid())
0083 continue;
0084
0085 if (eh.sharesInput(&th, TrackingRecHit::all)) {
0086 hitShared++;
0087 break;
0088 }
0089 }
0090 }
0091 if (hitShared == hitSeed) {
0092 tSeedsMatched[it] = true;
0093 newSeed.setCtfTrack(tSeed.ctfTrack());
0094 break;
0095 }
0096 if (hitShared == (hitSeed - 1)) {
0097 newSeed.setCtfTrack(tSeed.ctfTrack());
0098 } else if ((hitShared > 0 || tSeed.nHits() == 0) && !newSeed.isTrackerDriven()) {
0099
0100 unsigned int hitSharedOnTrack = 0;
0101 for (auto const& eh : newSeed.recHits()) {
0102 if (!eh.isValid())
0103 continue;
0104 for (auto const* th : tSeed.ctfTrack()->recHits()) {
0105 if (!th->isValid())
0106 continue;
0107
0108 if (eh.sharesInput(th, TrackingRecHit::some)) {
0109 hitSharedOnTrack++;
0110 break;
0111 }
0112 }
0113 }
0114 if (hitSharedOnTrack == hitSeed) {
0115 tSeedsMatched[it] = true;
0116 newSeed.setCtfTrack(tSeed.ctfTrack());
0117 break;
0118 }
0119 if (hitSharedOnTrack == (hitSeed - 1)) {
0120 newSeed.setCtfTrack(tSeed.ctfTrack());
0121 }
0122 }
0123 }
0124
0125 output->push_back(newSeed);
0126 }
0127
0128
0129 for (unsigned int it = 0; it < tSeeds.size(); it++) {
0130 if (!tSeedsMatched[it])
0131 output->push_back(tSeeds[it]);
0132 }
0133
0134
0135 iEvent.put(std::move(output));
0136 }
0137
0138 DEFINE_FWK_MODULE(ElectronSeedMerger);