CMS 3D CMS Logo

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