00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <iostream>
00011 #include <algorithm>
00012
00013 #include "FWCore/ServiceRegistry/interface/Service.h"
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/Framework/interface/EventSetup.h"
00017
00018 #include "DQMServices/Core/interface/MonitorElement.h"
00019 #include "DQMServices/Core/interface/DQMStore.h"
00020
00021 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00022
00023 #include "DQM/EcalCommon/interface/Numbers.h"
00024 #include "DQM/EcalCommon/interface/UtilsClient.h"
00025
00026 #include "DQM/EcalBarrelMonitorTasks/interface/EBDataCertificationTask.h"
00027
00028 EBDataCertificationTask::EBDataCertificationTask(const edm::ParameterSet& ps) {
00029
00030 dqmStore_ = edm::Service<DQMStore>().operator->();
00031
00032
00033 cloneME_ = ps.getUntrackedParameter<bool>("cloneME", true);
00034
00035 prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");
00036
00037 enableCleanup_ = ps.getUntrackedParameter<bool>("enableCleanup", false);
00038
00039 mergeRuns_ = ps.getUntrackedParameter<bool>("mergeRuns", false);
00040
00041 meEBDataCertificationSummary_ = 0;
00042 meEBDataCertificationSummaryMap_ = 0;
00043 for (int i = 0; i < 36; i++) {
00044 meEBDataCertification_[i] = 0;
00045 }
00046
00047 hDQM_ = 0;
00048 hDAQ_ = 0;
00049 hDCS_ = 0;
00050 hIntegrityByLumi_ = 0;
00051 hFrontendByLumi_ = 0;
00052 hSynchronizationByLumi_ = 0;
00053
00054 }
00055
00056 EBDataCertificationTask::~EBDataCertificationTask() {
00057
00058 }
00059
00060 void EBDataCertificationTask::beginJob(void){
00061
00062 char histo[200];
00063
00064 if ( dqmStore_ ) {
00065
00066 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
00067
00068 sprintf(histo, "CertificationSummary");
00069 meEBDataCertificationSummary_ = dqmStore_->bookFloat(histo);
00070 meEBDataCertificationSummary_->Fill(-1.0);
00071
00072 sprintf(histo, "CertificationSummaryMap");
00073 meEBDataCertificationSummaryMap_ = dqmStore_->book2D(histo,histo, 72, 0., 72., 34, 0., 34.);
00074 meEBDataCertificationSummaryMap_->setAxisTitle("jphi", 1);
00075 meEBDataCertificationSummaryMap_->setAxisTitle("jeta", 2);
00076
00077 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/CertificationContents");
00078
00079 for (int i = 0; i < 36; i++) {
00080 sprintf(histo, "EcalBarrel_%s", Numbers::sEB(i+1).c_str());
00081 meEBDataCertification_[i] = dqmStore_->bookFloat(histo);
00082 meEBDataCertification_[i]->Fill(-1.0);
00083 }
00084
00085 }
00086
00087 }
00088
00089 void EBDataCertificationTask::endJob(void) {
00090
00091 if ( enableCleanup_ ) this->cleanup();
00092
00093 }
00094
00095 void EBDataCertificationTask::beginLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const edm::EventSetup& iSetup){
00096
00097 }
00098
00099 void EBDataCertificationTask::endLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const edm::EventSetup& iSetup) {
00100
00101 this->reset();
00102
00103 char histo[200];
00104
00105 MonitorElement* me;
00106
00107
00108 float DQMVal[36];
00109 for (int i = 0; i < 36; i++) {
00110 DQMVal[i] = -1.;
00111 }
00112
00113 sprintf(histo, (prefixME_ + "/EBIntegrityTask/EBIT weighted integrity errors by lumi").c_str());
00114 me = dqmStore_->get(histo);
00115 hIntegrityByLumi_ = UtilsClient::getHisto<TH1F*>( me, cloneME_, hIntegrityByLumi_ );
00116
00117 sprintf(histo, (prefixME_ + "/EBStatusFlagsTask/FEStatus/EBSFT weighted frontend errors by lumi").c_str());
00118 me = dqmStore_->get(histo);
00119 hFrontendByLumi_ = UtilsClient::getHisto<TH1F*>( me, cloneME_, hFrontendByLumi_ );
00120
00121 sprintf(histo, (prefixME_ + "/EBRawDataTask/EBRDT FE synchronization errors by lumi").c_str());
00122 me = dqmStore_->get(histo);
00123 hSynchronizationByLumi_ = UtilsClient::getHisto<TH1F*>( me, cloneME_, hSynchronizationByLumi_ );
00124
00125 if( hIntegrityByLumi_ && hFrontendByLumi_ && hSynchronizationByLumi_ ) {
00126
00127 float integrityErrSum = 0.;
00128 float integrityQual = 1.0;
00129 float frontendErrSum = 0.;
00130 float frontendQual = 1.0;
00131 float synchronizationErrSum = 0.;
00132 float synchronizationQual = 1.0;
00133
00134 for ( int i=0; i<36; i++) {
00135 float ismIntegrityQual = 1.0;
00136 if( hIntegrityByLumi_->GetBinContent(0) > 0 ) {
00137 float errors = hIntegrityByLumi_->GetBinContent(i+1);
00138 ismIntegrityQual = 1.0 - errors/hIntegrityByLumi_->GetBinContent(0);
00139 integrityErrSum += errors;
00140 }
00141 float ismFrontendQual = 1.0;
00142 if( hFrontendByLumi_->GetBinContent(0) > 0 ) {
00143 float errors = hFrontendByLumi_->GetBinContent(i+1);
00144 ismFrontendQual = 1.0 - errors/hFrontendByLumi_->GetBinContent(0);
00145 frontendErrSum += errors;
00146 }
00147 float ismSynchronizationQual = 1.0;
00148 if( hSynchronizationByLumi_->GetBinContent(0) > 0 ) {
00149 float errors = hSynchronizationByLumi_->GetBinContent(i+1);
00150 ismSynchronizationQual = 1.0 - errors/hSynchronizationByLumi_->GetBinContent(0);
00151 synchronizationErrSum += errors;
00152 }
00153 float minVal= std::min(ismIntegrityQual,ismFrontendQual);
00154 DQMVal[i] = std::min(minVal,ismSynchronizationQual);
00155 }
00156
00157 if( hIntegrityByLumi_->GetBinContent(0) > 0 ) integrityQual = 1.0 - integrityErrSum/hIntegrityByLumi_->GetBinContent(0)/36.;
00158 if( hFrontendByLumi_->GetBinContent(0) > 0 ) frontendQual = 1.0 - frontendErrSum/hFrontendByLumi_->GetBinContent(0)/36.;
00159 if( hSynchronizationByLumi_->GetBinContent(0) > 0 ) synchronizationQual = 1.0 - synchronizationErrSum/hSynchronizationByLumi_->GetBinContent(0)/36.;
00160 float minVal = std::min(integrityQual,frontendQual);
00161 float totDQMVal = std::min(minVal,synchronizationQual);
00162
00163 sprintf(histo, (prefixME_ + "/EventInfo/reportSummary").c_str());
00164 me = dqmStore_->get(histo);
00165 if( me ) me->Fill(totDQMVal);
00166
00167 for ( int i=0; i<36; i++) {
00168 sprintf(histo, "EcalBarrel_%s", Numbers::sEB(i+1).c_str());
00169 me = dqmStore_->get(prefixME_ + "/EventInfo/reportSummaryContents/" + histo);
00170 if( me ) me->Fill(DQMVal[i]);
00171
00172 sprintf(histo, "reportSummaryMap");
00173 me = dqmStore_->get(prefixME_ + "/EventInfo/" + histo );
00174 if( me ) {
00175 for ( int iett = 0; iett < 34; iett++ ) {
00176 for ( int iptt = 0; iptt < 72; iptt++ ) {
00177 int ism = ( iett<17 ) ? iptt/4+1 : 18+iptt/4+1;
00178 if( i == (ism-1) ) me->setBinContent(iptt+1, iett+1, DQMVal[ism-1]);
00179 }
00180 }
00181 }
00182 }
00183
00184 }
00185
00186
00187 sprintf(histo, (prefixME_ + "/EventInfo/DAQSummaryMap").c_str());
00188 me = dqmStore_->get(histo);
00189 hDAQ_ = UtilsClient::getHisto<TH2F*>( me, cloneME_, hDAQ_ );
00190
00191 sprintf(histo, (prefixME_ + "/EventInfo/DCSSummaryMap").c_str());
00192 me = dqmStore_->get(histo);
00193 hDCS_ = UtilsClient::getHisto<TH2F*>( me, cloneME_, hDCS_ );
00194
00195 float sumCert = 0.;
00196 float sumCertEB[36];
00197 int nValidChannels = 0;
00198 int nValidChannelsEB[36];
00199
00200 for (int i = 0; i < 36; i++) {
00201 sumCertEB[i] = 0.;
00202 nValidChannelsEB[i] = 0;
00203 }
00204
00205 for ( int iett = 0; iett < 34; iett++ ) {
00206 for ( int iptt = 0; iptt < 72; iptt++ ) {
00207
00208 int ism = ( iett<17 ) ? iptt/4+1 : 18+iptt/4+1;
00209
00210 float xvalDQM = DQMVal[ism-1];
00211
00212 float xvalDAQ, xvalDCS;
00213 xvalDAQ = xvalDCS = -1.;
00214 float xcert = -1.;
00215
00216 if ( hDAQ_ ) xvalDAQ = hDAQ_->GetBinContent( iptt+1, iett+1 );
00217 if ( hDCS_ ) xvalDCS = hDCS_->GetBinContent( iptt+1, iett+1 );
00218
00219 if ( xvalDQM == -1 || ( xvalDAQ == -1 && xvalDCS == -1 ) ) {
00220
00221 xcert = 0.0;
00222 } else {
00223
00224 xcert = std::abs(xvalDQM) * std::abs(xvalDAQ) * std::abs(xvalDCS);
00225 }
00226
00227 if ( meEBDataCertificationSummaryMap_ ) meEBDataCertificationSummaryMap_->setBinContent( iptt+1, iett+1, xcert );
00228
00229 sumCertEB[ism-1] += xcert;
00230 nValidChannelsEB[ism-1]++;
00231
00232 sumCert += xcert;
00233 nValidChannels++;
00234
00235 }
00236 }
00237
00238 if( meEBDataCertificationSummary_ ) {
00239 if( nValidChannels>0 ) meEBDataCertificationSummary_->Fill( sumCert/nValidChannels );
00240 else meEBDataCertificationSummary_->Fill( 0.0 );
00241 }
00242
00243 for (int i = 0; i < 36; i++) {
00244 if( meEBDataCertification_[i] ) {
00245 if( nValidChannelsEB[i]>0 ) meEBDataCertification_[i]->Fill( sumCertEB[i]/nValidChannelsEB[i] );
00246 else meEBDataCertification_[i]->Fill( 0.0 );
00247 }
00248 }
00249
00250 }
00251
00252 void EBDataCertificationTask::beginRun(const edm::Run& r, const edm::EventSetup& c) {
00253
00254 if ( ! mergeRuns_ ) this->reset();
00255
00256 }
00257
00258 void EBDataCertificationTask::endRun(const edm::Run& r, const edm::EventSetup& c) {
00259
00260 this->reset();
00261
00262 char histo[200];
00263
00264 MonitorElement* me;
00265
00266 sprintf(histo, (prefixME_ + "/EventInfo/reportSummaryMap").c_str());
00267 me = dqmStore_->get(histo);
00268 hDQM_ = UtilsClient::getHisto<TH2F*>( me, cloneME_, hDQM_ );
00269
00270 sprintf(histo, (prefixME_ + "/EventInfo/DAQSummaryMap").c_str());
00271 me = dqmStore_->get(histo);
00272 hDAQ_ = UtilsClient::getHisto<TH2F*>( me, cloneME_, hDAQ_ );
00273
00274 sprintf(histo, (prefixME_ + "/EventInfo/DCSSummaryMap").c_str());
00275 me = dqmStore_->get(histo);
00276 hDCS_ = UtilsClient::getHisto<TH2F*>( me, cloneME_, hDCS_ );
00277
00278 float sumCert = 0.;
00279 float sumCertEB[36];
00280 int nValidChannels = 0;
00281 int nValidChannelsEB[36];
00282
00283 for (int i = 0; i < 36; i++) {
00284 sumCertEB[i] = 0.;
00285 nValidChannelsEB[i] = 0;
00286 }
00287
00288 for ( int iett = 0; iett < 34; iett++ ) {
00289 for ( int iptt = 0; iptt < 72; iptt++ ) {
00290
00291 float xvalDQM, xvalDAQ, xvalDCS;
00292 xvalDQM = xvalDAQ = xvalDCS = -1.;
00293 float xcert = -1;
00294
00295 if ( hDQM_ ) xvalDQM = hDQM_->GetBinContent( iptt+1, iett+1 );
00296 if ( hDAQ_ ) xvalDAQ = hDAQ_->GetBinContent( iptt+1, iett+1 );
00297 if ( hDCS_ ) xvalDCS = hDCS_->GetBinContent( iptt+1, iett+1 );
00298
00299 if ( xvalDQM == -1 || ( xvalDAQ == -1 && xvalDCS == -1 ) ) {
00300
00301 xcert = 0.0;
00302 } else {
00303
00304 xcert = std::abs(xvalDQM) * std::abs(xvalDAQ) * std::abs(xvalDCS);
00305 }
00306
00307 if ( meEBDataCertificationSummaryMap_ ) meEBDataCertificationSummaryMap_->setBinContent( iptt+1, iett+1, xcert );
00308
00309 int ism = ( iett<17 ) ? iptt/4+1 : 18+iptt/4+1;
00310
00311 sumCertEB[ism-1] += xcert;
00312 nValidChannelsEB[ism-1]++;
00313
00314 sumCert += xcert;
00315 nValidChannels++;
00316
00317 }
00318 }
00319
00320 if( meEBDataCertificationSummary_ ) {
00321 if( nValidChannels>0 ) {
00322 meEBDataCertificationSummary_->Fill( sumCert/nValidChannels );
00323 } else {
00324 meEBDataCertificationSummary_->Fill( 0.0 );
00325 }
00326 }
00327
00328 for (int i = 0; i < 36; i++) {
00329 if( meEBDataCertification_[i] ) {
00330 if( nValidChannelsEB[i]>0 ) {
00331 meEBDataCertification_[i]->Fill( sumCertEB[i]/nValidChannelsEB[i] );
00332 } else {
00333 meEBDataCertification_[i]->Fill( 0.0 );
00334 }
00335 }
00336 }
00337
00338 }
00339
00340 void EBDataCertificationTask::reset(void) {
00341
00342 if ( meEBDataCertificationSummary_ ) meEBDataCertificationSummary_->Reset();
00343
00344 for (int i = 0; i < 36; i++) {
00345 if ( meEBDataCertification_[i] ) meEBDataCertification_[i]->Reset();
00346 }
00347
00348 if ( meEBDataCertificationSummaryMap_ ) meEBDataCertificationSummaryMap_->Reset();
00349
00350 }
00351
00352
00353 void EBDataCertificationTask::cleanup(void){
00354
00355 if ( cloneME_ ) {
00356 if( hDQM_ ) delete hDQM_;
00357 if( hDAQ_ ) delete hDAQ_;
00358 if( hDCS_ ) delete hDCS_;
00359 if( hIntegrityByLumi_ ) delete hIntegrityByLumi_;
00360 if( hFrontendByLumi_ ) delete hFrontendByLumi_;
00361 if( hSynchronizationByLumi_ ) delete hSynchronizationByLumi_;
00362 }
00363 hDQM_ = 0;
00364 hDAQ_ = 0;
00365 hDCS_ = 0;
00366 hIntegrityByLumi_ = 0;
00367 hFrontendByLumi_ = 0;
00368 hSynchronizationByLumi_ = 0;
00369
00370 if ( dqmStore_ ) {
00371 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo");
00372 if ( meEBDataCertificationSummary_ ) dqmStore_->removeElement( meEBDataCertificationSummary_->getName() );
00373 if ( meEBDataCertificationSummaryMap_ ) dqmStore_->removeElement( meEBDataCertificationSummaryMap_->getName() );
00374
00375 dqmStore_->setCurrentFolder(prefixME_ + "/EventInfo/CertificationContents");
00376 for (int i = 0; i < 36; i++) {
00377 if ( meEBDataCertification_[i] ) dqmStore_->removeElement( meEBDataCertification_[i]->getName() );
00378 }
00379 }
00380
00381 }
00382
00383 void EBDataCertificationTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00384
00385 }