CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
cms::SiPixelCondObjForHLTReader Class Reference
Inheritance diagram for cms::SiPixelCondObjForHLTReader:
edm::one::EDAnalyzer< edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void endJob () override
 
 SiPixelCondObjForHLTReader (const edm::ParameterSet &iConfig)
 
 ~SiPixelCondObjForHLTReader () override=default
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::SharedResources >
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Attributes

std::map< uint32_t, double > _deadfrac_m
 
std::map< uint32_t, double > _noisyfrac_m
 
TH1F * _TH1F_Dead_all
 
TH1F * _TH1F_Dead_sum
 
TH1F * _TH1F_Gains_all
 
TH1F * _TH1F_Gains_bpix
 
TH1F * _TH1F_Gains_fpix
 
std::map< uint32_t, TH1F * > _TH1F_Gains_m
 
TH1F * _TH1F_Gains_sum
 
TH1F * _TH1F_Noisy_all
 
TH1F * _TH1F_Noisy_sum
 
TH1F * _TH1F_Pedestals_all
 
TH1F * _TH1F_Pedestals_bpix
 
TH1F * _TH1F_Pedestals_fpix
 
std::map< uint32_t, TH1F * > _TH1F_Pedestals_m
 
TH1F * _TH1F_Pedestals_sum
 
edm::ParameterSet conf_
 
std::unique_ptr< SiPixelGainCalibrationServiceBaseSiPixelGainCalibrationService_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 50 of file SiPixelCondObjForHLTReader.cc.

Constructor & Destructor Documentation

◆ SiPixelCondObjForHLTReader()

SiPixelCondObjForHLTReader::SiPixelCondObjForHLTReader ( const edm::ParameterSet iConfig)
explicit

Definition at line 87 of file SiPixelCondObjForHLTReader.cc.

References conf_, edm::EDConsumerBase::consumesCollector(), edm::ParameterSet::getParameter(), and SiPixelGainCalibrationService_.

