File indexing completed on 2024-04-06 11:58:43
0001
0002 #include <memory>
0003 #include <string>
0004
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/global/EDProducer.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013
0014 #include "DataFormats/DetId/interface/DetId.h"
0015 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0016 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0017 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0018 #include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
0019 #include "DataFormats/HcalRecHit/interface/HFRecHit.h"
0020 #include "DataFormats/HcalRecHit/interface/HORecHit.h"
0021 #include "DataFormats/Common/interface/ValueMap.h"
0022 #include "DataFormats/Common/interface/TriggerResults.h"
0023 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0024 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0025 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0026 #include "DataFormats/METReco/interface/PFMET.h"
0027 #include "DataFormats/METReco/interface/PFMETCollection.h"
0028 #include "DataFormats/VertexReco/interface/Vertex.h"
0029 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0030 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0031
0032 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0033 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0034 #include <iostream>
0035
0036 namespace alCaHcalDiJetsProducer {
0037 struct Counters {
0038 Counters() : nAll_(0), nSelect_(0) {}
0039 mutable std::atomic<unsigned int> nAll_, nSelect_;
0040 };
0041 }
0042
0043
0044
0045
0046
0047 class AlCaDiJetsProducer : public edm::global::EDProducer<edm::RunCache<alCaHcalDiJetsProducer::Counters>> {
0048 public:
0049 explicit AlCaDiJetsProducer(const edm::ParameterSet&);
0050 ~AlCaDiJetsProducer() override = default;
0051
0052 std::shared_ptr<alCaHcalDiJetsProducer::Counters> globalBeginRun(edm::Run const&,
0053 edm::EventSetup const&) const override {
0054 return std::make_shared<alCaHcalDiJetsProducer::Counters>();
0055 }
0056 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0057 void globalEndRun(edm::Run const& iRun, edm::EventSetup const&) const override;
0058 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0059
0060 private:
0061 bool select(reco::PFJetCollection) const;
0062
0063
0064
0065 edm::InputTag labelPFJet_, labelHBHE_, labelHF_, labelHO_, labelPFCandidate_, labelVertex_;
0066 double minPtJet_;
0067
0068 edm::EDGetTokenT<reco::PFJetCollection> tok_PFJet_;
0069 edm::EDGetTokenT<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>> tok_HBHE_;
0070 edm::EDGetTokenT<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>> tok_HF_;
0071 edm::EDGetTokenT<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>> tok_HO_;
0072
0073 edm::EDGetTokenT<reco::PFCandidateCollection> tok_PFCand_;
0074 edm::EDGetTokenT<reco::VertexCollection> tok_Vertex_;
0075 };
0076
0077 AlCaDiJetsProducer::AlCaDiJetsProducer(const edm::ParameterSet& iConfig) {
0078
0079 labelPFJet_ = iConfig.getParameter<edm::InputTag>("PFjetInput");
0080 labelHBHE_ = iConfig.getParameter<edm::InputTag>("HBHEInput");
0081 labelHF_ = iConfig.getParameter<edm::InputTag>("HFInput");
0082 labelHO_ = iConfig.getParameter<edm::InputTag>("HOInput");
0083
0084 labelPFCandidate_ = iConfig.getParameter<edm::InputTag>("particleFlowInput");
0085 labelVertex_ = iConfig.getParameter<edm::InputTag>("VertexInput");
0086 minPtJet_ = iConfig.getParameter<double>("MinPtJet");
0087
0088 tok_PFJet_ = consumes<reco::PFJetCollection>(labelPFJet_);
0089 tok_HBHE_ = consumes<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>>(labelHBHE_);
0090 tok_HF_ = consumes<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>>(labelHF_);
0091 tok_HO_ = consumes<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>>(labelHO_);
0092
0093 tok_PFCand_ = consumes<reco::PFCandidateCollection>(labelPFCandidate_);
0094 tok_Vertex_ = consumes<reco::VertexCollection>(labelVertex_);
0095
0096
0097 produces<reco::PFJetCollection>(labelPFJet_.encode());
0098 produces<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>>(labelHBHE_.encode());
0099 produces<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>>(labelHF_.encode());
0100 produces<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>>(labelHO_.encode());
0101
0102 produces<reco::PFCandidateCollection>(labelPFCandidate_.encode());
0103 produces<reco::VertexCollection>(labelVertex_.encode());
0104 }
0105
0106 void AlCaDiJetsProducer::globalEndRun(edm::Run const& iRun, edm::EventSetup const&) const {
0107 edm::LogVerbatim("AlcaDiJets") << "Accepts " << runCache(iRun.index())->nSelect_ << " events from a total of "
0108 << runCache(iRun.index())->nAll_ << " events";
0109 }
0110
0111 bool AlCaDiJetsProducer::select(reco::PFJetCollection jt) const {
0112 if (jt.size() < 2)
0113 return false;
0114 if (((jt.at(0)).pt()) < minPtJet_)
0115 return false;
0116 return true;
0117 }
0118
0119
0120 void AlCaDiJetsProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0121 ++(runCache(iEvent.getRun().index())->nAll_);
0122
0123
0124 edm::Handle<reco::PFJetCollection> pfjet = iEvent.getHandle(tok_PFJet_);
0125 if (!pfjet.isValid()) {
0126 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelPFJet_;
0127 return;
0128 }
0129 const reco::PFJetCollection pfjets = *(pfjet.product());
0130
0131 edm::Handle<reco::PFCandidateCollection> pfc = iEvent.getHandle(tok_PFCand_);
0132 if (!pfc.isValid()) {
0133 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelPFCandidate_;
0134 return;
0135 }
0136 const reco::PFCandidateCollection pfcand = *(pfc.product());
0137
0138 edm::Handle<reco::VertexCollection> vt = iEvent.getHandle(tok_Vertex_);
0139 if (!vt.isValid()) {
0140 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelVertex_;
0141 return;
0142 }
0143 const reco::VertexCollection vtx = *(vt.product());
0144
0145 auto hbhe = iEvent.getHandle(tok_HBHE_);
0146 if (!hbhe.isValid()) {
0147 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHBHE_;
0148 return;
0149 }
0150 const edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>> Hithbhe = *(hbhe.product());
0151
0152 auto ho = iEvent.getHandle(tok_HO_);
0153 if (!ho.isValid()) {
0154 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHO_;
0155 return;
0156 }
0157 const edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>> Hitho = *(ho.product());
0158
0159 auto hf = iEvent.getHandle(tok_HF_);
0160 if (!hf.isValid()) {
0161 edm::LogWarning("AlCaDiJets") << "AlCaDiJetsProducer: Error! can't get product " << labelHF_;
0162 return;
0163 }
0164 const edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>> Hithf = *(hf.product());
0165
0166
0167 bool accept = select(pfjets);
0168 if (accept) {
0169 ++(runCache(iEvent.getRun().index())->nSelect_);
0170
0171
0172 auto miniPFjetCollection = std::make_unique<reco::PFJetCollection>();
0173 for (reco::PFJetCollection::const_iterator pfjetItr = pfjets.begin(); pfjetItr != pfjets.end(); pfjetItr++) {
0174 miniPFjetCollection->push_back(*pfjetItr);
0175 }
0176
0177 auto miniPFCandCollection = std::make_unique<reco::PFCandidateCollection>();
0178 for (reco::PFCandidateCollection::const_iterator pfcItr = pfcand.begin(); pfcItr != pfcand.end(); pfcItr++) {
0179 miniPFCandCollection->push_back(*pfcItr);
0180 }
0181
0182 auto miniVtxCollection = std::make_unique<reco::VertexCollection>();
0183 for (reco::VertexCollection::const_iterator vtxItr = vtx.begin(); vtxItr != vtx.end(); vtxItr++) {
0184 miniVtxCollection->push_back(*vtxItr);
0185 }
0186
0187 auto miniHBHECollection =
0188 std::make_unique<edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>>();
0189 for (edm::SortedCollection<HBHERecHit, edm::StrictWeakOrdering<HBHERecHit>>::const_iterator hbheItr =
0190 Hithbhe.begin();
0191 hbheItr != Hithbhe.end();
0192 hbheItr++) {
0193 miniHBHECollection->push_back(*hbheItr);
0194 }
0195
0196 auto miniHOCollection = std::make_unique<edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>>();
0197 for (edm::SortedCollection<HORecHit, edm::StrictWeakOrdering<HORecHit>>::const_iterator hoItr = Hitho.begin();
0198 hoItr != Hitho.end();
0199 hoItr++) {
0200 miniHOCollection->push_back(*hoItr);
0201 }
0202
0203 auto miniHFCollection = std::make_unique<edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>>();
0204 for (edm::SortedCollection<HFRecHit, edm::StrictWeakOrdering<HFRecHit>>::const_iterator hfItr = Hithf.begin();
0205 hfItr != Hithf.end();
0206 hfItr++) {
0207 miniHFCollection->push_back(*hfItr);
0208 }
0209
0210
0211 iEvent.put(std::move(miniPFjetCollection), labelPFJet_.encode());
0212 iEvent.put(std::move(miniHBHECollection), labelHBHE_.encode());
0213 iEvent.put(std::move(miniHFCollection), labelHF_.encode());
0214 iEvent.put(std::move(miniHOCollection), labelHO_.encode());
0215
0216 iEvent.put(std::move(miniPFCandCollection), labelPFCandidate_.encode());
0217 iEvent.put(std::move(miniVtxCollection), labelVertex_.encode());
0218 }
0219 return;
0220 }
0221
0222 void AlCaDiJetsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0223 edm::ParameterSetDescription desc;
0224 desc.add<edm::InputTag>("PFjetInput", edm::InputTag("ak4PFJetsCHS")),
0225 desc.add<edm::InputTag>("HBHEInput", edm::InputTag("hbhereco"));
0226 desc.add<edm::InputTag>("HFInput", edm::InputTag("hfreco"));
0227 desc.add<edm::InputTag>("HOInput", edm::InputTag("horeco"));
0228 desc.add<edm::InputTag>("particleFlowInput", edm::InputTag("particleFlow"));
0229 desc.add<edm::InputTag>("VertexInput", edm::InputTag("offlinePrimaryVertices"));
0230 desc.add<double>("MinPtJet", 20.0);
0231 descriptions.add("alcaDiJetsProducer", desc);
0232 }
0233
0234 DEFINE_FWK_MODULE(AlCaDiJetsProducer);