Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:14

0001 // Producer for particle flow candidates. Plots Eta, Phi, Charge, Pt (log freq, bin)
0002 // for different types of particles described in python/defaults_cfi.py
0003 // It actually uses packedCandidates so that we need only MINIAOD contents to run this DQMAnalyzer.
0004 // note: for pt, log freq is done in this producer, but log freq is done by running
0005 // compare.py
0006 // author: Chosila Sutantawibul, April 23, 2020
0007 
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ServiceRegistry/interface/Service.h"
0011 
0012 #include "DataFormats/Common/interface/Handle.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/Framework/interface/MakerMacros.h"
0020 #include "FWCore/Utilities/interface/InputTag.h"
0021 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0022 
0023 #include "TH1F.h"
0024 
0025 #include <algorithm>
0026 #include <cmath>
0027 #include <fstream>
0028 #include <iostream>
0029 #include <memory>
0030 #include <ostream>
0031 #include <map>
0032 #include <string>
0033 #include <cstring>
0034 
0035 class PFCandidateAnalyzerDQM : public DQMEDAnalyzer {
0036 public:
0037   explicit PFCandidateAnalyzerDQM(const edm::ParameterSet&);
0038   void analyze(const edm::Event&, const edm::EventSetup&) override;
0039 
0040 protected:
0041   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0042 
0043 private:
0044   //from config file
0045   edm::InputTag PFCandTag;
0046   edm::EDGetTokenT<edm::View<pat::PackedCandidate>> PFCandToken;
0047   std::vector<double> etabins;
0048   std::map<std::string, MonitorElement*> me;
0049 
0050   std::map<uint32_t, std::string> pdgMap;
0051 };
0052 
0053 // constructor
0054 PFCandidateAnalyzerDQM::PFCandidateAnalyzerDQM(const edm::ParameterSet& iConfig) {
0055   PFCandTag = iConfig.getParameter<edm::InputTag>("PFCandType");
0056   PFCandToken = consumes<edm::View<pat::PackedCandidate>>(PFCandTag);
0057   etabins = iConfig.getParameter<std::vector<double>>("etabins");
0058 
0059   //create map of pdgId
0060   std::vector<uint32_t> pdgKeys = iConfig.getParameter<std::vector<uint32_t>>("pdgKeys");
0061   std::vector<std::string> pdgStrs = iConfig.getParameter<std::vector<std::string>>("pdgStrs");
0062   for (int i = 0, n = pdgKeys.size(); i < n; i++)
0063     pdgMap[pdgKeys[i]] = pdgStrs[i];
0064 }
0065 
0066 void PFCandidateAnalyzerDQM::bookHistograms(DQMStore::IBooker& booker, edm::Run const&, edm::EventSetup const&) {
0067   // all candidate
0068   booker.setCurrentFolder("ParticleFlow/PackedCandidates/AllCandidate");
0069 
0070   // for eta binning
0071   int n = etabins.size() - 1;
0072   float etabinArray[etabins.size()];
0073   std::copy(etabins.begin(), etabins.end(), etabinArray);
0074 
0075   //eta has variable bin sizes, use 4th def of TH1F constructor
0076   TH1F* etaHist = new TH1F("AllCandidateEta", "AllCandidateEta", n, etabinArray);
0077   me["AllCandidateEta"] = booker.book1D("AllCandidateEta", etaHist);
0078 
0079   me["AllCandidateLog10Pt"] = booker.book1D("AllCandidateLog10Pt", "AllCandidateLog10Pt", 120, -2, 4);
0080 
0081   //for phi binnings
0082   double nPhiBins = 73;
0083   double phiBinWidth = M_PI / (nPhiBins - 1) * 2.;
0084   me["AllCandidatePhi"] = booker.book1D(
0085       "AllCandidatePhi", "AllCandidatePhi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);
0086 
0087   me["AllCandidateCharge"] = booker.book1D("AllCandidateCharge", "AllCandidateCharge", 3, -1.5, 1.5);
0088   me["AllCandidatePtLow"] = booker.book1D("AllCandidatePtLow", "AllCandidatePtLow", 100, 0., 5.);
0089   me["AllCandidatePtMid"] = booker.book1D("AllCandidatePtMid", "AllCandidatePtMid", 100, 0., 200.);
0090   me["AllCandidatePtHigh"] = booker.book1D("AllCandidatePtHigh", "AllCandidatePtHigh", 100, 0., 1000.);
0091 
0092   std::string etaHistName;
0093   for (auto& pair : pdgMap) {
0094     booker.setCurrentFolder("ParticleFlow/PackedCandidates/" + pair.second);
0095 
0096     //TH1F only takes char*, so have to do conversions for histogram name
0097     etaHistName = pair.second + "Eta";
0098     TH1F* etaHist = new TH1F(etaHistName.c_str(), etaHistName.c_str(), n, etabinArray);
0099     me[pair.second + "Eta"] = booker.book1D(pair.second + "Eta", etaHist);
0100 
0101     me[pair.second + "Log10Pt"] = booker.book1D(pair.second + "Log10Pt", pair.second + "Log10Pt", 120, -2, 4);
0102     me[pair.second + "Phi"] = booker.book1D(
0103         pair.second + "Phi", pair.second + "Phi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);
0104     me[pair.second + "Charge"] = booker.book1D(pair.second + "Charge", pair.second + "Charge", 3, -1.5, 1.5);
0105     me[pair.second + "PtLow"] = booker.book1D(pair.second + "PtLow", pair.second + "PtLow", 100, 0., 5.);
0106     me[pair.second + "PtMid"] = booker.book1D(pair.second + "PtMid", pair.second + "PtMid", 100, 0., 200.);
0107     me[pair.second + "PtHigh"] = booker.book1D(pair.second + "PtHigh", pair.second + "PtHigh", 100, 0., 1000.);
0108   }
0109 }
0110 
0111 void PFCandidateAnalyzerDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0112   //retrieve
0113   edm::Handle<edm::View<pat::PackedCandidate>> pfHandle;
0114   iEvent.getByToken(PFCandToken, pfHandle);
0115 
0116   if (!pfHandle.isValid()) {
0117     edm::LogInfo("OutputInfo") << " failed to retrieve data required by ParticleFlow Task";
0118     edm::LogInfo("OutputInfo") << " ParticleFlow Task cannot continue...!";
0119     return;
0120   } else {
0121     //Analyze
0122     // Loop Over Particle Flow Candidates
0123 
0124     for (unsigned int i = 0; i < pfHandle->size(); i++) {
0125       // Fill Histograms for Candidate Methods
0126       // all candidates
0127       me["AllCandidateLog10Pt"]->Fill(log10(pfHandle->at(i).pt()));
0128       me["AllCandidateEta"]->Fill(pfHandle->at(i).eta());
0129       me["AllCandidatePhi"]->Fill(pfHandle->at(i).phi());
0130       me["AllCandidateCharge"]->Fill(pfHandle->at(i).charge());
0131       me["AllCandidatePtLow"]->Fill(pfHandle->at(i).pt());
0132       me["AllCandidatePtMid"]->Fill(pfHandle->at(i).pt());
0133       me["AllCandidatePtHigh"]->Fill(pfHandle->at(i).pt());
0134 
0135       int pdgId = abs(pfHandle->at(i).pdgId());
0136       if (pdgMap.find(pdgId) != pdgMap.end()) {
0137         me[pdgMap[pdgId] + "Log10Pt"]->Fill(log10(pfHandle->at(i).pt()));
0138         me[pdgMap[pdgId] + "Eta"]->Fill(pfHandle->at(i).eta());
0139         me[pdgMap[pdgId] + "Phi"]->Fill(pfHandle->at(i).phi());
0140         me[pdgMap[pdgId] + "Charge"]->Fill(pfHandle->at(i).charge());
0141         me[pdgMap[pdgId] + "PtLow"]->Fill(pfHandle->at(i).pt());
0142         me[pdgMap[pdgId] + "PtMid"]->Fill(pfHandle->at(i).pt());
0143         me[pdgMap[pdgId] + "PtHigh"]->Fill(pfHandle->at(i).pt());
0144       }
0145     }
0146   }
0147 }
0148 
0149 #include "FWCore/Framework/interface/MakerMacros.h"
0150 DEFINE_FWK_MODULE(PFCandidateAnalyzerDQM);