CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
SiPixelPhase1EfficiencyExtras Class Reference
Inheritance diagram for SiPixelPhase1EfficiencyExtras:

Public Member Functions

 SiPixelPhase1EfficiencyExtras (const edm::ParameterSet &conf)
 
 ~SiPixelPhase1EfficiencyExtras () override
 

Protected Member Functions

void beginRun (edm::Run const &run, edm::EventSetup const &eSetup) override
 
void dqmEndJob (DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter) override
 

Protected Attributes

const std::string effFolderName_
 
const std::string instLumiFolderName_
 
const std::string vtxFolderName_
 

Detailed Description

Definition at line 35 of file SiPixelPhase1EfficiencyExtras.cc.

Constructor & Destructor Documentation

SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras ( const edm::ParameterSet conf)
explicit

Definition at line 50 of file SiPixelPhase1EfficiencyExtras.cc.

51  : effFolderName_(iConfig.getParameter<std::string>("EffFolderName")),
52  vtxFolderName_(iConfig.getParameter<std::string>("VtxFolderName")),
53  instLumiFolderName_(iConfig.getParameter<std::string>("InstLumiFolderName")) {
54  LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras: Hello!" << endl;
55 }
Log< level::Info, false > LogInfo
SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras ( )
override

Definition at line 57 of file SiPixelPhase1EfficiencyExtras.cc.

57  {
58  LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras: Destructor" << endl;
59 }
Log< level::Info, false > LogInfo

Member Function Documentation

void SiPixelPhase1EfficiencyExtras::beginRun ( edm::Run const &  run,
edm::EventSetup const &  eSetup 
)
overrideprotected

Definition at line 61 of file SiPixelPhase1EfficiencyExtras.cc.

61 {}
void SiPixelPhase1EfficiencyExtras::dqmEndJob ( DQMStore::IBooker iBooker,
DQMStore::IGetter iGetter 
)
overrideprotected

Definition at line 66 of file SiPixelPhase1EfficiencyExtras.cc.

References dqm::implementation::IBooker::book2D(), effFolderName_, dqm::implementation::IGetter::get(), dqm::impl::MonitorElement::getBinContent(), dqm::impl::MonitorElement::getEntries(), dqm::impl::MonitorElement::getNbinsX(), dqm::legacy::MonitorElement::getTH2F(), dqm::legacy::MonitorElement::getTProfile(), instLumiFolderName_, dqm::impl::MonitorElement::setBinContent(), dqm::implementation::NavigatorBase::setCurrentFolder(), and vtxFolderName_.

