Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:00:06

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