CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DQM/HcalMonitorTasks/src/HcalCoarsePedestalMonitor.cc

Go to the documentation of this file.
00001 #include "DQM/HcalMonitorTasks/interface/HcalCoarsePedestalMonitor.h"
00002 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
00003 #include <cmath>
00004 
00005 #include "FWCore/Common/interface/TriggerNames.h" 
00006 #include "DataFormats/Common/interface/TriggerResults.h"
00007 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00008 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00009 #include "Geometry/HcalTowerAlgo/src/HcalHardcodeGeometryData.h" // for eta bounds
00010 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00011 #include "FWCore/Framework/interface/LuminosityBlock.h"
00012 
00013 // constructor
00014 HcalCoarsePedestalMonitor::HcalCoarsePedestalMonitor(const edm::ParameterSet& ps) 
00015 {
00016   Online_                = ps.getUntrackedParameter<bool>("online",false);
00017   mergeRuns_             = ps.getUntrackedParameter<bool>("mergeRuns",false);
00018   enableCleanup_         = ps.getUntrackedParameter<bool>("enableCleanup",false);
00019   debug_                 = ps.getUntrackedParameter<int>("debug",0);
00020   prefixME_              = ps.getUntrackedParameter<std::string>("subSystemFolder","Hcal/");
00021   if (prefixME_.substr(prefixME_.size()-1,prefixME_.size())!="/")
00022     prefixME_.append("/");
00023   subdir_                = ps.getUntrackedParameter<std::string>("TaskFolder","CoarsePedestalMonitor_Hcal"); 
00024   if (subdir_.size()>0 && subdir_.substr(subdir_.size()-1,subdir_.size())!="/")
00025     subdir_.append("/");
00026   subdir_=prefixME_+subdir_;
00027   AllowedCalibTypes_     = ps.getUntrackedParameter<std::vector<int> > ("AllowedCalibTypes");
00028   skipOutOfOrderLS_      = ps.getUntrackedParameter<bool>("skipOutOfOrderLS",true);
00029   NLumiBlocks_           = ps.getUntrackedParameter<int>("NLumiBlocks",4000);
00030   makeDiagnostics_       = ps.getUntrackedParameter<bool>("makeDiagnostics",false);
00031   digiLabel_             = ps.getUntrackedParameter<edm::InputTag>("digiLabel");
00032   ADCDiffThresh_         = ps.getUntrackedParameter<double>("ADCDiffThresh",1.);
00033   minEvents_             = ps.getUntrackedParameter<int>("minEvents",100); // minimum number of events needed before histograms are filled
00034   excludeHORing2_       = ps.getUntrackedParameter<bool>("excludeHORing2",false);
00035 
00036 }
00037 
00038 
00039 // destructor
00040 HcalCoarsePedestalMonitor::~HcalCoarsePedestalMonitor() {}
00041 
00042 
00043 
00044 void HcalCoarsePedestalMonitor::cleanup()
00045 {
00046   // Need to add code to clear out subfolders as well?
00047   if (debug_>0) std::cout <<"HcalCoarsePedestalMonitor::cleanup()"<<std::endl;
00048   if (!enableCleanup_) return;
00049   if (dbe_)
00050     {
00051       // removeContents doesn't remove subdirectories
00052       dbe_->setCurrentFolder(subdir_);
00053       dbe_->removeContents();
00054       dbe_->setCurrentFolder(subdir_+"CoarsePedestal_parameters");
00055       dbe_->removeContents();
00056       dbe_->setCurrentFolder(subdir_+"LSvalues");
00057       dbe_->removeContents();
00058     } // if(dbe_)
00059 
00060 } // void HcalCoarsePedestalMonitor::cleanup();
00061 
00062 
00063 void HcalCoarsePedestalMonitor::endRun(const edm::Run& run, const edm::EventSetup& c)
00064 {
00065   // Anything to do here?
00066 }
00067 
00068 void HcalCoarsePedestalMonitor::endJob()
00069 {
00070   if (debug_>0) std::cout <<"HcalCoarsePedestalMonitor::endJob()"<<std::endl;
00071   if (enableCleanup_) cleanup(); // when do we force cleanup?
00072 }
00073 
00074 
00075 void HcalCoarsePedestalMonitor::setup()
00076 {
00077   // Call base class setup
00078   HcalBaseDQMonitor::setup();
00079   if (!dbe_) return;
00080 
00081   /******* Set up all histograms  ********/
00082   if (debug_>1)
00083     std::cout <<"<HcalCoarsePedestalMonitor::beginRun>  Setting up histograms"<<std::endl;
00084 
00085   std::ostringstream name;
00086   dbe_->setCurrentFolder(subdir_ +"CoarsePedestalSumPlots");
00087   SetupEtaPhiHists(CoarsePedestalsSumByDepth,"Coarse Pedestal Summed Map","");
00088   SetupEtaPhiHists(CoarsePedestalsOccByDepth,"Coarse Pedestal Occupancy Map","");
00089   for (unsigned int i=0;i<CoarsePedestalsSumByDepth.depth.size();++i)
00090     (CoarsePedestalsSumByDepth.depth[i]->getTH2F())->SetOption("colz");
00091   for (unsigned int i=0;i<CoarsePedestalsOccByDepth.depth.size();++i)
00092     (CoarsePedestalsOccByDepth.depth[i]->getTH2F())->SetOption("colz");
00093 
00094 
00095   dbe_->setCurrentFolder(subdir_+"CoarsePedestal_parameters");
00096   MonitorElement* ADCDiffThresh = dbe_->bookFloat("ADCdiff_Problem_Threshold");
00097   ADCDiffThresh->Fill(ADCDiffThresh_);
00098   MonitorElement* minevents = dbe_->bookInt("minEventsNeededForPedestalCalculation");
00099   minevents->Fill(minEvents_);
00100   MonitorElement* excludeHORing2 = dbe_->bookInt("excludeHORing2");
00101   if (excludeHORing2_==true)
00102     excludeHORing2->Fill(1);
00103   else
00104     excludeHORing2->Fill(0);
00105 
00106   this->reset();
00107   return;
00108 } // void HcalCoarsePedestalMonitor::setup()
00109 
00110 void HcalCoarsePedestalMonitor::beginRun(const edm::Run& run, const edm::EventSetup& c)
00111 {
00112   HcalBaseDQMonitor::beginRun(run,c);
00113   if (mergeRuns_ && tevt_>0) return; // don't reset counters if merging runs
00114 
00115   if (tevt_==0) this->setup(); // create all histograms; not necessary if merging runs together
00116   if (mergeRuns_==false) this->reset(); // call reset at start of all runs
00117 } // void HcalCoarsePedestalMonitor::beginRun()
00118 
00119 
00120 void HcalCoarsePedestalMonitor::analyze(edm::Event const&e, edm::EventSetup const&s)
00121 {
00122   if (!IsAllowedCalibType()) return;
00123   if (LumiInOrder(e.luminosityBlock())==false) return;
00124 
00125 
00126   // try to get digis
00127   edm::Handle<HBHEDigiCollection> hbhe_digi;
00128   edm::Handle<HODigiCollection> ho_digi;
00129   edm::Handle<HFDigiCollection> hf_digi;
00130   edm::Handle<HcalUnpackerReport> report;
00131 
00132   if (!(e.getByLabel(digiLabel_,hbhe_digi)))
00133     {
00134       edm::LogWarning("HcalCoarsePedestalMonitor")<< digiLabel_<<" hbhe_digi not available";
00135       return;
00136     }
00137   
00138   if (!(e.getByLabel(digiLabel_,hf_digi)))
00139     {
00140       edm::LogWarning("HcalCoarsePedestalMonitor")<< digiLabel_<<" hf_digi not available";
00141       return;
00142     }
00143   if (!(e.getByLabel(digiLabel_,ho_digi)))
00144     {
00145       edm::LogWarning("HcalCoarsePedestalMonitor")<< digiLabel_<<" ho_digi not available";
00146       return;
00147     }
00148   if (!(e.getByLabel(digiLabel_,report)))
00149     {
00150       edm::LogWarning("HcalCoarsePedestalMonitor")<< digiLabel_<<" unpacker report not available";
00151       return;
00152     }
00153 
00154   // all objects grabbed; event is good
00155   if (debug_>1) std::cout <<"\t<HcalCoarsePedestalMonitor::analyze>  Processing good event! event # = "<<ievt_<<std::endl;
00156 
00157   HcalBaseDQMonitor::analyze(e,s); // base class increments ievt_, etc. counters
00158 
00159   // Digi collection was grabbed successfully; process the Event
00160   processEvent(*hbhe_digi, *ho_digi, *hf_digi, *report);
00161 
00162 } //void HcalCoarsePedestalMonitor::analyze(...)
00163 
00164 
00165 void HcalCoarsePedestalMonitor::processEvent(const HBHEDigiCollection& hbhe,
00166                                              const HODigiCollection& ho,
00167                                              const HFDigiCollection& hf,
00168                                              const HcalUnpackerReport& report)
00169 { 
00170   if(!dbe_) 
00171     { 
00172       if(debug_) 
00173         std::cout <<"HcalCoarsePedestalMonitor::processEvent   DQMStore not instantiated!!!"<<std::endl; 
00174       return; 
00175     }
00176 
00177   // Skip events in which minimal good digis found -- still getting some strange (calib?) events through DQM
00178   
00179   unsigned int allgooddigis= hbhe.size()+ho.size()+hf.size();
00180   // bad threshold:  ignore events in which bad outnumber good by more than 100:1
00181   // (one RBX in HBHE seems to send valid data occasionally even on QIE resets, which is why we can't just require allgooddigis==0 when looking for events to skip)
00182   if ((allgooddigis==0) ||
00183       (1.*report.badQualityDigis()>100*allgooddigis))
00184     {
00185       return;
00186     }
00187 
00189 
00190   unsigned int digisize=0;
00191   int depth=0, iphi=0, ieta=0, binEta=-9999;
00192 
00193   double value=0;
00194 
00195   for (HBHEDigiCollection::const_iterator j=hbhe.begin(); j!=hbhe.end(); ++j)
00196     {
00197         const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00198         if (digi.id().subdet()==HcalBarrel)
00199           {if (!HBpresent_) continue;}
00200         else if (digi.id().subdet()==HcalEndcap)
00201           {if (!HEpresent_) continue;}
00202         else continue;
00203         digisize=digi.size();
00204         if (digisize<8) 
00205           continue;
00206         
00207         depth=digi.id().depth();
00208         iphi=digi.id().iphi();
00209         ieta=digi.id().ieta();
00210         
00211         digi.id().subdet()==HcalBarrel ? 
00212           binEta=CalcEtaBin(HcalBarrel, ieta, depth) :
00213           binEta=CalcEtaBin(HcalEndcap, ieta, depth);
00214           
00215         // 'value' is the average pedestal over the 8 time slices.
00216         // In the CoarsePedestal client, we will divide the summed value by Nevents (in underflow bin)
00217         // in order to calculate average pedestals.
00218         value=0;
00219         for (unsigned int i=0;i<8;++i)
00220           value+=digi.sample(i).adc()/8.;
00221         pedestalsum_[binEta][iphi-1][depth-1]+=value;
00222         ++pedestalocc_[binEta][iphi-1][depth-1];
00223     }
00224 
00225 
00227   if (HOpresent_)
00228     {
00229       for (HODigiCollection::const_iterator j=ho.begin(); j!=ho.end(); ++j)
00230         {
00231           const HODataFrame digi = (const HODataFrame)(*j);
00232           digisize=digi.size();
00233           if (digisize<8) 
00234             continue;
00235           
00236           depth=digi.id().depth();
00237           iphi=digi.id().iphi();
00238           ieta=digi.id().ieta();
00239           binEta=CalcEtaBin(HcalOuter, ieta, depth);
00240 
00241           // Don't fill cells that are part of HO ring 2 if an exclusion is applied
00242           if (excludeHORing2_==true && abs(ieta)>10 && isSiPM(ieta,iphi,4)==false)
00243             continue;
00244 
00245           // 'value' is the average pedestal over the 8 time slices.
00246           // In the CoarsePedestal client, we will divide the summed value by Nevents (in underflow bin)
00247           // in order to calculate average pedestals.
00248           value=0;
00249           for (unsigned int i=0;i<8;++i)
00250             value+=digi.sample(i).adc()/8.;
00251           pedestalsum_[binEta][iphi-1][depth-1]+=value;
00252           ++pedestalocc_[binEta][iphi-1][depth-1];
00253         } // for (HODigiCollection)
00254     }
00255   
00257   if (HFpresent_)
00258     {
00259       for (HFDigiCollection::const_iterator j=hf.begin(); j!=hf.end(); ++j)
00260         {
00261           const HFDataFrame digi = (const HFDataFrame)(*j);
00262           digisize=digi.size();
00263           if (digisize<8) 
00264             continue;
00265           digisize=digi.size();
00266           if (digisize<8) 
00267             continue;
00268           
00269           depth=digi.id().depth();
00270           iphi=digi.id().iphi();
00271           ieta=digi.id().ieta();
00272           binEta=CalcEtaBin(HcalForward, ieta, depth);
00273           // 'value' is the average pedestal over the 8 time slices.
00274           // In the CoarsePedestal client, we will divide the summed value by Nevents (in underflow bin)
00275           // in order to calculate average pedestals.
00276           value=0;
00277           for (unsigned int i=0;i<8;++i)
00278             value+=digi.sample(i).adc()/8.;
00279           pedestalsum_[binEta][iphi-1][depth-1]+=value;
00280           ++pedestalocc_[binEta][iphi-1][depth-1];
00281         } // for (HFDigiCollection)
00282     } // if (HFpresent_)
00283 
00284   return;
00285 } // void HcalCoarsePedestalMonitor::processEvent(...)
00286 
00287 
00288 void HcalCoarsePedestalMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00289                                              const edm::EventSetup& c) 
00290 {
00291   HcalBaseDQMonitor::beginLuminosityBlock(lumiSeg,c);
00292   ProblemsCurrentLB->Reset();
00293 }
00294 
00295 void HcalCoarsePedestalMonitor::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00296                                            const edm::EventSetup& c)
00297 {
00298   if (LumiInOrder(lumiSeg.luminosityBlock())==false) return;
00299   fill_Nevents();
00300   return;
00301 }
00302 void HcalCoarsePedestalMonitor::fill_Nevents()
00303 {
00304 
00305   // require minimum number of events before processing
00306   if (ievt_<minEvents_)
00307     return;
00308 
00309 
00310   // Set underflow bin to the number of pedestal events processed
00311   // (Assume that number of ped events is the same for all channels/subdetectors)
00312   for (unsigned int i=0;i<CoarsePedestalsSumByDepth.depth.size();++i)
00313     CoarsePedestalsSumByDepth.depth[i]->setBinContent(0,0,ievt_);
00314   for (unsigned int i=0;i<CoarsePedestalsOccByDepth.depth.size();++i)
00315     CoarsePedestalsOccByDepth.depth[i]->setBinContent(0,0,ievt_);
00316 
00317 
00318   int iphi=-1, ieta=-99, idepth=0, calcEta=-99;
00319   // Loop over all depths, eta, phi
00320   for (int d=0;d<4;++d)
00321     {
00322       idepth=d+1;
00323       for (int phi=0;phi<72;++phi)
00324         {
00325           iphi=phi+1; // actual iphi value
00326           for (int eta=0;eta<83;++eta)
00327             {
00328               ieta=eta-41; // actual ieta value;
00329               if (validDetId(HcalBarrel, ieta, iphi, idepth))
00330                 {
00331                   calcEta = CalcEtaBin(HcalBarrel,ieta,idepth);
00332                   CoarsePedestalsSumByDepth.depth[d]->Fill(ieta,iphi,pedestalsum_[calcEta][phi][d]);
00333                   CoarsePedestalsOccByDepth.depth[d]->Fill(ieta,iphi,pedestalocc_[calcEta][phi][d]);
00334                 }
00335               if (validDetId(HcalEndcap, ieta, iphi, idepth))
00336                 {
00337                   calcEta = CalcEtaBin(HcalEndcap,ieta,idepth);
00338                   CoarsePedestalsSumByDepth.depth[d]->Fill(ieta,iphi,pedestalsum_[calcEta][phi][d]);
00339                   CoarsePedestalsOccByDepth.depth[d]->Fill(ieta,iphi,pedestalocc_[calcEta][phi][d]);
00340                 }
00341               if (validDetId(HcalOuter, ieta, iphi, idepth))
00342                 {
00343                   calcEta = CalcEtaBin(HcalOuter,ieta,idepth);
00344                   CoarsePedestalsSumByDepth.depth[d]->Fill(ieta,iphi,pedestalsum_[calcEta][phi][d]);
00345                   CoarsePedestalsOccByDepth.depth[d]->Fill(ieta,iphi,pedestalocc_[calcEta][phi][d]);
00346                 }
00347               if (validDetId(HcalForward, ieta, iphi, idepth))
00348                 {
00349                   calcEta = CalcEtaBin(HcalBarrel,ieta,idepth);
00350                   int zside=ieta/abs(ieta);
00351                   CoarsePedestalsSumByDepth.depth[d]->Fill(ieta+zside,iphi,pedestalsum_[calcEta][phi][d]);
00352                   CoarsePedestalsOccByDepth.depth[d]->Fill(ieta+zside,iphi,pedestalocc_[calcEta][phi][d]);
00353                 }
00354             }
00355         }
00356     }
00357   // Fill unphysical bins
00358   FillUnphysicalHEHFBins(CoarsePedestalsSumByDepth);
00359   FillUnphysicalHEHFBins(CoarsePedestalsOccByDepth);
00360 
00361   return;
00362 } // void HcalCoarsePedestalMonitor::fill_Nevents()
00363 
00364 
00365 void HcalCoarsePedestalMonitor::reset()
00366 {
00367   // then reset the MonitorElements
00368   zeroCounters();
00369   CoarsePedestalsSumByDepth.Reset();
00370   CoarsePedestalsOccByDepth.Reset();
00371   return;
00372 }
00373 
00374 void HcalCoarsePedestalMonitor::zeroCounters()
00375 {
00376   for (int i=0;i<85;++i)
00377     for (int j=0;j<72;++j)
00378       for (int k=0;k<4;++k)
00379         {
00380           pedestalsum_[i][j][k]=0;
00381           pedestalocc_[i][j][k]=0;
00382         }
00383 }
00384 DEFINE_FWK_MODULE(HcalCoarsePedestalMonitor);
00385