CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DQM/HcalMonitorClient/src/HcalMonitorClient.cc

Go to the documentation of this file.
00001 /*
00002  * \file HcalMonitorClient.cc
00003  * 
00004  * $Date: 2011/04/18 20:12:14 $
00005  * $Revision: 1.102 $
00006  * \author J. Temple
00007  * 
00008  */
00009 
00010 #include "DQM/HcalMonitorClient/interface/HcalMonitorClient.h"
00011 #include "DQM/HcalMonitorClient/interface/HcalDeadCellClient.h"
00012 #include "DQM/HcalMonitorClient/interface/HcalHotCellClient.h"
00013 #include "DQM/HcalMonitorClient/interface/HcalRecHitClient.h"
00014 #include "DQM/HcalMonitorClient/interface/HcalRawDataClient.h"
00015 #include "DQM/HcalMonitorClient/interface/HcalDigiClient.h"
00016 #include "DQM/HcalMonitorClient/interface/HcalTrigPrimClient.h"
00017 #include "DQM/HcalMonitorClient/interface/HcalBeamClient.h"
00018 #include "DQM/HcalMonitorClient/interface/HcalNZSClient.h"
00019 #include "DQM/HcalMonitorClient/interface/HcalSummaryClient.h"
00020 #include "DQM/HcalMonitorClient/interface/HcalDetDiagPedestalClient.h"
00021 #include "DQM/HcalMonitorClient/interface/HcalDetDiagLaserClient.h"
00022 #include "DQM/HcalMonitorClient/interface/HcalDetDiagLEDClient.h"
00023 #include "DQM/HcalMonitorClient/interface/HcalDetDiagNoiseMonitorClient.h"
00024 #include "DQM/HcalMonitorClient/interface/HcalDetDiagTimingClient.h"
00025 #include "DQM/HcalMonitorClient/interface/HcalCoarsePedestalClient.h"
00026 
00027 #include "FWCore/Framework/interface/Run.h"
00028 #include "FWCore/Framework/interface/LuminosityBlock.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032 #include "FWCore/ServiceRegistry/interface/Service.h"
00033 #include "FWCore/Framework/interface/ESHandle.h"
00034 
00035 #include "DQMServices/Core/interface/MonitorElement.h"
00036 #include "DQMServices/Core/interface/DQMStore.h"
00037 
00038 #include "CondFormats/HcalObjects/interface/HcalChannelStatus.h"
00039 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00040 #include "CondFormats/HcalObjects/interface/HcalCondObjectContainer.h"
00041 #include "CondFormats/DataRecord/interface/HcalChannelQualityRcd.h"
00042 
00043 #include "CalibCalorimetry/HcalAlgos/interface/HcalDbASCIIIO.h"
00044 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00045 
00046 #include <iostream>
00047 #include <iomanip>
00048 #include <fstream>
00049 
00050 #include "TROOT.h"
00051 #include "TH1.h"
00052 
00053 //'using' declarations should only be used within classes/functions, and 'using namespace std;' should not be used,
00054 // according to Bill Tanenbaum  (DQM development hypernews, 25 March 2010)
00055 
00056 HcalMonitorClient::HcalMonitorClient(const edm::ParameterSet& ps)
00057 {
00058   debug_ = ps.getUntrackedParameter<int>("debug",0);
00059   inputFile_ = ps.getUntrackedParameter<std::string>("inputFile","");
00060   mergeRuns_ = ps.getUntrackedParameter<bool>("mergeRuns", false);
00061   cloneME_ = ps.getUntrackedParameter<bool>("cloneME", true);
00062   prescaleFactor_ = ps.getUntrackedParameter<int>("prescaleFactor", -1);
00063   prefixME_ = ps.getUntrackedParameter<std::string>("subSystemFolder", "Hcal/");
00064   if (prefixME_.substr(prefixME_.size()-1,prefixME_.size())!="/")
00065     prefixME_.append("/");
00066   enableCleanup_ = ps.getUntrackedParameter<bool>("enableCleanup", false);
00067   enabledClients_ = ps.getUntrackedParameter<std::vector<std::string> >("enabledClients", enabledClients_);
00068 
00069   updateTime_ = ps.getUntrackedParameter<int>("UpdateTime",0);
00070   baseHtmlDir_ = ps.getUntrackedParameter<std::string>("baseHtmlDir", "");
00071   htmlUpdateTime_ = ps.getUntrackedParameter<int>("htmlUpdateTime", 0);
00072   htmlFirstUpdate_ = ps.getUntrackedParameter<int>("htmlFirstUpdate",20);
00073   databasedir_   = ps.getUntrackedParameter<std::string>("databaseDir","");
00074   databaseUpdateTime_ = ps.getUntrackedParameter<int>("databaseUpdateTime",0);
00075   databaseFirstUpdate_ = ps.getUntrackedParameter<int>("databaseFirstUpdate",10);
00076 
00077   saveByLumiSection_  = ps.getUntrackedParameter<bool>("saveByLumiSection",false);
00078   Online_                = ps.getUntrackedParameter<bool>("online",false);
00079 
00080 
00081   if (debug_>0)
00082     {
00083       std::cout <<"HcalMonitorClient:: The following clients are enabled:"<<std::endl;
00084       for (unsigned int i=0;i<enabledClients_.size();++i)
00085           std::cout <<enabledClients_[i]<<std::endl;
00086     } // if (debug_>0)
00087 
00088   // Set all EtaPhiHists pointers to 0 to start
00089   ChannelStatus=0; 
00090   ADC_PedestalFromDBByDepth=0;
00091   ADC_WidthFromDBByDepth=0;
00092   fC_PedestalFromDBByDepth=0;
00093   fC_WidthFromDBByDepth=0;
00094 
00095   // Add all relevant clients
00096   clients_.clear();
00097   clients_.reserve(14); // any reason to reserve ahead of time?
00098   summaryClient_=0;
00099 
00100   clients_.push_back(new HcalBaseDQClient((std::string)"HcalMonitorModule",ps));
00101   if (find(enabledClients_.begin(), enabledClients_.end(),"DeadCellMonitor")!=enabledClients_.end())
00102     clients_.push_back(new HcalDeadCellClient((std::string)"DeadCellMonitor",ps));
00103   if (find(enabledClients_.begin(), enabledClients_.end(),"HotCellMonitor")!=enabledClients_.end())
00104     clients_.push_back(new HcalHotCellClient((std::string)"HotCellMonitor",ps));
00105   if (find(enabledClients_.begin(), enabledClients_.end(),"RecHitMonitor")!=enabledClients_.end())
00106     clients_.push_back(new HcalRecHitClient((std::string)"RecHitMonitor",ps));
00107   if (find(enabledClients_.begin(), enabledClients_.end(),"DigiMonitor")!=enabledClients_.end())
00108     clients_.push_back(new HcalDigiClient((std::string)"DigiMonitor",ps));
00109   if (find(enabledClients_.begin(), enabledClients_.end(),"RawDataMonitor")!=enabledClients_.end())
00110     clients_.push_back(new HcalRawDataClient((std::string)"RawDataMonitor",ps));
00111   if (find(enabledClients_.begin(), enabledClients_.end(),"TrigPrimMonitor")!=enabledClients_.end())
00112     clients_.push_back(new HcalTrigPrimClient((std::string)"TrigPrimMonitor",ps));
00113   if (find(enabledClients_.begin(), enabledClients_.end(),"NZSMonitor")!=enabledClients_.end())
00114     clients_.push_back(new HcalNZSClient((std::string)"NZSMonitor",ps));
00115   if (find(enabledClients_.begin(), enabledClients_.end(),"BeamMonitor")!=enabledClients_.end())
00116     clients_.push_back(new HcalBeamClient((std::string)"BeamMonitor",ps));
00117   if (find(enabledClients_.begin(), enabledClients_.end(),"DetDiagPedestalMonitor")!=enabledClients_.end())
00118     clients_.push_back(new HcalDetDiagPedestalClient((std::string)"DetDiagPedestalMonitor",ps));
00119   if (find(enabledClients_.begin(), enabledClients_.end(),"DetDiagLaserMonitor")!=enabledClients_.end())
00120     clients_.push_back(new HcalDetDiagLaserClient((std::string)"DetDiagLaserMonitor",ps));
00121   if (find(enabledClients_.begin(), enabledClients_.end(),"DetDiagLEDMonitor")!=enabledClients_.end())
00122     clients_.push_back(new HcalDetDiagLEDClient((std::string)"DetDiagLEDMonitor",ps));
00123   if (find(enabledClients_.begin(), enabledClients_.end(),"DetDiagNoiseMonitor")!=enabledClients_.end())
00124     clients_.push_back(new HcalDetDiagNoiseMonitorClient((std::string)"DetDiagNoiseMonitor",ps));
00125   if (find(enabledClients_.begin(), enabledClients_.end(),"DetDiagTimingMonitor")!=enabledClients_.end())
00126     clients_.push_back(new HcalDetDiagTimingClient((std::string)"DetDiagTimingMonitor",ps));
00127  if (find(enabledClients_.begin(), enabledClients_.end(),"CoarsePedestalMonitor")!=enabledClients_.end())
00128     clients_.push_back(new HcalCoarsePedestalClient((std::string)"CoarsePedestalMonitor",ps));
00129 
00130   if (find(enabledClients_.begin(), enabledClients_.end(),"Summary")!=enabledClients_.end())
00131     summaryClient_ = new HcalSummaryClient((std::string)"ReportSummaryClient",ps);
00132   
00133 } // HcalMonitorClient constructor
00134 
00135 
00136 HcalMonitorClient::~HcalMonitorClient()
00137 {
00138   if (debug_>0) std::cout <<"<HcalMonitorClient>  Exiting..."<<std::endl;
00139   for (unsigned int i=0;i<clients_.size();++i)
00140     delete clients_[i];
00141   //if (summaryClient_) delete summaryClient_;
00142 
00143 }
00144 
00145 void HcalMonitorClient::beginJob(void)
00146 {
00147 
00148   begin_run_ = false;
00149   end_run_   = false;
00150 
00151   run_=-1;
00152   evt_=-1;
00153   ievt_=0;
00154   jevt_=0;
00155 
00156   current_time_ = time(NULL);
00157   last_time_html_ = 0; 
00158   last_time_db_ = 0;   
00159 
00160   // get hold of back-end interface
00161 
00162   dqmStore_ = edm::Service<DQMStore>().operator->();
00163 
00164   if ( inputFile_.size() != 0 ) 
00165     {
00166       if ( dqmStore_ )    dqmStore_->open(inputFile_);
00167     }
00168 
00169   for ( unsigned int i=0; i<clients_.size();++i ) 
00170     clients_[i]->beginJob();
00171 
00172   if ( summaryClient_ ) summaryClient_->beginJob();
00173   
00174 
00175 } // void HcalMonitorClient::beginJob(void)
00176 
00177 
00178 void HcalMonitorClient::beginRun(const edm::Run& r, const edm::EventSetup& c) 
00179 {
00180   if (debug_>0) std::cout <<"<HcalMonitorClient::beginRun(r,c)>"<<std::endl;
00181   begin_run_ = true;
00182   end_run_   = false;
00183 
00184   run_=r.id().run();
00185   evt_=0;
00186   jevt_=0;
00187   htmlcounter_=0;
00188 
00189   // Store list of bad channels and their values
00190   std::map <HcalDetId, unsigned int> badchannelmap; 
00191   badchannelmap.clear();
00192 
00193   // Let's get the channel status quality
00194   edm::ESHandle<HcalChannelQuality> p;
00195   c.get<HcalChannelQualityRcd>().get(p);
00196   chanquality_= new HcalChannelQuality(*p.product());
00197  
00198   if (dqmStore_ && ChannelStatus==0)
00199     {
00200       dqmStore_->setCurrentFolder(prefixME_+"HcalInfo/ChannelStatus");
00201       ChannelStatus=new EtaPhiHists;
00202       ChannelStatus->setup(dqmStore_,"ChannelStatus");
00203       std::stringstream x;
00204       for (unsigned int d=0;d<ChannelStatus->depth.size();++d)
00205         {
00206           ChannelStatus->depth[d]->Reset();
00207           x<<"1+log2(status) for HCAL depth "<<d+1;
00208           if (ChannelStatus->depth[d]) ChannelStatus->depth[d]->setTitle(x.str().c_str());
00209           x.str("");
00210         }
00211     }
00212 
00213   edm::ESHandle<HcalDbService> conditions;
00214   c.get<HcalDbRecord>().get(conditions);
00215   // Now let's setup pedestals
00216   if (dqmStore_ )
00217     {
00218       dqmStore_->setCurrentFolder(prefixME_+"HcalInfo/PedestalsFromCondDB");
00219       if (ADC_PedestalFromDBByDepth==0)
00220         {
00221           ADC_PedestalFromDBByDepth = new EtaPhiHists;
00222           ADC_PedestalFromDBByDepth->setup(dqmStore_,"ADC Pedestals From Conditions DB");
00223         }
00224       if (ADC_WidthFromDBByDepth==0)
00225         {
00226           ADC_WidthFromDBByDepth = new EtaPhiHists;
00227           ADC_WidthFromDBByDepth->setup(dqmStore_,"ADC Widths From Conditions DB");
00228         }
00229       if (fC_PedestalFromDBByDepth==0)
00230         {
00231           fC_PedestalFromDBByDepth = new EtaPhiHists;
00232           fC_PedestalFromDBByDepth->setup(dqmStore_,"fC Pedestals From Conditions DB");
00233         }
00234       if (fC_WidthFromDBByDepth==0)
00235         {
00236           fC_WidthFromDBByDepth = new EtaPhiHists;
00237           fC_WidthFromDBByDepth->setup(dqmStore_,"fC Widths From Conditions DB");
00238         }
00239       PlotPedestalValues(*conditions);
00240     }
00241 
00242   // Find only channels with non-zero quality, and add them to badchannelmap
00243   std::vector<DetId> mydetids = chanquality_->getAllChannels();
00244   for (std::vector<DetId>::const_iterator i = mydetids.begin();i!=mydetids.end();++i)
00245     {
00246       if (i->det()!=DetId::Hcal) continue; // not an hcal cell
00247       HcalDetId id=HcalDetId(*i);
00248       int status=(chanquality_->getValues(id))->getValue();
00249       //if (status!=status) status=-1;  // protects against NaN values
00250       // The above line doesn't seem to work in identifying NaNs;  ints for bad values come back as negative numbers (at least in run 146501)
00251       if (status==0) continue;
00252       badchannelmap[id]=status;
00253 
00254       // Fill Channel Status histogram
00255       if (dqmStore_==0) continue;
00256       int depth=id.depth();
00257       if (depth<1 || depth>4) continue;
00258       int ieta=id.ieta();
00259       int iphi=id.iphi();
00260       if (id.subdet()==HcalForward)
00261         ieta>0 ? ++ieta: --ieta;
00262 
00263       double logstatus = 0;
00264       // Fill ChannelStatus value with '-1' when a 'NaN' occurs
00265       if (status<0)
00266         logstatus=-1*(log2(-1.*status)+1);
00267       else
00268         logstatus=log2(1.*status)+1;
00269       if (ChannelStatus->depth[depth-1]) ChannelStatus->depth[depth-1]->Fill(ieta,iphi,logstatus);
00270     }
00271     
00272   for (unsigned int i=0;i<clients_.size();++i)
00273     {
00274       if (clients_[i]->name()=="RawDataMonitor") clients_[i]->setEventSetup(c);
00275       clients_[i]->beginRun();
00276       clients_[i]->setStatusMap(badchannelmap);
00277     }
00278   
00279   if (summaryClient_!=0)
00280     {
00281       summaryClient_->getFriends(clients_);
00282       summaryClient_->beginRun();
00283     }
00284 
00285 } // void HcalMonitorClient::beginRun(const Run& r, const EventSetup& c)
00286 
00287 void HcalMonitorClient::beginRun()
00288 {
00289   // What is the difference between this and beginRun above?
00290   // When would this be called?
00291   begin_run_ = true;
00292   end_run_   = false;
00293   jevt_ = 0;
00294   htmlcounter_=0;
00295 
00296   if (dqmStore_==0 || ChannelStatus!=0) return;
00297   dqmStore_->setCurrentFolder(prefixME_+"HcalInfo");
00298   ChannelStatus=new EtaPhiHists;
00299   ChannelStatus->setup(dqmStore_,"ChannelStatus");
00300   std::stringstream x;
00301   for (unsigned int d=0;d<ChannelStatus->depth.size();++d)
00302     {
00303       x<<"1+log2(status) for HCAL depth "<<d+1;
00304       if (ChannelStatus->depth[d]) ChannelStatus->depth[d]->setTitle(x.str().c_str());
00305       x.str("");
00306     }
00307 } // void HcalMonitorClient::beginRun()
00308 
00309 void HcalMonitorClient::setup(void)
00310 {
00311   // no setup required
00312 }
00313 
00314 void HcalMonitorClient::beginLuminosityBlock(const edm::LuminosityBlock &l, const edm::EventSetup &c) 
00315 {
00316   if (debug_>0) std::cout <<"<HcalMonitorClient::beginLuminosityBlock>"<<std::endl;
00317 } // void HcalMonitorClient::beginLuminosityBlock
00318 
00319 void HcalMonitorClient::analyze(const edm::Event & e, const edm::EventSetup & c)
00320 {
00321   if (debug_>4) 
00322     std::cout <<"HcalMonitorClient::analyze(const edm::Event&, const edm::EventSetup&) ievt_ = "<<ievt_<<std::endl;
00323   ievt_++;
00324   jevt_++;
00325 
00326   run_=e.id().run();
00327   evt_=e.id().event();
00328   if (prescaleFactor_>0 && jevt_%prescaleFactor_==0) this->analyze(e.luminosityBlock());
00329 
00330 } // void HcalMonitorClient::analyze(const edm::Event & e, const edm::EventSetup & c)
00331 
00332 void HcalMonitorClient::analyze(int LS)
00333 {
00334   if (debug_>0)
00335     std::cout <<"HcalMonitorClient::analyze() "<<std::endl;
00336   current_time_ = time(NULL);
00337   // no ievt_, jevt_ counters needed here:  this function gets called at endlumiblock, after default analyze function runs
00338   for (unsigned int i=0;i<clients_.size();++i)
00339     clients_[i]->analyze();
00340   if (summaryClient_!=0)
00341     {
00342       // Always call basic analyze to form histograms for each task
00343       summaryClient_->analyze(LS);
00344       // Call this if LS-by-LS enabling is set to true
00345       if (saveByLumiSection_==true)
00346         summaryClient_->fillReportSummaryLSbyLS(LS);
00347     }
00348 } // void HcalMonitorClient::analyze()
00349 
00350 
00351 void HcalMonitorClient::endLuminosityBlock(const edm::LuminosityBlock &l, const edm::EventSetup &c) 
00352 {
00353   if (debug_>0) std::cout <<"<HcalMonitorClient::endLuminosityBlock>"<<std::endl;
00354   current_time_ = time(NULL);
00355   if (updateTime_>0)
00356     {
00357       if ((current_time_-last_time_update_)<60*updateTime_)
00358         return;
00359       last_time_update_ = current_time_;
00360     }
00361   this->analyze(l.luminosityBlock());
00362 
00363   if (databaseUpdateTime_>0)
00364     {
00365       if (
00366           // first update occurs at after databaseFirstUpdate_ minutes
00367           (last_time_db_==0 && (current_time_-last_time_db_)>=60*databaseFirstUpdate_)
00368           ||
00369           // following updates follow once every databaseUpdateTime_ minutes
00370           ((current_time_-last_time_db_)>=60*databaseUpdateTime_)
00371           )
00372         {
00373           this->writeChannelStatus();
00374           last_time_db_=current_time_;
00375         }
00376     }
00377 
00378   if (htmlUpdateTime_>0)
00379     {
00380       if (
00381           (last_time_html_==0 && (current_time_-last_time_html_)>=60*htmlFirstUpdate_)
00382           // 
00383           ||((current_time_-last_time_html_)>=60*htmlUpdateTime_)
00384           ) // htmlUpdateTime_ in minutes
00385         {
00386           this->writeHtml();
00387           last_time_html_=current_time_;
00388         }
00389     }
00390 
00391 } // void HcalMonitorClient::endLuminosityBlock
00392 
00393 void HcalMonitorClient::endRun(void)
00394 {
00395   begin_run_ = false;
00396   end_run_   = true;
00397 
00398   // Always fill summaryClient at end of run (as opposed to the end-lumi fills, which may just contain info for a single LS)
00399   // At the end of this run, set LS=-1  (LS-based plotting in doesn't work yet anyway)
00400   if (summaryClient_)
00401     summaryClient_->analyze(-1);
00402 
00403   if (databasedir_.size()>0)
00404     this->writeChannelStatus();
00405   // writeHtml takes longer; run it last 
00406   // Also, don't run it if htmlUpdateTime_>0 -- it should have already been run
00407   if (baseHtmlDir_.size()>0 && htmlUpdateTime_==0)
00408     this->writeHtml();
00409 }
00410 
00411 void HcalMonitorClient::endRun(const edm::Run& r, const edm::EventSetup& c) 
00412 {
00413   // Set values here, because the "analyze" method occasionally times out, 
00414   // which keeps the endRun() call from being made.  This causes endJob to
00415   // crash, since end_run_ is still set to false at that point.
00416   begin_run_ = false;
00417   end_run_   = true;
00418 
00419   this->analyze();
00420   this->endRun();
00421 }
00422 
00423 void HcalMonitorClient::endJob(void)
00424 {
00425   // Temporary fix for crash of April 2011 in online DQM
00426   if (Online_==true)
00427     return;
00428 
00429   if (! end_run_)
00430     {
00431       this->analyze();
00432       this->endRun();
00433     }
00434   this->cleanup(); // currently does nothing
00435 
00436   for ( unsigned int i=0; i<clients_.size(); i++ ) 
00437     clients_[i]->endJob();
00438   //if ( summaryClient_ ) summaryClient_->endJob();
00439 
00440 } // void HcalMonitorClient::endJob(void)
00441 
00442 void HcalMonitorClient::cleanup(void)
00443 {
00444   if (!enableCleanup_) return;
00445   // other cleanup?
00446 } // void HcalMonitorClient::cleanup(void)
00447 
00448 
00449 void HcalMonitorClient::writeHtml()
00450 {
00451   if (debug_>0) std::cout << "Preparing HcalMonitorClient html output ..." << std::endl;
00452   
00453 
00454   // global ROOT style
00455   gStyle->Reset("Default");
00456   gStyle->SetCanvasColor(0);
00457   gStyle->SetPadColor(0);
00458   gStyle->SetFillColor(0);
00459   gStyle->SetTitleFillColor(10);
00460   //  gStyle->SetOptStat(0);
00461   gStyle->SetOptStat("ouemr");
00462   gStyle->SetPalette(1);
00463 
00464   char tmp[20];
00465 
00466   if(run_!=-1) sprintf(tmp, "DQM_%s_R%09d_%i", prefixME_.substr(0,prefixME_.size()-1).c_str(),run_,htmlcounter_);
00467   else sprintf(tmp, "DQM_%s_R%09d_%i", prefixME_.substr(0,prefixME_.size()-1).c_str(),0,htmlcounter_);
00468   std::string htmlDir = baseHtmlDir_ + "/" + tmp + "/";
00469   system(("/bin/mkdir -p " + htmlDir).c_str());
00470 
00471   ++htmlcounter_;
00472 
00473   ofstream htmlFile;
00474   htmlFile.open((htmlDir + "index.html").c_str());
00475 
00476   // html page header
00477   htmlFile << "<!DOCTYPE html PUBLIC \"-//W3C//DTD HTML 4.01 Transitional//EN\">  " << std::endl;
00478   htmlFile << "<html>  " << std::endl;
00479   htmlFile << "<head>  " << std::endl;
00480   htmlFile << "  <meta content=\"text/html; charset=ISO-8859-1\"  " << std::endl;
00481   htmlFile << " https-equiv=\"content-type\">  " << std::endl;
00482   htmlFile << "  <title>Hcal Data Quality Monitor</title> " << std::endl;
00483   htmlFile << "</head>  " << std::endl;
00484   htmlFile << "<body>  " << std::endl;
00485   htmlFile << "<br>  " << std::endl;
00486  htmlFile << "<center><h1>Hcal Data Quality Monitor</h1></center>" << std::endl;
00487   htmlFile << "<h2>Run Number:&nbsp;&nbsp;&nbsp;" << std::endl;
00488   htmlFile << "<span style=\"color: rgb(0, 0, 153);\">" << run_ <<"</span></h2> " << std::endl;
00489   htmlFile << "<h2>Events processed:&nbsp;&nbsp;&nbsp;" << std::endl;
00490   htmlFile << "<span style=\"color: rgb(0, 0, 153);\">" << ievt_ <<"</span></h2> " << std::endl;
00491   htmlFile << "<hr>" << std::endl;
00492   htmlFile << "<ul>" << std::endl;
00493 
00494   for (unsigned int i=0;i<clients_.size();++i)
00495     {
00496       if (clients_[i]->validHtmlOutput()==true)
00497         {
00498           clients_[i]->htmlOutput(htmlDir);
00499           // Always print this out?  Or only when validHtmlOutput is true? 
00500           htmlFile << "<table border=0 WIDTH=\"50%\"><tr>" << std::endl;
00501           htmlFile << "<td WIDTH=\"35%\"><a href=\"" << clients_[i]->name_ << ".html"<<"\">"<<clients_[i]->name_<<"</a></td>" << std::endl;
00502           if(clients_[i]->hasErrors_Temp()) htmlFile << "<td bgcolor=red align=center>This monitor task has errors.</td>" << std::endl;
00503           else if(clients_[i]->hasWarnings_Temp()) htmlFile << "<td bgcolor=yellow align=center>This monitor task has warnings.</td>" << std::endl;
00504           else if(clients_[i]->hasOther_Temp()) htmlFile << "<td bgcolor=aqua align=center>This monitor task has messages.</td>" << std::endl;
00505           else htmlFile << "<td bgcolor=lime align=center>This monitor task has no problems</td>" << std::endl;
00506           htmlFile << "</tr></table>" << std::endl;
00507         }
00508     }
00509 
00510   // Add call to reportSummary html output
00511   if (summaryClient_)
00512     {
00513       summaryClient_->htmlOutput(htmlDir);
00514       htmlFile << "<table border=0 WIDTH=\"50%\"><tr>" << std::endl;
00515       htmlFile << "<td WIDTH=\"35%\"><a href=\"" << summaryClient_->name_ << ".html"<<"\">"<<summaryClient_->name_<<"</a></td>" << std::endl;
00516       if(summaryClient_->hasErrors_Temp()) htmlFile << "<td bgcolor=red align=center>This monitor task has errors.</td>" << std::endl;
00517       else if(summaryClient_->hasWarnings_Temp()) htmlFile << "<td bgcolor=yellow align=center>This monitor task has warnings.</td>" << std::endl;
00518       else if(summaryClient_->hasOther_Temp()) htmlFile << "<td bgcolor=aqua align=center>This monitor task has messages.</td>" << std::endl;
00519       else htmlFile << "<td bgcolor=lime align=center>This monitor task has no problems</td>" << std::endl;
00520       htmlFile << "</tr></table>" << std::endl;
00521     }
00522 
00523   htmlFile << "</ul>" << std::endl;
00524 
00525   // html page footer
00526   htmlFile << "</body> " << std::endl;
00527   htmlFile << "</html> " << std::endl;
00528 
00529   htmlFile.close();
00530   if (debug_>0) std::cout << "HcalMonitorClient html output done..." << std::endl;
00531   
00532 } // void HcalMonitorClient::writeHtml()
00533 
00534 void HcalMonitorClient::writeChannelStatus()
00535 {
00536   if (databasedir_.size()==0) return;
00537   if (debug_>0) std::cout <<"<HcalMonitorClient::writeDBfile>  Writing file for database"<<std::endl;
00538 
00539   std::map<HcalDetId, unsigned int> myquality; //map of quality flags as reported by each client
00540   // Get status from all channels (we need to store all channels in case a bad channel suddenly becomes good)
00541   for (std::vector<HcalBaseDQClient*>::size_type i=0;i<clients_.size();++i)
00542     clients_[i]->updateChannelStatus(myquality);
00543 
00544   if (debug_>0) std::cout <<"<HcalMonitorClient::writeChannelStatus()>  myquality size = "<<myquality.size()<<std::endl;
00545 
00546   std::vector<DetId> mydetids = chanquality_->getAllChannels();
00547   HcalChannelQuality* newChanQual = new HcalChannelQuality();
00548 
00549   for (unsigned int i=0;i<mydetids.size();++i)
00550     {
00551       if (mydetids[i].det()!=DetId::Hcal) continue; // not hcal
00552       
00553       HcalDetId id=mydetids[i];
00554       // get original channel status item
00555       const HcalChannelStatus* origstatus=chanquality_->getValues(mydetids[i]);
00556       // make copy of status
00557       HcalChannelStatus* mystatus=new HcalChannelStatus(origstatus->rawId(),origstatus->getValue());
00558       // loop over myquality flags
00559       if (myquality.find(id)!=myquality.end())
00560         {
00561           
00562           // check dead cells
00563           if ((myquality[id]>>HcalChannelStatus::HcalCellDead)&0x1)
00564             mystatus->setBit(HcalChannelStatus::HcalCellDead);
00565           else
00566             mystatus->unsetBit(HcalChannelStatus::HcalCellDead);
00567           // check hot cells
00568           if ((myquality[id]>>HcalChannelStatus::HcalCellHot)&0x1)
00569             mystatus->setBit(HcalChannelStatus::HcalCellHot);
00570           else
00571             mystatus->unsetBit(HcalChannelStatus::HcalCellHot);
00572         } // if (myquality.find_...)
00573       newChanQual->addValues(*mystatus);
00574     } // for (unsigned int i=0;...)
00575   
00576   //Now dump out to text file
00577   std::ostringstream file;
00578   databasedir_=databasedir_+"/"; // add extra slash, just in case
00579   //file <<databasedir_<<"HcalDQMstatus_"<<run_<<".txt";
00580   file <<databasedir_<<"HcalDQMstatus.txt";
00581   std::ofstream outStream(file.str().c_str());
00582   outStream<<"###  Run # "<<run_<<std::endl;
00583   HcalDbASCIIIO::dumpObject (outStream, (*newChanQual));
00584   return;
00585 } // void HcalMonitorClient::writeChannelStatus()
00586 
00587 
00588 void HcalMonitorClient::PlotPedestalValues(const HcalDbService& cond)
00589 {
00590   const HcalQIEShape* shape_ = cond.getHcalShape(); // this one is generic
00591 
00592   double ADC_ped=0;
00593   double ADC_width=0;
00594   double fC_ped=0;
00595   double fC_width=0;
00596   double temp_ADC=0;
00597   double temp_fC=0;
00598 
00599   int ieta=-9999;
00600   int iphi=-9999;
00601   HcalCalibrations calibs_;
00602 
00603   ADC_PedestalFromDBByDepth->Reset();
00604   ADC_WidthFromDBByDepth->Reset();
00605   fC_PedestalFromDBByDepth->Reset();
00606   fC_WidthFromDBByDepth->Reset();
00607 
00608 
00609   for (int subdet=1; subdet<=4;++subdet)
00610     {
00611       for (int depth=0;depth<4;++depth)
00612         {
00613           int etabins= ADC_PedestalFromDBByDepth->depth[depth]->getNbinsX();
00614           int phibins = ADC_PedestalFromDBByDepth->depth[depth]->getNbinsY();
00615           for (int eta=0;eta<etabins;++eta)
00616             {
00617               ieta=CalcIeta(subdet,eta,depth+1);
00618               if (ieta==-9999) continue;
00619               for (int phi=0;phi<phibins;++phi)
00620                 {
00621                   iphi=phi+1;
00622                   if (!validDetId((HcalSubdetector)(subdet), ieta, iphi, depth+1)) continue;
00623                   HcalDetId detid((HcalSubdetector)(subdet), ieta, iphi, depth+1);
00624                   ADC_ped=0;
00625                   ADC_width=0;
00626                   fC_ped=0;
00627                   fC_width=0;
00628                   calibs_= cond.getHcalCalibrations(detid);  
00629                   const HcalPedestalWidth* pedw = cond.getPedestalWidth(detid);
00630                   const HcalQIECoder* channelCoder_ = cond.getHcalCoder(detid);
00631 
00632                   // Loop over capIDs
00633                   for (unsigned int capid=0;capid<4;++capid)
00634                     {
00635                       // Still need to determine how to convert widths to ADC or fC
00636                       // calibs_.pedestal value is always in fC, according to Radek
00637                       temp_fC = calibs_.pedestal(capid);
00638                       fC_ped+= temp_fC;
00639                       // convert to ADC from fC
00640                       temp_ADC=channelCoder_->adc(*shape_,
00641                                                   (float)calibs_.pedestal(capid),
00642                                                   capid);
00643                       ADC_ped+=temp_ADC;
00644                       // Pedestals assumed to be read out in fC
00645                       temp_fC=pedw->getSigma(capid,capid);
00646                       fC_width+=temp_fC;
00647                       temp_ADC=pedw->getSigma(capid,capid)*pow(1.*channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid),capid)/calibs_.pedestal(capid),2);
00648                       ADC_width+=temp_ADC;
00649                     }//capid loop
00650 
00651                   // Pedestal values are average over four cap IDs
00652                   // widths are sqrt(SUM [sigma_ii^2])/4.
00653                   fC_ped/=4.;
00654                   ADC_ped/=4.;
00655 
00656                   // Divide width by 2, or by four?
00657                   // Dividing by 2 gives subtracted results closer to zero -- estimate of variance?
00658                   fC_width=pow(fC_width,0.5)/2.;
00659                   ADC_width=pow(ADC_width,0.5)/2.;
00660 
00661                   if (debug_>1)
00662                     {
00663                       std::cout <<"<HcalMonitorClient::PlotPedestalValues> HcalDet ID = "<<(HcalSubdetector)subdet<<": ("<<ieta<<", "<<iphi<<", "<<depth<<")"<<std::endl;
00664                       std::cout <<"\tADC pedestal = "<<ADC_ped<<" +/- "<<ADC_width<<std::endl;
00665                       std::cout <<"\tfC pedestal = "<<fC_ped<<" +/- "<<fC_width<<std::endl;
00666                     }
00667                   // Shift HF by -/+1 when filling eta-phi histograms
00668                   int zside=0;
00669                   if (subdet==4)
00670                     {
00671                       if (ieta<0) zside=-1;
00672                       else zside=1;
00673                     }
00674                   ADC_PedestalFromDBByDepth->depth[depth]->Fill(ieta+zside,iphi,ADC_ped);
00675                   ADC_WidthFromDBByDepth->depth[depth]->Fill(ieta+zside, iphi, ADC_width);
00676                   fC_PedestalFromDBByDepth->depth[depth]->Fill(ieta+zside,iphi,fC_ped);
00677                   fC_WidthFromDBByDepth->depth[depth]->Fill(ieta+zside, iphi, fC_width);
00678                 } // phi loop
00679             } // eta loop
00680         } //depth loop
00681 
00682     } // subdet loop
00683   FillUnphysicalHEHFBins(*ADC_PedestalFromDBByDepth);
00684   FillUnphysicalHEHFBins(*ADC_WidthFromDBByDepth);
00685   FillUnphysicalHEHFBins(*fC_PedestalFromDBByDepth);
00686   FillUnphysicalHEHFBins(*fC_WidthFromDBByDepth);
00687 
00688   // Center ADC pedestal values near 3 +/- 1
00689   for (unsigned int i=0;i<ADC_PedestalFromDBByDepth->depth.size();++i)
00690   {
00691     ADC_PedestalFromDBByDepth->depth[i]->getTH2F()->SetMinimum(0);
00692     if (ADC_PedestalFromDBByDepth->depth[i]->getTH2F()->GetMaximum()<6)
00693       ADC_PedestalFromDBByDepth->depth[i]->getTH2F()->SetMaximum(6);
00694   }
00695 
00696   for (unsigned int i=0;i<ADC_WidthFromDBByDepth->depth.size();++i)
00697   {
00698     ADC_WidthFromDBByDepth->depth[i]->getTH2F()->SetMinimum(0);
00699     if (ADC_WidthFromDBByDepth->depth[i]->getTH2F()->GetMaximum()<2)
00700       ADC_WidthFromDBByDepth->depth[i]->getTH2F()->SetMaximum(2);
00701   }
00702 
00703 }
00704 
00705 DEFINE_FWK_MODULE(HcalMonitorClient);