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 42 of file RPixDetClusterizer.h.

Constructor & Destructor Documentation

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

Definition at line 13 of file RPixDetClusterizer.cc.

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

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

Definition at line 23 of file RPixDetClusterizer.cc.

23 {}

Member Function Documentation

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

Definition at line 25 of file RPixDetClusterizer.cc.

References ecalMGPA::adc(), CTPPSPixelAnalysisMask::analysisMask, calib_rpix_digi_map_, calibrate(), ElectronADCGain_, pwdgSkimBPark_cfi::electrons, Exception, make_cluster(), muonClassificationByHits_cfi::pixel, rpix_digi_set_, SeedADCThreshold_, SeedVector_, and verbosity_.

Referenced by CTPPSPixelClusterProducer::run().

29  {
30  std::map<uint32_t, CTPPSPixelROCAnalysisMask> const &mask = maskera->analysisMask;
31  std::map<uint32_t, CTPPSPixelROCAnalysisMask>::const_iterator mask_it = mask.find(detId);
32 
33  std::set<std::pair<unsigned char, unsigned char> > maskedPixels;
34  if (mask_it != mask.end())
35  maskedPixels = mask_it->second.maskedPixels;
36 
37  if (verbosity_)
38  edm::LogInfo("RPixDetClusterizer") << detId << " received digi.size()=" << digi.size();
39 
40  // creating a set of CTPPSPixelDigi's and fill it
41 
42  rpix_digi_set_.clear();
43  rpix_digi_set_.insert(digi.begin(), digi.end());
44  SeedVector_.clear();
45 
46  // calibrate digis here
47  calib_rpix_digi_map_.clear();
48 
49  for (auto const &RPdit : rpix_digi_set_) {
50  uint8_t row = RPdit.row();
51  uint8_t column = RPdit.column();
52  if (row > maxRow || column > maxCol)
53  throw cms::Exception("CTPPSDigiOutofRange") << " row = " << row << " column = " << column;
54 
55  std::pair<unsigned char, unsigned char> pixel = std::make_pair(row, column);
56  unsigned short adc = RPdit.adc();
57  unsigned short electrons = 0;
58 
59  // check if pixel is noisy or dead (i.e. in the mask)
60  const bool is_in = maskedPixels.find(pixel) != maskedPixels.end();
61  if (!is_in) {
62  //calibrate digi and store the new ones
63  electrons = calibrate(detId, adc, row, column, pcalibrations);
64  if (electrons < SeedADCThreshold_ * ElectronADCGain_)
65  electrons = SeedADCThreshold_ * ElectronADCGain_;
66  RPixCalibDigi calibDigi(row, column, adc, electrons);
67  unsigned int index = column * maxRow + row;
68  calib_rpix_digi_map_.insert(std::pair<unsigned int, RPixCalibDigi>(index, calibDigi));
69  SeedVector_.push_back(calibDigi);
70  }
71  }
72  if (verbosity_)
73  edm::LogInfo("RPixDetClusterizer") << " RPix set size = " << calib_rpix_digi_map_.size();
74 
75  for (auto SeedIt : SeedVector_) {
76  make_cluster(SeedIt, clusters);
77  }
78 }
std::map< unsigned int, RPixCalibDigi > calib_rpix_digi_map_
std::map< uint32_t, CTPPSPixelROCAnalysisMask > analysisMask
std::vector< RPixCalibDigi > SeedVector_
void make_cluster(RPixCalibDigi const &aSeed, std::vector< CTPPSPixelCluster > &clusters)
std::set< CTPPSPixelDigi > rpix_digi_set_
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
unsigned short SeedADCThreshold_
int calibrate(unsigned int, int, int, int, const CTPPSPixelGainCalibrations *pcalibration)
int RPixDetClusterizer::calibrate ( unsigned int  detId,
int  adc,
int  row,
int  col,
const CTPPSPixelGainCalibrations pcalibration 
)

Definition at line 121 of file RPixDetClusterizer.cc.

References 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().

