CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/DQM/HcalMonitorModule/src/HcalTimingMonitorModule.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HcalTimingMonitorModule
00004 // Class:      HcalTimingMonitorModule
00005 // 
00013 //
00014 // Original Author:  Dmitry Vishnevskiy
00015 //         Created:  Thu Mar 27 08:12:02 CET 2008
00016 // $Id: HcalTimingMonitorModule.cc,v 1.10 2010/04/04 16:00:38 temple Exp $
00017 //
00018 //
00019 
00020 static const int MAXGEN =10;
00021 static const int MAXRPC =20;
00022 static const int MAXDTBX=20;
00023 static const int MAXCSC =20;    
00024 static const int MAXGMT =20;
00025 static const int TRIG_DT =1;
00026 static const int TRIG_RPC=2;
00027 static const int TRIG_GCT=4;
00028 static const int TRIG_CSC=8;
00029 
00030 // system include files
00031 #include <iostream>
00032 #include <fstream>
00033 #include <cmath>
00034 #include <iosfwd>
00035 #include <bitset>
00036 #include <memory>
00037 
00038 // user include files
00039 #include "FWCore/Framework/interface/Frameworkfwd.h"
00040 #include "FWCore/Framework/interface/EDAnalyzer.h"
00041 #include "FWCore/Framework/interface/Event.h"
00042 #include "FWCore/Framework/interface/ESHandle.h"
00043 #include "FWCore/Framework/interface/MakerMacros.h"
00044 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00045 
00046 // this is to retrieve HCAL digi's
00047 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00048 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00049 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00050 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00051 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00052 
00053 // this is to retrieve GT digi's 
00054 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuRegionalCand.h"
00055 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
00056 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00057 #include "DataFormats/L1GlobalTrigger/interface/L1GtPsbWord.h"
00058 #include "DataFormats/L1GlobalTrigger/interface/L1GtFdlWord.h"
00059 
00060 // DQM stuff
00061 #include "DQMServices/Core/interface/DQMStore.h"
00062 #include "DQMServices/Core/interface/MonitorElement.h"
00063 #include "FWCore/ServiceRegistry/interface/Service.h"
00064  
00065 static const 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,
00066                    13.5,15.,17.,19.,21.,23.,25.,27.,29.5,32.5,35.5,38.5,42.,46.,50.,54.5,59.5,
00067                    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,
00068                    124.5,129.5,137.,147.,157.,167.,177.,187.,197.,209.5,224.5,239.5,254.5,272.,
00069                    292.,312.,334.5,359.5,384.5,359.5,384.5,409.5,434.5,459.5,484.5,509.5,534.5,
00070                    559.5,584.5,609.5,634.5,659.5,684.5,709.5,747.,797.,847.,897.,947.,997.,
00071                    1047.,1109.5,1184.5,1259.5,1334.5,1422.,1522.,1622.,1734.5,1859.5,1984.5,
00072                    1859.5,1984.5,2109.5,2234.5,2359.5,2484.5,2609.5,2734.5,2859.5,2984.5,
00073                    3109.5,3234.5,3359.5,3484.5,3609.5,3797.,4047.,4297.,4547.,4797.,5047.,
00074                    5297.,5609.5,5984.5,6359.5,6734.5,7172.,7672.,8172.,8734.5,9359.5,9984.5};
00075 
00076 
00077 class HcalTimingMonitorModule : public edm::EDAnalyzer {
00078    public:
00079       explicit HcalTimingMonitorModule(const edm::ParameterSet&);
00080       ~HcalTimingMonitorModule();
00081       void   initialize();
00082    
00083    private:
00084       virtual void beginJob() ;
00085       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00086       virtual void endJob() ;
00087       
00088       double GetTime(double *data,int n){
00089              int MaxI=-100; double Time=0,SumT=0,MaxT=-10;
00090              for(int j=0;j<n;++j) if(MaxT<data[j]){ MaxT=data[j]; MaxI=j; }
00091              if (MaxI>=0)
00092                {
00093                  Time=MaxI*data[MaxI];
00094                  SumT=data[MaxI];
00095                  if(MaxI>0){ Time+=(MaxI-1)*data[MaxI-1]; SumT+=data[MaxI-1]; }
00096                  if(MaxI<(n-1)){ Time+=(MaxI+1)*data[MaxI+1]; SumT+=data[MaxI+1]; }
00097                  Time=Time/SumT;
00098                }
00099              return Time;
00100       }
00101       bool isSignal(double *data,int n){
00102         int Imax=-1; double max=-100;
00103         for(int i=0;i<n;i++) if(data[i]>max){max=data[i]; Imax=i;}
00104         if(Imax==0 && Imax==(n-1)) return false;
00105         float sum=data[Imax-1]+data[Imax+1];
00106         if(data[Imax]>5.5 && sum>(data[Imax]*0.20)) return true;
00107         return false;
00108      }
00111      // to be undependent from th DB we want to calculate pedestals using first 100 events
00112      void set_hbhe(int eta,int phi,int depth,int cap,float val){
00113        HBHE[eta+50][phi][depth][cap]+=val;
00114        nHBHE[eta+50][phi][depth][cap]+=1.0;
00115      }   
00116      void set_ho(int eta,int phi,int depth,int cap,float val){
00117        HO[eta+50][phi][depth][cap]+=val;
00118        nHO[eta+50][phi][depth][cap]+=1.0;
00119      }   
00120      void set_hf(int eta,int phi,int depth,int cap,float val){
00121        HF[eta+50][phi][depth][cap]+=val;
00122        nHF[eta+50][phi][depth][cap]+=1.0;
00123      }
00124      double get_ped_hbhe(int eta,int phi,int depth,int cup){
00125       if(nHBHE[eta+50][phi][depth][cup]<10) return 2.5; 
00126       if(nHBHE[eta+50][phi][depth][cup]!=0){
00127          double ped=HBHE[eta+50][phi][depth][cup]/nHBHE[eta+50][phi][depth][cup];
00128          if(ped>1.5 && ped<4.5) return ped;
00129       } 
00130       return 99999; 
00131      }   
00132      double get_ped_ho(int eta,int phi,int depth,int cup){
00133       if(nHO[eta+50][phi][depth][cup]<10) return 2.5; 
00134       if(nHO[eta+50][phi][depth][cup]!=0){
00135          double ped=HO[eta+50][phi][depth][cup]/nHO[eta+50][phi][depth][cup];
00136          if(ped>1.5 && ped<4.5) return ped;
00137       }
00138       return 99999; 
00139      }   
00140      double get_ped_hf(int eta,int phi,int depth,int cup){
00141       if(nHF[eta+50][phi][depth][cup]<10) return 2.5; 
00142       if(nHF[eta+50][phi][depth][cup]!=0){
00143          double ped=HF[eta+50][phi][depth][cup]/nHF[eta+50][phi][depth][cup];
00144          if(ped>1.5 && ped<4.5) return ped;
00145       }
00146       return 99999; 
00147      }   
00148      double HBHE[100][73][5][4];
00149      double nHBHE[100][73][5][4];
00150      double HO[100][73][5][4];
00151      double nHO[100][73][5][4];   
00152      double HF[100][73][5][4];
00153      double nHF[100][73][5][4];
00156     int counterEvt_;
00157     int run_number; 
00158    
00159     int TrigCSC,TrigDT,TrigRPC,TrigGCT;
00160     
00161     edm::ParameterSet parameters_;
00162     DQMStore          *dbe_;
00163     std::string       monitorName_;
00164     int               prescaleLS_,prescaleEvt_;
00165     int               GCTTriggerBit1_;
00166     int               GCTTriggerBit2_;
00167     int               GCTTriggerBit3_;
00168     int               GCTTriggerBit4_;
00169     int               GCTTriggerBit5_;
00170     bool              CosmicsCorr_;
00171     bool              Debug_;
00172     
00173     MonitorElement *HBEnergy,*HEEnergy,*HOEnergy,*HFEnergy;
00174     MonitorElement *DTcand,*RPCbcand,*RPCfcand,*CSCcand,*OR; 
00175   
00176     MonitorElement *HBShapeDT; 
00177     MonitorElement *HBShapeRPC; 
00178     MonitorElement *HBShapeGCT; 
00179     MonitorElement *HOShapeDT; 
00180     MonitorElement *HOShapeRPC; 
00181     MonitorElement *HOShapeGCT; 
00182     MonitorElement *HEShapeCSCp; 
00183     MonitorElement *HEShapeCSCm; 
00184     MonitorElement *HFShapeCSCp; 
00185     MonitorElement *HFShapeCSCm; 
00186     
00187     MonitorElement *HBTimeDT; 
00188     MonitorElement *HBTimeRPC; 
00189     MonitorElement *HBTimeGCT; 
00190     MonitorElement *HOTimeDT; 
00191     MonitorElement *HOTimeRPC; 
00192     MonitorElement *HOTimeGCT; 
00193     MonitorElement *HETimeCSCp; 
00194     MonitorElement *HETimeCSCm;
00195     MonitorElement *HFTimeCSCp; 
00196     MonitorElement *HFTimeCSCm;
00197     
00198     std::string L1ADataLabel;
00199 };
00200 
00201 HcalTimingMonitorModule::HcalTimingMonitorModule(const edm::ParameterSet& iConfig){
00202   std::string str;   
00203    parameters_ = iConfig;
00204    dbe_ = edm::Service<DQMStore>().operator->();
00205    // Base folder for the contents of this job
00206    std::string subsystemname = parameters_.getUntrackedParameter<std::string>("subSystemFolder", "HcalTiming") ;
00207    
00208    monitorName_ = parameters_.getUntrackedParameter<std::string>("monitorName","HcalTiming");
00209    if (monitorName_ != "" ) monitorName_ =subsystemname+"/"+monitorName_+"/" ;
00210    counterEvt_=0;
00211    
00212    // some currently dummy things for compartability with GUI
00213    dbe_->setCurrentFolder(subsystemname+"/EventInfo/");
00214    str="reportSummary";
00215    dbe_->bookFloat(str)->Fill(1);     // Unknown status by default
00216    str="reportSummaryMap";
00217    MonitorElement* me=dbe_->book2D(str,str,5,0,5,1,0,1); // Unknown status by default
00218    TH2F* myhist=me->getTH2F();
00219    myhist->GetXaxis()->SetBinLabel(1,"HB");
00220    myhist->GetXaxis()->SetBinLabel(2,"HE");
00221    myhist->GetXaxis()->SetBinLabel(3,"HO");
00222    myhist->GetXaxis()->SetBinLabel(4,"HF");
00223    myhist->GetYaxis()->SetBinLabel(1,"Status");
00224    // Unknown status by default
00225    myhist->SetBinContent(1,1,-1);
00226    myhist->SetBinContent(2,1,-1);
00227    myhist->SetBinContent(3,1,-1);
00228    myhist->SetBinContent(4,1,-1);
00229    // Add ZDC at some point
00230    myhist->GetXaxis()->SetBinLabel(5,"ZDC");
00231    myhist->SetBinContent(5,1,-1); // no ZDC info known
00232    myhist->SetOption("textcolz");
00233      
00234    run_number=0;
00235    TrigCSC=TrigDT=TrigRPC=TrigGCT=0;
00236    L1ADataLabel   = iConfig.getUntrackedParameter<std::string>("L1ADataLabel" , "l1GtUnpack");
00237    prescaleLS_    = parameters_.getUntrackedParameter<int>("prescaleLS",  1);
00238    prescaleEvt_   = parameters_.getUntrackedParameter<int>("prescaleEvt", 1);
00239    GCTTriggerBit1_= parameters_.getUntrackedParameter<int>("GCTTriggerBit1", -1);         
00240    GCTTriggerBit2_= parameters_.getUntrackedParameter<int>("GCTTriggerBit2", -1);         
00241    GCTTriggerBit3_= parameters_.getUntrackedParameter<int>("GCTTriggerBit3", -1);         
00242    GCTTriggerBit4_= parameters_.getUntrackedParameter<int>("GCTTriggerBit4", -1);         
00243    GCTTriggerBit5_= parameters_.getUntrackedParameter<int>("GCTTriggerBit5", -1);         
00244    CosmicsCorr_   = parameters_.getUntrackedParameter<bool>("CosmicsCorr", true); 
00245    Debug_         = parameters_.getUntrackedParameter<bool>("Debug", true);    
00246    initialize();
00247 }
00248 
00249 HcalTimingMonitorModule::~HcalTimingMonitorModule(){}
00250 
00251 // ------------ method called once each job just before starting event loop  ------------
00252 void HcalTimingMonitorModule::beginJob(){}
00253 // ------------ method called once each job just after ending the event loop  ------------
00254 void HcalTimingMonitorModule::endJob(){}
00255 
00256 void HcalTimingMonitorModule::initialize(){
00257   std::string str;
00258   dbe_->setCurrentFolder(monitorName_+"DebugPlots");
00259   str="L1MuGMTReadoutRecord_getDTBXCands";   DTcand  =dbe_->book1D(str,str,5,-0.5,4.5);
00260   str="L1MuGMTReadoutRecord_getBrlRPCCands"; RPCbcand=dbe_->book1D(str,str,5,-0.5,4.5);
00261   str="L1MuGMTReadoutRecord_getFwdRPCCands"; RPCfcand=dbe_->book1D(str,str,5,-0.5,4.5);
00262   str="L1MuGMTReadoutRecord_getCSCCands";    CSCcand =dbe_->book1D(str,str,5,-0.5,4.5);
00263   str="DT_OR_RPCb_OR_RPCf_OR_CSC";           OR      =dbe_->book1D(str,str,5,-0.5,4.5);
00264   
00265   str="HB Tower Energy (LinADC-PED)"; HBEnergy=dbe_->book1D(str,str,1000,-10,90);
00266   str="HE Tower Energy (LinADC-PED)"; HEEnergy=dbe_->book1D(str,str,1000,-10,90);
00267   str="HO Tower Energy (LinADC-PED)"; HOEnergy=dbe_->book1D(str,str,1000,-10,90);
00268   str="HF Tower Energy (LinADC-PED)"; HFEnergy=dbe_->book1D(str,str,1000,-10,90);
00269   
00270   dbe_->setCurrentFolder(monitorName_+"ShapePlots");
00271   str="HB Shape (DT Trigger)";        HBShapeDT  =dbe_->book1D(str,str,10,-0.5,9.5); 
00272   str="HB Shape (RPC Trigger)";       HBShapeRPC =dbe_->book1D(str,str,10,-0.5,9.5); 
00273   str="HB Shape (GCT Trigger)";       HBShapeGCT =dbe_->book1D(str,str,10,-0.5,9.5); 
00274   str="HO Shape (DT Trigger)";        HOShapeDT  =dbe_->book1D(str,str,10,-0.5,9.5); 
00275   str="HO Shape (RPC Trigger)";       HOShapeRPC =dbe_->book1D(str,str,10,-0.5,9.5); 
00276   str="HO Shape (GCT Trigger)";       HOShapeGCT =dbe_->book1D(str,str,10,-0.5,9.5); 
00277   str="HE+ Shape (CSC Trigger)";      HEShapeCSCp=dbe_->book1D(str,str,10,-0.5,9.5); 
00278   str="HE- Shape (CSC Trigger)";      HEShapeCSCm=dbe_->book1D(str,str,10,-0.5,9.5); 
00279   str="HF+ Shape (CSC Trigger)";      HFShapeCSCp=dbe_->book1D(str,str,10,-0.5,9.5); 
00280   str="HF- Shape (CSC Trigger)";      HFShapeCSCm=dbe_->book1D(str,str,10,-0.5,9.5); 
00281   
00282   dbe_->setCurrentFolder(monitorName_+"TimingPlots");
00283   str="HB Timing (DT Trigger)";       HBTimeDT   =dbe_->book1D(str,str,100,0,10);
00284   str="HB Timing (RPC Trigger)";      HBTimeRPC  =dbe_->book1D(str,str,100,0,10);
00285   str="HB Timing (GCT Trigger)";      HBTimeGCT  =dbe_->book1D(str,str,100,0,10);
00286   str="HO Timing (DT Trigger)";       HOTimeDT   =dbe_->book1D(str,str,100,0,10);
00287   str="HO Timing (RPC Trigger)";      HOTimeRPC  =dbe_->book1D(str,str,100,0,10);
00288   str="HO Timing (GCT Trigger)";      HOTimeGCT  =dbe_->book1D(str,str,100,0,10);
00289   str="HE+ Timing (CSC Trigger)";     HETimeCSCp =dbe_->book1D(str,str,100,0,10);
00290   str="HE- Timing (CSC Trigger)";     HETimeCSCm =dbe_->book1D(str,str,100,0,10);
00291   str="HF+ Timing (CSC Trigger)";     HFTimeCSCp =dbe_->book1D(str,str,100,0,10);
00292   str="HF- Timing (CSC Trigger)";     HFTimeCSCm =dbe_->book1D(str,str,100,0,10);
00293 }
00294 
00295 void HcalTimingMonitorModule::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup){
00296 int HBcnt=0,HEcnt=0,HOcnt=0,HFcnt=0,eta,phi,depth,nTS;
00297 int TRIGGER=0;
00298    counterEvt_++;
00299    if (prescaleEvt_<1)  return;
00300    if (counterEvt_%prescaleEvt_!=0)  return;
00301 
00302    run_number=iEvent.id().run();
00303    // Check GCT trigger bits
00304    edm::Handle< L1GlobalTriggerReadoutRecord > gtRecord;
00305    
00306    if (!iEvent.getByLabel( L1ADataLabel, gtRecord))
00307      return;
00308    const TechnicalTriggerWord tWord = gtRecord->technicalTriggerWord();
00309    const DecisionWord         dWord = gtRecord->decisionWord();
00310    //bool HFselfTrigger   = tWord.at(9);
00311    //bool HOselfTrigger   = tWord.at(11);
00312    bool GCTTrigger1      = dWord.at(GCTTriggerBit1_);     
00313    bool GCTTrigger2      = dWord.at(GCTTriggerBit2_);     
00314    bool GCTTrigger3      = dWord.at(GCTTriggerBit3_);     
00315    bool GCTTrigger4      = dWord.at(GCTTriggerBit4_);     
00316    bool GCTTrigger5      = dWord.at(GCTTriggerBit5_);     
00317    if(GCTTrigger1 || GCTTrigger2 || GCTTrigger3 || GCTTrigger4 || GCTTrigger5){ TrigGCT++; TRIGGER=+TRIG_GCT; }
00318    
00321    // define trigger trigger source (example from GMT group)
00322    edm::Handle<L1MuGMTReadoutCollection> gmtrc_handle; 
00323    if (!iEvent.getByLabel(L1ADataLabel,gmtrc_handle)) return;
00324    L1MuGMTReadoutCollection const* gmtrc = gmtrc_handle.product();
00325    
00326         int idt   =0;
00327         int icsc  =0;
00328         int irpcb =0;
00329         int irpcf =0;
00330         int ndt[5]   = {0,0,0,0,0};
00331         int ncsc[5]  = {0,0,0,0,0};
00332         int nrpcb[5] = {0,0,0,0,0};
00333         int nrpcf[5] = {0,0,0,0,0};
00334         int N;
00335                 
00336         std::vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
00337         std::vector<L1MuGMTReadoutRecord>::const_iterator igmtrr;
00338         N=0;
00339         for(igmtrr=gmt_records.begin(); igmtrr!=gmt_records.end(); igmtrr++) {
00340                 std::vector<L1MuRegionalCand>::const_iterator iter1;
00341                 std::vector<L1MuRegionalCand> rmc;
00342                 // DTBX Trigger
00343                 rmc = igmtrr->getDTBXCands(); 
00344                 for(iter1=rmc.begin(); iter1!=rmc.end(); iter1++) {
00345                         if ( idt < MAXDTBX && !(*iter1).empty() ) {
00346                                 idt++; 
00347                                 if(N<5) ndt[N]++; 
00348                                  
00349                         }        
00350                 }
00351                 // CSC Trigger
00352                 rmc = igmtrr->getCSCCands(); 
00353                 for(iter1=rmc.begin(); iter1!=rmc.end(); iter1++) {
00354                         if ( icsc < MAXCSC && !(*iter1).empty() ) {
00355                                 icsc++; 
00356                                 if(N<5) ncsc[N]++; 
00357                         } 
00358                 }
00359                 // RPCb Trigger
00360                 rmc = igmtrr->getBrlRPCCands();
00361                 for(iter1=rmc.begin(); iter1!=rmc.end(); iter1++) {
00362                         if ( irpcb < MAXRPC && !(*iter1).empty() ) {
00363                                 irpcb++;
00364                                 if(N<5) nrpcb[N]++;
00365                                 
00366                         }  
00367                 }
00368                 // RPCfwd Trigger
00369                 rmc = igmtrr->getFwdRPCCands();
00370                 for(iter1=rmc.begin(); iter1!=rmc.end(); iter1++) {
00371                         if ( irpcf < MAXRPC && !(*iter1).empty() ) {
00372                                 irpcf++;
00373                                 if(N<5) nrpcf[N]++;
00374                                 
00375                         }  
00376                 }
00377                 
00378                 N++;
00379         }
00380         if(ndt[0]) DTcand->Fill(0);
00381         if(ndt[1]) DTcand->Fill(1);
00382         if(ndt[2]) DTcand->Fill(2);
00383         if(ndt[3]) DTcand->Fill(3);
00384         if(ndt[4]) DTcand->Fill(4);
00385         if(ncsc[0]) CSCcand->Fill(0);
00386         if(ncsc[1]) CSCcand->Fill(1);
00387         if(ncsc[2]) CSCcand->Fill(2);
00388         if(ncsc[3]) CSCcand->Fill(3);
00389         if(ncsc[4]) CSCcand->Fill(4);
00390         if(nrpcb[0]) RPCbcand->Fill(0);
00391         if(nrpcb[1]) RPCbcand->Fill(1);
00392         if(nrpcb[2]) RPCbcand->Fill(2);
00393         if(nrpcb[3]) RPCbcand->Fill(3);
00394         if(nrpcb[4]) RPCbcand->Fill(4);
00395         if(nrpcf[0]) RPCfcand->Fill(0);
00396         if(nrpcf[1]) RPCfcand->Fill(1);
00397         if(nrpcf[2]) RPCfcand->Fill(2);
00398         if(nrpcf[3]) RPCfcand->Fill(3);
00399         if(nrpcf[4]) RPCfcand->Fill(4);
00400         if(ndt[0]||nrpcb[0]||nrpcf[0]||ncsc[0]) OR->Fill(0);
00401         if(ndt[1]||nrpcb[1]||nrpcf[1]||ncsc[1]) OR->Fill(1);
00402         if(ndt[2]||nrpcb[2]||nrpcf[2]||ncsc[2]) OR->Fill(2);
00403         if(ndt[3]||nrpcb[3]||nrpcf[3]||ncsc[3]) OR->Fill(3);
00404         if(ndt[4]||nrpcb[4]||nrpcf[4]||ncsc[4]) OR->Fill(4);
00405         
00406         if(ncsc[1]>0 ) { TrigCSC++;   TRIGGER=+TRIG_CSC;  }
00407         if(ndt[1]>0  ) { TrigDT++;    TRIGGER=+TRIG_DT;   }
00408         if(nrpcb[1]>0) { TrigRPC++;   TRIGGER=+TRIG_RPC;  }
00409 
00412    if(counterEvt_<100){
00413      edm::Handle<HBHEDigiCollection> hbhe; 
00414      iEvent.getByType(hbhe);
00415      if (hbhe.isValid())
00416        {
00417          for(HBHEDigiCollection::const_iterator digi=hbhe->begin();digi!=hbhe->end();digi++){
00418            eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00419            if(digi->id().subdet()==HcalBarrel) HBcnt++;
00420            if(digi->id().subdet()==HcalEndcap) HEcnt++;
00421            for(int i=0;i<nTS;i++)
00422              if(digi->sample(i).adc()<20) set_hbhe(eta,phi,depth,digi->sample(i).capid(),adc2fC[digi->sample(i).adc()]);
00423          } 
00424        }  
00425      edm::Handle<HODigiCollection> ho; 
00426      iEvent.getByType(ho);
00427      if (ho.isValid())
00428      {
00429        for(HODigiCollection::const_iterator digi=ho->begin();digi!=ho->end();digi++){
00430          eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00431          HOcnt++;
00432          for(int i=0;i<nTS;i++)
00433                 if(digi->sample(i).adc()<20) set_ho(eta,phi,depth,digi->sample(i).capid(),adc2fC[digi->sample(i).adc()]);
00434        }   
00435      } // if
00436 
00437      edm::Handle<HFDigiCollection> hf;
00438      iEvent.getByType(hf);
00439      if (hf.isValid())
00440        {
00441          for(HFDigiCollection::const_iterator digi=hf->begin();digi!=hf->end();digi++){
00442            eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00443            HFcnt++;
00444            for(int i=0;i<nTS;i++) 
00445              if(digi->sample(i).adc()<20) set_hf(eta,phi,depth,digi->sample(i).capid(),adc2fC[digi->sample(i).adc()]);
00446          }   
00447        }
00448    } // if (counterEvt<100)
00449    else{
00451       double data[10];
00452       
00453       edm::Handle<HBHEDigiCollection> hbhe; 
00454       iEvent.getByType(hbhe);
00455       if (hbhe.isValid())
00456         {
00457           for(HBHEDigiCollection::const_iterator digi=hbhe->begin();digi!=hbhe->end();digi++){
00458             eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00459             if(nTS>10) nTS=10;
00460             if(digi->id().subdet()==HcalBarrel) HBcnt++;
00461             if(digi->id().subdet()==HcalEndcap) HEcnt++;
00462             double energy=0;
00463             for(int i=0;i<nTS;i++){
00464               data[i]=adc2fC[digi->sample(i).adc()]-get_ped_hbhe(eta,phi,depth,digi->sample(i).capid());
00465               energy+=data[i];
00466             }
00467             if(digi->id().subdet()==HcalBarrel) HBEnergy->Fill(energy); 
00468             if(digi->id().subdet()==HcalEndcap) HEEnergy->Fill(energy); 
00469             if(!isSignal(data,nTS)) continue;
00470             for(int i=0;i<nTS;i++){
00471               if(data[i]>-1.0){
00472                 if(digi->id().subdet()==HcalBarrel && (TRIGGER|TRIG_DT)==TRIG_DT)             HBShapeDT->Fill(i,data[i]);
00473                 if(digi->id().subdet()==HcalBarrel && (TRIGGER|TRIG_RPC)==TRIG_RPC)           HBShapeRPC->Fill(i,data[i]);
00474                 if(digi->id().subdet()==HcalBarrel && (TRIGGER|TRIG_GCT)==TRIG_GCT)           HBShapeGCT->Fill(i,data[i]); 
00475                 if(digi->id().subdet()==HcalEndcap && (TRIGGER|TRIG_CSC)==TRIG_CSC && eta>0)  HEShapeCSCp->Fill(i,data[i]);
00476                 if(digi->id().subdet()==HcalEndcap && (TRIGGER|TRIG_CSC)==TRIG_CSC && eta<0)  HEShapeCSCm->Fill(i,data[i]);   
00477               }  
00478             }
00479             double Time=GetTime(data,nTS);
00480             if(digi->id().subdet()==HcalBarrel){
00481               if(CosmicsCorr_) Time+=(7.5*sin((phi*5.0)/180.0*3.14159))/25.0;
00482               if((TRIGGER&TRIG_DT)==TRIG_DT)  HBTimeDT ->Fill(GetTime(data,nTS));
00483               if((TRIGGER&TRIG_RPC)==TRIG_RPC) HBTimeRPC->Fill(GetTime(data,nTS));
00484               if((TRIGGER&TRIG_GCT)==TRIG_GCT) HBTimeGCT->Fill(GetTime(data,nTS));
00485             }else{
00486               if(CosmicsCorr_) Time+=(3.5*sin((phi*5.0)/180.0*3.14159))/25.0; 
00487               if(digi->id().subdet()==HcalEndcap && (TRIGGER&TRIG_CSC)==TRIG_CSC && eta>0) HETimeCSCp->Fill(Time);
00488               if(digi->id().subdet()==HcalEndcap && (TRIGGER&TRIG_CSC)==TRIG_CSC && eta<0) HETimeCSCm->Fill(Time);  
00489             }
00490           }    
00491         } // if (...)
00492 
00493       edm::Handle<HODigiCollection> ho; 
00494       iEvent.getByType(ho);
00495       if (ho.isValid())
00496         {
00497           for(HODigiCollection::const_iterator digi=ho->begin();digi!=ho->end();digi++){
00498             eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00499             if(nTS>10) nTS=10;
00500             HOcnt++; 
00501             double energy=0;
00502             for(int i=0;i<nTS;i++){
00503               data[i]=adc2fC[digi->sample(i).adc()]-get_ped_ho(eta,phi,depth,digi->sample(i).capid());
00504               energy+=data[i];
00505             }     
00506             HOEnergy->Fill(energy); 
00507             if(!isSignal(data,nTS)) continue;
00508             for(int i=0;i<nTS;i++){
00509               if(data[i]>-1.0){
00510                 if((TRIGGER&TRIG_DT)==TRIG_DT)    HOShapeDT->Fill(i,data[i]);
00511                 if((TRIGGER&TRIG_RPC)==TRIG_RPC)  HOShapeRPC->Fill(i,data[i]);
00512                 if((TRIGGER&TRIG_GCT)==TRIG_GCT)  HOShapeGCT->Fill(i,data[i]);       
00513               }
00514             }  
00515             double Time=GetTime(data,nTS);
00516             if(CosmicsCorr_) Time+=(12.0*sin((phi*5.0)/180.0*3.14159))/25.0;    
00517             if((TRIGGER&TRIG_DT)==TRIG_DT)   HOTimeDT->Fill(Time);
00518             if((TRIGGER&TRIG_RPC)==TRIG_RPC) HOTimeRPC->Fill(Time);
00519             if((TRIGGER&TRIG_GCT)==TRIG_GCT) HOTimeGCT->Fill(Time);
00520           }   
00521         }// if (ho)
00522 
00523       edm::Handle<HFDigiCollection> hf; 
00524       iEvent.getByType(hf);
00525       if (hf.isValid())
00526         {
00527           for(HFDigiCollection::const_iterator digi=hf->begin();digi!=hf->end();digi++){
00528             eta=digi->id().ieta(); phi=digi->id().iphi(); depth=digi->id().depth(); nTS=digi->size();
00529             if(nTS>10) nTS=10;
00530             HFcnt++; 
00531             double energy=0;
00532             for(int i=0;i<nTS;i++){
00533                data[i]=adc2fC[digi->sample(i).adc()]-get_ped_hf(eta,phi,depth,digi->sample(i).capid());
00534                energy+=data[i]; 
00535             }
00536             HFEnergy->Fill(energy); 
00537             if(energy<15.0) continue;
00538             for(int i=0;i<nTS;i++){
00539               if(data[i]>-1.0){
00540                 if((TRIGGER&TRIG_CSC)==TRIG_CSC && eta>0)  HFShapeCSCp->Fill(i,data[i]); 
00541                 if((TRIGGER&TRIG_CSC)==TRIG_CSC && eta<0)  HFShapeCSCm->Fill(i,data[i]);    
00542               } 
00543             }
00544             if((TRIGGER&TRIG_CSC)==TRIG_CSC && eta>0) HFTimeCSCp->Fill(GetTime(data,nTS)); 
00545             if((TRIGGER&TRIG_CSC)==TRIG_CSC && eta<0) HFTimeCSCm->Fill(GetTime(data,nTS)); 
00546         }   
00547         } // if (hf)
00550    }
00551    if(Debug_) if((counterEvt_%100)==0) printf("Run: %i,Events processed: %i (HB: %i towers,HE: %i towers,HO: %i towers,HF: %i towers)"
00552                                         " CSC: %i DT: %i RPC: %i GCT: %i\n",
00553                                         run_number,counterEvt_,HBcnt,HEcnt,HOcnt,HFcnt,TrigCSC,TrigDT,TrigRPC,TrigGCT);
00556 }
00557 //define this as a plug-in
00558 DEFINE_FWK_MODULE(HcalTimingMonitorModule);
00559 
00560