CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/DQM/HcalMonitorTasks/src/HcalDetDiagNoiseMonitor.cc

Go to the documentation of this file.
00001 #include "DQM/HcalMonitorTasks/interface/HcalDetDiagNoiseMonitor.h"
00002 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00003 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00004 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00005 #include "DataFormats/HcalDigi/interface/HcalCalibrationEventTypes.h"
00006 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
00007 #include "FWCore/Common/interface/TriggerNames.h"
00008 
00009 #include "TFile.h"
00010 #include "TTree.h"
00011 #include <TVector2.h>
00012 #include <TVector3.h>
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 
00016 #include "DataFormats/Math/interface/deltaR.h"
00017 #include "DataFormats/Common/interface/TriggerResults.h"
00018 
00019 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
00020 #include "DataFormats/L1Trigger/interface/L1EmParticleFwd.h"
00021 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
00022 #include "DataFormats/L1Trigger/interface/L1MuonParticleFwd.h"
00023 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00024 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00025 #include "DataFormats/L1Trigger/interface/L1EtMissParticle.h"
00026 #include "DataFormats/L1Trigger/interface/L1EtMissParticleFwd.h"
00027 #include "DataFormats/L1GlobalCaloTrigger/interface/L1GctHFRingEtSums.h"
00028 #include "DataFormats/L1GlobalCaloTrigger/interface/L1GctHFBitCounts.h"
00029 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTExtendedCand.h"
00030 #include "DataFormats/Candidate/interface/Candidate.h"
00031 #include "L1Trigger/RegionalCaloTrigger/interface/L1RCTProducer.h" 
00032 #include "DataFormats/L1CaloTrigger/interface/L1CaloCollections.h"
00033 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00034 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00035 #include "HLTrigger/HLTanalyzers/interface/JetUtil.h"
00036 #include "DataFormats/METReco/interface/CaloMET.h"
00037 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00038 // #include "DataFormats/L1Trigger/interface/L1ParticleMap.h"
00039 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00040 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
00041 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMapRecord.h"
00042 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMapFwd.h"
00043 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerObjectMap.h"
00044 //#include "DataFormats/L1GlobalTrigger/interface/L1GtLogicParser.h"
00045 #include "DataFormats/TrackReco/interface/Track.h"
00046 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00047 #include "DataFormats/TrackReco/interface/TrackBase.h"
00048 #include "RecoMET/METAlgorithms/interface/HcalNoiseRBXArray.h"
00049 #include "DataFormats/VertexReco/interface/Vertex.h"
00050 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00051 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00052 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00053 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00054 #include "DQMServices/Core/interface/DQMStore.h"
00055 #include "DQMServices/Core/interface/MonitorElement.h"
00056 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00057 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00058 
00059 // this is to retrieve HCAL LogicalMap
00060 #include "CondFormats/HcalObjects/interface/HcalLogicalMap.h"
00061 // to retrive trigger information (local runs only)
00062 #include "TBDataFormats/HcalTBObjects/interface/HcalTBTriggerData.h"
00063 
00064 
00065 #include <math.h>
00066 
00067 using namespace reco;
00068 
00070 constexpr float adc2fC[128]={-0.5,0.5,1.5,2.5,3.5,4.5,5.5,6.5,7.5,8.5,9.5, 10.5,11.5,12.5,
00071                    13.5,15.,17.,19.,21.,23.,25.,27.,29.5,32.5,35.5,38.5,42.,46.,50.,54.5,59.5,
00072                    64.5,59.5,64.5,69.5,74.5,79.5,84.5,89.5,94.5,99.5,104.5,109.5,114.5,119.5,
00073                    124.5,129.5,137.,147.,157.,167.,177.,187.,197.,209.5,224.5,239.5,254.5,272.,
00074                    292.,312.,334.5,359.5,384.5,359.5,384.5,409.5,434.5,459.5,484.5,509.5,534.5,
00075                    559.5,584.5,609.5,634.5,659.5,684.5,709.5,747.,797.,847.,897.,947.,997.,
00076                    1047.,1109.5,1184.5,1259.5,1334.5,1422.,1522.,1622.,1734.5,1859.5,1984.5,
00077                    1859.5,1984.5,2109.5,2234.5,2359.5,2484.5,2609.5,2734.5,2859.5,2984.5,
00078                    3109.5,3234.5,3359.5,3484.5,3609.5,3797.,4047.,4297.,4547.,4797.,5047.,
00079                    5297.,5609.5,5984.5,6359.5,6734.5,7172.,7672.,8172.,8734.5,9359.5,9984.5};
00081 constexpr char const *subdets[11]={"HBM","HBP","HEM","HEP","HO1M","HO0","HO1P","HO2M","HO2P","HFM","HFP"};
00082 constexpr char const *HB_RBX[36]={
00083 "HBM01","HBM02","HBM03","HBM04","HBM05","HBM06","HBM07","HBM08","HBM09","HBM10","HBM11","HBM12","HBM13","HBM14","HBM15","HBM16","HBM17","HBM18",
00084 "HBP01","HBP02","HBP03","HBP04","HBP05","HBP06","HBP07","HBP08","HBP09","HBP10","HBP11","HBP12","HBP13","HBP14","HBP15","HBP16","HBP17","HBP18"};
00085 constexpr char const *HE_RBX[36]={
00086 "HEM01","HEM02","HEM03","HEM04","HEM05","HEM06","HEM07","HEM08","HEM09","HEM10","HEM11","HEM12","HEM13","HEM14","HEM15","HEM16","HEM17","HEM18",
00087 "HEP01","HEP02","HEP03","HEP04","HEP05","HEP06","HEP07","HEP08","HEP09","HEP10","HEP11","HEP12","HEP13","HEP14","HEP15","HEP16","HEP17","HEP18"};
00088 constexpr char const *HO_RBX[36]={
00089 "HO2M02","HO2M04","HO2M06","HO2M08","HO2M10","HO2M12","HO1M02","HO1M04","HO1M06","HO1M08","HO1M10","HO1M12",
00090 "HO001","HO002","HO003","HO004","HO005","HO006","HO007","HO008","HO009","HO010","HO011","HO012",
00091 "HO1P02","HO1P04","HO1P06","HO1P08","HO1P10","HO1P12","HO2P02","HO2P04","HO2P06","HO2P08","HO2P10","HO2P12",
00092 };
00093 
00094 
00095 class HcalDetDiagNoiseRMData{
00096 public:
00097   HcalDetDiagNoiseRMData(){
00098      reset();
00099      reset_LS();
00100   };
00101   void reset(){
00102     n_th_hi=n_th_300=n_pix_1=n_pix_8=pix=n_pix=0;
00103   }  
00104   void reset_LS(){
00105     n_th_hi_LS=n_th_300_LS=0;
00106   }
00107   int    n_th_hi;
00108   int    n_th_300;
00109   int    n_pix_1;
00110   int    n_pix_8;
00111   int    pix;
00112   int    n_pix; 
00113   int    n_th_hi_LS;
00114   int    n_th_300_LS;
00115 };
00116 class HcalDetDiagNoiseRMEvent{
00117 public:
00118   HcalDetDiagNoiseRMEvent(){
00119     reset();
00120   }
00121   void reset(){
00122    n_pix_hi=n_pix_lo=energy=n_zero=0;
00123   }
00124   int n_pix_hi;
00125   int n_pix_lo;
00126   int n_zero;
00127   float energy;
00128 };
00129 
00130 class HcalDetDiagNoiseRMSummary{
00131 public:
00132   HcalDetDiagNoiseRMSummary(){ 
00133      reset();
00134   }
00135   void reset(){
00136      for(int i=0;i<HcalFrontEndId::maxRmIndex;i++) rm[i].reset(); 
00137   }
00138   void reset_LS(){
00139      for(int i=0;i<HcalFrontEndId::maxRmIndex;i++) rm[i].reset_LS(); 
00140   }
00141   int GetRMindex(const std::string& rbx,int rm){
00142       if(rbx.substr(0,3)=="HO0"){
00143          int sect=atoi(rbx.substr(3,2).c_str());
00144          if(sect>12) return -1;
00145          if(rm==1 && (sect==2  || sect==3 || sect==6 || sect==7 || sect==10 || sect==11)) return -1;
00146          if(rm==4 && (sect==12 || sect==1 || sect==4 || sect==5 || sect==8  || sect==9 )) return -1;
00147       }
00148       if(rbx.substr(0,3)=="HO1" || rbx.substr(0,3)=="HO2"){ 
00149          int sect=atoi(rbx.substr(4,2).c_str());
00150          if(sect>12) return -1;
00151          if(sect==1 || sect==3 || sect==5 || sect==7 || sect==9 || sect==11) return -1;
00152       }
00153       HcalFrontEndId id(rbx,rm,1,1,1,1,1);
00154       if(id.rawId()==0) return -1;
00155       return id.rmIndex(); 
00156   }
00157   HcalDetDiagNoiseRMData rm[HcalFrontEndId::maxRmIndex];
00158 };
00159 
00160 HcalDetDiagNoiseMonitor::HcalDetDiagNoiseMonitor(const edm::ParameterSet& ps) :
00161   hcalTBTriggerDataTag_(ps.getParameter<edm::InputTag>("hcalTBTriggerDataTag"))
00162 {
00163   ievt_=0;
00164   run_number=-1;
00165   NoisyEvents=0;
00166   LocalRun=false; 
00167   dataset_seq_number=1;
00168   FirstOrbit=FirstOrbitLS=0xFFFFFFFF;
00169   LastOrbit=LastOrbitLS=0;
00170 
00171   Online_                = ps.getUntrackedParameter<bool>("online",false);
00172   mergeRuns_             = ps.getUntrackedParameter<bool>("mergeRuns",false);
00173   enableCleanup_         = ps.getUntrackedParameter<bool>("enableCleanup",false);
00174   debug_                 = ps.getUntrackedParameter<int>("debug",0);
00175   prefixME_              = ps.getUntrackedParameter<std::string>("subSystemFolder","Hcal/");
00176   if (prefixME_.substr(prefixME_.size()-1,prefixME_.size())!="/")
00177     prefixME_.append("/");
00178   subdir_                = ps.getUntrackedParameter<std::string>("TaskFolder","DetDiagNoiseMonitor_Hcal"); 
00179   if (subdir_.size()>0 && subdir_.substr(subdir_.size()-1,subdir_.size())!="/")
00180     subdir_.append("/");
00181   subdir_=prefixME_+subdir_;
00182   AllowedCalibTypes_     = ps.getUntrackedParameter<std::vector<int> > ("AllowedCalibTypes");
00183   skipOutOfOrderLS_      = ps.getUntrackedParameter<bool>("skipOutOfOrderLS",false);
00184   NLumiBlocks_           = ps.getUntrackedParameter<int>("NLumiBlocks",4000);
00185   makeDiagnostics_       = ps.getUntrackedParameter<bool>("makeDiagnostics",false);
00186 
00187   UseDB            = ps.getUntrackedParameter<bool>  ("UseDB"  , false);
00188   OutputFilePath   = ps.getUntrackedParameter<std::string>("OutputFilePath", "");
00189   HPDthresholdHi   = ps.getUntrackedParameter<double>("NoiseThresholdHPDhi",49.0);
00190   HPDthresholdLo   = ps.getUntrackedParameter<double>("NoiseThresholdHPDlo",10.0);
00191   SpikeThreshold   = ps.getUntrackedParameter<double>("NoiseSpikeThreshold",0.5);
00192   Overwrite        = ps.getUntrackedParameter<bool>  ("Overwrite",true);
00193 
00194   rawDataLabel_ =  ps.getUntrackedParameter<edm::InputTag>("RawDataLabel",edm::InputTag("source",""));
00195   digiLabel_     = ps.getUntrackedParameter<edm::InputTag>("digiLabel",edm::InputTag("hcalDigis"));
00196   L1ADataLabel_  = ps.getUntrackedParameter<edm::InputTag>("gtLabel");
00197   
00198   RMSummary = 0;
00199   needLogicalMap_=true;
00200   setupDone_ = false;
00201 }
00202 
00203 void HcalDetDiagNoiseMonitor::cleanup(){
00204   if(dbe_){
00205     dbe_->setCurrentFolder(subdir_);
00206     dbe_->removeContents();
00207     dbe_ = 0;
00208   }
00209 } 
00210 void HcalDetDiagNoiseMonitor::reset(){}
00211 
00212 
00213 void HcalDetDiagNoiseMonitor::beginRun(const edm::Run& run, const edm::EventSetup& c)
00214 {
00215   if (debug_>1) std::cout <<"HcalDetDiagNoiseMonitor::beginRun"<<std::endl;
00216   HcalBaseDQMonitor::beginRun(run,c);
00217 
00218   if (tevt_==0) this->setup(); // set up histograms if they have not been created before
00219   if (mergeRuns_==false)
00220     this->reset();
00221 
00222   return;
00223 
00224 } 
00225 
00226 void HcalDetDiagNoiseMonitor::setup(){
00227   if (setupDone_)
00228     return;
00229   setupDone_ = true;
00230   // Call base class setup
00231   HcalBaseDQMonitor::setup();
00232   if (!dbe_) return;
00233   RMSummary = new HcalDetDiagNoiseRMSummary();
00234 
00235   std::string name;
00236   if(dbe_!=NULL){    
00237      dbe_->setCurrentFolder(subdir_);   
00238      meEVT_ = dbe_->bookInt("HcalNoiseMonitor Event Number");
00239      dbe_->setCurrentFolder(subdir_+"Common Plots");
00240      
00241      name="RBX Pixel multiplicity";     PixelMult        = dbe_->book1D(name,name,73,0,73);
00242      name="HPD energy";                 HPDEnergy        = dbe_->book1D(name,name,200,0,2500);
00243      name="RBX energy";                 RBXEnergy        = dbe_->book1D(name,name,200,0,3500);
00244      name="Number of zero TS per RBX";  NZeroes          = dbe_->book1D(name,name,100,0,100);
00245      name="Trigger BX Tbit11";          TriggerBx11      = dbe_->book1D(name,name,4000,0,4000);
00246      name="Trigger BX Tbit12";          TriggerBx12      = dbe_->book1D(name,name,4000,0,4000);
00247 
00248      dbe_->setCurrentFolder(subdir_+"HBHE Plots");
00249      name="HBP HPD Noise Rate Pixel above 50fC"; HBP_Rate50    = dbe_->book1D(name,name,73,0,73);
00250      name="HBM HPD Noise Rate Pixel above 50fC"; HBM_Rate50    = dbe_->book1D(name,name,73,0,73);
00251      name="HEP HPD Noise Rate Pixel above 50fC"; HEP_Rate50    = dbe_->book1D(name,name,73,0,73);
00252      name="HEM HPD Noise Rate Pixel above 50fC"; HEM_Rate50    = dbe_->book1D(name,name,73,0,73);
00253      name="HBP HPD Noise Rate HPD above 300fC";  HBP_Rate300   = dbe_->book1D(name,name,73,0,73);
00254      name="HBM HPD Noise Rate HPD above 300fC";  HBM_Rate300   = dbe_->book1D(name,name,73,0,73);
00255      name="HEP HPD Noise Rate HPD above 300fC";  HEP_Rate300   = dbe_->book1D(name,name,73,0,73);
00256      name="HEM HPD Noise Rate HPD above 300fC";  HEM_Rate300   = dbe_->book1D(name,name,73,0,73);
00257 
00258      dbe_->setCurrentFolder(subdir_+"HO Plots");
00259      name="HO0  HPD Noise Rate Pixel above 50fC"; HO0_Rate50   = dbe_->book1D(name,name,49,0,49);
00260      name="HO1P HPD Noise Rate Pixel above 50fC"; HO1P_Rate50   = dbe_->book1D(name,name,48,0,48);
00261      name="HO1M HPD Noise Rate Pixel above 50fC"; HO1M_Rate50   = dbe_->book1D(name,name,48,0,48);
00262      name="HO0 HPD Noise Rate HPD above 300fC";   HO0_Rate300  = dbe_->book1D(name,name,48,0,48);
00263      name="HO1P HPD Noise Rate HPD abGetRMindexove 300fC";  HO1P_Rate300 = dbe_->book1D(name,name,48,0,48);
00264      name="HO1M HPD Noise Rate HPD above 300fC";  HO1M_Rate300 = dbe_->book1D(name,name,48,0,48);
00265       
00266 
00267      dbe_->setCurrentFolder(subdir_+"Noise Spike Plots");
00268 
00269      name="HB RM Spike Map";          HB_RBXmapSpikeCnt= dbe_->book2D(name,name,4,0.5,4.5,36,0.5,36.5);
00270      name="HE RM Spike Map";          HE_RBXmapSpikeCnt= dbe_->book2D(name,name,4,0.5,4.5,36,0.5,36.5);
00271      name="HO RM Spike Map";          HO_RBXmapSpikeCnt= dbe_->book2D(name,name,4,0.5,4.5,36,0.5,36.5);
00272 
00273      std::string title="RM";
00274      HB_RBXmapSpikeCnt->setAxisTitle(title);
00275      HE_RBXmapSpikeCnt->setAxisTitle(title);
00276      HO_RBXmapSpikeCnt->setAxisTitle(title);
00277  
00278      for(int i=0;i<36;i++){
00279         HB_RBXmapSpikeCnt->setBinLabel(i+1,HB_RBX[i],2);
00280         HE_RBXmapSpikeCnt->setBinLabel(i+1,HE_RBX[i],2);
00281         HO_RBXmapSpikeCnt->setBinLabel(i+1,HO_RBX[i],2);
00282      }
00283   } 
00284 
00285 
00286   return;
00287 } 
00288 
00289 void HcalDetDiagNoiseMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup){
00290   getLogicalMap(iSetup);
00291   if (!IsAllowedCalibType()) return;
00292   if (LumiInOrder(iEvent.luminosityBlock())==false) return;
00293   HcalBaseDQMonitor::analyze(iEvent, iSetup);
00294   bool isNoiseEvent=false;  
00295   if(!dbe_) return;
00296   int orbit=-1111;
00297   int bx=-1111;
00298 
00299   // for local runs 
00300   edm::Handle<HcalTBTriggerData> trigger_data;
00301   iEvent.getByLabel(hcalTBTriggerDataTag_, trigger_data);
00302   if(trigger_data.isValid()){
00303       if(trigger_data->triggerWord()>1000) isNoiseEvent=true;
00304       LocalRun=true;
00305   }
00306 
00307   // We do not want to look at Abort Gap events
00308   edm::Handle<FEDRawDataCollection> rawdata;
00309   iEvent.getByLabel(rawDataLabel_,rawdata);
00310   //checking FEDs for calibration information
00311   for(int i=FEDNumbering::MINHCALFEDID;i<=FEDNumbering::MAXHCALFEDID; i++) {
00312       const FEDRawData& fedData = rawdata->FEDData(i) ;
00313       if ( fedData.size() < 24 ) continue ;
00314       orbit= ((const HcalDCCHeader*)(fedData.data()))->getOrbitNumber();
00315       bx=((const HcalDCCHeader*)(fedData.data()))->getBunchId();
00316       if(((const HcalDCCHeader*)(fedData.data()))->getCalibType()!=hc_Null) return;
00317   }
00318 
00319   // Check GCT trigger bits
00320   edm::Handle< L1GlobalTriggerReadoutRecord > gtRecord;
00321   iEvent.getByLabel(L1ADataLabel_, gtRecord);
00322   if(gtRecord.isValid()){
00323     const TechnicalTriggerWord tWord = gtRecord->technicalTriggerWord();
00324     if(tWord.at(11) || tWord.at(12)) isNoiseEvent=true;
00325     if(tWord.at(11)){ TriggerBx11->Fill(bx);}
00326     if(tWord.at(12)){ TriggerBx12->Fill(bx);}
00327   }
00328  
00329   if(!isNoiseEvent) return;
00330   if(ievt_==0){ FirstOrbit=orbit; FirstOrbitLS=orbit; newLS=true;}
00331   if(LastOrbit <orbit) LastOrbit=orbit; 
00332   if(FirstOrbit>orbit) FirstOrbit=orbit;
00333   if(LastOrbitLS <orbit) LastOrbitLS=orbit; 
00334   if(FirstOrbitLS>orbit) FirstOrbitLS=orbit;
00335   if(newLS){ 
00336      FirstOrbitLS=orbit; 
00337      newLS=false;
00338   }
00339 
00340   if(!LocalRun){
00341      double TIME=(double)(LastOrbit-FirstOrbit)/11223.0;
00342      if(TIME>1800.0){
00343         UpdateHistos();
00344         SaveRates();
00345         RMSummary->reset();
00346         FirstOrbit=orbit; 
00347      }
00348   }
00349 
00350   meEVT_->Fill(++ievt_);
00351 
00352   run_number=iEvent.id().run();
00353 
00354   HcalDetDiagNoiseRMEvent RMs[HcalFrontEndId::maxRmIndex];
00355    
00356    edm::Handle<HBHEDigiCollection> hbhe; 
00357    iEvent.getByLabel(digiLabel_,hbhe);
00358    for(HBHEDigiCollection::const_iterator digi=hbhe->begin();digi!=hbhe->end();digi++){
00359      double max=-100/*,sum*/,energy=0; int n_zero=0;
00360      for(int i=0;i<digi->size();i++){
00361        //       sum=adc2fC[digi->sample(i).adc()&0xff]; 
00362        if(max<adc2fC[digi->sample(i).adc()&0xff]) max=adc2fC[digi->sample(i).adc()&0xff];
00363        if(adc2fC[digi->sample(i).adc()&0xff]==0) n_zero++;
00364      }
00365      HcalFrontEndId lmap_entry=logicalMap_->getHcalFrontEndId(digi->id());
00366      int index=lmap_entry.rmIndex(); if(index>=HcalFrontEndId::maxRmIndex) continue;
00367      RMs[index].n_zero++;
00368      if(max>HPDthresholdLo){
00369        for(int i=0;i<digi->size();i++) energy+=adc2fC[digi->sample(i).adc()&0xff]-2.5;
00370        RMs[index].n_pix_lo++;
00371        if(max>HPDthresholdHi){ RMs[index].n_pix_hi++; isNoiseEvent=true;}
00372        RMs[index].energy+=energy;
00373      }
00374    }
00375 
00376    edm::Handle<HODigiCollection> ho; 
00377    iEvent.getByLabel(digiLabel_,ho);
00378    for(HODigiCollection::const_iterator digi=ho->begin();digi!=ho->end();digi++){
00379      double max=-100,energy=0; int Eta=digi->id().ieta(); int Phi=digi->id().iphi(); int n_zero=0;
00380      for(int i=0;i<digi->size()-1;i++){
00381        if(max<adc2fC[digi->sample(i).adc()&0xff]) max=adc2fC[digi->sample(i).adc()&0xff];
00382        if(adc2fC[digi->sample(i).adc()&0xff]==0) n_zero++;
00383      }
00384      if((Eta>=11 && Eta<=15 && Phi>=59 && Phi<=70) || (Eta>=5 && Eta<=10 && Phi>=47 && Phi<=58)){
00385        continue; // ignory SiPMs
00386      }else{
00387        HcalFrontEndId lmap_entry=logicalMap_->getHcalFrontEndId(digi->id());
00388        int index=lmap_entry.rmIndex(); if(index>=HcalFrontEndId::maxRmIndex) continue;
00389        RMs[index].n_zero++;
00390        if(max>HPDthresholdLo){
00391          for(int i=0;i<digi->size();i++) energy+=adc2fC[digi->sample(i).adc()&0xff]-2.5;
00392          RMs[index].n_pix_lo++;
00393          if(max>HPDthresholdHi){ RMs[index].n_pix_hi++; isNoiseEvent=true;}
00394          RMs[index].energy+=energy;
00395        }
00396      }                    
00397    }   
00398 
00399    NoisyEvents++;
00400       
00401    // RMs loop
00402    for(int i=0;i<HcalFrontEndId::maxRmIndex;i++){
00403       if(RMs[i].n_pix_hi>0){
00404          HPDEnergy->Fill(RMs[i].energy);
00405          RMSummary->rm[i].n_th_hi++;
00406          RMSummary->rm[i].n_th_hi_LS++;
00407          if(RMs[i].energy>300) RMSummary->rm[i].n_th_300++;
00408          if(RMs[i].energy>300) RMSummary->rm[i].n_th_300_LS++;
00409          if(RMs[i].n_pix_lo>1) RMSummary->rm[i].n_pix_1++;
00410          if(RMs[i].n_pix_lo>8) RMSummary->rm[i].n_pix_8++;
00411          RMSummary->rm[i].pix+=RMs[i].n_pix_lo;
00412          RMSummary->rm[i].n_pix++;
00413       }
00414    }
00415 
00416    // RBX loop
00417    for(int sd=0;sd<7;sd++) for(int sect=1;sect<=18;sect++){
00418      std::stringstream tempss;
00419      tempss << std::setw(2) << std::setfill('0') << sect;
00420      std::string rbx= subdets[sd]+tempss.str();
00421      
00422      double rbx_energy=0;int pix_mult=0; int n_zero=0; bool isValidRBX=false;
00423      for(int rm=1;rm<=4;rm++){
00424        int index=RMSummary->GetRMindex(rbx,rm);
00425        if(index>0 && index<HcalFrontEndId::maxRmIndex){
00426          rbx_energy+=RMs[index].energy;
00427          pix_mult+=RMs[index].n_pix_lo; 
00428          n_zero+=RMs[index].n_zero; 
00429          isValidRBX=true;
00430        }
00431      }
00432      if(isValidRBX){
00433        PixelMult->Fill(pix_mult);
00434        RBXEnergy->Fill(rbx_energy);
00435        NZeroes->Fill(n_zero);
00436      }
00437    }
00438 
00439    if((ievt_%100)==0 && debug_>0)
00440      std::cout <<ievt_<<"\t"<<NoisyEvents<<std::endl;
00441    return;
00442 }
00443 
00444 void HcalDetDiagNoiseMonitor::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,const edm::EventSetup& c){
00445 int first_rbx=0,last_rbx=0;  
00446   //double TIME=(double)(LastOrbitLS-FirstOrbitLS)/11223.0;
00447   double TIME=23.0;
00448   newLS=true;
00449   if(TIME==0) return;
00450  
00451   for(int sd=0;sd<9;sd++){
00452           if(sd==0){ first_rbx=0;  last_rbx=18;} //HBM
00453           if(sd==1){ first_rbx=18; last_rbx=36;} //HBP
00454           if(sd==0 || sd==1){  // update HB plots
00455             for(int rbx=first_rbx;rbx<last_rbx;rbx++)for(int rm=1;rm<=4;rm++){
00456                int index=RMSummary->GetRMindex(HB_RBX[rbx],rm);
00457                if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00458                double val=RMSummary->rm[index].n_th_hi_LS/TIME;
00459                if(val>SpikeThreshold){
00460                   HB_RBXmapSpikeCnt->Fill(rm,rbx+1,1);
00461                   //printf("%s %i %f (%f)\n",HO_RBX[rbx].c_str(),rm,RMSummary->rm[index].n_th_hi_LS/TIME,TIME);
00462                }
00463                RMSummary->rm[index].reset_LS();
00464             }
00465           }
00466           if(sd==2){ first_rbx=0;  last_rbx=18;} //HEM
00467           if(sd==3){ first_rbx=18; last_rbx=36;} //HEP
00468           if(sd==2 || sd==3){  // update HB plots
00469             for(int rbx=first_rbx;rbx<last_rbx;rbx++)for(int rm=1;rm<=4;rm++){
00470               int index=RMSummary->GetRMindex(HE_RBX[rbx],rm);
00471               if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00472               double val=RMSummary->rm[index].n_th_hi_LS/TIME;
00473               if(val>SpikeThreshold){
00474                   HE_RBXmapSpikeCnt->Fill(rm,rbx+1,1);
00475                   //printf("%s %i %f (%f)\n",HO_RBX[rbx].c_str(),rm,RMSummary->rm[index].n_th_hi_LS/TIME,TIME);
00476               }
00477               RMSummary->rm[index].reset_LS();
00478             }
00479           }
00480           if(sd==4){ first_rbx=6;  last_rbx=12;}  //HO1M
00481           if(sd==5){ first_rbx=12;  last_rbx=24;} //HO0
00482           if(sd==6){ first_rbx=24;  last_rbx=30;} //HO1P
00483           if(sd>3){
00484             for(int rbx=first_rbx;rbx<last_rbx;rbx++)for(int rm=1;rm<=4;rm++){
00485               int index=RMSummary->GetRMindex(HO_RBX[rbx],rm);
00486               if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00487               double val=RMSummary->rm[index].n_th_hi_LS/TIME;
00488               if(val>SpikeThreshold){
00489                   HO_RBXmapSpikeCnt->Fill(rm,rbx+1,1);
00490                   //printf("%s %i %f (%f)\n",HO_RBX[rbx].c_str(),rm,RMSummary->rm[index].n_th_hi_LS/TIME,TIME);
00491               }
00492               RMSummary->rm[index].reset_LS();
00493            }
00494          }
00495     } //sd=0;sd<9
00496 }
00497 
00498 void HcalDetDiagNoiseMonitor::UpdateHistos(){ 
00499    int first_rbx=0; 
00500    double TIME=(double)(LastOrbit-FirstOrbit)/11223.0;
00501    if(TIME==0) return;
00502    for(int sd=0;sd<9;sd++){
00503       if(sd==0){ first_rbx=0; } //HBM
00504       if(sd==1){ first_rbx=18;} //HBP
00505       if(sd==0 || sd==1){  // update HB plots
00506           for(int rbx=0;rbx<18;rbx++)for(int rm=1;rm<=4;rm++){
00507              int index=RMSummary->GetRMindex(HB_RBX[rbx+first_rbx],rm);
00508              if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00509              if(sd==0){
00510                  HBM_Rate50->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_hi/TIME);
00511                  HBM_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00512              }
00513              if(sd==1){
00514                  HBP_Rate50->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_hi/TIME);
00515                  HBP_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00516              }
00517           }     
00518       }
00519       if(sd==2){ first_rbx=0;} //HEM
00520       if(sd==3){ first_rbx=18;} //HEP
00521       if(sd==2 || sd==3){  // update HB plots
00522           for(int rbx=0;rbx<18;rbx++)for(int rm=1;rm<=4;rm++){
00523              int index=RMSummary->GetRMindex(HE_RBX[rbx+first_rbx],rm);
00524              if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00525              if(sd==2){
00526                  HEM_Rate50->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_hi/TIME);
00527                  HEM_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00528              }
00529              if(sd==3){
00530                  HEP_Rate50->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_hi/TIME);
00531                  HEP_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00532              }
00533           }     
00534       }
00535       int n=0;
00536       if(sd==4){ first_rbx=6; n=6;}  //HO1M
00537       if(sd==5){ first_rbx=12;n=12;} //HO0
00538       if(sd==6){ first_rbx=24;n=6;} //HO1P
00539       if(sd>3){ // update HO plots
00540           for(int rbx=0;rbx<n;rbx++)for(int rm=1;rm<=4;rm++){
00541              int index=RMSummary->GetRMindex(HO_RBX[rbx+first_rbx],rm);
00542              if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00543              if(sd==4){
00544                  HO1M_Rate50->setBinContent(rbx*4*2+rm,RMSummary->rm[index].n_th_hi/TIME);
00545                  HO1M_Rate300->setBinContent(rbx*4*2+rm,RMSummary->rm[index].n_th_300/TIME);
00546              }
00547              if(sd==5){
00548                  HO0_Rate50->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_hi/TIME);
00549                  HO0_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00550              } 
00551              if(sd==5){
00552                  HO1P_Rate50->setBinContent(rbx*4*2+rm,RMSummary->rm[index].n_th_hi/TIME);
00553                  HO1P_Rate300->setBinContent(rbx*4+rm,RMSummary->rm[index].n_th_300/TIME);
00554              } 
00555          }
00556       }
00557    } //sd=0;sd<9
00558 } 
00559 
00560 void HcalDetDiagNoiseMonitor::SaveRates(){
00561 char   RBX[20];
00562 int    RM;
00563 float VAL1,VAL2,VAL3,VAL4,VAL5;
00564 char str[500]; 
00565     double TIME=(double)(LastOrbit-FirstOrbit)/11223.0;
00566     if(TIME==0) return;
00567     if(OutputFilePath.size()>0){
00568        if(!Overwrite){
00569           sprintf(str,"%sHcalDetDiagNoiseData_run%06i_%i.root",OutputFilePath.c_str(),run_number,dataset_seq_number);
00570        }else{
00571           sprintf(str,"%sHcalDetDiagNoiseData.root",OutputFilePath.c_str());
00572        }
00573        TFile *theFile = new TFile(str, "RECREATE");
00574        if(!theFile->IsOpen()) return;
00575        theFile->cd();
00576        sprintf(str,"%d",run_number); TObjString run(str);    run.Write("run number");
00577        sprintf(str,"%d",ievt_);      TObjString events(str); events.Write("Total events processed");
00578        sprintf(str,"%d",dataset_seq_number);      TObjString dsnum(str);  dsnum.Write("Dataset number");
00579        Long_t t; t=time(0); strftime(str,30,"%F %T",localtime(&t)); TObjString tm(str);  tm.Write("Dataset creation time");
00580 
00581        TTree *tree   =new TTree("HCAL Noise data","HCAL Noise data");
00582        if(tree==0)   return;
00583        tree->Branch("RBX",            &RBX,     "RBX/C");
00584        tree->Branch("rm",             &RM,      "rm/I");
00585        tree->Branch("RATE_50",        &VAL1,    "RATE_50");
00586        tree->Branch("RATE_300",       &VAL2,    "RATE_300");
00587        tree->Branch("RATE_PIX1",      &VAL3,    "RATE_PIX1");
00588        tree->Branch("RATE_PIX8",      &VAL4,    "RATE_PIX8");
00589        tree->Branch("RATE_PIXMEAN",   &VAL5,    "RATE_PIXMEAN");
00590        for(int rbx=0;rbx<36;rbx++) for(int rm=1;rm<=4;rm++){
00591            int index=RMSummary->GetRMindex(HB_RBX[rbx],rm);
00592            if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00593                sprintf(RBX,"%s",HB_RBX[rbx]);
00594                RM=rm;
00595                VAL1=RMSummary->rm[index].n_th_hi/TIME;
00596                VAL2=RMSummary->rm[index].n_th_300/TIME;
00597                VAL3=RMSummary->rm[index].n_pix_1/TIME;
00598                VAL4=RMSummary->rm[index].n_pix_8/TIME;
00599                if(RMSummary->rm[index].n_pix>0)VAL5=RMSummary->rm[index].pix/RMSummary->rm[index].n_pix; else VAL5=0;
00600                tree->Fill();
00601        }
00602        for(int rbx=0;rbx<36;rbx++) for(int rm=1;rm<=4;rm++){
00603            int index=RMSummary->GetRMindex(HE_RBX[rbx],rm);
00604            if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00605                sprintf(RBX,"%s",HE_RBX[rbx]);
00606                RM=rm;
00607                VAL1=RMSummary->rm[index].n_th_hi/TIME;
00608                VAL2=RMSummary->rm[index].n_th_300/TIME;
00609                VAL3=RMSummary->rm[index].n_pix_1/TIME;
00610                VAL4=RMSummary->rm[index].n_pix_8/TIME;
00611                if(RMSummary->rm[index].n_pix>0)VAL5=RMSummary->rm[index].pix/RMSummary->rm[index].n_pix; else VAL5=0;
00612                tree->Fill();
00613        }
00614        for(int rbx=0;rbx<36;rbx++) for(int rm=1;rm<=4;rm++){
00615            int index=RMSummary->GetRMindex(HO_RBX[rbx],rm);
00616            if(index<0 || index>=HcalFrontEndId::maxRmIndex) continue;
00617                sprintf(RBX,"%s",HO_RBX[rbx]);
00618                RM=rm;
00619                VAL1=RMSummary->rm[index].n_th_hi/TIME;
00620                VAL2=RMSummary->rm[index].n_th_300/TIME;
00621                VAL3=RMSummary->rm[index].n_pix_1/TIME;
00622                VAL4=RMSummary->rm[index].n_pix_8/TIME;
00623                if(RMSummary->rm[index].n_pix>0)VAL5=RMSummary->rm[index].pix/RMSummary->rm[index].n_pix; else VAL5=0;
00624                tree->Fill();
00625        }
00626        theFile->Write();
00627        theFile->Close();
00628        dataset_seq_number++;
00629 
00630    }
00631 }
00632 
00633 
00634 void HcalDetDiagNoiseMonitor::done(){}
00635  
00636 HcalDetDiagNoiseMonitor::~HcalDetDiagNoiseMonitor(){if(LocalRun) UpdateHistos(); SaveRates(); }
00637 
00638 DEFINE_FWK_MODULE (HcalDetDiagNoiseMonitor);