CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
PFClient_JetRes Class Reference

#include <PFClient_JetRes.h>

Inheritance diagram for PFClient_JetRes:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 PFClient_JetRes (const edm::ParameterSet &parameterSet)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void analyze (edm::Event const &, edm::EventSetup const &)
 
void beginJob ()
 
void createEfficiencyPlots (std::string &folder, std::string &name)
 
void createResolutionPlots (std::string &folder, std::string &name)
 
void doEfficiency ()
 
void doSummaries ()
 
void endJob ()
 
void endRun (edm::Run const &run, edm::EventSetup const &eSetup)
 
void getHistogramParameters (MonitorElement *me_slice, double &avarage, double &rms, double &mean, double &sigma)
 

Private Attributes

DQMStoredqmStore_
 
std::vector< std::string > effHistogramNames_
 
bool efficiencyFlag_
 
std::vector< std::string > folderNames_
 
std::vector< std::string > histogramNames_
 
std::vector< int > PtBins_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 12 of file PFClient_JetRes.h.

Constructor & Destructor Documentation

PFClient_JetRes::PFClient_JetRes ( const edm::ParameterSet parameterSet)

Definition at line 18 of file PFClient_JetRes.cc.

References effHistogramNames_, efficiencyFlag_, folderNames_, edm::ParameterSet::getParameter(), histogramNames_, and PtBins_.

19 {
20  folderNames_ = parameterSet.getParameter< std::vector<std::string> >( "FolderNames" );
21  histogramNames_ = parameterSet.getParameter< std::vector<std::string> >( "HistogramNames" );
22  efficiencyFlag_ = parameterSet.getParameter< bool> ("CreateEfficiencyPlots" );
23  effHistogramNames_ = parameterSet.getParameter< std::vector<std::string> >( "HistogramNamesForEfficiencyPlots" );
24  PtBins_ = parameterSet.getParameter< std::vector<int> >( "VariablePtBins" ) ;
25 
26 }
T getParameter(std::string const &) const
std::vector< std::string > effHistogramNames_
std::vector< int > PtBins_
std::vector< std::string > folderNames_
std::vector< std::string > histogramNames_

Member Function Documentation

void PFClient_JetRes::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
inlineprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 19 of file PFClient_JetRes.h.

19 {;}
void PFClient_JetRes::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 30 of file PFClient_JetRes.cc.

References dqmStore_, and cppFunctionSkipper::operator.

void PFClient_JetRes::createEfficiencyPlots ( std::string &  folder,
std::string &  name 
)
private

Definition at line 199 of file PFClient_JetRes.cc.

References DQMStore::book1D(), MonitorElement::DQM_KIND_TH1F, dqmStore_, postValidation_cfi::efficiency, DQMStore::get(), MonitorElement::getAxisTitle(), MonitorElement::getBinContent(), MonitorElement::getNbinsX(), MonitorElement::getTH1F(), MonitorElement::kind(), MonitorElement::Reset(), MonitorElement::setBinContent(), DQMStore::setCurrentFolder(), MonitorElement::setEfficiencyFlag(), AlCaHLTBitMon_QueryRunRegistry::string, SiStripMonitorClusterAlca_cfi::xmax, and SiStripMonitorClusterAlca_cfi::xmin.

Referenced by doEfficiency().

199  {
200  MonitorElement* me1 = dqmStore_->get(folder+"/"+name);
201  MonitorElement* me2 = dqmStore_->get(folder+"/"+name+"ref_");
202  if (!me1 || !me2) return;
203  MonitorElement* me_eff;
204  if ( (me1->kind() == MonitorElement::DQM_KIND_TH1F) &&
205  (me1->kind() == MonitorElement::DQM_KIND_TH1F) ) {
206  TH1* th1 = me1->getTH1F();
207  size_t nbinx = me1->getNbinsX();
208 
209  float xmin = th1->GetXaxis()->GetXmin();
210  float xmax = th1->GetXaxis()->GetXmax();
211  std::string xtit = me1->getAxisTitle(1);
212  std::string tit_new;
213  tit_new = ";"+xtit+";Efficiency";
214 
215  dqmStore_->setCurrentFolder(folder);
216  me_eff = dqmStore_->book1D("efficiency_"+name,tit_new, nbinx, xmin, xmax);
217 
218  double efficiency;
219  me_eff->Reset(); me_eff->setEfficiencyFlag();
220  for (size_t ix = 1; ix < nbinx+1; ++ix) {
221  float val1 = me1->getBinContent(ix);
222  float val2 = me2->getBinContent(ix);
223  if (val2 > 0.0) efficiency = val1/val2;
224  else efficiency = 0;
225  me_eff->setBinContent(ix,efficiency);
226  }
227  }
228 }
void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:964
std::string getAxisTitle(int axis=1) const
get x-, y- or z-axis title (axis=1, 2, 3 respectively)
DQMStore * dqmStore_
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:1718
TH1F * getTH1F(void) const
double getBinContent(int binx) const
get content of bin (1-D)
int getNbinsX(void) const
get # of bins in X-axis
void setEfficiencyFlag(void)
void Reset(void)
reset ME (ie. contents, errors, etc)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:677
void PFClient_JetRes::createResolutionPlots ( std::string &  folder,
std::string &  name 
)
private

