CMS 3D CMS Logo

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
42 //
43 #include <cstdlib>
44 #include <string>
45 
46 using namespace std;
47 using namespace edm;
48 
50  : conf_(iConfig),
51  rawin_(consumes<FEDRawDataCollection>(conf_.getParameter<edm::InputTag>("RawInput"))),
52  errin_(consumes<edm::DetSetVector<SiPixelRawDataError>>(conf_.getParameter<edm::InputTag>("ErrorInput"))),
53  saveFile(conf_.getUntrackedParameter<bool>("saveFile", false)),
54  slowDown(conf_.getUntrackedParameter<bool>("slowDown", false)),
55  dirName_(conf_.getUntrackedParameter<string>("DirName", "Pixel/FEDIntegrity/")) {
56  firstRun = true;
57  LogInfo("PixelDQM") << "SiPixelHLTSource::SiPixelHLTSource: Got DQM BackEnd interface" << endl;
58 }
59 
61  // do anything here that needs to be done at desctruction time
62  // (e.g. close files, deallocate resources etc.)
63  LogInfo("PixelDQM") << "SiPixelHLTSource::~SiPixelHLTSource: Destructor" << endl;
64 }
65 
67  LogInfo("PixelDQM") << " SiPixelHLTSource::beginJob - Initialisation ... " << std::endl;
68  iSetup.get<TrackerDigiGeometryRecord>().get(pDD);
69  if (firstRun) {
70  eventNo = 0;
71 
72  firstRun = false;
73  }
74 }
75 
77  // Book Monitoring Elements
78  bookMEs(iBooker);
79 }
80 
81 //------------------------------------------------------------------
82 // Method called for every event
83 //------------------------------------------------------------------
85  eventNo++;
86  // get raw input data
88  iEvent.getByToken(rawin_, rawinput);
89  // get error input data
91  iEvent.getByToken(errin_, errorinput);
92  if (!errorinput.isValid())
93  return;
94 
95  int fedId;
96 
97  for (fedId = 0; fedId <= 39; fedId++) {
98  // get event data for this fed
99  const FEDRawData &fedRawData = rawinput->FEDData(fedId);
100  if (fedRawData.size() != 0)
101  (meRawWords_)->Fill(fedId);
102  } // end for
103 
105 
106  for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++) {
107  if (GeomDetEnumerators::isTrackerPixel((*it)->subDetector())) {
108  uint32_t detId = (*it)->geographicalId();
109  edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(detId);
110  if (isearch != errorinput->end()) {
111  for (di = isearch->data.begin(); di != isearch->data.end(); di++) {
112  fedId = di->getFedId(); // FED the error came from
113  int errorType = di->getType(); // type of error
114  switch (errorType) {
115  case (35):
116  (meNErrors_)->Fill(fedId);
117  break;
118  case (36):
119  (meNErrors_)->Fill(fedId);
120  break;
121  case (37):
122  (meNErrors_)->Fill(fedId);
123  break;
124  case (38):
125  (meNErrors_)->Fill(fedId);
126  break;
127  default:
128  break;
129  }; // end switch
130  } // end for(di
131  } // end if( isearch
132  } // end if( ((*it)->subDetector()
133  } // for(TrackerGeometry
134 
135  edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(0xffffffff);
136 
137  if (isearch != errorinput->end()) { // Not at empty iterator
138  for (di = isearch->data.begin(); di != isearch->data.end(); di++) {
139  fedId = di->getFedId(); // FED the error came from
140  int errorType = di->getType(); // type of error
141  switch (errorType) {
142  case (35):
143  (meNErrors_)->Fill(fedId);
144  break;
145  case (36):
146  (meNErrors_)->Fill(fedId);
147  break;
148  case (37):
149  (meNErrors_)->Fill(fedId);
150  break;
151  case (38):
152  (meNErrors_)->Fill(fedId);
153  break;
154  case (39):
155  (meNCRCs_)->Fill(fedId);
156  break;
157  default:
158  break;
159  }; // end switch
160  } // end for(di
161  } // end if( isearch
162  // slow down...
163  if (slowDown)
164  usleep(100000);
165 }
166 
167 //------------------------------------------------------------------
168 // Book MEs
169 //------------------------------------------------------------------
171  iBooker.cd();
172  iBooker.setCurrentFolder(dirName_);
173 
174  std::string rawhid;
175  std::string errhid;
176  // Get collection name and instantiate Histo Id builder
177  edm::InputTag rawin = conf_.getParameter<edm::InputTag>("RawInput");
178  SiPixelHistogramId *RawHistogramId = new SiPixelHistogramId(rawin.label());
179  edm::InputTag errin = conf_.getParameter<edm::InputTag>("ErrorInput");
180  SiPixelHistogramId *ErrorHistogramId = new SiPixelHistogramId(errin.label());
181 
182  // Is a FED sending raw data
183  meRawWords_ = iBooker.book1D("FEDEntries", "Number of raw words", 40, -0.5, 39.5);
184  meRawWords_->setAxisTitle("Number of raw words", 1);
185 
186  // Number of CRC errors
187  meNCRCs_ = iBooker.book1D("FEDFatal", "Number of fatal errors", 40, -0.5, 39.5);
188  meNCRCs_->setAxisTitle("Number of fatal errors", 1);
189 
190  // Number of translation error words
191  meNErrors_ = iBooker.book1D("FEDNonFatal", "Number of non-fatal errors", 40, -0.5, 39.5);
192  meNErrors_->setAxisTitle("Number of non-fatal errors", 1);
193 
194  delete RawHistogramId;
195  delete ErrorHistogramId;
196 }
197 
T getParameter(std::string const &) const
edm::ParameterSet conf_
SiPixelHLTSource(const edm::ParameterSet &conf)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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
virtual void bookMEs(DQMStore::IBooker &)
std::string dirName_
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< FEDRawDataCollection > rawin_
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
MonitorElement * meRawWords_
void dqmBeginRun(const edm::Run &, edm::EventSetup const &) override
std::string const & label() const
Definition: InputTag.h:36
HLT enums.
T get() const
Definition: EventSetup.h:71
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
~SiPixelHLTSource() override
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
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)
void analyze(const edm::Event &, const edm::EventSetup &) override
bool isTrackerPixel(GeomDetEnumerators::SubDetector m)
Pixel error – collection of errors and error information.
Definition: Run.h:45