Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Implimentation of Dual-Slope signal scaling for the Inner Pixel #24325

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
Expand Up @@ -82,7 +82,11 @@ Phase2TrackerDigitizerAlgorithm::Phase2TrackerDigitizerAlgorithm(const edm::Para

ClusterWidth(conf_specific.getParameter<double>("ClusterWidth")), // Charge integration spread on the collection plane

doDigitalReadout(conf_specific.getParameter<bool>("DigitalReadout")), // Flag to decide analog or digital readout
// Allowed modes of readout which has following values :
// 0 ---> Digital or binary readout
// -1 ---> Analog readout, current digitizer (Inner Pixel) (TDR version) with no threshold subtraction
// Analog readout with dual slope with the "second" slope being 1/2^(n-1) and threshold subtraction (n = 1, 2, 3,4)
thePhase2ReadoutMode(conf_specific.getParameter<int>("Phase2ReadoutMode")),

// ADC calibration 1adc count(135e.
// Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
Expand Down Expand Up @@ -161,7 +165,8 @@ Phase2TrackerDigitizerAlgorithm::Phase2TrackerDigitizerAlgorithm(const edm::Para
<< theThresholdInE_Endcap
<< "\nthreshold in electron Barrel = "
<< theThresholdInE_Barrel
<< " " << theElectronPerADC << " " << theAdcFullScale
<< " ElectronPerADC " << theElectronPerADC
<< " ADC Scale (in bits) " << theAdcFullScale
<< " The delta cut-off is set to " << tMax
<< " pix-inefficiency " << AddPixelInefficiency;
}
Expand Down Expand Up @@ -972,10 +977,9 @@ void Phase2TrackerDigitizerAlgorithm::digitize(const Phase2TrackerGeomDetUnit* p
// DigitizerUtility::Amplitude sig_data = s.second;
const DigitizerUtility::Amplitude& sig_data = s.second;
float signalInElectrons = sig_data.ampl();
int adc;
unsigned short adc;
if (signalInElectrons >= theThresholdInE) { // check threshold
if (doDigitalReadout) adc = theAdcFullScale;
else adc = std::min( int(signalInElectrons / theElectronPerADC), theAdcFullScale );
adc = convertSignalToAdc(detID, signalInElectrons, theThresholdInE);
DigitizerUtility::DigiSimInfo info;
info.sig_tot = adc;
info.ot_bit = ( signalInElectrons > theHIPThresholdInE ? true : false);
Expand All @@ -989,3 +993,33 @@ void Phase2TrackerDigitizerAlgorithm::digitize(const Phase2TrackerGeomDetUnit* p
}
}
}
//
// Scale the Signal using Dual Slope option
//
int Phase2TrackerDigitizerAlgorithm::convertSignalToAdc(uint32_t detID, float signal_in_elec,float threshold) {
int signal_in_adc;
int temp_signal;
const int max_limit = 10;
if (thePhase2ReadoutMode == 0) signal_in_adc = theAdcFullScale;
else {

if (thePhase2ReadoutMode == -1) {
temp_signal = std::min( int(signal_in_elec / theElectronPerADC), theAdcFullScale );
} else {
// calculate the kink point and the slope
const int dualslope_param = std::min(abs(thePhase2ReadoutMode), max_limit);
const int kink_point = int(theAdcFullScale/2) +1;
temp_signal = std::floor((signal_in_elec-threshold) / theElectronPerADC) + 1;
if ( temp_signal > kink_point) temp_signal = std::floor((temp_signal - kink_point)/(pow(2, dualslope_param-1))) + kink_point;
}
signal_in_adc = std::min(temp_signal, theAdcFullScale );
LogInfo("Phase2TrackerDigitizerAlgorithm") << " DetId " << detID
<< " signal_in_elec " << signal_in_elec
<< " threshold " << threshold << " signal_above_thr "
<< (signal_in_elec-threshold)
<< " temp conversion " << std::floor((signal_in_elec-threshold)/theElectronPerADC)+1
<< " signal after slope correction " << temp_signal
<< " signal_in_adc " << signal_in_adc;
}
return signal_in_adc;
}
Expand Up @@ -116,7 +116,7 @@ class Phase2TrackerDigitizerAlgorithm {
const float ClusterWidth; // Gaussian charge cutoff width in sigma units

//-- make_digis
const bool doDigitalReadout; // Flag to decide analog or digital readout
const int thePhase2ReadoutMode; // Flag to decide readout mode (digital/Analog dual slope etc.)
const float theElectronPerADC; // Gain, number of electrons per adc count.
const int theAdcFullScale; // Saturation count, 255=8bit.
const float theNoiseInElectrons; // Noise (RMS) in units of electrons.
Expand Down Expand Up @@ -209,6 +209,9 @@ class Phase2TrackerDigitizerAlgorithm {
//for engine passed into the constructor from Digitizer
CLHEP::HepRandomEngine* rengine_;

// convert signal in electrons to ADC counts
int convertSignalToAdc(uint32_t detID,float signal_in_elec,float threshold);

double calcQ(float x) const {
auto xx = std::min(0.5f * x * x,12.5f);
return 0.5 * (1.0-std::copysign(std::sqrt(1.f- unsafe_expf<4>(-xx * (1.f + 0.2733f/(1.f + 0.147f * xx)))), x));
Expand Down
Expand Up @@ -9,9 +9,9 @@
OuterTrackerDigiSource = cms.InputTag("mix", "Tracker"),
GeometryType = cms.string('idealForDigi'),
NumberOfDigisPerDetH = cms.PSet(
Nbins = cms.int32(200),
Nbins = cms.int32(100),
xmin = cms.double(-0.5),
xmax = cms.double(200.5)
xmax = cms.double(99.5)
),
DigiOccupancySH = cms.PSet(
Nbins = cms.int32(51),
Expand Down Expand Up @@ -52,9 +52,9 @@
xmax = cms.double(2000.5)
),
NumberOfClustersPerDetH = cms.PSet(
Nbins = cms.int32(200),
Nbins = cms.int32(100),
xmin = cms.double(-0.5),
xmax = cms.double(200.5)
xmax = cms.double(99.5)
),
ClusterWidthH = cms.PSet(
Nbins = cms.int32(16),
Expand Down
Expand Up @@ -30,8 +30,8 @@
HIPThresholdInElectrons_Barrel = cms.double(1.0e10), # very high value to avoid Over threshold bit
HIPThresholdInElectrons_Endcap = cms.double(1.0e10), # very high value to avoid Over threshold bit
NoiseInElectrons = cms.double(0.0),
DigitalReadout = cms.bool(False), # Flag to decide analog or digital readout
AdcFullScale = cms.int32(16),
Phase2ReadoutMode = cms.int32(-1), # Flag to decide Readout Mode :Digital(0) or Analog (linear TDR (-1), dual slope with slope parameters (+1,+2,+3,+4) with threshold subtraction
AdcFullScale = cms.int32(15),
TofUpperCut = cms.double(12.5),
TofLowerCut = cms.double(-12.5),
AddNoisyPixels = cms.bool(False),
Expand Down Expand Up @@ -67,7 +67,7 @@
HIPThresholdInElectrons_Barrel = cms.double(1.0e10), # very high value to avoid Over threshold bit
HIPThresholdInElectrons_Endcap = cms.double(1.0e10), # very high value to avoid Over threshold bit
NoiseInElectrons = cms.double(200), # 30% of the readout noise (should be changed in future)
DigitalReadout = cms.bool(True), # Flag to decide analog or digital readout
Phase2ReadoutMode = cms.int32(0), # Flag to decide Readout Mode :Digital(0) or Analog (linear TDR (-1), dual slope with slope parameters (+1,+2,+3,+4) with threshold subtraction
AdcFullScale = cms.int32(255),
TofUpperCut = cms.double(12.5),
TofLowerCut = cms.double(-12.5),
Expand Down Expand Up @@ -105,7 +105,7 @@
HIPThresholdInElectrons_Barrel = cms.double(21000.), # 1.4 MIP considered as HIP
HIPThresholdInElectrons_Endcap = cms.double(21000.), # 1.4 MIP considered as HIP
NoiseInElectrons = cms.double(700), # 30% of the readout noise (should be changed in future)
DigitalReadout = cms.bool(True), # Flag to decide analog or digital readout
Phase2ReadoutMode = cms.int32(0), # Flag to decide Readout Mode :Digital(0) or Analog (linear TDR (-1), dual slope with slope parameters (+1,+2,+3,+4) with threshold subtraction
AdcFullScale = cms.int32(255),
TofUpperCut = cms.double(12.5),
TofLowerCut = cms.double(-12.5),
Expand Down Expand Up @@ -143,7 +143,7 @@
HIPThresholdInElectrons_Barrel = cms.double(1.0e10), # very high value to avoid Over threshold bit
HIPThresholdInElectrons_Endcap = cms.double(1.0e10), # very high value to avoid Over threshold bit
NoiseInElectrons = cms.double(1000), # 30% of the readout noise (should be changed in future)
DigitalReadout = cms.bool(True), # Flag to decide analog or digital readout
Phase2ReadoutMode = cms.int32(0), # Flag to decide Readout Mode :Digital(0) or Analog (linear TDR (-1), dual slope with slope parameters (+1,+2,+3,+4) with threshold subtraction
AdcFullScale = cms.int32(255),
TofUpperCut = cms.double(12.5),
TofLowerCut = cms.double(-12.5),
Expand Down
26 changes: 22 additions & 4 deletions SimTracker/SiPhase2Digitizer/test/Phase2TrackerMonitorDigi.cc
Expand Up @@ -131,6 +131,7 @@ void Phase2TrackerMonitorDigi::fillITPixelDigiHistos(const edm::Handle<edm::DetS
int nclus = 0;
int width = 1;
int position = 0;
std::vector<int> charges;
for (typename edm::DetSet< PixelDigi >::const_iterator di = DSViter->begin(); di != DSViter->end(); di++) {
int col = di->column(); // column
int row = di->row(); // row
Expand All @@ -149,14 +150,19 @@ void Phase2TrackerMonitorDigi::fillITPixelDigiHistos(const edm::Handle<edm::DetS
if (row_last == -1 ) {
position = row+1;
nclus++;
charges.push_back(adc);
} else {
if (abs(row - row_last) == 1 && col == col_last) {
position += row+1;
width++;
charges.push_back(adc);
} else {
position /= width;
local_mes.ClusterWidth->Fill(width);
local_mes.ClusterPosition->Fill(position);
for (auto v: charges) local_mes.ChargeOfDigisVsWidth->Fill(v, width);
charges.clear();
charges.push_back(adc);
width = 1;
position = row+1;
nclus++;
Expand Down Expand Up @@ -407,7 +413,6 @@ void Phase2TrackerMonitorDigi::bookLayerHistos(DQMStore::IBooker & ibooker, unsi
Parameters.getParameter<int32_t>("Nbins"),
Parameters.getParameter<double>("xmin"),
Parameters.getParameter<double>("xmax"));

Parameters = config_.getParameter<edm::ParameterSet>("DigiOccupancyPH");
HistoName.str("");
HistoName << "DigiOccupancyP_" << fname2.str();
Expand Down Expand Up @@ -440,7 +445,6 @@ void Phase2TrackerMonitorDigi::bookLayerHistos(DQMStore::IBooker & ibooker, unsi
Parameters.getParameter<double>("xmin"),
Parameters.getParameter<double>("xmax"));


Parameters = config_.getParameter<edm::ParameterSet>("NumberOfHitDetsPerLayerH");
HistoName.str("");
HistoName << "NumberOfHitDetectorsPerLayer_" << fname2.str();
Expand All @@ -456,7 +460,6 @@ void Phase2TrackerMonitorDigi::bookLayerHistos(DQMStore::IBooker & ibooker, unsi
Parameters.getParameter<int32_t>("Nbins"),
Parameters.getParameter<double>("xmin"),
Parameters.getParameter<double>("xmax"));

Parameters = config_.getParameter<edm::ParameterSet>("ClusterWidthH");
HistoName.str("");
HistoName << "ClusterWidth_" << fname2.str();
Expand Down Expand Up @@ -497,7 +500,7 @@ void Phase2TrackerMonitorDigi::bookLayerHistos(DQMStore::IBooker & ibooker, unsi
HistoName << "FractionOfOverThresholdDigisVaEta_" << fname2.str();
local_mes.FractionOfOvTBitsVsEta= ibooker.bookProfile(HistoName.str(), HistoName.str(),
EtaParameters.getParameter<int32_t>("Nbins"),EtaParameters.getParameter<double>("xmin"),EtaParameters.getParameter<double>("xmax"),
Parameters.getParameter<double>("xmin"),Parameters.getParameter<double>("xmax"),"");
Parameters.getParameter<double>("xmin"),Parameters.getParameter<double>("xmax"),"");
} else {

Parameters = config_.getParameter<edm::ParameterSet>("DigiChargeH");
Expand All @@ -507,10 +510,25 @@ void Phase2TrackerMonitorDigi::bookLayerHistos(DQMStore::IBooker & ibooker, unsi
Parameters.getParameter<int32_t>("Nbins"),
Parameters.getParameter<double>("xmin"),
Parameters.getParameter<double>("xmax"));

edm::ParameterSet WidthParameters = config_.getParameter<edm::ParameterSet>("ClusterWidthH");
HistoName.str("");
HistoName << "ChargeOfDigisVsWidth_" << fname2.str();
local_mes.ChargeOfDigisVsWidth = ibooker.book2D(HistoName.str(),HistoName.str(),
Parameters.getParameter<int32_t>("Nbins"),
Parameters.getParameter<double>("xmin"),
Parameters.getParameter<double>("xmax"),
WidthParameters.getParameter<int32_t>("Nbins"),
WidthParameters.getParameter<double>("xmin"),
WidthParameters.getParameter<double>("xmax"));


}

layerMEs.insert(std::make_pair(layer, local_mes));
}
}
void Phase2TrackerMonitorDigi::endLuminosityBlock(edm::LuminosityBlock const& lumiBlock, edm::EventSetup const& eSetup){
}
//define this as a plug-in
DEFINE_FWK_MODULE(Phase2TrackerMonitorDigi);
2 changes: 2 additions & 0 deletions SimTracker/SiPhase2Digitizer/test/Phase2TrackerMonitorDigi.h
Expand Up @@ -26,6 +26,7 @@ class Phase2TrackerMonitorDigi : public DQMEDAnalyzer{
edm::EventSetup const & iSetup ) override;
void dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) override;
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
void endLuminosityBlock(edm::LuminosityBlock const& lumiBlock, edm::EventSetup const& iSetup) override;


struct DigiMEs{
Expand All @@ -34,6 +35,7 @@ class Phase2TrackerMonitorDigi : public DQMEDAnalyzer{
MonitorElement* DigiOccupancyS;
MonitorElement* PositionOfDigis;
MonitorElement* ChargeOfDigis;
MonitorElement* ChargeOfDigisVsWidth;
MonitorElement* TotalNumberOfDigisPerLayer;
MonitorElement* NumberOfHitDetectorsPerLayer;
MonitorElement* NumberOfClustersPerDet;
Expand Down