CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/DQMServices/Components/src/DQMDaqInfo.cc

Go to the documentation of this file.
00001 #include "DQMServices/Components/src/DQMDaqInfo.h"
00002 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00003 
00004 DQMDaqInfo::DQMDaqInfo(const edm::ParameterSet& iConfig)  
00005 {   
00006 }
00007 
00008 DQMDaqInfo::~DQMDaqInfo()
00009 {  
00010 }
00011 
00012 void DQMDaqInfo::beginLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const  edm::EventSetup& iSetup){
00013   
00014   edm::eventsetup::EventSetupRecordKey recordKey(edm::eventsetup::EventSetupRecordKey::TypeTag::findType("RunInfoRcd"));
00015   
00016   if( 0 != iSetup.find( recordKey ) ) {
00017     
00018     edm::ESHandle<RunInfo> sumFED;
00019     iSetup.get<RunInfoRcd>().get(sumFED);    
00020    
00021     //const RunInfo* summaryFED=sumFED.product();
00022   
00023     std::vector<int> FedsInIds= sumFED->m_fed_in;   
00024 
00025     float  FedCount[9]={0., 0., 0., 0., 0., 0., 0., 0., 0.};
00026     
00027     for(unsigned int fedItr=0;fedItr<FedsInIds.size(); ++fedItr) {
00028       int fedID=FedsInIds[fedItr];     
00029 
00030       if(fedID>=PixelRange.first             &&  fedID<=PixelRange.second)        ++FedCount[Pixel];
00031       if(fedID>=TrackerRange.first           &&  fedID<=TrackerRange.second)      ++FedCount[SiStrip];
00032       if(fedID>=CSCRange.first               &&  fedID<=CSCRange.second)          ++FedCount[CSC];
00033       if(fedID>=RPCRange.first               &&  fedID<=RPCRange.second)          ++FedCount[RPC];
00034       if(fedID>=DTRange.first                &&  fedID<=DTRange.second)           ++FedCount[DT];
00035       if(fedID>=HcalRange.first              &&  fedID<=HcalRange.second)         ++FedCount[Hcal];       
00036       if(fedID>=ECALBarrRange.first          &&  fedID<=ECALBarrRange.second)     ++FedCount[EcalBarrel];      
00037       if((fedID>=ECALEndcapRangeLow.first    && fedID<=ECALEndcapRangeLow.second) ||
00038          (fedID>=ECALEndcapRangeHigh.first && fedID<=ECALEndcapRangeHigh.second)) ++FedCount[EcalEndcap];
00039       if(fedID>=L1TRange.first               &&  fedID<=L1TRange.second)          ++FedCount[L1T];
00040     
00041     }   
00042     
00043     for(int detIndex=0; detIndex<9; ++detIndex) { 
00044       DaqFraction[detIndex]->Fill( FedCount[detIndex]/NumberOfFeds[detIndex]);
00045     }
00046 
00047   }else{    
00048   
00049     for(int detIndex=0; detIndex<9; ++detIndex)  DaqFraction[detIndex]->Fill(-1);               
00050     return; 
00051   }
00052   
00053  
00054 }
00055 
00056 
00057 void DQMDaqInfo::endLuminosityBlock(const edm::LuminosityBlock&  lumiBlock, const  edm::EventSetup& iSetup){
00058 }
00059 
00060 
00061 void 
00062 DQMDaqInfo::beginJob()
00063 {
00064   dbe_ = 0;
00065   dbe_ = edm::Service<DQMStore>().operator->();
00066   
00067   std::string commonFolder = "/EventInfo/DAQContents";  
00068   std::string subsystFolder;
00069   std::string curentFolder;
00070   
00071   subsystFolder="Pixel";  
00072   curentFolder= subsystFolder+commonFolder;
00073   dbe_->setCurrentFolder(curentFolder.c_str());
00074   DaqFraction[Pixel]   = dbe_->bookFloat("PixelDaqFraction");
00075   
00076 
00077   subsystFolder="SiStrip";  
00078   curentFolder=subsystFolder+commonFolder;
00079   dbe_->setCurrentFolder(curentFolder.c_str());
00080   DaqFraction[SiStrip]    = dbe_->bookFloat("SiStripDaqFraction");
00081   
00082   subsystFolder="RPC";  
00083   curentFolder=subsystFolder+commonFolder;
00084   dbe_->setCurrentFolder(curentFolder.c_str());
00085   DaqFraction[RPC]        = dbe_->bookFloat("RPCDaqFraction");
00086   
00087   subsystFolder="CSC";  
00088   curentFolder=subsystFolder+commonFolder;
00089   dbe_->setCurrentFolder(curentFolder.c_str());
00090   DaqFraction[CSC]       = dbe_->bookFloat("CSCDaqFraction");
00091 
00092   subsystFolder="DT";  
00093   curentFolder=subsystFolder+commonFolder;
00094   dbe_->setCurrentFolder(curentFolder.c_str());
00095   DaqFraction[DT]         = dbe_->bookFloat("DTDaqFraction");
00096 
00097   subsystFolder="Hcal";  
00098   curentFolder=subsystFolder+commonFolder;
00099   dbe_->setCurrentFolder(curentFolder.c_str());
00100   DaqFraction[Hcal]       = dbe_->bookFloat("HcalDaqFraction");
00101 
00102   subsystFolder="EcalBarrel";  
00103   curentFolder=subsystFolder+commonFolder;
00104   dbe_->setCurrentFolder(curentFolder.c_str());
00105   DaqFraction[EcalBarrel]       = dbe_->bookFloat("EcalBarrDaqFraction");
00106 
00107   subsystFolder="EcalEndcap";  
00108   curentFolder=subsystFolder+commonFolder;
00109   dbe_->setCurrentFolder(curentFolder.c_str());
00110   DaqFraction[EcalEndcap]       = dbe_->bookFloat("EcalEndDaqFraction");
00111 
00112   subsystFolder="L1T";  
00113   curentFolder=subsystFolder+commonFolder;
00114   dbe_->setCurrentFolder(curentFolder.c_str());
00115   DaqFraction[L1T]       = dbe_->bookFloat("L1TDaqFraction");
00116 
00117 
00118   PixelRange.first  = FEDNumbering::MINSiPixelFEDID;
00119   PixelRange.second = FEDNumbering::MAXSiPixelFEDID;
00120   TrackerRange.first = FEDNumbering::MINSiStripFEDID;
00121   TrackerRange.second = FEDNumbering::MAXSiStripFEDID;
00122   CSCRange.first  = FEDNumbering::MINCSCFEDID;
00123   CSCRange.second = FEDNumbering::MAXCSCFEDID;
00124   RPCRange.first  = 790;
00125   RPCRange.second = 792;
00126   DTRange.first   = 770;
00127   DTRange.second  = 774;
00128   HcalRange.first  = FEDNumbering::MINHCALFEDID;
00129   HcalRange.second = FEDNumbering::MAXHCALFEDID;
00130   L1TRange.first  = FEDNumbering::MINTriggerGTPFEDID;
00131   L1TRange.second = FEDNumbering::MAXTriggerGTPFEDID;
00132   ECALBarrRange.first  = 610;    
00133   ECALBarrRange.second = 645;
00134   ECALEndcapRangeLow.first   = 601;
00135   ECALEndcapRangeLow.second  = 609;
00136   ECALEndcapRangeHigh.first  = 646;
00137   ECALEndcapRangeHigh.second = 654;
00138 
00139   NumberOfFeds[Pixel]   = PixelRange.second-PixelRange.first +1;
00140   NumberOfFeds[SiStrip] = TrackerRange.second-TrackerRange.first +1;
00141   NumberOfFeds[CSC]     = CSCRange.second-CSCRange.first  +1;
00142   NumberOfFeds[RPC]     = RPCRange.second-RPCRange.first  +1;
00143   NumberOfFeds[DT]      = DTRange.second-DTRange.first +1;
00144   NumberOfFeds[Hcal]    = HcalRange.second-HcalRange.first +1;  
00145   NumberOfFeds[EcalBarrel]    = ECALBarrRange.second-ECALBarrRange.first +1 ;
00146   NumberOfFeds[EcalEndcap]    = (ECALEndcapRangeLow.second-ECALEndcapRangeLow.first +1)+(ECALEndcapRangeHigh.second-ECALEndcapRangeHigh.first +1) ;
00147   NumberOfFeds[L1T]    = L1TRange.second-L1TRange.first +1;
00148 
00149 }
00150 
00151 
00152 void 
00153 DQMDaqInfo::endJob() {
00154 }
00155 
00156 
00157 
00158 void
00159 DQMDaqInfo::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00160 { 
00161  
00162 
00163 }
00164