CMS 3D CMS Logo

SiPixelCondObjAllPayloadsReader.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelCondObjAllPayloadsReader
4 // Class: SiPixelCondObjAllPayloadsReader
5 //
13 //
14 // Original Author: Vincenzo CHIOCHIA
15 // Created: Tue Oct 17 17:40:56 CEST 2006
16 // $Id: SiPixelCondObjAllPayloadsReader.h,v 1.4 2009/05/28 22:12:54 dlange Exp $
17 //
18 //
19 
20 // system includes
21 #include <memory>
22 
23 // user includes
41 
42 // ROOT includes
43 #include "TROOT.h"
44 #include "TFile.h"
45 #include "TTree.h"
46 #include "TBranch.h"
47 #include "TH1F.h"
48 #include <string>
49 
50 namespace cms {
51  class SiPixelCondObjAllPayloadsReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
52  public:
53  explicit SiPixelCondObjAllPayloadsReader(const edm::ParameterSet& iConfig);
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56  void analyze(const edm::Event&, const edm::EventSetup&) override;
57  void endJob() override;
58 
59  private:
61  std::unique_ptr<SiPixelGainCalibrationServiceBase> SiPixelGainCalibrationService_;
62 
63  std::map<uint32_t, TH1F*> _TH1F_Pedestals_m;
64  std::map<uint32_t, TH1F*> _TH1F_Gains_m;
69  };
70 } // namespace cms
71 
72 namespace cms {
73 
76  desc.setComment("EDAnalyzer to read per-module SiPixelGainCalibration payloads in the EventSetup, for any type");
77  desc.add<std::string>("payloadType", "HLT");
78  descriptions.addWithDefaultLabel(desc);
79  }
80 
82  : tkGeomToken_(esConsumes()) {
83  usesResource(TFileService::kSharedResource);
84  std::string payloadType = conf.getParameter<std::string>("payloadType");
85  if (strcmp(payloadType.c_str(), "HLT") == 0) {
86  SiPixelGainCalibrationService_ = std::make_unique<SiPixelGainCalibrationForHLTService>(conf, consumesCollector());
87  } else if (strcmp(payloadType.c_str(), "Offline") == 0) {
89  std::make_unique<SiPixelGainCalibrationOfflineService>(conf, consumesCollector());
90  } else if (strcmp(payloadType.c_str(), "Full") == 0) {
91  SiPixelGainCalibrationService_ = std::make_unique<SiPixelGainCalibrationService>(conf, consumesCollector());
92  }
93  }
94 
96  //Create Subdirectories
98  TFileDirectory subDirPed = fs->mkdir("Pedestals");
99  TFileDirectory subDirGain = fs->mkdir("Gains");
100  char name[128];
101 
102  unsigned int nmodules = 0;
103  uint32_t nchannels = 0;
104 
105  // Get the calibration data
106  SiPixelGainCalibrationService_->setESObjects(iSetup);
107  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
108  << "[SiPixelCondObjAllPayloadsReader::beginJob] End Reading CondObjects" << std::endl;
109 
110  // Get the Geometry
111  const TrackerGeometry* tkgeom = &iSetup.getData(tkGeomToken_);
112  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
113  << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
114 
115  //Get list of DetIDs
116  std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_->getDetIds();
117 
118  //Create histograms
119  _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
121  fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
122  _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
123  _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
124 
125  // Loop over DetId's
126  int ibin = 1;
127  for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
128  detid_iter++) {
129  uint32_t detid = *detid_iter;
130 
131  sprintf(name, "Pedestals_%d", detid);
132  _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 250, 0., 250.);
133  sprintf(name, "Gains_%d", detid);
134  _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
135 
136  DetId detIdObject(detid);
137  const PixelGeomDetUnit* _PixelGeomDetUnit =
138  dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
139  if (_PixelGeomDetUnit == nullptr) {
140  edm::LogError("SiPixelCondObjDisplay") << "[SiPixelCondObjAllPayloadsReader::beginJob] the detID " << detid
141  << " doesn't seem to belong to Tracker" << std::endl;
142  continue;
143  }
144 
145  nmodules++;
146 
147  const GeomDetUnit* geoUnit = tkgeom->idToDetUnit(detIdObject);
148  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
149  const PixelTopology& topol = pixDet->specificTopology();
150 
151  // Get the module sizes.
152  int nrows = topol.nrows(); // rows in x
153  int ncols = topol.ncolumns(); // cols in y
154 
155  for (int col_iter = 0; col_iter < ncols; col_iter++) {
156  for (int row_iter = 0; row_iter < nrows; row_iter++) {
157  nchannels++;
158 
159  float gain = SiPixelGainCalibrationService_->getGain(detid, col_iter, row_iter);
160  _TH1F_Gains_m[detid]->Fill(gain);
161  _TH1F_Gains_all->Fill(gain);
162 
163  float ped = SiPixelGainCalibrationService_->getPedestal(detid, col_iter, row_iter);
164  _TH1F_Pedestals_m[detid]->Fill(ped);
165  _TH1F_Pedestals_all->Fill(ped);
166 
167  //edm::LogPrint("SiPixelCondObjAllPayloadsReader") << " Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
168  }
169  }
170 
171  _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
172  _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
173  _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
174  _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
175 
176  ibin++;
177  }
178 
179  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
180  << "[SiPixelCondObjAllPayloadsReader::analyze] ---> PIXEL Modules " << nmodules << std::endl;
181  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
182  << "[SiPixelCondObjAllPayloadsReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
183  }
184 
185  // ------------ method called once each job just after ending the event loop ------------
187  edm::LogPrint("SiPixelCondObjAllPayloadsReader") << " ---> End job " << std::endl;
188  }
189 } // namespace cms
190 
191 using namespace cms;
static const std::string kSharedResource
Definition: TFileService.h:76
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual int ncolumns() const =0
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
SiPixelCondObjAllPayloadsReader(const edm::ParameterSet &iConfig)
virtual int nrows() const =0
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Log< level::Error, false > LogError
T * make(const Args &...args) const
make new ROOT object
int iEvent
Definition: GenABIO.cc:224
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
void analyze(const edm::Event &, const edm::EventSetup &) override
std::unique_ptr< SiPixelGainCalibrationServiceBase > SiPixelGainCalibrationService_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Log< level::Warning, true > LogPrint
Namespace of DDCMS conversion namespace.
Log< level::Info, false > LogInfo
Definition: DetId.h:17
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_