File indexing completed on 2021-06-16 03:19:28
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <iostream>
0022 #include <fstream>
0023
0024
0025 #include "FWCore/Framework/interface/Frameworkfwd.h"
0026 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Framework/interface/EventSetup.h"
0031
0032 #include "CalibCalorimetry/HcalTPGAlgos/interface/XMLProcessor.h"
0033 #include "CalibCalorimetry/HcalTPGAlgos/interface/LutXml.h"
0034 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0035 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0036 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0037 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0038 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0039
0040 #include "TString.h"
0041 #include "TH1D.h"
0042 #include "TH2D.h"
0043 #include "TProfile.h"
0044 #include "TCanvas.h"
0045 #include "TROOT.h"
0046 #include "TStyle.h"
0047 #include "TSystem.h"
0048
0049 class HcalLutAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0050 public:
0051 explicit HcalLutAnalyzer(const edm::ParameterSet&);
0052 ~HcalLutAnalyzer() override{};
0053 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0054
0055 private:
0056 void analyze(const edm::Event&, const edm::EventSetup&) override;
0057
0058 std::string inputDir;
0059 std::string plotsDir;
0060 std::vector<std::string> tags_;
0061 std::vector<std::string> quality_;
0062 std::vector<std::string> pedestals_;
0063 std::vector<std::string> gains_;
0064 std::vector<std::string> respcorrs_;
0065
0066 double Zmin;
0067 double Zmax;
0068 double Ymin;
0069 double Ymax;
0070 double Pmin;
0071 double Pmax;
0072
0073 edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tok_htopo_;
0074 };
0075
0076 HcalLutAnalyzer::HcalLutAnalyzer(const edm::ParameterSet& iConfig) {
0077 inputDir = iConfig.getParameter<std::string>("inputDir");
0078 plotsDir = iConfig.getParameter<std::string>("plotsDir");
0079 tags_ = iConfig.getParameter<std::vector<std::string> >("tags");
0080 quality_ = iConfig.getParameter<std::vector<std::string> >("quality");
0081 pedestals_ = iConfig.getParameter<std::vector<std::string> >("pedestals");
0082 gains_ = iConfig.getParameter<std::vector<std::string> >("gains");
0083 respcorrs_ = iConfig.getParameter<std::vector<std::string> >("respcorrs");
0084
0085 Zmin = iConfig.getParameter<double>("Zmin");
0086 Zmax = iConfig.getParameter<double>("Zmax");
0087 Ymin = iConfig.getParameter<double>("Ymin");
0088 Ymax = iConfig.getParameter<double>("Ymax");
0089 Pmin = iConfig.getParameter<double>("Pmin");
0090 Pmax = iConfig.getParameter<double>("Pmax");
0091
0092 tok_htopo_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
0093 }
0094
0095 void HcalLutAnalyzer::analyze(const edm::Event&, const edm::EventSetup& iSetup) {
0096 using namespace std;
0097
0098 const HcalTopology* topology = &iSetup.getData(tok_htopo_);
0099
0100 typedef std::vector<std::string> vstring;
0101 typedef std::map<unsigned long int, float> LUTINPUT;
0102
0103 static const int NVAR = 5;
0104 static const int NDET = 2;
0105 static const int NDEP = 7;
0106 static const int NLEV = 3;
0107
0108 const bool doRatio[NVAR] = {false, true, true, false, true};
0109 const char* titleVar[NVAR] = {"Pedestals", "RespCorrs", "Gains", "Threshold", "LUTs"};
0110 const char* titleHisR[NLEV] = {"Old", "New", "Ratio"};
0111 const char* titleHisD[NLEV] = {"Old", "New", "Difference"};
0112 const char* titleDet[4] = {"HBHE", "HF", "HEP17", "HO"};
0113 const int DEP[NDET] = {7, 4};
0114 const char* titleDep[NDEP] = {"depth1", "depth2", "depth3", "depth4", "depth5", "depth6", "depth7"};
0115
0116 TH2D* r[NVAR][NDET];
0117 TProfile* p[NVAR][NDET];
0118
0119 TH2D* h[NVAR][NLEV][NDET][NDEP];
0120 TH2D* hlut[4][NLEV];
0121 TH2D* hslope[2];
0122 TH2D* houtput[4][2];
0123
0124 for (int d = 0; d < 4; ++d) {
0125 for (int i = 0; i < 3; ++i) {
0126 hlut[d][i] = new TH2D(Form("Lut_%s_%s", titleDet[d], titleHisR[i]),
0127 Form("Input LUT, %s", (i == 2 ? "Ratio" : tags_[i].c_str())),
0128 260,
0129 0,
0130 260,
0131 240,
0132 0,
0133 i == NLEV - 1 ? 3 : 2400);
0134 hlut[d][i]->SetMarkerColor(d == 0 ? kBlue : d == 1 ? kGreen + 2 : d == 2 ? kRed : kCyan);
0135 hlut[d][i]->SetXTitle("raw adc");
0136 hlut[d][i]->SetYTitle("lin adc");
0137 }
0138 }
0139
0140 for (int d = 0; d < NDET; ++d) {
0141 hslope[d] = new TH2D(
0142 Form("GainLutScatter_%s", titleDet[d]), Form("Gain-Lutslope scatter, %s", titleDet[d]), 200, 0, 2, 200, 0, 2);
0143 hslope[d]->SetXTitle("Gain x RespCorr ratio");
0144 hslope[d]->SetYTitle("Lut ratio");
0145
0146 for (int j = 0; j < NVAR; ++j) {
0147 double rmin = doRatio[j] ? Ymin : -6;
0148 double rmax = doRatio[j] ? Ymax : 6;
0149 r[j][d] = new TH2D(Form("r%s_%s", titleVar[j], titleDet[d]),
0150 Form("%s, %s", titleVar[j], titleDet[d]),
0151 83,
0152 -41.5,
0153 41.5,
0154 250,
0155 rmin,
0156 rmax);
0157 r[j][d]->SetXTitle("iEta");
0158 r[j][d]->SetYTitle(doRatio[j] ? "New / Old" : "New - Old");
0159 p[j][d] = new TProfile(
0160 Form("p%s_%s", titleVar[j], titleDet[d]), Form("%s, %s", titleVar[j], titleDet[d]), 83, -41.5, 41.5);
0161 p[j][d]->SetXTitle("iEta");
0162 p[j][d]->SetYTitle(doRatio[j] ? "New / Old" : "New - Old");
0163 p[j][d]->SetMarkerStyle(20);
0164 p[j][d]->SetMarkerSize(0.9);
0165 p[j][d]->SetMarkerColor(kBlue);
0166
0167 for (int p = 0; p < DEP[d]; ++p) {
0168 for (int i = 0; i < NLEV; ++i) {
0169 const char* titHist = doRatio[j] ? titleHisR[i] : titleHisD[i];
0170 h[j][i][d][p] = new TH2D(Form("h%s_%s_%s_%s", titleVar[j], titHist, titleDet[d], titleDep[p]),
0171 Form("%s, %s, %s, %s", titleVar[j], titHist, titleDet[d], titleDep[p]),
0172 83,
0173 -41.5,
0174 41.5,
0175 72,
0176 0.5,
0177 72.5);
0178 h[j][i][d][p]->SetXTitle("iEta");
0179 h[j][i][d][p]->SetYTitle("iPhi");
0180 }
0181 }
0182 }
0183 }
0184
0185 for (int i = 0; i < 4; ++i) {
0186 int color = i == 0 ? kBlue : i == 1 ? kViolet : i == 2 ? kGreen + 2 : kRed;
0187 houtput[i][0] =
0188 new TH2D(Form("houtlut0_%d", i), Form("Output LUT, %s", tags_[0].c_str()), 2100, 0, 2100, 260, 0, 260);
0189 houtput[i][1] =
0190 new TH2D(Form("houtlut1_%d", i), Form("Output LUT, %s", tags_[1].c_str()), 2100, 0, 2100, 260, 0, 260);
0191 for (int j = 0; j < 2; ++j) {
0192 houtput[i][j]->SetMarkerColor(color);
0193 houtput[i][j]->SetLineColor(color);
0194 }
0195 }
0196
0197
0198 LUTINPUT lutgain[2];
0199 LUTINPUT lutresp[2];
0200 LUTINPUT lutpede[2];
0201
0202 assert(tags_.size() == 2);
0203 assert(gains_.size() == 2);
0204 assert(respcorrs_.size() == 2);
0205 assert(pedestals_.size() == 2);
0206
0207 unsigned long int iraw;
0208 int ieta, iphi, idep;
0209 string det, base;
0210 float val1, val2, val3, val4;
0211 float wid1, wid2, wid3, wid4;
0212 char buffer[1024];
0213
0214 std::vector<HcalDetId> BadChans[2];
0215 std::vector<HcalDetId> ZeroLuts[2];
0216
0217
0218 for (int ii = 0; ii < 2; ++ii) {
0219
0220 std::ifstream infile(
0221 edm::FileInPath(Form("%s/Gains/Gains_Run%s.txt", inputDir.c_str(), gains_[ii].c_str())).fullPath().c_str());
0222 assert(!infile.fail());
0223 while (!infile.eof()) {
0224 infile.getline(buffer, 1024);
0225 if (buffer[0] == '#')
0226 continue;
0227 std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> val2 >> val3 >> val4 >> iraw;
0228 if (det != "HB" && det != "HE" && det != "HF")
0229 continue;
0230
0231 float theval = (val1 + val2 + val3 + val4) / 4.0;
0232
0233 HcalSubdetector subdet = det == "HB" ? HcalBarrel
0234 : det == "HE" ? HcalEndcap
0235 : det == "HF" ? HcalForward
0236 : HcalOther;
0237
0238 HcalDetId id(subdet, ieta, iphi, idep);
0239 lutgain[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
0240 }
0241
0242
0243 std::ifstream infped(
0244 edm::FileInPath(Form("%s/Pedestals/Pedestals_Run%s.txt", inputDir.c_str(), pedestals_[ii].c_str()))
0245 .fullPath()
0246 .c_str());
0247 assert(!infped.fail());
0248 while (!infped.eof()) {
0249 infped.getline(buffer, 1024);
0250 if (buffer[0] == '#')
0251 continue;
0252 std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> val2 >> val3 >> val4 >> wid1 >> wid2 >>
0253 wid3 >> wid4 >> iraw;
0254 if (det != "HB" && det != "HE" && det != "HF")
0255 continue;
0256
0257 float theval = (val1 + val2 + val3 + val4) / 4.0;
0258
0259 HcalSubdetector subdet = det == "HB" ? HcalBarrel
0260 : det == "HE" ? HcalEndcap
0261 : det == "HF" ? HcalForward
0262 : HcalOther;
0263
0264 HcalDetId id(subdet, ieta, iphi, idep);
0265 lutpede[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
0266 }
0267
0268
0269 std::ifstream inresp(
0270 edm::FileInPath(Form("%s/RespCorrs/RespCorrs_Run%s.txt", inputDir.c_str(), respcorrs_[ii].c_str()))
0271 .fullPath()
0272 .c_str());
0273 assert(!inresp.fail());
0274 while (!inresp.eof()) {
0275 inresp.getline(buffer, 1024);
0276 if (buffer[0] == '#')
0277 continue;
0278 std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> val1 >> iraw;
0279 if (det != "HB" && det != "HE" && det != "HF")
0280 continue;
0281
0282 float theval = val1;
0283
0284 HcalSubdetector subdet = det == "HB" ? HcalBarrel
0285 : det == "HE" ? HcalEndcap
0286 : det == "HF" ? HcalForward
0287 : HcalOther;
0288
0289 HcalDetId id(subdet, ieta, iphi, idep);
0290 lutresp[ii].insert(LUTINPUT::value_type(id.rawId(), theval));
0291 }
0292
0293
0294 std::ifstream inchan(
0295 edm::FileInPath(Form("%s/ChannelQuality/ChannelQuality_Run%s.txt", inputDir.c_str(), quality_[ii].c_str()))
0296 .fullPath()
0297 .c_str());
0298 assert(!inchan.fail());
0299 while (!inchan.eof()) {
0300 inchan.getline(buffer, 1024);
0301 if (buffer[0] == '#')
0302 continue;
0303 std::istringstream(buffer) >> ieta >> iphi >> idep >> det >> base >> val1 >> iraw;
0304
0305 float theval = val1;
0306
0307 HcalSubdetector subdet = det == "HB" ? HcalBarrel
0308 : det == "HE" ? HcalEndcap
0309 : det == "HF" ? HcalForward
0310 : det == "HO" ? HcalOuter
0311 : HcalOther;
0312
0313 HcalDetId id(subdet, ieta, iphi, idep);
0314 if (theval != 0)
0315 BadChans[ii].push_back(id);
0316 }
0317 }
0318
0319 LutXml xmls1(edm::FileInPath(Form("%s/%s/%s.xml", inputDir.c_str(), tags_[0].c_str(), tags_[0].c_str())).fullPath());
0320 LutXml xmls2(edm::FileInPath(Form("%s/%s/%s.xml", inputDir.c_str(), tags_[1].c_str(), tags_[1].c_str())).fullPath());
0321
0322 xmls1.create_lut_map();
0323 xmls2.create_lut_map();
0324
0325 for (const auto& xml2 : xmls2) {
0326 HcalGenericDetId detid(xml2.first);
0327
0328 if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower) {
0329 HcalTrigTowerDetId tid(detid.rawId());
0330 if (!topology->validHT(tid))
0331 continue;
0332 const auto& lut2 = xml2.second;
0333
0334 int D = abs(tid.ieta()) < 29 ? (lut2.size() == 1024 ? 0 : 3) : tid.version() == 0 ? 1 : 2;
0335 for (size_t i = 0; i < lut2.size(); ++i) {
0336 if (int(i) % 4 == D)
0337 houtput[D][1]->Fill(i, lut2[i]);
0338 }
0339 } else if (topology->valid(detid)) {
0340 HcalDetId id(detid);
0341 HcalSubdetector subdet = id.subdet();
0342 int idet = int(subdet);
0343 const auto& lut2 = xml2.second;
0344 int hbhe = idet == HcalForward ? 1 : idet == HcalOuter ? 3 : lut2.size() == 128 ? 0 : 2;
0345 for (size_t i = 0; i < lut2.size(); ++i) {
0346 hlut[hbhe][1]->Fill(i, lut2[i]);
0347 if (hbhe == 2)
0348 hlut[hbhe][1]->Fill(i, lut2[i] & 0x3FF);
0349 }
0350 }
0351 }
0352
0353 for (const auto& xml1 : xmls1) {
0354 HcalGenericDetId detid(xml1.first);
0355 const auto& lut1 = xml1.second;
0356
0357 if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower) {
0358 HcalTrigTowerDetId tid(detid.rawId());
0359 if (!topology->validHT(tid))
0360 continue;
0361 int D = abs(tid.ieta()) < 29 ? (lut1.size() == 1024 ? 0 : 3) : tid.version() == 0 ? 1 : 2;
0362 for (size_t i = 0; i < lut1.size(); ++i) {
0363 if (int(i) % 4 == D)
0364 houtput[D][0]->Fill(i, lut1[i]);
0365 }
0366 } else if (topology->valid(detid)) {
0367 HcalDetId id(detid);
0368 HcalSubdetector subdet = id.subdet();
0369 int idet = int(subdet);
0370 const auto& lut1 = xml1.second;
0371 int hbhe = idet == HcalForward ? 1 : idet == HcalOuter ? 3 : lut1.size() == 128 ? 0 : 2;
0372 for (size_t i = 0; i < lut1.size(); ++i) {
0373 hlut[hbhe][0]->Fill(i, lut1[i]);
0374 if (hbhe == 2)
0375 hlut[hbhe][0]->Fill(i, lut1[i] & 0x3FF);
0376 }
0377 }
0378
0379 auto xml2 = xmls2.find(detid.rawId());
0380 if (xml2 == xmls2.end())
0381 continue;
0382
0383 if (detid.genericSubdet() == HcalGenericDetId::HcalGenTriggerTower)
0384 continue;
0385
0386 HcalDetId id(detid);
0387
0388 HcalSubdetector subdet = id.subdet();
0389 int idet = int(subdet);
0390 int ieta = id.ieta();
0391 int iphi = id.iphi();
0392 int idep = id.depth() - 1;
0393 unsigned long int iraw = id.rawId();
0394
0395 if (!topology->valid(detid))
0396 continue;
0397
0398 int hbhe = idet == HcalForward;
0399
0400 const auto& lut2 = xml2->second;
0401
0402 size_t size = lut1.size();
0403 if (size != lut2.size())
0404 continue;
0405
0406 std::vector<unsigned int> llut1(size);
0407 std::vector<unsigned int> llut2(size);
0408 for (size_t i = 0; i < size; ++i) {
0409 llut1[i] = hbhe == 0 ? lut1[i] & 0x3FF : lut1[i];
0410 llut2[i] = hbhe == 0 ? lut2[i] & 0x3FF : lut2[i];
0411 }
0412
0413 int threshold[2] = {0, 0};
0414
0415 for (size_t i = 0; i < size; ++i) {
0416 if (llut1[i] > 0) {
0417 threshold[0] = i;
0418 break;
0419 }
0420 if (i == size - 1) {
0421 ZeroLuts[0].push_back(id);
0422 }
0423 }
0424 for (size_t i = 0; i < size; ++i) {
0425 if (llut2[i] > 0) {
0426 threshold[1] = i;
0427 break;
0428 }
0429 if (i == size - 1) {
0430 ZeroLuts[1].push_back(id);
0431 }
0432 }
0433
0434 if (subdet != HcalBarrel && subdet != HcalEndcap && subdet != HcalForward)
0435 continue;
0436
0437
0438
0439 double xfill = 0;
0440 h[0][0][hbhe][idep]->Fill(ieta, iphi, lutpede[0][iraw]);
0441 h[0][1][hbhe][idep]->Fill(ieta, iphi, lutpede[1][iraw]);
0442 xfill = lutpede[1][iraw] - lutpede[0][iraw];
0443 h[0][2][hbhe][idep]->Fill(ieta, iphi, xfill);
0444 r[0][hbhe]->Fill(ieta, xfill);
0445 p[0][hbhe]->Fill(ieta, xfill);
0446
0447 h[1][0][hbhe][idep]->Fill(ieta, iphi, lutresp[0][iraw]);
0448 h[1][1][hbhe][idep]->Fill(ieta, iphi, lutresp[1][iraw]);
0449 xfill = lutresp[1][iraw] / lutresp[0][iraw];
0450 h[1][2][hbhe][idep]->Fill(ieta, iphi, xfill);
0451 r[1][hbhe]->Fill(ieta, xfill);
0452 p[1][hbhe]->Fill(ieta, xfill);
0453
0454 h[2][0][hbhe][idep]->Fill(ieta, iphi, lutgain[0][iraw]);
0455 h[2][1][hbhe][idep]->Fill(ieta, iphi, lutgain[1][iraw]);
0456 xfill = lutgain[1][iraw] / lutgain[0][iraw];
0457 h[2][2][hbhe][idep]->Fill(ieta, iphi, xfill);
0458 r[2][hbhe]->Fill(ieta, xfill);
0459 p[2][hbhe]->Fill(ieta, xfill);
0460
0461 h[3][0][hbhe][idep]->Fill(ieta, iphi, threshold[0]);
0462 h[3][1][hbhe][idep]->Fill(ieta, iphi, threshold[1]);
0463 xfill = threshold[1] - threshold[0];
0464 h[3][2][hbhe][idep]->Fill(ieta, iphi, xfill);
0465 r[3][hbhe]->Fill(ieta, xfill);
0466 p[3][hbhe]->Fill(ieta, xfill);
0467
0468 size_t maxvalue = hbhe == 0 ? 1023 : 2047;
0469
0470
0471 for (size_t i = 0; i < size; ++i) {
0472 hlut[hbhe][2]->Fill(i, llut1[i] == 0 ? 0 : (double)llut2[i] / llut1[i]);
0473
0474 if (i == size - 1 ||
0475 (llut1[i] == maxvalue || llut2[i] == maxvalue)) {
0476 if (llut1[i - 1] == 0 || llut2[i - 1] == 0) {
0477 break;
0478 }
0479 double condratio = lutgain[1][iraw] / lutgain[0][iraw] * lutresp[1][iraw] / lutresp[0][iraw];
0480 xfill = (double)llut2[i - 1] / llut1[i - 1];
0481 hslope[hbhe]->Fill(condratio, xfill);
0482
0483 h[4][0][hbhe][idep]->Fill(ieta, iphi, (double)llut1[i - 1] / (i - 1));
0484 h[4][1][hbhe][idep]->Fill(ieta, iphi, (double)llut2[i - 1] / (i - 1));
0485 h[4][2][hbhe][idep]->Fill(ieta, iphi, xfill);
0486 r[4][hbhe]->Fill(ieta, xfill);
0487 p[4][hbhe]->Fill(ieta, xfill);
0488
0489 break;
0490 }
0491 }
0492 }
0493
0494 gROOT->SetStyle("Plain");
0495 gStyle->SetPalette(1);
0496 gStyle->SetStatW(0.2);
0497 gStyle->SetStatH(0.1);
0498 gStyle->SetStatY(1.0);
0499 gStyle->SetStatX(0.9);
0500 gStyle->SetOptStat(110010);
0501 gStyle->SetOptFit(111111);
0502
0503
0504 TCanvas* cc = new TCanvas("cc", "cc", 0, 0, 1600, 1200);
0505 cc->SetGridy();
0506 for (int j = 0; j < NVAR; ++j) {
0507 gSystem->mkdir(TString(plotsDir) + "/_" + titleVar[j]);
0508 for (int d = 0; d < NDET; ++d) {
0509 cc->Clear();
0510 r[j][d]->Draw("colz");
0511 cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(r[j][d]->GetName()) + ".pdf");
0512
0513 cc->Clear();
0514 p[j][d]->Draw();
0515 if (doRatio[j]) {
0516 p[j][d]->SetMinimum(Pmin);
0517 p[j][d]->SetMaximum(Pmax);
0518 }
0519 cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(p[j][d]->GetName()) + ".pdf");
0520
0521 for (int i = 0; i < NLEV; ++i) {
0522 for (int p = 0; p < DEP[d]; ++p) {
0523 cc->Clear();
0524 h[j][i][d][p]->Draw("colz");
0525
0526 if (i == NLEV - 1) {
0527 if (doRatio[j]) {
0528 h[j][2][d][p]->SetMinimum(Zmin);
0529 h[j][2][d][p]->SetMaximum(Zmax);
0530 } else {
0531 h[j][2][d][p]->SetMinimum(-3);
0532 h[j][2][d][p]->SetMaximum(3);
0533 }
0534 }
0535 cc->Print(TString(plotsDir) + "/_" + titleVar[j] + "/" + TString(h[j][i][d][p]->GetName()) + ".pdf");
0536 }
0537 }
0538 }
0539 }
0540
0541 for (int i = 0; i < NLEV; ++i) {
0542 cc->Clear();
0543 hlut[0][i]->Draw();
0544 hlut[1][i]->Draw("sames");
0545 hlut[2][i]->Draw("sames");
0546 hlut[3][i]->Draw("sames");
0547 cc->Print(TString(plotsDir) + Form("LUT_%d.gif", i));
0548 }
0549 cc->Clear();
0550 hslope[0]->Draw("colz");
0551 cc->SetGridx();
0552 cc->SetGridy();
0553 cc->Print(TString(plotsDir) + "GainLutScatterHBHE.pdf");
0554 cc->Clear();
0555 hslope[1]->Draw("colz");
0556 cc->SetGridx();
0557 cc->SetGridy();
0558 cc->Print(TString(plotsDir) + "GainLutScatterLutHF.pdf");
0559
0560 for (int i = 0; i < 2; ++i) {
0561 cc->Clear();
0562 houtput[0][i]->Draw("box");
0563 houtput[1][i]->Draw("samebox");
0564 houtput[2][i]->Draw("samebox");
0565 houtput[3][i]->Draw("samebox");
0566 cc->Print(TString(plotsDir) + Form("OUT_%d.gif", i));
0567 }
0568 }
0569
0570 void HcalLutAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0571 edm::ParameterSetDescription desc;
0572 desc.setUnknown();
0573 descriptions.addDefault(desc);
0574 }
0575
0576 DEFINE_FWK_MODULE(HcalLutAnalyzer);