#include <DQM/SiPixelMonitorRawData/interface/SiPixelHLTSource.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob (edm::EventSetup const &) |
virtual void | bookMEs () |
virtual void | endJob () |
SiPixelHLTSource (const edm::ParameterSet &conf) | |
~SiPixelHLTSource () | |
Private Attributes | |
edm::ParameterSet | conf_ |
edm::InputTag | errin_ |
int | eventNo |
MonitorElement * | meNCRCs_ |
MonitorElement * | meNErrors_ |
MonitorElement * | meRawWords_ |
edm::ESHandle< TrackerGeometry > | pDD |
edm::InputTag | rawin_ |
bool | saveFile |
bool | slowDown |
DQMStore * | theDMBE |
Implementation: Takes raw data and error data as input, and uses it to populate three histograms indexed by FED id.
Definition at line 49 of file SiPixelHLTSource.h.
SiPixelHLTSource::SiPixelHLTSource | ( | const edm::ParameterSet & | conf | ) | [explicit] |
Definition at line 47 of file SiPixelHLTSource.cc.
References lat::endl(), and theDMBE.
00047 : 00048 conf_(iConfig), 00049 rawin_( conf_.getParameter<edm::InputTag>( "RawInput" ) ), 00050 errin_( conf_.getParameter<edm::InputTag>( "ErrorInput" ) ), 00051 saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ), 00052 slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ) 00053 { 00054 theDMBE = edm::Service<DQMStore>().operator->(); 00055 LogInfo ("PixelDQM") << "SiPixelHLTSource::SiPixelHLTSource: Got DQM BackEnd interface"<<endl; 00056 }
SiPixelHLTSource::~SiPixelHLTSource | ( | ) |
Definition at line 59 of file SiPixelHLTSource.cc.
References lat::endl().
00060 { 00061 // do anything here that needs to be done at desctruction time 00062 // (e.g. close files, deallocate resources etc.) 00063 LogInfo ("PixelDQM") << "SiPixelHLTSource::~SiPixelHLTSource: Destructor"<<endl; 00064 }
void SiPixelHLTSource::analyze | ( | const edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 91 of file SiPixelHLTSource.cc.
References edm::DetSetVector< T >::begin(), detId, edm::DetSetVector< T >::end(), errin_, eventNo, edm::Event::getByLabel(), edm::Handle< T >::isValid(), it, meNCRCs_, meNErrors_, meRawWords_, pDD, rawin_, FEDRawData::size(), and slowDown.
00092 { 00093 eventNo++; 00094 // get raw input data 00095 edm::Handle< FEDRawDataCollection > rawinput; 00096 iEvent.getByLabel( rawin_, rawinput ); 00097 // get error input data 00098 edm::Handle< edm::DetSetVector<SiPixelRawDataError> > errorinput; 00099 iEvent.getByLabel( errin_, errorinput ); 00100 if (!errorinput.isValid()) return; 00101 00102 int fedId; 00103 00104 for(fedId = 0; fedId <= 39; fedId++) { 00105 //get event data for this fed 00106 const FEDRawData& fedRawData = rawinput->FEDData( fedId ); 00107 if (fedRawData.size() != 0) (meRawWords_)->Fill(fedId); 00108 } // end for 00109 00110 edm::DetSet<SiPixelRawDataError>::const_iterator di; 00111 00112 for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){ 00113 if(dynamic_cast<PixelGeomDetUnit*>((*it))!=0){ 00114 uint32_t detId = (*it)->geographicalId(); 00115 edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(detId); 00116 if( isearch != errorinput->end() ) { 00117 for(di = isearch->data.begin(); di != isearch->data.end(); di++) { 00118 fedId = di->getFedId(); // FED the error came from 00119 int errorType = di->getType(); // type of error 00120 switch(errorType) { 00121 case(35) : (meNErrors_)->Fill(fedId); break; 00122 case(36) : (meNErrors_)->Fill(fedId); break; 00123 case(37) : (meNErrors_)->Fill(fedId); break; 00124 case(38) : (meNErrors_)->Fill(fedId); break; 00125 default : break; 00126 }; // end switch 00127 } // end for(di 00128 } // end if( isearch 00129 } // end if(dynamic_cast 00130 } // for(TrackerGeometry 00131 00132 edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(0xffffffff); 00133 00134 if( isearch != errorinput->end() ) { // Not at empty iterator 00135 for(di = isearch->data.begin(); di != isearch->data.end(); di++) { 00136 fedId = di->getFedId(); // FED the error came from 00137 int errorType = di->getType(); // type of error 00138 switch(errorType) { 00139 case(35) : (meNErrors_)->Fill(fedId); break; 00140 case(36) : (meNErrors_)->Fill(fedId); break; 00141 case(37) : (meNErrors_)->Fill(fedId); break; 00142 case(38) : (meNErrors_)->Fill(fedId); break; 00143 case(39) : (meNCRCs_)->Fill(fedId); break; 00144 default : break; 00145 }; // end switch 00146 } // end for(di 00147 } // end if( isearch 00148 // slow down... 00149 if(slowDown) usleep(100000); 00150 00151 }
void SiPixelHLTSource::beginJob | ( | edm::EventSetup const & | iSetup | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 67 of file SiPixelHLTSource.cc.
References bookMEs(), lat::endl(), eventNo, edm::EventSetup::get(), and pDD.
00067 { 00068 LogInfo ("PixelDQM") << " SiPixelHLTSource::beginJob - Initialisation ... " << std::endl; 00069 eventNo = 0; 00070 // Build map 00071 iSetup.get<TrackerDigiGeometryRecord>().get( pDD ); 00072 // Book Monitoring Elements 00073 bookMEs(); 00074 00075 }
void SiPixelHLTSource::bookMEs | ( | ) | [virtual] |
Definition at line 156 of file SiPixelHLTSource.cc.
References DQMStore::cd(), conf_, edm::ParameterSet::getParameter(), edm::InputTag::label(), meNCRCs_, meNErrors_, meRawWords_, MonitorElement::setAxisTitle(), DQMStore::setCurrentFolder(), DQMStore::setVerbose(), and theDMBE.
Referenced by beginJob().
00156 { 00157 00158 theDMBE->setVerbose(0); 00159 theDMBE->cd(); 00160 theDMBE->setCurrentFolder("Pixel/FEDIntegrity/"); 00161 00162 std::string rawhid; 00163 std::string errhid; 00164 // Get collection name and instantiate Histo Id builder 00165 edm::InputTag rawin = conf_.getParameter<edm::InputTag>( "RawInput" ); 00166 SiPixelHistogramId* RawHistogramId = new SiPixelHistogramId( rawin.label() ); 00167 edm::InputTag errin = conf_.getParameter<edm::InputTag>( "ErrorInput" ); 00168 SiPixelHistogramId* ErrorHistogramId = new SiPixelHistogramId( errin.label() ); 00169 // Get DQM interface 00170 DQMStore* theDMBE = edm::Service<DQMStore>().operator->(); 00171 00172 // Is a FED sending raw data 00173 meRawWords_ = theDMBE->book1D("FEDEntries","Number of raw words",40,0.,39.); 00174 meRawWords_->setAxisTitle("Number of raw words",1); 00175 00176 // Number of CRC errors 00177 meNCRCs_ = theDMBE->book1D("FEDFatal","Number of fatal errors",40,0.,39.); 00178 meNCRCs_->setAxisTitle("Number of fatal errors",1); 00179 00180 // Number of translation error words 00181 meNErrors_ = theDMBE->book1D("FEDNonFatal","Number of non-fatal errors",40,0.,39.); 00182 meNErrors_->setAxisTitle("Number of non-fatal errors",1); 00183 00184 delete RawHistogramId; 00185 delete ErrorHistogramId; 00186 00187 }
Reimplemented from edm::EDAnalyzer.
Definition at line 78 of file SiPixelHLTSource.cc.
References conf_, lat::endl(), edm::ParameterSet::getParameter(), DQMStore::save(), saveFile, and theDMBE.
00078 { 00079 00080 if(saveFile) { 00081 LogInfo ("PixelDQM") << " SiPixelHLTSource::endJob - Saving Root File " << std::endl; 00082 std::string outputFile = conf_.getParameter<std::string>("outputFile"); 00083 theDMBE->save( outputFile.c_str() ); 00084 } 00085 00086 }
edm::ParameterSet SiPixelHLTSource::conf_ [private] |
edm::InputTag SiPixelHLTSource::errin_ [private] |
int SiPixelHLTSource::eventNo [private] |
MonitorElement* SiPixelHLTSource::meNCRCs_ [private] |
MonitorElement* SiPixelHLTSource::meNErrors_ [private] |
MonitorElement* SiPixelHLTSource::meRawWords_ [private] |
edm::ESHandle<TrackerGeometry> SiPixelHLTSource::pDD [private] |
edm::InputTag SiPixelHLTSource::rawin_ [private] |
bool SiPixelHLTSource::saveFile [private] |
bool SiPixelHLTSource::slowDown [private] |
DQMStore* SiPixelHLTSource::theDMBE [private] |
Definition at line 67 of file SiPixelHLTSource.h.
Referenced by bookMEs(), endJob(), and SiPixelHLTSource().