CMS 3D CMS Logo

SiPixelCondObjForHLTReader.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelCondObjForHLTReader
4 // Class: SiPixelCondObjForHLTReader
5 //
13 //
14 // Original Author: Vincenzo CHIOCHIA
15 // Created: Tue Oct 17 17:40:56 CEST 2006
16 // $Id: SiPixelCondObjForHLTReader.h,v 1.7 2009/05/28 22:12:55 dlange Exp $
17 //
18 //
19 
20 // system includes
21 #include <memory>
22 
23 // user includes
41 
42 // ROOT includes
43 #include "TROOT.h"
44 #include "TFile.h"
45 #include "TTree.h"
46 #include "TBranch.h"
47 #include "TH1F.h"
48 
49 namespace cms {
50  class SiPixelCondObjForHLTReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
51  public:
52  explicit SiPixelCondObjForHLTReader(const edm::ParameterSet &iConfig);
53 
54  ~SiPixelCondObjForHLTReader() override = default;
55 
56  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
57  void analyze(const edm::Event &, const edm::EventSetup &) override;
58  void endJob() override;
59 
60  private:
62  //edm::ESHandle<SiPixelGainCalibration> SiPixelGainCalibration_;
64  std::unique_ptr<SiPixelGainCalibrationServiceBase> SiPixelGainCalibrationService_;
65 
66  std::map<uint32_t, TH1F *> _TH1F_Pedestals_m;
67  std::map<uint32_t, TH1F *> _TH1F_Gains_m;
68  std::map<uint32_t, double> _deadfrac_m;
69  std::map<uint32_t, double> _noisyfrac_m;
70 
83  };
84 } // namespace cms
85 
86 namespace cms {
88  : conf_(conf), tkGeomToken_(esConsumes()) {
89  if (conf_.getParameter<bool>("useSimRcd"))
91  std::make_unique<SiPixelGainCalibrationForHLTSimService>(conf_, consumesCollector());
92  else
94  std::make_unique<SiPixelGainCalibrationForHLTService>(conf_, consumesCollector());
95  }
96 
98  //Create Subdirectories
100  TFileDirectory subDirPed = fs->mkdir("Pedestals");
101  TFileDirectory subDirGain = fs->mkdir("Gains");
102  char name[128];
103 
104  unsigned int nmodules = 0;
105  uint32_t nchannels = 0;
106  uint32_t ndead = 0;
107  uint32_t nnoisy = 0;
108 
109  // Get the calibration data
110  SiPixelGainCalibrationService_->setESObjects(iSetup);
111  edm::LogInfo("SiPixelCondObjForHLTReader")
112  << "[SiPixelCondObjForHLTReader::beginJob] End Reading CondObjForHLTects" << std::endl;
113 
114  // Get the Geometry
115  const TrackerGeometry *tkgeom = &iSetup.getData(tkGeomToken_);
116  edm::LogInfo("SiPixelCondObjForHLTReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
117 
118  // Get the list of DetId's
119  std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_->getDetIds();
120 
121  //Create histograms
122  _TH1F_Dead_sum = fs->make<TH1F>(
123  "Summary_dead", "Dead pixel fraction (0=dead, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
124  _TH1F_Dead_all = fs->make<TH1F>("DeadAll",
125  "Dead pixel fraction (0=dead, 1=alive)",
126  50,
127  0.,
128  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
129  _TH1F_Noisy_sum = fs->make<TH1F>(
130  "Summary_noisy", "Noisy pixel fraction (0=noisy, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
131  _TH1F_Noisy_all = fs->make<TH1F>("NoisyAll",
132  "Noisy pixel fraction (0=noisy, 1=alive)",
133  50,
134  0.,
135  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
136  _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
138  fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
139  _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
140  _TH1F_Pedestals_bpix = fs->make<TH1F>("PedestalsBpix", "bpix Pedestals", 350, -100, 250);
141  _TH1F_Pedestals_fpix = fs->make<TH1F>("PedestalsFpix", "fpix Pedestals", 350, -100, 250);
142  _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
143  _TH1F_Gains_bpix = fs->make<TH1F>("GainsBpix", "bpix Gains", 100, 0, 10);
144  _TH1F_Gains_fpix = fs->make<TH1F>("GainsFpix", "fpix Gains", 100, 0, 10);
145 
146  TTree *tree = new TTree("tree", "tree");
147  uint32_t detid;
148  double gainmeanfortree, gainrmsfortree, pedmeanfortree, pedrmsfortree;
149  tree->Branch("detid", &detid, "detid/I");
150  tree->Branch("ped_mean", &pedmeanfortree, "ped_mean/D");
151  tree->Branch("ped_rms", &pedrmsfortree, "ped_rms/D");
152  tree->Branch("gain_mean", &gainmeanfortree, "gain_mean/D");
153  tree->Branch("gain_rms", &gainrmsfortree, "gain_rms/D");
154 
155  // Loop over DetId's
156  int ibin = 1;
157  for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
158  detid_iter++) {
159  detid = *detid_iter;
160 
161  sprintf(name, "Pedestals_%d", detid);
162  _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 350, -100., 250.);
163  sprintf(name, "Gains_%d", detid);
164  _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
165 
166  DetId detIdObject(detid);
167  const PixelGeomDetUnit *_PixelGeomDetUnit =
168  dynamic_cast<const PixelGeomDetUnit *>(tkgeom->idToDetUnit(DetId(detid)));
169  if (_PixelGeomDetUnit == nullptr) {
170  edm::LogError("SiPixelCondObjHLTDisplay") << "[SiPixelCondObjHLTReader::beginJob] the detID " << detid
171  << " doesn't seem to belong to Tracker" << std::endl;
172  continue;
173  }
174 
175  _deadfrac_m[detid] = 0.;
176  _noisyfrac_m[detid] = 0.;
177 
178  nmodules++;
179 
180  const GeomDetUnit *geoUnit = tkgeom->idToDetUnit(detIdObject);
181  const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
182  const PixelTopology &topol = pixDet->specificTopology();
183 
184  // Get the module sizes.
185  int nrows = topol.nrows(); // rows in x
186  int ncols = topol.ncolumns(); // cols in y
187  float nchannelspermod = 0;
188 
189  for (int col_iter = 0; col_iter < ncols; col_iter++) {
190  for (int row_iter = 0; row_iter < nrows; row_iter++) {
191  nchannelspermod++;
192  nchannels++;
193 
194  if (SiPixelGainCalibrationService_->isDead(detid, col_iter, row_iter)) {
195  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found dead pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
196  ndead++;
197  _deadfrac_m[detid]++;
198  continue;
199  } else if (SiPixelGainCalibrationService_->isNoisy(detid, col_iter, row_iter)) {
200  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found noisy pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
201  nnoisy++;
202  _noisyfrac_m[detid]++;
203  continue;
204  }
205 
206  float gain = SiPixelGainCalibrationService_->getGain(detid, col_iter, row_iter);
207  _TH1F_Gains_m[detid]->Fill(gain);
208  _TH1F_Gains_all->Fill(gain);
209 
210  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
211  _TH1F_Gains_bpix->Fill(gain);
212  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
213  _TH1F_Gains_fpix->Fill(gain);
214 
215  float ped = SiPixelGainCalibrationService_->getPedestal(detid, col_iter, row_iter);
216  _TH1F_Pedestals_m[detid]->Fill(ped);
217  _TH1F_Pedestals_all->Fill(ped);
218  // edm::LogPrint("SiPixelCondObjForHLTReader")<<"detid "<<detid<<" ped "<<ped<<std::endl;
219 
220  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
221  _TH1F_Pedestals_bpix->Fill(ped);
222  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
223  _TH1F_Pedestals_fpix->Fill(ped);
224 
225  // edm::LogPrint("SiPixelCondObjForHLTReader") <<" DetId "<<detid<<" Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
226  }
227  }
228 
229  _deadfrac_m[detid] /= nchannelspermod;
230  _noisyfrac_m[detid] /= nchannelspermod;
231  _TH1F_Dead_sum->SetBinContent(ibin, _deadfrac_m[detid]);
232  _TH1F_Dead_all->Fill(_deadfrac_m[detid]);
233  _TH1F_Noisy_sum->SetBinContent(ibin, _noisyfrac_m[detid]);
234  _TH1F_Noisy_all->Fill(_noisyfrac_m[detid]);
235  _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
236  _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
237  _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
238  _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
239 
240  gainmeanfortree = _TH1F_Gains_m[detid]->GetMean();
241  gainrmsfortree = _TH1F_Gains_m[detid]->GetRMS();
242  pedmeanfortree = _TH1F_Pedestals_m[detid]->GetMean();
243  pedrmsfortree = _TH1F_Pedestals_m[detid]->GetRMS();
244  //edm::LogPrint("SiPixelCondObjForHLTReader")<<"DetId "<<detid<<" GainMean "<<gainmeanfortree<<" RMS "<<gainrmsfortree<<" PedMean "<<pedmeanfortree<<" RMS "<<pedrmsfortree<<std::endl;
245  tree->Fill();
246 
247  ibin++;
248  }
249 
250  edm::LogInfo("SiPixelCondObjForHLTReader")
251  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Modules " << nmodules << std::endl;
252  edm::LogInfo("SiPixelCondObjForHLTReader")
253  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Channels (i.e. Number of Columns)" << nchannels
254  << std::endl;
255 
256  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> SUMMARY :" << std::endl;
257  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << ndead << " dead pixels" << std::endl;
258  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << nnoisy << " noisy pixels" << std::endl;
259  edm::LogPrint("SiPixelCondObjForHLTReader")
260  << "The Gain Mean is " << _TH1F_Gains_all->GetMean() << " with rms " << _TH1F_Gains_all->GetRMS() << std::endl;
261  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Gains_bpix->GetMean() << " with rms "
262  << _TH1F_Gains_bpix->GetRMS() << std::endl;
263  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Gains_fpix->GetMean() << " with rms "
264  << _TH1F_Gains_fpix->GetRMS() << std::endl;
265  edm::LogPrint("SiPixelCondObjForHLTReader") << "The Ped Mean is " << _TH1F_Pedestals_all->GetMean() << " with rms "
266  << _TH1F_Pedestals_all->GetRMS() << std::endl;
267  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Pedestals_bpix->GetMean()
268  << " with rms " << _TH1F_Pedestals_bpix->GetRMS() << std::endl;
269  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Pedestals_fpix->GetMean()
270  << " with rms " << _TH1F_Pedestals_fpix->GetRMS() << std::endl;
271  }
272 
273  // ------------ method called once each job just after ending the event loop ------------
275  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> End job " << std::endl;
276  }
277 
280  desc.setComment("EDAnalyzer to read per-module SiPixelGainCalibrationForHLT payloads in the EventSetup");
281  desc.add<bool>("useSimRcd", false);
282  desc.addUntracked<double>("maxRangeDeadPixHist", 0.001);
283  descriptions.addWithDefaultLabel(desc);
284  }
285 
286 } // namespace cms
287 
288 using namespace cms;
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual int ncolumns() const =0
std::map< uint32_t, double > _deadfrac_m
virtual int nrows() const =0
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Log< level::Error, false > LogError
T getUntrackedParameter(std::string const &, T const &) const
T * make(const Args &...args) const
make new ROOT object
int iEvent
Definition: GenABIO.cc:224
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
SiPixelCondObjForHLTReader(const edm::ParameterSet &iConfig)
std::map< uint32_t, TH1F * > _TH1F_Gains_m
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Log< level::Warning, true > LogPrint
void analyze(const edm::Event &, const edm::EventSetup &) override
Namespace of DDCMS conversion namespace.
Log< level::Info, false > LogInfo
Definition: DetId.h:17
std::unique_ptr< SiPixelGainCalibrationServiceBase > SiPixelGainCalibrationService_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
~SiPixelCondObjForHLTReader() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: tree.py:1
std::map< uint32_t, TH1F * > _TH1F_Pedestals_m
std::map< uint32_t, double > _noisyfrac_m