Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:34

0001 #include <memory>
0002 
0003 #include "CalibTracker/SiPixelESProducers/test/SiPixelFakeGainForHLTReader.h"
0004 
0005 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0006 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0007 namespace cms {
0008   SiPixelFakeGainForHLTReader::SiPixelFakeGainForHLTReader(const edm::ParameterSet& conf)
0009       : conf_(conf),
0010         trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0011         trackerGeomTokenBeginRun_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0012         SiPixelGainCalibrationForHLTService_(conf, consumesCollector()),
0013         filename_(conf.getParameter<std::string>("fileName")) {}
0014 
0015   void SiPixelFakeGainForHLTReader::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0016     unsigned int nmodules = 0;
0017     uint32_t nchannels = 0;
0018     fFile->cd();
0019 
0020     // Get the Geometry
0021     edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomToken_);
0022     edm::LogInfo("SiPixelFakeGainForHLTReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
0023 
0024     //  for(TrackerGeometry::DetContainer::const_iterator it = tkgeom->dets().begin(); it != tkgeom->dets().end(); it++){
0025     //   if( dynamic_cast<PixelGeomDetUnit*>((*it))!=0){
0026     //     uint32_t detid=((*it)->geographicalId()).rawId();
0027     // Get the list of DetId's
0028 
0029     std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationForHLTService_.getDetIds();
0030     // Loop over DetId's
0031     for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0032          detid_iter++) {
0033       uint32_t detid = *detid_iter;
0034 
0035       DetId detIdObject(detid);
0036       nmodules++;
0037       //if(nmodules>3) break;
0038 
0039       std::map<uint32_t, TH1F*>::iterator p_iter = _TH1F_Pedestals_m.find(detid);
0040       std::map<uint32_t, TH1F*>::iterator g_iter = _TH1F_Gains_m.find(detid);
0041 
0042       const GeomDetUnit* geoUnit = tkgeom->idToDetUnit(detIdObject);
0043       const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
0044       const PixelTopology& topol = pixDet->specificTopology();
0045 
0046       // Get the module sizes.
0047       int nrows = topol.nrows();     // rows in x
0048       int ncols = topol.ncolumns();  // cols in y
0049 
0050       for (int col_iter = 0; col_iter < ncols; col_iter++) {
0051         for (int row_iter = 0; row_iter < nrows; row_iter++) {
0052           nchannels++;
0053           float ped = SiPixelGainCalibrationForHLTService_.getPedestal(detid, col_iter, row_iter);
0054           float gain = SiPixelGainCalibrationForHLTService_.getGain(detid, col_iter, row_iter);
0055           p_iter->second->Fill(ped);
0056           g_iter->second->Fill(gain);
0057         }
0058       }
0059     }
0060 
0061     edm::LogInfo("SiPixelFakeGainForHLTReader")
0062         << "[SiPixelFakeGainForHLTReader::analyze] ---> PIXEL Modules  " << nmodules << std::endl;
0063     edm::LogInfo("SiPixelFakeGainForHLTReader")
0064         << "[SiPixelFakeGainForHLTReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
0065     fFile->ls();
0066     fFile->Write();
0067     fFile->Close();
0068   }
0069 
0070   void SiPixelFakeGainForHLTReader::endRun(const edm::Run& run, const edm::EventSetup& iSetup) {}
0071 
0072   // ----------------------------------------------------------------------
0073   void SiPixelFakeGainForHLTReader::beginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0074     static int first(1);
0075     if (1 == first) {
0076       first = 0;
0077       edm::LogInfo("SiPixelFakeGainForHLTReader")
0078           << "[SiPixelFakeGainForHLTReader::beginJob] Opening ROOT file  " << std::endl;
0079       fFile = new TFile(filename_.c_str(), "RECREATE");
0080       fFile->mkdir("Pedestals");
0081       fFile->mkdir("Gains");
0082       fFile->cd();
0083       char name[128];
0084 
0085       // Get Geometry
0086       edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomTokenBeginRun_);
0087 
0088       // Get the calibrationForHLT data
0089       SiPixelGainCalibrationForHLTService_.setESObjects(iSetup);
0090       edm::LogInfo("SiPixelFakeGainForHLTReader")
0091           << "[SiPixelFakeGainForHLTReader::beginJob] End Reading FakeGainForHLTects" << std::endl;
0092       // Get the list of DetId's
0093       std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationForHLTService_.getDetIds();
0094       // Loop over DetId's
0095       for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0096            detid_iter++) {
0097         uint32_t detid = *detid_iter;
0098 
0099         const PixelGeomDetUnit* _PixelGeomDetUnit =
0100             dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
0101         if (_PixelGeomDetUnit == 0) {
0102           edm::LogError("SiPixelFakeGainForHLTDisplay") << "[SiPixelFakeGainForHLTReader::beginJob] the detID " << detid
0103                                                         << " doesn't seem to belong to Tracker" << std::endl;
0104           continue;
0105         }
0106         // Book histograms
0107         sprintf(name, "Pedestals_%d", detid);
0108         fFile->cd();
0109         fFile->cd("Pedestals");
0110         _TH1F_Pedestals_m[detid] = new TH1F(name, name, 50, 0., 50.);
0111         sprintf(name, "Gains_%d", detid);
0112         fFile->cd();
0113         fFile->cd("Gains");
0114         _TH1F_Gains_m[detid] = new TH1F(name, name, 100, 0., 10.);
0115       }
0116     }
0117   }
0118 
0119   // ------------ method called once each job just after ending the event loop  ------------
0120   void SiPixelFakeGainForHLTReader::endJob() { std::cout << " ---> End job " << std::endl; }
0121 }  // namespace cms