File indexing completed on 2024-04-06 12:21:14
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include <memory>
0015
0016
0017 #include "FWCore/Framework/interface/Frameworkfwd.h"
0018 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 #include "FWCore/Framework/interface/MakerMacros.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "FWCore/Framework/interface/ESHandle.h"
0025
0026
0027 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0028
0029
0030 #include "FWCore/ServiceRegistry/interface/Service.h"
0031 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0032 #include "TH1.h"
0033 #include "TTree.h"
0034 #include "TF1.h"
0035
0036
0037 #include "L1Trigger/L1TNtuples/interface/L1AnalysisRecoPhoton.h"
0038
0039
0040
0041
0042
0043 class L1PhotonRecoTreeProducer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0044 public:
0045 explicit L1PhotonRecoTreeProducer(const edm::ParameterSet&);
0046 ~L1PhotonRecoTreeProducer() override;
0047
0048 private:
0049 void beginJob(void) override;
0050 void analyze(const edm::Event&, const edm::EventSetup&) override;
0051 void endJob() override;
0052
0053 public:
0054 L1Analysis::L1AnalysisRecoPhoton* photon;
0055
0056 L1Analysis::L1AnalysisRecoPhotonDataFormat* photon_data;
0057
0058 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0059
0060 private:
0061
0062 edm::Service<TFileService> fs_;
0063
0064
0065 TTree* tree_;
0066
0067
0068
0069 edm::EDGetTokenT<reco::PhotonCollection> PhotonToken_;
0070 edm::EDGetTokenT<edm::ValueMap<bool> > PhotonWP80MapToken_;
0071 edm::EDGetTokenT<edm::ValueMap<bool> > PhotonWP90MapToken_;
0072
0073
0074 bool photonsMissing_;
0075 unsigned int maxPhoton_;
0076 };
0077
0078 L1PhotonRecoTreeProducer::L1PhotonRecoTreeProducer(const edm::ParameterSet& iConfig) : photonsMissing_(false) {
0079 maxPhoton_ = iConfig.getParameter<unsigned int>("maxPhoton");
0080 PhotonToken_ =
0081 consumes<reco::PhotonCollection>(iConfig.getUntrackedParameter("PhotonToken", edm::InputTag("photons")));
0082
0083 PhotonWP80MapToken_ = consumes<edm::ValueMap<bool> >(
0084 iConfig.getUntrackedParameter("phoWP80MapToken", edm::InputTag("egmPhotonIDs:mvaPhoID-RunIIIWinter22-v1-wp80")));
0085 PhotonWP90MapToken_ = consumes<edm::ValueMap<bool> >(
0086 iConfig.getUntrackedParameter("phoWP90MapToken", edm::InputTag("egmPhotonIDs:mvaPhoID-RunIIIWinter22-v1-wp90")));
0087
0088 photon = new L1Analysis::L1AnalysisRecoPhoton();
0089 photon_data = photon->getData();
0090
0091 usesResource(TFileService::kSharedResource);
0092 tree_ = fs_->make<TTree>("PhotonRecoTree", "PhotonRecoTree");
0093 tree_->Branch("Photon", "L1Analysis::L1AnalysisRecoPhotonDataFormat", &photon_data, 32000, 3);
0094 }
0095
0096 L1PhotonRecoTreeProducer::~L1PhotonRecoTreeProducer() {
0097
0098
0099 }
0100
0101
0102
0103
0104
0105
0106 void L1PhotonRecoTreeProducer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0107 photon->Reset();
0108 edm::Handle<reco::PhotonCollection> recoPhotons;
0109 iEvent.getByToken(PhotonToken_, recoPhotons);
0110
0111 std::vector<edm::Handle<edm::ValueMap<bool> > > phoVIDDecisionHandles(2);
0112
0113 iEvent.getByToken(PhotonWP80MapToken_, phoVIDDecisionHandles[0]);
0114 iEvent.getByToken(PhotonWP90MapToken_, phoVIDDecisionHandles[1]);
0115
0116 if (recoPhotons.isValid() && phoVIDDecisionHandles[0].isValid() && phoVIDDecisionHandles[1].isValid()) {
0117 photon->SetPhoton(iEvent, iSetup, recoPhotons, phoVIDDecisionHandles, maxPhoton_);
0118 } else {
0119 if (!photonsMissing_) {
0120 edm::LogWarning("MissingProduct") << "Photons or photon ID not found. Branch will not be filled" << std::endl;
0121 }
0122 photonsMissing_ = true;
0123 }
0124
0125 tree_->Fill();
0126 }
0127
0128
0129 void L1PhotonRecoTreeProducer::beginJob(void) {}
0130
0131
0132 void L1PhotonRecoTreeProducer::endJob() {}
0133
0134 void L1PhotonRecoTreeProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0135 edm::ParameterSetDescription desc;
0136 desc.add<unsigned int>("maxPhoton", 20);
0137 descriptions.addWithDefaultLabel(desc);
0138 }
0139
0140
0141 DEFINE_FWK_MODULE(L1PhotonRecoTreeProducer);