CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HcalLutAnalyzer Class Reference

#include <Test/HcalLutAnalyzer/plugins/HcalLutAnalyzer.cc>

Inheritance diagram for HcalLutAnalyzer:
edm::one::EDAnalyzer< edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 HcalLutAnalyzer (const edm::ParameterSet &)
 
 ~HcalLutAnalyzer () override
 
- 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 Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::vector< std::string > effpedestals_
 
std::vector< std::string > gains_
 
std::string inputDir
 
std::vector< std::string > pedestals_
 
std::string plotsDir
 
double Pmax
 
double Pmin
 
std::vector< std::string > quality_
 
std::vector< std::string > respcorrs_
 
std::vector< std::string > tags_
 
edm::ESGetToken< HcalTopology, HcalRecNumberingRecordtok_htopo_
 
double Ymax
 
double Ymin
 
double Zmax
 
double Zmin
 

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

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 49 of file HcalLutAnalyzer.cc.

Constructor & Destructor Documentation

◆ HcalLutAnalyzer()

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

Definition at line 77 of file HcalLutAnalyzer.cc.

References effpedestals_, gains_, edm::ParameterSet::getParameter(), inputDir, pedestals_, plotsDir, Pmax, Pmin, quality_, respcorrs_, AlCaHLTBitMon_QueryRunRegistry::string, tags_, tok_htopo_, Ymax, Ymin, Zmax, and Zmin.

