CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DQM/HcalMonitorTasks/src/HcalDigiMonitor.cc

Go to the documentation of this file.
00001 #include "DQM/HcalMonitorTasks/interface/HcalDigiMonitor.h"
00002 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
00003 #include <cmath>
00004 
00005 #include "FWCore/Common/interface/TriggerNames.h" 
00006 #include "FWCore/Framework/interface/LuminosityBlock.h"
00007 #include "DataFormats/Common/interface/TriggerResults.h"
00008 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00009 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00010 #include "Geometry/HcalTowerAlgo/src/HcalHardcodeGeometryData.h" // for eta bounds
00011 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00012 
00013 // constructor
00014 HcalDigiMonitor::HcalDigiMonitor(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","DigiMonitor_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   hfRechitLabel_        = ps.getUntrackedParameter<edm::InputTag>("hfRechitLabel");
00033   shapeThresh_           = ps.getUntrackedParameter<int>("shapeThresh",20);
00034   //shapeThresh_ is used for plotting pulse shapes for all digis with pedestal-subtracted ADC sum > shapeThresh_;
00035   shapeThreshHB_ = ps.getUntrackedParameter<int>("shapeThreshHB",shapeThresh_);
00036   shapeThreshHE_ = ps.getUntrackedParameter<int>("shapeThreshHE",shapeThresh_);
00037   shapeThreshHF_ = ps.getUntrackedParameter<int>("shapeThreshHF",shapeThresh_);
00038   shapeThreshHO_ = ps.getUntrackedParameter<int>("shapeThreshHO",shapeThresh_);
00039   
00040   hltresultsLabel_       = ps.getUntrackedParameter<edm::InputTag>("HLTResultsLabel");
00041   MinBiasHLTBits_        = ps.getUntrackedParameter<std::vector<std::string> >("MinBiasHLTBits");
00042   excludeHORing2_       = ps.getUntrackedParameter<bool>("excludeHORing2",false);
00043   excludeHO1P02_        = ps.getUntrackedParameter<bool>("excludeHO1P02",false);
00044 
00045   if (debug_>0)
00046     std::cout <<"<HcalDigiMonitor> Digi shape ADC threshold set to: >" << shapeThresh_ <<" counts above nominal pedestal (3*10)"<< std::endl;
00047   
00048   // Specify which tests to run when looking for problem digis
00049   digi_checkoccupancy_ = ps.getUntrackedParameter<bool>("checkForMissingDigis",false); // off by default -- checked by dead cell monitor
00050   digi_checkcapid_     = ps.getUntrackedParameter<bool>("checkCapID",true);
00051   digi_checkdigisize_  = ps.getUntrackedParameter<bool>("checkDigiSize",true);
00052   digi_checkadcsum_    = ps.getUntrackedParameter<bool>("checkADCsum",true);
00053   digi_checkdverr_     = ps.getUntrackedParameter<bool>("checkDVerr",true);
00054   mindigisizeHBHE_     = ps.getUntrackedParameter<int>("minDigiSizeHBHE",1);
00055   maxdigisizeHBHE_     = ps.getUntrackedParameter<int>("maxDigiSizeHBHE",10);
00056   mindigisizeHO_       = ps.getUntrackedParameter<int>("minDigiSizeHO",1);
00057   maxdigisizeHO_       = ps.getUntrackedParameter<int>("maxDigiSizeHO",10);
00058   mindigisizeHF_       = ps.getUntrackedParameter<int>("minDigiSizeHF",1);
00059   maxdigisizeHF_       = ps.getUntrackedParameter<int>("maxDigiSizeHF",10);
00060 
00061 
00062   badChannelStatusMask_   = ps.getUntrackedParameter<int>("BadChannelStatusMask",
00063                                                           ps.getUntrackedParameter<int>("BadChannelStatusMask",
00064                                                                                         (1<<HcalChannelStatus::HcalCellDead)));  // identify channel status values to mask
00065   if (debug_>1)
00066     {
00067       std::cout <<"<HcalDigiMonitor> Checking for the following problems:"<<std::endl; 
00068       if (digi_checkcapid_) std::cout <<"\tChecking that cap ID rotation is correct;"<<std::endl;
00069       if (digi_checkdigisize_)
00070         {
00071           std::cout <<"\tChecking that HBHE digi size is between ["<<mindigisizeHBHE_<<" - "<<maxdigisizeHBHE_<<"];"<<std::endl;
00072           std::cout <<"\tChecking that HO digi size is between ["<<mindigisizeHO_<<" - "<<maxdigisizeHO_<<"];"<<std::endl;
00073           std::cout <<"\tChecking that HF digi size is between ["<<mindigisizeHF_<<" - "<<maxdigisizeHF_<<"];"<<std::endl;
00074         }
00075       if (digi_checkadcsum_) std::cout <<"\tChecking that ADC sum of digi is greater than 0;"<<std::endl; 
00076       if (digi_checkdverr_) std::cout <<"\tChecking that data valid bit is true and digi error bit is false;\n"<<std::endl;
00077     }
00078   
00079   shutOffOrbitTest_ = ps.getUntrackedParameter<bool>("shutOffOrbitTest",false);
00080   DigiMonitor_ExpectedOrbitMessageTime_=ps.getUntrackedParameter<int>("ExpectedOrbitMessageTime",3559); // -1 means that orbit mismatches won't be checked
00081 
00082   HFtiming_totaltime2D=0;
00083   HFtiming_occupancy2D=0;
00084   HFtiming_etaProfile=0;
00085   HFP_shape=0;
00086   HFM_shape=0;
00087 }
00088 
00089 // destructor
00090 HcalDigiMonitor::~HcalDigiMonitor() {}
00091 
00092 // Checks capid rotation; returns false if no problems with rotation
00093 static bool bitUpset(int last, int now){
00094   if(last ==-1) return false;
00095   int v = last+1; 
00096   if(v==4) v=0;
00097   if(v==now) return false;
00098   return true;
00099 } // static bool bitUpset(...)
00100 
00101 void HcalDigiMonitor::cleanup()
00102 {
00103   // Need to add code to clear out subfolders as well?
00104   if (debug_>0) std::cout <<"HcalDigiMonitor::cleanup()"<<std::endl;
00105   if (!enableCleanup_) return;
00106   if (dbe_)
00107     {
00108       // removeContents doesn't remove subdirectories
00109       dbe_->setCurrentFolder(subdir_);
00110       dbe_->removeContents();
00111       dbe_->setCurrentFolder(subdir_+"digi_parameters");  dbe_->removeContents();
00112       dbe_->setCurrentFolder(subdir_+"bad_digis/bad_digi_occupancy");  dbe_->removeContents();
00113       dbe_->setCurrentFolder(subdir_+"bad_digis/1D_digi_plots");  dbe_->removeContents();
00114       dbe_->setCurrentFolder(subdir_+"bad_digis/badcapID");  dbe_->removeContents();
00115       dbe_->setCurrentFolder(subdir_+"bad_digis/data_invalid_error");  dbe_->removeContents();
00116       dbe_->setCurrentFolder(subdir_+"bad_digis/bad_reportUnpackerErrors");  dbe_->removeContents();
00117       dbe_->setCurrentFolder(subdir_+"bad_digis/baddigisize");  dbe_->removeContents();
00118       dbe_->setCurrentFolder(subdir_+"digi_info");  dbe_->removeContents();
00119       dbe_->setCurrentFolder(subdir_+"bad_digis/badfibBCNoff");  dbe_->removeContents();
00120       dbe_->setCurrentFolder(subdir_+"good_digis/1D_digi_plots");  dbe_->removeContents();
00121       dbe_->setCurrentFolder(subdir_+"good_digis/digi_occupancy");  dbe_->removeContents();
00122       dbe_->setCurrentFolder(subdir_+"bad_digis/bad_digi_occupancy");  dbe_->removeContents();
00123       dbe_->setCurrentFolder(subdir_+"bad_digis");  dbe_->removeContents();
00124       dbe_->setCurrentFolder(subdir_+"good_digis/");  dbe_->removeContents();
00125       dbe_->setCurrentFolder(subdir_+"digi_info/HB");  dbe_->removeContents();
00126       dbe_->setCurrentFolder(subdir_+"digi_info/HE");  dbe_->removeContents();
00127       dbe_->setCurrentFolder(subdir_+"digi_info/HO");  dbe_->removeContents();
00128       dbe_->setCurrentFolder(subdir_+"digi_info/HF");  dbe_->removeContents();
00129       dbe_->setCurrentFolder(subdir_+"LSvalues");
00130       dbe_->removeContents();
00131     } // if(dbe_)
00132 
00133 } // void HcalDigiMonitor::cleanup();
00134 
00135 
00136 void HcalDigiMonitor::endRun(const edm::Run& run, const edm::EventSetup& c)
00137 {
00138   // Anything to do here?
00139 }
00140 
00141 void HcalDigiMonitor::endJob()
00142 {
00143   if (debug_>0) std::cout <<"HcalDigiMonitor::endJob()"<<std::endl;
00144   if (enableCleanup_) cleanup(); // when do we force cleanup?
00145 }
00146 
00147 
00148 void HcalDigiMonitor::setup()
00149 {
00150   // Call base class setup
00151   HcalBaseDQMonitor::setup();
00152   if (!dbe_) return;
00153 
00154   /******* Set up all histograms  ********/
00155   if (debug_>1)
00156     std::cout <<"<HcalDigiMonitor::beginRun>  Setting up histograms"<<std::endl;
00157 
00158   std::ostringstream name;
00159   dbe_->setCurrentFolder(subdir_);
00160   
00161   dbe_->setCurrentFolder(subdir_+"digi_parameters");
00162   MonitorElement* ExpectedOrbit = dbe_->bookInt("ExpectedOrbitMessageTime");
00163   ExpectedOrbit->Fill(DigiMonitor_ExpectedOrbitMessageTime_);
00164 
00165   MonitorElement* shapeT = dbe_->bookInt("DigiShapeThresh");
00166   shapeT->Fill(shapeThresh_);
00167   MonitorElement* shapeTHB = dbe_->bookInt("DigiShapeThreshHB");
00168   shapeTHB->Fill(shapeThreshHB_);
00169   MonitorElement* shapeTHE = dbe_->bookInt("DigiShapeThreshHE");
00170   shapeTHE->Fill(shapeThreshHE_);
00171   MonitorElement* shapeTHO = dbe_->bookInt("DigiShapeThreshHO");
00172   shapeTHO->Fill(shapeThreshHO_);
00173   MonitorElement* shapeTHF = dbe_->bookInt("DigiShapeThreshHF");
00174   shapeTHF->Fill(shapeThreshHF_);
00175   
00176   dbe_->setCurrentFolder(subdir_+"bad_digis/bad_digi_occupancy");
00177   SetupEtaPhiHists(DigiErrorsByDepth,"Bad Digi Map","");
00178   dbe_->setCurrentFolder(subdir_+"bad_digis/1D_digi_plots");
00179   ProblemsVsLB=dbe_->bookProfile("BadDigisVsLB","Number Bad Digis vs Luminosity block;Lumi block;# of Bad digis",
00180                                  NLumiBlocks_,0.5,NLumiBlocks_+0.5,100,0,10000);
00181   ProblemsVsLB_HB=dbe_->bookProfile("HB Bad Quality Digis vs LB","HB Bad Quality Digis vs Luminosity Block",
00182                                      NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00183                                      100,0,10000);   
00184   ProblemsVsLB_HE=dbe_->bookProfile("HE Bad Quality Digis vs LB","HE Bad Quality Digis vs Luminosity Block",
00185                                      NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00186                                      100,0,10000);
00187   ProblemsVsLB_HO=dbe_->bookProfile("HO Bad Quality Digis vs LB","HO Bad Quality Digis vs Luminosity Block",
00188                                      NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00189                                      100,0,10000);
00190   ProblemsVsLB_HF=dbe_->bookProfile("HF Bad Quality Digis vs LB","HF Bad Quality Digis vs Luminosity Block",
00191                                      NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00192                                      100,0,10000);
00193   ProblemsVsLB_HBHEHF=dbe_->bookProfile("HBHEHF Bad Quality Digis vs LB","HBHEHF Bad Quality Digis vs Luminosity Block",
00194                                      NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00195                                      100,0,10000);
00196 
00197   ProblemDigisInLastNLB_HBHEHF_alarm=dbe_->book1D("ProblemDigisInLastNLB_HBHEHF_alarm",
00198                                               "Total Number of ProblemDigis HBHEHF in last 10 LS. Last bin contains OverFlow",
00199                                               100,0,100);
00200 
00201 
00202   if (makeDiagnostics_) 
00203     {
00204       // by default, unpacked digis won't have these errors
00205       dbe_->setCurrentFolder(subdir_+"diagnostics/bad_digis/badcapID");
00206       SetupEtaPhiHists(DigiErrorsBadCapID," Digis with Bad Cap ID Rotation", "");
00207       dbe_->setCurrentFolder(subdir_+"diagnostics/bad_digis/data_invalid_error");
00208       SetupEtaPhiHists(DigiErrorsDVErr," Digis with Data Invalid or Error Bit Set", "");
00209     }
00210   
00211   if (Online_)
00212     {
00213       // Special histograms for Pawel's timing study
00214       dbe_->setCurrentFolder(subdir_+"HFTimingStudy");
00215       HFtiming_etaProfile=dbe_->bookProfile("HFTiming_etaProfile","HFTiming Eta Profile;ieta;average time (time slice)",83,-41.5,41.5,200,0,10);
00216       HFP_shape=dbe_->book1D("HFP_signal_shape","HFP signal shape",10,-0.5,9.5);
00217       HFM_shape=dbe_->book1D("HFM_signal_shape","HFM signal shape",10,-0.5,9.5);
00218       dbe_->setCurrentFolder(subdir_+"HFTimingStudy/sumplots");
00219       HFtiming_totaltime2D=dbe_->book2D("HFTiming_Total_Time","HFTiming Total Time",83,-41.5,41.5,72,0.5,72.5);
00220       HFtiming_occupancy2D=dbe_->book2D("HFTiming_Occupancy","HFTiming Occupancy",83,-41.5,41.5,72,0.5,72.5);
00221     }
00222   
00223   dbe_->setCurrentFolder(subdir_+"bad_digis/bad_reportUnpackerErrors");
00224   SetupEtaPhiHists(DigiErrorsUnpacker," Bad Unpacker Digis", "");
00225   
00226   dbe_->setCurrentFolder(subdir_+"bad_digis/baddigisize");
00227   SetupEtaPhiHists(DigiErrorsBadDigiSize," Digis with Bad Size", "");
00228   
00229   dbe_->setCurrentFolder(subdir_+"digi_info");
00230   
00231   h_valid_digis=dbe_->book1D("ValidEvents","Events with minimum number of valid digis",2,-0.5,1.5);
00232   h_valid_digis->setBinLabel(1,"Valid");
00233   h_valid_digis->setBinLabel(2,"Invalid");
00234   
00235   h_invalid_orbitnumMod103=dbe_->book1D("InvalidDigiEvents_ORN","Orbit Number (mod 103) for Events with Many Unpacker Errors",103,-0.5,102.5);
00236   h_invalid_bcn=dbe_->book1D("InvalidDigiEvents_BCN","Bunch Crossing Number fo Events with Many Unpacker Errors",3464,-0.5,3563.5);
00237 
00238   DigiSize = dbe_->book2D("Digi Size", "Digi Size",4,0,4,20,-0.5,19.5);
00239   DigiSize->setBinLabel(1,"HB",1);
00240   DigiSize->setBinLabel(2,"HE",1);
00241   DigiSize->setBinLabel(3,"HO",1);
00242   DigiSize->setBinLabel(4,"HF",1);
00243   DigiSize->setAxisTitle("Subdetector",1);
00244   DigiSize->setAxisTitle("Digi Size",2);
00245   
00246   dbe_->setCurrentFolder(subdir_+"bad_digis/badfibBCNoff");
00247   SetupEtaPhiHists(DigiErrorsBadFibBCNOff," Digis with non-zero Fiber Orbit Msg Idle BCN Offsets", "");
00248   
00249   dbe_->setCurrentFolder(subdir_+"good_digis/1D_digi_plots");
00250   HBocc_vs_LB=dbe_->bookProfile("HBoccVsLB","HB digi occupancy vs Luminosity Block;Lumi block;# of Good digis",
00251                                  NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00252                                  0,2600);
00253   HEocc_vs_LB=dbe_->bookProfile("HEoccVsLB","HE digi occupancy vs Luminosity Block;Lumi block;# of Good digis",
00254                                  NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00255                                  0,2600);
00256   HOocc_vs_LB=dbe_->bookProfile("HOoccVsLB","HO digi occupancy vs Luminosity Block;Lumi block;# of Good digis",
00257                                  NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00258                                  0,2200);
00259   HFocc_vs_LB=dbe_->bookProfile("HFoccVsLB","HF digi occupancy vs Luminosity Block;Lumi block;# of Good digis",
00260                                  NLumiBlocks_,0.5,NLumiBlocks_+0.5,
00261                                  0,1800);
00262   
00263   dbe_->setCurrentFolder(subdir_+"good_digis/digi_occupancy");
00264   SetupEtaPhiHists(DigiOccupancyByDepth," Digi Eta-Phi Occupancy Map","");
00265   DigiOccupancyPhi= dbe_->book1D("Digi Phi Occupancy Map",
00266                                   "Digi Phi Occupancy Map;i#phi;# of Events",
00267                                   72,0.5,72.5);
00268   DigiOccupancyEta= dbe_->book1D("Digi Eta Occupancy Map",
00269                                   "Digi Eta Occupancy Map;i#eta;# of Events",
00270                                   83,-41.5,41.5);
00271   DigiOccupancyVME = dbe_->book2D("Digi VME Occupancy Map",
00272                                    "Digi VME Occupancy Map;HTR Slot;VME Crate Id",
00273                                    40,-0.25,19.75,18,-0.5,17.5);
00274   
00275   DigiOccupancySpigot = dbe_->book2D("Digi Spigot Occupancy Map",
00276                                       "Digi Spigot Occupancy Map;Spigot;DCC Id",
00277                                       HcalDCCHeader::SPIGOT_COUNT,-0.5,HcalDCCHeader::SPIGOT_COUNT-0.5,
00278                                       36,-0.5,35.5);
00279   
00280   dbe_->setCurrentFolder(subdir_+"bad_digis/bad_digi_occupancy");
00281   DigiErrorVME = dbe_->book2D("Digi VME Error Map",
00282                                "Digi VME Error Map;HTR Slot;VME Crate Id",
00283                                40,-0.25,19.75,18,-0.5,17.5);
00284   
00285   DigiErrorSpigot = dbe_->book2D("Digi Spigot Error Map",
00286                                   "Digi Spigot Error Map;Spigot;DCC Id",
00287                                   HcalDCCHeader::SPIGOT_COUNT,-0.5,HcalDCCHeader::SPIGOT_COUNT-0.5,
00288                                   36,-0.5,35.5);
00289   
00290   dbe_->setCurrentFolder(subdir_+"bad_digis");
00291   int nbins = sizeof(bins_cellcount_new)/sizeof(float)-1;
00292 
00293   DigiBQ = dbe_->book1D("NumBadQualDigis","Number Bad Qual Digis within Digi Collection",nbins, bins_cellcount_new);
00294 
00295   nbins=sizeof(bins_fraccount_new)/sizeof(float)-1;
00296 
00297   DigiBQFrac =  dbe_->book1D("Bad Digi Fraction","Bad Digi Fraction;Bad Quality Digi Fraction for digis in collection; # of Events",
00298                              nbins, bins_fraccount_new);
00299   
00300   nbins = sizeof(bins_cellcount_new)/sizeof(float)-1;
00301   DigiUnpackerErrorCount = dbe_->book1D("Unpacker Error Count", "Number of Bad Digis from Unpacker; Bad Unpacker Digis; # of Events",nbins, bins_cellcount_new);
00302   
00303   nbins=sizeof(bins_fraccount_new)/sizeof(float)-1;
00304   DigiUnpackerErrorFrac = dbe_->book1D("Unpacker Bad Digi Fraction", 
00305                                        "Bad Digis From Unpacker/ (Bad Digis From Unpacker + Good Digis); Bad Unpacker Fraction; # of Events",
00306                                        nbins,bins_fraccount_new);
00307   
00308   dbe_->setCurrentFolder(subdir_+"good_digis/");
00309   DigiNum = dbe_->book1D("NumGoodDigis","Number of Digis;# of Good Digis;# of Events",DIGI_NUM+1,-0.5,DIGI_NUM+1-0.5);
00310     
00311   setupSubdetHists(hbHists,"HB");
00312   setupSubdetHists(heHists,"HE");
00313   setupSubdetHists(hoHists,"HO");
00314   setupSubdetHists(hfHists,"HF");
00315 
00316   this->reset();
00317   return;
00318 } // void HcalDigiMonitor::setup()
00319 
00320 void HcalDigiMonitor::beginRun(const edm::Run& run, const edm::EventSetup& c)
00321 {
00322   HcalBaseDQMonitor::beginRun(run,c);
00323   if (mergeRuns_ && tevt_>0) return; // don't reset counters if merging runs
00324 
00325   if (debug_>1) std::cout <<"\t<HcalDigiMonitor::setup> Getting conditions from DB!"<<std::endl;
00326   c.get<HcalDbRecord>().get(conditions_);
00327 
00328   // Get all pedestals by Cap ID
00329   edm::ESHandle<HcalChannelQuality> p;
00330   c.get<HcalChannelQualityRcd>().get(p);
00331   HcalChannelQuality *chanquality= new HcalChannelQuality(*p.product());
00332   std::vector<DetId> mydetids = chanquality->getAllChannels();
00333   PedestalsByCapId_.clear();
00334 
00335   const HcalQIEShape* shape = conditions_->getHcalShape();
00336   for (std::vector<DetId>::const_iterator chan = mydetids.begin();chan!=mydetids.end();++chan)
00337     {
00338       if (chan->det()!=DetId::Hcal) continue; // not hcal
00339       std::vector <double> peds;  // could be ints, right?
00340       peds.clear();
00341       HcalCalibrations calibs=conditions_->getHcalCalibrations(*chan);
00342       const HcalQIECoder* channelCoder = conditions_->getHcalCoder(*chan);
00343       //double total=0; // use this is we want to calculate average pedestal value
00344       for (int capid=0;capid<4;++capid)
00345         {
00346           // temp_ADC should be an int, right?
00347           double temp_ADC=channelCoder->adc(*shape,(float)calibs.pedestal(capid),capid);
00348           peds.push_back(temp_ADC);
00349           //total=total+temp_ADC;
00350         }
00351       //for (int capid=0;capid<4;++capid) peds.push_back(total/4.); // use this if we just want to use average value
00352       PedestalsByCapId_[*chan]=peds;
00353     } // loop on DetIds
00354 
00355   if (tevt_==0) this->setup(); // create all histograms; not necessary if merging runs together
00356   if (mergeRuns_==false) this->reset(); // call reset at start of all runs
00357   delete chanquality;
00358 
00359   // Get known dead cells for this run
00360   KnownBadCells_.clear();
00361   if (badChannelStatusMask_>0)
00362     {
00363       edm::ESHandle<HcalChannelQuality> p;
00364       c.get<HcalChannelQualityRcd>().get(p);
00365       HcalChannelQuality* chanquality= new HcalChannelQuality(*p.product());
00366       std::vector<DetId> mydetids = chanquality->getAllChannels();
00367       for (std::vector<DetId>::const_iterator i = mydetids.begin();
00368            i!=mydetids.end();
00369            ++i)
00370         {
00371           if (i->det()!=DetId::Hcal) continue; // not an hcal cell
00372           HcalDetId id=HcalDetId(*i);
00373           int status=(chanquality->getValues(id))->getValue();
00374           if ((status & badChannelStatusMask_))
00375             {
00376               KnownBadCells_[id.rawId()]=status;
00377             }
00378         } 
00379         delete chanquality;
00380     } // if (badChannelStatusMask_>0)
00381 
00382 } // void HcalDigiMonitor::beginRun()
00383 
00384 
00385 void HcalDigiMonitor::setupSubdetHists(DigiHists& hist, std::string subdet)
00386 {
00387   if (!dbe_) return;
00388   std::stringstream name;
00389   int nChan=0;
00390   if (subdet=="HB" || subdet=="HE") nChan=2592;
00391   else if (subdet == "HO") nChan=2160;
00392   else if (subdet == "HF") nChan=1728;
00393 
00394   dbe_->setCurrentFolder(subdir_+"digi_info/"+subdet);
00395   hist.shape = dbe_->book1D(subdet+" Digi Shape",subdet+" Digi Shape;Time Slice",10,-0.5,9.5);
00396   hist.shapeThresh = dbe_->book1D(subdet+" Digi Shape - over thresh",
00397                                   subdet+" Digi Shape - over thresh passing trigger and HF HT cuts;Time slice",
00398                                   10,-0.5,9.5);
00399   hist.ThreshCount = dbe_->book1D(subdet+" Total Digis Over Threshold",
00400                                   subdet+" Total Digis Over Threshold",
00401                                   1,-0.5,0.5);
00402   // Create plots of sums of adjacent time slices
00403   for (int ts=0;ts<9;++ts)
00404     {
00405       name<<subdet<<" Plus Time Slices "<<ts<<" and "<<ts+1;
00406       hist.TS_sum_plus.push_back(dbe_->book1D(name.str().c_str(),name.str().c_str(),50, 0., 50.));
00407       name.str("");
00408       name<<subdet<<" Minus Time Slices "<<ts<<" and "<<ts+1;
00409       hist.TS_sum_minus.push_back(dbe_->book1D(name.str().c_str(),name.str().c_str(),50, 0., 50.));
00410       name.str("");
00411     }
00412   hist.presample= dbe_->book1D(subdet+" Digi Presamples",subdet+" Digi Presamples",50,-0.5,49.5);
00413   hist.BQ = dbe_->book1D(subdet+" Bad Quality Digis",subdet+" Bad Quality Digis",nChan+1,-0.5,nChan+0.5);
00414   //(hist.BQ->getTH1F())->LabelsOption("v");
00415   hist.BQFrac = dbe_->book1D(subdet+" Bad Quality Digi Fraction",subdet+" Bad Quality Digi Fraction",DIGI_BQ_FRAC_NBINS,(0-0.5/(DIGI_BQ_FRAC_NBINS-1)),1+0.5/(DIGI_BQ_FRAC_NBINS-1));
00416   hist.DigiFirstCapID = dbe_->book1D(subdet+" Capid 1st Time Slice",subdet+" Capid for 1st Time Slice;CapId (T0)- 1st CapId (T0);# of Events",7,-3.5,3.5);
00417 
00418   hist.DVerr = dbe_->book1D(subdet+" Data Valid Err Bits",subdet+" QIE Data Valid Err Bits",4,-0.5,3.5);
00419   hist.DVerr ->setBinLabel(1,"Err=0, DV=0",1);
00420   hist.DVerr ->setBinLabel(2,"Err=0, DV=1",1);
00421   hist.DVerr ->setBinLabel(3,"Err=1, DV=0",1);
00422   hist.DVerr ->setBinLabel(4,"Err=1, DV=1",1);
00423   hist.CapID = dbe_->book1D(subdet+" CapID",subdet+" CapID",4,-0.5,3.5);
00424   hist.ADC = dbe_->book1D(subdet+" ADC count per time slice",subdet+" ADC count per time slice",200,-0.5,199.5);
00425   hist.ADCsum = dbe_->book1D(subdet+" ADC sum", subdet+" ADC sum",200,-0.5,199.5);
00426   hist.fibBCNOff = dbe_->book1D(subdet+" Fiber Orbit Message Idle BCN Offset", subdet+" Fiber Orbit Message Idle BCN Offset;Offset from Expected",
00427                                  15, -7.5, 7.5);
00428 }
00429 
00430 void HcalDigiMonitor::analyze(edm::Event const&e, edm::EventSetup const&s)
00431 {
00432   if (!IsAllowedCalibType()) return;
00433   if (LumiInOrder(e.luminosityBlock())==false) return;
00434 
00435   // Get HLT trigger information for HF timing study
00436   passedMinBiasHLT_=false;
00437 
00439   // check if detectors whether they were ON
00440   edm::Handle<DcsStatusCollection> dcsStatus;
00441   e.getByLabel("scalersRawToDigi", dcsStatus);
00442   
00443   if (dcsStatus.isValid() && dcsStatus->size() != 0) 
00444     {      
00445       if ((*dcsStatus)[0].ready(DcsStatus::HBHEa) &&
00446           (*dcsStatus)[0].ready(DcsStatus::HBHEb) &&   
00447           (*dcsStatus)[0].ready(DcsStatus::HBHEc))
00448         {       
00449           hbhedcsON = true;
00450           if (debug_) std::cout << "hbhe on" << std::endl;
00451         } 
00452       else hbhedcsON = false;
00453 
00454       if ((*dcsStatus)[0].ready(DcsStatus::HF))
00455         {
00456           hfdcsON = true;
00457           if (debug_) std::cout << "hf on" << std::endl;
00458         } 
00459       else hfdcsON = false;
00460     }
00462 
00463   edm::Handle<edm::TriggerResults> hltRes;
00464   if (!(e.getByLabel(hltresultsLabel_,hltRes)))
00465     {
00466       if (debug_>0) edm::LogWarning("HcalDigiMonitor")<<" Could not get HLT results with tag "<<hltresultsLabel_<<std::endl;
00467     }
00468   else
00469     {
00470       const edm::TriggerNames & triggerNames = e.triggerNames(*hltRes);
00471       const unsigned int nTrig(triggerNames.size());
00472       for (unsigned int i=0;i<nTrig;++i){
00473           // repeat for minbias triggers
00474           for (unsigned int k=0;k<MinBiasHLTBits_.size();++k)
00475             {
00476               // if (triggerNames.triggerName(i)==MinBiasHLTBits_[k] && hltRes->accept(i))
00477               if (triggerNames.triggerName(i).find(MinBiasHLTBits_[k])!=std::string::npos && hltRes->accept(i))
00478                 { 
00479                   passedMinBiasHLT_=true;
00480                   break;
00481                 }
00482             }
00483         }
00484     } //else
00485   
00486   // Now get collections we need
00487   HT_HFP_=0;
00488   HT_HFM_=0;
00489   //  bool rechitsFound=false;
00490   edm::Handle<HFRecHitCollection> hf_rechit;
00491   if (e.getByLabel(hfRechitLabel_,hf_rechit))
00492     {
00493       //      rechitsFound=true;
00494       for (HFRecHitCollection::const_iterator HF=hf_rechit->begin();HF!=hf_rechit->end();++HF)
00495         {
00496           float en=HF->energy();
00497           int ieta=HF->id().ieta();
00498           // ieta for HF starts at 29, so subtract away 29 when computing fEta
00499           double fEta=fabs(0.5*(theHFEtaBounds[abs(ieta)-28]+theHFEtaBounds[abs(ieta)-29]));
00500           ieta>0 ?  HT_HFP_+=en/cosh(fEta) : HT_HFM_+=en/cosh(fEta);
00501         }
00502     }
00503   else
00504     {
00505       // if no rechits found, form above-threshold plots based only on digi comparison to ADC threshold 
00506       HT_HFP_=999;
00507       HT_HFM_=999;
00508     }
00509 
00510   // try to get digis
00511   edm::Handle<HBHEDigiCollection> hbhe_digi;
00512   edm::Handle<HODigiCollection> ho_digi;
00513   edm::Handle<HFDigiCollection> hf_digi;
00514 
00515   if (!(e.getByLabel(digiLabel_,hbhe_digi)))
00516     {
00517       edm::LogWarning("HcalDigiMonitor")<< digiLabel_<<" hbhe_digi not available";
00518       return;
00519     }
00520   
00521   if (!(e.getByLabel(digiLabel_,hf_digi)))
00522     {
00523       edm::LogWarning("HcalDigiMonitor")<< digiLabel_<<" hf_digi not available";
00524       return;
00525     }
00526   if (!(e.getByLabel(digiLabel_,ho_digi)))
00527     {
00528       edm::LogWarning("HcalDigiMonitor")<< digiLabel_<<" ho_digi not available";
00529       return;
00530     }
00531   edm::Handle<HcalUnpackerReport> report;  
00532   if (!(e.getByLabel(digiLabel_,report)))
00533     {
00534       edm::LogWarning("HcalDigiMonitor")<< digiLabel_<<" unpacker report not available";
00535       return;
00536     }
00537 
00538   // all objects grabbed; event is good
00539   if (debug_>1) std::cout <<"\t<HcalDigiMonitor::analyze>  Processing good event! event # = "<<ievt_<<std::endl;
00540 
00541   HcalBaseDQMonitor::analyze(e,s); // base class increments ievt_, etc. counters
00542 
00543   // Digi collection was grabbed successfully; process the Event
00544   processEvent(*hbhe_digi, *ho_digi, *hf_digi, *conditions_,
00545                *report, e.orbitNumber(),e.bunchCrossing());
00546   
00547 } //void HcalDigiMonitor::analyze(...)
00548 
00549 void HcalDigiMonitor::processEvent(const HBHEDigiCollection& hbhe,
00550                                    const HODigiCollection& ho,
00551                                    const HFDigiCollection& hf,
00552                                    const HcalDbService& cond,
00553                                    const HcalUnpackerReport& report,
00554                                    int orN, int bcN)
00555 { 
00556   if(!dbe_) 
00557     { 
00558       if(debug_) 
00559         std::cout <<"HcalDigiMonitor::processEvent   DQMStore not instantiated!!!"<<std::endl; 
00560       return; 
00561     }
00562 
00563   // Skip events in which minimal good digis found -- still getting some strange (calib?) events through DQM
00564 
00565   DigiUnpackerErrorCount->Fill(report.badQualityDigis());
00566   
00567   unsigned int allgooddigis= hbhe.size()+ho.size()+hf.size();
00568   // bad threshold:  ignore events in which bad outnumber good by more than 100:1
00569   // (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)
00570   if ((allgooddigis==0) ||
00571       (1.*report.badQualityDigis()>100*allgooddigis))
00572     {
00573       h_valid_digis->Fill(1);
00574       if (bcN>-1)
00575         h_invalid_bcn->Fill(bcN);
00576       if (orN>-1)
00577         h_invalid_orbitnumMod103->Fill(orN%103);
00578 
00579       return;
00580     }
00581   
00582   h_valid_digis->Fill(0);
00583 
00584   // hbHists.count_bad=0;
00585   // hbHists.count_good=0;
00586   // heHists.count_bad=0;
00587   // heHists.count_good=0;
00588   // hoHists.count_bad=0;
00589   // hoHists.count_good=0;
00590   // hfHists.count_bad=0;
00591   // hfHists.count_good=0;
00592 
00593   // int HO0bad=0;
00594   // int HO12bad=0;
00595   // int HFlumibad=0;
00596 
00597   // Check unpacker report for bad digis
00598 
00599   typedef std::vector<DetId> DetIdVector;
00600 
00601   for ( DetIdVector::const_iterator baddigi_iter=report.bad_quality_begin(); 
00602         baddigi_iter != report.bad_quality_end();
00603         ++baddigi_iter)
00604     {
00605       HcalDetId id(baddigi_iter->rawId());
00606       int rDepth = id.depth();
00607       int rPhi   = id.iphi();
00608       int rEta   = id.ieta();
00609       int binEta = CalcEtaBin(id.subdet(), rEta, rDepth); // why is this here?
00610       if (id.subdet()==HcalBarrel) ++hbHists.count_bad;
00611       else if (id.subdet()==HcalEndcap) ++heHists.count_bad;
00612       else if (id.subdet()==HcalForward) 
00613         {
00614           ++hfHists.count_bad;
00615           if (rDepth==1 && (abs(rEta)==33 || abs(rEta)==34)) ++HFlumibad;
00616           else if (rDepth==2 && (abs(rEta)==35 || abs(rEta)==36)) ++HFlumibad;
00617         }
00618       else if (id.subdet()==HcalOuter) 
00619         {
00620           // Mark HORing+/-2 channels as present, HO/YB+/-2 has HV off (at 100V).
00621           if (excludeHORing2_==true && rDepth==4)
00622             if (abs(rEta)>=11 && abs(rEta)<=15 && !isSiPM(rEta,rPhi,rDepth)) continue;
00623 
00624           if (excludeHO1P02_==true)
00625             if( (rEta>4 && rEta<10) && (rPhi<=10 || rPhi>70) ) continue;
00626           
00627           if (KnownBadCells_.find(id)!=KnownBadCells_.end()) continue;
00628 
00629           ++hoHists.count_bad;
00630           if (abs(rEta)<5) ++HO0bad;
00631           else ++HO12bad;
00632         }
00633       else 
00634         continue; // skip anything that isn't HB, HE, HO, HF
00635       // extra protection against nonsensical values -- prevents occasional crashes
00636       if (binEta < 85 && binEta >= 0 
00637           && (rPhi-1) >= 0 && (rPhi-1)<72 
00638           && (rDepth-1) >= 0 && (rDepth-1)<4)
00639         {
00640           ++badunpackerreport[binEta][rPhi-1][rDepth-1];
00641           ++baddigis[binEta][rPhi-1][rDepth-1];  
00642         }
00643     }
00644 
00646 
00647   int firsthbcap=-1;
00648   int firsthecap=-1;
00649   int firsthocap=-1;
00650   int firsthfcap=-1;
00651 
00652   for (HBHEDigiCollection::const_iterator j=hbhe.begin(); j!=hbhe.end(); ++j)
00653     {
00654         const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00655           
00656         if (digi.id().subdet()==HcalBarrel)
00657           {
00658             if (!HBpresent_) continue;
00659             if (KnownBadCells_.find(digi.id())!=KnownBadCells_.end()) continue;
00660           
00661             process_Digi(digi, hbHists, firsthbcap);
00662           }
00663         else if (digi.id().subdet()==HcalEndcap)
00664           {
00665             if (!HEpresent_) continue;
00666             process_Digi(digi, heHists,firsthecap);
00667           }     
00668     }
00669 
00670   // // Fill good digis vs lumi block; also fill bad errors?
00671   HBocc_vs_LB->Fill(currentLS,hbHists.count_good);
00672   HEocc_vs_LB->Fill(currentLS,heHists.count_good);
00673 
00674   // Calculate number of bad quality cells and bad quality fraction
00675   if (HBpresent_ && (hbHists.count_good>0 || hbHists.count_bad>0))
00676     {
00677       int counter=hbHists.count_bad;
00678       if (counter<DIGI_SUBDET_NUM)
00679         ++hbHists.count_BQ[counter];
00680       float counter2 = (1.*hbHists.count_bad)/(hbHists.count_bad+hbHists.count_good)*(DIGI_BQ_FRAC_NBINS-1);
00681       if (counter2<DIGI_SUBDET_NUM) ++hbHists.count_BQFrac[(int)counter2];
00682     }
00683 
00684   if (HEpresent_ && (heHists.count_good>0 || heHists.count_bad>0))
00685     {
00686       int counter=heHists.count_bad;
00687       if (counter<DIGI_SUBDET_NUM)
00688         ++heHists.count_BQ[counter];
00689       float counter2 = (1.*heHists.count_bad)/(heHists.count_bad+heHists.count_good)*(DIGI_BQ_FRAC_NBINS-1);
00690       if (counter2<DIGI_SUBDET_NUM) ++heHists.count_BQFrac[int(counter2)];
00691     }
00692 
00694   if (HOpresent_)
00695     {
00696       for (HODigiCollection::const_iterator j=ho.begin(); j!=ho.end(); ++j)
00697         {
00698           const HODataFrame digi = (const HODataFrame)(*j);
00699           // Mark HORing+/-2 channels as present, HO/YB+/-2 has HV off (at 100V).
00700           if (excludeHORing2_==true && digi.id().depth()==4)
00701             if (abs(digi.id().ieta())>=11 && abs(digi.id().ieta())<=15 && 
00702                 !isSiPM(digi.id().ieta(),digi.id().iphi(),digi.id().depth())) continue;
00703           
00704           if (excludeHO1P02_==true)       
00705             if( (digi.id().ieta()>4 && digi.id().ieta()<10) 
00706                 && (digi.id().iphi()<=10 || digi.id().iphi()>70) ) continue;
00707 
00708           if (KnownBadCells_.find(digi.id())!=KnownBadCells_.end()) continue;
00709           
00710           process_Digi(digi, hoHists, firsthocap);
00711         } // for (HODigiCollection)
00712       
00713       if (hoHists.count_bad>0 || hoHists.count_good>0)
00714         {
00715           int counter=hoHists.count_bad;
00716           if (counter<DIGI_SUBDET_NUM)
00717             ++hoHists.count_BQ[counter];
00718 
00719           float counter2 = (1.*hoHists.count_bad)/(hoHists.count_bad+hoHists.count_good)*(DIGI_BQ_FRAC_NBINS-1);
00720           if (counter2<DIGI_SUBDET_NUM) ++hoHists.count_BQFrac[int(counter2)];
00721         }
00722       HOocc_vs_LB->Fill(currentLS,hoHists.count_good);
00723     } // if (HOpresent_)
00724   
00726   if (HFpresent_)
00727     {
00728       for (HFDigiCollection::const_iterator j=hf.begin(); j!=hf.end(); ++j)
00729         {
00730           const HFDataFrame digi = (const HFDataFrame)(*j);
00731           process_Digi(digi, hfHists, firsthfcap);
00732         } // for (HFDigiCollection)
00733 
00734       if (hfHists.count_bad>0 || hfHists.count_good>0)
00735         {
00736           int counter=hfHists.count_bad;
00737           if (counter<DIGI_SUBDET_NUM)
00738             ++hfHists.count_BQ[counter];
00739           float counter2 = (1.*hfHists.count_bad)/(hfHists.count_bad+hfHists.count_good)*(DIGI_BQ_FRAC_NBINS-1);
00740           if (counter2<DIGI_SUBDET_NUM) ++hfHists.count_BQFrac[int(counter2)];
00741         }
00742       HFocc_vs_LB->Fill(currentLS,hfHists.count_good);
00743     } // if (HFpresent_)
00744   
00745   // This only counts digis that are present but bad somehow; it does not count digis that are missing
00746   int count_good=hbHists.count_good+heHists.count_good+hoHists.count_good+hfHists.count_good;
00747   int count_bad=hbHists.count_bad+heHists.count_bad+hoHists.count_bad+hfHists.count_bad;
00748 
00749   if (count_good<DIGI_NUM)
00750     ++diginum[count_good];
00751 
00752   // Fill bad quality histograms
00753   DigiUnpackerErrorFrac->Fill(1.*report.badQualityDigis()/(report.badQualityDigis()+count_good));
00754   DigiBQ->Fill(count_bad);
00755   if (count_bad>0 || count_good>0)
00756     DigiBQFrac->Fill(1.*count_bad/(count_bad+count_good));
00757 
00758   // Call 'update' on all histograms so that they update in online DQM  
00759   UpdateHists(hbHists);
00760   UpdateHists(heHists);
00761   UpdateHists(hoHists);
00762   UpdateHists(hfHists);
00763 
00764   // Now update global (non-subdetector-specific) histograms
00765   DigiNum->update();
00766   DigiErrorVME->update();
00767   DigiErrorSpigot->update();
00768   DigiBQ->update();
00769   DigiBQFrac->update();
00770   DigiUnpackerErrorCount->update();
00771   DigiUnpackerErrorFrac->update();
00772 
00773   // Update eta-phi hists
00774   for (unsigned int zz=0;zz<DigiErrorOccupancyByDepth.depth.size();++zz)
00775       DigiErrorOccupancyByDepth.depth[zz]->update();
00776   for (unsigned int zz=0;zz<DigiErrorsByDepth.depth.size();++zz)
00777       DigiErrorsByDepth.depth[zz]->update();
00778   for (unsigned int zz=0;zz<DigiErrorsBadCapID.depth.size();++zz)
00779       DigiErrorsBadCapID.depth[zz]->update();
00780   for (unsigned int zz=0;zz<DigiErrorsDVErr.depth.size();++zz)
00781       DigiErrorsDVErr.depth[zz]->update();
00782   for (unsigned int zz=0;zz<DigiErrorsBadDigiSize.depth.size();++zz)
00783       DigiErrorsBadDigiSize.depth[zz]->update();
00784   for (unsigned int zz=0;zz<DigiErrorsBadADCSum.depth.size();++zz)
00785       DigiErrorsBadADCSum.depth[zz]->update();
00786   for (unsigned int zz=0;zz<DigiErrorsUnpacker.depth.size();++zz)
00787       DigiErrorsUnpacker.depth[zz]->update();
00788   for (unsigned int zz=0;zz<DigiErrorsBadFibBCNOff.depth.size();++zz)
00789     DigiErrorsBadFibBCNOff.depth[zz]->update();
00790 
00791   DigiOccupancyEta->update();
00792   DigiOccupancyPhi->update();
00793   DigiOccupancyVME->update();
00794   DigiOccupancySpigot->update();
00795   DigiSize->update();
00796 
00797   // // Fill problems vs. lumi block plots
00798   // ProblemsVsLB->Fill(currentLS,count_bad);
00799   // ProblemsVsLB_HB->Fill(currentLS,hbHists.count_bad);
00800   // ProblemsVsLB_HE->Fill(currentLS,heHists.count_bad);
00801   // ProblemsVsLB_HO->Fill(currentLS,hoHists.count_bad);
00802   // ProblemsVsLB_HF->Fill(currentLS,hfHists.count_bad);
00803   // ProblemsVsLB_HBHEHF->Fill(currentLS,hbHists.count_bad+heHists.count_bad+hfHists.count_bad);
00804 
00805   // // Fill the number of problem digis in each channel
00806   // ProblemsCurrentLB->Fill(-1,-1,1);  // event counter
00807   // ProblemsCurrentLB->Fill(0,0,hbHists.count_bad);
00808   // ProblemsCurrentLB->Fill(1,0,heHists.count_bad);
00809   // ProblemsCurrentLB->Fill(2,0,hoHists.count_bad);
00810   // ProblemsCurrentLB->Fill(3,0,hfHists.count_bad);
00811   // ProblemsCurrentLB->Fill(4,0,HO0bad);
00812   // ProblemsCurrentLB->Fill(5,0,HO12bad);
00813   // ProblemsCurrentLB->Fill(6,0,HFlumibad);
00814 
00815   // Call fill method every checkNevents
00816   //fill_Nevents();
00817   
00818   return;
00819 } // void HcalDigiMonitor::processEvent(...)
00820 
00821 
00822 
00823 template <class DIGI>
00824 int HcalDigiMonitor::process_Digi(DIGI& digi, DigiHists& h, int& firstcap)
00825 {
00826   int err=0x0;
00827   bool bitUp = false;
00828   int ADCcount=0;
00829 
00830   int shapeThresh=0;
00831   
00832   int mindigisize=1;
00833   int maxdigisize=10;
00834 
00835   if (digi.id().subdet()==HcalBarrel)
00836     {
00837       shapeThresh=shapeThreshHB_;
00838       mindigisize=mindigisizeHBHE_;
00839       maxdigisize=maxdigisizeHBHE_;
00840     }
00841   else if (digi.id().subdet()==HcalEndcap)
00842     {
00843       shapeThresh=shapeThreshHE_;
00844       mindigisize=mindigisizeHBHE_;
00845       maxdigisize=maxdigisizeHBHE_;
00846     }
00847   else if (digi.id().subdet()==HcalOuter)
00848     {
00849       shapeThresh=shapeThreshHO_;
00850       mindigisize=mindigisizeHO_;
00851       maxdigisize=maxdigisizeHO_;
00852     }
00853   else if (digi.id().subdet()==HcalForward)
00854     {
00855       shapeThresh=shapeThreshHF_;
00856       mindigisize=mindigisizeHF_;
00857       maxdigisize=maxdigisizeHF_;
00858     }
00859   int iEta = digi.id().ieta();
00860   int iPhi = digi.id().iphi();
00861   int iDepth = digi.id().depth();
00862   int calcEta = CalcEtaBin(digi.id().subdet(),iEta,iDepth);
00863 
00864   // Check that digi size is correct
00865   if (digi.size()<mindigisize || digi.size()>maxdigisize)
00866     {
00867       if (digi_checkdigisize_) err|=0x1;
00868       ++baddigisize[calcEta][iPhi-1][iDepth-1];
00869     }
00870   // Check digi size; if > 20, increment highest bin of digisize array
00871   if (digi.size()<20)
00872     ++digisize[static_cast<int>(digi.size())][digi.id().subdet()-1];
00873   else
00874     ++digisize[19][digi.id().subdet()-1];
00875 
00876   // loop over time slices of digi to check capID and errors
00877   ++h.count_presample[digi.presamples()];
00878 
00879   // Check CapID rotation
00880   if (firstcap==-1) firstcap = digi.sample(0).capid();
00881   int capdif = digi.sample(0).capid() - firstcap;
00882   //capdif = capdif%3 - capdif/3; // unnecessary?
00883   // capdif should run from -3 to +3
00884   if (capdif >-4 && capdif<4)
00885     ++h.capIDdiff[capdif+3];
00886   else
00887       ++h.capIDdiff[7];
00888 
00889   int last=-1;
00890   
00891   int offset = digi.fiberIdleOffset();
00892 
00893   // Only count BCN offset errors if ExpectedOrbitMessage Time is >-1
00894   // For offline (and thus cfg default), this won't be checked, since
00895   // we can't keep up to date with changes.
00896   if (offset != -1000 && DigiMonitor_ExpectedOrbitMessageTime_>-1) 
00897     {
00898       // increment counters only for non-zero offsets?
00899       ++h.fibbcnoff[offset + 7];
00900       if (offset != 0)
00901         {
00902           ++badFibBCNOff[calcEta][iPhi-1][iDepth-1];
00903           if (shutOffOrbitTest_ == false) err |= 0xF; // not an error if test turned off
00904         }
00905     }
00906 
00907   int tssum=0;
00908 
00909   bool digi_error=false;
00910 
00911   const int DigiSize=digi.size();
00912   for (int i=0;i<10;++i) pedSubtractedADC_[i]=0;
00913   const int pedSubADCsize=sizeof(pedSubtractedADC_)/sizeof(double);
00914 
00915   std::map<HcalDetId, std::vector<double> >::iterator  foundID = PedestalsByCapId_.find(digi.id());
00916   for (int i=0;i<DigiSize;++i)
00917     {
00918       int thisCapid = digi.sample(i).capid();
00919       if (thisCapid>=0 && thisCapid<4) ++h.capid[thisCapid];
00920 
00921       if (makeDiagnostics_)
00922         {
00923           if(bitUpset(last,thisCapid)) bitUp=true; // checking capID rotation
00924           last = thisCapid;
00925           // Check for digi error bits
00926           if (digi_checkdverr_)
00927             {
00928               if(digi.sample(i).er()) err=(err|0x2);
00929               if(!digi.sample(i).dv()) err=(err|0x2);
00930             }
00931           if ((digi_error==false) && (digi.sample(i).er() || !digi.sample(i).dv()))
00932             {
00933               ++digierrorsdverr[calcEta][iPhi-1][iDepth-1];
00934               digi_error=true; // only count 1 error per digi in this plot
00935             }
00936           ++h.dverr[static_cast<int>(2*digi.sample(i).er()+digi.sample(i).dv())];
00937         } // if (makeDiagnostics_)
00938 
00939       h.count_shape[i]+=digi.sample(i).adc();
00940       
00941       // Calculate ADC sum of adjacent samples -- still necessary?
00942       if (i==digi.size()-1) continue;
00943       tssum= digi.sample(i).adc()+digi.sample(i+1).adc();
00944       if (tssum<50 && tssum>=0)
00945         {
00946           if (iEta>0)
00947             ++h.tssumplus[tssum][i];
00948           else
00949             ++h.tssumminus[tssum][i];
00950         }
00951 
00952       if (digi.sample(i).adc()<0) ++h.adc[0];
00953       else if (digi.sample(i).adc()<200) ++h.adc[digi.sample(i).adc()];
00954       else ++h.adc[199];
00955 
00956       if (i>=pedSubADCsize) continue; // don't exceed maximum array length when checking digis
00957       
00958       if (foundID!=PedestalsByCapId_.end())
00959         {
00960           pedSubtractedADC_[i]=digi.sample(i).adc()-(foundID->second)[thisCapid];
00961           ADCcount+=(int)(digi.sample(i).adc()-(foundID->second)[thisCapid]);
00962         }
00963       else
00964         {
00965           pedSubtractedADC_[i]=digi.sample(i).adc()-3;
00966           ADCcount+=digi.sample(i).adc()-3; // default pedestal subtraction of 3 ADC counts
00967         }
00968     } // for (int i=0;i<digi.size();++i)
00969   
00970   // capid error found
00971   if(bitUp) 
00972     {
00973       if (digi_checkcapid_) err=(err|0x4);
00974       ++badcapID[calcEta][iPhi-1][iDepth-1];
00975     }
00976 
00977   // These plots generally don't get filled, unless we turn off the suppression of bad digis
00978   if (err>0)
00979     {
00980       if(uniqcounter[calcEta][iPhi-1][iDepth-1]<1)  
00981         {
00982           ++h.count_bad;
00983           ++baddigis[calcEta][iPhi-1][iDepth-1];
00984           ++errorVME[static_cast<int>(2*(digi.elecId().htrSlot()+0.5*digi.elecId().htrTopBottom()))][static_cast<int>(digi.elecId().readoutVMECrateId())];
00985           ++errorSpigot[static_cast<int>(digi.elecId().spigot())][static_cast<int>(digi.elecId().dccid())];
00986         }
00987       uniqcounter[calcEta][iPhi-1][iDepth-1]++; 
00988 
00989       return err;
00990     }
00991 
00992   if (ADCcount<0) ADCcount=0;
00993   if (ADCcount<199)
00994     ++h.adcsum[ADCcount];
00995   else
00996     ++h.adcsum[199]; // effective overflow bin
00997 
00998   // require larger threshold to look at pulse shapes
00999   
01000   if (ADCcount>shapeThresh && passedMinBiasHLT_  && HT_HFP_>1 && HT_HFM_>1)
01001     {
01002       h.ThreshCount->Fill(0,1);
01003       if (digi.id().subdet()!=HcalOuter || isSiPM(iEta,iPhi, iDepth)==false)
01004         {
01005           for (int i=0;i<pedSubADCsize;++i)
01006             h.count_shapeThresh[i]+=pedSubtractedADC_[i];
01007         }
01008     }
01009 
01010   // occupancy plots are only filled for good histograms
01011   ++h.count_good;
01012   ++occupancyEtaPhi[calcEta][iPhi-1][iDepth-1];
01013   ++occupancyEta[iEta+41];
01014   ++occupancyPhi[iPhi-1];
01015   // htr Slots run from 0-20, incremented by 0.5 for top/bottom
01016   ++occupancyVME[static_cast<int>(2*(digi.elecId().htrSlot()+0.5*digi.elecId().htrTopBottom()))][static_cast<int>(digi.elecId().readoutVMECrateId())];
01017   ++occupancySpigot[static_cast<int>(digi.elecId().spigot())][static_cast<int>(digi.elecId().dccid())];
01018 
01019   // Pawel's code for HF timing checks -- run only in online mode for non-calib events
01020   if (digi.id().subdet()==HcalForward 
01021       && Online_  //only run online 
01022       && currenttype_==0  // require non-calibration event
01023       && passedMinBiasHLT_ // require min bias trigger
01024       )
01025     {
01026       int maxtime=-1;
01027       double maxenergy=-1, fullenergy=0;
01028       int digisize=digi.size();
01029       for (int ff=0;ff<digisize;++ff)
01030         {
01031           fullenergy+=digi.sample(ff).nominal_fC()-2.5;
01032           if (digi.sample(ff).nominal_fC()-2.5>maxenergy)
01033             {
01034               maxenergy=digi.sample(ff).nominal_fC()-2.5;
01035               maxtime=ff;
01036             }
01037         }
01038       
01039       if (maxtime>=2 && maxtime<=5 && maxenergy>20 && maxenergy<100)  // only look between time slices 2-5; anything else should be nonsense
01040         {
01041           for (int ff=0;ff<digisize;++ff){
01042             if(fullenergy>0){
01043               if(digi.id().ieta()>0)HFP_shape->Fill(ff,(digi.sample(ff).nominal_fC()-2.5)/fullenergy);
01044               if(digi.id().ieta()<0)HFM_shape->Fill(ff,(digi.sample(ff).nominal_fC()-2.5)/fullenergy);
01045             }
01046           }
01047 
01048           double time_den=0, time_num=0;
01049           // form weighted time sum
01050           int startslice=std::max(0,maxtime-1);
01051           int endslice=std::min(digisize-1,maxtime+1);
01052           for (int ss=startslice;ss<=endslice;++ss)
01053             {
01054               // subtract 'default' pedestal of 2.5 fC
01055               time_num+=ss*(digi.sample(ss).nominal_fC()-2.5);
01056               time_den+=digi.sample(ss).nominal_fC()-2.5;
01057             }
01058 
01059           int myiphi=iPhi;
01060           if (iDepth==2) ++myiphi;
01061           if (HFtiming_etaProfile!=0 && time_den!=0)
01062             HFtiming_etaProfile->Fill(iEta,time_num/time_den);
01063           if (HFtiming_totaltime2D!=0 && time_den!=0)
01064             HFtiming_totaltime2D->Fill(iEta,myiphi,time_num/time_den);
01065           if (HFtiming_occupancy2D!=0 && time_den!=0)
01066             HFtiming_occupancy2D->Fill(iEta,myiphi,1);
01067         } //maxtime>-1
01068     } // if HcalForward
01069 
01070   return err;
01071 } // template <class DIGI> int HcalDigiMonitor::process_Digi
01072 
01073 void HcalDigiMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
01074                                              const edm::EventSetup& c) 
01075 {
01076   HcalBaseDQMonitor::beginLuminosityBlock(lumiSeg,c);
01077   ProblemsCurrentLB->Reset();
01078 }
01079 
01080 void HcalDigiMonitor::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
01081                                            const edm::EventSetup& c)
01082 {
01083   if (LumiInOrder(lumiSeg.luminosityBlock())==false) return;
01084 
01085   // Reset current LS histogram
01086   if (ProblemsCurrentLB)
01087     ProblemsCurrentLB->Reset();
01088   
01089   ProblemDigisInLastNLB_HBHEHF_alarm->Reset();
01090 
01091   //increase the number of LS counting, for alarmer. Only make alarms for HBHE
01092   if(hbhedcsON == true && hfdcsON == true && HBpresent_ == 1 && HEpresent_ == 1 && HFpresent_ == 1)
01093     ++alarmer_counter_; 
01094   else 
01095     alarmer_counter_ = 0;
01096 
01097   fill_Nevents();
01098 
01099   zeroCounters(); // reset counters of good/bad digis
01100  
01101   return;
01102 }
01103 
01104 void HcalDigiMonitor::fill_Nevents()
01105 {
01106   if (debug_>0)
01107     std::cout <<"<HcalDigiMonitor> Calling fill_Nevents for event  "<<tevt_<< " (processed events = "<<ievt_<<")"<<std::endl;
01108   int iPhi, iEta, iDepth;
01109   //  bool valid=false;
01110 
01111   // Fill problems vs. lumi block plots
01112   ProblemsVsLB->Fill(currentLS,hbHists.count_bad+heHists.count_bad+hoHists.count_bad+hfHists.count_bad);
01113   ProblemsVsLB_HB->Fill(currentLS,hbHists.count_bad);
01114   ProblemsVsLB_HE->Fill(currentLS,heHists.count_bad);
01115   ProblemsVsLB_HO->Fill(currentLS,hoHists.count_bad);
01116   ProblemsVsLB_HF->Fill(currentLS,hfHists.count_bad);
01117   ProblemsVsLB_HBHEHF->Fill(currentLS,hbHists.count_bad+heHists.count_bad+hfHists.count_bad);
01118 
01119   if( hbHists.count_bad+heHists.count_bad+hfHists.count_bad < 50 )
01120     alarmer_counter_ = 0;
01121     
01122   if( alarmer_counter_ >= 5 )
01123     ProblemDigisInLastNLB_HBHEHF_alarm->Fill( std::min(int(hbHists.count_bad+heHists.count_bad+hfHists.count_bad), 99) );
01124 
01125   // Fill the number of problem digis in each channel
01126   if (ProblemsCurrentLB)
01127     {      
01128       ProblemsCurrentLB->Fill(-1,-1,1);  // event counter
01129       ProblemsCurrentLB->Fill(0,0,hbHists.count_bad);
01130       ProblemsCurrentLB->Fill(1,0,heHists.count_bad);
01131       ProblemsCurrentLB->Fill(2,0,hoHists.count_bad);
01132       ProblemsCurrentLB->Fill(3,0,hfHists.count_bad);
01133       ProblemsCurrentLB->Fill(4,0,HO0bad);
01134       ProblemsCurrentLB->Fill(5,0,HO12bad);
01135       ProblemsCurrentLB->Fill(6,0,HFlumibad);
01136     }
01137 
01138   // Fill plots of sums of adjacent digi samples
01139   for (int i=0;i<10;++i)
01140     {
01141       for (int j=0;j<50;++j)
01142         {
01143           if (hbHists.tssumplus[j][i]>0) hbHists.TS_sum_plus[i]->Fill(j, hbHists.tssumplus[j][i]);
01144           if (hbHists.tssumminus[j][i]>0) hbHists.TS_sum_minus[i]->Fill(j, hbHists.tssumminus[j][i]);     
01145           if (heHists.tssumplus[j][i]>0) heHists.TS_sum_plus[i]->Fill(j, heHists.tssumplus[j][i]); 
01146           if (heHists.tssumminus[j][i]>0) heHists.TS_sum_minus[i]->Fill(j, heHists.tssumminus[j][i]); 
01147           if (hoHists.tssumplus[j][i]>0) hoHists.TS_sum_plus[i]->Fill(j, hoHists.tssumplus[j][i]); 
01148           if (hoHists.tssumminus[j][i]>0) hoHists.TS_sum_minus[i]->Fill(j, hoHists.tssumminus[j][i]); 
01149           if (hfHists.tssumplus[j][i]>0) hfHists.TS_sum_plus[i]->Fill(j, hfHists.tssumplus[j][i]); 
01150           if (hfHists.tssumminus[j][i]>0) hfHists.TS_sum_minus[i]->Fill(j, hfHists.tssumminus[j][i]); 
01151         }
01152     } // for (int i=0;i<10;++i)
01153 
01154   for (int i=0;i<DIGI_NUM;++i)
01155     {
01156       if (diginum[i]>0) DigiNum->Fill(i, diginum[i]);
01157       if (i>=DIGI_SUBDET_NUM) continue;
01158 
01159       if (hbHists.count_BQ[i]>0) hbHists.BQ->Fill(i, hbHists.count_BQ[i]);
01160       if (heHists.count_BQ[i]>0) heHists.BQ->Fill(i, heHists.count_BQ[i]);
01161       if (hoHists.count_BQ[i]>0) hoHists.BQ->Fill(i, hoHists.count_BQ[i]);
01162       if (hfHists.count_BQ[i]>0) hfHists.BQ->Fill(i, hfHists.count_BQ[i]);
01163     }//for int i=0;i<DIGI_NUM;++i)
01164 
01165 
01166   // Fill data-valid/error plots and capid plots
01167   for (int i=0;i<4;++i)
01168     {
01169       if (hbHists.dverr[i]>0) hbHists.DVerr->Fill(i, hbHists.dverr[i]);
01170       if (heHists.dverr[i]>0) heHists.DVerr->Fill(i, heHists.dverr[i]);
01171       if (hoHists.dverr[i]>0) hoHists.DVerr->Fill(i, hoHists.dverr[i]);
01172       if (hfHists.dverr[i]>0) hfHists.DVerr->Fill(i, hfHists.dverr[i]);
01173 
01174       if (hbHists.capid[i]>0) hbHists.CapID->Fill(i, hbHists.capid[i]);
01175       if (heHists.capid[i]>0) heHists.CapID->Fill(i, heHists.capid[i]);
01176       if (hoHists.capid[i]>0) hoHists.CapID->Fill(i, hoHists.capid[i]);
01177       if (hfHists.capid[i]>0) hfHists.CapID->Fill(i, hfHists.capid[i]);
01178     }
01179 
01180   for (int i=0;i<200;++i)
01181     {
01182       if (hbHists.adc[i]>0) hbHists.ADC->Fill(i, hbHists.adc[i]);
01183 
01184       if (heHists.adc[i]>0) heHists.ADC->Fill(i, heHists.adc[i]);
01185       if (hoHists.adc[i]>0) hoHists.ADC->Fill(i, hoHists.adc[i]);
01186       if (hfHists.adc[i]>0) hfHists.ADC->Fill(i, hfHists.adc[i]);
01187 
01188       if (hbHists.adcsum[i]>0) hbHists.ADCsum->Fill(i, hbHists.adcsum[i]);
01189       if (heHists.adcsum[i]>0) heHists.ADCsum->Fill(i, heHists.adcsum[i]);
01190       if (hoHists.adcsum[i]>0) hoHists.ADCsum->Fill(i, hoHists.adcsum[i]);
01191       if (hfHists.adcsum[i]>0) hfHists.ADCsum->Fill(i, hfHists.adcsum[i]);
01192     }
01193 
01194   for (int i = 0; i < 15; ++i)
01195     {
01196       if (hbHists.fibbcnoff[i]>0) hbHists.fibBCNOff->Fill(i-7, 
01197                                                           hbHists.fibbcnoff[i]);
01198       if (heHists.fibbcnoff[i]>0) heHists.fibBCNOff->Fill(i-7, 
01199                                                           heHists.fibbcnoff[i]);
01200       if (hfHists.fibbcnoff[i]>0) hfHists.fibBCNOff->Fill(i-7, 
01201                                                           hfHists.fibbcnoff[i]);
01202       if (hoHists.fibbcnoff[i]>0) hoHists.fibBCNOff->Fill(i-7,
01203                                                           hoHists.fibbcnoff[i]);
01204     }
01205 
01206   // Fill plots of bad fraction of digis found
01207   for (int i=0;i<DIGI_BQ_FRAC_NBINS;++i)
01208     {
01209       if (DIGI_BQ_FRAC_NBINS==1) break;
01210       if (hbHists.count_BQFrac[i]>0) hbHists.BQFrac->Fill(1.*i/(DIGI_BQ_FRAC_NBINS-1), hbHists.count_BQFrac[i]);
01211       if (heHists.count_BQFrac[i]>0) heHists.BQFrac->Fill(1.*i/(DIGI_BQ_FRAC_NBINS-1), heHists.count_BQFrac[i]);
01212       if (hoHists.count_BQFrac[i]>0) 
01213         {
01214           hoHists.BQFrac->Fill(1.*i/(DIGI_BQ_FRAC_NBINS), hoHists.count_BQFrac[i]);
01215         }
01216       if (hfHists.count_BQFrac[i]>0) hfHists.BQFrac->Fill(1.*i/(DIGI_BQ_FRAC_NBINS-1), hfHists.count_BQFrac[i]);
01217     }//for (int i=0;i<DIGI_BQ_FRAC_NBINS;++i)
01218 
01219   // Fill presample plots
01220   for (int i=0;i<50;++i)
01221     {
01222       if (hbHists.count_presample[i]>0) hbHists.presample->Fill(i, hbHists.count_presample[i]);
01223       if (heHists.count_presample[i]>0) heHists.presample->Fill(i, heHists.count_presample[i]);
01224       if (hoHists.count_presample[i]>0) hoHists.presample->Fill(i, hoHists.count_presample[i]);
01225       if (hfHists.count_presample[i]>0) hfHists.presample->Fill(i, hfHists.count_presample[i]);
01226     } //for (int i=0;i<50;++i)
01227 
01228   // Fill shape plots
01229   for (int i=0;i<10;++i)
01230     {
01231       if (hbHists.count_shape[i]>0) hbHists.shape->Fill(i, hbHists.count_shape[i]);
01232       if (hbHists.count_shapeThresh[i]>0) hbHists.shapeThresh->Fill(i, hbHists.count_shapeThresh[i]);
01233       if (heHists.count_shape[i]>0) heHists.shape->Fill(i, heHists.count_shape[i]);
01234       if (heHists.count_shapeThresh[i]>0) heHists.shapeThresh->Fill(i, heHists.count_shapeThresh[i]);
01235       if (hoHists.count_shape[i]>0) hoHists.shape->Fill(i, hoHists.count_shape[i]);
01236       if (hoHists.count_shapeThresh[i]>0) hoHists.shapeThresh->Fill(i, hoHists.count_shapeThresh[i]);
01237       if (hfHists.count_shape[i]>0) hfHists.shape->Fill(i, hfHists.count_shape[i]);
01238       if (hfHists.count_shapeThresh[i]>0) hfHists.shapeThresh->Fill(i, hfHists.count_shapeThresh[i]);
01239     }//  for (int i=0;i<10;++i)
01240 
01241   // Fill capID difference plots
01242   for (int i=0;i<8;++i)
01243     {
01244       if (hbHists.capIDdiff[i]>0) hbHists.DigiFirstCapID->Fill(i, hbHists.capIDdiff[i]);
01245       if (heHists.capIDdiff[i]>0) heHists.DigiFirstCapID->Fill(i, heHists.capIDdiff[i]);
01246       if (hoHists.capIDdiff[i]>0) hoHists.DigiFirstCapID->Fill(i, hoHists.capIDdiff[i]);
01247       if (hfHists.capIDdiff[i]>0) hfHists.DigiFirstCapID->Fill(i, hfHists.capIDdiff[i]);
01248     }
01249 
01250   // Fill VME plots
01251   for (int i=0;i<40;++i)
01252     {
01253       for (int j=0;j<18;++j)
01254         {
01255           if (errorVME[i][j]>0) DigiErrorVME->Fill(i, j,errorVME[i][j]);
01256           if (occupancyVME[i][j]>0) DigiOccupancyVME->Fill(i, j,occupancyVME[i][j]);
01257         }
01258     } //for (int i=0;i<40;++i)
01259 
01260   // Fill SPIGOT plots
01261   for (int i=0;i<HcalDCCHeader::SPIGOT_COUNT;++i)
01262     {
01263       for (int j=0;j<36;++j)
01264         {
01265           if (errorSpigot[i][j]>0) DigiErrorSpigot->Fill(i, j,errorSpigot[i][j]);
01266           if (occupancySpigot[i][j]>0) DigiOccupancySpigot->Fill(i, j,occupancySpigot[i][j]);
01267         }
01268     } //for (int i=0;i<HcalDCCHeader::SPIGOT_COUNT;++i)
01269 
01270   // Loop over subdetectors
01271   for (int sub=0;sub<4;++sub)
01272     {
01273       for (int dsize=0;dsize<20;++dsize)
01274         {
01275           if (digisize[dsize][sub]>0)
01276             DigiSize->Fill(sub,dsize,digisize[dsize][sub]);
01277         }
01278     } // for (int sub=0;sub<4;++sub)
01279 
01280   // Loop over eta, phi, depth
01281   for (int d=0;d<4;++d)
01282     {
01283       iDepth=d+1;
01284       DigiErrorsByDepth.depth[d]->setBinContent(0,0,ievt_); // underflow bin contains event counter
01285       DigiOccupancyByDepth.depth[d]->setBinContent(0,0,ievt_);
01286       DigiErrorsBadDigiSize.depth[d]->setBinContent(0,0,ievt_);
01287       DigiErrorsUnpacker.depth[d]->setBinContent(0,0,ievt_);
01288       DigiErrorsBadFibBCNOff.depth[d]->setBinContent(0,0,ievt_);
01289 
01290       for (int phi=0;phi<72;++phi)
01291         {
01292           iPhi=phi+1;
01293           DigiOccupancyPhi->Fill(iPhi,occupancyPhi[phi]);
01294           for (int eta=0;eta<83;++eta)
01295             {
01296               // DigiOccupanyEta uses 'true' ieta (included the overlap at +/- 29)
01297               iEta=eta-41;
01298               if (phi==0)
01299                 DigiOccupancyEta->Fill(iEta,occupancyEta[eta]);
01300               //              valid=false;
01301         
01302               // HB
01303               if (validDetId(HcalBarrel, iEta, iPhi, iDepth))
01304                 {
01305                   //              valid=true;
01306                   if (HBpresent_)
01307                     {
01308                       int calcEta = CalcEtaBin(HcalBarrel,iEta,iDepth);
01309 
01310                       DigiOccupancyByDepth.depth[d]->Fill(iEta, iPhi,
01311                                                     occupancyEtaPhi[calcEta][phi][d]);
01312                       
01313                       if (makeDiagnostics_)
01314                         {
01315                           DigiErrorsBadCapID.depth[d]->Fill(iEta, iPhi,
01316                                                             badcapID[calcEta][phi][d]);
01317                           DigiErrorsDVErr.depth[d]->Fill(iEta, iPhi,
01318                                                          digierrorsdverr[calcEta][phi][d]);
01319                         }
01320                       DigiErrorsBadDigiSize.depth[d]->Fill(iEta, iPhi,
01321                                                            baddigisize[calcEta][phi][d]);
01322                       DigiErrorsBadFibBCNOff.depth[d]->Fill(iEta, iPhi,
01323                                                             badFibBCNOff[calcEta][phi][d]);
01324                       DigiErrorsUnpacker.depth[d]->Fill(iEta, iPhi,
01325                                                         badunpackerreport[calcEta][phi][d]);
01326                       DigiErrorsByDepth.depth[d]->Fill(iEta, iPhi,
01327                                                        baddigis[calcEta][phi][d]);
01328                       // Use this for testing purposes only
01329                       //DigiErrorsByDepth[d]->Fill(iEta, iPhi, ievt_);
01330                     } // if (HBpresent_)
01331                 } // validDetId(HB)
01332               // HE
01333               if (validDetId(HcalEndcap, iEta, iPhi, iDepth))
01334                 {
01335                   //              valid=true;
01336                   if (HEpresent_)
01337                     {
01338                       int calcEta = CalcEtaBin(HcalEndcap,iEta,iDepth);
01339 
01340                       DigiOccupancyByDepth.depth[d]->Fill(iEta, iPhi,
01341                                                     occupancyEtaPhi[calcEta][phi][d]);
01342                       
01343                       if (makeDiagnostics_)
01344                         {
01345                           DigiErrorsBadCapID.depth[d]->Fill(iEta, iPhi,
01346                                                             badcapID[calcEta][phi][d]);
01347                           DigiErrorsDVErr.depth[d]->Fill(iEta, iPhi,
01348                                                          digierrorsdverr[calcEta][phi][d]);
01349                         }
01350                       DigiErrorsBadDigiSize.depth[d]->Fill(iEta, iPhi,
01351                                                            baddigisize[calcEta][phi][d]);
01352                       DigiErrorsBadFibBCNOff.depth[d]->Fill(iEta, iPhi,
01353                                                             badFibBCNOff[calcEta][phi][d]);
01354                       DigiErrorsUnpacker.depth[d]->Fill(iEta, iPhi,
01355                                                         badunpackerreport[calcEta][phi][d]);
01356                       DigiErrorsByDepth.depth[d]->Fill(iEta, iPhi,
01357                                                        baddigis[calcEta][phi][d]);
01358                     } // if (HEpresent_)
01359                 } // valid HE found
01360               // HO
01361               if (validDetId(HcalOuter,iEta,iPhi,iDepth))
01362                 {
01363                   //              valid=true;
01364                   if (HOpresent_)
01365                     {
01366                       int calcEta = CalcEtaBin(HcalOuter,iEta,iDepth);
01367                       DigiOccupancyByDepth.depth[d]->Fill(iEta, iPhi,
01368                                                           occupancyEtaPhi[calcEta][phi][d]);
01369                       if (makeDiagnostics_)
01370                         {
01371                           DigiErrorsBadCapID.depth[d]->Fill(iEta, iPhi,
01372                                                             badcapID[calcEta][phi][d]);
01373                           DigiErrorsDVErr.depth[d]->Fill(iEta, iPhi,
01374                                                          digierrorsdverr[calcEta][phi][d]);
01375                         }
01376                       DigiErrorsBadDigiSize.depth[d]->Fill(iEta, iPhi,
01377                                                            baddigisize[calcEta][phi][d]);
01378                       DigiErrorsBadFibBCNOff.depth[d]->Fill(iEta, iPhi,
01379                                                             badFibBCNOff[calcEta][phi][d]);
01380                       DigiErrorsUnpacker.depth[d]->Fill(iEta, iPhi,
01381                                                         badunpackerreport[calcEta][phi][d]);
01382                       
01383                       DigiErrorsByDepth.depth[d]->Fill(iEta,iPhi,
01384                                                        baddigis[calcEta][phi][d]);
01385                     } // if (HOpresent_)
01386                 }//validDetId(HO)
01387               // HF
01388               if (validDetId(HcalForward,iEta,iPhi,iDepth))
01389                 {
01390                   //              valid=true;
01391                   if (HFpresent_)
01392                     {
01393                       int calcEta = CalcEtaBin(HcalForward,iEta,iDepth);
01394                       int zside = iEta/abs(iEta);
01395                       DigiOccupancyByDepth.depth[d]->Fill(iEta+zside, iPhi,
01396                                                     occupancyEtaPhi[calcEta][phi][d]);
01397                       
01398                       if (makeDiagnostics_)
01399                         {
01400                           DigiErrorsBadCapID.depth[d]->Fill(iEta+zside, iPhi,
01401                                                             badcapID[calcEta][phi][d]);
01402                           DigiErrorsDVErr.depth[d]->Fill(iEta+zside, iPhi,
01403                                                          digierrorsdverr[calcEta][phi][d]);
01404                         }
01405                       DigiErrorsBadDigiSize.depth[d]->Fill(iEta+zside, iPhi,
01406                                                      baddigisize[calcEta][phi][d]);
01407                       DigiErrorsBadFibBCNOff.depth[d]->Fill(iEta+zside, iPhi,
01408                                                          badFibBCNOff[calcEta][phi][d]);
01409                       DigiErrorsUnpacker.depth[d]->Fill(iEta+zside, iPhi,
01410                                                         badunpackerreport[calcEta][phi][d]);
01411                       DigiErrorsByDepth.depth[d]->Fill(iEta+zside, iPhi,
01412                                                        baddigis[calcEta][phi][d]);
01413                       
01414                     } // if (HFpresent_)
01415                 }
01416             } // for (int eta=0;...)
01417         } // for (int phi=0;...)
01418     } // for (int d=0;...)
01419 
01420   // Now fill all the unphysical cell values
01421   FillUnphysicalHEHFBins(DigiErrorsByDepth);
01422   if (makeDiagnostics_)
01423     {
01424       FillUnphysicalHEHFBins(DigiErrorsBadCapID);
01425       FillUnphysicalHEHFBins(DigiErrorsDVErr);
01426     }
01427   FillUnphysicalHEHFBins(DigiErrorsBadDigiSize);
01428   FillUnphysicalHEHFBins(DigiOccupancyByDepth);
01429   FillUnphysicalHEHFBins(DigiErrorsBadFibBCNOff);
01430   FillUnphysicalHEHFBins(DigiErrorsUnpacker);
01431 
01432   //  zeroCounters(); // reset counters of good/bad digis
01433  
01434   return;
01435 } // void HcalDigiMonitor::fill_Nevents()
01436 
01437 
01438 void HcalDigiMonitor::zeroCounters()
01439 {
01440   // Set all histogram counters back to 0
01441   // Call this after all every N evnets
01442 
01443   /******** Zero all counters *******/
01444   for (int d=0; d<DEPTHBINS; d++) {
01445     for (int eta=0; eta<ETABINS; eta++) {
01446       for (int phi=0; phi<PHIBINS; phi++){
01447         uniqcounter[eta][phi][d] = 0.0;
01448       }
01449     }
01450   }
01451 
01452   hbHists.count_bad=0;
01453   hbHists.count_good=0;
01454   heHists.count_bad=0;
01455   heHists.count_good=0;
01456   hoHists.count_bad=0;
01457   hoHists.count_good=0;
01458   hfHists.count_bad=0;
01459   hfHists.count_good=0;
01460 
01461   HO0bad=0;
01462   HO12bad=0;
01463   HFlumibad=0;
01464 
01465   for (int i=0;i<85;++i)
01466     {
01467       occupancyEta[i]=0;
01468       if (i<72)
01469         occupancyPhi[i]=0;
01470       for (int j=0;j<72;++j)
01471         {
01472           for (int k=0;k<4;++k)
01473             {
01474               baddigis[i][j][k]=0;
01475               badcapID[i][j][k]=0;
01476               baddigisize[i][j][k]=0;
01477               occupancyEtaPhi[i][j][k]=0;
01478               digierrorsdverr[i][j][k]=0;
01479               badFibBCNOff[i][j][k]=0;
01480               badunpackerreport[i][j][k]=0;
01481             }
01482         } // for (int j=0;j<72;++i)
01483     } // for (int i=0;i<85;++i)
01484 
01485   for (int i=0;i<40;++i)
01486     {
01487       for (int j=0;j<18;++j)
01488         {
01489           occupancyVME[i][j]=0;
01490           errorVME[i][j]=0;
01491         }
01492     }
01493 
01494   for (int i=0;i<HcalDCCHeader::SPIGOT_COUNT;++i)
01495     {
01496       for (int j=0;j<36;++j)
01497         {
01498           occupancySpigot[i][j]=0;
01499           errorSpigot[i][j]=0;
01500         }
01501     }
01502 
01503 
01504   for (int i=0;i<20;++i)
01505     {
01506       for (int j=0;j<4;++j)
01507         digisize[i][j]=0;
01508     }
01509 
01510   for (int i=0;i<DIGI_NUM;++i)
01511     {
01512       diginum[i]=0;
01513       
01514       // set all DigiHists counters to 0
01515       if (i<4)
01516         {
01517           hbHists.dverr[i]=0;
01518           heHists.dverr[i]=0;
01519           hoHists.dverr[i]=0;
01520           hfHists.dverr[i]=0;
01521           hbHists.capid[i]=0;
01522           heHists.capid[i]=0;
01523           hoHists.capid[i]=0;
01524           hfHists.capid[i]=0;
01525         }
01526       if (i<8)
01527         {
01528           hbHists.capIDdiff[i]=0;
01529           heHists.capIDdiff[i]=0;
01530           hoHists.capIDdiff[i]=0;
01531           hfHists.capIDdiff[i]=0;
01532         }
01533 
01534       if (i<10)
01535         {
01536           hbHists.count_shape[i]=0;
01537           heHists.count_shape[i]=0;
01538           hoHists.count_shape[i]=0;
01539           hfHists.count_shape[i]=0;
01540          
01541           hbHists.count_shapeThresh[i]=0;
01542           heHists.count_shapeThresh[i]=0;
01543           hoHists.count_shapeThresh[i]=0;
01544           hfHists.count_shapeThresh[i]=0;
01545         }
01546       if (i<50)
01547         {
01548           hbHists.count_presample[i]=0;
01549           heHists.count_presample[i]=0;
01550           hoHists.count_presample[i]=0;
01551           hfHists.count_presample[i]=0;
01552           for (int j=0;j<10;++j)
01553             {
01554               hbHists.tssumplus[i][j]=0;
01555               heHists.tssumplus[i][j]=0;
01556               hoHists.tssumplus[i][j]=0;
01557               hfHists.tssumplus[i][j]=0;
01558               hbHists.tssumminus[i][j]=0;
01559               heHists.tssumminus[i][j]=0;
01560               hoHists.tssumminus[i][j]=0;
01561               hfHists.tssumminus[i][j]=0;
01562             }
01563         }
01564 
01565       if (i<15)
01566         {
01567           hbHists.fibbcnoff[i]=0;
01568           heHists.fibbcnoff[i]=0;
01569           hoHists.fibbcnoff[i]=0;
01570           hfHists.fibbcnoff[i]=0;
01571         }
01572 
01573       if (i<200)
01574         {
01575           hbHists.adc[i]=0;
01576           heHists.adc[i]=0;
01577           hoHists.adc[i]=0;
01578           hfHists.adc[i]=0;
01579           hbHists.adcsum[i]=0;
01580           heHists.adcsum[i]=0;
01581           hoHists.adcsum[i]=0;
01582           hfHists.adcsum[i]=0;
01583         }
01584       if (i<DIGI_SUBDET_NUM)
01585         {
01586           hbHists.count_BQ[i]=0;
01587           heHists.count_BQ[i]=0;
01588           hoHists.count_BQ[i]=0;
01589           hfHists.count_BQ[i]=0;
01590         }
01591       if (i<DIGI_BQ_FRAC_NBINS)
01592         {
01593           hbHists.count_BQFrac[i]=0;
01594           heHists.count_BQFrac[i]=0;
01595           hoHists.count_BQFrac[i]=0;
01596           hfHists.count_BQFrac[i]=0;
01597         }
01598     } // for (int i=0;i<DIGI_NUM;++i)
01599 
01600 
01601   return;
01602 }
01603 
01604 void HcalDigiMonitor::UpdateHists(DigiHists& h)
01605 {
01606   // call update command for all histograms (should make them update when running in online DQM?)
01607   h.shape->update();
01608   h.shapeThresh->update();
01609   h.presample->update();
01610   h.BQ->update();
01611   h.BQFrac->update();
01612   h.DigiFirstCapID->update();
01613   h.DVerr->update();
01614   h.CapID->update();
01615   h.ADC->update();
01616   h.ADCsum->update();
01617   h.fibBCNOff->update();
01618 
01619   for (unsigned int i=0;i<h.TS_sum_plus.size();++i)
01620     h.TS_sum_plus[i]->update();
01621   for (unsigned int i=0;i<h.TS_sum_minus.size();++i)
01622     h.TS_sum_minus[i]->update();
01623 } //void HcalDigiMonitor::UpdateHists(DigiHists& h)
01624 
01625 
01626 void HcalDigiMonitor::reset()
01627 {
01628   // reset the temporary histograms
01629   zeroCounters();
01630   
01631   // then reset the MonitorElements
01632 
01633   ProblemDigisInLastNLB_HBHEHF_alarm->Reset();
01634   alarmer_counter_ = 0;
01635  
01636   hbhedcsON = true; hfdcsON = true;
01637 
01638   DigiErrorsByDepth.Reset();
01639   DigiErrorsBadCapID.Reset();
01640   DigiErrorsDVErr.Reset();
01641   DigiErrorsBadDigiSize.Reset();
01642   DigiErrorsBadADCSum.Reset();
01643   DigiErrorsUnpacker.Reset();
01644   DigiErrorsBadFibBCNOff.Reset();
01645   DigiOccupancyByDepth.Reset();
01646   DigiErrorOccupancyByDepth.Reset();
01647 
01648   DigiOccupancyEta->Reset();
01649   DigiOccupancyPhi->Reset();
01650   DigiOccupancyVME->Reset();
01651   DigiOccupancySpigot->Reset();
01652   DigiErrorVME->Reset();
01653   DigiErrorSpigot->Reset();
01654   
01655   DigiBQ->Reset();
01656   DigiBQFrac->Reset();
01657   DigiUnpackerErrorCount->Reset();
01658   DigiUnpackerErrorFrac->Reset();
01659 
01660   DigiNum->Reset();
01661 
01662   hbHists.shape->Reset();
01663   hbHists.shapeThresh->Reset();
01664   hbHists.presample->Reset();
01665   hbHists.BQ->Reset();
01666   hbHists.BQFrac->Reset();
01667   hbHists.DigiFirstCapID->Reset();
01668   hbHists.DVerr->Reset();
01669   hbHists.CapID->Reset();
01670   hbHists.ADC->Reset();
01671   hbHists.ADCsum->Reset();
01672   hbHists.fibBCNOff->Reset();
01673   for (unsigned int i=0;i<hbHists.TS_sum_plus.size();++i)
01674     hbHists.TS_sum_plus[i]->Reset();
01675   for (unsigned int i=0;i<hbHists.TS_sum_minus.size();++i)
01676     hbHists.TS_sum_minus[i]->Reset();
01677 
01678   heHists.shape->Reset();
01679   heHists.shapeThresh->Reset();
01680   heHists.presample->Reset();
01681   heHists.BQ->Reset();
01682   heHists.BQFrac->Reset();
01683   heHists.DigiFirstCapID->Reset();
01684   heHists.DVerr->Reset();
01685   heHists.CapID->Reset();
01686   heHists.ADC->Reset();
01687   heHists.ADCsum->Reset();
01688   heHists.fibBCNOff->Reset();
01689   for (unsigned int i=0;i<heHists.TS_sum_plus.size();++i)
01690     heHists.TS_sum_plus[i]->Reset();
01691   for (unsigned int i=0;i<heHists.TS_sum_minus.size();++i)
01692     heHists.TS_sum_minus[i]->Reset();
01693 
01694   hoHists.shape->Reset();
01695   hoHists.shapeThresh->Reset();
01696   hoHists.presample->Reset();
01697   hoHists.BQ->Reset();
01698   hoHists.BQFrac->Reset();
01699   hoHists.DigiFirstCapID->Reset();
01700   hoHists.DVerr->Reset();
01701   hoHists.CapID->Reset();
01702   hoHists.ADC->Reset();
01703   hoHists.ADCsum->Reset();
01704   hoHists.fibBCNOff->Reset();
01705   for (unsigned int i=0;i<hoHists.TS_sum_plus.size();++i)
01706     hoHists.TS_sum_plus[i]->Reset();
01707   for (unsigned int i=0;i<hoHists.TS_sum_minus.size();++i)
01708     hoHists.TS_sum_minus[i]->Reset();
01709 
01710   hfHists.shape->Reset();
01711   hfHists.shapeThresh->Reset();
01712   hfHists.presample->Reset();
01713   hfHists.BQ->Reset();
01714   hfHists.BQFrac->Reset();
01715   hfHists.DigiFirstCapID->Reset();
01716   hfHists.DVerr->Reset();
01717   hfHists.CapID->Reset();
01718   hfHists.ADC->Reset();
01719   hfHists.ADCsum->Reset();
01720   hfHists.fibBCNOff->Reset();
01721   for (unsigned int i=0;i<hfHists.TS_sum_plus.size();++i)
01722     hfHists.TS_sum_plus[i]->Reset();
01723   for (unsigned int i=0;i<hfHists.TS_sum_minus.size();++i)
01724     hfHists.TS_sum_minus[i]->Reset();
01725 
01726   return;
01727 }
01728 DEFINE_FWK_MODULE(HcalDigiMonitor);
01729