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
PhotonPostprocessing Class Reference

#include <PhotonPostprocessing.h>

Inheritance diagram for PhotonPostprocessing:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
virtual void endLuminosityBlock (const edm::LuminosityBlock &, const edm::EventSetup &)
 
virtual void endRun (const edm::Run &, const edm::EventSetup &)
 
 PhotonPostprocessing (const edm::ParameterSet &pset)
 
virtual ~PhotonPostprocessing ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

void dividePlots (MonitorElement *dividend, MonitorElement *numerator, MonitorElement *denominator, std::string type)
 
void dividePlots (MonitorElement *dividend, MonitorElement *numerator, double denominator)
 
virtual void runPostprocessing ()
 

Private Attributes

bool batch_
 
MonitorElementbkgDeadChEt_
 
MonitorElementbkgDeadChEta_
 
MonitorElementbkgDeadChPhi_
 
MonitorElementbkgRecoEffEt_
 
MonitorElementbkgRecoEffEta_
 
MonitorElementbkgRecoEffPhi_
 
MonitorElementconvEffEtaOneTrack_
 
MonitorElementconvEffEtaTwoTracks_
 
MonitorElementconvEffEtaTwoTracksAndVtxProbGT0005_
 
MonitorElementconvEffEtaTwoTracksAndVtxProbGT0_
 
MonitorElementconvEffEtOneTrack_
 
MonitorElementconvEffEtTwoTracks_
 
MonitorElementconvEffPhiTwoTracks_
 
MonitorElementconvEffROneTrack_
 
MonitorElementconvEffRTwoTracks_
 
MonitorElementconvEffRTwoTracksAndVtxProbGT0005_
 
MonitorElementconvEffRTwoTracksAndVtxProbGT0_
 
MonitorElementconvEffZTwoTracks_
 
MonitorElementconvFakeRateEtaTwoTracks_
 
MonitorElementconvFakeRateEtTwoTracks_
 
MonitorElementconvFakeRatePhiTwoTracks_
 
MonitorElementconvFakeRateRTwoTracks_
 
MonitorElementconvFakeRateZTwoTracks_
 
MonitorElementconvVsEt_ [2]
 
std::stringstream currentFolder_
 
DQMStoredbe_
 
int etaBin
 
int etaBin2
 
double etaMax
 
double etaMin
 
int etBin
 
double etMax
 
double etMin
 
bool fastSim_
 
std::string inputFileName_
 
bool isRunCentrally_
 
std::string outputFileName_
 
edm::ParameterSet parameters_
 
int phiBin
 
double phiMax
 
double phiMin
 
MonitorElementphoDeadChEt_
 
MonitorElementphoDeadChEta_
 
MonitorElementphoDeadChPhi_
 
MonitorElementphoRecoEffEt_
 
MonitorElementphoRecoEffEta_
 
MonitorElementphoRecoEffPhi_
 
int rBin
 
double rMax
 
double rMin
 
bool standAlone_
 
int verbosity_
 
int zBin
 
double zMax
 
double zMin
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

$Id: PhotonPostprocessing

Date:
2012/03/22 15:21:21

author: Nancy Marinelli, U. of Notre Dame, US

$Id: PhotonPostprocessing

Date:
2012/03/22 15:21:24

authors: Nancy Marinelli, U. of Notre Dame, US

Definition at line 53 of file PhotonPostprocessing.h.

Constructor & Destructor Documentation

PhotonPostprocessing::PhotonPostprocessing ( const edm::ParameterSet pset)
explicit

Definition at line 26 of file PhotonPostprocessing.cc.

References dbe_, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, jptDQMConfig_cff::etMax, reco::tau::qcuts::etMin(), edm::ParameterSet::getParameter(), cppFunctionSkipper::operator, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, and DQMStore::setVerbose().

