File indexing completed on 2024-04-06 12:33:18
0001 #include "FWCore/Framework/interface/global/EDProducer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Utilities/interface/InputTag.h"
0006
0007 #include "DataFormats/Common/interface/View.h"
0008 #include "DataFormats/Math/interface/deltaR.h"
0009 #include "DataFormats/MuonReco/interface/Muon.h"
0010 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0011 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0012 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0013 #include "DataFormats/JetReco/interface/Jet.h"
0014 #include "DataFormats/TauReco/interface/PFTau.h"
0015 #include "DataFormats/TauReco/interface/PFTauDiscriminator.h"
0016 #include "DataFormats/TrackReco/interface/Track.h"
0017
0018 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0019
0020 #include <algorithm>
0021 #include <memory>
0022 #include <vector>
0023
0024
0025
0026
0027 class ZllArbitrator : public edm::global::EDProducer<> {
0028 public:
0029 explicit ZllArbitrator(edm::ParameterSet const&);
0030 void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
0031
0032 private:
0033 edm::EDGetTokenT<std::vector<reco::CompositeCandidate>> srcZCand_;
0034 };
0035
0036
0037
0038
0039
0040 ZllArbitrator::ZllArbitrator(edm::ParameterSet const& iConfig)
0041 : srcZCand_{consumes<std::vector<reco::CompositeCandidate>>(
0042 iConfig.getParameter<edm::InputTag>("ZCandidateCollection"))} {
0043 produces<std::vector<reco::CompositeCandidate>>();
0044 }
0045
0046
0047
0048
0049
0050 void ZllArbitrator::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0051 edm::Handle<std::vector<reco::CompositeCandidate>> zCandidates;
0052 iEvent.getByToken(srcZCand_, zCandidates);
0053
0054 auto bestZ = std::make_unique<std::vector<reco::CompositeCandidate>>();
0055 if (!zCandidates->empty()) {
0056
0057 double constexpr ZmassPDG{91.18};
0058
0059 auto bestZCand = std::min_element(
0060 std::cbegin(*zCandidates), std::cend(*zCandidates), [ZmassPDG](auto const& firstCand, auto const& secondCand) {
0061 return std::abs(firstCand.mass() - ZmassPDG) < std::abs(secondCand.mass() - ZmassPDG);
0062 });
0063 bestZ->push_back(*bestZCand);
0064 }
0065
0066 iEvent.put(std::move(bestZ));
0067 }
0068
0069 using BestMassZArbitrationProducer = ZllArbitrator;
0070
0071 DEFINE_FWK_MODULE(BestMassZArbitrationProducer);