CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
RPixDetClusterizer Class Reference

#include <RPixDetClusterizer.h>

Public Member Functions

void buildClusters (unsigned int detId, const std::vector< CTPPSPixelDigi > &digi, std::vector< CTPPSPixelCluster > &clusters, const CTPPSPixelGainCalibrations *pcalibration, const CTPPSPixelAnalysisMask *mask)
 
int calibrate (unsigned int, int, int, int, const CTPPSPixelGainCalibrations *pcalibration)
 
void make_cluster (RPixCalibDigi const &aSeed, std::vector< CTPPSPixelCluster > &clusters)
 
 RPixDetClusterizer (edm::ParameterSet const &conf)
 
 ~RPixDetClusterizer ()
 

Private Attributes

unsigned short ADCThreshold_
 
std::map< unsigned int, RPixCalibDigicalib_rpix_digi_map_
 
std::string CalibrationFile_
 
bool doSingleCalibration_
 
double ElectronADCGain_
 
const edm::ParameterSetparams_
 
std::set< CTPPSPixelDigirpix_digi_set_
 
unsigned short SeedADCThreshold_
 
std::vector< RPixCalibDigiSeedVector_
 
int VcaltoElectronGain_
 
int VcaltoElectronOffset_
 
int verbosity_
 

Detailed Description

Definition at line 40 of file RPixDetClusterizer.h.

Constructor & Destructor Documentation

◆ RPixDetClusterizer()

RPixDetClusterizer::RPixDetClusterizer ( edm::ParameterSet const &  conf)

Definition at line 16 of file RPixDetClusterizer.cc.

References ADCThreshold_, doSingleCalibration_, ElectronADCGain_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), SeedADCThreshold_, VcaltoElectronGain_, VcaltoElectronOffset_, and verbosity_.

16  : params_(conf), SeedVector_(0) {
17  verbosity_ = conf.getUntrackedParameter<int>("RPixVerbosity");
18  SeedADCThreshold_ = conf.getParameter<int>("SeedADCThreshold");
19  ADCThreshold_ = conf.getParameter<int>("ADCThreshold");
20  ElectronADCGain_ = conf.getParameter<double>("ElectronADCGain");
21  VcaltoElectronGain_ = conf.getParameter<int>("VCaltoElectronGain");
22  VcaltoElectronOffset_ = conf.getParameter<int>("VCaltoElectronOffset");
23  doSingleCalibration_ = conf.getParameter<bool>("doSingleCalibration");
24 }
std::vector< RPixCalibDigi > SeedVector_
const edm::ParameterSet & params_
unsigned short SeedADCThreshold_
unsigned short ADCThreshold_

◆ ~RPixDetClusterizer()

RPixDetClusterizer::~RPixDetClusterizer ( )

Definition at line 26 of file RPixDetClusterizer.cc.

26 {}

Member Function Documentation

◆ buildClusters()

void RPixDetClusterizer::buildClusters ( unsigned int  detId,
const std::vector< CTPPSPixelDigi > &  digi,
std::vector< CTPPSPixelCluster > &  clusters,
const CTPPSPixelGainCalibrations pcalibration,
const CTPPSPixelAnalysisMask mask 
)

Definition at line 28 of file RPixDetClusterizer.cc.

References gpuClustering::adc, CTPPSPixelAnalysisMask::analysisMask, calib_rpix_digi_map_, calibrate(), bsc_activity_cfg::clusters, ElectronADCGain_, pwdgSkimBPark_cfi::electrons, Exception, make_cluster(), gpuClustering::pixelStatus::mask, muonClassificationByHits_cfi::pixel, rpix_digi_set_, SeedADCThreshold_, SeedVector_, CTPPSPixelIndices::transformToModule(), and verbosity_.

Referenced by CTPPSPixelClusterProducer::run().

