CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
DeDxHitCalibrator Class Reference
Inheritance diagram for DeDxHitCalibrator:
edm::stream::EDProducer<>

Public Types

typedef std::pair< uint32_t, unsigned char > ChipId
 
- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

 DeDxHitCalibrator (const edm::ParameterSet &)
 
 ~DeDxHitCalibrator () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &)
 

Static Public Attributes

static constexpr int kIsBelow = 1
 
static constexpr int kIsNormal = 0
 
static constexpr int kIsOver = 2
 
static constexpr int PXB = 0
 
static constexpr int PXF = 1
 
static constexpr int TECThick = 6
 
static constexpr int TECThin = 5
 
static constexpr int TIB = 2
 
static constexpr int TID = 3
 
static constexpr int TOB = 4
 

Private Member Functions

void beginRun (edm::Run const &, const edm::EventSetup &) override
 
float correctEnergy (const float &, const ChipId &)
 
std::pair< double, double > fitStripCluster (const std::vector< std::pair< double, int > > &, const double &, const double &)
 
void getAlphaBeta (const std::vector< double > &, const std::vector< std::pair< double, int > > &, CLHEP::HepMatrix &, CLHEP::HepVector &, const std::vector< bool > &, const double &, const double &)
 
double getChi2 (const std::vector< double > &, const std::vector< std::pair< double, int > > &, const double &, const double &)
 
int getDetId (const DetId &, const float &)
 
