CMS 3D CMS Logo

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