Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:23

0001 #include <memory>
0002 #include <fstream>
0003 #include <iostream>
0004 
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/ParameterSet/interface/FileInPath.h"
0011 
0012 #include "DataFormats/DetId/interface/DetId.h"
0013 #include "DataFormats/EcalDetId/interface/ESDetId.h"
0014 #include "DataFormats/EcalDigi/interface/ESDataFrame.h"
0015 #include "DataFormats/EcalRawData/interface/EcalRawDataCollections.h"
0016 
0017 #include "DQMServices/Core/interface/DQMStore.h"
0018 
0019 #include "DQM/EcalPreshowerMonitorModule/interface/ESPedestalTask.h"
0020 
0021 using namespace cms;
0022 using namespace edm;
0023 using namespace std;
0024 
0025 ESPedestalTask::ESPedestalTask(const edm::ParameterSet& ps) {
0026   digitoken_ = consumes<ESDigiCollection>(ps.getParameter<InputTag>("DigiLabel"));
0027   lookup_ = ps.getUntrackedParameter<FileInPath>("LookupTable");
0028   outputFile_ = ps.getUntrackedParameter<string>("OutputFile", "");
0029   prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
0030 
0031   for (int i = 0; i < 2; ++i)
0032     for (int j = 0; j < 2; ++j)
0033       for (int k = 0; k < 40; ++k)
0034         for (int l = 0; l < 40; ++l)
0035           senCount_[i][j][k][l] = -1;
0036 
0037   ievt_ = 0;
0038 }
0039 
0040 void ESPedestalTask::bookHistograms(DQMStore::IBooker& iBooker, Run const&, EventSetup const&) {
0041   int iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
0042   int senZ_[4288], senP_[4288], senX_[4288], senY_[4288];
0043 
0044   // read in look-up table
0045   ifstream file(lookup_.fullPath().c_str());
0046   if (!file.is_open())
0047     throw cms::Exception("FileNotFound") << lookup_.fullPath();
0048 
0049   file >> nLines_;
0050 
0051   for (int i = 0; i < nLines_; ++i) {
0052     file >> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
0053 
0054     senZ_[i] = iz;
0055     senP_[i] = ip;
0056     senX_[i] = ix;
0057     senY_[i] = iy;
0058 
0059     iz = (senZ_[i] == 1) ? 0 : 1;
0060     senCount_[iz][senP_[i] - 1][senX_[i] - 1][senY_[i] - 1] = i;
0061   }
0062 
0063   char hname[300];
0064 
0065   iBooker.setCurrentFolder(prefixME_ + "/ESPedestalTask");
0066 
0067   for (int i = 0; i < nLines_; ++i) {
0068     for (int is = 0; is < 32; ++is) {
0069       sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is + 1);
0070       meADC_[i][is] = iBooker.book1D(hname, hname, 1000, 899.5, 1899.5);
0071     }
0072   }
0073 }
0074 
0075 void ESPedestalTask::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
0076   ievt_++;
0077   runNum_ = e.id().run();
0078 
0079   Handle<ESDigiCollection> digis;
0080   e.getByToken(digitoken_, digis);
0081 
0082   runtype_ = 1;  // Let runtype_ = 1
0083 
0084   // Digis
0085   int zside, plane, ix, iy, strip, iz;
0086   if (digis.isValid()) {
0087     for (ESDigiCollection::const_iterator digiItr = digis->begin(); digiItr != digis->end(); ++digiItr) {
0088       ESDataFrame dataframe = (*digiItr);
0089       ESDetId id = dataframe.id();
0090 
0091       zside = id.zside();
0092       plane = id.plane();
0093       ix = id.six();
0094       iy = id.siy();
0095       strip = id.strip();
0096       iz = (zside == 1) ? 0 : 1;
0097 
0098       if (meADC_[senCount_[iz][plane - 1][ix - 1][iy - 1]][strip - 1]) {
0099         if (runtype_ == 1) {
0100           meADC_[senCount_[iz][plane - 1][ix - 1][iy - 1]][strip - 1]->Fill(dataframe.sample(0).adc());
0101           meADC_[senCount_[iz][plane - 1][ix - 1][iy - 1]][strip - 1]->Fill(dataframe.sample(1).adc());
0102           meADC_[senCount_[iz][plane - 1][ix - 1][iy - 1]][strip - 1]->Fill(dataframe.sample(2).adc());
0103         } else if (runtype_ == 3) {
0104           meADC_[senCount_[iz][plane - 1][ix - 1][iy - 1]][strip - 1]->Fill(dataframe.sample(1).adc());
0105         }
0106       }
0107     }
0108   }
0109 }
0110 
0111 DEFINE_FWK_MODULE(ESPedestalTask);