32  {
33  std::map<uint32_t, CTPPSPixelROCAnalysisMask> const &mask = maskera->analysisMask;
34  std::set<std::pair<unsigned char, unsigned char> > maskedPixels;
35 
36  // read and store masked pixels after converting ROC numbering into module numbering
38  for (auto const &det : mask) {
39  uint32_t planeId = det.first & ~rocMask;
40 
41  if (planeId == detId) {
42  unsigned int rocNum = (det.first & rocMask) >> rocOffset;
43  if (rocNum > 5)
44  throw cms::Exception("InvalidRocNumber") << "roc number from mask: " << rocNum;
45  for (auto const &paio : det.second.maskedPixels) {
46  std::pair<unsigned char, unsigned char> pa = paio;
47  int modCol, modRow;
48  pI.transformToModule(pa.second, pa.first, rocNum, modCol, modRow);
49  std::pair<int, int> modPix(modRow, modCol);
50  maskedPixels.insert(modPix);
51  }
52  }
53  }
54 
55  if (verbosity_)
56  edm::LogInfo("RPixDetClusterizer") << detId << " received digi.size()=" << digi.size();
57 
58  // creating a set of CTPPSPixelDigi's and fill it
59 
60  rpix_digi_set_.clear();
61  rpix_digi_set_.insert(digi.begin(), digi.end());
62  SeedVector_.clear();
63 
64  // calibrate digis here
65  calib_rpix_digi_map_.clear();
66 
67  for (auto const &RPdit : rpix_digi_set_) {
68  uint8_t row = RPdit.row();
69  uint8_t column = RPdit.column();
70  if (row > maxRow || column > maxCol)
71  throw cms::Exception("CTPPSDigiOutofRange") << " row = " << row << " column = " << column;
72 
73  std::pair<unsigned char, unsigned char> pixel = std::make_pair(row, column);
74  unsigned short adc = RPdit.adc();
75  unsigned short electrons = 0;
76 
77  // check if pixel is noisy or dead (i.e. in the mask)
78  const bool is_in = maskedPixels.find(pixel) != maskedPixels.end();
79  if (!is_in) {
80  //calibrate digi and store the new ones
81  electrons = calibrate(detId, adc, row, column, pcalibrations);
84  RPixCalibDigi calibDigi(row, column, adc, electrons);
85  unsigned int index = column * maxRow + row;
86  calib_rpix_digi_map_.insert(std::pair<unsigned int, RPixCalibDigi>(index, calibDigi));
87  SeedVector_.push_back(calibDigi);
88  }
89  }
90  if (verbosity_)
91  edm::LogInfo("RPixDetClusterizer") << " RPix set size = " << calib_rpix_digi_map_.size();
92 
93  for (auto SeedIt : SeedVector_) {
94  make_cluster(SeedIt, clusters);
95  }
96 }
std::map< unsigned int, RPixCalibDigi > calib_rpix_digi_map_
constexpr uint32_t mask
Definition: gpuClustering.h:24
std::vector< RPixCalibDigi > SeedVector_
void make_cluster(RPixCalibDigi const &aSeed, std::vector< CTPPSPixelCluster > &clusters)
std::set< CTPPSPixelDigi > rpix_digi_set_
unsigned short SeedADCThreshold_
Log< level::Info, false > LogInfo
int calibrate(unsigned int, int, int, int, const CTPPSPixelGainCalibrations *pcalibration)
int transformToModule(const int colROC, const int rowROC, const int rocId, int &col, int &row) const
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ calibrate()

int RPixDetClusterizer::calibrate ( unsigned int  detId,
int  adc,
int  row,
int  col,
const CTPPSPixelGainCalibrations pcalibration 
)

Definition at line 137 of file RPixDetClusterizer.cc.

References gpuClustering::adc, cuy::col, doSingleCalibration_, ElectronADCGain_, pwdgSkimBPark_cfi::electrons, PedestalClient_cfi::gain, CTPPSPixelGainCalibration::getDetId(), CTPPSPixelGainCalibration::getGain(), CTPPSPixelGainCalibrations::getGainCalibration(), CTPPSPixelGainCalibration::getPed(), createfilelist::int, EcalCondDBWriter_cfi::pedestal, VcaltoElectronGain_, VcaltoElectronOffset_, and verbosity_.

Referenced by buildClusters().

138  {
139  float gain = 0;
140  float pedestal = 0;
141  int electrons = 0;
142 
143  if (!doSingleCalibration_) {
144  const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(detId);
145 
146  if (DetCalibs.getDetId() != 0) {
147  gain = DetCalibs.getGain(col, row) * highRangeCal / lowRangeCal;
148  pedestal = DetCalibs.getPed(col, row);
149  float vcal = (adc - pedestal) * gain;
151  } else {
153  pedestal = 0;
154  electrons = int(adc * gain - pedestal);
155  }
156 
157  if (verbosity_)
158  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
159  << " gain = " << gain << " pedestal = " << pedestal;
160  } else {
162  pedestal = 0;
163  electrons = int(adc * gain - pedestal);
164  if (verbosity_)
165  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
166  << " ElectronADCGain = " << gain << " pedestal = " << pedestal;
167  }
168  if (electrons < 0) {
169  edm::LogInfo("RPixCalibration") << "RPixDetClusterizer::calibrate: *** electrons < 0 *** --> " << electrons
170  << " --> electrons = 0";
171  electrons = 0;
172  }
173 
174  return electrons;
175 }
float getPed(const int &col, const int &row) const
float getGain(const int &col, const int &row) const
Log< level::Info, false > LogInfo
col
Definition: cuy.py:1009
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ make_cluster()

void RPixDetClusterizer::make_cluster ( RPixCalibDigi const &  aSeed,
std::vector< CTPPSPixelCluster > &  clusters 
)

