Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:23

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelCondObjReader
0004 // Class:      SiPixelCondObjReader
0005 //
0006 /**\class SiPixelCondObjReader SiPixelCondObjReader.h SiPixel/plugins/SiPixelCondObjReader.cc
0007 
0008  Description: Test analyzer for reading pixel calibration from the DB
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Vincenzo CHIOCHIA
0015 //         Created:  Tue Oct 17 17:40:56 CEST 2006
0016 // $Id: SiPixelCondObjReader.h,v 1.8 2009/05/28 22:12:55 dlange Exp $
0017 //
0018 //
0019 
0020 // system includes
0021 #include <memory>
0022 
0023 // user includes
0024 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0025 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0026 #include "FWCore/ServiceRegistry/interface/Service.h"
0027 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0028 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 #include "FWCore/Framework/interface/Frameworkfwd.h"
0031 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/EventSetup.h"
0034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0037 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0038 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationService.h"
0039 
0040 // ROOT includes
0041 #include "TROOT.h"
0042 #include "TFile.h"
0043 #include "TTree.h"
0044 #include "TBranch.h"
0045 #include "TH1F.h"
0046 
0047 namespace cms {
0048   class SiPixelCondObjReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0049   public:
0050     explicit SiPixelCondObjReader(const edm::ParameterSet &iConfig);
0051 
0052     ~SiPixelCondObjReader() override = default;
0053 
0054     static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0055     void analyze(const edm::Event &, const edm::EventSetup &) override;
0056     void endJob() override;
0057 
0058   private:
0059     edm::ParameterSet conf_;
0060     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
0061     SiPixelGainCalibrationService SiPixelGainCalibrationService_;
0062 
0063     std::map<uint32_t, TH1F *> _TH1F_Pedestals_m;
0064     std::map<uint32_t, TH1F *> _TH1F_Gains_m;
0065     std::map<uint32_t, double> _deadfrac_m;
0066     std::map<uint32_t, double> _noisyfrac_m;
0067 
0068     TH1F *_TH1F_Dead_sum;
0069     TH1F *_TH1F_Noisy_sum;
0070     TH1F *_TH1F_Gains_sum;
0071     TH1F *_TH1F_Pedestals_sum;
0072     TH1F *_TH1F_Dead_all;
0073     TH1F *_TH1F_Noisy_all;
0074     TH1F *_TH1F_Gains_all;
0075     TH1F *_TH1F_Pedestals_all;
0076     TH1F *_TH1F_Gains_bpix;
0077     TH1F *_TH1F_Gains_fpix;
0078     TH1F *_TH1F_Pedestals_bpix;
0079     TH1F *_TH1F_Pedestals_fpix;
0080   };
0081 }  // namespace cms
0082 
0083 namespace cms {
0084   SiPixelCondObjReader::SiPixelCondObjReader(const edm::ParameterSet &conf)
0085       : conf_(conf), tkGeomToken_(esConsumes()), SiPixelGainCalibrationService_(conf, consumesCollector()) {
0086     usesResource(TFileService::kSharedResource);
0087   }
0088 
0089   void SiPixelCondObjReader::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0090     //Create Subdirectories
0091     edm::Service<TFileService> fs;
0092     TFileDirectory subDirPed = fs->mkdir("Pedestals");
0093     TFileDirectory subDirGain = fs->mkdir("Gains");
0094     char name[128];
0095 
0096     unsigned int nmodules = 0;
0097     uint32_t nchannels = 0;
0098     uint32_t ndead = 0;
0099     uint32_t nnoisy = 0;
0100 
0101     // Get the calibration data
0102     SiPixelGainCalibrationService_.setESObjects(iSetup);
0103     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::beginJob] End Reading CondObjects" << std::endl;
0104 
0105     // Get the Geometry
0106     const TrackerGeometry *tkgeom = &iSetup.getData(tkGeomToken_);
0107     edm::LogInfo("SiPixelCondObjReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
0108 
0109     // Get the list of DetId's
0110     std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_.getDetIds();
0111 
0112     //Create histograms
0113     _TH1F_Dead_sum = fs->make<TH1F>(
0114         "Summary_dead", "Dead pixel fraction (0=dead, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0115     _TH1F_Dead_all = fs->make<TH1F>("DeadAll",
0116                                     "Dead pixel fraction (0=dead, 1=alive)",
0117                                     50,
0118                                     0.,
0119                                     conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
0120     _TH1F_Noisy_sum = fs->make<TH1F>(
0121         "Summary_noisy", "Noisy pixel fraction (0=noisy, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0122     _TH1F_Noisy_all = fs->make<TH1F>("NoisyAll",
0123                                      "Noisy pixel fraction (0=noisy, 1=alive)",
0124                                      50,
0125                                      0.,
0126                                      conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
0127     _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0128     _TH1F_Pedestals_sum =
0129         fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0130     _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
0131     _TH1F_Pedestals_bpix = fs->make<TH1F>("PedestalsBpix", "bpix Pedestals", 350, -100, 250);
0132     _TH1F_Pedestals_fpix = fs->make<TH1F>("PedestalsFpix", "fpix Pedestals", 350, -100, 250);
0133     _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
0134     _TH1F_Gains_bpix = fs->make<TH1F>("GainsBpix", "bpix Gains", 100, 0, 10);
0135     _TH1F_Gains_fpix = fs->make<TH1F>("GainsFpix", "fpix Gains", 100, 0, 10);
0136 
0137     TTree *tree = new TTree("tree", "tree");
0138     uint32_t detid;
0139     double gainmeanfortree, gainrmsfortree, pedmeanfortree, pedrmsfortree;
0140     tree->Branch("detid", &detid, "detid/I");
0141     tree->Branch("ped_mean", &pedmeanfortree, "ped_mean/D");
0142     tree->Branch("ped_rms", &pedrmsfortree, "ped_rms/D");
0143     tree->Branch("gain_mean", &gainmeanfortree, "gain_mean/D");
0144     tree->Branch("gain_rms", &gainrmsfortree, "gain_rms/D");
0145 
0146     // Loop over DetId's
0147     int ibin = 1;
0148     for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0149          detid_iter++) {
0150       detid = *detid_iter;
0151 
0152       sprintf(name, "Pedestals_%d", detid);
0153       _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 350, -100., 250.);
0154       sprintf(name, "Gains_%d", detid);
0155       _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
0156 
0157       DetId detIdObject(detid);
0158       const PixelGeomDetUnit *_PixelGeomDetUnit =
0159           dynamic_cast<const PixelGeomDetUnit *>(tkgeom->idToDetUnit(DetId(detid)));
0160       if (_PixelGeomDetUnit == nullptr) {
0161         edm::LogError("SiPixelCondObjDisplay") << "[SiPixelCondObjReader::beginJob] the detID " << detid
0162                                                << " doesn't seem to belong to Tracker" << std::endl;
0163         continue;
0164       }
0165 
0166       _deadfrac_m[detid] = 0.;
0167       _noisyfrac_m[detid] = 0.;
0168 
0169       nmodules++;
0170 
0171       const GeomDetUnit *geoUnit = tkgeom->idToDetUnit(detIdObject);
0172       const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
0173       const PixelTopology &topol = pixDet->specificTopology();
0174 
0175       // Get the module sizes.
0176       int nrows = topol.nrows();     // rows in x
0177       int ncols = topol.ncolumns();  // cols in y
0178       float nchannelspermod = 0;
0179 
0180       for (int col_iter = 0; col_iter < ncols; col_iter++) {
0181         for (int row_iter = 0; row_iter < nrows; row_iter++) {
0182           nchannelspermod++;
0183           nchannels++;
0184 
0185           if (SiPixelGainCalibrationService_.isDead(detid, col_iter, row_iter)) {
0186             //      edm::LogPrint("SiPixelCondObjReader") << "found dead pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
0187             ndead++;
0188             _deadfrac_m[detid]++;
0189             continue;
0190           } else if (SiPixelGainCalibrationService_.isNoisy(detid, col_iter, row_iter)) {
0191             //      edm::LogPrint("SiPixelCondObjReader") << "found noisy pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
0192             nnoisy++;
0193             _noisyfrac_m[detid]++;
0194             continue;
0195           }
0196 
0197           float gain = SiPixelGainCalibrationService_.getGain(detid, col_iter, row_iter);
0198           _TH1F_Gains_m[detid]->Fill(gain);
0199           _TH1F_Gains_all->Fill(gain);
0200 
0201           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
0202             _TH1F_Gains_bpix->Fill(gain);
0203           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
0204             _TH1F_Gains_fpix->Fill(gain);
0205 
0206           float ped = SiPixelGainCalibrationService_.getPedestal(detid, col_iter, row_iter);
0207           _TH1F_Pedestals_m[detid]->Fill(ped);
0208           _TH1F_Pedestals_all->Fill(ped);
0209           //     edm::LogPrint("SiPixelCondObjReader")<<"detid  "<<detid<<"     ped "<<ped<<std::endl;
0210 
0211           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
0212             _TH1F_Pedestals_bpix->Fill(ped);
0213           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
0214             _TH1F_Pedestals_fpix->Fill(ped);
0215 
0216           //     edm::LogPrint("SiPixelCondObjReader") <<"    DetId "<<detid<<"       Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
0217         }
0218       }
0219 
0220       _deadfrac_m[detid] /= nchannelspermod;
0221       _noisyfrac_m[detid] /= nchannelspermod;
0222       _TH1F_Dead_sum->SetBinContent(ibin, _deadfrac_m[detid]);
0223       _TH1F_Dead_all->Fill(_deadfrac_m[detid]);
0224       _TH1F_Noisy_sum->SetBinContent(ibin, _noisyfrac_m[detid]);
0225       _TH1F_Noisy_all->Fill(_noisyfrac_m[detid]);
0226       _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
0227       _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
0228       _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
0229       _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
0230 
0231       gainmeanfortree = _TH1F_Gains_m[detid]->GetMean();
0232       gainrmsfortree = _TH1F_Gains_m[detid]->GetRMS();
0233       pedmeanfortree = _TH1F_Pedestals_m[detid]->GetMean();
0234       pedrmsfortree = _TH1F_Pedestals_m[detid]->GetRMS();
0235       edm::LogPrint("SiPixelCondObjReader")
0236           << "DetId " << detid << "       GainMean " << gainmeanfortree << " RMS " << gainrmsfortree << "      PedMean "
0237           << pedmeanfortree << " RMS " << pedrmsfortree << std::endl;
0238       tree->Fill();
0239 
0240       ibin++;
0241     }
0242 
0243     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::analyze] ---> PIXEL Modules  " << nmodules
0244                                          << std::endl;
0245     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::analyze] ---> PIXEL Channels " << nchannels
0246                                          << std::endl;
0247 
0248     edm::LogPrint("SiPixelCondObjReader") << " ---> SUMMARY :" << std::endl;
0249     edm::LogPrint("SiPixelCondObjReader") << "Encounted " << ndead << " dead pixels" << std::endl;
0250     edm::LogPrint("SiPixelCondObjReader") << "Encounted " << nnoisy << " noisy pixels" << std::endl;
0251     edm::LogPrint("SiPixelCondObjReader")
0252         << "The Gain Mean is " << _TH1F_Gains_all->GetMean() << " with rms " << _TH1F_Gains_all->GetRMS() << std::endl;
0253     edm::LogPrint("SiPixelCondObjReader") << "         in BPIX " << _TH1F_Gains_bpix->GetMean() << " with rms "
0254                                           << _TH1F_Gains_bpix->GetRMS() << std::endl;
0255     edm::LogPrint("SiPixelCondObjReader") << "         in FPIX " << _TH1F_Gains_fpix->GetMean() << " with rms "
0256                                           << _TH1F_Gains_fpix->GetRMS() << std::endl;
0257     edm::LogPrint("SiPixelCondObjReader") << "The Ped Mean is " << _TH1F_Pedestals_all->GetMean() << " with rms "
0258                                           << _TH1F_Pedestals_all->GetRMS() << std::endl;
0259     edm::LogPrint("SiPixelCondObjReader") << "         in BPIX " << _TH1F_Pedestals_bpix->GetMean() << " with rms "
0260                                           << _TH1F_Pedestals_bpix->GetRMS() << std::endl;
0261     edm::LogPrint("SiPixelCondObjReader") << "         in FPIX " << _TH1F_Pedestals_fpix->GetMean() << " with rms "
0262                                           << _TH1F_Pedestals_fpix->GetRMS() << std::endl;
0263   }
0264 
0265   // ------------ method called once each job just after ending the event loop  ------------
0266   void SiPixelCondObjReader::endJob() { edm::LogPrint("SiPixelCondObjReader") << " ---> End job " << std::endl; }
0267 
0268   void SiPixelCondObjReader::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0269     edm::ParameterSetDescription desc;
0270     desc.setComment("EDAnalyzer to read per-module SiPixelGainCalibration payloads in the EventSetup");
0271     desc.addUntracked<double>("maxRangeDeadPixHist", 0.001);
0272     descriptions.addWithDefaultLabel(desc);
0273   }
0274 
0275 }  // namespace cms
0276 using namespace cms;
0277 DEFINE_FWK_MODULE(SiPixelCondObjReader);