CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes | Static Private Attributes
RPixDummyROCSimulator Class Reference

#include <RPixDummyROCSimulator.h>

Public Member Functions

void ConvertChargeToHits (const std::map< unsigned short, double > &signals, std::map< unsigned short, std::vector< std::pair< int, double > > > &theSignalProvenance, std::vector< CTPPSPixelDigi > &output_digi, std::vector< std::vector< std::pair< int, double > > > &output_digi_links, const CTPPSPixelGainCalibrations *pcalibration)
 
 RPixDummyROCSimulator (const edm::ParameterSet &params, uint32_t det_id)
 

Private Types

typedef std::set< unsigned short > dead_pixel_set
 

Private Attributes

double dead_pixel_probability_
 
dead_pixel_set dead_pixels_
 
bool dead_pixels_simulation_on_
 
uint32_t det_id_
 
bool doSingleCalibration_
 
double electron_per_adc_
 
bool links_persistence_
 
double threshold_
 
int VcaltoElectronGain_
 
int VcaltoElectronOffset_
 
int verbosity_
 

Static Private Attributes

static constexpr double highRangeCal_ = 1800.
 
static constexpr double lowRangeCal_ = 260.
 
static constexpr int maxADC_ = 255
 

Detailed Description

Definition at line 13 of file RPixDummyROCSimulator.h.

Member Typedef Documentation

◆ dead_pixel_set

typedef std::set<unsigned short> RPixDummyROCSimulator::dead_pixel_set
private

Definition at line 24 of file RPixDummyROCSimulator.h.

Constructor & Destructor Documentation

◆ RPixDummyROCSimulator()

RPixDummyROCSimulator::RPixDummyROCSimulator ( const edm::ParameterSet params,
uint32_t  det_id 
)

Definition at line 6 of file RPixDummyROCSimulator.cc.

References dead_pixel_probability_, dead_pixels_simulation_on_, doSingleCalibration_, electron_per_adc_, links_persistence_, submitPVValidationJobs::params, threshold_, VcaltoElectronGain_, VcaltoElectronOffset_, and verbosity_.

6  : det_id_(det_id) {
7  threshold_ = params.getParameter<double>("RPixDummyROCThreshold");
8  electron_per_adc_ = params.getParameter<double>("RPixDummyROCElectronPerADC");
9  VcaltoElectronGain_ = params.getParameter<int>("VCaltoElectronGain");
10  VcaltoElectronOffset_ = params.getParameter<int>("VCaltoElectronOffset");
11  doSingleCalibration_ = params.getParameter<bool>("doSingleCalibration");
12  dead_pixel_probability_ = params.getParameter<double>("RPixDeadPixelProbability");
13  dead_pixels_simulation_on_ = params.getParameter<bool>("RPixDeadPixelSimulationOn");
14  verbosity_ = params.getParameter<int>("RPixVerbosity");
15  links_persistence_ = params.getParameter<bool>("CTPPSPixelDigiSimHitRelationsPersistence");
16 }

Member Function Documentation

◆ ConvertChargeToHits()

void RPixDummyROCSimulator::ConvertChargeToHits ( const std::map< unsigned short, double > &  signals,
std::map< unsigned short, std::vector< std::pair< int, double > > > &  theSignalProvenance,
std::vector< CTPPSPixelDigi > &  output_digi,
std::vector< std::vector< std::pair< int, double > > > &  output_digi_links,
const CTPPSPixelGainCalibrations pcalibration 
)

set maximum for 8 bits adc

Definition at line 18 of file RPixDummyROCSimulator.cc.

References gpuClustering::adc, cuy::col, dead_pixels_, dead_pixels_simulation_on_, det_id_, doSingleCalibration_, electron_per_adc_, PedestalClient_cfi::gain, CTPPSPixelGainCalibration::getDetId(), CTPPSPixelGainCalibration::getGain(), CTPPSPixelGainCalibrations::getGainCalibration(), CTPPSPixelGainCalibration::getNCols(), CTPPSPixelGainCalibration::getPed(), highRangeCal_, mps_fire::i, createfilelist::int, links_persistence_, lowRangeCal_, maxADC_, EcalCondDBWriter_cfi::pedestal, threshold_, VcaltoElectronGain_, VcaltoElectronOffset_, and verbosity_.

