CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 // system include files
21 
22 // user include files
24 
25 
26 
27 
28 
29 
31 
38 
40 #include <sstream>
41 
42 //
43 // constants, enums and typedefs
44 //
45 
46 
47 //
48 // static data member definitions
49 //
50 
51 //
52 // constructors and destructor
53 //
55  src_(iConfig.getParameter<edm::InputTag>("src")),
56  iEventCounter_(0),
57  ignore_non_pattern_(iConfig.getParameter<bool>("ignoreNonPattern")),
58  control_pattern_size_(iConfig.getParameter<bool>("checkPatternEachEvent")),
59  includeErrors_(iConfig.getUntrackedParameter<bool>("includeErrors",false)),
60  errorType(iConfig.getUntrackedParameter<int>("errorTypeNumber",1)),
61  conf_(iConfig),
62  number_of_pixels_per_pattern_(0),
63  use_realeventnumber_(iConfig.getParameter<bool>("useRealEventNumber"))
64 
65 {
66  tPixelDigi = consumes<edm::DetSetVector<PixelDigi>> (src_);
67  //register your products
68  produces< edm::DetSetVector<SiPixelCalibDigi> >();
69  if(includeErrors_)
70  produces< edm::DetSetVector<SiPixelCalibDigiError> > ();
71 
72 }
73 
74 
76 {
77 
78  // do anything here that needs to be done at desctruction time
79  // (e.g. close files, deallocate resources etc.)
80 
81 }
82 
83 
84 //
85 // member functions
86 //
87 
89 // function description:
90 // this function checks where/when the pattern changes so that one can fill the data from the temporary container into the event
91 bool
93 {
94  // std::cout << "in store() " << std::endl;
96  // std::cout << "now at event " << iEventCounter_ <<" where we save the calibration information into the CMSSW digi";
97  return 1;
98  }
99  else if(iEventCounter_==calib_->expectedTotalEvents())
100  return 1;
101  else
102  return 0;
103  return 1;
104 }
106 // function description:
107 // 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
108 void
110 {
111 
112  // figure out which calibration point we're on now..
113  short icalibpoint = calib_->vcalIndexForEvent(iEventCounter_);
115  iEvent.getByToken( tPixelDigi, pixelDigis );
116 
117  edm::LogInfo("SiPixelCalibProducer") << "in fill(), calibpoint " << icalibpoint <<" ndigis " << pixelDigis->size() << std::endl;
118  // loop over the data and store things
120  for(digiIter=pixelDigis->begin(); digiIter!=pixelDigis->end(); ++digiIter){// ITERATOR OVER DET IDs
121  uint32_t detid = digiIter->id;
122  edm::DetSet<PixelDigi>::const_iterator ipix; // ITERATOR OVER DIGI DATA
123 
124  for(ipix = digiIter->data.begin(); ipix!=digiIter->end(); ++ipix){
125  // fill in the appropriate location of the temporary data container
126  fillPixel(detid,ipix->row(),ipix->column(),icalibpoint,ipix->adc());
127  }
128  }
129 }
131 // function description:
132 // this is the function where we check the cabling map and see if we can assign a fed id to the det ID.
133 // returns false if no fed <-> detid association was found
135  // edm::LogInfo("SiPixelCalibProducer") << "in checkFED" << std::endl;
136 
137  if(detid_to_fedid_[detid])
138  return true;
139  for(int fedid=0; fedid<=40; ++fedid){
140  // edm::LogInfo("SiPixelCalibProducer") << " looking at fedid " << fedid << std::endl;
141  SiPixelFrameConverter converter(theCablingMap_.product(),fedid);
142  if(converter.hasDetUnit(detid)){
143  detid_to_fedid_[detid]=fedid;
144  edm::LogInfo("SiPixelCalibDigiProducer") << "matched detid " << detid << " to fed " << detid_to_fedid_[detid] << std::endl;
145  return true;
146  }
147  }
148  return false;
149 }
150 
152 // function description:
153 // this is the function where we look in the maps to find the correct calibration digi container, after which the data is filled.
154 void SiPixelCalibDigiProducer::fillPixel(uint32_t detid, short row, short col, short ipoint, short adc){
155  // edm::LogInfo("SiPixelCalibProducer") << " in fillpixel()" << std::endl;
156 
157  // edm::LogInfo("SiPixelCalibProducer") << "in fillPixel " << detid << " " << row << " " << col << " " << ipoint << " " << adc << std::endl;
158  if(!checkFED(detid)){
159  edm::LogError("SiPixelCalibDigiProducer") << " was unable to match detid " << detid << " to a FED!" << std::endl;
160  return;
161  }
162  if(!checkPixel(detid,row,col)){
163  return;
164  }
165  // now the check if the pixel exists and fill
166  //
167  pixelstruct temppixelworker;
168  temppixelworker.first=detid;
169  temppixelworker.second.first=row;
170  temppixelworker.second.second=col;
171  std::map<pixelstruct,SiPixelCalibDigi>::const_iterator ipix = intermediate_data_.find(temppixelworker);
172 
173  if(ipix == intermediate_data_.end()){
174  SiPixelCalibDigi tempdigi(calib_->nVCal());
175  tempdigi.setrowcol(row,col);
176  intermediate_data_[temppixelworker]=tempdigi;
177  }
178 
179  intermediate_data_[temppixelworker].fill(ipoint,adc);
180  return;
181 
182 }
184 // function description:
185 // 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
186 void
188  // edm::LogInfo("SiPixelCalibProducer") << "in clear() " << std::endl;
189  // this is where we empty the containers so they can be re-filled
190  // the idea: the detPixelMap_ container shrinks/expands as a function
191  // of the number of pixels looked at at one particular time...
192  // unlike the intermediate_data_ container which only expands when
193  // detPixelMap_ becomes bigger than intermedate_data_
194 
195  // shrink the detPixelMap_
196  uint32_t tempsize = intermediate_data_.size();
197  if(tempsize>number_of_pixels_per_pattern_){
198  edm::LogError("SiPixelCalibDigiProducer") << "Number of pixels in pattern is now: " << tempsize<< ", size is was " << number_of_pixels_per_pattern_ << std::endl;
200  }
201 
203  intermediate_data_.clear();
204 
205  // and erase the error bits
206  error_data_.erase(error_data_.begin(), error_data_.end());
207  error_data_.clear();
208 }
209 
211 // function description:
212 // This method gets the pattern from the calib_ (SiPixelCalibConfiguration) object and fills a vector of pairs that is easier to check
213 void
215  // edm::LogInfo("SiPixelCalibProducer") << "in setPattern()" << std::endl;
216  uint32_t patternnumber = abs(iEventCounter_-1)/pattern_repeat_;
217  uint32_t rowpatternnumber = patternnumber/calib_->nColumnPatterns();
218  uint32_t colpatternnumber = patternnumber%calib_->nColumnPatterns();
219  edm::LogInfo("SiPixelCalibDigiProducer") << " rowpatternnumbers = " << rowpatternnumber << " " << colpatternnumber << " " << patternnumber << std::endl;
220  // update currentpattern_
221  std::vector<short> calibcols = calib_->getColumnPattern();
222  std::vector<short> calibrows = calib_->getRowPattern();
223  std::vector<short> temprowvals(0);
224  std::vector<short> tempcolvals(0);
225  uint32_t nminuscol=0;
226  uint32_t nminusrow=0;
227  uint32_t npatterns=0;
228  for(uint32_t icol=0; icol<calibcols.size(); icol++){
229  if(calibcols[icol]==-1){
230  nminuscol++;
231  }
232  else if(nminuscol==colpatternnumber){
233  //edm::LogInfo("SiPixelCalibProducer") << "col " << calibcols[icol] << std::endl;
234  short val=calibcols[icol];
235  tempcolvals.push_back(val);
236  }
237  else if (nminuscol> colpatternnumber)
238  break;
239  }
240  for(uint32_t irow=0; irow<calibrows.size(); irow++){
241  // edm::LogInfo("SiPixelCalibProducer") << "row " << irow <<" "<< nminusrow<<" " << calibrows[irow] << std::endl;
242  if(calibrows[irow]==-1)
243  nminusrow++;
244  else if(nminusrow==rowpatternnumber){
245  short val=calibrows[irow];
246  temprowvals.push_back(val);
247  }
248  else if(nminusrow>rowpatternnumber)
249  break;
250  }
251  //now clean up the currentpattern_;
252  while(currentpattern_.size()>temprowvals.size()*tempcolvals.size()){
253  currentpattern_.erase(currentpattern_.end());
254  }
255  for(uint32_t irow=0; irow<temprowvals.size(); irow++){
256  for(uint32_t icol=0; icol<tempcolvals.size(); icol++){
257  std::pair<short,short> pattern(temprowvals[irow],tempcolvals[icol]);
258  npatterns++;
259  if(npatterns>currentpattern_.size())
260  currentpattern_.push_back(pattern);
261  else
262  currentpattern_[npatterns-1]=pattern;
263  }
264  }
265 }
266 
268 // function description:
269 // produce method. This is the main loop method
270 void
272 {
273  // edm::LogInfo("SiPixelCalibDigiProducer") <<"in produce() " << std::endl;
274  using namespace edm;
275  iSetup.get<SiPixelCalibConfigurationRcd>().get(calib_);
276  iSetup.get<TrackerDigiGeometryRecord>().get( theGeometry_ );
278  pattern_repeat_=calib_->getNTriggers()*calib_->nVCal();
280  iEventCounter_= iEvent.id().event()-1;
281  }
282  else
283  iEventCounter_++;
285  setPattern();
286 
287  // edm::LogInfo("SiPixelCalibDigiProducer") << "now starting fill..." << std::endl;
288  fill(iEvent,iSetup); // fill method where the actual looping over the digis is done.
289  // edm::LogInfo("SiPixelCalibDigiProducer") << "done filling..." << std::endl;
290  std::auto_ptr<edm::DetSetVector<SiPixelCalibDigi> > pOut(new edm::DetSetVector<SiPixelCalibDigi>);
291  std::auto_ptr<edm::DetSetVector<SiPixelCalibDigiError> > pErr (new edm::DetSetVector<SiPixelCalibDigiError> );
292 
293  // copy things over into pOut if necessary (this is only once per pattern)
294  if(store()){
295  // edm::LogInfo("SiPixelCalibDigiProducer") << "in loop" << std::endl;
296  for(std::map<pixelstruct,SiPixelCalibDigi>::const_iterator idet=intermediate_data_.begin(); idet!=intermediate_data_.end();++idet){
297  uint32_t detid=idet->first.first;
299  if(! checkPixel(idet->first.first,idet->first.second.first,idet->first.second.second))
300  continue;
301  }
302 
303 
304  SiPixelCalibDigi tempdigi=idet->second;
305  edm::DetSet<SiPixelCalibDigi> & detSet = pOut->find_or_insert(detid);
306  detSet.data.push_back(tempdigi);
307  }
308  if(includeErrors_){
309  for(std::map<pixelstruct,SiPixelCalibDigiError>::const_iterator ierr=error_data_.begin(); ierr!=error_data_.end();++ierr){
310  uint32_t detid=ierr->first.first;
311  SiPixelCalibDigiError temperror = ierr->second;
312  edm::DetSet<SiPixelCalibDigiError> & errSet = pErr->find_or_insert(detid);
313  errSet.data.push_back(temperror);
314  }
315  }
316  edm::LogInfo("INFO") << "now filling event " << iEventCounter_ << " as pixel pattern changes every " << pattern_repeat_ << " events..." << std::endl;
317  clear();
318  }
319  iEvent.put(pOut);
320  if(includeErrors_)
321  iEvent.put(pErr);
322 }
323 //-----------------------------------------------
324 // method to check that the pixels are actually valid...
325 bool SiPixelCalibDigiProducer::checkPixel(uint32_t detid, short row, short col){
326 
327  if(!control_pattern_size_ && !store())
328  return true;
329 
330  if( !ignore_non_pattern_ )
331  return true;
332 
333 
334  edm::LogInfo("SiPixelCalibDigiProducer") << "Event" << iEventCounter_ << ",now in checkpixel() " << std::endl;
335  if(currentpattern_.size()==0)
336  setPattern();
337  // uint32_t iroc;
338  uint32_t fedid = detid_to_fedid_[detid];
339 
340  SiPixelFrameConverter formatter(theCablingMap_.product(),fedid);
341  sipixelobjects::DetectorIndex detector ={detid, row, col};
343 
344  formatter.toCabling(cabling,detector);
345  // cabling should now contain cabling.roc and cabling.dcol and cabling.pxid
346 
347  // however, the coordinates now need to be converted from dcl, pxid to the row,col coordinates used in the calibration info
349  loc.dcol = cabling.dcol;
350  loc.pxid = cabling.pxid;
351  sipixelobjects::LocalPixel locpixel(loc);
352  currentpair_.first = locpixel.rocRow();
353  currentpair_.second = locpixel.rocCol();
354 
355  for(uint32_t i=0; i<currentpattern_.size(); ++i){
356  // edm::LogInfo("SiPixelCalibDigiProducer") << "found pair " << currentpair_.first << "," << currentpair_.second << " calib " << currentpattern_[i].first << ","<< currentpattern_[i].second << " input " << row << "," << col << std::endl;
358  return true;
359  }
360  }
361  std::ostringstream errorlog;
362  errorlog << "DETID "<<detid<<", row, col (offline)="<<row<<","<<col<<" row, col (ROC) ="<<currentpair_.first<<","<< currentpair_.second<< " found no match in list of patterns: " ;
363  for(uint32_t i=0; i<currentpattern_.size(); ++i){
364  if(i!=0 && i!=currentpattern_.size()-1)
365  errorlog<<" ";
366  errorlog<<"(";
367  errorlog<<currentpattern_[i].first;
368  errorlog<<",";
369  errorlog<<currentpattern_[i].second;
370  errorlog<<")";
371  }
372  edm::LogError("ERROR") << errorlog.str() << std::endl;
373  if(includeErrors_){// book the error
374 
375  pixelstruct temppixelworker;
376  temppixelworker.first=detid;
377  temppixelworker.second.first=row;
378  temppixelworker.second.second=col;
379  std::map<pixelstruct,SiPixelCalibDigiError>::const_iterator ierr = error_data_.find(temppixelworker);
380  if(ierr== error_data_.end()){
381  SiPixelCalibDigiError temperr(row,col,1);
382  error_data_[temppixelworker]=temperr;
383  }
384  }
385 
386  return false;
387 }
388 
int adc(sample_type sample)
get the ADC sample (12 bits)
EventNumber_t event() const
Definition: EventID.h:41
int i
Definition: DBlmapReader.cc:9
std::pair< uint32_t, std::pair< short, short > > pixelstruct
virtual void fill(edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual bool checkFED(uint32_t detid)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
list pattern
Definition: chain.py:104
std::map< uint32_t, uint32_t > detid_to_fedid_
identify pixel inside single ROC
Definition: LocalPixel.h:7
int iEvent
Definition: GenABIO.cc:230
std::map< pixelstruct, SiPixelCalibDigi > intermediate_data_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
edm::ESHandle< SiPixelFedCablingMap > theCablingMap_
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void setrowcol(uint16_t row, uint16_t col)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual bool checkPixel(uint32_t detid, short row, short col)
std::map< pixelstruct, SiPixelCalibDigiError > error_data_
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
double collumn and pixel ID in double collumn representation
Definition: LocalPixel.h:22
std::vector< std::pair< short, short > > currentpattern_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
virtual void fillPixel(uint32_t detid, short row, short col, short ipoint, short adc)
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > tPixelDigi
edm::EventID id() const
Definition: EventBase.h:60
collection_type data
Definition: DetSet.h:78
SiPixelCalibDigiProducer(const edm::ParameterSet &iConfig)
edm::ESHandle< SiPixelCalibConfiguration > calib_
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:350
volatile std::atomic< bool > shutdown_flag false
std::pair< short, short > currentpair_
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
int col
Definition: cuy.py:1008
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
edm::ESHandle< TrackerGeometry > theGeometry_