Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:16

0001 /*
0002  * \file DQMHcalPhiSymAlCaReco.cc
0003  *
0004  * \author Olga Kodolova
0005  *
0006  *
0007  *
0008  * Description: Monitoring of Phi Symmetry Calibration Stream
0009  */
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 
0019 // DQM include files
0020 
0021 #include "DQMServices/Core/interface/DQMStore.h"
0022 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0023 
0024 // work on collections
0025 
0026 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0027 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0028 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0029 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0030 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0031 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0032 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0033 
0034 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
0035 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0036 
0037 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
0038 #include "EventFilter/HcalRawToDigi/interface/HcalHTRData.h"
0039 
0040 class DQMHcalPhiSymAlCaReco : public DQMOneEDAnalyzer<> {
0041 public:
0042   DQMHcalPhiSymAlCaReco(const edm::ParameterSet &);
0043   ~DQMHcalPhiSymAlCaReco() override;
0044 
0045 protected:
0046   //  void beginRun(const edm::Run& r, const edm::EventSetup& c);
0047   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0048   void analyze(const edm::Event &e, const edm::EventSetup &c) override;
0049 
0050   void dqmEndRun(const edm::Run &r, const edm::EventSetup &c) override;
0051 
0052 private:
0053   int eventCounter_;
0054 
0055   //
0056   // Monitor elements
0057   //
0058   MonitorElement *hiDistrMBPl2D_;
0059   MonitorElement *hiDistrNoisePl2D_;
0060   MonitorElement *hiDistrMBMin2D_;
0061   MonitorElement *hiDistrNoiseMin2D_;
0062 
0063   MonitorElement *hiDistrMB2Pl2D_;
0064   MonitorElement *hiDistrNoise2Pl2D_;
0065   MonitorElement *hiDistrMB2Min2D_;
0066   MonitorElement *hiDistrNoise2Min2D_;
0067 
0068   MonitorElement *hiDistrVarMBPl2D_;
0069   MonitorElement *hiDistrVarNoisePl2D_;
0070   MonitorElement *hiDistrVarMBMin2D_;
0071   MonitorElement *hiDistrVarNoiseMin2D_;
0072 
0073   MonitorElement *hiDistrHBHEsize1D_;
0074   MonitorElement *hiDistrHFsize1D_;
0075 
0076   MonitorElement *hFEDsize;
0077   MonitorElement *hHcalIsZS;
0078   MonitorElement *hL1Id;
0079 
0080   int hiDistr_y_nbin_;
0081   int hiDistr_x_nbin_;
0082   double hiDistr_y_min_;
0083   double hiDistr_y_max_;
0084   double hiDistr_x_min_;
0085   double hiDistr_x_max_;
0086 
0087   int hiDistr_r_nbin_;
0088   double ihbhe_size_;
0089   double ihf_size_;
0090 
0091   bool perLSsaving_;  //to avoid nanoDQMIO crashing, driven by  DQMServices/Core/python/DQMStore_cfi.py
0092 
0093   /// object to monitor
0094 
0095   edm::EDGetTokenT<HBHERecHitCollection> hbherecoMB;
0096   edm::InputTag horecoMB;
0097   edm::EDGetTokenT<HFRecHitCollection> hfrecoMB;
0098 
0099   edm::EDGetTokenT<HBHERecHitCollection> hbherecoNoise;
0100   edm::InputTag horecoNoise;
0101   edm::EDGetTokenT<HFRecHitCollection> hfrecoNoise;
0102 
0103   edm::EDGetTokenT<FEDRawDataCollection> rawInLabel_;
0104 
0105   /// DQM folder name
0106   std::string folderName_;
0107 
0108   /// Write to file
0109   bool saveToFile_;
0110 
0111   // period of ZS
0112   unsigned int period_;
0113 
0114   /// Output file name if required
0115   std::string fileName_;
0116 };
0117 
0118 // ******************************************
0119 // constructors
0120 // *****************************************
0121 
0122 DQMHcalPhiSymAlCaReco::DQMHcalPhiSymAlCaReco(const edm::ParameterSet &ps) : eventCounter_(0) {
0123   //
0124   // Input from configurator file
0125   //
0126   folderName_ = ps.getUntrackedParameter<std::string>("FolderName", "ALCAStreamHcalPhiSym");
0127 
0128   hbherecoMB = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInputMB"));
0129   horecoMB = ps.getParameter<edm::InputTag>("hoInputMB");
0130   hfrecoMB = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInputMB"));
0131 
0132   hbherecoNoise = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInputNoise"));
0133   horecoNoise = ps.getParameter<edm::InputTag>("hoInputNoise");
0134   hfrecoNoise = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInputNoise"));
0135 
0136   rawInLabel_ = consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("rawInputLabel"));
0137 
0138   period_ = ps.getParameter<unsigned int>("period");
0139 
0140   saveToFile_ = ps.getUntrackedParameter<bool>("SaveToFile", false);
0141   fileName_ = ps.getUntrackedParameter<std::string>("FileName", "MonitorAlCaHcalPhiSym.root");
0142 
0143   perLSsaving_ = (ps.getUntrackedParameter<bool>("perLSsaving", false));
0144 
0145   // histogram parameters
0146 
0147   // Distribution of rechits in iPhi, iEta
0148   hiDistr_y_nbin_ = ps.getUntrackedParameter<int>("hiDistr_y_nbin", 72);
0149   hiDistr_y_min_ = ps.getUntrackedParameter<double>("hiDistr_y_min", 0.5);
0150   hiDistr_y_max_ = ps.getUntrackedParameter<double>("hiDistr_y_max", 72.5);
0151   hiDistr_x_nbin_ = ps.getUntrackedParameter<int>("hiDistr_x_nbin", 41);
0152   hiDistr_x_min_ = ps.getUntrackedParameter<double>("hiDistr_x_min", 0.5);
0153   hiDistr_x_max_ = ps.getUntrackedParameter<double>("hiDistr_x_max", 41.5);
0154   // Check for NZS
0155   hiDistr_r_nbin_ = ps.getUntrackedParameter<int>("hiDistr_r_nbin", 100);
0156   ihbhe_size_ = ps.getUntrackedParameter<double>("ihbhe_size_", 5184.);
0157   ihf_size_ = ps.getUntrackedParameter<double>("ihf_size_", 1728.);
0158 }
0159 
0160 DQMHcalPhiSymAlCaReco::~DQMHcalPhiSymAlCaReco() {}
0161 
0162 //--------------------------------------------------------
0163 void DQMHcalPhiSymAlCaReco::bookHistograms(DQMStore::IBooker &ibooker,
0164                                            edm::Run const &irun,
0165                                            edm::EventSetup const &isetup) {
0166   // create and cd into new folder
0167   ibooker.setCurrentFolder(folderName_);
0168 
0169   eventCounter_ = 0;
0170 
0171   hFEDsize = ibooker.book1D("hFEDsize", "HCAL FED size (kB)", 200, -0.5, 20.5);
0172   hFEDsize->setAxisTitle("kB", 1);
0173 
0174   hHcalIsZS = ibooker.book1D("hHcalIsZS", "Hcal Is ZS", 4, -1.5, 2.5);
0175   hHcalIsZS->setBinLabel(2, "NZS");
0176   hHcalIsZS->setBinLabel(3, "ZS");
0177 
0178   char hname[50];
0179   sprintf(hname, "L1 Event Number %% %i", period_);
0180   hL1Id = ibooker.book1D("hL1Id", hname, 4200, -99.5, 4099.5);
0181   hL1Id->setAxisTitle(hname);
0182 
0183   // book some histograms 1D
0184   double xmin = 0.1;
0185   double xmax = 1.1;
0186   hiDistrHBHEsize1D_ = ibooker.book1D("DistrHBHEsize", "Size of HBHE Collection", hiDistr_r_nbin_, xmin, xmax);
0187   hiDistrHFsize1D_ = ibooker.book1D("DistrHFsize", "Size of HF Collection", hiDistr_r_nbin_, xmin, xmax);
0188 
0189   // First moment
0190   hiDistrMBPl2D_ = ibooker.book2D("MBdepthPl1",
0191                                   "iphi- +ieta signal distribution at depth1",
0192                                   hiDistr_x_nbin_,
0193                                   hiDistr_x_min_,
0194                                   hiDistr_x_max_,
0195                                   hiDistr_y_nbin_,
0196                                   hiDistr_y_min_,
0197                                   hiDistr_y_max_);
0198 
0199   hiDistrMBPl2D_->setAxisTitle("i#phi ", 2);
0200   hiDistrMBPl2D_->setAxisTitle("i#eta ", 1);
0201 
0202   hiDistrNoisePl2D_ = ibooker.book2D("NoisedepthPl1",
0203                                      "iphi-ieta noise distribution at depth1",
0204                                      hiDistr_x_nbin_ + 1,
0205                                      hiDistr_x_min_ - 1.,
0206                                      hiDistr_x_max_,
0207                                      hiDistr_y_nbin_ + 1,
0208                                      hiDistr_y_min_ - 1.,
0209                                      hiDistr_y_max_);
0210 
0211   hiDistrNoisePl2D_->setAxisTitle("i#phi ", 2);
0212   hiDistrNoisePl2D_->setAxisTitle("i#eta ", 1);
0213   // Second moment
0214   hiDistrMB2Pl2D_ = ibooker.book2D("MB2depthPl1",
0215                                    "iphi- +ieta signal distribution at depth1",
0216                                    hiDistr_x_nbin_,
0217                                    hiDistr_x_min_,
0218                                    hiDistr_x_max_,
0219                                    hiDistr_y_nbin_,
0220                                    hiDistr_y_min_,
0221                                    hiDistr_y_max_);
0222 
0223   hiDistrMB2Pl2D_->setAxisTitle("i#phi ", 2);
0224   hiDistrMB2Pl2D_->setAxisTitle("i#eta ", 1);
0225 
0226   hiDistrNoise2Pl2D_ = ibooker.book2D("Noise2depthPl1",
0227                                       "iphi-ieta noise distribution at depth1",
0228                                       hiDistr_x_nbin_,
0229                                       hiDistr_x_min_,
0230                                       hiDistr_x_max_,
0231                                       hiDistr_y_nbin_,
0232                                       hiDistr_y_min_,
0233                                       hiDistr_y_max_);
0234 
0235   hiDistrNoise2Pl2D_->setAxisTitle("i#phi ", 2);
0236   hiDistrNoise2Pl2D_->setAxisTitle("i#eta ", 1);
0237 
0238   // Variance
0239   hiDistrVarMBPl2D_ = ibooker.book2D("VarMBdepthPl1",
0240                                      "iphi- +ieta signal distribution at depth1",
0241                                      hiDistr_x_nbin_,
0242                                      hiDistr_x_min_,
0243                                      hiDistr_x_max_,
0244                                      hiDistr_y_nbin_,
0245                                      hiDistr_y_min_,
0246                                      hiDistr_y_max_);
0247 
0248   hiDistrVarMBPl2D_->setAxisTitle("i#phi ", 2);
0249   hiDistrVarMBPl2D_->setAxisTitle("i#eta ", 1);
0250 
0251   hiDistrVarNoisePl2D_ = ibooker.book2D("VarNoisedepthPl1",
0252                                         "iphi-ieta noise distribution at depth1",
0253                                         hiDistr_x_nbin_,
0254                                         hiDistr_x_min_,
0255                                         hiDistr_x_max_,
0256                                         hiDistr_y_nbin_,
0257                                         hiDistr_y_min_,
0258                                         hiDistr_y_max_);
0259 
0260   hiDistrVarNoisePl2D_->setAxisTitle("i#phi ", 2);
0261   hiDistrVarNoisePl2D_->setAxisTitle("i#eta ", 1);
0262 
0263   //==================================================================================
0264   // First moment
0265   hiDistrMBMin2D_ = ibooker.book2D("MBdepthMin1",
0266                                    "iphi- +ieta signal distribution at depth1",
0267                                    hiDistr_x_nbin_,
0268                                    hiDistr_x_min_,
0269                                    hiDistr_x_max_,
0270                                    hiDistr_y_nbin_,
0271                                    hiDistr_y_min_,
0272                                    hiDistr_y_max_);
0273 
0274   hiDistrMBMin2D_->setAxisTitle("i#phi ", 2);
0275   hiDistrMBMin2D_->setAxisTitle("i#eta ", 1);
0276 
0277   hiDistrNoiseMin2D_ = ibooker.book2D("NoisedepthMin1",
0278                                       "iphi-ieta noise distribution at depth1",
0279                                       hiDistr_x_nbin_,
0280                                       hiDistr_x_min_,
0281                                       hiDistr_x_max_,
0282                                       hiDistr_y_nbin_,
0283                                       hiDistr_y_min_,
0284                                       hiDistr_y_max_);
0285 
0286   hiDistrNoiseMin2D_->setAxisTitle("i#phi ", 2);
0287   hiDistrNoiseMin2D_->setAxisTitle("i#eta ", 1);
0288   // Second moment
0289   hiDistrMB2Min2D_ = ibooker.book2D("MB2depthMin1",
0290                                     "iphi- +ieta signal distribution at depth1",
0291                                     hiDistr_x_nbin_,
0292                                     hiDistr_x_min_,
0293                                     hiDistr_x_max_,
0294                                     hiDistr_y_nbin_,
0295                                     hiDistr_y_min_,
0296                                     hiDistr_y_max_);
0297 
0298   hiDistrMB2Min2D_->setAxisTitle("i#phi ", 2);
0299   hiDistrMB2Min2D_->setAxisTitle("i#eta ", 1);
0300 
0301   hiDistrNoise2Min2D_ = ibooker.book2D("Noise2depthMin1",
0302                                        "iphi-ieta noise distribution at depth1",
0303                                        hiDistr_x_nbin_,
0304                                        hiDistr_x_min_,
0305                                        hiDistr_x_max_,
0306                                        hiDistr_y_nbin_,
0307                                        hiDistr_y_min_,
0308                                        hiDistr_y_max_);
0309 
0310   hiDistrNoise2Min2D_->setAxisTitle("i#phi ", 2);
0311   hiDistrNoise2Min2D_->setAxisTitle("i#eta ", 1);
0312 
0313   // Variance
0314   hiDistrVarMBMin2D_ = ibooker.book2D("VarMBdepthMin1",
0315                                       "iphi- +ieta signal distribution at depth1",
0316                                       hiDistr_x_nbin_,
0317                                       hiDistr_x_min_,
0318                                       hiDistr_x_max_,
0319                                       hiDistr_y_nbin_,
0320                                       hiDistr_y_min_,
0321                                       hiDistr_y_max_);
0322 
0323   hiDistrVarMBMin2D_->setAxisTitle("i#phi ", 2);
0324   hiDistrVarMBMin2D_->setAxisTitle("i#eta ", 1);
0325 
0326   hiDistrVarNoiseMin2D_ = ibooker.book2D("VarNoisedepthMin1",
0327                                          "iphi-ieta noise distribution at depth1",
0328                                          hiDistr_x_nbin_,
0329                                          hiDistr_x_min_,
0330                                          hiDistr_x_max_,
0331                                          hiDistr_y_nbin_,
0332                                          hiDistr_y_min_,
0333                                          hiDistr_y_max_);
0334 
0335   hiDistrVarNoiseMin2D_->setAxisTitle("i#phi ", 2);
0336   hiDistrVarNoiseMin2D_->setAxisTitle("i#eta ", 1);
0337 }
0338 
0339 //--------------------------------------------------------
0340 // void DQMHcalPhiSymAlCaReco::beginRun(const edm::Run& r, const EventSetup&
0341 // context) {
0342 ////   eventCounter_ = 0;
0343 //}
0344 
0345 //--------------------------------------------------------
0346 
0347 //-------------------------------------------------------------
0348 
0349 void DQMHcalPhiSymAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0350   eventCounter_++;
0351 
0352   edm::Handle<FEDRawDataCollection> rawIn;
0353   iEvent.getByToken(rawInLabel_, rawIn);
0354 
0355   if (!rawIn.isValid()) {
0356     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
0357     return;
0358   }
0359 
0360   // get HCAL FEDs:
0361   std::vector<int> selFEDs;
0362   for (int i = FEDNumbering::MINHCALFEDID; i <= FEDNumbering::MAXHCALFEDID; i++) {
0363     selFEDs.push_back(i);
0364   }
0365 
0366   //  std::cout<<" Size of FED "<<selFEDs.size()<<std::endl;
0367 
0368   const FEDRawDataCollection *rdc = rawIn.product();
0369 
0370   bool hcalIsZS = false;
0371   int lvl1ID = 0;
0372   bool lvl1IDFound = false;
0373   for (unsigned int k = 0; k < selFEDs.size(); k++) {
0374     const FEDRawData &fedData = rdc->FEDData(selFEDs[k]);
0375     // std::cout<<fedData.size()*std::pow(1024.,-1)<<std::endl;
0376     hFEDsize->Fill(fedData.size() * std::pow(1024., -1), 1);
0377 
0378     // get HCAL DCC Header for each FEDRawData
0379     const HcalDCCHeader *dccHeader = (const HcalDCCHeader *)(fedData.data());
0380     if (dccHeader) {
0381       // walk through the HTR data...
0382       HcalHTRData htr;
0383 
0384       int nspigot = 0;
0385       for (int spigot = 0; spigot < HcalDCCHeader::SPIGOT_COUNT; spigot++) {
0386         nspigot++;
0387 
0388         if (!dccHeader->getSpigotPresent(spigot))
0389           continue;
0390 
0391         // Load the given decoder with the pointer and length from this spigot.
0392         dccHeader->getSpigotData(spigot, htr, fedData.size());
0393 
0394         if (k != 20 && nspigot != 14) {
0395           if (!htr.isUnsuppressed()) {
0396             hcalIsZS = true;
0397           }
0398         }
0399       }
0400     }
0401 
0402     // try to get the lvl1ID from the HCAL fed
0403     if (!lvl1IDFound && (fedData.size() > 0)) {
0404       // get FED Header for FEDRawData
0405       FEDHeader fedHeader(fedData.data());
0406       lvl1ID = fedHeader.lvl1ID();
0407       lvl1IDFound = true;
0408     }
0409   }  // loop over HcalFEDs
0410 
0411   hHcalIsZS->Fill(hcalIsZS);
0412   hL1Id->Fill(lvl1ID % period_);
0413 
0414   edm::Handle<HBHERecHitCollection> hbheNS;
0415   iEvent.getByToken(hbherecoNoise, hbheNS);
0416 
0417   if (!hbheNS.isValid()) {
0418     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
0419     return;
0420   }
0421 
0422   edm::Handle<HBHERecHitCollection> hbheMB;
0423   iEvent.getByToken(hbherecoMB, hbheMB);
0424 
0425   if (!hbheMB.isValid()) {
0426     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
0427     return;
0428   }
0429 
0430   edm::Handle<HFRecHitCollection> hfNS;
0431   iEvent.getByToken(hfrecoNoise, hfNS);
0432 
0433   if (!hfNS.isValid()) {
0434     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
0435     return;
0436   }
0437 
0438   edm::Handle<HFRecHitCollection> hfMB;
0439   iEvent.getByToken(hfrecoMB, hfMB);
0440 
0441   if (!hfMB.isValid()) {
0442     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
0443     return;
0444   }
0445 
0446   const HBHERecHitCollection HithbheNS = *(hbheNS.product());
0447 
0448   hiDistrHBHEsize1D_->Fill(HithbheNS.size() / ihbhe_size_);
0449 
0450   for (HBHERecHitCollection::const_iterator hbheItr = HithbheNS.begin(); hbheItr != HithbheNS.end(); hbheItr++) {
0451     DetId id = (*hbheItr).detid();
0452     HcalDetId hid = HcalDetId(id);
0453 
0454     if (hid.depth() == 1) {
0455       if (hid.ieta() > 0) {
0456         hiDistrNoisePl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
0457         hiDistrNoise2Pl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0458       } else {
0459         hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy());
0460         hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0461       }
0462     }
0463   }
0464 
0465   const HBHERecHitCollection HithbheMB = *(hbheMB.product());
0466 
0467   for (HBHERecHitCollection::const_iterator hbheItr = HithbheMB.begin(); hbheItr != HithbheMB.end(); hbheItr++) {
0468     DetId id = (*hbheItr).detid();
0469     HcalDetId hid = HcalDetId(id);
0470 
0471     if (hid.depth() == 1) {
0472       if (hid.ieta() > 0) {
0473         hiDistrMBPl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
0474         hiDistrMB2Pl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0475       } else {
0476         hiDistrMBMin2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy());
0477         hiDistrMB2Min2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0478       }
0479     }
0480   }
0481 
0482   const HFRecHitCollection HithfNS = *(hfNS.product());
0483 
0484   hiDistrHFsize1D_->Fill(HithfNS.size() / ihf_size_);
0485 
0486   for (HFRecHitCollection::const_iterator hbheItr = HithfNS.begin(); hbheItr != HithfNS.end(); hbheItr++) {
0487     DetId id = (*hbheItr).detid();
0488     HcalDetId hid = HcalDetId(id);
0489 
0490     if (hid.depth() == 1) {
0491       if (hid.ieta() > 0) {
0492         hiDistrNoisePl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
0493         hiDistrNoise2Pl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0494       } else {
0495         hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy());
0496         hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0497       }
0498     }
0499   }
0500 
0501   const HFRecHitCollection HithfMB = *(hfMB.product());
0502 
0503   for (HFRecHitCollection::const_iterator hbheItr = HithfMB.begin(); hbheItr != HithfMB.end(); hbheItr++) {
0504     DetId id = (*hbheItr).detid();
0505     HcalDetId hid = HcalDetId(id);
0506 
0507     if (hid.depth() == 1) {
0508       if (hid.ieta() > 0) {
0509         hiDistrMBPl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
0510         hiDistrMB2Pl2D_->Fill(hid.ieta(), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0511       } else {
0512         hiDistrMBMin2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy());
0513         hiDistrMB2Min2D_->Fill(fabs(hid.ieta()), hid.iphi(), hbheItr->energy() * hbheItr->energy());
0514       }
0515     }
0516   }
0517 
0518 }  // analyze
0519 
0520 //--------------------------------------------------------
0521 //--------------------------------------------------------
0522 void DQMHcalPhiSymAlCaReco::dqmEndRun(const edm::Run &r, const edm::EventSetup &context) {
0523   // Keep Variances
0524   if (eventCounter_ > 0 && !perLSsaving_) {
0525     for (int k = 0; k <= hiDistr_x_nbin_; k++) {
0526       for (int j = 0; j <= hiDistr_y_nbin_; j++) {
0527         // First moment
0528         float cc1 = hiDistrMBPl2D_->getBinContent(k, j);
0529         cc1 = cc1 * 1. / eventCounter_;
0530         float cc2 = hiDistrNoisePl2D_->getBinContent(k, j);
0531         cc2 = cc2 * 1. / eventCounter_;
0532         float cc3 = hiDistrMBMin2D_->getBinContent(k, j);
0533         cc3 = cc3 * 1. / eventCounter_;
0534         float cc4 = hiDistrNoiseMin2D_->getBinContent(k, j);
0535         cc4 = cc4 * 1. / eventCounter_;
0536         // Second moment
0537         float cc11 = hiDistrMB2Pl2D_->getBinContent(k, j);
0538         cc11 = cc11 * 1. / eventCounter_;
0539         hiDistrVarMBPl2D_->setBinContent(k, j, cc11 - cc1 * cc1);
0540         float cc22 = hiDistrNoise2Pl2D_->getBinContent(k, j);
0541         cc22 = cc22 * 1. / eventCounter_;
0542         hiDistrVarNoisePl2D_->setBinContent(k, j, cc22 - cc2 * cc2);
0543         float cc33 = hiDistrMB2Min2D_->getBinContent(k, j);
0544         cc33 = cc33 * 1. / eventCounter_;
0545         hiDistrVarMBMin2D_->setBinContent(k, j, cc33 - cc3 * cc3);
0546         float cc44 = hiDistrNoise2Min2D_->getBinContent(k, j);
0547         cc44 = cc44 * 1. / eventCounter_;
0548         hiDistrVarNoiseMin2D_->setBinContent(k, j, cc44 - cc4 * cc4);
0549       }
0550     }
0551   }
0552 }
0553 
0554 DEFINE_FWK_MODULE(DQMHcalPhiSymAlCaReco);