CMS 3D CMS Logo

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