Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <iostream>
00011
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016
00017 #include "CondFormats/EcalObjects/interface/EcalDCSTowerStatus.h"
00018 #include "CondFormats/DataRecord/interface/EcalDCSTowerStatusRcd.h"
00019
00020 #include "DQMServices/Core/interface/MonitorElement.h"
00021 #include "DQMServices/Core/interface/DQMStore.h"
00022
00023 #include "DQM/EcalCommon/interface/Numbers.h"
00024
00025 #include "DQM/EcalBarrelMonitorTasks/interface/EBDcsInfoTask.h"
00026
00027 EBDcsInfoTask::EBDcsInfoTask(const edm::ParameterSet& ps) {
00028
00029 dqmStore_ = edm::Service<DQMStore>().operator->();
00030
00031 prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");
00032
00033 enableCleanup_ = ps.getUntrackedParameter<bool>("enableCleanup", false);
00034
00035 mergeRuns_ = ps.getUntrackedParameter<bool>("mergeRuns", false);
00036
00037 meEBDcsFraction_ = 0;
00038 meEBDcsActiveMap_ = 0;
00039 for (int i = 0; i < 36; i++) {
00040 meEBDcsActive_[i] = 0;
00041 }
00042
00043 }
00044
00045 EBDcsInfoTask::~EBDcsInfoTask() {
00046
00047 }
00048
00049 void EBDcsInfoTask::beginJob(void){
00050
00051 std::string name;
00052
00053 if ( dqmStore_ ) {
00054
00055 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
00056
00057 name = "DCSSummary";
00058 meEBDcsFraction_ = dqmStore_->bookFloat(name);
00059 meEBDcsFraction_->Fill(0.0);
00060
00061 name = "DCSSummaryMap";
00062 meEBDcsActiveMap_ = dqmStore_->book2D(name, name, 72, 0., 72., 34, 0., 34.);
00063 meEBDcsActiveMap_->setAxisTitle("jphi", 1);
00064 meEBDcsActiveMap_->setAxisTitle("jeta", 2);
00065
00066 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/DCSContents");
00067
00068 for (int i = 0; i < 36; i++) {
00069 name = "EcalBarrel_" + Numbers::sEB(i+1);
00070 meEBDcsActive_[i] = dqmStore_->bookFloat(name);
00071 meEBDcsActive_[i]->Fill(-1.0);
00072 }
00073
00074 }
00075
00076 }
00077
00078 void EBDcsInfoTask::endJob(void) {
00079
00080 if ( enableCleanup_ ) this->cleanup();
00081
00082 }
00083
00084 void EBDcsInfoTask::beginLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const edm::EventSetup& iSetup){
00085
00086
00087 for ( int iett = 0; iett < 34; iett++ ) {
00088 for ( int iptt = 0; iptt < 72; iptt++ ) {
00089 readyLumi[iptt][iett] = 1;
00090 }
00091 }
00092
00093 if ( !iSetup.find( edm::eventsetup::EventSetupRecordKey::makeKey<EcalDCSTowerStatusRcd>() ) ) {
00094 edm::LogWarning("EBDcsInfoTask") << "EcalDCSTowerStatus record not found";
00095 return;
00096 }
00097
00098 edm::ESHandle<EcalDCSTowerStatus> pDCSStatus;
00099 iSetup.get<EcalDCSTowerStatusRcd>().get(pDCSStatus);
00100 if ( !pDCSStatus.isValid() ) {
00101 edm::LogWarning("EBDcsInfoTask") << "EcalDCSTowerStatus record not valid";
00102 return;
00103 }
00104 const EcalDCSTowerStatus* dcsStatus = pDCSStatus.product();
00105
00106 for(int iz=-1; iz<=1; iz+=2) {
00107 for(int iptt=0 ; iptt<72; iptt++) {
00108 for(int iett=0 ; iett<17; iett++) {
00109 if (EcalTrigTowerDetId::validDetId(iz,EcalBarrel,iett+1,iptt+1 )){
00110
00111 EcalTrigTowerDetId ebid(iz,EcalBarrel,iett+1,iptt+1);
00112
00113 uint16_t dbStatus = 0;
00114 EcalDCSTowerStatus::const_iterator dcsStatusIt = dcsStatus->find( ebid.rawId() );
00115 if ( dcsStatusIt != dcsStatus->end() ) dbStatus = dcsStatusIt->getStatusCode();
00116
00117 if ( dbStatus > 0 ) {
00118 int ipttEB = iptt;
00119 int iettEB = (iz==-1) ? iett : 17+iett;
00120 readyRun[ipttEB][iettEB] = 0;
00121 readyLumi[ipttEB][iettEB] = 0;
00122 }
00123
00124 }
00125 }
00126 }
00127 }
00128
00129
00130 }
00131
00132 void EBDcsInfoTask::endLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const edm::EventSetup& iSetup) {
00133
00134 this->fillMonitorElements(readyLumi);
00135
00136 }
00137
00138 void EBDcsInfoTask::beginRun(const edm::Run& r, const edm::EventSetup& c) {
00139
00140 if ( ! mergeRuns_ ) this->reset();
00141
00142 for ( int iett = 0; iett < 34; iett++ ) {
00143 for ( int iptt = 0; iptt < 72; iptt++ ) {
00144 readyRun[iptt][iett] = 1;
00145 }
00146 }
00147
00148 }
00149
00150 void EBDcsInfoTask::endRun(const edm::Run& r, const edm::EventSetup& c) {
00151
00152 this->fillMonitorElements(readyRun);
00153
00154 }
00155
00156 void EBDcsInfoTask::reset(void) {
00157
00158 if ( meEBDcsFraction_ ) meEBDcsFraction_->Reset();
00159
00160 for (int i = 0; i < 36; i++) {
00161 if ( meEBDcsActive_[i] ) meEBDcsActive_[i]->Reset();
00162 }
00163
00164 if ( meEBDcsActiveMap_ ) meEBDcsActiveMap_->Reset();
00165
00166 }
00167
00168
00169 void EBDcsInfoTask::cleanup(void){
00170
00171 if ( dqmStore_ ) {
00172
00173 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
00174
00175 if ( meEBDcsFraction_ ) dqmStore_->removeElement( meEBDcsFraction_->getName() );
00176
00177 if ( meEBDcsActiveMap_ ) dqmStore_->removeElement( meEBDcsActiveMap_->getName() );
00178
00179 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/DCSContents");
00180
00181 for (int i = 0; i < 36; i++) {
00182 if ( meEBDcsActive_[i] ) dqmStore_->removeElement( meEBDcsActive_[i]->getName() );
00183 }
00184
00185 }
00186
00187 }
00188
00189 void EBDcsInfoTask::fillMonitorElements(int ready[72][34]) {
00190
00191 float readySum[36];
00192 for ( int ism = 0; ism < 36; ism++ ) readySum[ism] = 0;
00193 float readySumTot = 0.;
00194
00195 for ( int iett = 0; iett < 34; iett++ ) {
00196 for ( int iptt = 0; iptt < 72; iptt++ ) {
00197
00198 if(meEBDcsActiveMap_) meEBDcsActiveMap_->setBinContent( iptt+1, iett+1, ready[iptt][iett] );
00199
00200 int ism = ( iett<17 ) ? iptt/4 : 18+iptt/4;
00201 if(ready[iptt][iett]) {
00202 readySum[ism]++;
00203 readySumTot++;
00204 }
00205
00206 }
00207 }
00208
00209 for ( int ism = 0; ism < 36; ism++ ) {
00210 if( meEBDcsActive_[ism] ) meEBDcsActive_[ism]->Fill( readySum[ism]/68. );
00211 }
00212
00213 if( meEBDcsFraction_ ) meEBDcsFraction_->Fill(readySumTot/34./72.);
00214
00215 }
00216
00217 void EBDcsInfoTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00218
00219 }