Definition at line 81 of file PFClient_JetRes.cc.

References PDRates::average, DQMStore::book1D(), MonitorElement::DQM_KIND_TH2D, MonitorElement::DQM_KIND_TH2F, MonitorElement::DQM_KIND_TH2S, dqmStore_, DQMStore::get(), getHistogramParameters(), MonitorElement::getNbinsY(), MonitorElement::getTH2F(), MonitorElement::kind(), timingPdfMaker::mean, PtBins_, plotscripts::rms(), MonitorElement::setBinContent(), MonitorElement::setBinError(), DQMStore::setCurrentFolder(), MonitorElement::setEfficiencyFlag(), MonitorElement::setEntries(), AlCaHLTBitMon_QueryRunRegistry::string, fw3dlego::xbins, SiStripMonitorClusterAlca_cfi::ymax, and SiStripMonitorClusterAlca_cfi::ymin.

Referenced by doSummaries().

81  {
82  MonitorElement* me = dqmStore_->get(folder+"/"+name);
83  if (!me) return;
84 
85  MonitorElement* pT[PtBins_.size() -1];
86  std::vector<double> pTEntries(PtBins_.size()-1, 0) ;
87 
88  //std::vector<std::string> pTRange (PtBins_.size() -1) ;
89  //char* pTRange[PtBins_.size() -1] ;
90  TString pTRange[PtBins_.size() -1] ;
91  //float pTCenter[PtBins_.size() -1] ;
92 
93  MonitorElement* me_average;
94  MonitorElement* me_rms;
95  MonitorElement* me_mean;
96  MonitorElement* me_sigma;
97 
98  if ( (me->kind() == MonitorElement::DQM_KIND_TH2F) ||
100  (me->kind() == MonitorElement::DQM_KIND_TH2D) ) {
101  TH2* th = me->getTH2F();
102  //size_t nbinx = me->getNbinsX();
103  size_t nbinx = PtBins_.size() -1;
104  size_t nbiny = me->getNbinsY();
105 
106  float ymin = th->GetYaxis()->GetXmin();
107  float ymax = th->GetYaxis()->GetXmax();
108  std::string xtit = th->GetXaxis()->GetTitle();
109  //std::string ytit = th->GetYaxis()->GetTitle();
110  std::string ytit = "#Deltap_{T}/p_{T}";
111 
112  float* xbins = new float[nbinx+1];
113  for (size_t ix = 1; ix < nbinx+1; ++ix) {
114  //xbins[ix-1] = th->GetBinLowEdge(ix);
115  xbins[ix-1] = PtBins_[ix-1];
116  //snprintf(pTRange[ix-1].data(), 15, "Pt%d_%d", PtBins_[ix-1], PtBins_[ix]);
117  pTRange[ix-1] = TString::Format("Pt%d_%d", PtBins_[ix-1], PtBins_[ix]) ;
118  if (name == "BRdelta_et_Over_et_VS_et_") pTRange[ix-1] = TString::Format("BRPt%d_%d", PtBins_[ix-1], PtBins_[ix]) ;
119  else if (name == "ERdelta_et_Over_et_VS_et_") pTRange[ix-1] = TString::Format("ERPt%d_%d", PtBins_[ix-1], PtBins_[ix]) ;
120 
121  //pTCenter[ix-1] = (PtBins_[ix] - PtBins_[ix-1]) / 2. ;
122  if (ix == nbinx) {
123  //xbins[ix] = th->GetXaxis()->GetBinUpEdge(ix);
124  xbins[ix] = PtBins_[ix];
125  }
126  }
127 
128  std::string tit_new;
129  dqmStore_->setCurrentFolder(folder);
130  //MonitorElement* me_slice = dqmStore_->book1D("PFlowSlice", "PFlowSlice", nbiny, ymin, ymax);
131 
132  tit_new = "Average "+ytit+";"+xtit+";Average_"+ytit ;
133  me_average = dqmStore_->book1D("average_"+name,tit_new, nbinx, xbins);
134  me_average->setEfficiencyFlag();
135  tit_new = "RMS "+ytit+";"+xtit+";RMS_"+ytit ;
136  me_rms = dqmStore_->book1D("rms_"+name,tit_new, nbinx, xbins);
137  me_rms->setEfficiencyFlag();
138  tit_new = ";"+xtit+";Mean_"+ytit;
139  me_mean = dqmStore_->book1D("mean_"+name,tit_new, nbinx, xbins);
140  me_mean->setEfficiencyFlag();
141  tit_new = ";"+xtit+";Sigma_"+ytit;
142  me_sigma = dqmStore_->book1D("sigma_"+name,tit_new, nbinx, xbins);
143  me_sigma->setEfficiencyFlag();
144 
145  double average, rms, mean, sigma;
146  for (size_t ix = 1; ix < nbinx+1; ++ix) {
147  //me_slice->Reset();
148  if (name == "delta_et_Over_et_VS_et_") pT[ix-1] = dqmStore_->book1D(pTRange[ix-1], TString::Format("Total %s;%s;Events", ytit.data(), ytit.data() ), nbiny, ymin, ymax) ;
149  if (name == "BRdelta_et_Over_et_VS_et_") pT[ix-1] = dqmStore_->book1D(pTRange[ix-1], TString::Format("Barrel %s;%s;Events", ytit.data(), ytit.data()), nbiny, ymin, ymax) ;
150  else if (name == "ERdelta_et_Over_et_VS_et_") pT[ix-1] = dqmStore_->book1D(pTRange[ix-1], TString::Format("Endcap %s;%s;Events", ytit.data(), ytit.data() ), nbiny, ymin, ymax) ;
151 
152  for (size_t iy = 0; iy <= nbiny+1; ++iy) // add under and overflow
153  if (th->GetBinContent(ix,iy)) {
154  //me_slice->setBinContent(iy,th->GetBinContent(ix,iy));
155  pT[ix-1]->setBinContent(iy, th->GetBinContent(ix,iy));
156  pT[ix-1]->setBinError(iy, th->GetBinError(ix,iy));
157  pTEntries[ix-1] += th->GetBinContent(ix,iy); }
158 
159  pT[ix-1]->setEntries( pTEntries[ix-1] ) ;
160 
161  //getHistogramParameters(me_slice, average, rms, mean, sigma);
162  getHistogramParameters(pT[ix-1], average, rms, mean, sigma);
163  me_average->setBinContent(ix,average);
164  me_rms->setBinContent(ix,rms);
165  me_mean->setBinContent(ix,mean);
166  me_sigma->setBinContent(ix,sigma);
167  }
168  //if (me_slice) dqmStore_->removeElement(me_slice->getName());
169  delete [] xbins;
170  }
171 }
void setBinContent(int binx, double content)
set content of bin (1-D)
const double xbins[]
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:964
int getNbinsY(void) const
get # of bins in Y-axis
DQMStore * dqmStore_
void getHistogramParameters(MonitorElement *me_slice, double &avarage, double &rms, double &mean, double &sigma)
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
Kind kind(void) const
Get the type of the monitor element.
void setEntries(double nentries)
set # of entries
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1718
std::vector< int > PtBins_
int average
Definition: PDRates.py:137
void setEfficiencyFlag(void)
TH2F * getTH2F(void) const
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:677
void PFClient_JetRes::doEfficiency ( )
private

