CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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:
53 
54  void analyze(const edm::Event &, const edm::EventSetup &) override;
55  void endJob() override;
56 
57  private:
59  //edm::ESHandle<SiPixelGainCalibration> SiPixelGainCalibration_;
61  std::unique_ptr<SiPixelGainCalibrationServiceBase> SiPixelGainCalibrationService_;
62 
63  std::map<uint32_t, TH1F *> _TH1F_Pedestals_m;
64  std::map<uint32_t, TH1F *> _TH1F_Gains_m;
65  std::map<uint32_t, double> _deadfrac_m;
66  std::map<uint32_t, double> _noisyfrac_m;
67 
80  };
81 } // namespace cms
82 
83 namespace cms {
85  : conf_(conf), tkGeomToken_(esConsumes()) {
86  if (conf_.getParameter<bool>("useSimRcd"))
88  std::make_unique<SiPixelGainCalibrationForHLTSimService>(conf_, consumesCollector());
89  else
91  std::make_unique<SiPixelGainCalibrationForHLTService>(conf_, consumesCollector());
92  }
93 
95  //Create Subdirectories
97  TFileDirectory subDirPed = fs->mkdir("Pedestals");
98  TFileDirectory subDirGain = fs->mkdir("Gains");
99  char name[128];
100 
101  unsigned int nmodules = 0;
102  uint32_t nchannels = 0;
103  uint32_t ndead = 0;
104  uint32_t nnoisy = 0;
105 
106  // Get the calibration data
107  SiPixelGainCalibrationService_->setESObjects(iSetup);
108  edm::LogInfo("SiPixelCondObjForHLTReader")
109  << "[SiPixelCondObjForHLTReader::beginJob] End Reading CondObjForHLTects" << std::endl;
110 
111  // Get the Geometry
112  const TrackerGeometry *tkgeom = &iSetup.getData(tkGeomToken_);
113  edm::LogInfo("SiPixelCondObjForHLTReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
114 
115  // Get the list of DetId's
116  std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_->getDetIds();
117 
118  //Create histograms
119  _TH1F_Dead_sum = fs->make<TH1F>(
120  "Summary_dead", "Dead pixel fraction (0=dead, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
121  _TH1F_Dead_all = fs->make<TH1F>("DeadAll",
122  "Dead pixel fraction (0=dead, 1=alive)",
123  50,
124  0.,
125  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
126  _TH1F_Noisy_sum = fs->make<TH1F>(
127  "Summary_noisy", "Noisy pixel fraction (0=noisy, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
128  _TH1F_Noisy_all = fs->make<TH1F>("NoisyAll",
129  "Noisy pixel fraction (0=noisy, 1=alive)",
130  50,
131  0.,
132  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
133  _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
135  fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
136  _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
137  _TH1F_Pedestals_bpix = fs->make<TH1F>("PedestalsBpix", "bpix Pedestals", 350, -100, 250);
138  _TH1F_Pedestals_fpix = fs->make<TH1F>("PedestalsFpix", "fpix Pedestals", 350, -100, 250);
139  _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
140  _TH1F_Gains_bpix = fs->make<TH1F>("GainsBpix", "bpix Gains", 100, 0, 10);
141  _TH1F_Gains_fpix = fs->make<TH1F>("GainsFpix", "fpix Gains", 100, 0, 10);
142 
143  TTree *tree = new TTree("tree", "tree");
144  uint32_t detid;
145  double gainmeanfortree, gainrmsfortree, pedmeanfortree, pedrmsfortree;
146  tree->Branch("detid", &detid, "detid/I");
147  tree->Branch("ped_mean", &pedmeanfortree, "ped_mean/D");
148  tree->Branch("ped_rms", &pedrmsfortree, "ped_rms/D");
149  tree->Branch("gain_mean", &gainmeanfortree, "gain_mean/D");
150  tree->Branch("gain_rms", &gainrmsfortree, "gain_rms/D");
151 
152  // Loop over DetId's
153  int ibin = 1;
154  for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
155  detid_iter++) {
156  detid = *detid_iter;
157 
158  sprintf(name, "Pedestals_%d", detid);
159  _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 350, -100., 250.);
160  sprintf(name, "Gains_%d", detid);
161  _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
162 
163  DetId detIdObject(detid);
164  const PixelGeomDetUnit *_PixelGeomDetUnit =
165  dynamic_cast<const PixelGeomDetUnit *>(tkgeom->idToDetUnit(DetId(detid)));
166  if (_PixelGeomDetUnit == nullptr) {
167  edm::LogError("SiPixelCondObjHLTDisplay") << "[SiPixelCondObjHLTReader::beginJob] the detID " << detid
168  << " doesn't seem to belong to Tracker" << std::endl;
169  continue;
170  }
171 
172  _deadfrac_m[detid] = 0.;
173  _noisyfrac_m[detid] = 0.;
174 
175  nmodules++;
176 
177  const GeomDetUnit *geoUnit = tkgeom->idToDetUnit(detIdObject);
178  const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
179  const PixelTopology &topol = pixDet->specificTopology();
180 
181  // Get the module sizes.
182  int nrows = topol.nrows(); // rows in x
183  int ncols = topol.ncolumns(); // cols in y
184  float nchannelspermod = 0;
185 
186  for (int col_iter = 0; col_iter < ncols; col_iter++) {
187  for (int row_iter = 0; row_iter < nrows; row_iter++) {
188  nchannelspermod++;
189  nchannels++;
190 
191  if (SiPixelGainCalibrationService_->isDead(detid, col_iter, row_iter)) {
192  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found dead pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
193  ndead++;
194  _deadfrac_m[detid]++;
195  continue;
196  } else if (SiPixelGainCalibrationService_->isNoisy(detid, col_iter, row_iter)) {
197  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found noisy pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
198  nnoisy++;
199  _noisyfrac_m[detid]++;
200  continue;
201  }
202 
203  float gain = SiPixelGainCalibrationService_->getGain(detid, col_iter, row_iter);
204  _TH1F_Gains_m[detid]->Fill(gain);
205  _TH1F_Gains_all->Fill(gain);
206 
207  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
208  _TH1F_Gains_bpix->Fill(gain);
209  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
210  _TH1F_Gains_fpix->Fill(gain);
211 
212  float ped = SiPixelGainCalibrationService_->getPedestal(detid, col_iter, row_iter);
213  _TH1F_Pedestals_m[detid]->Fill(ped);
214  _TH1F_Pedestals_all->Fill(ped);
215  // edm::LogPrint("SiPixelCondObjForHLTReader")<<"detid "<<detid<<" ped "<<ped<<std::endl;
216 
217  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
218  _TH1F_Pedestals_bpix->Fill(ped);
219  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
220  _TH1F_Pedestals_fpix->Fill(ped);
221 
222  // edm::LogPrint("SiPixelCondObjForHLTReader") <<" DetId "<<detid<<" Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
223  }
224  }
225 
226  _deadfrac_m[detid] /= nchannelspermod;
227  _noisyfrac_m[detid] /= nchannelspermod;
228  _TH1F_Dead_sum->SetBinContent(ibin, _deadfrac_m[detid]);
229  _TH1F_Dead_all->Fill(_deadfrac_m[detid]);
230  _TH1F_Noisy_sum->SetBinContent(ibin, _noisyfrac_m[detid]);
231  _TH1F_Noisy_all->Fill(_noisyfrac_m[detid]);
232  _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
233  _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
234  _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
235  _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
236 
237  gainmeanfortree = _TH1F_Gains_m[detid]->GetMean();
238  gainrmsfortree = _TH1F_Gains_m[detid]->GetRMS();
239  pedmeanfortree = _TH1F_Pedestals_m[detid]->GetMean();
240  pedrmsfortree = _TH1F_Pedestals_m[detid]->GetRMS();
241  //edm::LogPrint("SiPixelCondObjForHLTReader")<<"DetId "<<detid<<" GainMean "<<gainmeanfortree<<" RMS "<<gainrmsfortree<<" PedMean "<<pedmeanfortree<<" RMS "<<pedrmsfortree<<std::endl;
242  tree->Fill();
243 
244  ibin++;
245  }
246 
247  edm::LogInfo("SiPixelCondObjForHLTReader")
248  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Modules " << nmodules << std::endl;
249  edm::LogInfo("SiPixelCondObjForHLTReader")
250  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Channels (i.e. Number of Columns)" << nchannels
251  << std::endl;
252 
253  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> SUMMARY :" << std::endl;
254  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << ndead << " dead pixels" << std::endl;
255  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << nnoisy << " noisy pixels" << std::endl;
256  edm::LogPrint("SiPixelCondObjForHLTReader")
257  << "The Gain Mean is " << _TH1F_Gains_all->GetMean() << " with rms " << _TH1F_Gains_all->GetRMS() << std::endl;
258  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Gains_bpix->GetMean() << " with rms "
259  << _TH1F_Gains_bpix->GetRMS() << std::endl;
260  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Gains_fpix->GetMean() << " with rms "
261  << _TH1F_Gains_fpix->GetRMS() << std::endl;
262  edm::LogPrint("SiPixelCondObjForHLTReader") << "The Ped Mean is " << _TH1F_Pedestals_all->GetMean() << " with rms "
263  << _TH1F_Pedestals_all->GetRMS() << std::endl;
264  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Pedestals_bpix->GetMean()
265  << " with rms " << _TH1F_Pedestals_bpix->GetRMS() << std::endl;
266  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Pedestals_fpix->GetMean()
267  << " with rms " << _TH1F_Pedestals_fpix->GetRMS() << std::endl;
268  }
269 
270  // ------------ method called once each job just after ending the event loop ------------
272  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> End job " << std::endl;
273  }
274 } // namespace cms
275 
276 using namespace cms;
T getUntrackedParameter(std::string const &, T const &) const
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual int ncolumns() const =0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< uint32_t, double > _deadfrac_m
virtual int nrows() const =0
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Log< level::Error, false > LogError
bool getData(T &iHolder) const
Definition: EventSetup.h:128
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
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
T * make(const Args &...args) const
make new ROOT object
Log< level::Warning, true > LogPrint
void analyze(const edm::Event &, const edm::EventSetup &) override
Log< level::Info, false > LogInfo
Definition: DetId.h:17
std::unique_ptr< SiPixelGainCalibrationServiceBase > SiPixelGainCalibrationService_
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::map< uint32_t, TH1F * > _TH1F_Pedestals_m
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::map< uint32_t, double > _noisyfrac_m