File indexing completed on 2025-04-17 02:42:25
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018 #include <limits>
0019 #include <memory>
0020 #include <vector>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/Utilities/interface/StreamID.h"
0029
0030 #include "DataFormats/Common/interface/Wrapper.h"
0031 #include "DataFormats/Common/interface/ValueMap.h"
0032 #include "DataFormats/Math/interface/LorentzVector.h"
0033 #include "DataFormats/Math/interface/deltaPhi.h"
0034 #include "DataFormats/Scouting/interface/Run3ScoutingElectron.h"
0035
0036
0037
0038
0039
0040 class Run3ScoutingElectronBestTrackProducer : public edm::stream::EDProducer<> {
0041 public:
0042 explicit Run3ScoutingElectronBestTrackProducer(const edm::ParameterSet&);
0043 static void fillDescriptions(edm::ConfigurationDescriptions&);
0044
0045 private:
0046 void produce(edm::Event&, const edm::EventSetup&) override;
0047
0048 template <typename T>
0049 void putValueMap(edm::Event&, edm::Handle<Run3ScoutingElectronCollection>&, const std::vector<T>&, const std::string&);
0050
0051 const edm::EDGetTokenT<std::vector<Run3ScoutingElectron>> run3ScoutingElectronToken_;
0052 std::vector<double> trackPtMin_;
0053 std::vector<double> trackChi2OverNdofMax_;
0054 std::vector<double> relativeEnergyDifferenceMax_;
0055 std::vector<double> deltaPhiMax_;
0056 };
0057
0058
0059
0060
0061 Run3ScoutingElectronBestTrackProducer::Run3ScoutingElectronBestTrackProducer(const edm::ParameterSet& iConfig)
0062 : run3ScoutingElectronToken_(
0063 consumes<std::vector<Run3ScoutingElectron>>(iConfig.getParameter<edm::InputTag>("Run3ScoutingElectron"))) {
0064 trackPtMin_ = iConfig.getParameter<std::vector<double>>("TrackPtMin");
0065 trackChi2OverNdofMax_ = iConfig.getParameter<std::vector<double>>("TrackChi2OverNdofMax");
0066 relativeEnergyDifferenceMax_ = iConfig.getParameter<std::vector<double>>("RelativeEnergyDifferenceMax");
0067 deltaPhiMax_ = iConfig.getParameter<std::vector<double>>("DeltaPhiMax");
0068
0069 if (trackPtMin_.size() != 2) {
0070 throw cms::Exception("Run3ScoutingElectronBestTrackProducer")
0071 << "TrackPtMin must have exactly 2 elements for EB and EE respectively!";
0072 }
0073 if (trackChi2OverNdofMax_.size() != 2) {
0074 throw cms::Exception("Run3ScoutingElectronBestTrackProducer")
0075 << "TrackChi2OverNdofMax must have exactly 2 elements for EB and EE respectively!";
0076 }
0077 if (relativeEnergyDifferenceMax_.size() != 2) {
0078 throw cms::Exception("Run3ScoutingElectronBestTrackProducer")
0079 << "RelativeEnergyDifferenceMax must have exactly 2 elements for EB and EE respectively!";
0080 }
0081 if (deltaPhiMax_.size() != 2) {
0082 throw cms::Exception("Run3ScoutingElectronBestTrackProducer")
0083 << "DeltaPhiMax must have exactly 2 elements for EB and EE respectively!";
0084 }
0085
0086 produces<edm::ValueMap<int>>("Run3ScoutingElectronBestTrackIndex");
0087 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackd0");
0088 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackdz");
0089 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackpt");
0090 produces<edm::ValueMap<float>>("Run3ScoutingElectronTracketa");
0091 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackphi");
0092 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackpMode");
0093 produces<edm::ValueMap<float>>("Run3ScoutingElectronTracketaMode");
0094 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackphiMode");
0095 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackqoverpModeError");
0096 produces<edm::ValueMap<float>>("Run3ScoutingElectronTrackchi2overndf");
0097 produces<edm::ValueMap<int>>("Run3ScoutingElectronTrackcharge");
0098 }
0099
0100
0101
0102
0103
0104
0105 void Run3ScoutingElectronBestTrackProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0106 edm::Handle<std::vector<Run3ScoutingElectron>> run3ScoutingElectronHandle;
0107 iEvent.getByToken(run3ScoutingElectronToken_, run3ScoutingElectronHandle);
0108
0109 if (!run3ScoutingElectronHandle.isValid()) {
0110
0111 edm::LogWarning("Run3ScoutingElectronBestTrackProducer")
0112 << "No Run3ScoutingElectron collection found in the event!";
0113 return;
0114 }
0115
0116 const size_t num_electrons = run3ScoutingElectronHandle->size();
0117 std::vector<int> besttrk_idx(num_electrons, -1);
0118 std::vector<float> besttrk_d0s(num_electrons, std::numeric_limits<float>::max());
0119 std::vector<float> besttrk_dzs(num_electrons, std::numeric_limits<float>::max());
0120 std::vector<float> besttrk_pts(num_electrons, std::numeric_limits<float>::max());
0121 std::vector<float> besttrk_etas(num_electrons, std::numeric_limits<float>::max());
0122 std::vector<float> besttrk_phis(num_electrons, std::numeric_limits<float>::max());
0123 std::vector<float> besttrk_pModes(num_electrons, std::numeric_limits<float>::max());
0124 std::vector<float> besttrk_etaModes(num_electrons, std::numeric_limits<float>::max());
0125 std::vector<float> besttrk_phiModes(num_electrons, std::numeric_limits<float>::max());
0126 std::vector<float> besttrk_qoverpModeErrors(num_electrons, std::numeric_limits<float>::max());
0127 std::vector<float> besttrk_chi2overndfs(num_electrons, std::numeric_limits<float>::max());
0128 std::vector<int> besttrk_charges(num_electrons, std::numeric_limits<int>::max());
0129
0130 for (size_t iElectron = 0; iElectron < num_electrons; ++iElectron) {
0131 const Run3ScoutingElectron& electron = run3ScoutingElectronHandle->at(iElectron);
0132 const math::PtEtaPhiMLorentzVector cluster(electron.pt(), electron.eta(), electron.phi(), 0.0005);
0133
0134 double besttrack_ediff = std::numeric_limits<double>::max();
0135
0136 for (unsigned int i = 0; i < electron.trkpt().size(); ++i) {
0137 const unsigned int eta_idx = (std::abs(electron.trketa()[i]) < 1.479) ? 0 : 1;
0138 if (electron.trkpt()[i] < trackPtMin_[eta_idx])
0139 continue;
0140 if (electron.trkchi2overndf()[i] > trackChi2OverNdofMax_[eta_idx])
0141 continue;
0142
0143 const math::PtEtaPhiMLorentzVector gsftrack(
0144 electron.trkpt()[i], electron.trketa()[i], electron.trkphi()[i], 0.0005);
0145
0146 if (deltaPhi(cluster.phi(), gsftrack.phi()) > deltaPhiMax_[eta_idx])
0147 continue;
0148
0149 const double track_ediff = std::abs((cluster.energy() - gsftrack.energy()) / cluster.energy());
0150 if (track_ediff > relativeEnergyDifferenceMax_[eta_idx])
0151 continue;
0152
0153 if (track_ediff < besttrack_ediff) {
0154 besttrack_ediff = track_ediff;
0155 besttrk_idx[iElectron] = i;
0156 }
0157 }
0158
0159 if (besttrk_idx[iElectron] >= 0) {
0160 besttrk_d0s[iElectron] = electron.trkd0()[besttrk_idx[iElectron]];
0161 besttrk_dzs[iElectron] = electron.trkdz()[besttrk_idx[iElectron]];
0162 besttrk_pts[iElectron] = electron.trkpt()[besttrk_idx[iElectron]];
0163 besttrk_etas[iElectron] = electron.trketa()[besttrk_idx[iElectron]];
0164 besttrk_phis[iElectron] = electron.trkphi()[besttrk_idx[iElectron]];
0165 if (!electron.trkpMode().empty()) {
0166 besttrk_pModes[iElectron] = electron.trkpMode()[besttrk_idx[iElectron]];
0167 besttrk_etaModes[iElectron] = electron.trketaMode()[besttrk_idx[iElectron]];
0168 besttrk_phiModes[iElectron] = electron.trkphiMode()[besttrk_idx[iElectron]];
0169 besttrk_qoverpModeErrors[iElectron] = electron.trkqoverpModeError()[besttrk_idx[iElectron]];
0170 }
0171 besttrk_chi2overndfs[iElectron] = electron.trkchi2overndf()[besttrk_idx[iElectron]];
0172 besttrk_charges[iElectron] = electron.trkcharge()[besttrk_idx[iElectron]];
0173 }
0174 }
0175
0176 putValueMap<int>(iEvent, run3ScoutingElectronHandle, besttrk_idx, "Run3ScoutingElectronBestTrackIndex");
0177 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_d0s, "Run3ScoutingElectronTrackd0");
0178 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_dzs, "Run3ScoutingElectronTrackdz");
0179 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_pts, "Run3ScoutingElectronTrackpt");
0180 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_etas, "Run3ScoutingElectronTracketa");
0181 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_phis, "Run3ScoutingElectronTrackphi");
0182 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_pModes, "Run3ScoutingElectronTrackpMode");
0183 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_etaModes, "Run3ScoutingElectronTracketaMode");
0184 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_phiModes, "Run3ScoutingElectronTrackphiMode");
0185 putValueMap<float>(
0186 iEvent, run3ScoutingElectronHandle, besttrk_qoverpModeErrors, "Run3ScoutingElectronTrackqoverpModeError");
0187 putValueMap<float>(iEvent, run3ScoutingElectronHandle, besttrk_chi2overndfs, "Run3ScoutingElectronTrackchi2overndf");
0188 putValueMap<int>(iEvent, run3ScoutingElectronHandle, besttrk_charges, "Run3ScoutingElectronTrackcharge");
0189 }
0190
0191
0192 void Run3ScoutingElectronBestTrackProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0193 edm::ParameterSetDescription desc;
0194 desc.add<edm::InputTag>(("Run3ScoutingElectron"), edm::InputTag("hltScoutingEgammaPacker"));
0195 desc.add<std::vector<double>>(("TrackPtMin"), {0.0, 0.0});
0196 desc.add<std::vector<double>>(("TrackChi2OverNdofMax"), {9999.0, 9999.0});
0197 desc.add<std::vector<double>>(("RelativeEnergyDifferenceMax"), {9999.0, 9999.0});
0198 desc.add<std::vector<double>>(("DeltaPhiMax"), {9999.0, 9999.0});
0199 descriptions.add("Run3ScoutingElectronBestTrackProducer", desc);
0200 }
0201
0202
0203 template <typename T>
0204 void Run3ScoutingElectronBestTrackProducer::putValueMap(edm::Event& iEvent,
0205 edm::Handle<Run3ScoutingElectronCollection>& handle,
0206 const std::vector<T>& values,
0207 const std::string& label) {
0208 std::unique_ptr<edm::ValueMap<T>> valuemap(new edm::ValueMap<T>());
0209 typename edm::ValueMap<T>::Filler filler(*valuemap);
0210 filler.insert(handle, values.begin(), values.end());
0211 filler.fill();
0212 iEvent.put(std::move(valuemap), label);
0213 }
0214
0215
0216 DEFINE_FWK_MODULE(Run3ScoutingElectronBestTrackProducer);