CMS 3D CMS Logo

SiPixelPhase1EfficiencyExtras.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1EfficiencyExtras
4 // Class: SiPixelPhase1EfficiencyExtras
5 //
13 //
14 // Original Author: Jack Sisson, Julie Hogan
15 // Created: 7 July, 2021
16 //
17 //
18 
19 // Framework
27 
28 // DQM Framework
31 
32 #include <algorithm>
33 
35 public:
38 
39 protected:
40  void beginRun(edm::Run const& run, edm::EventSetup const& eSetup) override;
41 
42  void dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) override;
43 
47 };
48 
50  : effFolderName_(iConfig.getParameter<std::string>("EffFolderName")),
51  vtxFolderName_(iConfig.getParameter<std::string>("VtxFolderName")),
52  instLumiFolderName_(iConfig.getParameter<std::string>("InstLumiFolderName")) {
53  edm::LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras: Hello!";
54 }
55 
57  edm::LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras: Destructor";
58 }
59 
61 
62 //------------------------------------------------------------------
63 // Method called for every event
64 //------------------------------------------------------------------
67 
68  //Get the existing histos
69  MonitorElement* vtx_v_lumi = iGetter.get(vtxFolderName_ + "/NumberOfGoodPVtxVsLS_GenTk");
70 
71  MonitorElement* scalLumi_v_lumi = iGetter.get(instLumiFolderName_ + "/lumiVsLS");
72 
73  MonitorElement* eff_v_lumi_forward =
74  iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXDisk_PXForward");
75 
76  MonitorElement* eff_v_lumi_barrel =
77  iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXLayer_PXBarrel");
78 
79  //set up some booleans that will tell us which graphs to create
80  bool createNvtx = true;
81  bool createInstLumi = true;
82 
83  //check which of the MEs exist and respond appropriately
84  if (!eff_v_lumi_forward) {
85  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
86  << "no hitefficiency_per_Lumisection_per_PXDisk_PXForward ME is available in " << effFolderName_;
87  return;
88  }
89  if (!eff_v_lumi_barrel) {
90  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
91  << "no hitefficiency_per_Lumisection_per_PXLayer_PXBarrel ME is available in " << effFolderName_;
92  return;
93  }
94  if (!vtx_v_lumi) {
95  edm::LogWarning("SiPixelPhase1EfficiencyExtras")
96  << "no NumberOfGoodPVtxVsLS_GenTK ME is available in " << vtxFolderName_;
97  createNvtx = false;
98  }
99  if (!scalLumi_v_lumi) {
100  edm::LogWarning("SiPixelPhase1EfficiencyExtras") << "no lumiVsLS ME is available in " << instLumiFolderName_;
101  createInstLumi = false;
102  }
103 
104  //If the existing MEs are empty, set the boolean to skip booking
105  if (vtx_v_lumi && vtx_v_lumi->getEntries() == 0)
106  createNvtx = false;
107  if (scalLumi_v_lumi && scalLumi_v_lumi->getEntries() == 0)
108  createInstLumi = false;
109 
110  // if the max mean lumi is not higher than zero, do not create profiles with respect to lumi
111  if (createInstLumi and scalLumi_v_lumi->getTProfile()->GetMaximum() <= 0.)
112  createInstLumi = false;
113 
114  double eff = 0.0;
115 
116  //Will pass if nvtx ME exists and is not empty
117  if (createNvtx) {
118  //Book new histos
119  MonitorElement* eff_v_vtx_barrel =
120  iBooker.book2D("hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel",
121  "hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel; meanNvtx; PXLayer",
122  500,
123  0,
124  100,
125  3,
126  .5,
127  3.5);
128 
129  MonitorElement* eff_v_vtx_forward =
130  iBooker.book2D("hitefficiency_per_meanNvtx_per_PXDisk_PXForward",
131  "hitefficiency_per_meanNvtx_per_PXDisk_PXForward; meanNvtx; PXDisk",
132  500,
133  0,
134  100,
135  7,
136  -3.5,
137  3.5);
138 
139  //initialize variables
140  int numLumiNvtx = int(vtx_v_lumi->getNbinsX());
141  double nvtx = 0.0;
142  int binNumVtx = 0;
143 
144  //For loop to loop through lumisections
145  for (int iLumi = 1; iLumi < numLumiNvtx - 1; iLumi++) {
146  //get the meanNvtx for each lumi
147  nvtx = vtx_v_lumi->getBinContent(iLumi);
148 
149  //Filter out useless iterations
150  if (nvtx != 0) {
151  //Grab the bin number for the nvtx
152  binNumVtx = eff_v_vtx_barrel->getTH2F()->FindBin(nvtx);
153 
154  //loop through the layers
155  for (int iLayer = 1; iLayer < 8; iLayer++) {
156  //get the eff at the lumisection and layer
157  eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
158 
159  //set the efficiency in the new histo
160  eff_v_vtx_forward->setBinContent(binNumVtx, iLayer, eff);
161  }
162 
163  //loop through the layers
164  for (int iLayer = 1; iLayer < 5; iLayer++) {
165  //get the efficiency for each lumi at each layer
166  eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
167 
168  //set the efficiency
169  eff_v_vtx_barrel->setBinContent(binNumVtx, iLayer, eff);
170  }
171  }
172  }
173  }
174  // Will pass if InstLumi ME exists, is not empty, and max mean lumi is larger than zero
175  if (createInstLumi) {
176  //Get the max value of inst lumi for plot (ensuring yMax2 is larger than zero)
177  int yMax2 = std::max(1., scalLumi_v_lumi->getTProfile()->GetMaximum());
178  yMax2 *= 1.1;
179 
180  //Book new histos
181  MonitorElement* eff_v_scalLumi_barrel =
182  iBooker.book2D("hitefficiency_per_scalLumi_per_PXLayer_PXBarrel",
183  "hitefficiency_per_scalLumi_per_PXLayer_PXBarrel; scal inst lumi E30; PXLayer",
184  500,
185  0,
186  yMax2,
187  3,
188  .5,
189  3.5);
190 
191  MonitorElement* eff_v_scalLumi_forward =
192  iBooker.book2D("hitefficiency_per_scalLumi_per_PXDisk_PXForward",
193  "hitefficiency_per_scalLumi_per_PXDisk_PXForward; scal inst lumi E30; PXDisk",
194  500,
195  0,
196  yMax2,
197  7,
198  -3.5,
199  3.5);
200 
201  //initialize variables
202  int numLumiScal = int(scalLumi_v_lumi->getNbinsX());
203  double scalLumi = 0.0;
204  int binNumScal = 0;
205 
206  //For loop to loop through lumisections
207  for (int iLumi = 1; iLumi < numLumiScal - 1; iLumi++) {
208  //get the inst lumi for each lumi
209  scalLumi = scalLumi_v_lumi->getBinContent(iLumi);
210 
211  //Filter out useless iterations
212  if (scalLumi > 0) {
213  //Grab the bin number for the inst lumi
214  binNumScal = eff_v_scalLumi_barrel->getTH2F()->FindBin(scalLumi);
215 
216  //loop through the layers
217  for (int iLayer = 1; iLayer < 8; iLayer++) {
218  //get the eff at the lumisection and layer
219  eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
220 
221  //set the efficiency in the new histo
222  eff_v_scalLumi_forward->setBinContent(binNumScal, iLayer, eff);
223  }
224 
225  //loop through the layers
226  for (int iLayer = 1; iLayer < 5; iLayer++) {
227  //get the eff at the lumisection and layer
228  eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
229 
230  //set the efficiency in the new histo
231  eff_v_scalLumi_barrel->setBinContent(binNumScal, iLayer, eff);
232  }
233  }
234  }
235  } else
236  return;
237 }
238 
239 //define this as a plug-in
virtual TProfile * getTProfile() const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual TH2F * getTH2F() const
void dqmEndJob(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter) override
virtual double getEntries() const
get # of entries
Log< level::Info, false > LogInfo
SiPixelPhase1EfficiencyExtras(const edm::ParameterSet &conf)
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
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:212
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:680
virtual int getNbinsX() const
get # of bins in X-axis
Log< level::Warning, false > LogWarning
void beginRun(edm::Run const &run, edm::EventSetup const &eSetup) override
Definition: Run.h:45
virtual double getBinContent(int binx) const
get content of bin (1-D)