CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFClient.cc
Go to the documentation of this file.
2 
7 
8 
12 
13 //
14 // -- Constructor
15 //
17 {
18  folderNames_ = parameterSet.getParameter< std::vector<std::string> >( "FolderNames" );
19  histogramNames_ = parameterSet.getParameter< std::vector<std::string> >( "HistogramNames" );
20  efficiencyFlag_ = parameterSet.getParameter< bool> ("CreateEfficiencyPlots" );
21  effHistogramNames_ = parameterSet.getParameter< std::vector<std::string> >( "HistogramNamesForEfficiencyPlots" );
22 }
23 //
24 // -- BeginJob
25 //
27 
29 }
30 //
31 // -- EndJobBegin Run
32 //
33 void PFClient::endRun(edm::Run const& run, edm::EventSetup const& eSetup) {
34  doSummaries();
36 }
37 //
38 // -- EndJob
39 //
41 
42 }
43 //
44 // -- Create Summaries
45 //
47 
48  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin();
49  ifolder != folderNames_.end(); ifolder++) {
50  std::string path = "ParticleFlow/"+(*ifolder);
51 
52  for (std::vector<std::string>::const_iterator ihist = histogramNames_.begin();
53  ihist != histogramNames_.end(); ihist++) {
54  std::string hname = (*ihist);
55  createResolutionPlots(path, hname);
56  }
57  }
58 }
59 //
60 // -- Create Summaries
61 //
63  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin();
64  ifolder != folderNames_.end(); ifolder++) {
65  std::string path = "ParticleFlow/"+(*ifolder);
66 
67  for (std::vector<std::string>::const_iterator ihist = effHistogramNames_.begin();
68  ihist != effHistogramNames_.end(); ihist++) {
69  std::string hname = (*ihist);
70  createEfficiencyPlots(path, hname);
71  }
72  }
73 }
74 //
75 // -- Create Resolution Plots
76 //
77 void PFClient::createResolutionPlots(std::string& folder, std::string& name) {
78  MonitorElement* me = dqmStore_->get(folder+"/"+name);
79  if (!me) return;
80  MonitorElement* me_average;
81  MonitorElement* me_rms;
82  MonitorElement* me_mean;
83  MonitorElement* me_sigma;
84  if ( (me->kind() == MonitorElement::DQM_KIND_TH2F) ||
86  (me->kind() == MonitorElement::DQM_KIND_TH2D) ) {
87  TH2* th = me->getTH2F();
88  size_t nbinx = me->getNbinsX();
89  size_t nbiny = me->getNbinsY();
90 
91  float ymin = th->GetYaxis()->GetXmin();
92  float ymax = th->GetYaxis()->GetXmax();
93  std::string xtit = th->GetXaxis()->GetTitle();
94  std::string ytit = th->GetYaxis()->GetTitle();
95  float* xbins = new float[nbinx+1];
96  for (size_t ix = 1; ix < nbinx+1; ++ix) {
97  xbins[ix-1] = th->GetBinLowEdge(ix);
98  if (ix == nbinx) xbins[ix] = th->GetXaxis()->GetBinUpEdge(ix);
99  }
100 
101  std::string tit_new;
102  dqmStore_->setCurrentFolder(folder);
103  MonitorElement* me_slice = dqmStore_->book1D("PFlowSlice","PFlowSlice",nbiny,ymin,ymax);
104 
105  tit_new = ";"+xtit+";Average_"+ytit;
106  me_average = dqmStore_->book1D("average_"+name,tit_new, nbinx, xbins);
107  tit_new = ";"+xtit+";RMS_"+ytit;
108  me_rms = dqmStore_->book1D("rms_"+name,tit_new, nbinx, xbins);
109  tit_new = ";"+xtit+";Mean_"+ytit;
110  me_mean = dqmStore_->book1D("mean_"+name,tit_new, nbinx, xbins);
111  tit_new = ";"+xtit+";Sigma_"+ytit;
112  me_sigma = dqmStore_->book1D("sigma_"+name,tit_new, nbinx, xbins);
113 
114  double average, rms, mean, sigma;
115  for (size_t ix = 1; ix < nbinx+1; ++ix) {
116  me_slice->Reset();
117  for (size_t iy = 1; iy < nbiny+1; ++iy) {
118  me_slice->setBinContent(iy,th->GetBinContent(ix,iy));
119  }
120  getHistogramParameters(me_slice, average, rms, mean, sigma);
121  me_average->setBinContent(ix,average);
122  me_rms->setBinContent(ix,rms);
123  me_mean->setBinContent(ix,mean);
124  me_sigma->setBinContent(ix,sigma);
125  }
126  if (me_slice) dqmStore_->removeElement(me_slice->getName());
127  delete [] xbins;
128  }
129 }
130 //
131 // -- Get Histogram Parameters
132 //
134  double& rms,double& mean, double& sigma) {
135  average = 0.0;
136  rms = 0.0;
137  mean = 0.0;
138  sigma = 0.0;
139 
140  if (!me_slice) return;
141  if (me_slice->kind() == MonitorElement::DQM_KIND_TH1F) {
142  average = me_slice->getMean();
143  rms = me_slice->getRMS();
144  TH1F* th_slice = me_slice->getTH1F();
145  if (th_slice && th_slice->GetEntries() > 0) {
146  th_slice->Fit( "gaus","Q0");
147  TF1* gaus = th_slice->GetFunction( "gaus" );
148  if (gaus) {
149  sigma = gaus->GetParameter(2);
150  mean = gaus->GetParameter(1);
151  }
152  }
153  }
154 }
155 //
156 // -- Create Resolution Plots
157 //
158 void PFClient::createEfficiencyPlots(std::string& folder, std::string& name) {
159  MonitorElement* me1 = dqmStore_->get(folder+"/"+name);
160  MonitorElement* me2 = dqmStore_->get(folder+"/"+name+"ref_");
161  if (!me1 || !me2) return;
162  MonitorElement* me_eff;
163  if ( (me1->kind() == MonitorElement::DQM_KIND_TH1F) &&
164  (me1->kind() == MonitorElement::DQM_KIND_TH1F) ) {
165  TH1* th1 = me1->getTH1F();
166  size_t nbinx = me1->getNbinsX();
167 
168  float xmin = th1->GetXaxis()->GetXmin();
169  float xmax = th1->GetXaxis()->GetXmax();
170  std::string xtit = me1->getAxisTitle(1);
171  std::string tit_new;
172  tit_new = ";"+xtit+";Efficiency";
173 
174  dqmStore_->setCurrentFolder(folder);
175  me_eff = dqmStore_->book1D("efficiency_"+name,tit_new, nbinx, xmin, xmax);
176 
177  double efficiency;
178  me_eff->Reset();
179  for (size_t ix = 1; ix < nbinx+1; ++ix) {
180  float val1 = me1->getBinContent(ix);
181  float val2 = me2->getBinContent(ix);
182  if (val2 > 0.0) efficiency = val1/val2;
183  else efficiency = 0;
184  me_eff->setBinContent(ix,efficiency);
185  }
186  }
187 }
T getParameter(std::string const &) const
const std::string & getName(void) const
get name of ME
void setBinContent(int binx, double content)
set content of bin (1-D)
void createEfficiencyPlots(std::string &folder, std::string &name)
Definition: PFClient.cc:158
const double xbins[]
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string getAxisTitle(int axis=1) const
get x-, y- or z-axis title (axis=1, 2, 3 respectively)
void doSummaries()
Definition: PFClient.cc:46
std::vector< std::string > histogramNames_
Definition: PFClient.h:31
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
void doEfficiency()
Definition: PFClient.cc:62
void beginJob()
Definition: PFClient.cc:26
void endRun(edm::Run const &run, edm::EventSetup const &eSetup)
Definition: PFClient.cc:33
int getNbinsY(void) const
get # of bins in Y-axis
list path
Definition: scaleCards.py:51
void createResolutionPlots(std::string &folder, std::string &name)
Definition: PFClient.cc:77
void removeElement(const std::string &name)
Definition: DQMStore.cc:2572
Kind kind(void) const
Get the type of the monitor element.
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1468
PFClient(const edm::ParameterSet &parameterSet)
Definition: PFClient.cc:16
std::vector< std::string > folderNames_
Definition: PFClient.h:30
TH1F * getTH1F(void) const
int average
Definition: PDRates.py:137
double getBinContent(int binx) const
get content of bin (1-D)
void endJob()
Definition: PFClient.cc:40
double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
std::vector< std::string > effHistogramNames_
Definition: PFClient.h:32
int getNbinsX(void) const
get # of bins in X-axis
TH2F * getTH2F(void) const
void Reset(void)
reset ME (ie. contents, errors, etc)
DQMStore * dqmStore_
Definition: PFClient.h:35
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
Definition: Run.h:33
bool efficiencyFlag_
Definition: PFClient.h:33
void getHistogramParameters(MonitorElement *me_slice, double &avarage, double &rms, double &mean, double &sigma)
Definition: PFClient.cc:133