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  void analyze(const edm::Event&, const edm::EventSetup&) override;
56  void endJob() override;
57 
58  private:
60  std::unique_ptr<SiPixelGainCalibrationServiceBase> SiPixelGainCalibrationService_;
61 
62  std::map<uint32_t, TH1F*> _TH1F_Pedestals_m;
63  std::map<uint32_t, TH1F*> _TH1F_Gains_m;
68  };
69 } // namespace cms
70 
71 namespace cms {
73  : tkGeomToken_(esConsumes()) {
74  usesResource(TFileService::kSharedResource);
75  std::string payloadType = conf.getParameter<std::string>("payloadType");
76  if (strcmp(payloadType.c_str(), "HLT") == 0) {
77  SiPixelGainCalibrationService_ = std::make_unique<SiPixelGainCalibrationForHLTService>(conf, consumesCollector());
78  } else if (strcmp(payloadType.c_str(), "Offline") == 0) {
80  std::make_unique<SiPixelGainCalibrationOfflineService>(conf, consumesCollector());
81  } else if (strcmp(payloadType.c_str(), "Full") == 0) {
82  SiPixelGainCalibrationService_ = std::make_unique<SiPixelGainCalibrationService>(conf, consumesCollector());
83  }
84  }
85 
87  //Create Subdirectories
89  TFileDirectory subDirPed = fs->mkdir("Pedestals");
90  TFileDirectory subDirGain = fs->mkdir("Gains");
91  char name[128];
92 
93  unsigned int nmodules = 0;
94  uint32_t nchannels = 0;
95 
96  // Get the calibration data
97  SiPixelGainCalibrationService_->setESObjects(iSetup);
98  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
99  << "[SiPixelCondObjAllPayloadsReader::beginJob] End Reading CondObjects" << std::endl;
100 
101  // Get the Geometry
102  const TrackerGeometry* tkgeom = &iSetup.getData(tkGeomToken_);
103  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
104  << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
105 
106  //Get list of DetIDs
107  std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_->getDetIds();
108 
109  //Create histograms
110  _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
112  fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
113  _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
114  _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
115 
116  // Loop over DetId's
117  int ibin = 1;
118  for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
119  detid_iter++) {
120  uint32_t detid = *detid_iter;
121 
122  sprintf(name, "Pedestals_%d", detid);
123  _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 250, 0., 250.);
124  sprintf(name, "Gains_%d", detid);
125  _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
126 
127  DetId detIdObject(detid);
128  const PixelGeomDetUnit* _PixelGeomDetUnit =
129  dynamic_cast<const PixelGeomDetUnit*>(tkgeom->idToDetUnit(DetId(detid)));
130  if (_PixelGeomDetUnit == nullptr) {
131  edm::LogError("SiPixelCondObjDisplay") << "[SiPixelCondObjAllPayloadsReader::beginJob] the detID " << detid
132  << " doesn't seem to belong to Tracker" << std::endl;
133  continue;
134  }
135 
136  nmodules++;
137 
138  const GeomDetUnit* geoUnit = tkgeom->idToDetUnit(detIdObject);
139  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
140  const PixelTopology& topol = pixDet->specificTopology();
141 
142  // Get the module sizes.
143  int nrows = topol.nrows(); // rows in x
144  int ncols = topol.ncolumns(); // cols in y
145 
146  for (int col_iter = 0; col_iter < ncols; col_iter++) {
147  for (int row_iter = 0; row_iter < nrows; row_iter++) {
148  nchannels++;
149 
150  float gain = SiPixelGainCalibrationService_->getGain(detid, col_iter, row_iter);
151  _TH1F_Gains_m[detid]->Fill(gain);
152  _TH1F_Gains_all->Fill(gain);
153 
154  float ped = SiPixelGainCalibrationService_->getPedestal(detid, col_iter, row_iter);
155  _TH1F_Pedestals_m[detid]->Fill(ped);
156  _TH1F_Pedestals_all->Fill(ped);
157 
158  //edm::LogPrint("SiPixelCondObjAllPayloadsReader") << " Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
159  }
160  }
161 
162  _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
163  _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
164  _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
165  _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
166 
167  ibin++;
168  }
169 
170  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
171  << "[SiPixelCondObjAllPayloadsReader::analyze] ---> PIXEL Modules " << nmodules << std::endl;
172  edm::LogInfo("SiPixelCondObjAllPayloadsReader")
173  << "[SiPixelCondObjAllPayloadsReader::analyze] ---> PIXEL Channels " << nchannels << std::endl;
174  }
175 
176  // ------------ method called once each job just after ending the event loop ------------
178  edm::LogPrint("SiPixelCondObjAllPayloadsReader") << " ---> End job " << std::endl;
179  }
180 } // namespace cms
181 
182 using namespace cms;
static const std::string kSharedResource
Definition: TFileService.h:76
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual int ncolumns() const =0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
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_