CMS 3D CMS Logo

SiPixelCalibDigiProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiPixelCalibDigiProducer
00004 // Class:      SiPixelCalibDigiProducer
00005 // 
00013 //
00014 // Original Author:  Freya Blekman
00015 //         Created:  Wed Oct 31 15:28:52 CET 2007
00016 // $Id: SiPixelCalibDigiProducer.cc,v 1.16 2008/08/29 14:59:27 fblekman Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 
00023 // user include files
00024 #include "FWCore/Framework/interface/EventSetup.h"
00025 
00026 
00027 
00028 
00029 
00030 
00031 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h" 
00032 
00033 #include "DataFormats/SiPixelDigi/interface/SiPixelCalibDigiError.h"
00034 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameConverter.h"
00035 #include "CondFormats/SiPixelObjects/interface/ElectronicIndex.h"
00036 #include "CondFormats/SiPixelObjects/interface/DetectorIndex.h"
00037 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
00038 
00039 #include "SiPixelCalibDigiProducer.h"
00040 #include <sstream>
00041 
00042 //
00043 // constants, enums and typedefs
00044 //
00045 
00046 
00047 //
00048 // static data member definitions
00049 //
00050 
00051 //
00052 // constructors and destructor
00053 //
00054 SiPixelCalibDigiProducer::SiPixelCalibDigiProducer(const edm::ParameterSet& iConfig):
00055   src_(iConfig.getParameter<edm::InputTag>("src")),
00056   iEventCounter_(0),
00057   ignore_non_pattern_(iConfig.getParameter<bool>("ignoreNonPattern")),
00058   control_pattern_size_(iConfig.getParameter<bool>("checkPatternEachEvent")),
00059   includeErrors_(iConfig.getUntrackedParameter<bool>("includeErrors",false)),
00060   errorType(iConfig.getUntrackedParameter<int>("errorTypeNumber",1)),
00061   conf_(iConfig),
00062   number_of_pixels_per_pattern_(0),
00063   use_realeventnumber_(iConfig.getParameter<bool>("useRealEventNumber"))
00064 
00065 {
00066    //register your products
00067   produces< edm::DetSetVector<SiPixelCalibDigi> >();
00068   if(includeErrors_)
00069     produces< edm::DetSetVector<SiPixelCalibDigiError> > ();
00070 
00071 }
00072 
00073 
00074 SiPixelCalibDigiProducer::~SiPixelCalibDigiProducer()
00075 {
00076  
00077    // do anything here that needs to be done at desctruction time
00078    // (e.g. close files, deallocate resources etc.)
00079 
00080 }
00081 
00082 
00083 //
00084 // member functions
00085 //
00086 
00088 // function description:
00089 // this function checks where/when the pattern changes so that one can fill the data from the temporary container into the event 
00090 bool 
00091 SiPixelCalibDigiProducer::store()
00092 {
00093   //  std::cout << "in store() " << std::endl;
00094   if(iEventCounter_%pattern_repeat_==0){
00095     //    std::cout << "now at event " << iEventCounter_ <<" where we save the calibration information into the CMSSW digi";
00096     return 1;
00097   }
00098   else if(iEventCounter_==calib_->expectedTotalEvents())
00099     return 1;
00100   else
00101     return 0;
00102   return 1;
00103 }
00105 // function description:
00106 // 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
00107 void
00108 SiPixelCalibDigiProducer::fill(edm::Event& iEvent, const edm::EventSetup& iSetup)
00109 {
00110 
00111   // figure out which calibration point we're on now..
00112   short icalibpoint = calib_->vcalIndexForEvent(iEventCounter_);
00113   edm::Handle< edm::DetSetVector<PixelDigi> > pixelDigis;
00114   iEvent.getByLabel( src_, pixelDigis );
00115   
00116   edm::LogInfo("SiPixelCalibProducer") << "in fill(), calibpoint " << icalibpoint <<" ndigis " << pixelDigis->size() <<  std::endl;
00117     // loop over the data and store things
00118   edm::DetSetVector<PixelDigi>::const_iterator digiIter;
00119   for(digiIter=pixelDigis->begin(); digiIter!=pixelDigis->end(); ++digiIter){// ITERATOR OVER DET IDs
00120     uint32_t detid = digiIter->id;
00121     edm::DetSet<PixelDigi>::const_iterator ipix; // ITERATOR OVER DIGI DATA  
00122 
00123     for(ipix = digiIter->data.begin(); ipix!=digiIter->end(); ++ipix){
00124       // fill in the appropriate location of the temporary data container
00125       fillPixel(detid,ipix->row(),ipix->column(),icalibpoint,ipix->adc());
00126     }
00127   }
00128 }
00130 // function description:
00131 // this is the function where we check the cabling map and see if we can assign a fed id to the det ID. 
00132 // returns false if no fed <-> detid association was found
00133 bool SiPixelCalibDigiProducer::checkFED(uint32_t detid){
00134   //  edm::LogInfo("SiPixelCalibProducer") << "in checkFED" << std::endl;
00135 
00136   if(detid_to_fedid_[detid])
00137     return true;
00138   for(int fedid=0; fedid<=40; ++fedid){
00139     //    edm::LogInfo("SiPixelCalibProducer") << " looking at fedid " << fedid << std::endl;
00140     SiPixelFrameConverter converter(theCablingMap_.product(),fedid);
00141     if(converter.hasDetUnit(detid)){
00142       detid_to_fedid_[detid]=fedid;
00143       edm::LogInfo("SiPixelCalibDigiProducer") << "matched detid " << detid << " to fed " << detid_to_fedid_[detid] << std::endl;
00144       return true;
00145     }
00146   }
00147   return false;
00148 }
00149 
00151 // function description:
00152 // this is the function where we look in the maps to find the correct calibration digi container, after which the data is filled.
00153 void SiPixelCalibDigiProducer::fillPixel(uint32_t detid, short row, short col, short ipoint, short adc){
00154   //  edm::LogInfo("SiPixelCalibProducer") << " in fillpixel()" << std::endl;
00155  
00156   //  edm::LogInfo("SiPixelCalibProducer") << "in fillPixel " << detid << " " << row << " " << col << " " << ipoint << " " << adc << std::endl;
00157   if(!checkFED(detid)){
00158     edm::LogError("SiPixelCalibDigiProducer") << " was unable to match detid " << detid << " to a FED!" << std::endl;
00159     return;
00160   }
00161   if(!checkPixel(detid,row,col)){
00162     return;
00163   }
00164   // now the check if the pixel exists and fill
00165   //
00166   pixelstruct temppixelworker;
00167   temppixelworker.first=detid;
00168   temppixelworker.second.first=row;
00169   temppixelworker.second.second=col;  
00170   std::map<pixelstruct,SiPixelCalibDigi>::const_iterator ipix = intermediate_data_.find(temppixelworker);
00171   
00172   if(ipix == intermediate_data_.end()){
00173     SiPixelCalibDigi tempdigi(calib_->nVCal());
00174     tempdigi.setrowcol(row,col);
00175     intermediate_data_[temppixelworker]=tempdigi;
00176   }
00177 
00178   intermediate_data_[temppixelworker].fill(ipoint,adc);
00179   return;
00180 
00181 }
00183 // function description:
00184 // 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 
00185 void
00186 SiPixelCalibDigiProducer::clear(){
00187   //  edm::LogInfo("SiPixelCalibProducer") << "in clear() " << std::endl;
00188   // this is where we empty the containers so they can be re-filled
00189   // the idea: the detPixelMap_ container shrinks/expands as a function 
00190   // of the number of pixels looked at at one particular time... 
00191   // unlike the intermediate_data_ container which only expands when 
00192   // detPixelMap_ becomes bigger than intermedate_data_
00193   
00194   // shrink the detPixelMap_
00195   uint32_t tempsize = intermediate_data_.size();
00196   if(tempsize>number_of_pixels_per_pattern_){
00197     edm::LogError("SiPixelCalibDigiProducer") << "Number of pixels in pattern is now: " << tempsize<< ", size is was " << number_of_pixels_per_pattern_ << std::endl;
00198     number_of_pixels_per_pattern_=tempsize;
00199   }
00200 
00201   intermediate_data_.erase(intermediate_data_.begin(),intermediate_data_.end());
00202   intermediate_data_.clear();
00203 
00204   // and erase the error bits
00205   error_data_.erase(error_data_.begin(), error_data_.end());
00206   error_data_.clear();
00207 }
00208 
00210 // function description:
00211 // This method gets the pattern from the calib_ (SiPixelCalibConfiguration) object and fills a vector of pairs that is easier to check 
00212 void 
00213 SiPixelCalibDigiProducer::setPattern(){
00214   //  edm::LogInfo("SiPixelCalibProducer") << "in setPattern()" << std::endl;
00215   uint32_t patternnumber = abs(iEventCounter_-1)/pattern_repeat_;
00216   uint32_t rowpatternnumber = patternnumber/calib_->nColumnPatterns();
00217   uint32_t colpatternnumber = patternnumber%calib_->nColumnPatterns();
00218   edm::LogInfo("SiPixelCalibDigiProducer") << " rowpatternnumbers = " << rowpatternnumber << " " << colpatternnumber << " " << patternnumber << std::endl;
00219   // update currentpattern_
00220   std::vector<short> calibcols = calib_->getColumnPattern();
00221   std::vector<short> calibrows = calib_->getRowPattern();
00222   std::vector<short> temprowvals(0);
00223   std::vector<short> tempcolvals(0);
00224   uint32_t nminuscol=0;
00225   uint32_t nminusrow=0;
00226   uint32_t npatterns=0;
00227   for(uint32_t icol=0; icol<calibcols.size(); icol++){
00228     if(calibcols[icol]==-1){
00229       nminuscol++;
00230     }
00231     else if(nminuscol==colpatternnumber){
00232       //edm::LogInfo("SiPixelCalibProducer") << "col " << calibcols[icol] << std::endl;
00233       short val=calibcols[icol];
00234       tempcolvals.push_back(val);
00235     }
00236     else if (nminuscol> colpatternnumber)
00237       break;
00238   }
00239   for(uint32_t irow=0; irow<calibrows.size(); irow++){
00240     // edm::LogInfo("SiPixelCalibProducer") << "row " << irow <<" "<< nminusrow<<" "  << calibrows[irow] << std::endl;
00241     if(calibrows[irow]==-1)
00242       nminusrow++;
00243     else if(nminusrow==rowpatternnumber){
00244       short val=calibrows[irow];
00245       temprowvals.push_back(val);
00246     }
00247     else if(nminusrow>rowpatternnumber)
00248       break;
00249   }
00250   //now clean up the currentpattern_;
00251   while(currentpattern_.size()>temprowvals.size()*tempcolvals.size()){
00252     currentpattern_.erase(currentpattern_.end());
00253   }
00254   for(uint32_t irow=0; irow<temprowvals.size(); irow++){
00255     for(uint32_t icol=0; icol<tempcolvals.size(); icol++){
00256       std::pair<short,short> pattern(temprowvals[irow],tempcolvals[icol]);
00257       npatterns++;
00258       if(npatterns>currentpattern_.size())
00259         currentpattern_.push_back(pattern);
00260       else
00261         currentpattern_[npatterns-1]=pattern;
00262     }
00263   }
00264 }
00265 
00267 // function description:
00268 // produce method. This is the main loop method
00269 void
00270 SiPixelCalibDigiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00271 {
00272   //  edm::LogInfo("SiPixelCalibDigiProducer") <<"in produce() " << std::endl;
00273   using namespace edm;
00274   iSetup.get<SiPixelCalibConfigurationRcd>().get(calib_);
00275   iSetup.get<TrackerDigiGeometryRecord>().get( theGeometry_ );
00276   iSetup.get<SiPixelFedCablingMapRcd>().get(theCablingMap_);
00277   pattern_repeat_=calib_->getNTriggers()*calib_->nVCal();
00278   if(use_realeventnumber_){
00279     iEventCounter_= iEvent.id().event()-1;
00280   }
00281   else
00282     iEventCounter_++;
00283   if(iEventCounter_%pattern_repeat_==1)
00284     setPattern();
00285   
00286   //  edm::LogInfo("SiPixelCalibDigiProducer") << "now starting fill..." << std::endl;
00287   fill(iEvent,iSetup); // fill method where the actual looping over the digis is done.
00288   //  edm::LogInfo("SiPixelCalibDigiProducer") << "done filling..." << std::endl;
00289   std::auto_ptr<edm::DetSetVector<SiPixelCalibDigi> > pOut(new edm::DetSetVector<SiPixelCalibDigi>);
00290   std::auto_ptr<edm::DetSetVector<SiPixelCalibDigiError> > pErr (new edm::DetSetVector<SiPixelCalibDigiError> );
00291   
00292   // copy things over into pOut if necessary (this is only once per pattern)
00293   if(store()){
00294     //    edm::LogInfo("SiPixelCalibDigiProducer") << "in loop" << std::endl;
00295     for(std::map<pixelstruct,SiPixelCalibDigi>::const_iterator idet=intermediate_data_.begin(); idet!=intermediate_data_.end();++idet){
00296       uint32_t detid=idet->first.first;
00297       if(!control_pattern_size_){
00298         if(! checkPixel(idet->first.first,idet->first.second.first,idet->first.second.second))
00299           continue;
00300       }
00301 
00302       
00303       SiPixelCalibDigi tempdigi=idet->second;
00304       edm::DetSet<SiPixelCalibDigi> & detSet = pOut->find_or_insert(detid);
00305       detSet.data.push_back(tempdigi);
00306     }
00307     if(includeErrors_){
00308       for(std::map<pixelstruct,SiPixelCalibDigiError>::const_iterator ierr=error_data_.begin(); ierr!=error_data_.end();++ierr){
00309         uint32_t detid=ierr->first.first;
00310         SiPixelCalibDigiError temperror = ierr->second;
00311         edm::DetSet<SiPixelCalibDigiError> & errSet = pErr->find_or_insert(detid);
00312         errSet.data.push_back(temperror);
00313       }
00314     }
00315     edm::LogInfo("INFO") << "now filling event " << iEventCounter_ << " as pixel pattern changes every " <<  pattern_repeat_ << " events..." << std::endl;
00316     clear();
00317   }
00318   iEvent.put(pOut);
00319   if(includeErrors_)
00320     iEvent.put(pErr);
00321 }
00322 //-----------------------------------------------
00323 //  method to check that the pixels are actually valid...
00324 bool SiPixelCalibDigiProducer::checkPixel(uint32_t detid, short row, short col){
00325 
00326   if(!control_pattern_size_ && !store())
00327     return true;
00328   
00329   if( !ignore_non_pattern_ )
00330     return true;
00331   
00332   
00333   edm::LogInfo("SiPixelCalibDigiProducer") << "Event" << iEventCounter_ << ",now in checkpixel() " << std::endl;
00334   if(currentpattern_.size()==0)
00335     setPattern();
00336   //  uint32_t iroc;
00337   uint32_t fedid = detid_to_fedid_[detid];
00338 
00339   SiPixelFrameConverter formatter(theCablingMap_.product(),fedid);
00340   sipixelobjects::DetectorIndex detector ={detid, row, col};
00341   sipixelobjects::ElectronicIndex cabling;
00342   
00343   formatter.toCabling(cabling,detector);
00344   // cabling should now contain cabling.roc and cabling.dcol  and cabling.pxid
00345 
00346   // however, the coordinates now need to be converted from dcl, pxid to the row,col coordinates used in the calibration info
00347   sipixelobjects::LocalPixel::DcolPxid loc;
00348   loc.dcol = cabling.dcol;
00349   loc.pxid = cabling.pxid;
00350   sipixelobjects::LocalPixel locpixel(loc);
00351   currentpair_.first = locpixel.rocRow();
00352   currentpair_.second = locpixel.rocCol();
00353 
00354   for(uint32_t i=0; i<currentpattern_.size(); ++i){
00355     //    edm::LogInfo("SiPixelCalibDigiProducer") << "found pair " << currentpair_.first << "," << currentpair_.second << " calib " << currentpattern_[i].first << ","<< currentpattern_[i].second << " input " << row << "," << col << std::endl;
00356     if(currentpair_==currentpattern_[i]){
00357       return true;
00358     }
00359   }
00360   std::ostringstream errorlog;
00361   errorlog <<  "DETID "<<detid<<", row, col (offline)="<<row<<","<<col<<" row, col (ROC) ="<<currentpair_.first<<","<< currentpair_.second<< " found no match in list of patterns: " ;
00362   for(uint32_t i=0; i<currentpattern_.size(); ++i){
00363     if(i!=0 && i!=currentpattern_.size()-1)
00364       errorlog<<" ";
00365     errorlog<<"(";
00366     errorlog<<currentpattern_[i].first;
00367     errorlog<<",";
00368     errorlog<<currentpattern_[i].second;
00369     errorlog<<")";
00370   }
00371   edm::LogError("ERROR") << errorlog.str() << std::endl;
00372   if(includeErrors_){// book the error
00373     
00374     pixelstruct temppixelworker;
00375     temppixelworker.first=detid;
00376     temppixelworker.second.first=row;
00377     temppixelworker.second.second=col;
00378     std::map<pixelstruct,SiPixelCalibDigiError>::const_iterator ierr = error_data_.find(temppixelworker);
00379     if(ierr== error_data_.end()){
00380       SiPixelCalibDigiError temperr(row,col,1);
00381       error_data_[temppixelworker]=temperr;
00382     }
00383   }
00384     
00385   return false;
00386 }
00387 
00388 
00389 // ------------ method called once each job just before starting event loop  ------------
00390 void 
00391 SiPixelCalibDigiProducer::beginRun(const edm::Run & iRun, const edm::EventSetup& iSetup)
00392 {
00393 }
00394 
00395 // ------------ method called once each job just after ending the event loop  ------------
00396 void SiPixelCalibDigiProducer::endJob() {
00397 }
00398 DEFINE_FWK_MODULE(SiPixelCalibDigiProducer);

Generated on Tue Jun 9 17:25:46 2009 for CMSSW by  doxygen 1.5.4