Definition at line 66 of file PFClient_JetRes.cc.

References createEfficiencyPlots(), effHistogramNames_, folderNames_, cmsHarvester::path, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by endRun().

66  {
67  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin();
68  ifolder != folderNames_.end(); ifolder++) {
69  std::string path = "ParticleFlow/"+(*ifolder);
70 
71  for (std::vector<std::string>::const_iterator ihist = effHistogramNames_.begin();
72  ihist != effHistogramNames_.end(); ihist++) {
73  std::string hname = (*ihist);
74  createEfficiencyPlots(path, hname);
75  }
76  }
77 }
tuple path
else: Piece not in the list, fine.
void createEfficiencyPlots(std::string &folder, std::string &name)
std::vector< std::string > effHistogramNames_
std::vector< std::string > folderNames_
void PFClient_JetRes::doSummaries ( )
private

Definition at line 50 of file PFClient_JetRes.cc.

References createResolutionPlots(), folderNames_, histogramNames_, cmsHarvester::path, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by endRun().

50  {
51 
52  for (std::vector<std::string>::const_iterator ifolder = folderNames_.begin();
53  ifolder != folderNames_.end(); ifolder++) {
54  std::string path = "ParticleFlow/"+(*ifolder);
55 
56  for (std::vector<std::string>::const_iterator ihist = histogramNames_.begin();
57  ihist != histogramNames_.end(); ihist++) {
58  std::string hname = (*ihist);
59  createResolutionPlots(path, hname);
60  }
61  }
62 }
tuple path
else: Piece not in the list, fine.
void createResolutionPlots(std::string &folder, std::string &name)
std::vector< std::string > folderNames_
std::vector< std::string > histogramNames_
void PFClient_JetRes::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 44 of file PFClient_JetRes.cc.