27 {
28 
29  dbe_ = 0;
31  dbe_->setVerbose(0);
32  parameters_ = pset;
33 
34 
35  standAlone_ = pset.getParameter<bool>("standAlone");
36  batch_ = pset.getParameter<bool>("batch");
37  outputFileName_ = pset.getParameter<string>("OutputFileName");
38  inputFileName_ = pset.getParameter<std::string>("InputFileName");
39  isRunCentrally_= pset.getParameter<bool>("isRunCentrally");
40  fastSim_ = pset.getParameter<bool>("fastSim");
41 
42  etMin = parameters_.getParameter<double>("etMin");
43  etMax = parameters_.getParameter<double>("etMax");
44  etBin = parameters_.getParameter<int>("etBin");
45 
46 
47  etaMin = parameters_.getParameter<double>("etaMin");
48  etaMax = parameters_.getParameter<double>("etaMax");
49  etaBin = parameters_.getParameter<int>("etaBin");
50  etaBin2 = parameters_.getParameter<int>("etaBin2");
51 
52  phiMin = parameters_.getParameter<double>("phiMin");
53  phiMax = parameters_.getParameter<double>("phiMax");
54  phiBin = parameters_.getParameter<int>("phiBin");
55 
56  rMin = parameters_.getParameter<double>("rMin");
57  rMax = parameters_.getParameter<double>("rMax");
58  rBin = parameters_.getParameter<int>("rBin");
59 
60  zMin = parameters_.getParameter<double>("zMin");
61  zMax = parameters_.getParameter<double>("zMax");
62  zBin = parameters_.getParameter<int>("zBin");
63 
64 
65 
66 }
T getParameter(std::string const &) const
edm::ParameterSet parameters_
void setVerbose(unsigned level)
Definition: DQMStore.cc:393
PhotonPostprocessing::~PhotonPostprocessing ( )
virtual

Definition at line 70 of file PhotonPostprocessing.cc.

71 {}

Member Function Documentation

