Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-17 04:57:59

0001 // -*- C++ -*-
0002 //
0003 // Package:    CalibTracker/SiPixelESProducers
0004 // Class:      SiPixelFakeGainForHLTReader
0005 //
0006 /**\class SiPixelFakeGainForHLTReader SiPixelFakeGainForHLTReader.cc SiPixelESProducers/test/SiPixelFakeGainForHLTReader.h
0007 
0008  Description: Test analyzer for fake pixel calibrationForHLT
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 //
0017 //
0018 
0019 // system includes
0020 #include <memory>
0021 
0022 // user includes
0023 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationForHLTService.h"
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/MakerMacros.h"
0026 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0028 #include "FWCore/PluginManager/interface/ModuleDef.h"
0029 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0030 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0031 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0032 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0033 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0035 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0036 
0037 // ROOT file
0038 #include "TROOT.h"
0039 #include "TFile.h"
0040 #include "TTree.h"
0041 #include "TBranch.h"
0042 #include "TH1F.h"
0043 
0044 namespace cms {
0045   class SiPixelFakeGainForHLTReader : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0046   public:
0047     explicit SiPixelFakeGainForHLTReader(const edm::ParameterSet& iConfig);
0048     ~SiPixelFakeGainForHLTReader() override = default;
0049 
0050     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0051 
0052     virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
0053     virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
0054     virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
0055     virtual void endJob() override;
0056 
0057   private:
0058     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
0059     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomTokenBeginRun_;
0060     SiPixelGainCalibrationForHLTService SiPixelGainCalibrationForHLTService_;
0061 
0062     std::map<uint32_t, TH1F*> _TH1F_Pedestals_m;
0063     std::map<uint32_t, TH1F*> _TH1F_Gains_m;
0064     std::string filename_;
0065     TFile* fFile;
0066   };
0067 }  // namespace cms
0068 
0069 namespace cms {
0070   SiPixelFakeGainForHLTReader::SiPixelFakeGainForHLTReader(const edm::ParameterSet& conf)
0071       : trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0072         trackerGeomTokenBeginRun_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0073         SiPixelGainCalibrationForHLTService_(conf, consumesCollector()),
0074         filename_(conf.getParameter<std::string>("fileName")) {}
0075 
0076   void SiPixelFakeGainForHLTReader::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0077     unsigned int nmodules = 0;
0078     uint32_t nchannels = 0;
0079     fFile->cd();
0080 
0081     // Get the Geometry
0082     edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomToken_);
0083     edm::LogInfo("SiPixelFakeGainForHLTReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
0084 
0085     //  for(TrackerGeometry::DetContainer::const_iterator it = tkgeom->dets().begin(); it != tkgeom->dets().end(); it++){
0086     //   if( dynamic_cast<PixelGeomDetUnit*>((*it))!=0){
0087     //     uint32_t detid=((*it)->geographicalId()).rawId();
0088     // Get the list of DetId's
0089 
0090     std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationForHLTService_.getDetIds();
0091     // Loop over DetId's
0092     for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0093          detid_iter++) {
0094       uint32_t detid = *detid_iter;
0095 
0096       DetId detIdObject(detid);
0097       nmodules++;
0098       //if(nmodules>3) break;
0099 
0100       std::map<uint32_t, TH1F*>::iterator p_iter = _TH1F_Pedestals_m.find(detid);
0101       std::map<uint32_t, TH1F*>::iterator g_iter = _TH1F_Gains_m.find(detid);
0102 
0103       const GeomDetUnit* geoUnit = tkgeom->idToDetUnit(detIdObject);
0104       const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
0105       const PixelTopology& topol = pixDet->specificTopology();
0106 
0107       // Get the module sizes.
0108       int nrows = topol.nrows();     // rows in x
0109       int ncols = topol.ncolumns();  // cols in y
0110 
0111       for (int col_iter = 0; col_iter < ncols; col_iter++) {
0112         for (int row_iter = 0; row_iter < nrows; row_iter++) {
0113           nchannels++;
0114           float ped = SiPixelGainCalibrationForHLTService_.getPedestal(detid, col_iter, row_iter);
0115           float gain = SiPixelGainCalibrationForHLTService_.getGain(detid, col_iter, row_iter);
0116           p_iter->second->Fill(ped);
0117           g_iter->second->Fill(gain);
0118         }
0119       }
0120     }
0121 
0122     edm::LogInfo("SiPixelFakeGainForHLTReader")
0123         << "[SiPixelFakeGainForHLTReader::analyze] ---> PIXEL Modules  " << nmodules << std::endl;
0124     edm::LogInfo("SiPixelFakeGainForHLTReader")
0125         << "[SiPixelFakeGainForHLTReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
0126     fFile->ls();
0127     fFile->Write();
0128     fFile->Close();
0129   }
0130 
0131   void SiPixelFakeGainForHLTReader::endRun(const edm::Run& run, const edm::EventSetup& iSetup) {}
0132 
0133   // ----------------------------------------------------------------------
0134   void SiPixelFakeGainForHLTReader::beginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0135     static int first(1);
0136     if (1 == first) {
0137       first = 0;
0138       edm::LogInfo("SiPixelFakeGainForHLTReader")
0139           << "[SiPixelFakeGainForHLTReader::beginJob] Opening ROOT file  " << std::endl;
0140       fFile = new TFile(filename_.c_str(), "RECREATE");
0141       fFile->mkdir("Pedestals");
0142       fFile->mkdir("Gains");
0143       fFile->cd();
0144       char name[128];
0145 
0146       // Get Geometry
0147       edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomTokenBeginRun_);
0148 
0149       // Get the calibrationForHLT data
0150       SiPixelGainCalibrationForHLTService_.setESObjects(iSetup);
0151       edm::LogInfo("SiPixelFakeGainForHLTReader")
0152           << "[SiPixelFakeGainForHLTReader::beginJob] End Reading FakeGainForHLTects" << std::endl;
0153       // Get the list of DetId's
0154       std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationForHLTService_.getDetIds();
0155       // Loop over DetId's
0156       for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0157            detid_iter++) {
0158         uint32_t detid = *detid_iter;
0159 
0160         const PixelGeomDetUnit* _PixelGeomDetUnit =
0161             dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
0162         if (_PixelGeomDetUnit == 0) {
0163           edm::LogError("SiPixelFakeGainForHLTDisplay") << "[SiPixelFakeGainForHLTReader::beginJob] the detID " << detid
0164                                                         << " doesn't seem to belong to Tracker" << std::endl;
0165           continue;
0166         }
0167         // Book histograms
0168         sprintf(name, "Pedestals_%d", detid);
0169         fFile->cd();
0170         fFile->cd("Pedestals");
0171         _TH1F_Pedestals_m[detid] = new TH1F(name, name, 50, 0., 50.);
0172         sprintf(name, "Gains_%d", detid);
0173         fFile->cd();
0174         fFile->cd("Gains");
0175         _TH1F_Gains_m[detid] = new TH1F(name, name, 100, 0., 10.);
0176       }
0177     }
0178   }
0179 
0180   // ------------ method called once each job just after ending the event loop  ------------
0181   void SiPixelFakeGainForHLTReader::endJob() { std::cout << " ---> End job " << std::endl; }
0182 
0183   void SiPixelFakeGainForHLTReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0184     edm::ParameterSetDescription desc;
0185     desc.add<std::string>("fileName", "out.root");
0186     descriptions.addWithDefaultLabel(desc);
0187   }
0188 }  // namespace cms
0189 
0190 using cms::SiPixelFakeGainForHLTReader;
0191 DEFINE_FWK_MODULE(SiPixelFakeGainForHLTReader);