CMS 3D CMS Logo

PFCandidateAnalyzerDQM.cc
Go to the documentation of this file.
1 // Producer for particle flow candidates. Plots Eta, Phi, Charge, Pt (log freq, bin)
2 // for different types of particles described in python/defaults_cfi.py
3 // It actually uses packedCandidates so that we need only MINIAOD contents to run this DQMAnalyzer.
4 // note: for pt, log freq is done in this producer, but log freq is done by running
5 // compare.py
6 // author: Chosila Sutantawibul, April 23, 2020
7 
11 
15 
22 
23 #include "TH1F.h"
24 
25 #include <algorithm>
26 #include <cmath>
27 #include <fstream>
28 #include <iostream>
29 #include <memory>
30 #include <ostream>
31 #include <map>
32 #include <string>
33 #include <cstring>
34 
36 public:
38  void analyze(const edm::Event&, const edm::EventSetup&) override;
39 
40 protected:
41  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
42 
43 private:
44  //from config file
47  std::vector<double> etabins;
48  std::map<std::string, MonitorElement*> me;
49 
50  std::map<uint32_t, std::string> pdgMap;
51 };
52 
53 // constructor
55  PFCandTag = iConfig.getParameter<edm::InputTag>("PFCandType");
56  PFCandToken = consumes<edm::View<pat::PackedCandidate>>(PFCandTag);
57  etabins = iConfig.getParameter<std::vector<double>>("etabins");
58 
59  //create map of pdgId
60  std::vector<uint32_t> pdgKeys = iConfig.getParameter<std::vector<uint32_t>>("pdgKeys");
61  std::vector<std::string> pdgStrs = iConfig.getParameter<std::vector<std::string>>("pdgStrs");
62  for (int i = 0, n = pdgKeys.size(); i < n; i++)
63  pdgMap[pdgKeys[i]] = pdgStrs[i];
64 }
65 
67  // all candidate
68  booker.setCurrentFolder("ParticleFlow/PackedCandidates/AllCandidate");
69 
70  // for eta binning
71  int n = etabins.size() - 1;
72  float etabinArray[etabins.size()];
73  std::copy(etabins.begin(), etabins.end(), etabinArray);
74 
75  //eta has variable bin sizes, use 4th def of TH1F constructor
76  TH1F* etaHist = new TH1F("AllCandidateEta", "AllCandidateEta", n, etabinArray);
77  me["AllCandidateEta"] = booker.book1D("AllCandidateEta", etaHist);
78 
79  me["AllCandidateLog10Pt"] = booker.book1D("AllCandidateLog10Pt", "AllCandidateLog10Pt", 120, -2, 4);
80 
81  //for phi binnings
82  double nPhiBins = 73;
83  double phiBinWidth = M_PI / (nPhiBins - 1) * 2.;
84  me["AllCandidatePhi"] = booker.book1D(
85  "AllCandidatePhi", "AllCandidatePhi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);
86 
87  me["AllCandidateCharge"] = booker.book1D("AllCandidateCharge", "AllCandidateCharge", 3, -1.5, 1.5);
88  me["AllCandidatePtLow"] = booker.book1D("AllCandidatePtLow", "AllCandidatePtLow", 100, 0., 5.);
89  me["AllCandidatePtMid"] = booker.book1D("AllCandidatePtMid", "AllCandidatePtMid", 100, 0., 200.);
90  me["AllCandidatePtHigh"] = booker.book1D("AllCandidatePtHigh", "AllCandidatePtHigh", 100, 0., 1000.);
91 
92  std::string etaHistName;
93  for (auto& pair : pdgMap) {
94  booker.setCurrentFolder("ParticleFlow/PackedCandidates/" + pair.second);
95 
96  //TH1F only takes char*, so have to do conversions for histogram name
97  etaHistName = pair.second + "Eta";
98  TH1F* etaHist = new TH1F(etaHistName.c_str(), etaHistName.c_str(), n, etabinArray);
99  me[pair.second + "Eta"] = booker.book1D(pair.second + "Eta", etaHist);
100 
101  me[pair.second + "Log10Pt"] = booker.book1D(pair.second + "Log10Pt", pair.second + "Log10Pt", 120, -2, 4);
102  me[pair.second + "Phi"] = booker.book1D(
103  pair.second + "Phi", pair.second + "Phi", nPhiBins, -M_PI - 0.25 * phiBinWidth, +M_PI + 0.75 * phiBinWidth);
104  me[pair.second + "Charge"] = booker.book1D(pair.second + "Charge", pair.second + "Charge", 3, -1.5, 1.5);
105  me[pair.second + "PtLow"] = booker.book1D(pair.second + "PtLow", pair.second + "PtLow", 100, 0., 5.);
106  me[pair.second + "PtMid"] = booker.book1D(pair.second + "PtMid", pair.second + "PtMid", 100, 0., 200.);
107  me[pair.second + "PtHigh"] = booker.book1D(pair.second + "PtHigh", pair.second + "PtHigh", 100, 0., 1000.);
108  }
109 }
110 
112  //retrieve
114  iEvent.getByToken(PFCandToken, pfHandle);
115 
116  if (!pfHandle.isValid()) {
117  edm::LogInfo("OutputInfo") << " failed to retrieve data required by ParticleFlow Task";
118  edm::LogInfo("OutputInfo") << " ParticleFlow Task cannot continue...!";
119  return;
120  } else {
121  //Analyze
122  // Loop Over Particle Flow Candidates
123 
124  for (unsigned int i = 0; i < pfHandle->size(); i++) {
125  // Fill Histograms for Candidate Methods
126  // all candidates
127  me["AllCandidateLog10Pt"]->Fill(log10(pfHandle->at(i).pt()));
128  me["AllCandidateEta"]->Fill(pfHandle->at(i).eta());
129  me["AllCandidatePhi"]->Fill(pfHandle->at(i).phi());
130  me["AllCandidateCharge"]->Fill(pfHandle->at(i).charge());
131  me["AllCandidatePtLow"]->Fill(pfHandle->at(i).pt());
132  me["AllCandidatePtMid"]->Fill(pfHandle->at(i).pt());
133  me["AllCandidatePtHigh"]->Fill(pfHandle->at(i).pt());
134 
135  int pdgId = abs(pfHandle->at(i).pdgId());
136  if (pdgMap.find(pdgId) != pdgMap.end()) {
137  me[pdgMap[pdgId] + "Log10Pt"]->Fill(log10(pfHandle->at(i).pt()));
138  me[pdgMap[pdgId] + "Eta"]->Fill(pfHandle->at(i).eta());
139  me[pdgMap[pdgId] + "Phi"]->Fill(pfHandle->at(i).phi());
140  me[pdgMap[pdgId] + "Charge"]->Fill(pfHandle->at(i).charge());
141  me[pdgMap[pdgId] + "PtLow"]->Fill(pfHandle->at(i).pt());
142  me[pdgMap[pdgId] + "PtMid"]->Fill(pfHandle->at(i).pt());
143  me[pdgMap[pdgId] + "PtHigh"]->Fill(pfHandle->at(i).pt());
144  }
145  }
146  }
147 }
148 
Handle.h
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
filterCSVwithJSON.copy
copy
Definition: filterCSVwithJSON.py:36
ESHandle.h
edm::Run
Definition: Run.h:45
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm::LogInfo
Definition: MessageLogger.h:254
offsetAnalyzerDQM_cff.pdgKeys
pdgKeys
Definition: offsetAnalyzerDQM_cff.py:70
dqm::implementation::NavigatorBase::setCurrentFolder
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
PFCandidateAnalyzerDQM::pdgMap
std::map< uint32_t, std::string > pdgMap
Definition: PFCandidateAnalyzerDQM.cc:50
DQMStore.h
PFCandidateAnalyzerDQM::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: PFCandidateAnalyzerDQM.cc:111
edm::Handle
Definition: AssociativeIterator.h:50
MakerMacros.h
PFCandidateAnalyzerDQM::me
std::map< std::string, MonitorElement * > me
Definition: PFCandidateAnalyzerDQM.cc:48
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PFCandidateAnalyzerDQM::PFCandidateAnalyzerDQM
PFCandidateAnalyzerDQM(const edm::ParameterSet &)
Definition: PFCandidateAnalyzerDQM.cc:54
PFCandidateAnalyzerDQM::PFCandToken
edm::EDGetTokenT< edm::View< pat::PackedCandidate > > PFCandToken
Definition: PFCandidateAnalyzerDQM.cc:46
Service.h
DQMEDAnalyzer.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
DQMEDAnalyzer
Definition: DQMEDAnalyzer.py:1
edm::ParameterSet
Definition: ParameterSet.h:36
Event.h
PackedCandidate.h
PFCandidateAnalyzerDQM
Definition: PFCandidateAnalyzerDQM.cc:35
iEvent
int iEvent
Definition: GenABIO.cc:224
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
edm::EventSetup
Definition: EventSetup.h:57
PFCandidateAnalyzerDQM::bookHistograms
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: PFCandidateAnalyzerDQM.cc:66
InputTag.h
offsetAnalyzerDQM_cff.pdgStrs
pdgStrs
Definition: offsetAnalyzerDQM_cff.py:71
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
PFCandidateAnalyzerDQM::PFCandTag
edm::InputTag PFCandTag
Definition: PFCandidateAnalyzerDQM.cc:45
Frameworkfwd.h
dqm::implementation::IBooker
Definition: DQMStore.h:43
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
ecaldqm::binning::nPhiBins
Definition: MESetBinningUtils.h:69
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
edm::Event
Definition: Event.h:73
PFCandidateAnalyzerDQM::etabins
std::vector< double > etabins
Definition: PFCandidateAnalyzerDQM.cc:47
edm::InputTag
Definition: InputTag.h:15
dqm::implementation::IBooker::book1D
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98