44  {
45 
46 }
void PFClient_JetRes::endRun ( edm::Run const &  run,
edm::EventSetup const &  eSetup 
)
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 37 of file PFClient_JetRes.cc.

References doEfficiency(), doSummaries(), and efficiencyFlag_.

37  {
38  doSummaries();
40 }
void PFClient_JetRes::getHistogramParameters ( MonitorElement me_slice,
double &  avarage,
double &  rms,
double &  mean,
double &  sigma 
)
private

Definition at line 175 of file PFClient_JetRes.cc.

References MonitorElement::DQM_KIND_TH1F, MonitorElement::getMean(), MonitorElement::getRMS(), MonitorElement::getTH1F(), and MonitorElement::kind().

Referenced by createResolutionPlots().

176  {
177  average = 0.0;
178  rms = 0.0;
179  mean = 0.0;
180  sigma = 0.0;
181 
182  if (!me_slice) return;
183  if (me_slice->kind() == MonitorElement::DQM_KIND_TH1F) {
184  average = me_slice->getMean();
185  rms = me_slice->getRMS();
186  TH1F* th_slice = me_slice->getTH1F();
187  if (th_slice && th_slice->GetEntries() > 0) {
188  //need our own copy for thread safety
189  TF1 gaus("mygaus","gaus");
190  th_slice->Fit( &gaus,"Q0");
191  sigma = gaus.GetParameter(2);
192  mean = gaus.GetParameter(1);
193  }
194  }
195 }
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
Kind kind(void) const
Get the type of the monitor element.
TH1F * getTH1F(void) const
int average
Definition: PDRates.py:137
double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)

Member Data Documentation

DQMStore* PFClient_JetRes::dqmStore_
private

Definition at line 37 of file PFClient_JetRes.h.

Referenced by beginJob(), createEfficiencyPlots(), and createResolutionPlots().

std::vector<std::string> PFClient_JetRes::effHistogramNames_
private

Definition at line 32 of file PFClient_JetRes.h.

Referenced by doEfficiency(), and PFClient_JetRes().

bool PFClient_JetRes::efficiencyFlag_
private

Definition at line 35 of file PFClient_JetRes.h.

Referenced by endRun(), and PFClient_JetRes().

std::vector<std::string> PFClient_JetRes::folderNames_
private

Definition at line 30 of file PFClient_JetRes.h.

Referenced by doEfficiency(), doSummaries(), and PFClient_JetRes().

std::vector<std::string> PFClient_JetRes::histogramNames_
private

Definition at line 31 of file PFClient_JetRes.h.

Referenced by doSummaries(), and PFClient_JetRes().

std::vector<int> PFClient_JetRes::PtBins_
private

Definition at line 33 of file PFClient_JetRes.h.

Referenced by createResolutionPlots(), and PFClient_JetRes().