void processHitInfo (const reco::DeDxHitInfo &, const float &trackMomentum, reco::DeDxHitCollection &, reco::DeDxHitCollection &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const bool applyGain_
 
edm::ESHandle< DeDxCalibrationdedxCalib_
 
const edm::ESGetToken< DeDxCalibration, DeDxCalibrationRcddedxCalibToken_
 
const edm::EDGetTokenT< reco::DeDxHitInfoAssdedxHitInfoToken_
 
const double MeVPerElectron_
 
SiPixelGainCalibrationOfflineService pixelCalib_
 
const int pixelSaturationThr_
 
edm::ESHandle< TrackerGeometrytkGeom_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 
edm::ESHandle< TrackerTopologytkTopo_
 
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtkTopoToken_
 
const edm::EDGetTokenT< reco::TrackCollectiontracksToken_
 
const int VCaltoElectronGain_
 
const int VCaltoElectronGain_L1_
 
const int VCaltoElectronOffset_
 
const int VCaltoElectronOffset_L1_
 

Detailed Description

Definition at line 27 of file DeDxHitCalibrator.cc.

Member Typedef Documentation

◆ ChipId

typedef std::pair<uint32_t, unsigned char> DeDxHitCalibrator::ChipId

Definition at line 31 of file DeDxHitCalibrator.cc.

Constructor & Destructor Documentation

◆ DeDxHitCalibrator()

DeDxHitCalibrator::DeDxHitCalibrator ( const edm::ParameterSet iConfig)
explicit

Definition at line 74 of file DeDxHitCalibrator.cc.

75  : applyGain_(iConfig.getParameter<bool>("applyGain")),
76  MeVPerElectron_(iConfig.getParameter<double>("MeVPerElectron")),
77  VCaltoElectronGain_(iConfig.getParameter<int>("VCaltoElectronGain")),
78  VCaltoElectronGain_L1_(iConfig.getParameter<int>("VCaltoElectronGain_L1")),
79  VCaltoElectronOffset_(iConfig.getParameter<int>("VCaltoElectronOffset")),
80  VCaltoElectronOffset_L1_(iConfig.getParameter<int>("VCaltoElectronOffset_L1")),
81  pixelSaturationThr_(iConfig.getParameter<int>("pixelSaturationThr")),
82  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("trackProducer"))),
83  dedxHitInfoToken_(consumes<reco::DeDxHitInfoAss>(iConfig.getParameter<edm::InputTag>("dedxHitInfo"))),
84  dedxCalibToken_(esConsumes<DeDxCalibration, DeDxCalibrationRcd, edm::Transition::BeginRun>()),
85  tkGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
86  tkTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
87  pixelCalib_(iConfig, consumesCollector()) {
88  produces<reco::TrackDeDxHitsCollection>("PixelHits");
89  produces<reco::TrackDeDxHitsCollection>("StripHits");
90 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const edm::EDGetTokenT< reco::DeDxHitInfoAss > dedxHitInfoToken_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_
const edm::ESGetToken< DeDxCalibration, DeDxCalibrationRcd > dedxCalibToken_
const int VCaltoElectronGain_
const int VCaltoElectronOffset_
const int VCaltoElectronOffset_L1_
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
SiPixelGainCalibrationOfflineService pixelCalib_
const int VCaltoElectronGain_L1_
const double MeVPerElectron_
const int pixelSaturationThr_

◆ ~DeDxHitCalibrator()

DeDxHitCalibrator::~DeDxHitCalibrator ( )
overridedefault

Member Function Documentation

◆ beginRun()

void DeDxHitCalibrator::beginRun ( edm::Run const &  ,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 106 of file DeDxHitCalibrator.cc.

References dedxCalib_, dedxCalibToken_, edm::EventSetup::getHandle(), tkGeom_, tkGeomToken_, tkTopo_, and tkTopoToken_.

106  {
108  tkGeom_ = iSetup.getHandle(tkGeomToken_);
109  tkTopo_ = iSetup.getHandle(tkTopoToken_);
110 }
edm::ESHandle< TrackerGeometry > tkGeom_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
edm::ESHandle< TrackerTopology > tkTopo_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_
edm::ESHandle< DeDxCalibration > dedxCalib_
const edm::ESGetToken< DeDxCalibration, DeDxCalibrationRcd > dedxCalibToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130

◆ correctEnergy()

float DeDxHitCalibrator::correctEnergy ( const float &  energy,
const ChipId detId 
)
private

Definition at line 242 of file DeDxHitCalibrator.cc.

References applyGain_, dedxCalib_, hcalRecHitTable_cff::detId, HBHEDarkening_cff::energy, g, and DeDxCalibration::gain().

Referenced by processHitInfo().

242  {
243  const auto& g = dedxCalib_->gain().find(detId);
244  if (applyGain_ && g != dedxCalib_->gain().end())
245  return energy * g->second;
246  return energy;
247 }
const std::map< ChipId, float > & gain() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
edm::ESHandle< DeDxCalibration > dedxCalib_

◆ fillDescriptions()

void DeDxHitCalibrator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 92 of file DeDxHitCalibrator.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

92  {
94  desc.add<bool>("applyGain", true);
95  desc.add<double>("MeVPerElectron", 3.61e-06);
96  desc.add<int>("VCaltoElectronGain", 65);
97  desc.add<int>("VCaltoElectronGain_L1", 65);
98  desc.add<int>("VCaltoElectronOffset", -414);
99  desc.add<int>("VCaltoElectronOffset_L1", -414);
100  desc.add<int>("pixelSaturationThr", 254);
101  desc.add<edm::InputTag>("trackProducer", edm::InputTag("generalTracks"));
102  desc.add<edm::InputTag>("dedxHitInfo", edm::InputTag("dedxHitInfo"));
103  descriptions.add("dedxHitCalibrator", desc);
104 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ fitStripCluster()

std::pair< double, double > DeDxHitCalibrator::fitStripCluster ( const std::vector< std::pair< double, int > > &  b,
const double &  coupling,
const double &  iSigma 
)
private

Definition at line 365 of file DeDxHitCalibrator.cc.

References funct::abs(), isotrackApplyRegressor::alpha, b, HLT_2024v14_cff::beta, dumpMFGeometry_cfg::delta, change_name::diff, MillePedeFileConverter_cfg::e, RemoveAddSevLevel::flag, getAlphaBeta(), getChi2(), mps_fire::i, cuy::ib, dqmiolumiharvest::j, GetRecoTauVFromDQM_MC_cff::next, convertSQLiteXML::ok, mathSSE::sqrt(), and x.

Referenced by processHitInfo().

367  {
368  const auto& npar = b.size();
369  std::vector<bool> isFix(npar, false);
370  std::vector<double> x;
371  x.reserve(b.size());
372  for (const auto& ib : b)
373  x.emplace_back(ib.first);
374 
375  bool ok = false;
376  CLHEP::HepMatrix hessian(npar, npar);
377  int iter = 0;
378  do {
379  double diff(0);
380  auto old = getChi2(x, b, coupling, iSigma);
381  do {
382  CLHEP::HepMatrix alpha(npar, npar, 0.);
383  CLHEP::HepVector beta(npar, 0.);
384  getAlphaBeta(x, b, alpha, beta, isFix, coupling, iSigma);
385  const auto& delta = CLHEP::solve(alpha, -beta);
386  double lambda(1.);
387  std::vector<double> xn(npar);
388  for (size_t i = 0; i < npar; i++)
389  xn[i] = x[i] + lambda * delta[i];
390  const auto& next = getChi2(xn, b, coupling, iSigma);
391  diff = old - next;
392  if (diff > 0)
393  x = xn;
394  old = next;
395  hessian = alpha;
396  iter++;
397  } while (diff > 1e-6 && iter < 100);
398  // Check if we have negatives
399  const auto& mi = std::min_element(x.begin(), x.end()) - x.begin();
400  if (x[mi] < 0) {
401  x[mi] = 0.;
402  isFix[mi] = true;
403  } else
404  ok = true;
405  } while (!ok && iter < 100);
406  hessian /= 2.;
407  int flag;
408  hessian.invert(flag);
409 
410  double var2(0.);
411  for (size_t i = 0; i < npar; i++)
412  for (size_t j = 0; j < npar; j++)
413  if (!isFix[i] && !isFix[j])
414  var2 += hessian[i][j];
415  var2 = std::abs(var2);
416 
417  double sum(0.);
418  for (size_t i = 0; i < npar; i++)
419  if (!isFix[i])
420  sum += x[i];
421 
422  return {sum, std::sqrt(var2)};
423 }
T sqrt(T t)
Definition: SSEVec.h:23
double getChi2(const std::vector< double > &, const std::vector< std::pair< double, int > > &, const double &, const double &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void getAlphaBeta(const std::vector< double > &, const std::vector< std::pair< double, int > > &, CLHEP::HepMatrix &, CLHEP::HepVector &, const std::vector< bool > &, const double &, const double &)
double b
Definition: hdecay.h:120
ib
Definition: cuy.py:661

◆ getAlphaBeta()

void DeDxHitCalibrator::getAlphaBeta ( const std::vector< double > &  x,
const std::vector< std::pair< double, int > > &  b,
CLHEP::HepMatrix &  alpha,
CLHEP::HepVector &  beta,
const std::vector< bool > &  isFix,
const double &  coupling,
const double &  iSigma 
)
private

Definition at line 288 of file DeDxHitCalibrator.cc.

References a0, isotrackTrainRegressor::a1, isotrackApplyRegressor::alpha, b, HLT_2024v14_cff::beta, PVValHelper::dx, f, mps_fire::i, kIsBelow, kIsNormal, kIsOver, submitPVResolutionJobs::q, edm::second(), x, and y.

Referenced by fitStripCluster().

294  {
295  const auto& npar = b.size();
296  const auto a0 = coupling * iSigma;
297  const auto a1 = (1 - 2 * coupling) * iSigma;
298  const auto a00 = a0 * a0;
299  const auto a11 = a1 * a1;
300  const auto a01 = a0 * a1;
301  std::vector<double> y(npar, 0.);
302  for (size_t i = 0; i < npar; i++)
303  if (!isFix[i]) {
304  const auto dx = coupling * x[i];
305  if (i >= 1)
306  y[i - 1] += dx;
307  y[i] += x[i] - 2 * dx;
308  if (i + 1 < npar)
309  y[i + 1] += dx;
310  }
311  for (size_t i = 0; i < npar; i++) {
312  auto q = (y[i] - b[i].first) * iSigma;
313  int f(0);
314  if (b[i].second == kIsNormal)
315  f = 2;
316  else if (b[i].second == kIsBelow) {
317  q += 2;
318  if (q > 0)
319  f = 1;
320  } else if (b[i].second == kIsOver) {
321  q -= 2;
322  if (q < 0)
323  f = 1;
324  }
325  if (f > 0) {
326  if (i >= 1)
327  if (!isFix[i - 1]) {
328  alpha[i - 1][i - 1] += f * a00;
329  if (!isFix[i]) {
330  alpha[i - 1][i] += f * a01;
331  alpha[i][i - 1] += f * a01;
332  }
333  beta[i - 1] += f * q * a0;
334  }
335  if (!isFix[i]) {
336  alpha[i][i] += f * a11;
337  beta[i] += f * q * a1;
338  }
339  if (i + 1 < npar)
340  if (!isFix[i + 1]) {
341  alpha[i + 1][i + 1] += f * a00;
342  if (!isFix[i]) {
343  alpha[i + 1][i] += f * a01;
344  alpha[i][i + 1] += f * a01;
345  }
346  beta[i + 1] += f * q * a0;
347  }
348  }
349  // Penalty for negatives
350  if (!isFix[i])
351  if (x[i] < 0) {
352  alpha[i][i] += 2 * iSigma * iSigma;
353  q = x[i] * iSigma;
354  beta[i] += 2 * q * iSigma;
355  }
356  }
357  for (size_t i = 0; i < npar; i++)
358  if (isFix[i]) {
359  alpha[i][i] = 1.;
360  beta[i] = 0.;
361  }
362 }
U second(std::pair< T, U > const &p)
double f[11][100]
static constexpr float a0
double b
Definition: hdecay.h:120
static constexpr int kIsOver
static constexpr int kIsNormal
static constexpr int kIsBelow

◆ getChi2()

double DeDxHitCalibrator::getChi2 ( const std::vector< double > &  x,
const std::vector< std::pair< double, int > > &  b,
const double &  coupling,
const double &  iSigma 
)
private

Definition at line 250 of file DeDxHitCalibrator.cc.

References b, isoTrack_cff::chi2, PVValHelper::dx, mps_fire::i, kIsBelow, kIsNormal, kIsOver, submitPVResolutionJobs::q, edm::second(), x, and y.

Referenced by fitStripCluster().

253  {
254  const auto& npar = b.size();
255  std::vector<double> y(npar, 0.);
256  for (size_t i = 0; i < npar; i++) {
257  const auto dx = coupling * x[i];
258  if (i >= 1)
259  y[i - 1] += dx;
260  y[i] += x[i] - 2 * dx;
261  if (i + 1 < npar)
262  y[i + 1] += dx;
263  }
264  double chi2(0.);
265  for (size_t i = 0; i < npar; i++) {
266  auto q = (b[i].first - y[i]) * iSigma;
267  if (b[i].second == kIsNormal)
268  chi2 += q * q;
269  else if (b[i].second == kIsBelow) {
270  q -= 2;
271  if (q < 0)
272  chi2 += 0.5 * q * q;
273  } else if (b[i].second == kIsOver) {
274  q += 2;
275  if (q > 0)
276  chi2 += 0.5 * q * q;
277  }
278  // Penalty for negatives
279  if (x[i] < 0) {
280  q = x[i] * iSigma;
281  chi2 += q * q;
282  }
283  }
284  return chi2;
285 }
U second(std::pair< T, U > const &p)
double b
Definition: hdecay.h:120
static constexpr int kIsOver
static constexpr int kIsNormal
static constexpr int kIsBelow

◆ getDetId()

int DeDxHitCalibrator::getDetId ( const DetId id,
const float &  thickness 
)
private

Definition at line 234 of file DeDxHitCalibrator.cc.

References MillePedeFileConverter_cfg::e, TECThick, TECThin, and ppsPixelTopologyESSourceRun2_cfi::thickness.

Referenced by processHitInfo().

234  {
235  const auto& subdet = id.subdetId() - 1;
236  if (subdet == TECThin && thickness > 400e-4)
237  return TECThick;
238  return subdet;
239 }
static constexpr int TECThick
static constexpr int TECThin

◆ processHitInfo()

void DeDxHitCalibrator::processHitInfo ( const reco::DeDxHitInfo info,
const float &  trackMomentum,
reco::DeDxHitCollection pixelHits,
reco::DeDxHitCollection stripHits 
)
private

Definition at line 137 of file DeDxHitCalibrator.cc.

References a, gpuClustering::adc, DeDxCalibration::alpha(), b, ALCARECOTkAlJpsiMuMu_cff::charge, cuy::col, reco::DeDxHitInfo::Compatible, reco::DeDxHitInfo::Complete, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), correctEnergy(), dedxCalib_, dumpMFGeometry_cfg::delta, hcalRecHitTable_cff::detId, MillePedeFileConverter_cfg::e, HBHEDarkening_cff::energy, fitStripCluster(), getDetId(), SiPixelGainCalibrationOfflineService::getGain(), SiPixelGainCalibrationOfflineService::getPedestal(), mps_fire::i, TrackerGeometry::idToDet(), info(), createfilelist::int, CastorSimpleRecAlgoImpl::isSaturated(), dqmiolumiharvest::j, kIsBelow, kIsNormal, kIsOver, CrabHelper::log, WZElectronSkims53X_cff::max, sistrip::MeVperADCStrip, MeVPerElectron_, ecaldqm::nChannels, pixelCalib_, TrackingMonitor_cfi::pixelCluster, InitialStepPreSplitting_cff::pixelHits, pixelSaturationThr_, TrackerTopology::pxbLayer(), PXF, mps_fire::result, rpixValues::ROCSizeInX, rpixValues::ROCSizeInY, DeDxCalibration::sigma(), mathSSE::sqrt(), TrackingMonitor_cfi::stripCluster, InitialStep_cff::stripHits, sistrip::STRIPS_PER_APV, ppsPixelTopologyESSourceRun2_cfi::thickness, DeDxCalibration::thr(), tkGeom_, tkTopo_, reco::btau::trackMomentum, VCaltoElectronGain_, VCaltoElectronGain_L1_, VCaltoElectronOffset_, and VCaltoElectronOffset_L1_.

Referenced by produce().

140  {
141  // Energy loss parametrization from https://doi.org/10.1016/j.nima.2012.06.064
142  static constexpr double a = 0.07; //energy loss factor in silicon
143  static constexpr double l0 = 450e-4; //reference path length in um
144 
145  for (size_t i = 0; i < info.size(); i++) {
146  // Require hits to be complete and compatible
147  const auto& type = info.type(i);
149  continue;
150 
151  // Effective path length
152  const auto& pl = info.pathlength(i);
153  const auto pathLength = pl * (1. + a * std::log(pl / l0));
154  const DetId& detId = info.detId(i);
155 
156  // Strip
157  if (const auto& stripCluster = info.stripCluster(i)) {
158  const auto& thickness =
159  dynamic_cast<const StripGeomDetUnit*>(tkGeom_->idToDet(detId))->surface().bounds().thickness();
160  const auto& det = getDetId(detId, thickness);
161  const auto& chipId = ChipId(detId, stripCluster->barycenter() / sistrip::STRIPS_PER_APV);
162  // Fill measured strip deposits
163  const auto& thr = dedxCalib_->thr()[det];
164  std::vector<std::pair<double, int> > b;
165  b.reserve(stripCluster->amplitudes().size() + 2);
166  b.emplace_back(thr, kIsBelow);
167  for (const auto& adc : stripCluster->amplitudes()) {
168  if (adc > 253)
169  b.emplace_back(254., kIsOver);
170  else
171  b.emplace_back(adc + 0.5, kIsNormal);
172  }
173  b.emplace_back(thr, kIsBelow);
174  // Fit
175  const auto& result = fitStripCluster(b, dedxCalib_->alpha()[det], 1. / dedxCalib_->sigma()[det]);
176  // ADC -> e- -> MeV
177  const auto& energy = correctEnergy(result.first * sistrip::MeVperADCStrip, chipId);
178  const auto& esigma = correctEnergy(result.second * sistrip::MeVperADCStrip, chipId);
179  // Compute charge
180  const auto& charge = pathLength != 0. ? energy / pathLength : 0.;
181  stripHits.emplace_back(charge, trackMomentum, pathLength, 0, esigma);
182  }
183 
184  // Pixel
185  if (const auto& pixelCluster = info.pixelCluster(i)) {
186  const auto& thickness =
187  dynamic_cast<const PixelGeomDetUnit*>(tkGeom_->idToDet(detId))->surface().bounds().thickness();
188  const auto& det = getDetId(detId, thickness);
189  const auto& chipId =
190  ChipId(detId, (int(pixelCluster->x() / ROCSizeInX) << 3) + int(pixelCluster->y() / ROCSizeInY));
191  // Collect adc
192  double delta(0);
193  bool isSaturated(false);
194  for (size_t j = 0; j < pixelCluster->pixelADC().size(); j++) {
195  const auto& elec = pixelCluster->pixelADC()[j];
196  delta += elec * MeVPerElectron_;
197  if (isSaturated)
198  continue;
199  // ( pos , (x , y) )
200  const auto& row = pixelCluster->minPixelRow() + pixelCluster->pixelOffset()[2 * j];
201  const auto& col = pixelCluster->minPixelCol() + pixelCluster->pixelOffset()[2 * j + 1];
202  // Go back to adc
203  const auto& DBgain = pixelCalib_.getGain(detId, col, row);
204  const auto& DBpedestal = pixelCalib_.getPedestal(detId, col, row);
206  isSaturated = true;
207  else if (DBgain > 0.) {
208  double vcal;
209  const auto& theLayer = (detId.subdetId() == 1) ? tkTopo_->pxbLayer(detId) : 0;
210  if (theLayer == 1)
212  else
213  vcal = (elec - VCaltoElectronOffset_) / VCaltoElectronGain_;
214  const auto adc = std::round(DBpedestal + vcal / DBgain);
215  if (adc > pixelSaturationThr_)
216  isSaturated = true;
217  }
218  }
219  // Compute energy
220  const auto& energy = correctEnergy(delta, chipId);
221  if (det != PXF && energy <= 50e-3)
222  continue;
223  // Estimate sigma for cluster
224  const unsigned char& nChannels = pixelCluster->pixelADC().size();
225  const auto esigma = 10e-3 * std::sqrt(nChannels);
226  // Compute charge
227  const auto& charge = pathLength != 0. ? energy / pathLength : 0.;
228  pixelHits.emplace_back(charge, trackMomentum, pathLength, isSaturated, esigma);
229  }
230  }
231 }
std::pair< uint32_t, unsigned char > ChipId
unsigned int pxbLayer(const DetId &id) const
static const TGPicture * info(bool iBackgroundIsBlack)
std::pair< double, double > fitStripCluster(const std::vector< std::pair< double, int > > &, const double &, const double &)
edm::ESHandle< TrackerGeometry > tkGeom_
static constexpr int Compatible
Definition: DeDxHitInfo.h:40
const std::vector< double > & thr() const
const std::vector< double > & alpha() const
static constexpr int nChannels
edm::ESHandle< TrackerTopology > tkTopo_
edm::ESHandle< DeDxCalibration > dedxCalib_
T sqrt(T t)
Definition: SSEVec.h:23
int getDetId(const DetId &, const float &)
const int VCaltoElectronGain_
const int VCaltoElectronOffset_
const TrackerGeomDet * idToDet(DetId) const override
Definition: DetId.h:17
const int VCaltoElectronOffset_L1_
double b
Definition: hdecay.h:120
static constexpr int Complete
Definition: DeDxHitInfo.h:40
bool isSaturated(const Digi &digi, const int &maxADCvalue, int ifirst, int n)
static constexpr int kIsOver
SiPixelGainCalibrationOfflineService pixelCalib_
const int VCaltoElectronGain_L1_
static const uint16_t STRIPS_PER_APV
double a
Definition: hdecay.h:121
static constexpr float MeVperADCStrip
col
Definition: cuy.py:1009
float correctEnergy(const float &, const ChipId &)
const std::vector< double > & sigma() const
float getPedestal(const uint32_t &detID, const int &col, const int &row) override
const double MeVPerElectron_
static constexpr int PXF
static constexpr int kIsNormal
uint16_t *__restrict__ uint16_t const *__restrict__ adc
constexpr int ROCSizeInX
float getGain(const uint32_t &detID, const int &col, const int &row) override
static constexpr int kIsBelow
constexpr int ROCSizeInY
const int pixelSaturationThr_

◆ produce()

void DeDxHitCalibrator::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 112 of file DeDxHitCalibrator.cc.

References dedxEstimators_cff::dedxHitInfo, dedxHitInfoToken_, mps_fire::i, iEvent, eostools::move(), pixelCalib_, InitialStepPreSplitting_cff::pixelHits, processHitInfo(), SiPixelGainCalibrationServicePayloadGetter< thePayloadObject, theDBRecordType >::setESObjects(), InitialStep_cff::stripHits, HLT_2024v14_cff::track, DiMuonV_cfg::tracks, and tracksToken_.

112  {
113  const auto& tracks = iEvent.getHandle(tracksToken_);
114  const auto& dedxHitInfo = iEvent.get(dedxHitInfoToken_);
115  pixelCalib_.setESObjects(iSetup);
116 
117  // creates the output collection
118  auto pixelTrackDeDxHitAss = std::make_unique<reco::TrackDeDxHitsCollection>(reco::TrackRefProd(tracks));
119  auto stripTrackDeDxHitAss = std::make_unique<reco::TrackDeDxHitsCollection>(reco::TrackRefProd(tracks));
120 
121  for (size_t i = 0; i < tracks->size(); i++) {
122  const auto& track = reco::TrackRef(tracks, i);
123  const auto& dedxHits = dedxHitInfo[track];
124  if (dedxHits.isNull())
125  continue;
127  processHitInfo(*dedxHits, track->p(), pixelHits, stripHits);
128  pixelTrackDeDxHitAss->setValue(i, pixelHits);
129  stripTrackDeDxHitAss->setValue(i, stripHits);
130  }
131 
132  iEvent.put(std::move(pixelTrackDeDxHitAss), "PixelHits");
133  iEvent.put(std::move(stripTrackDeDxHitAss), "StripHits");
134 }
void processHitInfo(const reco::DeDxHitInfo &, const float &trackMomentum, reco::DeDxHitCollection &, reco::DeDxHitCollection &)
std::vector< DeDxHit > DeDxHitCollection
Definition: DeDxHit.h:45
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< reco::DeDxHitInfoAss > dedxHitInfoToken_
void setESObjects(const edm::EventSetup &es) override
edm::RefProd< TrackCollection > TrackRefProd
persistent reference to a Track collection
Definition: TrackFwd.h:26
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
SiPixelGainCalibrationOfflineService pixelCalib_
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ applyGain_

const bool DeDxHitCalibrator::applyGain_
private

Definition at line 58 of file DeDxHitCalibrator.cc.

Referenced by correctEnergy().

◆ dedxCalib_

edm::ESHandle<DeDxCalibration> DeDxHitCalibrator::dedxCalib_
private

Definition at line 69 of file DeDxHitCalibrator.cc.

Referenced by beginRun(), correctEnergy(), and processHitInfo().

◆ dedxCalibToken_

const edm::ESGetToken<DeDxCalibration, DeDxCalibrationRcd> DeDxHitCalibrator::dedxCalibToken_
private

Definition at line 64 of file DeDxHitCalibrator.cc.

Referenced by beginRun().

◆ dedxHitInfoToken_

const edm::EDGetTokenT<reco::DeDxHitInfoAss> DeDxHitCalibrator::dedxHitInfoToken_
private

Definition at line 63 of file DeDxHitCalibrator.cc.

Referenced by produce().

◆ kIsBelow

constexpr int DeDxHitCalibrator::kIsBelow = 1
static

Definition at line 29 of file DeDxHitCalibrator.cc.

Referenced by getAlphaBeta(), getChi2(), and processHitInfo().

◆ kIsNormal

constexpr int DeDxHitCalibrator::kIsNormal = 0
static

Definition at line 29 of file DeDxHitCalibrator.cc.

Referenced by getAlphaBeta(), getChi2(), and processHitInfo().

◆ kIsOver

constexpr int DeDxHitCalibrator::kIsOver = 2
static

Definition at line 29 of file DeDxHitCalibrator.cc.

Referenced by getAlphaBeta(), getChi2(), and processHitInfo().

◆ MeVPerElectron_

const double DeDxHitCalibrator::MeVPerElectron_
private

Definition at line 59 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ pixelCalib_

SiPixelGainCalibrationOfflineService DeDxHitCalibrator::pixelCalib_
private

Definition at line 68 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo(), and produce().

◆ pixelSaturationThr_

const int DeDxHitCalibrator::pixelSaturationThr_
private

Definition at line 61 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ PXB

constexpr int DeDxHitCalibrator::PXB = 0
static

Definition at line 30 of file DeDxHitCalibrator.cc.

◆ PXF

constexpr int DeDxHitCalibrator::PXF = 1
static

Definition at line 30 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ TECThick

constexpr int DeDxHitCalibrator::TECThick = 6
static

Definition at line 30 of file DeDxHitCalibrator.cc.

Referenced by getDetId().

◆ TECThin

constexpr int DeDxHitCalibrator::TECThin = 5
static

Definition at line 30 of file DeDxHitCalibrator.cc.

Referenced by getDetId().

◆ TIB

constexpr int DeDxHitCalibrator::TIB = 2
static

Definition at line 30 of file DeDxHitCalibrator.cc.

◆ TID

constexpr int DeDxHitCalibrator::TID = 3
static

Definition at line 30 of file DeDxHitCalibrator.cc.

◆ tkGeom_

edm::ESHandle<TrackerGeometry> DeDxHitCalibrator::tkGeom_
private

Definition at line 70 of file DeDxHitCalibrator.cc.

Referenced by beginRun(), and processHitInfo().

◆ tkGeomToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> DeDxHitCalibrator::tkGeomToken_
private

Definition at line 65 of file DeDxHitCalibrator.cc.

Referenced by beginRun().

◆ tkTopo_

edm::ESHandle<TrackerTopology> DeDxHitCalibrator::tkTopo_
private

Definition at line 71 of file DeDxHitCalibrator.cc.

Referenced by beginRun(), and processHitInfo().

◆ tkTopoToken_

const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> DeDxHitCalibrator::tkTopoToken_
private

Definition at line 66 of file DeDxHitCalibrator.cc.

Referenced by beginRun().

◆ TOB

constexpr int DeDxHitCalibrator::TOB = 4
static

Definition at line 30 of file DeDxHitCalibrator.cc.

◆ tracksToken_

const edm::EDGetTokenT<reco::TrackCollection> DeDxHitCalibrator::tracksToken_
private

Definition at line 62 of file DeDxHitCalibrator.cc.

Referenced by produce().

◆ VCaltoElectronGain_

const int DeDxHitCalibrator::VCaltoElectronGain_
private

Definition at line 60 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ VCaltoElectronGain_L1_

const int DeDxHitCalibrator::VCaltoElectronGain_L1_
private

Definition at line 60 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ VCaltoElectronOffset_

const int DeDxHitCalibrator::VCaltoElectronOffset_
private

Definition at line 60 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().

◆ VCaltoElectronOffset_L1_

const int DeDxHitCalibrator::VCaltoElectronOffset_L1_
private

Definition at line 60 of file DeDxHitCalibrator.cc.

Referenced by processHitInfo().