CMS 3D CMS Logo

OffsetAnalyzerDQM.cc
Go to the documentation of this file.
8 
12 
13 #include <vector>
14 #include <map>
15 
17 public:
19  void analyze(const edm::Event&, const edm::EventSetup&) override;
20 
21 protected:
22  //Book histograms
23  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
24  void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override {}
25  void dqmEndRun(const edm::Run&, const edm::EventSetup&) override {}
26  int getEtaIndex(float eta);
27 
28 private:
29  struct Plot1D {
31  int nxbins;
32  double xlow, xhigh;
34 
35  Plot1D() : name(""), title(""), dir(""), nxbins(0), xlow(0), xhigh(0), plot(nullptr) {}
36  Plot1D(const std::string& n, const std::string& t, const std::string& d, int nx, double x0, double x1)
37  : name(n), title(t), dir(d), nxbins(nx), xlow(x0), xhigh(x1), plot(nullptr) {}
38 
39  virtual void book(DQMStore::IBooker& booker) {
40  booker.setCurrentFolder(dir);
41  plot = booker.book1D(name, title, nxbins, xlow, xhigh);
42  }
43 
44  virtual void fill(float value) {
45  assert(plot != nullptr);
46  plot->Fill(value);
47  }
48 
49  virtual ~Plot1D() {}
50  };
51 
52  struct PlotProfile : public Plot1D {
53  std::vector<double> xbins;
54  int nybins;
55  double ylow, yhigh;
56  PlotProfile() : Plot1D(), xbins(0), nybins(0), ylow(0), yhigh(0) {}
58  const std::string& t,
59  const std::string& d,
60  int nx,
61  double x0,
62  double x1,
63  const std::vector<double>& vx,
64  int ny,
65  double y0,
66  double y1)
67  : Plot1D(n, t, d, nx, x0, x1), xbins(vx), nybins(ny), ylow(y0), yhigh(y1) {}
68 
69  void book(DQMStore::IBooker& booker) override {
70  booker.setCurrentFolder(dir);
71  plot = booker.bookProfile(name, title, xbins.size() - 1, &xbins[0], nybins, ylow, yhigh, " ");
72  }
73  //make other booker methods for uniform binning
74 
75  void fill2D(double value1, double value2) {
76  assert(plot != nullptr);
77  plot->Fill(value1, value2);
78  }
79  };
80 
81  std::map<std::string, Plot1D> th1dPlots;
82  std::map<std::string, PlotProfile> offsetPlots;
83  std::map<int, std::string> pdgMap;
84 
86  std::vector<std::string> pftypes;
87  std::vector<double> etabins;
88 
89  int muHigh;
90  int npvHigh;
91 
95 };
96 
98  offsetPlotBaseName = iConfig.getParameter<std::string>("offsetPlotBaseName");
99 
100  pvToken = consumes<edm::View<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("pvTag"));
101  muToken = consumes<edm::View<PileupSummaryInfo>>(iConfig.getParameter<edm::InputTag>("muTag"));
102  pfToken = consumes<edm::View<pat::PackedCandidate>>(iConfig.getParameter<edm::InputTag>("pfTag"));
103 
104  etabins = iConfig.getParameter<std::vector<double>>("etabins");
105  pftypes = iConfig.getParameter<std::vector<std::string>>("pftypes");
106 
107  muHigh = iConfig.getUntrackedParameter<int>("muHigh");
108  npvHigh = iConfig.getUntrackedParameter<int>("npvHigh");
109 
110  //initialize offset plots
111  const auto& offset_psets = iConfig.getParameter<std::vector<edm::ParameterSet>>("offsetPlots");
112  for (auto& pset : offset_psets) {
113  std::string name = pset.getParameter<std::string>("name");
114  std::string title = pset.getParameter<std::string>("title");
115  std::string dir = pset.getParameter<std::string>("dir");
116  std::vector<double> vx = pset.getParameter<std::vector<double>>("vx");
117  int ny = pset.getParameter<uint32_t>("ny");
118  double y0 = pset.getParameter<double>("y0");
119  double y1 = pset.getParameter<double>("y1");
120 
121  offsetPlots[name] = PlotProfile(name, title, dir, 0, 0, 0, vx, ny, y0, y1);
122  }
123 
124  //initialize th1d
125  const auto& th1d_psets = iConfig.getParameter<std::vector<edm::ParameterSet>>("th1dPlots");
126  for (auto& pset : th1d_psets) {
127  std::string name = pset.getParameter<std::string>("name");
128  std::string title = pset.getParameter<std::string>("title");
129  std::string dir = pset.getParameter<std::string>("dir");
130  int nx = pset.getParameter<uint32_t>("nx");
131  double x0 = pset.getParameter<double>("x0");
132  double x1 = pset.getParameter<double>("x1");
133 
134  th1dPlots[name] = Plot1D(name, title, dir, nx, x0, x1);
135  }
136 
137  //create pdg map
138  std::vector<uint32_t> pdgKeys = iConfig.getParameter<std::vector<uint32_t>>("pdgKeys");
139  std::vector<std::string> pdgStrs = iConfig.getParameter<std::vector<std::string>>("pdgStrs");
140  for (int i = 0, n = pdgKeys.size(); i < n; i++)
141  pdgMap[pdgKeys[i]] = pdgStrs[i];
142 }
143 
145  //std::cout << "OffsetAnalyzerDQM booking offset histograms" << std::endl;
146  for (auto& pair : offsetPlots) {
147  pair.second.book(booker);
148  }
149  for (auto& pair : th1dPlots) {
150  pair.second.book(booker);
151  }
152 }
153 
155  //npv//
157  iEvent.getByToken(pvToken, vertexHandle);
158 
159  unsigned int nPVall = vertexHandle->size();
160  bool isGoodPV[nPVall];
161  for (size_t i = 0; i < nPVall; ++i)
162  isGoodPV[i] = false;
163 
164  int npv = 0;
165  for (unsigned int i = 0; i < nPVall; i++) {
166  const auto& pv = vertexHandle->at(i);
167 
168  if (!pv.isFake() && pv.ndof() >= 4 && fabs(pv.z()) <= 24.0 && fabs(pv.position().rho()) <= 2.0) {
169  npv++;
170  isGoodPV[i] = true;
171  }
172  }
173  th1dPlots["npv"].fill(npv);
174  int npv_in_range = npv;
175  if (npv_in_range < 0)
176  npv_in_range = 0;
177  else if (npv_in_range >= npvHigh)
178  npv_in_range = npvHigh - 1; // make sure int_mu won't lead to non-existing ME
179 
180  //mu//
181  int int_mu = -1;
183  if (iEvent.getByToken(muToken, muHandle)) {
184  const auto& summary = *muHandle;
185  auto it = std::find_if(summary.begin(), summary.end(), [](const auto& s) { return s.getBunchCrossing() == 0; });
186 
187  if (it->getBunchCrossing() != 0) {
188  edm::LogError("OffsetAnalyzerDQM") << "Cannot find the in-time pileup info " << it->getBunchCrossing();
189  } else {
190  float mu = it->getTrueNumInteractions();
191  th1dPlots["mu"].fill(mu);
192  int_mu = mu + 0.5;
193  }
194  }
195  if (int_mu >= muHigh)
196  int_mu = muHigh - 1; // make sure int_mu won't lead to non-existing ME
197 
198  //create map of pftypes vs total energy / eta
199  std::map<std::string, std::vector<double>> m_pftype_etaE;
200  int nEta = etabins.size() - 1;
201  for (const auto& pftype : pftypes)
202  m_pftype_etaE[pftype].assign(nEta, 0.0);
203 
204  //pf particles//
206  iEvent.getByToken(pfToken, pfHandle);
207 
208  for (unsigned int i = 0, n = pfHandle->size(); i < n; i++) {
209  const auto& cand = pfHandle->at(i);
210 
211  int etaIndex = getEtaIndex(cand.eta());
212  std::string pftype = pdgMap[abs(cand.pdgId())];
213  if (etaIndex == -1 || pftype.empty())
214  continue;
215 
216  if (pftype == "chm") { //check charged hadrons ONLY
217  bool attached = false;
218 
219  for (unsigned int ipv = 0; ipv < nPVall && !attached; ipv++) {
220  if (isGoodPV[ipv] && cand.fromPV(ipv) == 3)
221  attached = true; //pv used in fit
222  }
223  if (!attached)
224  pftype = "chu"; //unmatched charged hadron
225  }
227  /*
228  reco::TrackRef candTrkRef( cand.trackRef() );
229  if ( pftype == "chm" && !candTrkRef.isNull() ) { //check charged hadrons ONLY
230  bool attached = false;
231 
232  for (auto ipv=vertexHandle->begin(), endpv=vertexHandle->end(); ipv != endpv && !attached; ++ipv) {
233  if ( !ipv->isFake() && ipv->ndof() >= 4 && fabs(ipv->z()) < 24 ) { //must be attached to a good pv
234 
235  for(auto ivtrk=ipv->tracks_begin(), endvtrk=ipv->tracks_end(); ivtrk != endvtrk && !attached; ++ivtrk) {
236  reco::TrackRef pvTrkRef(ivtrk->castTo<reco::TrackRef>());
237  if (pvTrkRef == candTrkRef) attached = true;
238  }
239  }
240  }
241  if (!attached) pftype = "chu"; //unmatched charged hadron
242  }
243 */
245  m_pftype_etaE[pftype][etaIndex] += cand.et();
246  }
247 
248  for (const auto& pair : m_pftype_etaE) {
249  std::string pftype = pair.first;
250  std::vector<double> etaE = pair.second;
251 
252  std::string offset_name_npv = offsetPlotBaseName + "_npv" + std::to_string(npv_in_range) + "_" + pftype;
253  if (offsetPlots.find(offset_name_npv) == offsetPlots.end())
254  return; //npv is out of range ()
255 
256  for (int i = 0; i < nEta; i++) {
257  double eta = 0.5 * (etabins[i] + etabins[i + 1]);
258  offsetPlots[offset_name_npv].fill2D(eta, etaE[i]);
259  }
260 
261  if (int_mu != -1) {
262  std::string offset_name_mu = offsetPlotBaseName + "_mu" + std::to_string(int_mu) + "_" + pftype;
263  if (offsetPlots.find(offset_name_mu) == offsetPlots.end())
264  return; //mu is out of range
265 
266  for (int i = 0; i < nEta; i++) {
267  double eta = 0.5 * (etabins[i] + etabins[i + 1]);
268  offsetPlots[offset_name_mu].fill2D(eta, etaE[i]);
269  }
270  }
271  }
272 }
273 
275  int nEta = etabins.size() - 1;
276 
277  for (int i = 0; i < nEta; i++) {
278  if (etabins[i] <= eta && eta < etabins[i + 1])
279  return i;
280  }
281  if (eta == etabins[nEta])
282  return nEta - 1;
283  else
284  return -1;
285 }
286 
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void book(DQMStore::IBooker &booker) override
virtual void fill(float value)
edm::EDGetTokenT< edm::View< PileupSummaryInfo > > muToken
std::map< int, std::string > pdgMap
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
#define nullptr
std::map< std::string, Plot1D > th1dPlots
OffsetAnalyzerDQM(const edm::ParameterSet &)
std::vector< double > etabins
void Fill(long long x)
edm::EDGetTokenT< edm::View< pat::PackedCandidate > > pfToken
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int getEtaIndex(float eta)
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, char const *option="s")
Definition: DQMStore.cc:333
void analyze(const edm::Event &, const edm::EventSetup &) override
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: value.py:1
std::string offsetPlotBaseName
d
Definition: ztail.py:151
std::map< std::string, PlotProfile > offsetPlots
Plot1D(const std::string &n, const std::string &t, const std::string &d, int nx, double x0, double x1)
PlotProfile(const std::string &n, const std::string &t, const std::string &d, int nx, double x0, double x1, const std::vector< double > &vx, int ny, double y0, double y1)
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
void fill2D(double value1, double value2)
edm::EDGetTokenT< edm::View< reco::Vertex > > pvToken
void dqmEndRun(const edm::Run &, const edm::EventSetup &) override
virtual void book(DQMStore::IBooker &booker)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: Run.h:45
std::vector< std::string > pftypes