66  {
68 
69  //Get the existing histos
70  MonitorElement* vtx_v_lumi = iGetter.get(vtxFolderName_ + "/NumberOfGoodPVtxVsLS_GenTk");
71 
72  MonitorElement* scalLumi_v_lumi = iGetter.get(instLumiFolderName_ + "/lumiVsLS");
73 
74  MonitorElement* eff_v_lumi_forward =
75  iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXDisk_PXForward");
76 
77  MonitorElement* eff_v_lumi_barrel =
78  iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXLayer_PXBarrel");
79 
80  //set up some booleans that will tell us which graphs to create
81  bool createNvtx = true;
82  bool createInstLumi = true;
83 
84  //check which of the MEs exist and respond appropriately
85  if (!eff_v_lumi_forward) {
86  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
87  << "no hitefficiency_per_Lumisection_per_PXDisk_PXForward ME is available in " << effFolderName_ << std::endl;
88  return;
89  }
90  if (!eff_v_lumi_barrel) {
91  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
92  << "no hitefficiency_per_Lumisection_per_PXLayer_PXBarrel ME is available in " << effFolderName_ << std::endl;
93  return;
94  }
95  if (!vtx_v_lumi) {
96  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
97  << "no NumberOfGoodPVtxVsLS_GenTK ME is available in " << vtxFolderName_ << std::endl;
98  createNvtx = false;
99  }
100  if (!scalLumi_v_lumi) {
101  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
102  << "no lumiVsLS ME is available in " << instLumiFolderName_ << std::endl;
103  createInstLumi = false;
104  }
105 
106  //If the existing MEs are empty, set the boolean to skip booking
107  if (vtx_v_lumi && vtx_v_lumi->getEntries() == 0)
108  createNvtx = false;
109  if (scalLumi_v_lumi && scalLumi_v_lumi->getEntries() == 0)
110  createInstLumi = false;
111 
112  double eff = 0.0;
113 
114  //Will pass if nvtx ME exists and is not empty
115  if (createNvtx) {
116  //Book new histos
117  MonitorElement* eff_v_vtx_barrel =
118  iBooker.book2D("hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel",
119  "hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel; meanNvtx; PXLayer",
120  500,
121  0,
122  100,
123  3,
124  .5,
125  3.5);
126 
127  MonitorElement* eff_v_vtx_forward =
128  iBooker.book2D("hitefficiency_per_meanNvtx_per_PXDisk_PXForward",
129  "hitefficiency_per_meanNvtx_per_PXDisk_PXForward; meanNvtx; PXDisk",
130  500,
131  0,
132  100,
133  7,
134  -3.5,
135  3.5);
136 
137  //initialize variables
138  int numLumiNvtx = int(vtx_v_lumi->getNbinsX());
139  double nvtx = 0.0;
140  int binNumVtx = 0;
141 
142  //For loop to loop through lumisections
143  for (int iLumi = 1; iLumi < numLumiNvtx - 1; iLumi++) {
144  //get the meanNvtx for each lumi
145  nvtx = vtx_v_lumi->getBinContent(iLumi);
146 
147  //Filter out useless iterations
148  if (nvtx != 0) {
149  //Grab the bin number for the nvtx
150  binNumVtx = eff_v_vtx_barrel->getTH2F()->FindBin(nvtx);
151 
152  //loop through the layers
153  for (int iLayer = 1; iLayer < 8; iLayer++) {
154  //get the eff at the lumisection and layer
155  eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
156 
157  //set the efficiency in the new histo
158  eff_v_vtx_forward->setBinContent(binNumVtx, iLayer, eff);
159  }
160 
161  //loop through the layers
162  for (int iLayer = 1; iLayer < 5; iLayer++) {
163  //get the efficiency for each lumi at each layer
164  eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
165 
166  //set the efficiency
167  eff_v_vtx_barrel->setBinContent(binNumVtx, iLayer, eff);
168  }
169  }
170  }
171  }
172  // Will pass if InstLumi ME exists and is not empty
173  if (createInstLumi) {
174  //Get the max value of inst lumi for plot
175  int yMax2 = scalLumi_v_lumi->getTProfile()->GetMaximum();
176  yMax2 = yMax2 + yMax2 * .1;
177 
178  //Book new histos
179  MonitorElement* eff_v_scalLumi_barrel =
180  iBooker.book2D("hitefficiency_per_scalLumi_per_PXLayer_PXBarrel",
181  "hitefficiency_per_scalLumi_per_PXLayer_PXBarrel; scal inst lumi E30; PXLayer",
182  500,
183  0,
184  yMax2,
185  3,
186  .5,
187  3.5);
188 
189  MonitorElement* eff_v_scalLumi_forward =
190  iBooker.book2D("hitefficiency_per_scalLumi_per_PXDisk_PXForward",
191  "hitefficiency_per_scalLumi_per_PXDisk_PXForward; scal inst lumi E30; PXDisk",
192  500,
193  0,
194  yMax2,
195  7,
196  -3.5,
197  3.5);
198 
199  //initialize variables
200  int numLumiScal = int(scalLumi_v_lumi->getNbinsX());
201  double scalLumi = 0.0;
202  int binNumScal = 0;
203 
204  //For loop to loop through lumisections
205  for (int iLumi = 1; iLumi < numLumiScal - 1; iLumi++) {
206  //get the inst lumi for each lumi
207  scalLumi = scalLumi_v_lumi->getBinContent(iLumi);
208 
209  //Filter out useless iterations
210  if (scalLumi != 0) {
211  //Grab the bin number for the inst lumi
212  binNumScal = eff_v_scalLumi_barrel->getTH2F()->FindBin(scalLumi);
213 
214  //loop through the layers
215  for (int iLayer = 1; iLayer < 8; iLayer++) {
216  //get the eff at the lumisection and layer
217  eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
218 
219  //set the efficiency in the new histo
220  eff_v_scalLumi_forward->setBinContent(binNumScal, iLayer, eff);
221  }
222 
223  //loop through the layers
224  for (int iLayer = 1; iLayer < 5; iLayer++) {
225  //get the eff at the lumisection and layer
226  eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
227 
228  //set the efficiency in the new histo
229  eff_v_scalLumi_barrel->setBinContent(binNumScal, iLayer, eff);
230  }
231  }
232  }
233  } else
234  return;
235 }
virtual TH2F * getTH2F() const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
virtual int getNbinsX() const
get # of bins in X-axis
virtual double getEntries() const
get # of entries
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
virtual double getBinContent(int binx) const
get content of bin (1-D)
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual TProfile * getTProfile() const
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
Log< level::Warning, false > LogWarning

Member Data Documentation

const std::string SiPixelPhase1EfficiencyExtras::effFolderName_
protected

Definition at line 45 of file SiPixelPhase1EfficiencyExtras.cc.

Referenced by dqmEndJob().

const std::string SiPixelPhase1EfficiencyExtras::instLumiFolderName_
protected

Definition at line 47 of file SiPixelPhase1EfficiencyExtras.cc.

Referenced by dqmEndJob().

const std::string SiPixelPhase1EfficiencyExtras::vtxFolderName_
protected

Definition at line 46 of file SiPixelPhase1EfficiencyExtras.cc.

Referenced by dqmEndJob().