CMS 3D CMS Logo

OffsetDQMPostProcessor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/RecoParticleFlow
4 // Class: OffsetDQMPostProcessor.cc
5 //
6 // Original Author: "Kenichi Hatakeyama"
7 //
8 
15 
18 
19 //
20 // class decleration
21 //
22 
24 public:
26  ~OffsetDQMPostProcessor() override;
27 
28 private:
30 
33  float offsetR;
34  std::vector<std::string> pftypes;
35  std::vector<std::string> offsetVariableTypes;
36  int muHigh;
37  int npvHigh;
38 
39  bool debug = false;
40 };
41 
42 //
43 // constructors and destructor
44 //
46  offsetPlotBaseName = iConfig.getParameter<std::string>("offsetPlotBaseName");
47  offsetDir = iConfig.getParameter<std::string>("offsetDir");
48  offsetVariableTypes = iConfig.getParameter<std::vector<std::string> >("offsetVariableTypes");
49  offsetR = iConfig.getUntrackedParameter<double>("offsetR");
50  pftypes = iConfig.getParameter<std::vector<std::string> >("pftypes");
51  muHigh = iConfig.getUntrackedParameter<int>("muHigh");
52  npvHigh = iConfig.getUntrackedParameter<int>("npvHigh");
53 };
54 
56 
57 // ------------ method called right after a run ends ------------
60 
61  std::string stitle;
62  std::vector<MonitorElement*> vME;
63  std::vector<std::string> MEStrings = iget_.getMEs();
64  std::for_each(MEStrings.begin(), MEStrings.end(), [&](auto& s) { s.insert(0, offsetDir); });
65 
66  // temporary ME and root objects
67  MonitorElement* mtmp;
68  TProfile* hproftmp;
69  TH1F* htmp;
70  TH1F* hscaled;
71 
72  //
73  // Offset plots vs eta
74  //
75  for (std::vector<std::string>::const_iterator i = offsetVariableTypes.begin(); i != offsetVariableTypes.end(); ++i) {
76  //
77  // getting the average value for Npv and mu
78  //
79  stitle = offsetDir + (*i);
80  std::vector<std::string>::const_iterator it = std::find(MEStrings.begin(), MEStrings.end(), stitle);
81  if (it == MEStrings.end())
82  continue;
83  mtmp = iget_.get(stitle);
84  float avg = mtmp->getMean();
85  int iavg = int(avg + 0.5); // integer version for identifying correcping ME, in order to get the rounding correctly
86 
87  if (avg < 1.)
88  avg = 1.; // protection against this value going too low
89 
90  if (iavg < 0)
91  iavg = 0; // checking lower bound (avoid division by zero)
92  if (*i == "npv" && iavg >= npvHigh)
93  iavg = npvHigh - 1; // checking upper bound
94  else if (*i == "mu" && iavg >= muHigh)
95  iavg = muHigh - 1; //
96 
97  //
98  // storing the value
99  //
100  stitle = (*i) + "_mean";
101  MonitorElement* MEmean = ibook_.bookFloat(stitle);
102  MEmean->Fill(avg);
103  vME.push_back(MEmean);
104 
105  //
106  // for each pf types
107  //
108  for (std::vector<std::string>::const_iterator j = pftypes.begin(); j != pftypes.end(); ++j) {
109  // accessing profiles
110  std::string str_base = *i + std::to_string(iavg);
111  if ((*i) == "npv")
112  stitle = offsetDir + "npvPlots/" + str_base + "/" + offsetPlotBaseName + "_" + str_base + "_" + (*j);
113  else if ((*i) == "mu")
114  stitle = offsetDir + "muPlots/" + str_base + "/" + offsetPlotBaseName + "_" + str_base + "_" + (*j);
115  else
116  return;
117 
118  // making scaled plot and ME
119  mtmp = iget_.get(stitle);
120  hproftmp = (TProfile*)mtmp->getTProfile();
121  htmp = (TH1F*)hproftmp->ProjectionX();
122  TAxis* xaxis = (TAxis*)htmp->GetXaxis();
123  stitle = offsetPlotBaseName + "_" + str_base + "_" + *j;
124  hscaled = new TH1F(stitle.c_str(), stitle.c_str(), xaxis->GetNbins(), xaxis->GetXbins()->GetArray());
125 
126  htmp->Scale(pow(offsetR, 2) / 2. / float(avg)); // pi*R^2 / (deltaEta*2pi) / <mu or NPV>
127  for (int ibin = 1; ibin <= hscaled->GetNbinsX(); ibin++) { // 1/deltaEta part
128  hscaled->SetBinContent(ibin, htmp->GetBinContent(ibin) / htmp->GetBinWidth(ibin));
129  hscaled->SetBinError(ibin, htmp->GetBinError(ibin) / htmp->GetBinWidth(ibin));
130  }
131 
132  // storing new ME
133  stitle = offsetPlotBaseName + "_" + *i + "_" + *j;
134  mtmp = ibook_.book1D(stitle.c_str(), hscaled);
135  vME.push_back(mtmp);
136  }
137  }
138 
139  //
140  // Checks
141  //
142  if (debug) {
143  for (std::vector<MonitorElement*>::const_iterator i = vME.begin(); i != vME.end(); ++i)
144  (*i)->getTH1F()->Print();
145  }
146 }
147 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * bookFloat(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:80
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
virtual TProfile * getTProfile() const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
virtual std::vector< std::string > getMEs() const
Definition: DQMStore.cc:737
std::string to_string(const V &value)
Definition: OMSAccess.h:71
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
T getUntrackedParameter(std::string const &, T const &) const
void Fill(long long x)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< std::string > offsetVariableTypes
std::vector< std::string > pftypes
virtual TH1F * getTH1F() const
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:690
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
OffsetDQMPostProcessor(const edm::ParameterSet &)
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29