File indexing completed on 2024-04-06 12:30:52
0001 #include "SimPPS/PPSPixelDigiProducer/interface/RPixDummyROCSimulator.h"
0002 #include <vector>
0003 #include "TRandom.h"
0004 #include <iostream>
0005
0006 RPixDummyROCSimulator::RPixDummyROCSimulator(const edm::ParameterSet ¶ms, uint32_t det_id) : det_id_(det_id) {
0007 threshold_ = params.getParameter<double>("RPixDummyROCThreshold");
0008 electron_per_adc_ = params.getParameter<double>("RPixDummyROCElectronPerADC");
0009 VcaltoElectronGain_ = params.getParameter<int>("VCaltoElectronGain");
0010 VcaltoElectronOffset_ = params.getParameter<int>("VCaltoElectronOffset");
0011 doSingleCalibration_ = params.getParameter<bool>("doSingleCalibration");
0012 dead_pixel_probability_ = params.getParameter<double>("RPixDeadPixelProbability");
0013 dead_pixels_simulation_on_ = params.getParameter<bool>("RPixDeadPixelSimulationOn");
0014 verbosity_ = params.getParameter<int>("RPixVerbosity");
0015 links_persistence_ = params.getParameter<bool>("CTPPSPixelDigiSimHitRelationsPersistence");
0016 }
0017
0018 void RPixDummyROCSimulator::ConvertChargeToHits(
0019 const std::map<unsigned short, double> &signals,
0020 std::map<unsigned short, std::vector<std::pair<int, double> > > &theSignalProvenance,
0021 std::vector<CTPPSPixelDigi> &output_digi,
0022 std::vector<std::vector<std::pair<int, double> > > &output_digi_links,
0023 const CTPPSPixelGainCalibrations *pcalibrations) {
0024 for (std::map<unsigned short, double>::const_iterator i = signals.begin(); i != signals.end(); ++i) {
0025
0026 unsigned short pixel_no = i->first;
0027 if (verbosity_)
0028 edm::LogInfo("PPS") << "RPixDummyROCSimulator "
0029 << "Dummy ROC adc and threshold : " << i->second << ", " << threshold_;
0030 if (i->second > threshold_ && (!dead_pixels_simulation_on_ || dead_pixels_.find(pixel_no) == dead_pixels_.end())) {
0031 float gain = 0;
0032 float pedestal = 0;
0033 int adc = 0;
0034 uint32_t col = pixel_no / 160;
0035 uint32_t row = pixel_no % 160;
0036
0037 const CTPPSPixelGainCalibration &DetCalibs = pcalibrations->getGainCalibration(det_id_);
0038
0039
0040 if (col >= DetCalibs.getNCols())
0041 continue;
0042
0043 if (doSingleCalibration_) {
0044 adc = int(round(i->second / electron_per_adc_));
0045 } else {
0046 if (DetCalibs.getDetId() != 0) {
0047 gain = DetCalibs.getGain(col, row) * highRangeCal_ / lowRangeCal_;
0048 pedestal = DetCalibs.getPed(col, row);
0049 adc = int(round((i->second - VcaltoElectronOffset_) / (gain * VcaltoElectronGain_) + pedestal));
0050 }
0051 }
0052
0053 if (adc >= maxADC_)
0054 adc = maxADC_;
0055 if (adc < 0)
0056 adc = 0;
0057 output_digi.push_back(CTPPSPixelDigi(row, col, adc));
0058 if (links_persistence_) {
0059 output_digi_links.push_back(theSignalProvenance[pixel_no]);
0060 if (verbosity_) {
0061 edm::LogInfo("PPS") << "RPixDummyROCSimulator "
0062 << "digi links size=" << theSignalProvenance[pixel_no].size();
0063 for (unsigned int u = 0; u < theSignalProvenance[pixel_no].size(); ++u) {
0064 edm::LogInfo("PPS") << "RPixDummyROCSimulator "
0065 << " digi: particle=" << theSignalProvenance[pixel_no][u].first
0066 << " energy [electrons]=" << theSignalProvenance[pixel_no][u].second;
0067 }
0068 }
0069 }
0070 }
0071 }
0072
0073 if (verbosity_) {
0074 for (unsigned int i = 0; i < output_digi.size(); ++i) {
0075 edm::LogInfo("RPixDummyROCSimulator")
0076 << "Dummy ROC Simulator " << det_id_ << " row= "
0077 << output_digi[i].row() << " col= " << output_digi[i].column() << " adc= " << output_digi[i].adc();
0078 }
0079 }
0080 }