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/SiPixelFakeGainReader.h"
0004 
0005 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0006 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0007 namespace cms {
0008   SiPixelFakeGainReader::SiPixelFakeGainReader(const edm::ParameterSet& conf)
0009       : conf_(conf),
0010         trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0011         trackerGeomTokenBeginRun_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0012         SiPixelGainCalibrationService_(conf, consumesCollector()),
0013         filename_(conf.getParameter<std::string>("fileName")) {}
0014 
0015   void SiPixelFakeGainReader::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("SiPixelFakeGainReader") << " 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_ = SiPixelGainCalibrationService_.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 
0054           float ped = SiPixelGainCalibrationService_.getPedestal(detid, col_iter, row_iter);
0055           p_iter->second->Fill(ped);
0056 
0057           float gain = SiPixelGainCalibrationService_.getGain(detid, col_iter, row_iter);
0058           g_iter->second->Fill(gain);
0059 
0060           //std::cout << "       Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
0061         }
0062       }
0063     }
0064 
0065     edm::LogInfo("SiPixelFakeGainReader")
0066         << "[SiPixelFakeGainReader::analyze] ---> PIXEL Modules  " << nmodules << std::endl;
0067     edm::LogInfo("SiPixelFakeGainReader")
0068         << "[SiPixelFakeGainReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
0069     fFile->ls();
0070     fFile->Write();
0071     fFile->Close();
0072   }
0073 
0074   void SiPixelFakeGainReader::endRun(const edm::Run& run, const edm::EventSetup& iSetup) {}
0075 
0076   // ------------ method called once each job just before starting event loop  ------------
0077   void SiPixelFakeGainReader::beginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0078     static int first(1);
0079     if (1 == first) {
0080       first = 0;
0081       edm::LogInfo("SiPixelFakeGainReader") << "[SiPixelFakeGainReader::beginJob] Opening ROOT file  " << std::endl;
0082       fFile = new TFile(filename_.c_str(), "RECREATE");
0083       fFile->mkdir("Pedestals");
0084       fFile->mkdir("Gains");
0085       fFile->cd();
0086       char name[128];
0087 
0088       // Get Geometry
0089       edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomTokenBeginRun_);
0090 
0091       // Get the calibration data
0092       SiPixelGainCalibrationService_.setESObjects(iSetup);
0093       edm::LogInfo("SiPixelFakeGainReader")
0094           << "[SiPixelFakeGainReader::beginJob] End Reading FakeGainects" << std::endl;
0095       // Get the list of DetId's
0096       std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_.getDetIds();
0097       // Loop over DetId's
0098       for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
0099            detid_iter++) {
0100         uint32_t detid = *detid_iter;
0101 
0102         const PixelGeomDetUnit* _PixelGeomDetUnit =
0103             dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
0104         if (_PixelGeomDetUnit == 0) {
0105           edm::LogError("SiPixelFakeGainDisplay") << "[SiPixelFakeGainReader::beginJob] the detID " << detid
0106                                                   << " doesn't seem to belong to Tracker" << std::endl;
0107           continue;
0108         }
0109         // Book histograms
0110         sprintf(name, "Pedestals_%d", detid);
0111         fFile->cd();
0112         fFile->cd("Pedestals");
0113         _TH1F_Pedestals_m[detid] = new TH1F(name, name, 50, 0., 50.);
0114         sprintf(name, "Gains_%d", detid);
0115         fFile->cd();
0116         fFile->cd("Gains");
0117         _TH1F_Gains_m[detid] = new TH1F(name, name, 100, 0., 10.);
0118       }
0119     }
0120   }
0121 
0122   // ------------ method called once each job just after ending the event loop  ------------
0123   void SiPixelFakeGainReader::endJob() { std::cout << " ---> End job " << std::endl; }
0124 }  // namespace cms