122  {
123  float gain = 0;
124  float pedestal = 0;
125  int electrons = 0;
126 
127  if (!doSingleCalibration_) {
128  const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(detId);
129 
130  if (DetCalibs.getDetId() != 0) {
131  gain = DetCalibs.getGain(col, row) * highRangeCal / lowRangeCal;
132  pedestal = DetCalibs.getPed(col, row);
133  float vcal = (adc - pedestal) * gain;
134  electrons = int(vcal * VcaltoElectronGain_ + VcaltoElectronOffset_);
135  } else {
136  gain = ElectronADCGain_;
137  pedestal = 0;
138  electrons = int(adc * gain - pedestal);
139  }
140 
141  if (verbosity_)
142  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
143  << " gain = " << gain << " pedestal = " << pedestal;
144  } else {
145  gain = ElectronADCGain_;
146  pedestal = 0;
147  electrons = int(adc * gain - pedestal);
148  if (verbosity_)
149  edm::LogInfo("RPixCalibration") << "calibrate: adc = " << adc << " electrons = " << electrons
150  << " ElectronADCGain = " << gain << " pedestal = " << pedestal;
151  }
152  if (electrons < 0) {
153  edm::LogInfo("RPixCalibration") << "RPixDetClusterizer::calibrate: *** electrons < 0 *** --> " << electrons
154  << " --> electrons = 0";
155  electrons = 0;
156  }
157 
158  return electrons;
159 }
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
float getPed(const int &col, const int &row) const
col
Definition: cuy.py:1010
float getGain(const int &col, const int &row) const
void RPixDetClusterizer::make_cluster ( RPixCalibDigi const &  aSeed,
std::vector< CTPPSPixelCluster > &  clusters 
)

check if seed already used

Definition at line 80 of file RPixDetClusterizer.cc.

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

Referenced by buildClusters().

80  {
82  unsigned int seedIndex = aSeed.column() * maxRow + aSeed.row();
83  if (calib_rpix_digi_map_.find(seedIndex) == calib_rpix_digi_map_.end()) {
84  return;
85  }
86  // creating a temporary cluster
87  RPixTempCluster atempCluster;
88 
89  // filling the cluster with the seed
90  atempCluster.addPixel(aSeed.row(), aSeed.column(), aSeed.electrons());
91  calib_rpix_digi_map_.erase(seedIndex);
92 
93  while (!atempCluster.empty()) {
94  //This is the standard algorithm to find and add a pixel
95  auto curInd = atempCluster.top();
96  atempCluster.pop();
97  for (auto c = std::max(0, int(atempCluster.col[curInd]) - 1);
98  c < std::min(int(atempCluster.col[curInd]) + 2, maxCol);
99  ++c) {
100  for (auto r = std::max(0, int(atempCluster.row[curInd]) - 1);
101  r < std::min(int(atempCluster.row[curInd]) + 2, maxRow);
102  ++r) {
103  unsigned int currIndex = c * maxRow + r;
104  if (calib_rpix_digi_map_.find(currIndex) != calib_rpix_digi_map_.end()) {
105  if (!atempCluster.addPixel(r, c, calib_rpix_digi_map_[currIndex].electrons())) {
106  CTPPSPixelCluster acluster(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
107  clusters.push_back(acluster);
108  return;
109  }
110  calib_rpix_digi_map_.erase(currIndex);
111  }
112  }
113  }
114 
115  } // while accretion
116 
117  CTPPSPixelCluster cluster(atempCluster.isize, atempCluster.adc, atempCluster.row, atempCluster.col);
118  clusters.push_back(cluster);
119 }
unsigned short adc[MAXSIZE]
std::map< unsigned int, RPixCalibDigi > calib_rpix_digi_map_
unsigned short isize
uint8_t row[MAXSIZE]
T min(T a, T b)
Definition: MathUtil.h:58
bool addPixel(unsigned char myrow, unsigned char mycol, unsigned short const iadc)
uint8_t col[MAXSIZE]
unsigned short top() const

Member Data Documentation

unsigned short RPixDetClusterizer::ADCThreshold_
private

Definition at line 62 of file RPixDetClusterizer.h.

Referenced by RPixDetClusterizer().

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

Definition at line 58 of file RPixDetClusterizer.h.

Referenced by buildClusters(), and make_cluster().

std::string RPixDetClusterizer::CalibrationFile_
private

Definition at line 67 of file RPixDetClusterizer.h.

bool RPixDetClusterizer::doSingleCalibration_
private

Definition at line 66 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

double RPixDetClusterizer::ElectronADCGain_
private

Definition at line 63 of file RPixDetClusterizer.h.

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

const edm::ParameterSet& RPixDetClusterizer::params_
private

Definition at line 59 of file RPixDetClusterizer.h.

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

Definition at line 57 of file RPixDetClusterizer.h.

Referenced by buildClusters().

unsigned short RPixDetClusterizer::SeedADCThreshold_
private

Definition at line 61 of file RPixDetClusterizer.h.

Referenced by buildClusters(), and RPixDetClusterizer().

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

Definition at line 68 of file RPixDetClusterizer.h.

Referenced by buildClusters().

int RPixDetClusterizer::VcaltoElectronGain_
private

Definition at line 64 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

int RPixDetClusterizer::VcaltoElectronOffset_
private

Definition at line 65 of file RPixDetClusterizer.h.

Referenced by calibrate(), and RPixDetClusterizer().

int RPixDetClusterizer::verbosity_
private

Definition at line 60 of file RPixDetClusterizer.h.

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