23  {
24  for (std::map<unsigned short, double>::const_iterator i = signals.begin(); i != signals.end(); ++i) {
25  //one threshold per hybrid
26  unsigned short pixel_no = i->first;
27  if (verbosity_)
28  edm::LogInfo("PPS") << "RPixDummyROCSimulator "
29  << "Dummy ROC adc and threshold : " << i->second << ", " << threshold_;
30  if (i->second > threshold_ && (!dead_pixels_simulation_on_ || dead_pixels_.find(pixel_no) == dead_pixels_.end())) {
31  float gain = 0;
32  float pedestal = 0;
33  int adc = 0;
34  uint32_t col = pixel_no / 160;
35  uint32_t row = pixel_no % 160;
36 
37  const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(det_id_);
38 
39  // Avoid exception due to col > 103 in case of 2x2 plane. To be removed
40  if (col >= DetCalibs.getNCols())
41  continue;
42 
44  adc = int(round(i->second / electron_per_adc_));
45  } else {
46  if (DetCalibs.getDetId() != 0) {
47  gain = DetCalibs.getGain(col, row) * highRangeCal_ / lowRangeCal_; // *highRangeCal/lowRangeCal
48  pedestal = DetCalibs.getPed(col, row);
49  adc = int(round((i->second - VcaltoElectronOffset_) / (gain * VcaltoElectronGain_) + pedestal));
50  }
51  }
53  if (adc >= maxADC_)
54  adc = maxADC_;
55  if (adc < 0)
56  adc = 0;
57  output_digi.push_back(CTPPSPixelDigi(row, col, adc));
58  if (links_persistence_) {
59  output_digi_links.push_back(theSignalProvenance[pixel_no]);
60  if (verbosity_) {
61  edm::LogInfo("PPS") << "RPixDummyROCSimulator "
62  << "digi links size=" << theSignalProvenance[pixel_no].size();
63  for (unsigned int u = 0; u < theSignalProvenance[pixel_no].size(); ++u) {
64  edm::LogInfo("PPS") << "RPixDummyROCSimulator "
65  << " digi: particle=" << theSignalProvenance[pixel_no][u].first
66  << " energy [electrons]=" << theSignalProvenance[pixel_no][u].second;
67  }
68  }
69  }
70  }
71  }
72 
73  if (verbosity_) {
74  for (unsigned int i = 0; i < output_digi.size(); ++i) {
75  edm::LogInfo("RPixDummyROCSimulator")
76  << "Dummy ROC Simulator " << det_id_ << " row= " //output_digi[i].GetDetId()<<" "
77  << output_digi[i].row() << " col= " << output_digi[i].column() << " adc= " << output_digi[i].adc();
78  }
79  }
80 }
float getPed(const int &col, const int &row) const
static constexpr int maxADC_
static constexpr double highRangeCal_
float getGain(const int &col, const int &row) const
Log< level::Info, false > LogInfo
col
Definition: cuy.py:1009
static constexpr double lowRangeCal_
uint16_t *__restrict__ uint16_t const *__restrict__ adc

Member Data Documentation

◆ dead_pixel_probability_

double RPixDummyROCSimulator::dead_pixel_probability_
private

Definition at line 30 of file RPixDummyROCSimulator.h.

Referenced by RPixDummyROCSimulator().

◆ dead_pixels_

dead_pixel_set RPixDummyROCSimulator::dead_pixels_
private

Definition at line 32 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits().

◆ dead_pixels_simulation_on_

bool RPixDummyROCSimulator::dead_pixels_simulation_on_
private

Definition at line 31 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ det_id_

uint32_t RPixDummyROCSimulator::det_id_
private

Definition at line 29 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits().

◆ doSingleCalibration_

bool RPixDummyROCSimulator::doSingleCalibration_
private

Definition at line 38 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ electron_per_adc_

double RPixDummyROCSimulator::electron_per_adc_
private

Definition at line 35 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ highRangeCal_

constexpr double RPixDummyROCSimulator::highRangeCal_ = 1800.
staticprivate

Definition at line 25 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits().

◆ links_persistence_

bool RPixDummyROCSimulator::links_persistence_
private

Definition at line 39 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ lowRangeCal_

constexpr double RPixDummyROCSimulator::lowRangeCal_ = 260.
staticprivate

Definition at line 26 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits().

◆ maxADC_

constexpr int RPixDummyROCSimulator::maxADC_ = 255
staticprivate

Definition at line 27 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits().

◆ threshold_

double RPixDummyROCSimulator::threshold_
private

Definition at line 34 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ VcaltoElectronGain_

int RPixDummyROCSimulator::VcaltoElectronGain_
private

Definition at line 36 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ VcaltoElectronOffset_

int RPixDummyROCSimulator::VcaltoElectronOffset_
private

Definition at line 37 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().

◆ verbosity_

int RPixDummyROCSimulator::verbosity_
private

Definition at line 33 of file RPixDummyROCSimulator.h.

Referenced by ConvertChargeToHits(), and RPixDummyROCSimulator().