Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:40

0001 #include <iostream>
0002 #include "FWCore/Utilities/interface/Exception.h"
0003 
0004 #include "RecoPPS/Local/interface/RPixDetClusterizer.h"
0005 #include "DataFormats/CTPPSDetId/interface/CTPPSPixelDetId.h"
0006 
0007 namespace {
0008   constexpr int maxCol = CTPPSPixelCluster::MAXCOL;
0009   constexpr int maxRow = CTPPSPixelCluster::MAXROW;
0010   constexpr double highRangeCal = 1800.;
0011   constexpr double lowRangeCal = 260.;
0012   constexpr int rocMask = 0xE000;
0013   constexpr int rocOffset = 13;
0014 }  // namespace
0015 
0016 RPixDetClusterizer::RPixDetClusterizer(edm::ParameterSet const &conf) : SeedVector_(0) {
0017   verbosity_ = conf.getUntrackedParameter<int>("RPixVerbosity");
0018   SeedADCThreshold_ = conf.getParameter<int>("SeedADCThreshold");
0019   ADCThreshold_ = conf.getParameter<int>("ADCThreshold");
0020   ElectronADCGain_ = conf.getParameter<double>("ElectronADCGain");
0021   VcaltoElectronGain_ = conf.getParameter<int>("VCaltoElectronGain");
0022   VcaltoElectronOffset_ = conf.getParameter<int>("VCaltoElectronOffset");
0023   doSingleCalibration_ = conf.getParameter<bool>("doSingleCalibration");
0024 }
0025 
0026 RPixDetClusterizer::~RPixDetClusterizer() {}
0027 
0028 void RPixDetClusterizer::buildClusters(unsigned int detId,
0029                                        const std::vector<CTPPSPixelDigi> &digi,
0030                                        std::vector<CTPPSPixelCluster> &clusters,
0031                                        const CTPPSPixelGainCalibrations *pcalibrations,
0032                                        const CTPPSPixelAnalysisMask *maskera) {
0033   std::map<uint32_t, CTPPSPixelROCAnalysisMask> const &mask = maskera->analysisMask;
0034   std::set<std::pair<unsigned char, unsigned char> > maskedPixels;
0035 
0036   bool fullyMaskedRoc[6][6] = {{false}};
0037   CTPPSPixelDetId currentId(detId);
0038 
0039   // read and store masked pixels after converting ROC numbering into module numbering
0040   CTPPSPixelIndices pI;
0041   for (auto const &det : mask) {
0042     uint32_t planeId = det.first & ~rocMask;
0043 
0044     if (planeId == detId) {
0045       unsigned int rocNum = (det.first & rocMask) >> rocOffset;
0046       if (rocNum > 5)
0047         throw cms::Exception("InvalidRocNumber") << "roc number from mask: " << rocNum;
0048       if (det.second.fullMask)
0049         fullyMaskedRoc[currentId.plane()][rocNum] = true;
0050       for (auto const &paio : det.second.maskedPixels) {
0051         std::pair<unsigned char, unsigned char> pa = paio;
0052         int modCol, modRow;
0053         pI.transformToModule(pa.second, pa.first, rocNum, modCol, modRow);
0054         std::pair<int, int> modPix(modRow, modCol);
0055         maskedPixels.insert(modPix);
0056       }
0057     }
0058   }
0059 
0060   if (verbosity_)
0061     edm::LogInfo("RPixDetClusterizer") << detId << " received digi.size()=" << digi.size();
0062 
0063   // creating a set of CTPPSPixelDigi's and fill it
0064 
0065   rpix_digi_set_.clear();
0066   rpix_digi_set_.insert(digi.begin(), digi.end());
0067   SeedVector_.clear();
0068 
0069   // calibrate digis here
0070   calib_rpix_digi_map_.clear();
0071 
0072   for (auto const &RPdit : rpix_digi_set_) {
0073     uint8_t row = RPdit.row();
0074     uint8_t column = RPdit.column();
0075     if (row > maxRow || column > maxCol)
0076       throw cms::Exception("CTPPSDigiOutofRange") << " row = " << row << "  column = " << column;
0077 
0078     std::pair<unsigned char, unsigned char> pixel = std::make_pair(row, column);
0079     unsigned short adc = RPdit.adc();
0080     unsigned short electrons = 0;
0081 
0082     // check if pixel is noisy or dead (i.e. in the mask)
0083     const bool is_in = maskedPixels.find(pixel) != maskedPixels.end();
0084 
0085     // check if pixel inside a fully masked roc
0086     int piano = currentId.plane();
0087     int rocId = pI.getROCId(column, row);
0088 
0089     if (!is_in && !fullyMaskedRoc[piano][rocId]) {
0090       //calibrate digi and store the new ones
0091       electrons = calibrate(detId, adc, row, column, pcalibrations);
0092       if (electrons < SeedADCThreshold_ * ElectronADCGain_)
0093         electrons = SeedADCThreshold_ * ElectronADCGain_;
0094       RPixCalibDigi calibDigi(row, column, adc, electrons);
0095       unsigned int index = column * maxRow + row;
0096       calib_rpix_digi_map_.insert(std::pair<unsigned int, RPixCalibDigi>(index, calibDigi));
0097       SeedVector_.push_back(calibDigi);
0098     }
0099   }
0100   if (verbosity_)
0101     edm::LogInfo("RPixDetClusterizer") << " RPix set size = " << calib_rpix_digi_map_.size();
0102 
0103   for (auto SeedIt : SeedVector_) {
0104     make_cluster(SeedIt, clusters);
0105   }
0106 }
0107 
0108 void RPixDetClusterizer::make_cluster(RPixCalibDigi const &aSeed, std::vector<CTPPSPixelCluster> &clusters) {
0109   /// check if seed already used
0110   unsigned int seedIndex = aSeed.column() * maxRow + aSeed.row();
0111   if (calib_rpix_digi_map_.find(seedIndex) == calib_rpix_digi_map_.end()) {
0112     return;
0113   }
0114   // creating a temporary cluster
0115   RPixTempCluster atempCluster;
0116 
0117   // filling the cluster with the seed
0118   atempCluster.addPixel(aSeed.row(), aSeed.column(), aSeed.electrons());
0119   calib_rpix_digi_map_.erase(seedIndex);
0120 
0121   while (!atempCluster.empty()) {
0122     //This is the standard algorithm to find and add a pixel
0123     auto curInd = atempCluster.top();
0124     atempCluster.pop();
0125     for (auto c = std::max(0, int(atempCluster.col[curInd]) - 1);
0126          c < std::min(int(atempCluster.col[curInd]) + 2, maxCol);
0127          ++c) {
0128       for (auto r = std::max(0, int(atempCluster.row[curInd]) - 1);
0129            r < std::min(int(atempCluster.row[curInd]) + 2, maxRow);
0130            ++r) {
0131         unsigned int currIndex = c * maxRow + r;
0132         if (calib_rpix_digi_map_.find(currIndex) != calib_rpix_digi_map_.end()) {
0133           if (!atempCluster.addPixel(r, c, calib_rpix_digi_map_[currIndex].electrons())) {
0134             clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
0135             return;
0136           }
0137           calib_rpix_digi_map_.erase(currIndex);
0138         }
0139       }
0140     }
0141 
0142   }  // while accretion
0143 
0144   clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
0145 }
0146 
0147 int RPixDetClusterizer::calibrate(
0148     unsigned int detId, int adc, int row, int col, const CTPPSPixelGainCalibrations *pcalibrations) {
0149   float gain = 0;
0150   float pedestal = 0;
0151   int electrons = 0;
0152 
0153   if (!doSingleCalibration_) {
0154     const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(detId);
0155 
0156     if (DetCalibs.getDetId() != 0) {
0157       gain = DetCalibs.getGain(col, row) * highRangeCal / lowRangeCal;
0158       pedestal = DetCalibs.getPed(col, row);
0159       float vcal = (adc - pedestal) * gain;
0160       electrons = int(vcal * VcaltoElectronGain_ + VcaltoElectronOffset_);
0161     } else {
0162       gain = ElectronADCGain_;
0163       pedestal = 0;
0164       electrons = int(adc * gain - pedestal);
0165     }
0166 
0167     if (verbosity_)
0168       edm::LogInfo("RPixCalibration") << "calibrate:  adc = " << adc << " electrons = " << electrons
0169                                       << " gain = " << gain << " pedestal = " << pedestal;
0170   } else {
0171     gain = ElectronADCGain_;
0172     pedestal = 0;
0173     electrons = int(adc * gain - pedestal);
0174     if (verbosity_)
0175       edm::LogInfo("RPixCalibration") << "calibrate:  adc = " << adc << " electrons = " << electrons
0176                                       << " ElectronADCGain = " << gain << " pedestal = " << pedestal;
0177   }
0178   if (electrons < 0) {
0179     edm::LogInfo("RPixCalibration") << "RPixDetClusterizer::calibrate: *** electrons < 0 *** --> " << electrons
0180                                     << "  --> electrons = 0";
0181     electrons = 0;
0182   }
0183 
0184   return electrons;
0185 }