77  {
78  inputDir = iConfig.getParameter<std::string>("inputDir");
79  plotsDir = iConfig.getParameter<std::string>("plotsDir");
80  tags_ = iConfig.getParameter<std::vector<std::string> >("tags");
81  quality_ = iConfig.getParameter<std::vector<std::string> >("quality");
82  pedestals_ = iConfig.getParameter<std::vector<std::string> >("pedestals");
83  effpedestals_ = iConfig.getParameter<std::vector<std::string> >("effpedestals");
84  gains_ = iConfig.getParameter<std::vector<std::string> >("gains");
85  respcorrs_ = iConfig.getParameter<std::vector<std::string> >("respcorrs");
86 
87  Zmin = iConfig.getParameter<double>("Zmin");
88  Zmax = iConfig.getParameter<double>("Zmax");
89  Ymin = iConfig.getParameter<double>("Ymin");
90  Ymax = iConfig.getParameter<double>("Ymax");
91  Pmin = iConfig.getParameter<double>("Pmin");
92  Pmax = iConfig.getParameter<double>("Pmax");
93 
94  tok_htopo_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
95 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::string plotsDir
std::vector< std::string > effpedestals_
std::vector< std::string > quality_
std::vector< std::string > gains_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > tok_htopo_
std::string inputDir
std::vector< std::string > pedestals_
std::vector< std::string > respcorrs_
std::vector< std::string > tags_

◆ ~HcalLutAnalyzer()

HcalLutAnalyzer::~HcalLutAnalyzer ( )
inlineoverride

Definition at line 52 of file HcalLutAnalyzer.cc.

52 {};

Member Function Documentation

◆ analyze()

void HcalLutAnalyzer::analyze ( const edm::Event ,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 97 of file HcalLutAnalyzer.cc.

References funct::abs(), cms::cuda::assert(), newFWLiteAna::base, edmScanValgrind::buffer, LutXml::create_lut_map(), ztail::d, histoStyle::doRatio, effpedestals_, HcalObjRepresent::Fill(), edm::FileInPath::fullPath(), gains_, HcalGenericDetId::genericSubdet(), edm::EventSetup::getData(), h, photonIsolationHIProducer_cfi::hbhe, HcalBarrel, HcalEndcap, HcalForward, HcalGenericDetId::HcalGenTriggerTower, HcalOther, HcalOuter, mps_fire::i, l1ctLayer2EG_cff::id, LEDCalibrationChannels::ieta, cuy::ii, timingPdfMaker::infile, inputDir, createfilelist::int, LEDCalibrationChannels::iphi, dqmiolumiharvest::j, AlCaHLTBitMon_ParallelJobs::p, pedestals_, plotsDir, Pmax, Pmin, quality_, alignCSCRings::r, DetId::rawId(), respcorrs_, findQualityFiles::size, tags_, remoteMonitoring_LASER_era2018_cfg::threshold, tok_htopo_, HcalTopology::valid(), HcalTopology::validHT(), Ymax, Ymin, Zmax, and Zmin.

97  {
98  using namespace std;
99 
100  const HcalTopology* topology = &iSetup.getData(tok_htopo_);
101 
102  typedef std::vector<std::string> vstring;
103  typedef std::map<unsigned long int, float> LUTINPUT;
104 
105  static const int NVAR = 6; //variables
106  static const int NDET = 2; //detectors
107  static const int NDEP = 7; //depths
108  static const int NLEV = 3; //old,new,ratio
109 
110  const bool doRatio[NVAR] = {false, true, true, false, true, false};
111  const char* titleVar[NVAR] = {"Pedestals", "RespCorrs", "Gains", "Threshold", "LUTs", "EffPedestals"};
112  const char* titleHisR[NLEV] = {"Old", "New", "Ratio"};
113  const char* titleHisD[NLEV] = {"Old", "New", "Difference"};
114  const char* titleDet[4] = {"HBHE", "HF", "HEP17", "HO"};
115  const int DEP[NDET] = {7, 4};
116  const char* titleDep[NDEP] = {"depth1", "depth2", "depth3", "depth4", "depth5", "depth6", "depth7"};
117 
118  TH2D* r[NVAR][NDET];
119  TProfile* p[NVAR][NDET];
120 
121  TH2D* h[NVAR][NLEV][NDET][NDEP];
122  TH2D* hlut[4][NLEV];
123  TH2D* hslope[2];
124  TH2D* houtput[4][2];
125 
126  for (int d = 0; d < 4; ++d) {
127  for (int i = 0; i < 3; ++i) {
128  hlut[d][i] = new TH2D(Form("Lut_%s_%s", titleDet[d], titleHisR[i]),
129  Form("Input LUT, %s", (i == 2 ? "Ratio" : tags_[i].c_str())),
130  260,
131  0,
132  260,
133  240,
134  0,
135  i == NLEV - 1 ? 3 : 2400);
136  hlut[d][i]->SetMarkerColor(d == 0 ? kBlue : d == 1 ? kGreen + 2 : d == 2 ? kRed : kCyan);
137  hlut[d][i]->SetXTitle("raw adc");
138  hlut[d][i]->SetYTitle("lin adc");
139  }
140  }
141 
142  for (int d = 0; d < NDET; ++d) {
143  hslope[d] = new TH2D(
144  Form("GainLutScatter_%s", titleDet[d]), Form("Gain-Lutslope scatter, %s", titleDet[d]), 200, 0, 2, 200, 0, 2);
145  hslope[d]->SetXTitle("Gain x RespCorr ratio");
146  hslope[d]->SetYTitle("Lut ratio");
147 
148  for (int j = 0; j < NVAR; ++j) {
149  double rmin = doRatio[j] ? Ymin : -4;
150  double rmax = doRatio[j] ? Ymax : 4;
151  r[j][d] = new TH2D(Form("r%s_%s", titleVar[j], titleDet[d]),
152  Form("%s, %s", titleVar[j], titleDet[d]),
153  83,
154  -41.5,
155  41.5,
156  250,
157  rmin,
158  rmax);
159  r[j][d]->SetXTitle("iEta");
160  r[j][d]->SetYTitle(doRatio[j] ? "New / Old" : "New - Old");
161  p[j][d] = new TProfile(
162  Form("p%s_%s", titleVar[j], titleDet[d]), Form("%s, %s", titleVar[j], titleDet[d]), 83, -41.5, 41.5);
163  p[j][d]->SetXTitle("iEta");
164  p[j][d]->SetYTitle(doRatio[j] ? "New / Old" : "New - Old");
165  p[j][d]->SetMarkerStyle(20);
166  p[j][d]->SetMarkerSize(0.9);
167  p[j][d]->SetMarkerColor(kBlue);
168 
169  for (int p = 0; p < DEP[d]; ++p) {
170  for (int i = 0; i < NLEV; ++i) {
171  const char* titHist = doRatio[j] ? titleHisR[i] : titleHisD[i];
172  h[j][i][d][p] = new TH2D(Form("h%s_%s_%s_%s", titleVar[j], titHist, titleDet[d], titleDep[p]),
173  Form("%s, %s, %s, %s", titleVar[j], titHist, titleDet[d], titleDep[p]),
174  83,
175  -41.5,
176  41.5,
177  72,
178  0.5,
179  72.5);
180  h[j][i][d][p]->SetXTitle("iEta");
181  h[j][i][d][p]->SetYTitle("iPhi");
182  }
183  }
184  }
185  }
186 
187  for (int i = 0; i < 4; ++i) {
188  int color = i == 0 ? kBlue : i == 1 ? kViolet : i == 2 ? kGreen + 2 : kRed;
189  houtput[i][0] =
190  new TH2D(Form("houtlut0_%d", i), Form("Output LUT, %s", tags_[0].c_str()), 2100, 0, 2100, 260, 0, 260);
191  houtput[i][1] =
192  new TH2D(Form("houtlut1_%d", i), Form("Output LUT, %s", tags_[1].c_str()), 2100, 0, 2100, 260, 0, 260);
193  for (int j = 0; j < 2; ++j) {
194  houtput[i][j]->SetMarkerColor(color);
195  houtput[i][j]->SetLineColor(color);
196  }
197  }
198 
199  //FILL LUT INPUT DATA
200  LUTINPUT lutgain[2];
201  LUTINPUT lutresp[2];
202  LUTINPUT lutpede[2];
203  LUTINPUT luteffpede[2];
204 
205  assert(tags_.size() == 2);
206  assert(gains_.size() == 2);
207  assert(respcorrs_.size() == 2);
208  assert(pedestals_.size() == 2);
209  assert(effpedestals_.size() == 2);
210 
211  unsigned long int iraw;
212  int ieta, iphi, idep;
213  string det, base;
214  float val1, val2, val3, val4;
215  float wid1, wid2, wid3, wid4;
216  char buffer[1024];
217 
218  std::vector<HcalDetId> BadChans[2];
219  std::vector<HcalDetId> ZeroLuts[2];
220 
221  //Read input condition files
222  for (int ii = 0; ii < 2; ++ii) {
223  //Gains
224  std::ifstream infile(
225  edm::FileInPath(Form("%s/Gains/Gains_Run%s.txt", inputDir.c_str(), gains_[ii].c_str())).fullPath().c_str());
226  assert(!infile.fail());
227  while (!infile.eof()) {
228  infile.getline(buffer, 1024);
229  if (buffer[0] == '#')
230  continue;
231  std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> val2 >> val3 >> val4 >> iraw;
232  if (det != "HB" && det != "HE" && det != "HF")
233  continue;
234 
235  float theval = (val1 + val2 + val3 + val4) / 4.0;
236 
237  HcalSubdetector subdet = det == "HB" ? HcalBarrel
238  : det == "HE" ? HcalEndcap
239  : det == "HF" ? HcalForward
240  : HcalOther;
241 
242  HcalDetId id(subdet, ieta, iphi, idep);
243  lutgain[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
244  }
245 
246  //Effective Pedestals
247  std::ifstream infeffped(
249  Form("%s/EffectivePedestals/EffectivePedestals_Run%s.txt", inputDir.c_str(), effpedestals_[ii].c_str()))
250  .fullPath()
251  .c_str());
252  assert(!infeffped.fail());
253  while (!infeffped.eof()) {
254  infeffped.getline(buffer, 1024);
255  if (buffer[0] == '#')
256  continue;
257  std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> val2 >> val3 >> val4 >> wid1 >> wid2 >>
258  wid3 >> wid4 >> iraw;
259  if (det != "HB" && det != "HE" && det != "HF")
260  continue;
261 
262  float theval = (val1 + val2 + val3 + val4) / 4.0;
263 
264  HcalSubdetector subdet = det == "HB" ? HcalBarrel
265  : det == "HE" ? HcalEndcap
266  : det == "HF" ? HcalForward
267  : HcalOther;
268 
269  HcalDetId id(subdet, ieta, iphi, idep);
270  luteffpede[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
271  }
272 
273  //Pedestals
274  std::ifstream infped(
275  edm::FileInPath(Form("%s/Pedestals/Pedestals_Run%s.txt", inputDir.c_str(), pedestals_[ii].c_str()))
276  .fullPath()
277  .c_str());
278  assert(!infped.fail());
279  while (!infped.eof()) {
280  infped.getline(buffer, 1024);
281  if (buffer[0] == '#')
282  continue;
283  std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> val2 >> val3 >> val4 >> wid1 >> wid2 >>
284  wid3 >> wid4 >> iraw;
285  if (det != "HB" && det != "HE" && det != "HF")
286  continue;
287 
288  float theval = (val1 + val2 + val3 + val4) / 4.0;
289 
290  HcalSubdetector subdet = det == "HB" ? HcalBarrel
291  : det == "HE" ? HcalEndcap
292  : det == "HF" ? HcalForward
293  : HcalOther;
294 
295  HcalDetId id(subdet, ieta, iphi, idep);
296  lutpede[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
297  }
298 
299  //Response corrections
300  std::ifstream inresp(
301  edm::FileInPath(Form("%s/RespCorrs/RespCorrs_Run%s.txt", inputDir.c_str(), respcorrs_[ii].c_str()))
302  .fullPath()
303  .c_str());
304  assert(!inresp.fail());
305  while (!inresp.eof()) {
306  inresp.getline(buffer, 1024);
307  if (buffer[0] == '#')
308  continue;
309  std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> iraw;
310  if (det != "HB" && det != "HE" && det != "HF")
311  continue;
312 
313  float theval = val1;
314 
315  HcalSubdetector subdet = det == "HB" ? HcalBarrel
316  : det == "HE" ? HcalEndcap
317  : det == "HF" ? HcalForward
318  : HcalOther;
319 
320  HcalDetId id(subdet, ieta, iphi, idep);
321  lutresp[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
322  }
323 
324  //ChannelQuality
325  std::ifstream inchan(
326  edm::FileInPath(Form("%s/ChannelQuality/ChannelQuality_Run%s.txt", inputDir.c_str(), quality_[ii].c_str()))
327  .fullPath()
328  .c_str());
329  assert(!inchan.fail());
330  while (!inchan.eof()) {
331  inchan.getline(buffer, 1024);
332  if (buffer[0] == '#')
333  continue;
334  std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> base >> val1 >> iraw;
335 
336  float theval = val1;
337 
338  HcalSubdetector subdet = det == "HB" ? HcalBarrel
339  : det == "HE" ? HcalEndcap
340  : det == "HF" ? HcalForward
341  : det == "HO" ? HcalOuter
342  : HcalOther;
343 
344  HcalDetId id(subdet, ieta, iphi, idep);
345  if (theval != 0)
346  BadChans[ii].push_back(id);
347  }
348  }
349 
350  LutXml xmls1(edm::FileInPath(Form("%s/%s/%s.xml", inputDir.c_str(), tags_[0].c_str(), tags_[0].c_str())).fullPath());
351  LutXml xmls2(edm::FileInPath(Form("%s/%s/%s.xml", inputDir.c_str(), tags_[1].c_str(), tags_[1].c_str())).fullPath());
352 
353  xmls1.create_lut_map();
354  xmls2.create_lut_map();
355 
356  for (const auto& xml2 : xmls2) {
357  HcalGenericDetId detid(xml2.first);
358 
359  if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower) {
360  HcalTrigTowerDetId tid(detid.rawId());
361  if (!topology->validHT(tid))
362  continue;
363  const auto& lut2 = xml2.second;
364 
365  int D = abs(tid.ieta()) < 29 ? (lut2.size() == 1024 ? 0 : 3) : tid.version() == 0 ? 1 : 2;
366  for (size_t i = 0; i < lut2.size(); ++i) {
367  if (int(i) % 4 == D)
368  houtput[D][1]->Fill(i, lut2[i]);
369  }
370  } else if (topology->valid(detid)) {
371  HcalDetId id(detid);
372  HcalSubdetector subdet = id.subdet();
373  int idet = int(subdet);
374  const auto& lut2 = xml2.second;
375  int hbhe = idet == HcalForward ? 1 : idet == HcalOuter ? 3 : lut2.size() == 128 ? 0 : 2;
376  for (size_t i = 0; i < lut2.size(); ++i) {
377  hlut[hbhe][1]->Fill(i, lut2[i]);
378  if (hbhe == 2)
379  hlut[hbhe][1]->Fill(i, lut2[i] & 0x3FF);
380  }
381  }
382  }
383 
384  for (const auto& xml1 : xmls1) {
385  HcalGenericDetId detid(xml1.first);
386  const auto& lut1 = xml1.second;
387 
388  if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower) {
389  HcalTrigTowerDetId tid(detid.rawId());
390  if (!topology->validHT(tid))
391  continue;
392  int D = abs(tid.ieta()) < 29 ? (lut1.size() == 1024 ? 0 : 3) : tid.version() == 0 ? 1 : 2;
393  for (size_t i = 0; i < lut1.size(); ++i) {
394  if (int(i) % 4 == D)
395  houtput[D][0]->Fill(i, lut1[i]);
396  }
397  } else if (topology->valid(detid)) {
398  HcalDetId id(detid);
399  HcalSubdetector subdet = id.subdet();
400  int idet = int(subdet);
401  const auto& lut1 = xml1.second;
402  int hbhe = idet == HcalForward ? 1 : idet == HcalOuter ? 3 : lut1.size() == 128 ? 0 : 2;
403  for (size_t i = 0; i < lut1.size(); ++i) {
404  hlut[hbhe][0]->Fill(i, lut1[i]);
405  if (hbhe == 2)
406  hlut[hbhe][0]->Fill(i, lut1[i] & 0x3FF);
407  }
408  }
409 
410  auto xml2 = xmls2.find(detid.rawId());
411  if (xml2 == xmls2.end())
412  continue;
413 
414  if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower)
415  continue;
416 
417  HcalDetId id(detid);
418 
419  HcalSubdetector subdet = id.subdet();
420  int idet = int(subdet);
421  int ieta = id.ieta();
422  int iphi = id.iphi();
423  int idep = id.depth() - 1;
424  unsigned long int iraw = id.rawId();
425 
426  if (!topology->valid(detid))
427  continue;
428 
429  int hbhe = idet == HcalForward;
430 
431  const auto& lut2 = xml2->second;
432 
433  size_t size = lut1.size();
434  if (size != lut2.size())
435  continue;
436 
437  std::vector<unsigned int> llut1(size);
438  std::vector<unsigned int> llut2(size);
439  for (size_t i = 0; i < size; ++i) {
440  llut1[i] = hbhe == 0 ? lut1[i] & 0x3FF : lut1[i];
441  llut2[i] = hbhe == 0 ? lut2[i] & 0x3FF : lut2[i];
442  }
443 
444  int threshold[2] = {0, 0};
445  //Thresholds
446  for (size_t i = 0; i < size; ++i) {
447  if (llut1[i] > 0) {
448  threshold[0] = i;
449  break;
450  }
451  if (i == size - 1) {
452  ZeroLuts[0].push_back(id);
453  }
454  }
455  for (size_t i = 0; i < size; ++i) {
456  if (llut2[i] > 0) {
457  threshold[1] = i;
458  break;
459  }
460  if (i == size - 1) {
461  ZeroLuts[1].push_back(id);
462  }
463  }
464 
465  if (subdet != HcalBarrel && subdet != HcalEndcap && subdet != HcalForward)
466  continue;
467 
468  //fill conditions
469 
470  double xfill = 0;
471  h[0][0][hbhe][idep]->Fill(ieta, iphi, lutpede[0][iraw]);
472  h[0][1][hbhe][idep]->Fill(ieta, iphi, lutpede[1][iraw]);
473  xfill = lutpede[1][iraw] - lutpede[0][iraw];
474  h[0][2][hbhe][idep]->Fill(ieta, iphi, xfill);
475  r[0][hbhe]->Fill(ieta, xfill);
476  p[0][hbhe]->Fill(ieta, xfill);
477 
478  h[1][0][hbhe][idep]->Fill(ieta, iphi, lutresp[0][iraw]);
479  h[1][1][hbhe][idep]->Fill(ieta, iphi, lutresp[1][iraw]);
480  xfill = lutresp[1][iraw] / lutresp[0][iraw];
481  h[1][2][hbhe][idep]->Fill(ieta, iphi, xfill);
482  r[1][hbhe]->Fill(ieta, xfill);
483  p[1][hbhe]->Fill(ieta, xfill);
484 
485  h[2][0][hbhe][idep]->Fill(ieta, iphi, lutgain[0][iraw]);
486  h[2][1][hbhe][idep]->Fill(ieta, iphi, lutgain[1][iraw]);
487  xfill = lutgain[1][iraw] / lutgain[0][iraw];
488  h[2][2][hbhe][idep]->Fill(ieta, iphi, xfill);
489  r[2][hbhe]->Fill(ieta, xfill);
490  p[2][hbhe]->Fill(ieta, xfill);
491 
492  h[3][0][hbhe][idep]->Fill(ieta, iphi, threshold[0]);
493  h[3][1][hbhe][idep]->Fill(ieta, iphi, threshold[1]);
494  xfill = threshold[1] - threshold[0];
495  h[3][2][hbhe][idep]->Fill(ieta, iphi, xfill);
496  r[3][hbhe]->Fill(ieta, xfill);
497  p[3][hbhe]->Fill(ieta, xfill);
498 
499  h[5][0][hbhe][idep]->Fill(ieta, iphi, luteffpede[0][iraw]);
500  h[5][1][hbhe][idep]->Fill(ieta, iphi, luteffpede[1][iraw]);
501  xfill = luteffpede[1][iraw] - luteffpede[0][iraw];
502  h[5][2][hbhe][idep]->Fill(ieta, iphi, xfill);
503  r[5][hbhe]->Fill(ieta, xfill);
504  p[5][hbhe]->Fill(ieta, xfill);
505 
506  size_t maxvalue = hbhe == 0 ? 1023 : 2047;
507 
508  //LutDifference
509  for (size_t i = 0; i < size; ++i) {
510  hlut[hbhe][2]->Fill(i, llut1[i] == 0 ? 0 : (double)llut2[i] / llut1[i]);
511 
512  if (i == size - 1 ||
513  (llut1[i] == maxvalue || llut2[i] == maxvalue)) { //Fill with only the last one before the maximum
514  if (llut1[i - 1] == 0 || llut2[i - 1] == 0) {
515  break;
516  }
517  double condratio = lutgain[1][iraw] / lutgain[0][iraw] * lutresp[1][iraw] / lutresp[0][iraw];
518  xfill = (double)llut2[i - 1] / llut1[i - 1];
519  hslope[hbhe]->Fill(condratio, xfill);
520 
521  h[4][0][hbhe][idep]->Fill(ieta, iphi, (double)llut1[i - 1] / (i - 1));
522  h[4][1][hbhe][idep]->Fill(ieta, iphi, (double)llut2[i - 1] / (i - 1));
523  h[4][2][hbhe][idep]->Fill(ieta, iphi, xfill);
524  r[4][hbhe]->Fill(ieta, xfill);
525  p[4][hbhe]->Fill(ieta, xfill);
526 
527  break;
528  }
529  }
530  }
531 
532  gROOT->SetStyle("Plain");
533  gStyle->SetPalette(1);
534  gStyle->SetStatW(0.2);
535  gStyle->SetStatH(0.1);
536  gStyle->SetStatY(1.0);
537  gStyle->SetStatX(0.9);
538  gStyle->SetOptStat(110010);
539  gStyle->SetOptFit(111111);
540 
541  //Draw and Print
542  TCanvas* cc = new TCanvas("cc", "cc", 0, 0, 1600, 1200);
543  cc->SetGridy();
544  for (int j = 0; j < NVAR; ++j) {
545  gSystem->mkdir(TString(plotsDir) + "/_" + titleVar[j]);
546  for (int d = 0; d < NDET; ++d) {
547  cc->Clear();
548  r[j][d]->Draw("colz");
549  cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(r[j][d]->GetName()) + ".pdf");
550 
551  cc->Clear();
552  p[j][d]->Draw();
553  if (doRatio[j]) {
554  p[j][d]->SetMinimum(Pmin);
555  p[j][d]->SetMaximum(Pmax);
556  }
557  cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(p[j][d]->GetName()) + ".pdf");
558 
559  for (int i = 0; i < NLEV; ++i) {
560  for (int p = 0; p < DEP[d]; ++p) {
561  cc->Clear();
562  h[j][i][d][p]->Draw("colz");
563 
564  if (i == NLEV - 1) {
565  if (doRatio[j]) {
566  h[j][2][d][p]->SetMinimum(Zmin);
567  h[j][2][d][p]->SetMaximum(Zmax);
568  } else {
569  h[j][2][d][p]->SetMinimum(-3);
570  h[j][2][d][p]->SetMaximum(3);
571  }
572  }
573  cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(h[j][i][d][p]->GetName()) + ".pdf");
574  }
575  }
576  }
577  }
578 
579  for (int i = 0; i < NLEV; ++i) {
580  cc->Clear();
581  hlut[0][i]->Draw();
582  hlut[1][i]->Draw("sames");
583  hlut[2][i]->Draw("sames");
584  hlut[3][i]->Draw("sames");
585  cc->Print(TString(plotsDir) + Form("LUT_%d.gif", i));
586  }
587  cc->Clear();
588  hslope[0]->Draw("colz");
589  cc->SetGridx();
590  cc->SetGridy();
591  cc->Print(TString(plotsDir) + "GainLutScatterHBHE.pdf");
592  cc->Clear();
593  hslope[1]->Draw("colz");
594  cc->SetGridx();
595  cc->SetGridy();
596  cc->Print(TString(plotsDir) + "GainLutScatterLutHF.pdf");
597 
598  for (int i = 0; i < 2; ++i) {
599  cc->Clear();
600  houtput[0][i]->Draw("box");
601  houtput[1][i]->Draw("samebox");
602  houtput[2][i]->Draw("samebox");
603  houtput[3][i]->Draw("samebox");
604  cc->Print(TString(plotsDir) + Form("OUT_%d.gif", i));
605  }
606 }
size
Write out results.
bool validHT(const HcalTrigTowerDetId &id) const
Definition: LutXml.h:27
std::string fullPath() const
Definition: FileInPath.cc:161
vector< string > vstring
Definition: ExoticaDQM.cc:8
bool valid(const DetId &id) const override
base
Main Program
Definition: newFWLiteAna.py:92
int create_lut_map(void)
Definition: LutXml.cc:338
std::string plotsDir
std::vector< std::string > effpedestals_
assert(be >=bs)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::string > quality_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::vector< std::string > gains_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > tok_htopo_
d
Definition: ztail.py:151
std::string inputDir
ii
Definition: cuy.py:589
std::vector< std::string > pedestals_
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
std::vector< std::string > respcorrs_
std::vector< std::string > tags_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4

