CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/CondFormats/SiPixelObjects/src/SiPixelGainCalibrationForHLT.cc

Go to the documentation of this file.
00001 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLT.h"
00002 #include "FWCore/Utilities/interface/Exception.h"
00003 #include <cstring>
00004 #include <algorithm>
00005 
00006 //
00007 // Constructors
00008 //
00009 SiPixelGainCalibrationForHLT::SiPixelGainCalibrationForHLT() :
00010   minPed_(0.),
00011   maxPed_(255.),
00012   minGain_(0.),
00013   maxGain_(255.),
00014   numberOfRowsToAverageOver_(80),
00015   nBinsToUseForEncoding_(253),
00016   deadFlag_(255),
00017   noisyFlag_(254)
00018 {
00019    if (deadFlag_ > 0xFF)
00020       throw cms::Exception("GainCalibration Payload configuration error")
00021          << "[SiPixelGainCalibrationHLT::SiPixelGainCalibrationHLT] Dead flag was set to " << deadFlag_ << ", and it must be set less than or equal to 255";
00022 }
00023 //
00024 SiPixelGainCalibrationForHLT::SiPixelGainCalibrationForHLT(float minPed, float maxPed, float minGain, float maxGain) :
00025   minPed_(minPed),
00026   maxPed_(maxPed),
00027   minGain_(minGain),
00028   maxGain_(maxGain),
00029   numberOfRowsToAverageOver_(80),
00030   nBinsToUseForEncoding_(253),
00031   deadFlag_(255),
00032   noisyFlag_(254)
00033 {
00034    if (deadFlag_ > 0xFF)
00035       throw cms::Exception("GainCalibration Payload configuration error")
00036          << "[SiPixelGainCalibrationHLT::SiPixelGainCalibrationHLT] Dead flag was set to " << deadFlag_ << ", and it must be set less than or equal to 255";
00037 }
00038 
00039 bool SiPixelGainCalibrationForHLT::put(const uint32_t& DetId, Range input, const int& nCols) {
00040   // put in SiPixelGainCalibrationForHLT of DetId
00041 
00042   Registry::iterator p = std::lower_bound(indexes.begin(),indexes.end(),DetId,SiPixelGainCalibrationForHLT::StrictWeakOrdering());
00043   if (p!=indexes.end() && p->detid==DetId)
00044     return false;
00045   
00046   size_t sd= input.second-input.first;
00047   DetRegistry detregistry;
00048   detregistry.detid=DetId;
00049   detregistry.ibegin=v_pedestals.size();
00050   detregistry.iend=v_pedestals.size()+sd;
00051   detregistry.ncols=nCols;
00052   indexes.insert(p,detregistry);
00053 
00054   v_pedestals.insert(v_pedestals.end(),input.first,input.second);
00055   return true;
00056 }
00057 
00058 const int SiPixelGainCalibrationForHLT::getNCols(const uint32_t& DetId) const {
00059   // get number of columns of DetId
00060   RegistryIterator p = std::lower_bound(indexes.begin(),indexes.end(),DetId,SiPixelGainCalibrationForHLT::StrictWeakOrdering());
00061   if (p==indexes.end()|| p->detid!=DetId) 
00062     return 0;
00063   else
00064   {
00065     return p->ncols;
00066   }
00067 }
00068 
00069 const SiPixelGainCalibrationForHLT::Range SiPixelGainCalibrationForHLT::getRange(const uint32_t& DetId) const {
00070   // get SiPixelGainCalibrationForHLT Range of DetId
00071   
00072   RegistryIterator p = std::lower_bound(indexes.begin(),indexes.end(),DetId,SiPixelGainCalibrationForHLT::StrictWeakOrdering());
00073   if (p==indexes.end()|| p->detid!=DetId) 
00074     return SiPixelGainCalibrationForHLT::Range(v_pedestals.end(),v_pedestals.end()); 
00075   else 
00076     return SiPixelGainCalibrationForHLT::Range(v_pedestals.begin()+p->ibegin,v_pedestals.begin()+p->iend);
00077 }
00078 
00079 const std::pair<const SiPixelGainCalibrationForHLT::Range, const int>
00080 SiPixelGainCalibrationForHLT::getRangeAndNCols(const uint32_t& DetId) const {
00081   RegistryIterator p = std::lower_bound(indexes.begin(),indexes.end(),DetId,SiPixelGainCalibrationForHLT::StrictWeakOrdering());
00082   if (p==indexes.end()|| p->detid!=DetId) 
00083     return std::make_pair(SiPixelGainCalibrationForHLT::Range(v_pedestals.end(),v_pedestals.end()), 0); 
00084   else 
00085     return std::make_pair(SiPixelGainCalibrationForHLT::Range(v_pedestals.begin()+p->ibegin,v_pedestals.begin()+p->iend), p->ncols);
00086 }
00087 
00088 void SiPixelGainCalibrationForHLT::getDetIds(std::vector<uint32_t>& DetIds_) const {
00089   // returns vector of DetIds in map
00090   SiPixelGainCalibrationForHLT::RegistryIterator begin = indexes.begin();
00091   SiPixelGainCalibrationForHLT::RegistryIterator end   = indexes.end();
00092   for (SiPixelGainCalibrationForHLT::RegistryIterator p=begin; p != end; ++p) {
00093     DetIds_.push_back(p->detid);
00094   }
00095 }
00096 
00097 void SiPixelGainCalibrationForHLT::setData(float ped, float gain, std::vector<char>& vped, bool thisColumnIsDead, bool thisColumnIsNoisy){
00098   
00099   float theEncodedGain=0;
00100   float theEncodedPed=0;
00101   if(!thisColumnIsDead && !thisColumnIsNoisy){
00102     theEncodedGain = encodeGain(gain);
00103     theEncodedPed = encodePed (ped);
00104   }
00105 
00106   unsigned int ped_   = (static_cast<unsigned int>(theEncodedPed))  & 0xFF; 
00107   unsigned int gain_  = (static_cast<unsigned int>(theEncodedGain)) & 0xFF;
00108 
00109   if (thisColumnIsDead)
00110   {
00111      ped_  = deadFlag_ & 0xFF;
00112      gain_ = deadFlag_ & 0xFF;
00113   }
00114   else if (thisColumnIsNoisy)
00115   {
00116      ped_  = noisyFlag_ & 0xFF;
00117      gain_ = noisyFlag_ & 0xFF;
00118   }
00119 
00120   unsigned int data = (ped_ << 8) | gain_ ;
00121   vped.resize(vped.size()+2);
00122   // insert in vector of char
00123   ::memcpy((void*)(&vped[vped.size()-2]),(void*)(&data),2);
00124 }
00125 
00126 float SiPixelGainCalibrationForHLT::getPed(const int& col, const int& row, const Range& range, const int& nCols, bool& isDeadColumn, bool& isNoisyColumn) const {
00127    // TODO MERGE THIS FUNCTION WITH GET GAIN, then provide wrappers
00128 
00129   // determine what averaged data block we are in (there should be 1 or 2 of these depending on if plaquette is 1 by X or 2 by X
00130   unsigned int lengthOfColumnData  = (range.second-range.first)/nCols;
00131   unsigned int lengthOfAveragedDataInEachColumn = 2;  // we always only have two values per column averaged block 
00132   unsigned int numberOfDataBlocksToSkip = row / numberOfRowsToAverageOver_;
00133 
00134   const DecodingStructure & s = (const DecodingStructure & ) *(range.first+col*lengthOfColumnData + lengthOfAveragedDataInEachColumn*numberOfDataBlocksToSkip);
00135 
00136   if ((s.ped & 0xFF) == deadFlag_)
00137      isDeadColumn = true;
00138   else if ((s.ped & 0xFF) == noisyFlag_)
00139      isNoisyColumn = true;
00140 
00141   int maxRow = (lengthOfColumnData/lengthOfAveragedDataInEachColumn)*numberOfRowsToAverageOver_ - 1;
00142   if (col >= nCols || row > maxRow){
00143     throw cms::Exception("CorruptedData")
00144       << "[SiPixelGainCalibrationForHLT::getPed] Pixel out of range: col " << col << " row: " << row;
00145   }  
00146   return decodePed(s.ped & 0xFF);  
00147 }
00148 
00149 float SiPixelGainCalibrationForHLT::getGain(const int& col, const int& row, const Range& range, const int& nCols, bool& isDeadColumn, bool& isNoisyColumn) const {
00150 
00151   // determine what averaged data block we are in (there should be 1 or 2 of these depending on if plaquette is 1 by X or 2 by X
00152   unsigned int lengthOfColumnData  = (range.second-range.first)/nCols;
00153   unsigned int lengthOfAveragedDataInEachColumn = 2;  // we always only have two values per column averaged block 
00154   unsigned int numberOfDataBlocksToSkip = row / numberOfRowsToAverageOver_;
00155 
00156   const DecodingStructure & s = (const DecodingStructure & ) *(range.first+col*lengthOfColumnData + lengthOfAveragedDataInEachColumn*numberOfDataBlocksToSkip);
00157 
00158   if ((s.gain & 0xFF) == deadFlag_)
00159      isDeadColumn = true;
00160   else if ((s.gain & 0xFF) == noisyFlag_)
00161      isNoisyColumn = true;
00162      
00163   int maxRow = (lengthOfColumnData/lengthOfAveragedDataInEachColumn)*numberOfRowsToAverageOver_ - 1;
00164   if (col >= nCols || row > maxRow){
00165     throw cms::Exception("CorruptedData")
00166       << "[SiPixelGainCalibrationForHLT::getGain] Pixel out of range: col " << col << " row: " << row;
00167   }  
00168   return decodeGain(s.gain & 0xFF);  
00169 
00170 }
00171 
00172 float SiPixelGainCalibrationForHLT::encodeGain( const float& gain ) {
00173   
00174   if(gain < minGain_ || gain > maxGain_ ) {
00175     throw cms::Exception("InsertFailure")
00176       << "[SiPixelGainCalibrationForHLT::encodeGain] Trying to encode gain (" << gain << ") out of range [" << minGain_ << "," << maxGain_ << "]\n";
00177   } else {
00178     double precision   = (maxGain_-minGain_)/static_cast<float>(nBinsToUseForEncoding_);
00179     float  encodedGain = (float)((gain-minGain_)/precision);
00180     return encodedGain;
00181   }
00182 
00183 }
00184 
00185 float SiPixelGainCalibrationForHLT::encodePed( const float& ped ) {
00186 
00187   if(ped < minPed_ || ped > maxPed_ ) {
00188     throw cms::Exception("InsertFailure")
00189       << "[SiPixelGainCalibrationForHLT::encodePed] Trying to encode pedestal (" << ped << ") out of range [" << minPed_ << "," << maxPed_ << "]\n";
00190   } else {
00191     double precision   = (maxPed_-minPed_)/static_cast<float>(nBinsToUseForEncoding_);
00192     float  encodedPed = (float)((ped-minPed_)/precision);
00193     return encodedPed;
00194   }
00195 
00196 }
00197 
00198 float SiPixelGainCalibrationForHLT::decodePed( unsigned int ped ) const {
00199 
00200   double precision = (maxPed_-minPed_)/static_cast<float>(nBinsToUseForEncoding_);
00201   float decodedPed = (float)(ped*precision + minPed_);
00202   return decodedPed;
00203 
00204 }
00205 
00206 float SiPixelGainCalibrationForHLT::decodeGain( unsigned int gain ) const {
00207 
00208   double precision = (maxGain_-minGain_)/static_cast<float>(nBinsToUseForEncoding_);
00209   float decodedGain = (float)(gain*precision + minGain_);
00210   return decodedGain;
00211 
00212 }
00213 
00214