SiPixelFakeGainReader

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194
// -*- C++ -*-
//
// Package:    CalibTracker/SiPixelESProducers
// Class:      SiPixelFakeGainReader
//
/**\class SiPixelFakeGainReader SiPixelFakeGainReader.cc SiPixelESProducers/test/SiPixelFakeGainReader.cc

 Description: Test analyzer for fake pixel calibration

 Implementation:
     <Notes on implementation>
*/
//
// Original Author:  Vincenzo CHIOCHIA
//         Created:  Tue Oct 17 17:40:56 CEST 2006
//
//

// system includes
#include <memory>

// user includes
#include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationService.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
#include "Geometry/CommonTopologies/interface/PixelTopology.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

// ROOT includes
#include "TROOT.h"
#include "TFile.h"
#include "TTree.h"
#include "TBranch.h"
#include "TH1F.h"

namespace cms {
  class SiPixelFakeGainReader : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
  public:
    explicit SiPixelFakeGainReader(const edm::ParameterSet& iConfig);
    ~SiPixelFakeGainReader() override = default;

    static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

    virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
    virtual void endRun(const edm::Run&, const edm::EventSetup&) override;
    virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
    virtual void endJob() override;

  private:
    edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
    edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomTokenBeginRun_;
    SiPixelGainCalibrationService SiPixelGainCalibrationService_;

    std::map<uint32_t, TH1F*> _TH1F_Pedestals_m;
    std::map<uint32_t, TH1F*> _TH1F_Gains_m;
    std::string filename_;
    TFile* fFile;
  };
}  // namespace cms

namespace cms {
  SiPixelFakeGainReader::SiPixelFakeGainReader(const edm::ParameterSet& conf)
      : trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
        trackerGeomTokenBeginRun_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
        SiPixelGainCalibrationService_(conf, consumesCollector()),
        filename_(conf.getParameter<std::string>("fileName")) {}

  void SiPixelFakeGainReader::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
    unsigned int nmodules = 0;
    uint32_t nchannels = 0;
    fFile->cd();

    // Get the Geometry
    edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomToken_);
    edm::LogInfo("SiPixelFakeGainReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;

    //  for(TrackerGeometry::DetContainer::const_iterator it = tkgeom->dets().begin(); it != tkgeom->dets().end(); it++){
    //   if( dynamic_cast<PixelGeomDetUnit*>((*it))!=0){
    //     uint32_t detid=((*it)->geographicalId()).rawId();
    // Get the list of DetId's

    std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_.getDetIds();
    // Loop over DetId's
    for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
         detid_iter++) {
      uint32_t detid = *detid_iter;

      DetId detIdObject(detid);
      nmodules++;
      //if(nmodules>3) break;

      std::map<uint32_t, TH1F*>::iterator p_iter = _TH1F_Pedestals_m.find(detid);
      std::map<uint32_t, TH1F*>::iterator g_iter = _TH1F_Gains_m.find(detid);

      const GeomDetUnit* geoUnit = tkgeom->idToDetUnit(detIdObject);
      const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
      const PixelTopology& topol = pixDet->specificTopology();

      // Get the module sizes.
      int nrows = topol.nrows();     // rows in x
      int ncols = topol.ncolumns();  // cols in y

      for (int col_iter = 0; col_iter < ncols; col_iter++) {
        for (int row_iter = 0; row_iter < nrows; row_iter++) {
          nchannels++;

          float ped = SiPixelGainCalibrationService_.getPedestal(detid, col_iter, row_iter);
          p_iter->second->Fill(ped);

          float gain = SiPixelGainCalibrationService_.getGain(detid, col_iter, row_iter);
          g_iter->second->Fill(gain);

          //std::cout << "       Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
        }
      }
    }

    edm::LogInfo("SiPixelFakeGainReader")
        << "[SiPixelFakeGainReader::analyze] ---> PIXEL Modules  " << nmodules << std::endl;
    edm::LogInfo("SiPixelFakeGainReader")
        << "[SiPixelFakeGainReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
    fFile->ls();
    fFile->Write();
    fFile->Close();
  }

  void SiPixelFakeGainReader::endRun(const edm::Run& run, const edm::EventSetup& iSetup) {}

  // ------------ method called once each job just before starting event loop  ------------
  void SiPixelFakeGainReader::beginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
    static int first(1);
    if (1 == first) {
      first = 0;
      edm::LogInfo("SiPixelFakeGainReader") << "[SiPixelFakeGainReader::beginJob] Opening ROOT file  " << std::endl;
      fFile = new TFile(filename_.c_str(), "RECREATE");
      fFile->mkdir("Pedestals");
      fFile->mkdir("Gains");
      fFile->cd();
      char name[128];

      // Get Geometry
      edm::ESHandle<TrackerGeometry> tkgeom = iSetup.getHandle(trackerGeomTokenBeginRun_);

      // Get the calibration data
      SiPixelGainCalibrationService_.setESObjects(iSetup);
      edm::LogInfo("SiPixelFakeGainReader")
          << "[SiPixelFakeGainReader::beginJob] End Reading FakeGainects" << std::endl;
      // Get the list of DetId's
      std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_.getDetIds();
      // Loop over DetId's
      for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
           detid_iter++) {
        uint32_t detid = *detid_iter;

        const PixelGeomDetUnit* _PixelGeomDetUnit =
            dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
        if (_PixelGeomDetUnit == 0) {
          edm::LogError("SiPixelFakeGainDisplay") << "[SiPixelFakeGainReader::beginJob] the detID " << detid
                                                  << " doesn't seem to belong to Tracker" << std::endl;
          continue;
        }
        // Book histograms
        sprintf(name, "Pedestals_%d", detid);
        fFile->cd();
        fFile->cd("Pedestals");
        _TH1F_Pedestals_m[detid] = new TH1F(name, name, 50, 0., 50.);
        sprintf(name, "Gains_%d", detid);
        fFile->cd();
        fFile->cd("Gains");
        _TH1F_Gains_m[detid] = new TH1F(name, name, 100, 0., 10.);
      }
    }
  }

  // ------------ method called once each job just after ending the event loop  ------------
  void SiPixelFakeGainReader::endJob() { std::cout << " ---> End job " << std::endl; }

  void SiPixelFakeGainReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
    edm::ParameterSetDescription desc;
    desc.add<std::string>("fileName", "out.root");
    descriptions.addWithDefaultLabel(desc);
  }
}  // namespace cms

using cms::SiPixelFakeGainReader;
DEFINE_FWK_MODULE(SiPixelFakeGainReader);