void PhotonPostprocessing::analyze ( const edm::Event e,
const edm::EventSetup esup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 78 of file PhotonPostprocessing.cc.

79 {}
void PhotonPostprocessing::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 73 of file PhotonPostprocessing.cc.

74 {
75 
76 }
void PhotonPostprocessing::dividePlots ( MonitorElement dividend,
MonitorElement numerator,
MonitorElement denominator,
std::string  type 
)
private

Definition at line 251 of file PhotonPostprocessing.cc.

References MonitorElement::getBinContent(), MonitorElement::getNbinsX(), j, MonitorElement::setBinContent(), MonitorElement::setBinError(), mathSSE::sqrt(), and relativeConstraints::value.

251  {
252  double value,err;
253  for (int j=1; j<=numerator->getNbinsX(); j++){
254 
255  if (denominator->getBinContent(j)!=0){
256  if (type=="effic")
257  value = ((double) numerator->getBinContent(j))/((double) denominator->getBinContent(j));
258  else if (type=="fakerate")
259  value = 1-((double) numerator->getBinContent(j))/((double) denominator->getBinContent(j));
260  else return;
261  err = sqrt( value*(1-value) / ((double) denominator->getBinContent(j)) );
262  dividend->setBinContent(j, value);
263  if ( err !=0 ) dividend->setBinError(j,err);
264  }
265  else {
266  dividend->setBinContent(j, 0);
267  dividend->setBinError(j,0);
268  }
269 
270  }
271 
272 
273 }
type
Definition: HCALResponse.h:22
void setBinContent(int binx, double content)
set content of bin (1-D)
T sqrt(T t)
Definition: SSEVec.h:46
int j
Definition: DBlmapReader.cc:9
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
double getBinContent(int binx) const
get content of bin (1-D)
int getNbinsX(void) const
get # of bins in X-axis
void PhotonPostprocessing::dividePlots ( MonitorElement dividend,
MonitorElement numerator,
double  denominator 
)
private

Definition at line 276 of file PhotonPostprocessing.cc.

References MonitorElement::getBinContent(), MonitorElement::getNbinsX(), j, MonitorElement::setBinContent(), MonitorElement::setBinError(), mathSSE::sqrt(), and relativeConstraints::value.

276  {
277  double value,err;
278 
279  for (int j=1; j<=numerator->getNbinsX(); j++){
280  if (denominator!=0){
281  value = ((double) numerator->getBinContent(j))/denominator;
282  err = sqrt( value*(1-value) / denominator);
283  dividend->setBinContent(j, value);
284  dividend->setBinError(j,err);
285  }
286  else {
287  dividend->setBinContent(j, 0);
288  }
289  }
290 
291 }
void setBinContent(int binx, double content)
set content of bin (1-D)
T sqrt(T t)
Definition: SSEVec.h:46
int j
Definition: DBlmapReader.cc:9
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
double getBinContent(int binx) const
get content of bin (1-D)
int getNbinsX(void) const
get # of bins in X-axis
void PhotonPostprocessing::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 82 of file PhotonPostprocessing.cc.

82  {
83 
85 
86 }
virtual void runPostprocessing()
void PhotonPostprocessing::endLuminosityBlock ( const edm::LuminosityBlock lumi,
const edm::EventSetup setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 243 of file PhotonPostprocessing.cc.

244 {
245 
246 
247 }
void PhotonPostprocessing::endRun ( const edm::Run run,
const edm::EventSetup setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 88 of file PhotonPostprocessing.cc.

88  {
89 
91 
92 }
virtual void runPostprocessing()
void PhotonPostprocessing::runPostprocessing ( )
privatevirtual

Definition at line 95 of file PhotonPostprocessing.cc.

References DQMStore::book1D(), dbe_, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, jptDQMConfig_cff::etMax, reco::tau::qcuts::etMin(), DQMStore::get(), DQMStore::open(), jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, DQMStore::save(), and DQMStore::setCurrentFolder().

96 {
97 
98 
99  std::string simInfoPathName = "EgammaV/PhotonValidator/SimulationInfo/";
100  std::string convPathName = "EgammaV/PhotonValidator/ConversionInfo/";
101  std::string effPathName = "EgammaV/PhotonValidator/Efficiencies/";
102  std::string photonPathName = "EgammaV/PhotonValidator/Photons/";
103 
105 
106  dbe_->setCurrentFolder(effPathName);
107  // Photon reconstruction efficiencies
108  string histname = "recoEffVsEta";
109  phoRecoEffEta_ = dbe_->book1D(histname,"Photon reconstruction efficiency vs simulated #eta",etaBin,etaMin, etaMax);
110  histname = "recoEffVsPhi";
111  phoRecoEffPhi_ = dbe_->book1D(histname,"Photon reconstruction efficiency vs simulated #phi",phiBin,phiMin, phiMax);
112  histname = "recoEffVsEt";
113  phoRecoEffEt_ = dbe_->book1D(histname,"Photon reconstruction efficiency vs simulated Et",etBin,etMin, etMax) ;
114  // Fraction of photons with at least one dead channel
115  histname = "deadChVsEta";
116  phoDeadChEta_ = dbe_->book1D(histname,"Fraction of photons with >=1 dead Xtal vs simulated #eta",etaBin,etaMin, etaMax);
117  histname = "deadChVsPhi";
118  phoDeadChPhi_ = dbe_->book1D(histname,"Fraction of photons with >=1 dead Xtal vs simulated #phi",phiBin,phiMin, phiMax);
119  histname = "deadChVsEt";
120  phoDeadChEt_ = dbe_->book1D(histname,"Fraction of photons with >=1 dead Xtal vs simulated Et",etBin,etMin, etMax) ;
121 
122  if ( ! isRunCentrally_ ) {
123  histname = "convVsEt";
124  convVsEt_[0] = dbe_->book1D(histname+"Barrel","Fraction of good conversions in R9<0.93 vs Et ",etBin,etMin, etMax) ;
125  convVsEt_[1] = dbe_->book1D(histname+"Endcap","Fraction of good conversions in R9<0.93 vs Et ",etBin,etMin, etMax) ;
126  }
127 
128 
129 
130  // Conversion reconstruction efficiency
131  histname = "convEffVsEtaTwoTracks";
132  convEffEtaTwoTracks_ = dbe_->book1D(histname,histname,etaBin2,etaMin, etaMax);
133 
134  histname = "convEffVsPhiTwoTracks";
135  convEffPhiTwoTracks_ = dbe_->book1D(histname,histname,phiBin,phiMin,phiMax);
136 
137  histname = "convEffVsRTwoTracks";
138  convEffRTwoTracks_ = dbe_->book1D(histname,histname,rBin,rMin, rMax);
139 
140  histname = "convEffVsZTwoTracks";
141  convEffZTwoTracks_ = dbe_->book1D(histname,histname,zBin,zMin,zMax);
142 
143  histname = "convEffVsEtTwoTracks";
144  convEffEtTwoTracks_ = dbe_->book1D(histname,histname,etBin,etMin, etMax);
145  //
146  histname = "convEffVsEtaTwoTracksAndVtxProbGT0";
148  histname = "convEffVsEtaTwoTracksAndVtxProbGT0005";
150  histname = "convEffVsRTwoTracksAndVtxProbGT0";
152  histname = "convEffVsRTwoTracksAndVtxProbGT0005";
154  //
155  histname = "convEffVsEtaOneTrack";
156  convEffEtaOneTrack_ = dbe_->book1D(histname,histname,etaBin2,etaMin, etaMax);
157  histname = "convEffVsROneTrack";
158  convEffROneTrack_ = dbe_->book1D(histname,histname,rBin,rMin, rMax);
159  histname = "convEffVsEtOneTrack";
160  convEffEtOneTrack_ = dbe_->book1D(histname,histname,etBin,etMin, etMax);
161  // Fake rate
162  histname = "convFakeRateVsEtaTwoTracks";
163  convFakeRateEtaTwoTracks_ = dbe_->book1D(histname,histname,etaBin2,etaMin, etaMax);
164  histname = "convFakeRateVsPhiTwoTracks";
166  histname = "convFakeRateVsRTwoTracks";
167  convFakeRateRTwoTracks_ = dbe_->book1D(histname,histname,rBin,rMin, rMax);
168  histname = "convFakeRateVsZTwoTracks";
169  convFakeRateZTwoTracks_ = dbe_->book1D(histname,histname,zBin,zMin,zMax);
170  histname = "convFakeRateVsEtTwoTracks";
171  convFakeRateEtTwoTracks_ = dbe_->book1D(histname,histname,etBin,etMin, etMax);
172 
173  histname = "bkgEffVsEta";
174  bkgRecoEffEta_ = dbe_->book1D(histname,"Bkg reconstruction efficiency vs simulated #eta",etaBin,etaMin, etaMax);
175  histname = "bkgEffVsPhi";
176  bkgRecoEffPhi_ = dbe_->book1D(histname,"Bkg reconstruction efficiency vs simulated #phi",phiBin,phiMin, phiMax);
177  histname = "bkgEffVsEt";
178  bkgRecoEffEt_ = dbe_->book1D(histname,"Bkg reconstruction efficiency vs simulated Et",etBin,etMin, etMax) ;
179  // Fraction of photons with at least one dead channel
180  histname = "deadChVsEtaBkg";
181  bkgDeadChEta_ = dbe_->book1D(histname,"Fraction of bkg with >=1 dead Xtal vs simulated #eta",etaBin,etaMin, etaMax);
182  histname = "deadChVsPhiBkg";
183  bkgDeadChPhi_ = dbe_->book1D(histname,"Fraction of bkg with >=1 dead Xtal vs simulated #phi",phiBin,phiMin, phiMax);
184  histname = "deadChVsEtBkg";
185  bkgDeadChEt_ = dbe_->book1D(histname,"Fraction of bkg with >=1 dead Xtal vs simulated Et",etBin,etMin, etMax) ;
186 
187 
188 
189  // efficiencies
190  if ( ! isRunCentrally_ ) {
191  dividePlots(dbe_->get(effPathName+"convVsEtBarrel"),dbe_->get(photonPathName+"EtR9Less093ConvBarrel"),dbe_->get(photonPathName+"EtR9Less093Barrel"), "effic");
192  dividePlots(dbe_->get(effPathName+"convVsEtEndcap"),dbe_->get(photonPathName+"EtR9Less093ConvEndcap"),dbe_->get(photonPathName+"EtR9Less093Endcap"), "effic");
193  }
194 
195  dividePlots(dbe_->get(effPathName+"recoEffVsEta"),dbe_->get(simInfoPathName+"h_MatchedSimPhoEta"),dbe_->get(simInfoPathName+"h_SimPhoEta"), "effic");
196  dividePlots(dbe_->get(effPathName+"recoEffVsPhi"),dbe_->get(simInfoPathName+"h_MatchedSimPhoPhi"),dbe_->get(simInfoPathName+"h_SimPhoPhi"),"effic");
197  dividePlots(dbe_->get(effPathName+"recoEffVsEt"),dbe_->get(simInfoPathName+"h_MatchedSimPhoEt"),dbe_->get(simInfoPathName+"h_SimPhoEt"),"effic");
198  // fraction of photons with at least one dead channel
199  dividePlots(dbe_->get(effPathName+"deadChVsEta"),dbe_->get(simInfoPathName+"h_MatchedSimPhoBadChEta"),dbe_->get(simInfoPathName+"h_MatchedSimPhoEta"), "effic");
200  dividePlots(dbe_->get(effPathName+"deadChVsPhi"),dbe_->get(simInfoPathName+"h_MatchedSimPhoBadChPhi"),dbe_->get(simInfoPathName+"h_MatchedSimPhoPhi"),"effic");
201  dividePlots(dbe_->get(effPathName+"deadChVsEt"), dbe_->get(simInfoPathName+"h_MatchedSimPhoBadChEt"),dbe_->get(simInfoPathName+"h_MatchedSimPhoEt"),"effic");
202  //
203  if ( ! fastSim_ ) {
204  dividePlots(dbe_->get(effPathName+"convEffVsEtaTwoTracks"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEta"),dbe_->get(simInfoPathName+"h_VisSimConvEta"),"effic");
205  dividePlots(dbe_->get(effPathName+"convEffVsPhiTwoTracks"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksPhi"),dbe_->get(simInfoPathName+"h_VisSimConvPhi"),"effic");
206  dividePlots(dbe_->get(effPathName+"convEffVsRTwoTracks"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksR"),dbe_->get(simInfoPathName+"h_VisSimConvR"),"effic");
207  dividePlots(dbe_->get(effPathName+"convEffVsZTwoTracks"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksZ"),dbe_->get(simInfoPathName+"h_VisSimConvZ"),"effic");
208  dividePlots(dbe_->get(effPathName+"convEffVsEtTwoTracks"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEt"),dbe_->get(simInfoPathName+"h_VisSimConvEt"),"effic");
209  dividePlots(dbe_->get(effPathName+"convEffVsEtaTwoTracksAndVtxProbGT0"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEtaAndVtxPGT0"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEta"),"effic");
210  dividePlots(dbe_->get(effPathName+"convEffVsEtaTwoTracksAndVtxProbGT0005"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEtaAndVtxPGT0005"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksEta"),"effic");
211  dividePlots(dbe_->get(effPathName+"convEffVsRTwoTracksAndVtxProbGT0"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksRAndVtxPGT0"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksR"),"effic");
212  dividePlots(dbe_->get(effPathName+"convEffVsRTwoTracksAndVtxProbGT0005"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksRAndVtxPGT0005"),dbe_->get(simInfoPathName+"h_SimConvTwoMTracksR"),"effic");
213  //
214  dividePlots(dbe_->get(effPathName+"convEffVsEtaOneTrack"),dbe_->get(simInfoPathName+"h_SimConvOneMTracksEta"),dbe_->get(simInfoPathName+"h_VisSimConvEta"),"effic");
215  dividePlots(dbe_->get(effPathName+"convEffVsROneTrack"),dbe_->get(simInfoPathName+"h_SimConvOneMTracksR"),dbe_->get(simInfoPathName+"h_VisSimConvR"),"effic");
216  dividePlots(dbe_->get(effPathName+"convEffVsEtOneTrack"),dbe_->get(simInfoPathName+"h_SimConvOneMTracksEt"),dbe_->get(simInfoPathName+"h_VisSimConvEt"),"effic");
217  // fake rate
218  dividePlots(dbe_->get(effPathName+"convFakeRateVsEtaTwoTracks"),dbe_->get(convPathName+"h_RecoConvTwoMTracksEta"),dbe_->get(convPathName+"h_RecoConvTwoTracksEta"),"fakerate");
219  dividePlots(dbe_->get(effPathName+"convFakeRateVsPhiTwoTracks"),dbe_->get(convPathName+"h_RecoConvTwoMTracksPhi"),dbe_->get(convPathName+"h_RecoConvTwoTracksPhi"),"fakerate");
220  dividePlots(dbe_->get(effPathName+"convFakeRateVsRTwoTracks"),dbe_->get(convPathName+"h_RecoConvTwoMTracksR"),dbe_->get(convPathName+"h_RecoConvTwoTracksR"),"fakerate");
221  dividePlots(dbe_->get(effPathName+"convFakeRateVsZTwoTracks"),dbe_->get(convPathName+"h_RecoConvTwoMTracksZ"),dbe_->get(convPathName+"h_RecoConvTwoTracksZ"),"fakerate");
222  dividePlots(dbe_->get(effPathName+"convFakeRateVsEtTwoTracks"),dbe_->get(convPathName+"h_RecoConvTwoMTracksEt"),dbe_->get(convPathName+"h_RecoConvTwoTracksEt"),"fakerate");
223  }
224  // Background efficiency
225  dividePlots(dbe_->get(effPathName+"bkgEffVsEta"),dbe_->get(simInfoPathName+"h_MatchedSimJetEta"),dbe_->get(simInfoPathName+"h_SimJetEta"), "effic");
226  dividePlots(dbe_->get(effPathName+"bkgEffVsPhi"),dbe_->get(simInfoPathName+"h_MatchedSimJetPhi"),dbe_->get(simInfoPathName+"h_SimJetPhi"),"effic");
227  dividePlots(dbe_->get(effPathName+"bkgEffVsEt"),dbe_->get(simInfoPathName+"h_MatchedSimJetEt"),dbe_->get(simInfoPathName+"h_SimJetEt"),"effic");
228  // fraction of photons with at least one dead channel
229  dividePlots(dbe_->get(effPathName+"deadChVsEtaBkg"),dbe_->get(simInfoPathName+"h_MatchedSimJetBadChEta"),dbe_->get(simInfoPathName+"h_MatchedSimJetEta"), "effic");
230  dividePlots(dbe_->get(effPathName+"deadChVsPhiBkg"),dbe_->get(simInfoPathName+"h_MatchedSimJetBadChPhi"),dbe_->get(simInfoPathName+"h_MatchedSimJetPhi"),"effic");
231  dividePlots(dbe_->get(effPathName+"deadChVsEtBkg"), dbe_->get(simInfoPathName+"h_MatchedSimJetBadChEt"),dbe_->get(simInfoPathName+"h_MatchedSimJetEt"),"effic");
232 
233 
234 
236  else if(batch_) dbe_->save(inputFileName_);
237 
238 
239 
240 }
MonitorElement * convEffEtaTwoTracksAndVtxProbGT0005_
MonitorElement * convEffEtOneTrack_
MonitorElement * convEffRTwoTracksAndVtxProbGT0_
MonitorElement * bkgRecoEffPhi_
MonitorElement * phoRecoEffEta_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
MonitorElement * phoRecoEffPhi_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2113
MonitorElement * phoDeadChEta_
MonitorElement * convEffEtaOneTrack_
MonitorElement * phoDeadChPhi_
MonitorElement * bkgDeadChEta_
MonitorElement * convEffRTwoTracks_
MonitorElement * bkgRecoEffEt_
MonitorElement * phoRecoEffEt_
MonitorElement * convVsEt_[2]
MonitorElement * bkgDeadChEt_
MonitorElement * phoDeadChEt_
MonitorElement * convFakeRateZTwoTracks_
MonitorElement * convFakeRateEtaTwoTracks_
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
MonitorElement * convEffPhiTwoTracks_
MonitorElement * convFakeRateEtTwoTracks_
MonitorElement * bkgDeadChPhi_
MonitorElement * convEffRTwoTracksAndVtxProbGT0005_
MonitorElement * convEffEtaTwoTracksAndVtxProbGT0_
MonitorElement * convFakeRatePhiTwoTracks_
MonitorElement * convEffROneTrack_
MonitorElement * convEffEtaTwoTracks_
MonitorElement * bkgRecoEffEta_
void dividePlots(MonitorElement *dividend, MonitorElement *numerator, MonitorElement *denominator, std::string type)
bool open(const std::string &filename, bool overwrite=false, const std::string &path="", const std::string &prepend="", OpenRunDirs stripdirs=KeepRunDirs, bool fileMustExist=true)
Definition: DQMStore.cc:2432
MonitorElement * convFakeRateRTwoTracks_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
MonitorElement * convEffZTwoTracks_
MonitorElement * convEffEtTwoTracks_

Member Data Documentation

bool PhotonPostprocessing::batch_
private

Definition at line 86 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgDeadChEt_
private

Definition at line 150 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgDeadChEta_
private

Definition at line 148 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgDeadChPhi_
private

Definition at line 149 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgRecoEffEt_
private

Definition at line 146 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgRecoEffEta_
private

Definition at line 144 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::bkgRecoEffPhi_
private

Definition at line 145 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtaOneTrack_
private

Definition at line 134 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtaTwoTracks_
private

Definition at line 123 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtaTwoTracksAndVtxProbGT0005_
private

Definition at line 130 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtaTwoTracksAndVtxProbGT0_
private

Definition at line 129 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtOneTrack_
private

Definition at line 136 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffEtTwoTracks_
private

Definition at line 127 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffPhiTwoTracks_
private

Definition at line 124 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffROneTrack_
private

Definition at line 135 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffRTwoTracks_
private

Definition at line 125 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffRTwoTracksAndVtxProbGT0005_
private

Definition at line 132 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffRTwoTracksAndVtxProbGT0_
private

Definition at line 131 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convEffZTwoTracks_
private

Definition at line 126 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convFakeRateEtaTwoTracks_
private

Definition at line 138 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convFakeRateEtTwoTracks_
private

Definition at line 142 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convFakeRatePhiTwoTracks_
private

Definition at line 139 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convFakeRateRTwoTracks_
private

Definition at line 140 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convFakeRateZTwoTracks_
private

Definition at line 141 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::convVsEt_[2]
private

Definition at line 152 of file PhotonPostprocessing.h.

std::stringstream PhotonPostprocessing::currentFolder_
private

Definition at line 92 of file PhotonPostprocessing.h.

DQMStore* PhotonPostprocessing::dbe_
private

Definition at line 79 of file PhotonPostprocessing.h.

int PhotonPostprocessing::etaBin
private

Definition at line 100 of file PhotonPostprocessing.h.

int PhotonPostprocessing::etaBin2
private

Definition at line 101 of file PhotonPostprocessing.h.

double PhotonPostprocessing::etaMax
private

Definition at line 99 of file PhotonPostprocessing.h.

double PhotonPostprocessing::etaMin
private

Definition at line 98 of file PhotonPostprocessing.h.

int PhotonPostprocessing::etBin
private

Definition at line 97 of file PhotonPostprocessing.h.

double PhotonPostprocessing::etMax
private

Definition at line 96 of file PhotonPostprocessing.h.

double PhotonPostprocessing::etMin
private

Definition at line 95 of file PhotonPostprocessing.h.

bool PhotonPostprocessing::fastSim_
private

Definition at line 88 of file PhotonPostprocessing.h.

std::string PhotonPostprocessing::inputFileName_
private

Definition at line 90 of file PhotonPostprocessing.h.

bool PhotonPostprocessing::isRunCentrally_
private

Definition at line 87 of file PhotonPostprocessing.h.

std::string PhotonPostprocessing::outputFileName_
private

Definition at line 89 of file PhotonPostprocessing.h.

edm::ParameterSet PhotonPostprocessing::parameters_
private
int PhotonPostprocessing::phiBin
private

Definition at line 104 of file PhotonPostprocessing.h.

double PhotonPostprocessing::phiMax
private

Definition at line 103 of file PhotonPostprocessing.h.

double PhotonPostprocessing::phiMin
private

Definition at line 102 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoDeadChEt_
private

Definition at line 120 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoDeadChEta_
private

Definition at line 118 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoDeadChPhi_
private

Definition at line 119 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoRecoEffEt_
private

Definition at line 116 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoRecoEffEta_
private

Definition at line 114 of file PhotonPostprocessing.h.

MonitorElement* PhotonPostprocessing::phoRecoEffPhi_
private

Definition at line 115 of file PhotonPostprocessing.h.

int PhotonPostprocessing::rBin
private

Definition at line 107 of file PhotonPostprocessing.h.

double PhotonPostprocessing::rMax
private

Definition at line 106 of file PhotonPostprocessing.h.

double PhotonPostprocessing::rMin
private

Definition at line 105 of file PhotonPostprocessing.h.

bool PhotonPostprocessing::standAlone_
private

Definition at line 85 of file PhotonPostprocessing.h.

int PhotonPostprocessing::verbosity_
private

Definition at line 80 of file PhotonPostprocessing.h.

int PhotonPostprocessing::zBin
private

Definition at line 110 of file PhotonPostprocessing.h.

double PhotonPostprocessing::zMax
private

Definition at line 109 of file PhotonPostprocessing.h.

double PhotonPostprocessing::zMin
private

Definition at line 108 of file PhotonPostprocessing.h.