Go to the documentation of this file.00001 #include <memory>
00002 #include <fstream>
00003 #include <iostream>
00004
00005 #include "FWCore/ServiceRegistry/interface/Service.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/ParameterSet/interface/FileInPath.h"
00011
00012 #include "DataFormats/DetId/interface/DetId.h"
00013 #include "DataFormats/EcalDetId/interface/ESDetId.h"
00014 #include "DataFormats/EcalDigi/interface/ESDataFrame.h"
00015 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00016 #include "DataFormats/EcalRawData/interface/EcalRawDataCollections.h"
00017
00018 #include "DQMServices/Core/interface/MonitorElement.h"
00019 #include "DQMServices/Core/interface/DQMStore.h"
00020
00021 #include "DQM/EcalPreshowerMonitorModule/interface/ESPedestalTask.h"
00022
00023 using namespace cms;
00024 using namespace edm;
00025 using namespace std;
00026
00027 ESPedestalTask::ESPedestalTask(const edm::ParameterSet& ps) {
00028
00029 init_ = false;
00030
00031 dqmStore_ = Service<DQMStore>().operator->();
00032
00033 digilabel_ = ps.getParameter<InputTag>("DigiLabel");
00034 lookup_ = ps.getUntrackedParameter<FileInPath>("LookupTable");
00035 outputFile_ = ps.getUntrackedParameter<string>("OutputFile","");
00036 prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
00037
00038 for (int i=0; i<2; ++i)
00039 for (int j=0; j<2; ++j)
00040 for (int k=0; k<40; ++k)
00041 for (int l=0; l<40; ++l)
00042 senCount_[i][j][k][l] = -1;
00043
00044 }
00045
00046 ESPedestalTask::~ESPedestalTask() {
00047 }
00048
00049 void ESPedestalTask::beginJob(void) {
00050 ievt_ = 0;
00051 }
00052
00053 void ESPedestalTask::beginRun(const Run& r, const EventSetup& c) {
00054
00055 if ( ! mergeRuns_ ) this->reset();
00056
00057 }
00058
00059 void ESPedestalTask::endRun(const Run& r, const EventSetup& c) {
00060 }
00061
00062 void ESPedestalTask::reset(void) {
00063
00064 }
00065
00066 void ESPedestalTask::setup(void) {
00067
00068 init_ = true;
00069
00070 int iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
00071 int senZ_[4288], senP_[4288], senX_[4288], senY_[4288];
00072
00073
00074 ifstream file;
00075 file.open(lookup_.fullPath().c_str());
00076 if( file.is_open() ) {
00077
00078 file >> nLines_;
00079
00080 for (int i=0; i<nLines_; ++i) {
00081 file>> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
00082
00083 senZ_[i] = iz;
00084 senP_[i] = ip;
00085 senX_[i] = ix;
00086 senY_[i] = iy;
00087
00088 iz = (senZ_[i]==1) ? 0:1;
00089 senCount_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1] = i;
00090 }
00091 }
00092 else {
00093 cout<<"ESPedestalTask : Look up table file can not be found in "<<lookup_.fullPath().c_str()<<endl;
00094 }
00095
00096 char hname[300];
00097
00098 if (dqmStore_) {
00099 dqmStore_->setCurrentFolder(prefixME_ + "/ESPedestalTask");
00100
00101 for (int i=0; i<nLines_; ++i) {
00102 for (int is=0; is<32; ++is) {
00103 sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is+1);
00104 meADC_[i][is] = dqmStore_->book1D(hname, hname, 1000, 899.5, 1899.5);
00105 }
00106 }
00107 }
00108
00109 }
00110
00111 void ESPedestalTask::cleanup(void){
00112
00113 if ( ! init_ ) return;
00114
00115 init_ = false;
00116 }
00117
00118 void ESPedestalTask::endJob(void) {
00119
00120 LogInfo("ESPedestalTask") << "analyzed " << ievt_ << " events";
00121
00122 if ( enableCleanup_ ) this->cleanup();
00123
00124 }
00125
00126 void ESPedestalTask::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
00127
00128 if ( ! init_ ) this->setup();
00129
00130 ievt_++;
00131 runNum_ = e.id().run();
00132
00133 Handle<ESDigiCollection> digis;
00134 e.getByLabel(digilabel_, digis);
00135
00136 runtype_ = 1;
00137
00138
00139 int zside, plane, ix, iy, strip, iz;
00140 if (digis.isValid()) {
00141 for (ESDigiCollection::const_iterator digiItr = digis->begin(); digiItr != digis->end(); ++digiItr) {
00142
00143 ESDataFrame dataframe = (*digiItr);
00144 ESDetId id = dataframe.id();
00145
00146 zside = id.zside();
00147 plane = id.plane();
00148 ix = id.six();
00149 iy = id.siy();
00150 strip = id.strip();
00151 iz = (zside==1) ? 0:1;
00152
00153 if (meADC_[senCount_[iz][plane-1][ix-1][iy-1]][strip-1]) {
00154 if(runtype_ == 1){
00155 meADC_[senCount_[iz][plane-1][ix-1][iy-1]][strip-1]->Fill(dataframe.sample(0).adc());
00156 meADC_[senCount_[iz][plane-1][ix-1][iy-1]][strip-1]->Fill(dataframe.sample(1).adc());
00157 meADC_[senCount_[iz][plane-1][ix-1][iy-1]][strip-1]->Fill(dataframe.sample(2).adc());
00158 } else if(runtype_ == 3) {
00159 meADC_[senCount_[iz][plane-1][ix-1][iy-1]][strip-1]->Fill(dataframe.sample(1).adc());
00160 }
00161 }
00162
00163 }
00164 }
00165
00166 }
00167
00168 DEFINE_FWK_MODULE(ESPedestalTask);