check if seed already used

Definition at line 98 of file RPixDetClusterizer.cc.

References RPixTempCluster::adc, RPixTempCluster::addPixel(), c, calib_rpix_digi_map_, bsc_activity_cfg::clusters, RPixTempCluster::col, CTPPSPixelDigi::column(), RPixCalibDigi::electrons(), RPixTempCluster::empty(), RPixTempCluster::isize, SiStripPI::max, SiStripPI::min, RPixTempCluster::pop(), alignCSCRings::r, CTPPSPixelDigi::row(), RPixTempCluster::row, and RPixTempCluster::top().

Referenced by buildClusters().

98  {
100  unsigned int seedIndex = aSeed.column() * maxRow + aSeed.row();
101  if (calib_rpix_digi_map_.find(seedIndex) == calib_rpix_digi_map_.end()) {
102  return;
103  }
104  // creating a temporary cluster
105  RPixTempCluster atempCluster;
106 
107  // filling the cluster with the seed
108  atempCluster.addPixel(aSeed.row(), aSeed.column(), aSeed.electrons());
109  calib_rpix_digi_map_.erase(seedIndex);
110 
111  while (!atempCluster.empty()) {
112  //This is the standard algorithm to find and add a pixel
113  auto curInd = atempCluster.top();
114  atempCluster.pop();
115  for (auto c = std::max(0, int(atempCluster.col[curInd]) - 1);
116  c < std::min(int(atempCluster.col[curInd]) + 2, maxCol);
117  ++c) {
118  for (auto r = std::max(0, int(atempCluster.row[curInd]) - 1);
119  r < std::min(int(atempCluster.row[curInd]) + 2, maxRow);
120  ++r) {
121  unsigned int currIndex = c * maxRow + r;
122  if (calib_rpix_digi_map_.find(currIndex) != calib_rpix_digi_map_.end()) {
123  if (!atempCluster.addPixel(r, c, calib_rpix_digi_map_[currIndex].electrons())) {
124  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
125  return;
126  }
127  calib_rpix_digi_map_.erase(currIndex);
128  }
129  }
130  }
131 
132  } // while accretion
133 
134  clusters.emplace_back(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
135 }
unsigned short adc[MAXSIZE]
std::map< unsigned int, RPixCalibDigi > calib_rpix_digi_map_
unsigned short isize
uint8_t row[MAXSIZE]
bool addPixel(unsigned char myrow, unsigned char mycol, unsigned short const iadc)
uint8_t col[MAXSIZE]
unsigned short top() const

Member Data Documentation

◆ ADCThreshold_

unsigned short RPixDetClusterizer::ADCThreshold_
private

Definition at line 60 of file RPixDetClusterizer.h.

Referenced by RPixDetClusterizer().

◆ calib_rpix_digi_map_

std::map<unsigned int, RPixCalibDigi> RPixDetClusterizer::calib_rpix_digi_map_
private

Definition at line 56 of file RPixDetClusterizer.h.

Referenced by buildClusters(), and make_cluster().

◆ CalibrationFile_

std::string RPixDetClusterizer::CalibrationFile_
private

Definition at line 65 of file RPixDetClusterizer.h.

◆ doSingleCalibration_

bool RPixDetClusterizer::doSingleCalibration_
private

Definition at line 64 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

◆ ElectronADCGain_

double RPixDetClusterizer::ElectronADCGain_
private

Definition at line 61 of file RPixDetClusterizer.h.

Referenced by buildClusters(), calibrate(), and RPixDetClusterizer().

◆ params_

const edm::ParameterSet& RPixDetClusterizer::params_
private

Definition at line 57 of file RPixDetClusterizer.h.

◆ rpix_digi_set_

std::set<CTPPSPixelDigi> RPixDetClusterizer::rpix_digi_set_
private

Definition at line 55 of file RPixDetClusterizer.h.

Referenced by buildClusters().

◆ SeedADCThreshold_

unsigned short RPixDetClusterizer::SeedADCThreshold_
private

Definition at line 59 of file RPixDetClusterizer.h.

Referenced by buildClusters(), and RPixDetClusterizer().

◆ SeedVector_

std::vector<RPixCalibDigi> RPixDetClusterizer::SeedVector_
private

Definition at line 66 of file RPixDetClusterizer.h.

Referenced by buildClusters().

◆ VcaltoElectronGain_

int RPixDetClusterizer::VcaltoElectronGain_
private

Definition at line 62 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

◆ VcaltoElectronOffset_

int RPixDetClusterizer::VcaltoElectronOffset_
private

Definition at line 63 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

◆ verbosity_

int RPixDetClusterizer::verbosity_
private

Definition at line 58 of file RPixDetClusterizer.h.

Referenced by buildClusters(), calibrate(), and RPixDetClusterizer().