Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-11-17 01:40:49

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     void analyze(const edm::Event &, const edm::EventSetup &) override;
0054     void endJob() override;
0055 
0056   private:
0057     edm::ParameterSet conf_;
0058     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
0059     SiPixelGainCalibrationService SiPixelGainCalibrationService_;
0060 
0061     std::map<uint32_t, TH1F *> _TH1F_Pedestals_m;
0062     std::map<uint32_t, TH1F *> _TH1F_Gains_m;
0063     std::map<uint32_t, double> _deadfrac_m;
0064     std::map<uint32_t, double> _noisyfrac_m;
0065 
0066     TH1F *_TH1F_Dead_sum;
0067     TH1F *_TH1F_Noisy_sum;
0068     TH1F *_TH1F_Gains_sum;
0069     TH1F *_TH1F_Pedestals_sum;
0070     TH1F *_TH1F_Dead_all;
0071     TH1F *_TH1F_Noisy_all;
0072     TH1F *_TH1F_Gains_all;
0073     TH1F *_TH1F_Pedestals_all;
0074     TH1F *_TH1F_Gains_bpix;
0075     TH1F *_TH1F_Gains_fpix;
0076     TH1F *_TH1F_Pedestals_bpix;
0077     TH1F *_TH1F_Pedestals_fpix;
0078   };
0079 }  // namespace cms
0080 
0081 namespace cms {
0082   SiPixelCondObjReader::SiPixelCondObjReader(const edm::ParameterSet &conf)
0083       : conf_(conf), tkGeomToken_(esConsumes()), SiPixelGainCalibrationService_(conf, consumesCollector()) {
0084     usesResource(TFileService::kSharedResource);
0085   }
0086 
0087   void SiPixelCondObjReader::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0088     //Create Subdirectories
0089     edm::Service<TFileService> fs;
0090     TFileDirectory subDirPed = fs->mkdir("Pedestals");
0091     TFileDirectory subDirGain = fs->mkdir("Gains");
0092     char name[128];
0093 
0094     unsigned int nmodules = 0;
0095     uint32_t nchannels = 0;
0096     uint32_t ndead = 0;
0097     uint32_t nnoisy = 0;
0098 
0099     // Get the calibration data
0100     SiPixelGainCalibrationService_.setESObjects(iSetup);
0101     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::beginJob] End Reading CondObjects" << std::endl;
0102 
0103     // Get the Geometry
0104     const TrackerGeometry *tkgeom = &iSetup.getData(tkGeomToken_);
0105     edm::LogInfo("SiPixelCondObjReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
0106 
0107     // Get the list of DetId's
0108     std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_.getDetIds();
0109 
0110     //Create histograms
0111     _TH1F_Dead_sum = fs->make<TH1F>(
0112         "Summary_dead", "Dead pixel fraction (0=dead, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0113     _TH1F_Dead_all = fs->make<TH1F>("DeadAll",
0114                                     "Dead pixel fraction (0=dead, 1=alive)",
0115                                     50,
0116                                     0.,
0117                                     conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
0118     _TH1F_Noisy_sum = fs->make<TH1F>(
0119         "Summary_noisy", "Noisy pixel fraction (0=noisy, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0120     _TH1F_Noisy_all = fs->make<TH1F>("NoisyAll",
0121                                      "Noisy pixel fraction (0=noisy, 1=alive)",
0122                                      50,
0123                                      0.,
0124                                      conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
0125     _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0126     _TH1F_Pedestals_sum =
0127         fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
0128     _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
0129     _TH1F_Pedestals_bpix = fs->make<TH1F>("PedestalsBpix", "bpix Pedestals", 350, -100, 250);
0130     _TH1F_Pedestals_fpix = fs->make<TH1F>("PedestalsFpix", "fpix Pedestals", 350, -100, 250);
0131     _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
0132     _TH1F_Gains_bpix = fs->make<TH1F>("GainsBpix", "bpix Gains", 100, 0, 10);
0133     _TH1F_Gains_fpix = fs->make<TH1F>("GainsFpix", "fpix Gains", 100, 0, 10);
0134 
0135     TTree *tree = new TTree("tree", "tree");
0136     uint32_t detid;
0137     double gainmeanfortree, gainrmsfortree, pedmeanfortree, pedrmsfortree;
0138     tree->Branch("detid", &detid, "detid/I");
0139     tree->Branch("ped_mean", &pedmeanfortree, "ped_mean/D");
0140     tree->Branch("ped_rms", &pedrmsfortree, "ped_rms/D");
0141     tree->Branch("gain_mean", &gainmeanfortree, "gain_mean/D");
0142     tree->Branch("gain_rms", &gainrmsfortree, "gain_rms/D");
0143 
0144     // Loop over DetId's
0145     int ibin = 1;
0146     for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0147          detid_iter++) {
0148       detid = *detid_iter;
0149 
0150       sprintf(name, "Pedestals_%d", detid);
0151       _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 350, -100., 250.);
0152       sprintf(name, "Gains_%d", detid);
0153       _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
0154 
0155       DetId detIdObject(detid);
0156       const PixelGeomDetUnit *_PixelGeomDetUnit =
0157           dynamic_cast<const PixelGeomDetUnit *>(tkgeom->idToDetUnit(DetId(detid)));
0158       if (_PixelGeomDetUnit == nullptr) {
0159         edm::LogError("SiPixelCondObjDisplay") << "[SiPixelCondObjReader::beginJob] the detID " << detid
0160                                                << " doesn't seem to belong to Tracker" << std::endl;
0161         continue;
0162       }
0163 
0164       _deadfrac_m[detid] = 0.;
0165       _noisyfrac_m[detid] = 0.;
0166 
0167       nmodules++;
0168 
0169       const GeomDetUnit *geoUnit = tkgeom->idToDetUnit(detIdObject);
0170       const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
0171       const PixelTopology &topol = pixDet->specificTopology();
0172 
0173       // Get the module sizes.
0174       int nrows = topol.nrows();     // rows in x
0175       int ncols = topol.ncolumns();  // cols in y
0176       float nchannelspermod = 0;
0177 
0178       for (int col_iter = 0; col_iter < ncols; col_iter++) {
0179         for (int row_iter = 0; row_iter < nrows; row_iter++) {
0180           nchannelspermod++;
0181           nchannels++;
0182 
0183           if (SiPixelGainCalibrationService_.isDead(detid, col_iter, row_iter)) {
0184             //      edm::LogPrint("SiPixelCondObjReader") << "found dead pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
0185             ndead++;
0186             _deadfrac_m[detid]++;
0187             continue;
0188           } else if (SiPixelGainCalibrationService_.isNoisy(detid, col_iter, row_iter)) {
0189             //      edm::LogPrint("SiPixelCondObjReader") << "found noisy pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
0190             nnoisy++;
0191             _noisyfrac_m[detid]++;
0192             continue;
0193           }
0194 
0195           float gain = SiPixelGainCalibrationService_.getGain(detid, col_iter, row_iter);
0196           _TH1F_Gains_m[detid]->Fill(gain);
0197           _TH1F_Gains_all->Fill(gain);
0198 
0199           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
0200             _TH1F_Gains_bpix->Fill(gain);
0201           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
0202             _TH1F_Gains_fpix->Fill(gain);
0203 
0204           float ped = SiPixelGainCalibrationService_.getPedestal(detid, col_iter, row_iter);
0205           _TH1F_Pedestals_m[detid]->Fill(ped);
0206           _TH1F_Pedestals_all->Fill(ped);
0207           //     edm::LogPrint("SiPixelCondObjReader")<<"detid  "<<detid<<"     ped "<<ped<<std::endl;
0208 
0209           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
0210             _TH1F_Pedestals_bpix->Fill(ped);
0211           if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
0212             _TH1F_Pedestals_fpix->Fill(ped);
0213 
0214           //     edm::LogPrint("SiPixelCondObjReader") <<"    DetId "<<detid<<"       Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
0215         }
0216       }
0217 
0218       _deadfrac_m[detid] /= nchannelspermod;
0219       _noisyfrac_m[detid] /= nchannelspermod;
0220       _TH1F_Dead_sum->SetBinContent(ibin, _deadfrac_m[detid]);
0221       _TH1F_Dead_all->Fill(_deadfrac_m[detid]);
0222       _TH1F_Noisy_sum->SetBinContent(ibin, _noisyfrac_m[detid]);
0223       _TH1F_Noisy_all->Fill(_noisyfrac_m[detid]);
0224       _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
0225       _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
0226       _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
0227       _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
0228 
0229       gainmeanfortree = _TH1F_Gains_m[detid]->GetMean();
0230       gainrmsfortree = _TH1F_Gains_m[detid]->GetRMS();
0231       pedmeanfortree = _TH1F_Pedestals_m[detid]->GetMean();
0232       pedrmsfortree = _TH1F_Pedestals_m[detid]->GetRMS();
0233       edm::LogPrint("SiPixelCondObjReader")
0234           << "DetId " << detid << "       GainMean " << gainmeanfortree << " RMS " << gainrmsfortree << "      PedMean "
0235           << pedmeanfortree << " RMS " << pedrmsfortree << std::endl;
0236       tree->Fill();
0237 
0238       ibin++;
0239     }
0240 
0241     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::analyze] ---> PIXEL Modules  " << nmodules
0242                                          << std::endl;
0243     edm::LogInfo("SiPixelCondObjReader") << "[SiPixelCondObjReader::analyze] ---> PIXEL Channels " << nchannels
0244                                          << std::endl;
0245 
0246     edm::LogPrint("SiPixelCondObjReader") << " ---> SUMMARY :" << std::endl;
0247     edm::LogPrint("SiPixelCondObjReader") << "Encounted " << ndead << " dead pixels" << std::endl;
0248     edm::LogPrint("SiPixelCondObjReader") << "Encounted " << nnoisy << " noisy pixels" << std::endl;
0249     edm::LogPrint("SiPixelCondObjReader")
0250         << "The Gain Mean is " << _TH1F_Gains_all->GetMean() << " with rms " << _TH1F_Gains_all->GetRMS() << std::endl;
0251     edm::LogPrint("SiPixelCondObjReader") << "         in BPIX " << _TH1F_Gains_bpix->GetMean() << " with rms "
0252                                           << _TH1F_Gains_bpix->GetRMS() << std::endl;
0253     edm::LogPrint("SiPixelCondObjReader") << "         in FPIX " << _TH1F_Gains_fpix->GetMean() << " with rms "
0254                                           << _TH1F_Gains_fpix->GetRMS() << std::endl;
0255     edm::LogPrint("SiPixelCondObjReader") << "The Ped Mean is " << _TH1F_Pedestals_all->GetMean() << " with rms "
0256                                           << _TH1F_Pedestals_all->GetRMS() << std::endl;
0257     edm::LogPrint("SiPixelCondObjReader") << "         in BPIX " << _TH1F_Pedestals_bpix->GetMean() << " with rms "
0258                                           << _TH1F_Pedestals_bpix->GetRMS() << std::endl;
0259     edm::LogPrint("SiPixelCondObjReader") << "         in FPIX " << _TH1F_Pedestals_fpix->GetMean() << " with rms "
0260                                           << _TH1F_Pedestals_fpix->GetRMS() << std::endl;
0261   }
0262 
0263   // ------------ method called once each job just after ending the event loop  ------------
0264   void SiPixelCondObjReader::endJob() { edm::LogPrint("SiPixelCondObjReader") << " ---> End job " << std::endl; }
0265 }  // namespace cms
0266 using namespace cms;
0267 DEFINE_FWK_MODULE(SiPixelCondObjReader);