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 > 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 76 of file HcalLutAnalyzer.cc.

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

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

References funct::abs(), cms::cuda::assert(), newFWLiteAna::base, edmScanValgrind::buffer, LutXml::create_lut_map(), ztail::d, histoStyle::doRatio, 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, triggerObjects_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.

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

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

570  {
572  desc.setUnknown();
573  descriptions.addDefault(desc);
574 }
void addDefault(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ gains_

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

Definition at line 63 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 71 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Pmin

double HcalLutAnalyzer::Pmin
private

Definition at line 70 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 64 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 73 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Ymax

double HcalLutAnalyzer::Ymax
private

Definition at line 69 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Ymin

double HcalLutAnalyzer::Ymin
private

Definition at line 68 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Zmax

double HcalLutAnalyzer::Zmax
private

Definition at line 67 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().

◆ Zmin

double HcalLutAnalyzer::Zmin
private

Definition at line 66 of file HcalLutAnalyzer.cc.

Referenced by analyze(), and HcalLutAnalyzer().