Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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