◆ fillDescriptions()

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

Definition at line 608 of file HcalLutAnalyzer.cc.

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

608  {
610  desc.setUnknown();
611  descriptions.addDefault(desc);
612 }
void addDefault(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ effpedestals_

std::vector<std::string> HcalLutAnalyzer::effpedestals_
private

Definition at line 63 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ gains_

std::vector<std::string> HcalLutAnalyzer::gains_
private

Definition at line 64 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ inputDir

std::string HcalLutAnalyzer::inputDir
private

Definition at line 58 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ pedestals_

std::vector<std::string> HcalLutAnalyzer::pedestals_
private

Definition at line 62 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ plotsDir

std::string HcalLutAnalyzer::plotsDir
private

Definition at line 59 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Pmax

double HcalLutAnalyzer::Pmax
private

Definition at line 72 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Pmin

double HcalLutAnalyzer::Pmin
private

Definition at line 71 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ quality_

std::vector<std::string> HcalLutAnalyzer::quality_
private

Definition at line 61 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ respcorrs_

std::vector<std::string> HcalLutAnalyzer::respcorrs_
private

Definition at line 65 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ tags_

std::vector<std::string> HcalLutAnalyzer::tags_
private

Definition at line 60 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ tok_htopo_

edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> HcalLutAnalyzer::tok_htopo_
private

Definition at line 74 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Ymax

double HcalLutAnalyzer::Ymax
private

Definition at line 70 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Ymin

double HcalLutAnalyzer::Ymin
private

Definition at line 69 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Zmax

double HcalLutAnalyzer::Zmax
private

Definition at line 68 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Zmin

double HcalLutAnalyzer::Zmin
private

Definition at line 67 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().