File indexing completed on 2024-04-06 12:23:59
0001 #include "FWCore/Framework/interface/MakerMacros.h"
0002
0003 #include "FWCore/Framework/interface/stream/EDProducer.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Utilities/interface/InputTag.h"
0007
0008 #include "DataFormats/PatCandidates/interface/EventHypothesis.h"
0009 #include "DataFormats/PatCandidates/interface/EventHypothesisLooper.h"
0010 #include "DataFormats/Common/interface/ValueMap.h"
0011 #include "DataFormats/Math/interface/deltaR.h"
0012 #include "DataFormats/MuonReco/interface/Muon.h"
0013 #include "DataFormats/TrackReco/interface/Track.h"
0014 #include "DataFormats/JetReco/interface/CaloJet.h"
0015
0016 class TestEventHypothesisWriter : public edm::stream::EDProducer<> {
0017 public:
0018 TestEventHypothesisWriter(const edm::ParameterSet &iConfig);
0019 void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0020 void runTests(const pat::EventHypothesis &h);
0021
0022 private:
0023 edm::EDGetTokenT<edm::View<reco::Candidate>> jetsToken_;
0024 edm::EDGetTokenT<edm::View<reco::Candidate>> muonsToken_;
0025 };
0026
0027 TestEventHypothesisWriter::TestEventHypothesisWriter(const edm::ParameterSet &iConfig)
0028 : jetsToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("jets"))),
0029 muonsToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("muons"))) {
0030 produces<std::vector<pat::EventHypothesis>>();
0031 produces<edm::ValueMap<double>>("deltaR");
0032 }
0033
0034 void TestEventHypothesisWriter::runTests(const pat::EventHypothesis &h) {
0035 using namespace std;
0036 using namespace pat::eventhypothesis;
0037 cout << "Test 1: Print the muon " << h["mu"]->pt() << endl;
0038
0039 for (size_t i = 0; i < h.count() - 2; ++i) {
0040 cout << "Test 2." << (i + 1) << ": Getting of the other jets: " << h.get("other jet", i)->et() << endl;
0041 }
0042
0043 cout << "Test 3: count: " << (h.count() - 2) << " vs " << h.count("other jet") << endl;
0044
0045 cout << "Test 4: regexp count: " << (h.count() - 1) << " vs " << h.count(".*jet") << endl;
0046
0047 cout << "Test 5.0: all with muon: " << h.all("mu").size() << endl;
0048 cout << "Test 5.1: all with muon: " << h.all("mu").front()->pt() << endl;
0049 cout << "Test 5.2: all with other jets: " << h.all("other jet").size() << endl;
0050 cout << "Test 5.3: all with regex: " << h.all(".*jet").size() << endl;
0051
0052 cout << "Test 6.0: get as : " << h.getAs<reco::CaloJet>("nearest jet")->maxEInHadTowers() << endl;
0053
0054 cout << "Loopers" << endl;
0055 cout << "Test 7.0: simple looper on all" << endl;
0056 for (CandLooper jet = h.loop(); jet; ++jet) {
0057 cout << "\titem " << jet.index() << ", role " << jet.role() << ": " << jet->et() << endl;
0058 }
0059 cout << "Test 7.1: simple looper on jets" << endl;
0060 for (CandLooper jet = h.loop(".*jet"); jet; ++jet) {
0061 cout << "\titem " << jet.index() << ", role " << jet.role() << ": " << jet->et() << endl;
0062 }
0063 cout << "Test 7.2: loopAs on jets" << endl;
0064 for (Looper<reco::CaloJet> jet = h.loopAs<reco::CaloJet>(".*jet"); jet; ++jet) {
0065 cout << "\titem " << jet.index() << ", role " << jet.role() << ": " << jet->maxEInHadTowers() << endl;
0066 }
0067 }
0068
0069 void TestEventHypothesisWriter::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0070 using namespace edm;
0071 using namespace std;
0072 using reco::Candidate;
0073 using reco::CandidatePtr;
0074
0075 auto hyps = std::make_unique<std::vector<pat::EventHypothesis>>();
0076 ;
0077 vector<double> deltaRs;
0078
0079 Handle<View<Candidate>> hMu;
0080 iEvent.getByToken(muonsToken_, hMu);
0081
0082 Handle<View<Candidate>> hJet;
0083 iEvent.getByToken(jetsToken_, hJet);
0084
0085
0086 for (size_t imu = 0, nmu = hMu->size(); imu < nmu; ++imu) {
0087 pat::EventHypothesis h;
0088 CandidatePtr mu = hMu->ptrAt(imu);
0089 h.add(mu, "mu");
0090
0091 int bestj = -1;
0092 double drmin = 99.0;
0093 for (size_t ij = 0, nj = hJet->size(); ij < nj; ++ij) {
0094 CandidatePtr jet = hJet->ptrAt(ij);
0095 if (jet->et() < 10)
0096 break;
0097 double dr = deltaR(*jet, *mu);
0098 if (dr < drmin) {
0099 bestj = ij;
0100 drmin = dr;
0101 }
0102 }
0103 if (bestj == -1)
0104 continue;
0105
0106 h.add(hJet->ptrAt(bestj), "nearest jet");
0107
0108 for (size_t ij = 0, nj = hJet->size(); ij < nj; ++ij) {
0109 if (ij == size_t(bestj))
0110 continue;
0111 CandidatePtr jet = hJet->ptrAt(ij);
0112 if (jet->et() < 10)
0113 break;
0114 h.add(jet, "other jet");
0115 }
0116
0117
0118 deltaRs.push_back(drmin);
0119
0120 runTests(h);
0121
0122 hyps->push_back(h);
0123 }
0124
0125 std::cout << "Found " << deltaRs.size() << " possible options" << std::endl;
0126
0127
0128 OrphanHandle<vector<pat::EventHypothesis>> handle = iEvent.put(std::move(hyps));
0129 auto deltaRMap = std::make_unique<ValueMap<double>>();
0130
0131 ValueMap<double>::Filler filler(*deltaRMap);
0132 filler.insert(handle, deltaRs.begin(), deltaRs.end());
0133 filler.fill();
0134
0135 iEvent.put(std::move(deltaRMap), "deltaR");
0136 }
0137
0138 DEFINE_FWK_MODULE(TestEventHypothesisWriter);