88  : conf_(conf), tkGeomToken_(esConsumes()) {
89  if (conf_.getParameter<bool>("useSimRcd"))
91  std::make_unique<SiPixelGainCalibrationForHLTSimService>(conf_, consumesCollector());
92  else
94  std::make_unique<SiPixelGainCalibrationForHLTService>(conf_, consumesCollector());
95  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::unique_ptr< SiPixelGainCalibrationServiceBase > SiPixelGainCalibrationService_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_

◆ ~SiPixelCondObjForHLTReader()

cms::SiPixelCondObjForHLTReader::~SiPixelCondObjForHLTReader ( )
overridedefault

Member Function Documentation

◆ analyze()

void SiPixelCondObjForHLTReader::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 97 of file SiPixelCondObjForHLTReader.cc.

References _deadfrac_m, _noisyfrac_m, _TH1F_Dead_all, _TH1F_Dead_sum, _TH1F_Gains_all, _TH1F_Gains_bpix, _TH1F_Gains_fpix, _TH1F_Gains_m, _TH1F_Gains_sum, _TH1F_Noisy_all, _TH1F_Noisy_sum, _TH1F_Pedestals_all, _TH1F_Pedestals_bpix, _TH1F_Pedestals_fpix, _TH1F_Pedestals_m, _TH1F_Pedestals_sum, conf_, TrackerGeometry::dets(), compareTotals::fs, PedestalClient_cfi::gain, edm::EventSetup::getData(), edm::ParameterSet::getUntrackedParameter(), TrackerGeometry::idToDetUnit(), TFileDirectory::make(), Skims_PA_cff::name, hgcalPlots::ncols, PixelTopology::ncolumns(), PixelTopology::nrows(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, SiPixelGainCalibrationService_, PixelGeomDetUnit::specificTopology(), DetId::subdetId(), and tkGeomToken_.

97  {
98  //Create Subdirectories
100  TFileDirectory subDirPed = fs->mkdir("Pedestals");
101  TFileDirectory subDirGain = fs->mkdir("Gains");
102  char name[128];
103 
104  unsigned int nmodules = 0;
105  uint32_t nchannels = 0;
106  uint32_t ndead = 0;
107  uint32_t nnoisy = 0;
108 
109  // Get the calibration data
110  SiPixelGainCalibrationService_->setESObjects(iSetup);
111  edm::LogInfo("SiPixelCondObjForHLTReader")
112  << "[SiPixelCondObjForHLTReader::beginJob] End Reading CondObjForHLTects" << std::endl;
113 
114  // Get the Geometry
115  const TrackerGeometry *tkgeom = &iSetup.getData(tkGeomToken_);
116  edm::LogInfo("SiPixelCondObjForHLTReader") << " There are " << tkgeom->dets().size() << " detectors" << std::endl;
117 
118  // Get the list of DetId's
119  std::vector<uint32_t> vdetId_ = SiPixelGainCalibrationService_->getDetIds();
120 
121  //Create histograms
122  _TH1F_Dead_sum = fs->make<TH1F>(
123  "Summary_dead", "Dead pixel fraction (0=dead, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
124  _TH1F_Dead_all = fs->make<TH1F>("DeadAll",
125  "Dead pixel fraction (0=dead, 1=alive)",
126  50,
127  0.,
128  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
129  _TH1F_Noisy_sum = fs->make<TH1F>(
130  "Summary_noisy", "Noisy pixel fraction (0=noisy, 1=alive)", vdetId_.size() + 1, 0, vdetId_.size() + 1);
131  _TH1F_Noisy_all = fs->make<TH1F>("NoisyAll",
132  "Noisy pixel fraction (0=noisy, 1=alive)",
133  50,
134  0.,
135  conf_.getUntrackedParameter<double>("maxRangeDeadPixHist", 0.001));
136  _TH1F_Gains_sum = fs->make<TH1F>("Summary_Gain", "Gain Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
138  fs->make<TH1F>("Summary_Pedestal", "Pedestal Summary", vdetId_.size() + 1, 0, vdetId_.size() + 1);
139  _TH1F_Pedestals_all = fs->make<TH1F>("PedestalsAll", "all Pedestals", 350, -100, 250);
140  _TH1F_Pedestals_bpix = fs->make<TH1F>("PedestalsBpix", "bpix Pedestals", 350, -100, 250);
141  _TH1F_Pedestals_fpix = fs->make<TH1F>("PedestalsFpix", "fpix Pedestals", 350, -100, 250);
142  _TH1F_Gains_all = fs->make<TH1F>("GainsAll", "all Gains", 100, 0, 10);
143  _TH1F_Gains_bpix = fs->make<TH1F>("GainsBpix", "bpix Gains", 100, 0, 10);
144  _TH1F_Gains_fpix = fs->make<TH1F>("GainsFpix", "fpix Gains", 100, 0, 10);
145 
146  TTree *tree = new TTree("tree", "tree");
147  uint32_t detid;
148  double gainmeanfortree, gainrmsfortree, pedmeanfortree, pedrmsfortree;
149  tree->Branch("detid", &detid, "detid/I");
150  tree->Branch("ped_mean", &pedmeanfortree, "ped_mean/D");
151  tree->Branch("ped_rms", &pedrmsfortree, "ped_rms/D");
152  tree->Branch("gain_mean", &gainmeanfortree, "gain_mean/D");
153  tree->Branch("gain_rms", &gainrmsfortree, "gain_rms/D");
154 
155  // Loop over DetId's
156  int ibin = 1;
157  for (std::vector<uint32_t>::const_iterator detid_iter = vdetId_.begin(); detid_iter != vdetId_.end();
158  detid_iter++) {
159  detid = *detid_iter;
160 
161  sprintf(name, "Pedestals_%d", detid);
162  _TH1F_Pedestals_m[detid] = subDirPed.make<TH1F>(name, name, 350, -100., 250.);
163  sprintf(name, "Gains_%d", detid);
164  _TH1F_Gains_m[detid] = subDirGain.make<TH1F>(name, name, 100, 0., 10.);
165 
166  DetId detIdObject(detid);
167  const PixelGeomDetUnit *_PixelGeomDetUnit =
168  dynamic_cast<const PixelGeomDetUnit *>(tkgeom->idToDetUnit(DetId(detid)));
169  if (_PixelGeomDetUnit == nullptr) {
170  edm::LogError("SiPixelCondObjHLTDisplay") << "[SiPixelCondObjHLTReader::beginJob] the detID " << detid
171  << " doesn't seem to belong to Tracker" << std::endl;
172  continue;
173  }
174 
175  _deadfrac_m[detid] = 0.;
176  _noisyfrac_m[detid] = 0.;
177 
178  nmodules++;
179 
180  const GeomDetUnit *geoUnit = tkgeom->idToDetUnit(detIdObject);
181  const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
182  const PixelTopology &topol = pixDet->specificTopology();
183 
184  // Get the module sizes.
185  int nrows = topol.nrows(); // rows in x
186  int ncols = topol.ncolumns(); // cols in y
187  float nchannelspermod = 0;
188 
189  for (int col_iter = 0; col_iter < ncols; col_iter++) {
190  for (int row_iter = 0; row_iter < nrows; row_iter++) {
191  nchannelspermod++;
192  nchannels++;
193 
194  if (SiPixelGainCalibrationService_->isDead(detid, col_iter, row_iter)) {
195  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found dead pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
196  ndead++;
197  _deadfrac_m[detid]++;
198  continue;
199  } else if (SiPixelGainCalibrationService_->isNoisy(detid, col_iter, row_iter)) {
200  // edm::LogPrint("SiPixelCondObjForHLTReader") << "found noisy pixel " << detid << " " <<col_iter << "," << row_iter << std::endl;
201  nnoisy++;
202  _noisyfrac_m[detid]++;
203  continue;
204  }
205 
206  float gain = SiPixelGainCalibrationService_->getGain(detid, col_iter, row_iter);
207  _TH1F_Gains_m[detid]->Fill(gain);
208  _TH1F_Gains_all->Fill(gain);
209 
210  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
211  _TH1F_Gains_bpix->Fill(gain);
212  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
213  _TH1F_Gains_fpix->Fill(gain);
214 
215  float ped = SiPixelGainCalibrationService_->getPedestal(detid, col_iter, row_iter);
216  _TH1F_Pedestals_m[detid]->Fill(ped);
217  _TH1F_Pedestals_all->Fill(ped);
218  // edm::LogPrint("SiPixelCondObjForHLTReader")<<"detid "<<detid<<" ped "<<ped<<std::endl;
219 
220  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
221  _TH1F_Pedestals_bpix->Fill(ped);
222  if (detIdObject.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
223  _TH1F_Pedestals_fpix->Fill(ped);
224 
225  // edm::LogPrint("SiPixelCondObjForHLTReader") <<" DetId "<<detid<<" Col "<<col_iter<<" Row "<<row_iter<<" Ped "<<ped<<" Gain "<<gain<<std::endl;
226  }
227  }
228 
229  _deadfrac_m[detid] /= nchannelspermod;
230  _noisyfrac_m[detid] /= nchannelspermod;
231  _TH1F_Dead_sum->SetBinContent(ibin, _deadfrac_m[detid]);
232  _TH1F_Dead_all->Fill(_deadfrac_m[detid]);
233  _TH1F_Noisy_sum->SetBinContent(ibin, _noisyfrac_m[detid]);
234  _TH1F_Noisy_all->Fill(_noisyfrac_m[detid]);
235  _TH1F_Gains_sum->SetBinContent(ibin, _TH1F_Gains_m[detid]->GetMean());
236  _TH1F_Gains_sum->SetBinError(ibin, _TH1F_Gains_m[detid]->GetRMS());
237  _TH1F_Pedestals_sum->SetBinContent(ibin, _TH1F_Pedestals_m[detid]->GetMean());
238  _TH1F_Pedestals_sum->SetBinError(ibin, _TH1F_Pedestals_m[detid]->GetRMS());
239 
240  gainmeanfortree = _TH1F_Gains_m[detid]->GetMean();
241  gainrmsfortree = _TH1F_Gains_m[detid]->GetRMS();
242  pedmeanfortree = _TH1F_Pedestals_m[detid]->GetMean();
243  pedrmsfortree = _TH1F_Pedestals_m[detid]->GetRMS();
244  //edm::LogPrint("SiPixelCondObjForHLTReader")<<"DetId "<<detid<<" GainMean "<<gainmeanfortree<<" RMS "<<gainrmsfortree<<" PedMean "<<pedmeanfortree<<" RMS "<<pedrmsfortree<<std::endl;
245  tree->Fill();
246 
247  ibin++;
248  }
249 
250  edm::LogInfo("SiPixelCondObjForHLTReader")
251  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Modules " << nmodules << std::endl;
252  edm::LogInfo("SiPixelCondObjForHLTReader")
253  << "[SiPixelCondObjForHLTReader::analyze] ---> PIXEL Channels (i.e. Number of Columns)" << nchannels
254  << std::endl;
255 
256  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> SUMMARY :" << std::endl;
257  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << ndead << " dead pixels" << std::endl;
258  edm::LogPrint("SiPixelCondObjForHLTReader") << "Encounted " << nnoisy << " noisy pixels" << std::endl;
259  edm::LogPrint("SiPixelCondObjForHLTReader")
260  << "The Gain Mean is " << _TH1F_Gains_all->GetMean() << " with rms " << _TH1F_Gains_all->GetRMS() << std::endl;
261  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Gains_bpix->GetMean() << " with rms "
262  << _TH1F_Gains_bpix->GetRMS() << std::endl;
263  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Gains_fpix->GetMean() << " with rms "
264  << _TH1F_Gains_fpix->GetRMS() << std::endl;
265  edm::LogPrint("SiPixelCondObjForHLTReader") << "The Ped Mean is " << _TH1F_Pedestals_all->GetMean() << " with rms "
266  << _TH1F_Pedestals_all->GetRMS() << std::endl;
267  edm::LogPrint("SiPixelCondObjForHLTReader") << " in BPIX " << _TH1F_Pedestals_bpix->GetMean()
268  << " with rms " << _TH1F_Pedestals_bpix->GetRMS() << std::endl;
269  edm::LogPrint("SiPixelCondObjForHLTReader") << " in FPIX " << _TH1F_Pedestals_fpix->GetMean()
270  << " with rms " << _TH1F_Pedestals_fpix->GetRMS() << std::endl;
271  }
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
std::map< uint32_t, double > _deadfrac_m
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 getUntrackedParameter(std::string const &, T const &) const
T * make(const Args &...args) const
make new ROOT object
std::map< uint32_t, TH1F * > _TH1F_Gains_m
Log< level::Warning, true > LogPrint
Log< level::Info, false > LogInfo
Definition: DetId.h:17
std::unique_ptr< SiPixelGainCalibrationServiceBase > SiPixelGainCalibrationService_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
Definition: tree.py:1
std::map< uint32_t, TH1F * > _TH1F_Pedestals_m
std::map< uint32_t, double > _noisyfrac_m

◆ endJob()

void SiPixelCondObjForHLTReader::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 274 of file SiPixelCondObjForHLTReader.cc.

274  {
275  edm::LogPrint("SiPixelCondObjForHLTReader") << " ---> End job " << std::endl;
276  }
Log< level::Warning, true > LogPrint

◆ fillDescriptions()

void SiPixelCondObjForHLTReader::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 278 of file SiPixelCondObjForHLTReader.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), and submitPVResolutionJobs::desc.

278  {
280  desc.setComment("EDAnalyzer to read per-module SiPixelGainCalibrationForHLT payloads in the EventSetup");
281  desc.add<bool>("useSimRcd", false);
282  desc.addUntracked<double>("maxRangeDeadPixHist", 0.001);
283  descriptions.addWithDefaultLabel(desc);
284  }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ _deadfrac_m

std::map<uint32_t, double> cms::SiPixelCondObjForHLTReader::_deadfrac_m
private

Definition at line 68 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _noisyfrac_m

std::map<uint32_t, double> cms::SiPixelCondObjForHLTReader::_noisyfrac_m
private

Definition at line 69 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Dead_all

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Dead_all
private

Definition at line 75 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Dead_sum

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Dead_sum
private

Definition at line 71 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Gains_all

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Gains_all
private

Definition at line 77 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Gains_bpix

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Gains_bpix
private

Definition at line 79 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Gains_fpix

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Gains_fpix
private

Definition at line 80 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Gains_m

std::map<uint32_t, TH1F *> cms::SiPixelCondObjForHLTReader::_TH1F_Gains_m
private

Definition at line 67 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Gains_sum

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Gains_sum
private

Definition at line 73 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Noisy_all

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Noisy_all
private

Definition at line 76 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Noisy_sum

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Noisy_sum
private

Definition at line 72 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Pedestals_all

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Pedestals_all
private

Definition at line 78 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Pedestals_bpix

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Pedestals_bpix
private

Definition at line 81 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Pedestals_fpix

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Pedestals_fpix
private

Definition at line 82 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Pedestals_m

std::map<uint32_t, TH1F *> cms::SiPixelCondObjForHLTReader::_TH1F_Pedestals_m
private

Definition at line 66 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ _TH1F_Pedestals_sum

TH1F* cms::SiPixelCondObjForHLTReader::_TH1F_Pedestals_sum
private

Definition at line 74 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().

◆ conf_

edm::ParameterSet cms::SiPixelCondObjForHLTReader::conf_
private

Definition at line 61 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze(), and SiPixelCondObjForHLTReader().

◆ SiPixelGainCalibrationService_

std::unique_ptr<SiPixelGainCalibrationServiceBase> cms::SiPixelCondObjForHLTReader::SiPixelGainCalibrationService_
private

Definition at line 64 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze(), and SiPixelCondObjForHLTReader().

◆ tkGeomToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> cms::SiPixelCondObjForHLTReader::tkGeomToken_
private

Definition at line 63 of file SiPixelCondObjForHLTReader.cc.

Referenced by analyze().