CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelHLTSource.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelMonitorRawData
4 // Class: SiPixelHLTSource
5 //
17 //
18 // Original Author: Andrew York
19 //
20 // Framework
23 // DQM Framework
28 // Geometry
34 // DataFormats
40 //
41 #include <string>
42 #include <stdlib.h>
43 
44 using namespace std;
45 using namespace edm;
46 
48  conf_(iConfig),
49  rawin_( consumes<FEDRawDataCollection>( conf_.getParameter<edm::InputTag>( "RawInput" ) ) ),
50  errin_( consumes<edm::DetSetVector<SiPixelRawDataError> >( conf_.getParameter<edm::InputTag>( "ErrorInput" ) ) ),
51  saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
52  slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
53  dirName_( conf_.getUntrackedParameter<std::string>("DirName","Pixel/FEDIntegrity/") )
54 {
56  LogInfo ("PixelDQM") << "SiPixelHLTSource::SiPixelHLTSource: Got DQM BackEnd interface"<<endl;
57 }
58 
59 
61 {
62  // do anything here that needs to be done at desctruction time
63  // (e.g. close files, deallocate resources etc.)
64  LogInfo ("PixelDQM") << "SiPixelHLTSource::~SiPixelHLTSource: Destructor"<<endl;
65 }
66 
67 
69  firstRun = true;
70 }
71 
73  LogInfo ("PixelDQM") << " SiPixelHLTSource::beginJob - Initialisation ... " << std::endl;
74  iSetup.get<TrackerDigiGeometryRecord>().get( pDD );
75  if(firstRun){
76  eventNo = 0;
77  // Build map
78  // Book Monitoring Elements
79  bookMEs();
80  firstRun = false;
81  }
82 }
83 
84 
86 
87  if(saveFile) {
88  LogInfo ("PixelDQM") << " SiPixelHLTSource::endJob - Saving Root File " << std::endl;
90  theDMBE->save( outputFile.c_str() );
91  }
92 
93 }
94 
95 //------------------------------------------------------------------
96 // Method called for every event
97 //------------------------------------------------------------------
99 {
100  eventNo++;
101  // get raw input data
103  iEvent.getByToken( rawin_, rawinput );
104  // get error input data
106  iEvent.getByToken( errin_, errorinput );
107  if (!errorinput.isValid()) return;
108 
109  int fedId;
110 
111  for(fedId = 0; fedId <= 39; fedId++) {
112  //get event data for this fed
113  const FEDRawData& fedRawData = rawinput->FEDData( fedId );
114  if (fedRawData.size() != 0) (meRawWords_)->Fill(fedId);
115  } // end for
116 
118 
119  for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
120  if( ((*it)->subDetector()==GeomDetEnumerators::PixelBarrel) || ((*it)->subDetector()==GeomDetEnumerators::PixelEndcap) ){
121  uint32_t detId = (*it)->geographicalId();
122  edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(detId);
123  if( isearch != errorinput->end() ) {
124  for(di = isearch->data.begin(); di != isearch->data.end(); di++) {
125  fedId = di->getFedId(); // FED the error came from
126  int errorType = di->getType(); // type of error
127  switch(errorType) {
128  case(35) : (meNErrors_)->Fill(fedId); break;
129  case(36) : (meNErrors_)->Fill(fedId); break;
130  case(37) : (meNErrors_)->Fill(fedId); break;
131  case(38) : (meNErrors_)->Fill(fedId); break;
132  default : break;
133  }; // end switch
134  } // end for(di
135  } // end if( isearch
136  } // end if( ((*it)->subDetector()
137  } // for(TrackerGeometry
138 
139  edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(0xffffffff);
140 
141  if( isearch != errorinput->end() ) { // Not at empty iterator
142  for(di = isearch->data.begin(); di != isearch->data.end(); di++) {
143  fedId = di->getFedId(); // FED the error came from
144  int errorType = di->getType(); // type of error
145  switch(errorType) {
146  case(35) : (meNErrors_)->Fill(fedId); break;
147  case(36) : (meNErrors_)->Fill(fedId); break;
148  case(37) : (meNErrors_)->Fill(fedId); break;
149  case(38) : (meNErrors_)->Fill(fedId); break;
150  case(39) : (meNCRCs_)->Fill(fedId); break;
151  default : break;
152  }; // end switch
153  } // end for(di
154  } // end if( isearch
155  // slow down...
156  if(slowDown) usleep(100000);
157 
158 }
159 
160 //------------------------------------------------------------------
161 // Book MEs
162 //------------------------------------------------------------------
164 
165  theDMBE->cd();
167 
168  std::string rawhid;
169  std::string errhid;
170  // Get collection name and instantiate Histo Id builder
171  edm::InputTag rawin = conf_.getParameter<edm::InputTag>( "RawInput" );
172  SiPixelHistogramId* RawHistogramId = new SiPixelHistogramId( rawin.label() );
173  edm::InputTag errin = conf_.getParameter<edm::InputTag>( "ErrorInput" );
174  SiPixelHistogramId* ErrorHistogramId = new SiPixelHistogramId( errin.label() );
175  // Get DQM interface
177 
178  // Is a FED sending raw data
179  meRawWords_ = theDMBE->book1D("FEDEntries","Number of raw words",40,-0.5,39.5);
180  meRawWords_->setAxisTitle("Number of raw words",1);
181 
182  // Number of CRC errors
183  meNCRCs_ = theDMBE->book1D("FEDFatal","Number of fatal errors",40,-0.5,39.5);
184  meNCRCs_->setAxisTitle("Number of fatal errors",1);
185 
186  // Number of translation error words
187  meNErrors_ = theDMBE->book1D("FEDNonFatal","Number of non-fatal errors",40,-0.5,39.5);
188  meNErrors_->setAxisTitle("Number of non-fatal errors",1);
189 
190  delete RawHistogramId;
191  delete ErrorHistogramId;
192 
193 }
194 
T getParameter(std::string const &) const
edm::ParameterSet conf_
virtual void endJob()
SiPixelHLTSource(const edm::ParameterSet &conf)
virtual void beginJob()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void cd(void)
go to top directory (ie. root)
Definition: DQMStore.cc:561
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * meNErrors_
edm::EDGetTokenT< edm::DetSetVector< SiPixelRawDataError > > errin_
MonitorElement * meNCRCs_
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
std::string dirName_
int iEvent
Definition: GenABIO.cc:243
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
virtual void beginRun(const edm::Run &, edm::EventSetup const &)
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2296
bool isValid() const
Definition: HandleBase.h:76
edm::EDGetTokenT< FEDRawDataCollection > rawin_
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:356
MonitorElement * meRawWords_
virtual void analyze(const edm::Event &, const edm::EventSetup &)
const T & get() const
Definition: EventSetup.h:55
std::string const & label() const
Definition: InputTag.h:42
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:341
volatile std::atomic< bool > shutdown_flag false
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
edm::ESHandle< TrackerGeometry > pDD
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
virtual void bookMEs()
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Pixel error – collection of errors and error information.
Definition: Run.h:41