CMS 3D CMS Logo

SiPixelCalibDigiProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelCalibDigiProducer
4 // Class: SiPixelCalibDigiProducer
5 //
13 //
14 // Original Author: Freya Blekman
15 // Created: Wed Oct 31 15:28:52 CET 2007
16 //
17 //
18 
19 // system include files
20 
21 // user include files
23 
25 
32 
34 #include <sstream>
35 
36 //
37 // constants, enums and typedefs
38 //
39 
40 //
41 // static data member definitions
42 //
43 
44 //
45 // constructors and destructor
46 //
48  : src_(iConfig.getParameter<edm::InputTag>("src")),
49  iEventCounter_(0),
50  ignore_non_pattern_(iConfig.getParameter<bool>("ignoreNonPattern")),
51  control_pattern_size_(iConfig.getParameter<bool>("checkPatternEachEvent")),
52  includeErrors_(iConfig.getUntrackedParameter<bool>("includeErrors", false)),
53  errorType(iConfig.getUntrackedParameter<int>("errorTypeNumber", 1)),
54  conf_(iConfig),
55  number_of_pixels_per_pattern_(0),
56  use_realeventnumber_(iConfig.getParameter<bool>("useRealEventNumber"))
57 
58 {
59  tPixelDigi = consumes<edm::DetSetVector<PixelDigi>>(src_);
60  //register your products
61  produces<edm::DetSetVector<SiPixelCalibDigi>>();
62  if (includeErrors_)
63  produces<edm::DetSetVector<SiPixelCalibDigiError>>();
64 
65  calibToken_ = esConsumes<SiPixelCalibConfiguration, SiPixelCalibConfigurationRcd>();
66  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
67  cablingMapToken_ = esConsumes<SiPixelFedCablingMap, SiPixelFedCablingMapRcd>();
68 }
69 
71  // do anything here that needs to be done at desctruction time
72  // (e.g. close files, deallocate resources etc.)
73 }
74 
75 //
76 // member functions
77 //
78 
80 // function description:
81 // this function checks where/when the pattern changes so that one can fill the data from the temporary container into the event
83  // std::cout << "in store() " << std::endl;
84  if (iEventCounter_ % pattern_repeat_ == 0) {
85  // std::cout << "now at event " << iEventCounter_ <<" where we save the calibration information into the CMSSW digi";
86  return true;
87  } else if (iEventCounter_ == calib_->expectedTotalEvents())
88  return true;
89  else
90  return false;
91  return true;
92 }
94 // function description:
95 // fill function, uses maps to keep track of which pixel is where in the local storage container. Called every event to fill and compress the digi data into calibdig format
97  // figure out which calibration point we're on now..
98  short icalibpoint = calib_->vcalIndexForEvent(iEventCounter_);
100  iEvent.getByToken(tPixelDigi, pixelDigis);
101 
102  edm::LogInfo("SiPixelCalibProducer") << "in fill(), calibpoint " << icalibpoint << " ndigis " << pixelDigis->size()
103  << std::endl;
104  // loop over the data and store things
106  for (digiIter = pixelDigis->begin(); digiIter != pixelDigis->end(); ++digiIter) { // ITERATOR OVER DET IDs
107  uint32_t detid = digiIter->id;
108  edm::DetSet<PixelDigi>::const_iterator ipix; // ITERATOR OVER DIGI DATA
109 
110  for (ipix = digiIter->data.begin(); ipix != digiIter->end(); ++ipix) {
111  // fill in the appropriate location of the temporary data container
112  fillPixel(detid, ipix->row(), ipix->column(), icalibpoint, ipix->adc());
113  }
114  }
115 }
117 // function description:
118 // this is the function where we check the cabling map and see if we can assign a fed id to the det ID.
119 // returns false if no fed <-> detid association was found
121  // edm::LogInfo("SiPixelCalibProducer") << "in checkFED" << std::endl;
122 
123  if (detid_to_fedid_[detid])
124  return true;
125  for (int fedid = 0; fedid <= 40; ++fedid) {
126  // edm::LogInfo("SiPixelCalibProducer") << " looking at fedid " << fedid << std::endl;
128  if (converter.hasDetUnit(detid)) {
129  detid_to_fedid_[detid] = fedid;
130  edm::LogInfo("SiPixelCalibDigiProducer")
131  << "matched detid " << detid << " to fed " << detid_to_fedid_[detid] << std::endl;
132  return true;
133  }
134  }
135  return false;
136 }
137 
139 // function description:
140 // this is the function where we look in the maps to find the correct calibration digi container, after which the data is filled.
141 void SiPixelCalibDigiProducer::fillPixel(uint32_t detid, short row, short col, short ipoint, short adc) {
142  // edm::LogInfo("SiPixelCalibProducer") << " in fillpixel()" << std::endl;
143 
144  // edm::LogInfo("SiPixelCalibProducer") << "in fillPixel " << detid << " " << row << " " << col << " " << ipoint << " " << adc << std::endl;
145  if (!checkFED(detid)) {
146  edm::LogError("SiPixelCalibDigiProducer") << " was unable to match detid " << detid << " to a FED!" << std::endl;
147  return;
148  }
149  if (!checkPixel(detid, row, col)) {
150  return;
151  }
152  // now the check if the pixel exists and fill
153  //
154  pixelstruct temppixelworker;
155  temppixelworker.first = detid;
156  temppixelworker.second.first = row;
157  temppixelworker.second.second = col;
158  std::map<pixelstruct, SiPixelCalibDigi>::const_iterator ipix = intermediate_data_.find(temppixelworker);
159 
160  if (ipix == intermediate_data_.end()) {
161  SiPixelCalibDigi tempdigi(calib_->nVCal());
162  tempdigi.setrowcol(row, col);
163  intermediate_data_[temppixelworker] = tempdigi;
164  }
165 
166  intermediate_data_[temppixelworker].fill(ipoint, adc);
167  return;
168 }
170 // function description:
171 // this function cleans up after everything in the pattern is filled. This involves setting the content of the intermediate_data_ containers to zero and completely emptying the map
173  // edm::LogInfo("SiPixelCalibProducer") << "in clear() " << std::endl;
174  // this is where we empty the containers so they can be re-filled
175  // the idea: the detPixelMap_ container shrinks/expands as a function
176  // of the number of pixels looked at at one particular time...
177  // unlike the intermediate_data_ container which only expands when
178  // detPixelMap_ becomes bigger than intermedate_data_
179 
180  // shrink the detPixelMap_
181  uint32_t tempsize = intermediate_data_.size();
182  if (tempsize > number_of_pixels_per_pattern_) {
183  edm::LogError("SiPixelCalibDigiProducer") << "Number of pixels in pattern is now: " << tempsize << ", size is was "
184  << number_of_pixels_per_pattern_ << std::endl;
186  }
187 
189  intermediate_data_.clear();
190 
191  // and erase the error bits
192  error_data_.erase(error_data_.begin(), error_data_.end());
193  error_data_.clear();
194 }
195 
197 // function description:
198 // This method gets the pattern from the calib_ (SiPixelCalibConfiguration) object and fills a vector of pairs that is easier to check
200  // edm::LogInfo("SiPixelCalibProducer") << "in setPattern()" << std::endl;
201  uint32_t patternnumber = (iEventCounter_ - 1) / pattern_repeat_;
202  uint32_t rowpatternnumber = patternnumber / calib_->nColumnPatterns();
203  uint32_t colpatternnumber = patternnumber % calib_->nColumnPatterns();
204  edm::LogInfo("SiPixelCalibDigiProducer")
205  << " rowpatternnumbers = " << rowpatternnumber << " " << colpatternnumber << " " << patternnumber << std::endl;
206  // update currentpattern_
207  std::vector<short> calibcols = calib_->getColumnPattern();
208  std::vector<short> calibrows = calib_->getRowPattern();
209  std::vector<short> temprowvals(0);
210  std::vector<short> tempcolvals(0);
211  uint32_t nminuscol = 0;
212  uint32_t nminusrow = 0;
213  uint32_t npatterns = 0;
214  for (uint32_t icol = 0; icol < calibcols.size(); icol++) {
215  if (calibcols[icol] == -1) {
216  nminuscol++;
217  } else if (nminuscol == colpatternnumber) {
218  //edm::LogInfo("SiPixelCalibProducer") << "col " << calibcols[icol] << std::endl;
219  short val = calibcols[icol];
220  tempcolvals.push_back(val);
221  } else if (nminuscol > colpatternnumber)
222  break;
223  }
224  for (uint32_t irow = 0; irow < calibrows.size(); irow++) {
225  // edm::LogInfo("SiPixelCalibProducer") << "row " << irow <<" "<< nminusrow<<" " << calibrows[irow] << std::endl;
226  if (calibrows[irow] == -1)
227  nminusrow++;
228  else if (nminusrow == rowpatternnumber) {
229  short val = calibrows[irow];
230  temprowvals.push_back(val);
231  } else if (nminusrow > rowpatternnumber)
232  break;
233  }
234  //now clean up the currentpattern_;
235  while (currentpattern_.size() > temprowvals.size() * tempcolvals.size()) {
236  currentpattern_.erase(currentpattern_.end());
237  }
238  for (uint32_t irow = 0; irow < temprowvals.size(); irow++) {
239  for (uint32_t icol = 0; icol < tempcolvals.size(); icol++) {
240  std::pair<short, short> pattern(temprowvals[irow], tempcolvals[icol]);
241  npatterns++;
242  if (npatterns > currentpattern_.size())
243  currentpattern_.push_back(pattern);
244  else
245  currentpattern_[npatterns - 1] = pattern;
246  }
247  }
248 }
249 
251 // function description:
252 // produce method. This is the main loop method
254  // edm::LogInfo("SiPixelCalibDigiProducer") <<"in produce() " << std::endl;
255  using namespace edm;
256  calib_ = iSetup.getHandle(calibToken_);
260  if (use_realeventnumber_) {
261  iEventCounter_ = iEvent.id().event() - 1;
262  } else
263  iEventCounter_++;
264  if (iEventCounter_ % pattern_repeat_ == 1)
265  setPattern();
266 
267  // edm::LogInfo("SiPixelCalibDigiProducer") << "now starting fill..." << std::endl;
268  fill(iEvent, iSetup); // fill method where the actual looping over the digis is done.
269  // edm::LogInfo("SiPixelCalibDigiProducer") << "done filling..." << std::endl;
270  auto pOut = std::make_unique<edm::DetSetVector<SiPixelCalibDigi>>();
271  auto pErr = std::make_unique<edm::DetSetVector<SiPixelCalibDigiError>>();
272 
273  // copy things over into pOut if necessary (this is only once per pattern)
274  if (store()) {
275  // edm::LogInfo("SiPixelCalibDigiProducer") << "in loop" << std::endl;
276  for (std::map<pixelstruct, SiPixelCalibDigi>::const_iterator idet = intermediate_data_.begin();
277  idet != intermediate_data_.end();
278  ++idet) {
279  uint32_t detid = idet->first.first;
280  if (!control_pattern_size_) {
281  if (!checkPixel(idet->first.first, idet->first.second.first, idet->first.second.second))
282  continue;
283  }
284 
285  SiPixelCalibDigi tempdigi = idet->second;
286  edm::DetSet<SiPixelCalibDigi>& detSet = pOut->find_or_insert(detid);
287  detSet.data.push_back(tempdigi);
288  }
289  if (includeErrors_) {
290  for (std::map<pixelstruct, SiPixelCalibDigiError>::const_iterator ierr = error_data_.begin();
291  ierr != error_data_.end();
292  ++ierr) {
293  uint32_t detid = ierr->first.first;
294  SiPixelCalibDigiError temperror = ierr->second;
295  edm::DetSet<SiPixelCalibDigiError>& errSet = pErr->find_or_insert(detid);
296  errSet.data.push_back(temperror);
297  }
298  }
299  edm::LogInfo("INFO") << "now filling event " << iEventCounter_ << " as pixel pattern changes every "
300  << pattern_repeat_ << " events..." << std::endl;
301  clear();
302  }
303  iEvent.put(std::move(pOut));
304  if (includeErrors_)
305  iEvent.put(std::move(pErr));
306 }
307 //-----------------------------------------------
308 // method to check that the pixels are actually valid...
309 bool SiPixelCalibDigiProducer::checkPixel(uint32_t detid, short row, short col) {
310  if (!control_pattern_size_ && !store())
311  return true;
312 
313  if (!ignore_non_pattern_)
314  return true;
315 
316  edm::LogInfo("SiPixelCalibDigiProducer") << "Event" << iEventCounter_ << ",now in checkpixel() " << std::endl;
317  if (currentpattern_.empty())
318  setPattern();
319  // uint32_t iroc;
320  uint32_t fedid = detid_to_fedid_[detid];
321 
325 
326  formatter.toCabling(cabling, detector);
327  // cabling should now contain cabling.roc and cabling.dcol and cabling.pxid
328 
329  // however, the coordinates now need to be converted from dcl, pxid to the row,col coordinates used in the calibration info
331  loc.dcol = cabling.dcol;
332  loc.pxid = cabling.pxid;
333  sipixelobjects::LocalPixel locpixel(loc);
334  currentpair_.first = locpixel.rocRow();
335  currentpair_.second = locpixel.rocCol();
336 
337  for (uint32_t i = 0; i < currentpattern_.size(); ++i) {
338  // edm::LogInfo("SiPixelCalibDigiProducer") << "found pair " << currentpair_.first << "," << currentpair_.second << " calib " << currentpattern_[i].first << ","<< currentpattern_[i].second << " input " << row << "," << col << std::endl;
339  if (currentpair_ == currentpattern_[i]) {
340  return true;
341  }
342  }
343  std::ostringstream errorlog;
344  errorlog << "DETID " << detid << ", row, col (offline)=" << row << "," << col
345  << " row, col (ROC) =" << currentpair_.first << "," << currentpair_.second
346  << " found no match in list of patterns: ";
347  for (uint32_t i = 0; i < currentpattern_.size(); ++i) {
348  if (i != 0 && i != currentpattern_.size() - 1)
349  errorlog << " ";
350  errorlog << "(";
351  errorlog << currentpattern_[i].first;
352  errorlog << ",";
353  errorlog << currentpattern_[i].second;
354  errorlog << ")";
355  }
356  edm::LogError("ERROR") << errorlog.str() << std::endl;
357  if (includeErrors_) { // book the error
358 
359  pixelstruct temppixelworker;
360  temppixelworker.first = detid;
361  temppixelworker.second.first = row;
362  temppixelworker.second.second = col;
363  std::map<pixelstruct, SiPixelCalibDigiError>::const_iterator ierr = error_data_.find(temppixelworker);
364  if (ierr == error_data_.end()) {
365  SiPixelCalibDigiError temperr(row, col, 1);
366  error_data_[temppixelworker] = temperr;
367  }
368  }
369 
370  return false;
371 }
372 
short vcalIndexForEvent(const uint32_t &eventnumber) const
virtual void fill(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual bool checkFED(uint32_t detid)
std::vector< short > getColumnPattern() const
std::vector< short > getRowPattern() const
std::map< pixelstruct, SiPixelCalibDigiError > error_data_
std::map< pixelstruct, SiPixelCalibDigi > intermediate_data_
Log< level::Error, false > LogError
identify pixel inside single ROC
Definition: LocalPixel.h:7
std::vector< std::pair< short, short > > currentpattern_
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
edm::ESHandle< SiPixelFedCablingMap > theCablingMap_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void setrowcol(uint16_t row, uint16_t col)
edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > cablingMapToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual bool checkPixel(uint32_t detid, short row, short col)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:325
double collumn and pixel ID in double collumn representation
Definition: LocalPixel.h:19
Log< level::Info, false > LogInfo
std::pair< uint32_t, std::pair< short, short > > pixelstruct
virtual void fillPixel(uint32_t detid, short row, short col, short ipoint, short adc)
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > tPixelDigi
collection_type data
Definition: DetSet.h:80
HLT enums.
edm::ESGetToken< SiPixelCalibConfiguration, SiPixelCalibConfigurationRcd > calibToken_
SiPixelCalibDigiProducer(const edm::ParameterSet &iConfig)
col
Definition: cuy.py:1009
edm::ESHandle< SiPixelCalibConfiguration > calib_
std::pair< short, short > currentpair_
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:314
std::map< uint32_t, uint32_t > detid_to_fedid_
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
def move(src, dest)
Definition: eostools.py:511
edm::ESHandle< TrackerGeometry > theGeometry_
uint16_t *__restrict__ uint16_t const *__restrict__ adc