CMS 3D CMS Logo

DTDAQInfo.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2008/12/12 18:04:17 $
00006  *  $Revision: 1.1 $
00007  *  \author G. Cerminara - INFN Torino
00008  */
00009 
00010 
00011 #include "DQM/DTMonitorClient/src/DTDAQInfo.h"
00012 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00013 
00014 #include "FWCore/ServiceRegistry/interface/Service.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/Framework/interface/EventSetup.h"
00017 
00018 #include "DQMServices/Core/interface/DQMStore.h"
00019 #include "DQMServices/Core/interface/MonitorElement.h"
00020 
00021 #include "CondFormats/RunInfo/interface/RunInfo.h"
00022 #include "CondFormats/RunInfo/interface/RunSummary.h"
00023 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
00024 
00025 
00026 
00027 
00028 
00029 using namespace std;
00030 using namespace edm;
00031 
00032 
00033 
00034 DTDAQInfo::DTDAQInfo(const ParameterSet& pset) {}
00035 
00036 
00037 
00038 
00039 DTDAQInfo::~DTDAQInfo() {}
00040 
00041 
00042 
00043 void DTDAQInfo::beginJob(const EventSetup& setup){
00044   // get the DQMStore
00045   theDbe = Service<DQMStore>().operator->();
00046   
00047   // book the ME
00048   theDbe->setCurrentFolder("DT/EventInfo/DAQContents");
00049   // global fraction
00050   totalDAQFraction = theDbe->bookFloat("DTDaqFraction");  
00051   totalDAQFraction->Fill(-1);
00052   // Wheel "fractions" -> will be 0 or 1
00053   for(int wheel = -2; wheel != 3; ++wheel) {
00054     stringstream streams;
00055     streams << "DT_Wheel" << wheel;
00056     daqFractions[wheel] = theDbe->bookFloat(streams.str());
00057     daqFractions[wheel]->Fill(-1);
00058   }
00059 
00060   //
00061 
00062 }
00063 
00064 
00065 
00066 void DTDAQInfo::beginLuminosityBlock(const LuminosityBlock& lumi, const  EventSetup& setup) {
00067   
00068   // create a record key for RunInfoRcd
00069   eventsetup::EventSetupRecordKey recordKey(eventsetup::EventSetupRecordKey::TypeTag::findType("RunInfoRcd"));
00070 
00071   if(setup.find(recordKey) != 0) { 
00072     //get fed summary information
00073     ESHandle<RunInfo> sumFED;
00074     setup.get<RunInfoRcd>().get(sumFED);    
00075     vector<int> fedInIDs = sumFED->m_fed_in;   
00076 
00077     int fedCount=0;
00078 
00079     // the range of DT feds
00080     static int FEDIDmin = FEDNumbering::getDTFEDIds().first;
00081     static int FEDIDMax = FEDNumbering::getDTFEDIds().second;
00082     static int nFeds = FEDIDMax-FEDIDmin-1;
00083 
00084     // loop on all active feds
00085     for(vector<int>::const_iterator fed = fedInIDs.begin();
00086         fed != fedInIDs.end();
00087         ++fed) {
00088       // check if the fed is in the DT range
00089       int wheelNumber = *fed - 772;; 
00090       if(wheelNumber >= -2 && wheelNumber < 3) { // this is a DT FED
00091         daqFractions[wheelNumber]->Fill(1);
00092         fedCount++;
00093       }
00094     }   
00095 
00096     //Fill total fraction of active feds
00097     if(nFeds > 0) totalDAQFraction->Fill(fedCount/nFeds);
00098     else {
00099       totalDAQFraction->Fill(-1);
00100       for(int wheel = -2; wheel != 3; ++wheel) {
00101         daqFractions[wheel]->Fill(-1);
00102       }
00103     }
00104   } else {      
00105     LogWarning("DQM|DTMonitorClient|DTDAQInfo") << "*** Warning: record key not found for RunInfoRcd" << endl;
00106     totalDAQFraction->Fill(-1);               
00107     for(int wheel = -2; wheel != 3; ++wheel) {
00108       daqFractions[wheel]->Fill(-1);
00109     }
00110     return; 
00111   }
00112 }
00113 
00114 
00115 
00116 
00117 void DTDAQInfo::endLuminosityBlock(const LuminosityBlock&  lumi, const  EventSetup& setup){}
00118 
00119 
00120 
00121 void DTDAQInfo::endJob() {}
00122 
00123 
00124 
00125 void DTDAQInfo::analyze(const Event& event, const EventSetup& setup){}
00126 
00127 
00128 

Generated on Tue Jun 9 17:32:34 2009 for CMSSW by  doxygen 1.5.4