|
| 1 | +// -*- C++ -*- |
| 2 | +// |
| 3 | +// Package: PhysicsTools/Scouting |
| 4 | +// Class: Run3ScoutingElectronBestTrackProducer |
| 5 | +// |
| 6 | +/** |
| 7 | + Description: Choose the most suitable track for a given scouting electron |
| 8 | + Implementation: |
| 9 | + Allows for ID selections on the tracks before associating them to the electrons |
| 10 | +*/ |
| 11 | +// |
| 12 | +// Original Author: Abanti Ranadhir Sahasransu and Patin Inkaew |
| 13 | +// Created: Fri, 31 Jan 2025 14:43:20 GMT |
| 14 | +// |
| 15 | +// |
| 16 | + |
| 17 | +// system include files |
| 18 | +#include <limits> |
| 19 | +#include <memory> |
| 20 | +#include <vector> |
| 21 | + |
| 22 | +// user include files |
| 23 | +#include "FWCore/Framework/interface/Frameworkfwd.h" |
| 24 | +#include "FWCore/Framework/interface/stream/EDProducer.h" |
| 25 | +#include "FWCore/Framework/interface/Event.h" |
| 26 | +#include "FWCore/Framework/interface/MakerMacros.h" |
| 27 | +#include "FWCore/ParameterSet/interface/ParameterSet.h" |
| 28 | +#include "FWCore/Utilities/interface/StreamID.h" |
| 29 | + |
| 30 | +#include "DataFormats/Common/interface/Wrapper.h" |
| 31 | +#include "DataFormats/Common/interface/ValueMap.h" |
| 32 | +#include "DataFormats/Math/interface/LorentzVector.h" |
| 33 | +#include "DataFormats/Math/interface/deltaPhi.h" |
| 34 | +#include "DataFormats/Scouting/interface/Run3ScoutingElectron.h" |
| 35 | + |
| 36 | +// |
| 37 | +// class declaration |
| 38 | +// |
| 39 | + |
| 40 | +class Run3ScoutingElectronBestTrackProducer : public edm::stream::EDProducer<> { |
| 41 | +public: |
| 42 | + explicit Run3ScoutingElectronBestTrackProducer(const edm::ParameterSet&); |
| 43 | + static void fillDescriptions(edm::ConfigurationDescriptions&); |
| 44 | + |
| 45 | +private: |
| 46 | + void produce(edm::Event&, const edm::EventSetup&) override; |
| 47 | + |
| 48 | + template <typename T> |
| 49 | + void putValueMap(edm::Event&, edm::Handle<Run3ScoutingElectronCollection>&, const std::vector<T>&, const std::string&); |
| 50 | + |
| 51 | + const edm::EDGetTokenT<std::vector<Run3ScoutingElectron>> run3ScoutingElectronToken_; |
| 52 | + std::vector<double> trackPtMin_; |
| 53 | + std::vector<double> trackChi2OverNdofMax_; |
| 54 | + std::vector<double> relativeEnergyDifferenceMax_; |
| 55 | + std::vector<double> deltaPhiMax_; |
| 56 | +}; |
| 57 | + |
| 58 | +// |
| 59 | +// constructors and destructor |
| 60 | +// |
| 61 | +Run3ScoutingElectronBestTrackProducer::Run3ScoutingElectronBestTrackProducer(const edm::ParameterSet& iConfig) |
| 62 | + : run3ScoutingElectronToken_( |
| 63 | + consumes<std::vector<Run3ScoutingElectron>>(iConfig.getParameter<edm::InputTag>("Run3ScoutingElectron"))) { |
| 64 | + trackPtMin_ = iConfig.getParameter<std::vector<double>>("TrackPtMin"); |
| 65 | + trackChi2OverNdofMax_ = iConfig.getParameter<std::vector<double>>("TrackChi2OverNdofMax"); |
| 66 | + relativeEnergyDifferenceMax_ = iConfig.getParameter<std::vector<double>>("RelativeEnergyDifferenceMax"); |
| 67 | + deltaPhiMax_ = iConfig.getParameter<std::vector<double>>("DeltaPhiMax"); |
| 68 | + |
| 69 | + if (trackPtMin_.size() != 2) { |
| 70 | + throw cms::Exception("Run3ScoutingElectronBestTrackProducer") |
| 71 | + << "TrackPtMin must have exactly 2 elements for EB and EE respectively!"; |
| 72 | + } |
| 73 | + if (trackChi2OverNdofMax_.size() != 2) { |
| 74 | + throw cms::Exception("Run3ScoutingElectronBestTrackProducer") |
| 75 | + << "TrackChi2OverNdofMax must have exactly 2 elements for EB and EE respectively!"; |
| 76 | + } |
| 77 | + if (relativeEnergyDifferenceMax_.size() != 2) { |
| 78 | + throw cms::Exception("Run3ScoutingElectronBestTrackProducer") |
| 79 | + << "RelativeEnergyDifferenceMax must have exactly 2 elements for EB and EE respectively!"; |
| 80 | + } |
| 81 | + if (deltaPhiMax_.size() != 2) { |
| 82 | + throw cms::Exception("Run3ScoutingElectronBestTrackProducer") |
| 83 | + << "DeltaPhiMax must have exactly 2 elements for EB and EE respectively!"; |
| 84 | + } |
| 85 | + |
| 86 | + produces<edm::ValueMap<int>>("Run3ScoutingElectronBestTrackIndex"); |
| 87 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackd0"); |
| 88 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackdz"); |
| 89 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackpt"); |
| 90 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTracketa"); |
| 91 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackphi"); |
| 92 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackpMode"); |
| 93 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTracketaMode"); |
| 94 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackphiMode"); |
| 95 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackqoverpModeError"); |
| 96 | + produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackchi2overndf"); |
| 97 | + produces<edm::ValueMap<int>>("Run3ScoutingElectronTrackcharge"); |
| 98 | +} |
| 99 | + |
| 100 | +// |
| 101 | +// member functions |
| 102 | +// |
| 103 | + |
| 104 | +// ------------ method called to produce the data ------------ |
| 105 | +void Run3ScoutingElectronBestTrackProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) { |
| 106 | + edm::Handle<std::vector<Run3ScoutingElectron>> run3ScoutingElectronHandle; |
| 107 | + iEvent.getByToken(run3ScoutingElectronToken_, run3ScoutingElectronHandle); |
| 108 | + |
| 109 | + if (!run3ScoutingElectronHandle.isValid()) { |
| 110 | + // Handle the absence as a warning |
| 111 | + edm::LogWarning("Run3ScoutingElectronBestTrackProducer") |
| 112 | + << "No Run3ScoutingElectron collection found in the event!"; |
| 113 | + return; |
| 114 | + } |
| 115 | + |
| 116 | + const size_t num_electrons = run3ScoutingElectronHandle->size(); |
| 117 | + std::vector<int> besttrk_idx(num_electrons, -1); |
| 118 | + std::vector<float> besttrk_d0s(num_electrons, std::numeric_limits<float>::max()); |
| 119 | + std::vector<float> besttrk_dzs(num_electrons, std::numeric_limits<float>::max()); |
| 120 | + std::vector<float> besttrk_pts(num_electrons, std::numeric_limits<float>::max()); |
| 121 | + std::vector<float> besttrk_etas(num_electrons, std::numeric_limits<float>::max()); |
| 122 | + std::vector<float> besttrk_phis(num_electrons, std::numeric_limits<float>::max()); |
| 123 | + std::vector<float> besttrk_pModes(num_electrons, std::numeric_limits<float>::max()); |
| 124 | + std::vector<float> besttrk_etaModes(num_electrons, std::numeric_limits<float>::max()); |
| 125 | + std::vector<float> besttrk_phiModes(num_electrons, std::numeric_limits<float>::max()); |
| 126 | + std::vector<float> besttrk_qoverpModeErrors(num_electrons, std::numeric_limits<float>::max()); |
| 127 | + std::vector<float> besttrk_chi2overndfs(num_electrons, std::numeric_limits<float>::max()); |
| 128 | + std::vector<int> besttrk_charges(num_electrons, std::numeric_limits<int>::max()); |
| 129 | + |
| 130 | + for (size_t iElectron = 0; iElectron < num_electrons; ++iElectron) { |
| 131 | + const Run3ScoutingElectron& electron = run3ScoutingElectronHandle->at(iElectron); |
| 132 | + const math::PtEtaPhiMLorentzVector cluster(electron.pt(), electron.eta(), electron.phi(), 0.0005); |
| 133 | + |
| 134 | + double besttrack_ediff = std::numeric_limits<double>::max(); |
| 135 | + |
| 136 | + for (unsigned int i = 0; i < electron.trkpt().size(); ++i) { |
| 137 | + const unsigned int eta_idx = (std::abs(electron.trketa()[i]) < 1.479) ? 0 : 1; |
| 138 | + if (electron.trkpt()[i] < trackPtMin_[eta_idx]) |
| 139 | + continue; |
| 140 | + if (electron.trkchi2overndf()[i] > trackChi2OverNdofMax_[eta_idx]) |
| 141 | + continue; |
| 142 | + |
| 143 | + const math::PtEtaPhiMLorentzVector gsftrack( |
| 144 | + electron.trkpt()[i], electron.trketa()[i], electron.trkphi()[i], 0.0005); |
| 145 | + |
| 146 | + if (deltaPhi(cluster.phi(), gsftrack.phi()) > deltaPhiMax_[eta_idx]) |
| 147 | + continue; |
| 148 | + |
| 149 | + const double track_ediff = std::abs((cluster.energy() - gsftrack.energy()) / cluster.energy()); |
| 150 | + if (track_ediff > relativeEnergyDifferenceMax_[eta_idx]) |
| 151 | + continue; |
| 152 | + |
| 153 | + if (track_ediff < besttrack_ediff) { |
| 154 | + besttrack_ediff = track_ediff; |
| 155 | + besttrk_idx[iElectron] = i; |
| 156 | + } |
| 157 | + } |
| 158 | + |
| 159 | + if (besttrk_idx[iElectron] >= 0) { |
| 160 | + besttrk_d0s[iElectron] = electron.trkd0()[besttrk_idx[iElectron]]; |
| 161 | + besttrk_dzs[iElectron] = electron.trkdz()[besttrk_idx[iElectron]]; |
| 162 | + besttrk_pts[iElectron] = electron.trkpt()[besttrk_idx[iElectron]]; |
| 163 | + besttrk_etas[iElectron] = electron.trketa()[besttrk_idx[iElectron]]; |
| 164 | + besttrk_phis[iElectron] = electron.trkphi()[besttrk_idx[iElectron]]; |
| 165 | + if (!electron.trkpMode().empty()) { |
| 166 | + besttrk_pModes[iElectron] = electron.trkpMode()[besttrk_idx[iElectron]]; |
| 167 | + besttrk_etaModes[iElectron] = electron.trketaMode()[besttrk_idx[iElectron]]; |
| 168 | + besttrk_phiModes[iElectron] = electron.trkphiMode()[besttrk_idx[iElectron]]; |
| 169 | + besttrk_qoverpModeErrors[iElectron] = electron.trkqoverpModeError()[besttrk_idx[iElectron]]; |
| 170 | + } |
| 171 | + besttrk_chi2overndfs[iElectron] = electron.trkchi2overndf()[besttrk_idx[iElectron]]; |
| 172 | + besttrk_charges[iElectron] = electron.trkcharge()[besttrk_idx[iElectron]]; |
| 173 | + } |
| 174 | + } |
| 175 | + |
| 176 | + putValueMap<int>(iEvent, run3ScoutingElectronHandle, besttrk_idx, "Run3ScoutingElectronBestTrackIndex"); |
| 177 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_d0s, "Run3ScoutingElectronTrackd0"); |
| 178 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_dzs, "Run3ScoutingElectronTrackdz"); |
| 179 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_pts, "Run3ScoutingElectronTrackpt"); |
| 180 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_etas, "Run3ScoutingElectronTracketa"); |
| 181 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_phis, "Run3ScoutingElectronTrackphi"); |
| 182 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_pModes, "Run3ScoutingElectronTrackpMode"); |
| 183 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_etaModes, "Run3ScoutingElectronTracketaMode"); |
| 184 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_phiModes, "Run3ScoutingElectronTrackphiMode"); |
| 185 | + putValueMap<float>( |
| 186 | + iEvent, run3ScoutingElectronHandle, besttrk_qoverpModeErrors, "Run3ScoutingElectronTrackqoverpModeError"); |
| 187 | + putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_chi2overndfs, "Run3ScoutingElectronTrackchi2overndf"); |
| 188 | + putValueMap<int>(iEvent, run3ScoutingElectronHandle, besttrk_charges, "Run3ScoutingElectronTrackcharge"); |
| 189 | +} |
| 190 | + |
| 191 | +// ------------ method fills 'descriptions' with the allowed parameters for the module ------------ |
| 192 | +void Run3ScoutingElectronBestTrackProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) { |
| 193 | + edm::ParameterSetDescription desc; |
| 194 | + desc.add<edm::InputTag>(("Run3ScoutingElectron"), edm::InputTag("hltScoutingEgammaPacker")); |
| 195 | + desc.add<std::vector<double>>(("TrackPtMin"), {0.0, 0.0}); |
| 196 | + desc.add<std::vector<double>>(("TrackChi2OverNdofMax"), {9999.0, 9999.0}); |
| 197 | + desc.add<std::vector<double>>(("RelativeEnergyDifferenceMax"), {9999.0, 9999.0}); |
| 198 | + desc.add<std::vector<double>>(("DeltaPhiMax"), {9999.0, 9999.0}); |
| 199 | + descriptions.add("Run3ScoutingElectronBestTrackProducer", desc); |
| 200 | +} |
| 201 | + |
| 202 | +// ------------ method template for putting value maps into the event ------------ |
| 203 | +template <typename T> |
| 204 | +void Run3ScoutingElectronBestTrackProducer::putValueMap(edm::Event& iEvent, |
| 205 | + edm::Handle<Run3ScoutingElectronCollection>& handle, |
| 206 | + const std::vector<T>& values, |
| 207 | + const std::string& label) { |
| 208 | + std::unique_ptr<edm::ValueMap<T>> valuemap(new edm::ValueMap<T>()); |
| 209 | + typename edm::ValueMap<T>::Filler filler(*valuemap); |
| 210 | + filler.insert(handle, values.begin(), values.end()); |
| 211 | + filler.fill(); |
| 212 | + iEvent.put(std::move(valuemap), label); |
| 213 | +} |
| 214 | + |
| 215 | +//define this as a plug-in |
| 216 | +DEFINE_FWK_MODULE(Run3ScoutingElectronBestTrackProducer); |
0 commit comments