CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelDigitizerAlgorithm.cc
Go to the documentation of this file.
1 //class SiPixelDigitizerAlgorithm SimTracker/SiPixelDigitizer/src/SiPixelDigitizerAlgoithm.cc
2 
3 // Original Author Danek Kotlinski
4 // Ported in CMSSW by Michele Pioppi-INFN perugia
5 // Added DB capabilities by F.Blekman, Cornell University
6 // Created: Mon Sep 26 11:08:32 CEST 2005
7 // Add tof, change AddNoise to tracked. 4/06
8 // Change drift direction. 6/06 d.k.
9 // Add the statuis (non-rate dependent) inefficiency.
10 // -1 - no ineffciency
11 // 0 - static inefficency only
12 // 1,2 - low-lumi rate dependent inefficency added
13 // 10 - high-lumi inefficiency added
14 // Adopt the correct drift sign convetion from Morris Swartz. d.k. 8/06
15 // Add more complex misscalinbration, change kev/e to 3.61, diff=3.7,d.k.9/06
16 // Add the readout channel electronic noise. d.k. 3/07
17 // Lower the pixel noise from 500 to 175elec.
18 // Change the input threshold from noise units to electrons.
19 // Lower the amount of static dead pixels from 0.01 to 0.001.
20 // Modify to the new random number services. d.k. 5/07
21 // Protect against sigma=0 (delta tracks on the surface). d.k.5/07
22 // Change the TOF cut to lower and upper limit. d.k. 7/07
23 //
24 // July 2008: Split Lorentz Angle configuration in BPix/FPix (V. Cuplov)
25 // tanLorentzAngleperTesla_FPix=0.0912 and tanLorentzAngleperTesla_BPix=0.106
26 // Sept. 2008: Disable Pixel modules which are declared dead in the configuration python file. (V. Cuplov)
27 // Oct. 2008: Accessing/Reading the Lorentz angle from the DataBase instead of the cfg file. (V. Cuplov)
28 // Accessing dead modules from the DB. Implementation done and tested on a test.db
29 // Do not use this option for now. The PixelQuality Objects are not in the official DB yet.
30 // Feb. 2009: Split Fpix and Bpix threshold and use official numbers (V. Cuplov)
31 // ThresholdInElectrons_FPix = 2870 and ThresholdInElectrons_BPix = 3700
32 // update the electron to VCAL conversion using: VCAL_electrons = VCAL * 65.5 - 414
33 // Feb. 2009: Threshold gaussian smearing (V. Cuplov)
34 // March, 2009: changed DB access to *SimRcd objects (to de-couple the DB objects from reco chain) (F. Blekman)
35 // May, 2009: Pixel charge VCAL smearing. (V. Cuplov)
36 // November, 2009: new parameterization of the pixel response. (V. Cuplov)
37 // December, 2009: Fix issue with different compilers.
38 // October, 2010: Improvement: Removing single dead ROC (V. Cuplov)
39 // November, 2010: Bug fix in removing TBMB/A half-modules (V. Cuplov)
40 // February, 2011: Time improvement in DriftDirection() (J. Bashir Butt)
41 // June, 2011: Bug Fix for pixels on ROC edges in module_killing_DB() (J. Bashir Butt)
42 // February, 2018: Implement cluster charge reweighting (P. Schuetze, with code from A. Hazi)
43 #include <iostream>
44 #include <iomanip>
45 
47 
54 
55 #include <gsl/gsl_sf_erf.h>
57 #include "CLHEP/Random/RandGaussQ.h"
58 #include "CLHEP/Random/RandFlat.h"
59 #include "CLHEP/Random/RandGeneral.h"
60 
61 //#include "PixelIndices.h"
65 
71 
72 // Accessing dead pixel modules from the DB:
74 
76 
90 
97 
98 // Geometry
101 
103 
104 using namespace edm;
105 using namespace sipixelobjects;
106 
108  if (use_ineff_from_db_) { // load gain calibration service fromdb...
109  theSiPixelGainCalibrationService_->setESObjects(es);
110  }
111  if (use_deadmodule_DB_) {
112  SiPixelBadModule_ = &es.getData(SiPixelBadModuleToken_);
113  }
114  if (use_LorentzAngle_DB_) {
115  // Get Lorentz angle from DB record
116  SiPixelLorentzAngle_ = &es.getData(SiPixelLorentzAngleToken_);
117  }
118  //gets the map and geometry from the DB (to kill ROCs)
119  map_ = &es.getData(mapToken_);
120  geom_ = &es.getData(geomToken_);
121 
122  if (KillBadFEDChannels) {
123  scenarioProbability_ = &es.getData(scenarioProbabilityToken_);
124  quality_map = &es.getData(PixelFEDChannelCollectionMapToken_);
125 
126  SiPixelQualityProbabilities::probabilityMap m_probabilities = scenarioProbability_->getProbability_Map();
127  std::vector<std::string> allScenarios;
128 
129  std::transform(quality_map->begin(),
130  quality_map->end(),
131  std::back_inserter(allScenarios),
132  [](const PixelFEDChannelCollectionMap::value_type& pair) { return pair.first; });
133 
134  std::vector<std::string> allScenariosInProb;
135 
136  for (auto it = m_probabilities.begin(); it != m_probabilities.end(); ++it) {
137  //int PUbin = it->first;
138  for (const auto& entry : it->second) {
139  auto scenario = entry.first;
140  auto probability = entry.second;
141  if (probability != 0) {
142  if (std::find(allScenariosInProb.begin(), allScenariosInProb.end(), scenario) == allScenariosInProb.end()) {
143  allScenariosInProb.push_back(scenario);
144  }
145  } // if prob!=0
146  } // loop on the scenarios for that PU bin
147  } // loop on PU bins
148 
149  std::vector<std::string> notFound;
150  std::copy_if(allScenariosInProb.begin(),
151  allScenariosInProb.end(),
152  std::back_inserter(notFound),
153  [&allScenarios](const std::string& arg) {
154  return (std::find(allScenarios.begin(), allScenarios.end(), arg) == allScenarios.end());
155  });
156 
157  if (!notFound.empty()) {
158  for (const auto& entry : notFound) {
159  LogError("SiPixelFEDChannelContainer")
160  << "The requested scenario: " << entry << " is not found in the map!! \n";
161  }
162  throw cms::Exception("SiPixelDigitizerAlgorithm") << "Found: " << notFound.size()
163  << " missing scenario(s) in SiPixelStatusScenariosRcd while "
164  "present in SiPixelStatusScenarioProbabilityRcd \n";
165  }
166  }
167 
168  TheNewSiPixelChargeReweightingAlgorithmClass->init(es);
169 }
170 
171 //=========================================================================
172 
174  : mapToken_(iC.esConsumes()),
175  geomToken_(iC.esConsumes()),
176 
177  _signal(),
178  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
179  use_ineff_from_db_(conf.getParameter<bool>("useDB")),
180  use_module_killing_(conf.getParameter<bool>("killModules")), // boolean to kill or not modules
181  use_deadmodule_DB_(conf.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
182  use_LorentzAngle_DB_(conf.getParameter<bool>("LorentzAngle_DB")), // boolean to access Lorentz angle from DB
183 
184  DeadModules(use_deadmodule_DB_ ? Parameters()
185  : conf.getParameter<Parameters>("DeadModules")), // get dead module from cfg file
186 
187  TheNewSiPixelChargeReweightingAlgorithmClass(),
188 
189  // Common pixel parameters
190  // These are parameters which are not likely to be changed
191  GeVperElectron(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
192  Sigma0(0.00037), // Charge diffusion constant 7->3.7
193  Dist300(0.0300), // normalized to 300micron Silicon
194  alpha2Order(conf.getParameter<bool>("Alpha2Order")), // switch on/off of E.B effect
195  ClusterWidth(3.), // Charge integration spread on the collection plane
196 
197  // get external parameters:
198  // To account for upgrade geometries do not assume the number
199  // of layers or disks.
200  NumberOfBarrelLayers(conf.exists("NumPixelBarrel") ? conf.getParameter<int>("NumPixelBarrel") : 3),
201  NumberOfEndcapDisks(conf.exists("NumPixelEndcap") ? conf.getParameter<int>("NumPixelEndcap") : 2),
202 
203  // ADC calibration 1adc count(135e.
204  // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
205  // Be carefull, this parameter is also used in SiPixelDet.cc to
206  // calculate the noise in adc counts from noise in electrons.
207  // Both defaults should be the same.
208  theElectronPerADC(conf.getParameter<double>("ElectronPerAdc")),
209 
210  // ADC saturation value, 255(8bit adc.
211  //theAdcFullScale(conf.getUntrackedParameter<int>("AdcFullScale",255)),
212  theAdcFullScale(conf.getParameter<int>("AdcFullScale")),
213 
214  // Noise in electrons:
215  // Pixel cell noise, relevant for generating noisy pixels
216  theNoiseInElectrons(conf.getParameter<double>("NoiseInElectrons")),
217 
218  // Fill readout noise, including all readout chain, relevant for smearing
219  //theReadoutNoise(conf.getUntrackedParameter<double>("ReadoutNoiseInElec",500.)),
220  theReadoutNoise(conf.getParameter<double>("ReadoutNoiseInElec")),
221 
222  // Pixel threshold in units of noise:
223  // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
224  // Pixel threshold in electron units.
225  theThresholdInE_FPix(conf.getParameter<double>("ThresholdInElectrons_FPix")),
226  theThresholdInE_BPix(conf.getParameter<double>("ThresholdInElectrons_BPix")),
227  theThresholdInE_BPix_L1(conf.exists("ThresholdInElectrons_BPix_L1")
228  ? conf.getParameter<double>("ThresholdInElectrons_BPix_L1")
229  : theThresholdInE_BPix),
230  theThresholdInE_BPix_L2(conf.exists("ThresholdInElectrons_BPix_L2")
231  ? conf.getParameter<double>("ThresholdInElectrons_BPix_L2")
232  : theThresholdInE_BPix),
233 
234  // Add threshold gaussian smearing:
235  theThresholdSmearing_FPix(conf.getParameter<double>("ThresholdSmearing_FPix")),
236  theThresholdSmearing_BPix(conf.getParameter<double>("ThresholdSmearing_BPix")),
237  theThresholdSmearing_BPix_L1(conf.exists("ThresholdSmearing_BPix_L1")
238  ? conf.getParameter<double>("ThresholdSmearing_BPix_L1")
239  : theThresholdSmearing_BPix),
240  theThresholdSmearing_BPix_L2(conf.exists("ThresholdSmearing_BPix_L2")
241  ? conf.getParameter<double>("ThresholdSmearing_BPix_L2")
242  : theThresholdSmearing_BPix),
243 
244  // electrons to VCAL conversion needed in misscalibrate()
245  electronsPerVCAL(conf.getParameter<double>("ElectronsPerVcal")),
246  electronsPerVCAL_Offset(conf.getParameter<double>("ElectronsPerVcal_Offset")),
247  electronsPerVCAL_L1(conf.exists("ElectronsPerVcal_L1") ? conf.getParameter<double>("ElectronsPerVcal_L1")
248  : electronsPerVCAL),
249  electronsPerVCAL_L1_Offset(conf.exists("ElectronsPerVcal_L1_Offset")
250  ? conf.getParameter<double>("ElectronsPerVcal_L1_Offset")
251  : electronsPerVCAL_Offset),
252 
253  //theTofCut 12.5, cut in particle TOD +/- 12.5ns
254  //theTofCut(conf.getUntrackedParameter<double>("TofCut",12.5)),
255  theTofLowerCut(conf.getParameter<double>("TofLowerCut")),
256  theTofUpperCut(conf.getParameter<double>("TofUpperCut")),
257 
258  // Get the Lorentz angle from the cfg file:
259  tanLorentzAnglePerTesla_FPix(use_LorentzAngle_DB_ ? 0.0
260  : conf.getParameter<double>("TanLorentzAnglePerTesla_FPix")),
261  tanLorentzAnglePerTesla_BPix(use_LorentzAngle_DB_ ? 0.0
262  : conf.getParameter<double>("TanLorentzAnglePerTesla_BPix")),
263 
264  // signal response new parameterization: split Fpix and BPix
265  FPix_p0(conf.getParameter<double>("FPix_SignalResponse_p0")),
266  FPix_p1(conf.getParameter<double>("FPix_SignalResponse_p1")),
267  FPix_p2(conf.getParameter<double>("FPix_SignalResponse_p2")),
268  FPix_p3(conf.getParameter<double>("FPix_SignalResponse_p3")),
269 
270  BPix_p0(conf.getParameter<double>("BPix_SignalResponse_p0")),
271  BPix_p1(conf.getParameter<double>("BPix_SignalResponse_p1")),
272  BPix_p2(conf.getParameter<double>("BPix_SignalResponse_p2")),
273  BPix_p3(conf.getParameter<double>("BPix_SignalResponse_p3")),
274 
275  // Add noise
276  addNoise(conf.getParameter<bool>("AddNoise")),
277 
278  // Smear the pixel charge with a gaussian which RMS is a function of the
279  // pixel charge (Danek's study)
280  addChargeVCALSmearing(conf.getParameter<bool>("ChargeVCALSmearing")),
281 
282  // Add noisy pixels
283  addNoisyPixels(conf.getParameter<bool>("AddNoisyPixels")),
284 
285  // Fluctuate charge in track subsegments
286  fluctuateCharge(conf.getUntrackedParameter<bool>("FluctuateCharge", true)),
287 
288  // Control the pixel inefficiency
289  AddPixelInefficiency(conf.getParameter<bool>("AddPixelInefficiency")),
290  KillBadFEDChannels(conf.getParameter<bool>("KillBadFEDChannels")),
291 
292  // Add threshold gaussian smearing:
293  addThresholdSmearing(conf.getParameter<bool>("AddThresholdSmearing")),
294 
295  // Get the constants for the miss-calibration studies
296  doMissCalibrate(conf.getParameter<bool>("MissCalibrate")), // Enable miss-calibration
297  theGainSmearing(conf.getParameter<double>("GainSmearing")), // sigma of the gain smearing
298  theOffsetSmearing(conf.getParameter<double>("OffsetSmearing")), //sigma of the offset smearing
299 
300  // Add pixel radiation damage for upgrade studies
301  AddPixelAging(conf.getParameter<bool>("DoPixelAging")),
302  UseReweighting(conf.getParameter<bool>("UseReweighting")),
303 
304  // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
305  //tMax(0.030), // In MeV.
306  //tMax(conf.getUntrackedParameter<double>("deltaProductionCut",0.030)),
307  tMax(conf.getParameter<double>("deltaProductionCut")),
308 
309  fluctuate(fluctuateCharge ? new SiG4UniversalFluctuation() : nullptr),
310  theNoiser(addNoise ? new GaussianTailNoiseGenerator() : nullptr),
311  calmap(doMissCalibrate ? initCal() : std::map<int, CalParameters, std::less<int> >()),
312  theSiPixelGainCalibrationService_(use_ineff_from_db_ ? new SiPixelGainCalibrationOfflineSimService(conf, iC)
313  : nullptr),
314  pixelEfficiencies_(conf, AddPixelInefficiency, NumberOfBarrelLayers, NumberOfEndcapDisks),
315  pixelAging_(conf, AddPixelAging, NumberOfBarrelLayers, NumberOfEndcapDisks) {
316  if (use_deadmodule_DB_) {
317  //string to specify SiPixelQuality label
318  SiPixelBadModuleToken_ = iC.esConsumes(edm::ESInputTag("", conf.getParameter<std::string>("SiPixelQualityLabel")));
319  }
320  if (use_LorentzAngle_DB_) {
322  }
324  // TODO: in practice the bunchspacing is known at MixingModule
325  // construction time, and thus we could declare the consumption of
326  // the actual product. In principle, however, MixingModule is
327  // capable of updating (parts of) its configuration from the
328  // EventSetup, so if that capability is really needed we'd need to
329  // invent something new (similar to mayConsume in the ESProducer
330  // side). So for now, let's consume both payloads.
332  }
333  if (KillBadFEDChannels) {
336  }
337 
338  LogInfo("PixelDigitizer ") << "SiPixelDigitizerAlgorithm constructed"
339  << "Configuration parameters:"
340  << "Threshold/Gain = "
341  << "threshold in electron FPix = " << theThresholdInE_FPix
342  << "threshold in electron BPix = " << theThresholdInE_BPix
343  << "threshold in electron BPix Layer1 = " << theThresholdInE_BPix_L1
344  << "threshold in electron BPix Layer2 = " << theThresholdInE_BPix_L2 << " "
345  << theElectronPerADC << " " << theAdcFullScale << " The delta cut-off is set to " << tMax
346  << " pix-inefficiency " << AddPixelInefficiency;
347 
348  TheNewSiPixelChargeReweightingAlgorithmClass = std::make_unique<SiPixelChargeReweightingAlgorithm>(conf, iC);
349 }
350 
351 std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> > SiPixelDigitizerAlgorithm::initCal() const {
352  std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> > calmap;
353  // Prepare for the analog amplitude miss-calibration
354  LogDebug("PixelDigitizer ") << " miss-calibrate the pixel amplitude \n";
355 
356  const bool ReadCalParameters = false;
357  if (ReadCalParameters) { // Read the calibration files from file
358  // read the calibration constants from a file (testing only)
359  std::ifstream in_file; // data file pointer
360  char filename[80] = "phCalibrationFit_C0.dat";
361 
362  in_file.open(filename, std::ios::in); // in C++
363  if (in_file.bad()) {
364  LogInfo("PixelDigitizer ") << " File not found \n ";
365  return calmap; // signal error
366  }
367  LogInfo("PixelDigitizer ") << " file opened : " << filename << "\n";
368 
369  char line[500];
370  for (int i = 0; i < 3; i++) {
371  in_file.getline(line, 500, '\n');
372  LogInfo("PixelDigitizer ") << line << "\n";
373  }
374 
375  LogInfo("PixelDigitizer ") << " test map"
376  << "\n";
377 
378  float par0, par1, par2, par3;
379  int colid, rowid;
381  // Read MC tracks
382  for (int i = 0; i < (52 * 80); i++) { // loop over tracks
383  in_file >> par0 >> par1 >> par2 >> par3 >> name >> colid >> rowid;
384  if (in_file.bad()) { // check for errors
385  LogError("PixelDigitizer") << "Cannot read data file for calmap"
386  << "\n";
387  return calmap;
388  }
389  if (in_file.eof() != 0) {
390  LogError("PixelDigitizer") << "calmap " << in_file.eof() << " " << in_file.gcount() << " " << in_file.fail()
391  << " " << in_file.good() << " end of file "
392  << "\n";
393  return calmap;
394  }
395 
396  LogDebug("PixelDigitizer ") << " line " << i << " " << par0 << " " << par1 << " " << par2 << " " << par3 << " "
397  << colid << " " << rowid << "\n";
398 
399  CalParameters onePix;
400  onePix.p0 = par0;
401  onePix.p1 = par1;
402  onePix.p2 = par2;
403  onePix.p3 = par3;
404 
405  // Convert ROC pixel index to channel
406  int chan = PixelIndices::pixelToChannelROC(rowid, colid);
407  calmap.insert(std::pair<int, CalParameters>(chan, onePix));
408 
409  // Testing the index conversion, can be skipped
410  std::pair<int, int> p = PixelIndices::channelToPixelROC(chan);
411  if (rowid != p.first)
412  LogInfo("PixelDigitizer ") << " wrong channel row " << rowid << " " << p.first << "\n";
413  if (colid != p.second)
414  LogInfo("PixelDigitizer ") << " wrong channel col " << colid << " " << p.second << "\n";
415 
416  } // pixel loop in a ROC
417 
418  LogInfo("PixelDigitizer ") << " map size " << calmap.size() << " max " << calmap.max_size() << " "
419  << calmap.empty() << "\n";
420 
421  // LogInfo("PixelDigitizer ") << " map size " << calmap.size() << "\n";
422  // map<int,CalParameters,std::less<int> >::iterator ix,it;
423  // map<int,CalParameters,std::less<int> >::const_iterator ip;
424  // for (ix = calmap.begin(); ix != calmap.end(); ++ix) {
425  // int i = (*ix).first;
426  // std::pair<int,int> p = channelToPixelROC(i);
427  // it = calmap.find(i);
428  // CalParameters y = (*it).second;
429  // CalParameters z = (*ix).second;
430  // LogInfo("PixelDigitizer ") << i <<" "<<p.first<<" "<<p.second<<" "<<y.p0<<" "<<z.p0<<" "<<calmap[i].p0<<"\n";
431 
432  // //int dummy=0;
433  // //cin>>dummy;
434  // }
435 
436  } // end if readparameters
437  return calmap;
438 } // end initCal()
439 
440 //=========================================================================
442  LogDebug("PixelDigitizer") << "SiPixelDigitizerAlgorithm deleted";
443 }
444 
445 // Read DynIneff Scale factors from Configuration
448  int NumberOfBarrelLayers,
449  int NumberOfEndcapDisks) {
450  // pixel inefficiency
451  // Don't use Hard coded values, read inefficiencies in from DB/python config or don't use any
452  int NumberOfTotLayers = NumberOfBarrelLayers + NumberOfEndcapDisks;
454  if (AddPixelInefficiency) {
455  FromConfig = conf.exists("thePixelColEfficiency_BPix1") && conf.exists("thePixelColEfficiency_BPix2") &&
456  conf.exists("thePixelColEfficiency_BPix3") && conf.exists("thePixelColEfficiency_FPix1") &&
457  conf.exists("thePixelColEfficiency_FPix2") && conf.exists("thePixelEfficiency_BPix1") &&
458  conf.exists("thePixelEfficiency_BPix2") && conf.exists("thePixelEfficiency_BPix3") &&
459  conf.exists("thePixelEfficiency_FPix1") && conf.exists("thePixelEfficiency_FPix2") &&
460  conf.exists("thePixelChipEfficiency_BPix1") && conf.exists("thePixelChipEfficiency_BPix2") &&
461  conf.exists("thePixelChipEfficiency_BPix3") && conf.exists("thePixelChipEfficiency_FPix1") &&
462  conf.exists("thePixelChipEfficiency_FPix2");
463  if (NumberOfBarrelLayers == 3)
464  FromConfig = FromConfig && conf.exists("theLadderEfficiency_BPix1") && conf.exists("theLadderEfficiency_BPix2") &&
465  conf.exists("theLadderEfficiency_BPix3") && conf.exists("theModuleEfficiency_BPix1") &&
466  conf.exists("theModuleEfficiency_BPix2") && conf.exists("theModuleEfficiency_BPix3") &&
467  conf.exists("thePUEfficiency_BPix1") && conf.exists("thePUEfficiency_BPix2") &&
468  conf.exists("thePUEfficiency_BPix3") && conf.exists("theInnerEfficiency_FPix1") &&
469  conf.exists("theInnerEfficiency_FPix2") && conf.exists("theOuterEfficiency_FPix1") &&
470  conf.exists("theOuterEfficiency_FPix2") && conf.exists("thePUEfficiency_FPix_Inner") &&
471  conf.exists("thePUEfficiency_FPix_Outer") && conf.exists("theInstLumiScaleFactor");
472  if (NumberOfBarrelLayers >= 4)
473  FromConfig = FromConfig && conf.exists("thePixelColEfficiency_BPix4") &&
474  conf.exists("thePixelEfficiency_BPix4") && conf.exists("thePixelChipEfficiency_BPix4");
475  if (NumberOfEndcapDisks >= 3)
476  FromConfig = FromConfig && conf.exists("thePixelColEfficiency_FPix3") &&
477  conf.exists("thePixelEfficiency_FPix3") && conf.exists("thePixelChipEfficiency_FPix3");
478  if (FromConfig) {
479  LogInfo("PixelDigitizer ") << "The PixelDigitizer inefficiency configuration is read from the config file.\n";
480  theInstLumiScaleFactor = conf.getParameter<double>("theInstLumiScaleFactor");
481  int i = 0;
482  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix1");
483  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix2");
484  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix3");
485  if (NumberOfBarrelLayers >= 4) {
486  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix4");
487  }
488  //
489  i = 0;
490  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix1");
491  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix2");
492  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix3");
493  if (NumberOfBarrelLayers >= 4) {
494  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix4");
495  }
496  //
497  i = 0;
498  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix1");
499  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix2");
500  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix3");
501  if (NumberOfBarrelLayers >= 4) {
502  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix4");
503  }
504  //
505  if (NumberOfBarrelLayers == 3) {
506  i = 0;
507  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix1");
508  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix2");
509  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix3");
510  if (((theLadderEfficiency_BPix[0].size() != 20) || (theLadderEfficiency_BPix[1].size() != 32) ||
511  (theLadderEfficiency_BPix[2].size() != 44)) &&
512  (NumberOfBarrelLayers == 3))
513  throw cms::Exception("Configuration") << "Wrong ladder number in efficiency config!";
514  //
515  i = 0;
516  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix1");
517  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix2");
518  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix3");
519  if (((theModuleEfficiency_BPix[0].size() != 4) || (theModuleEfficiency_BPix[1].size() != 4) ||
520  (theModuleEfficiency_BPix[2].size() != 4)) &&
521  (NumberOfBarrelLayers == 3))
522  throw cms::Exception("Configuration") << "Wrong module number in efficiency config!";
523  //
524  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix1"));
525  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix2"));
526  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix3"));
527  if (((thePUEfficiency[0].empty()) || (thePUEfficiency[1].empty()) || (thePUEfficiency[2].empty())) &&
528  (NumberOfBarrelLayers == 3))
529  throw cms::Exception("Configuration")
530  << "At least one PU efficiency (BPix) number is needed in efficiency config!";
531  }
532  // The next is needed for Phase2 Tracker studies
533  if (NumberOfBarrelLayers >= 5) {
534  if (NumberOfTotLayers > 20) {
535  throw cms::Exception("Configuration") << "SiPixelDigitizer was given more layers than it can handle";
536  }
537  // For Phase2 tracker layers just set the outermost BPix inefficiency to 99.9% THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
538  for (int j = 5; j <= NumberOfBarrelLayers; j++) {
539  thePixelColEfficiency[j - 1] = 0.999;
540  thePixelEfficiency[j - 1] = 0.999;
541  thePixelChipEfficiency[j - 1] = 0.999;
542  }
543  }
544  //
545  i = FPixIndex;
546  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix1");
547  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix2");
548  if (NumberOfEndcapDisks >= 3) {
549  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix3");
550  }
551  i = FPixIndex;
552  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix1");
553  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix2");
554  if (NumberOfEndcapDisks >= 3) {
555  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix3");
556  }
557  i = FPixIndex;
558  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix1");
559  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix2");
560  if (NumberOfEndcapDisks >= 3) {
561  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix3");
562  }
563  // The next is needed for Phase2 Tracker studies
564  if (NumberOfEndcapDisks >= 4) {
565  if (NumberOfTotLayers > 20) {
566  throw cms::Exception("Configuration") << "SiPixelDigitizer was given more layers than it can handle";
567  }
568  // For Phase2 tracker layers just set the extra FPix disk inefficiency to 99.9% THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
569  for (int j = 4 + FPixIndex; j <= NumberOfEndcapDisks + NumberOfBarrelLayers; j++) {
570  thePixelColEfficiency[j - 1] = 0.999;
571  thePixelEfficiency[j - 1] = 0.999;
572  thePixelChipEfficiency[j - 1] = 0.999;
573  }
574  }
575  //FPix Dynamic Inefficiency
576  if (NumberOfBarrelLayers == 3) {
577  i = FPixIndex;
578  theInnerEfficiency_FPix[i++] = conf.getParameter<double>("theInnerEfficiency_FPix1");
579  theInnerEfficiency_FPix[i++] = conf.getParameter<double>("theInnerEfficiency_FPix2");
580  i = FPixIndex;
581  theOuterEfficiency_FPix[i++] = conf.getParameter<double>("theOuterEfficiency_FPix1");
582  theOuterEfficiency_FPix[i++] = conf.getParameter<double>("theOuterEfficiency_FPix2");
583  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_FPix_Inner"));
584  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_FPix_Outer"));
585  if (((thePUEfficiency[3].empty()) || (thePUEfficiency[4].empty())) && (NumberOfEndcapDisks == 2))
586  throw cms::Exception("Configuration")
587  << "At least one (FPix) PU efficiency number is needed in efficiency config!";
588  pu_scale.resize(thePUEfficiency.size());
589  }
590  } else
591  LogInfo("PixelDigitizer ") << "The PixelDigitizer inefficiency configuration is read from the database.\n";
592  }
593  // the first "NumberOfBarrelLayers" settings [0],[1], ... , [NumberOfBarrelLayers-1] are for the barrel pixels
594  // the next "NumberOfEndcapDisks" settings [NumberOfBarrelLayers],[NumberOfBarrelLayers+1], ... [NumberOfEndcapDisks+NumberOfBarrelLayers-1]
595 }
596 
597 // Read DynIneff Scale factors from DB
599  LogDebug("PixelDigitizer ") << " In SiPixelDigitizerAlgorithm::init_DynIneffDB " << AddPixelInefficiency << " "
600  << pixelEfficiencies_.FromConfig << "\n";
604  }
605 }
606 
609  theInstLumiScaleFactor = SiPixelDynamicInefficiency->gettheInstLumiScaleFactor();
610  const std::map<uint32_t, double>& PixelGeomFactorsDBIn = SiPixelDynamicInefficiency->getPixelGeomFactors();
611  const std::map<uint32_t, double>& ColGeomFactorsDB = SiPixelDynamicInefficiency->getColGeomFactors();
612  const std::map<uint32_t, double>& ChipGeomFactorsDB = SiPixelDynamicInefficiency->getChipGeomFactors();
613  const std::map<uint32_t, std::vector<double> >& PUFactors = SiPixelDynamicInefficiency->getPUFactors();
614  std::vector<uint32_t> DetIdmasks = SiPixelDynamicInefficiency->getDetIdmasks();
615 
616  // Loop on all modules, initialize map for easy access
617  for (const auto& it_module : geom->detUnits()) {
618  if (dynamic_cast<PixelGeomDetUnit const*>(it_module) == nullptr)
619  continue;
620  const DetId detid = it_module->geographicalId();
621  uint32_t rawid = detid.rawId();
622  PixelGeomFactors[rawid] = 1;
623  ColGeomFactors[rawid] = 1;
624  ChipGeomFactors[rawid] = 1;
625  PixelGeomFactorsROCStdPixels[rawid] = std::vector<double>(16, 1);
626  PixelGeomFactorsROCBigPixels[rawid] = std::vector<double>(16, 1);
627  }
628 
629  // ROC level inefficiency for phase 1 (disentangle scale factors for big and std size pixels)
630  std::map<uint32_t, double> PixelGeomFactorsDB;
631 
632  LogDebug("PixelDigitizer ") << " Check PixelEfficiencies -- PixelGeomFactorsDBIn "
633  << "\n";
635  for (auto db_factor : PixelGeomFactorsDBIn) {
636  LogDebug("PixelDigitizer ") << " db_factor " << db_factor.first << " " << db_factor.second << "\n";
637 
638  int shift = DetId(db_factor.first).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel) ? BPixRocIdShift
639  : FPixRocIdShift;
640  unsigned int rocMask = rocIdMaskBits << shift;
641  unsigned int rocId = (((db_factor.first) & rocMask) >> shift);
642  if (rocId != 0) {
643  rocId--;
644  unsigned int rawid = db_factor.first & (~rocMask);
645  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(geom->idToDet(rawid));
646  PixelTopology const* topology = &(theGeomDet->specificTopology());
647  const int nPixelsInROC = topology->rowsperroc() * topology->colsperroc();
648  const int nBigPixelsInROC = 2 * topology->rowsperroc() + topology->colsperroc() - 2;
649  double factor = db_factor.second;
650  double badFraction = 1 - factor;
651  double bigPixelFraction = static_cast<double>(nBigPixelsInROC) / nPixelsInROC;
652  double stdPixelFraction = 1. - bigPixelFraction;
653 
654  double badFractionBig = std::min(bigPixelFraction, badFraction);
655  double badFractionStd = std::max(0., badFraction - badFractionBig);
656  double badFractionBigReNormalized = badFractionBig / bigPixelFraction;
657  double badFractionStdReNormalized = badFractionStd / stdPixelFraction;
658  PixelGeomFactorsROCStdPixels[rawid][rocId] *= (1. - badFractionStdReNormalized);
659  PixelGeomFactorsROCBigPixels[rawid][rocId] *= (1. - badFractionBigReNormalized);
660  } else {
661  PixelGeomFactorsDB[db_factor.first] = db_factor.second;
662  }
663  }
664  } // is Phase 1 geometry
665  else {
666  PixelGeomFactorsDB = PixelGeomFactorsDBIn;
667  }
668 
669  LogDebug("PixelDigitizer ")
670  << " Check PixelEfficiencies -- Loop on all modules and store module level geometrical scale factors "
671  << "\n";
672  // Loop on all modules, store module level geometrical scale factors
673  for (const auto& it_module : geom->detUnits()) {
674  if (dynamic_cast<PixelGeomDetUnit const*>(it_module) == nullptr)
675  continue;
676  const DetId detid = it_module->geographicalId();
677  uint32_t rawid = detid.rawId();
678  for (auto db_factor : PixelGeomFactorsDB) {
679  LogDebug("PixelDigitizer ") << " db_factor PixelGeomFactorsDB " << db_factor.first << " "
680  << db_factor.second << "\n";
681  if (matches(detid, DetId(db_factor.first), DetIdmasks))
682  PixelGeomFactors[rawid] *= db_factor.second;
683  }
684  for (auto db_factor : ColGeomFactorsDB) {
685  LogDebug("PixelDigitizer ") << " db_factor ColGeomFactorsDB " << db_factor.first << " " << db_factor.second
686  << "\n";
687  if (matches(detid, DetId(db_factor.first), DetIdmasks))
688  ColGeomFactors[rawid] *= db_factor.second;
689  }
690  for (auto db_factor : ChipGeomFactorsDB) {
691  LogDebug("PixelDigitizer ") << " db_factor ChipGeomFactorsDB " << db_factor.first << " "
692  << db_factor.second << "\n";
693  if (matches(detid, DetId(db_factor.first), DetIdmasks))
694  ChipGeomFactors[rawid] *= db_factor.second;
695  }
696  }
697 
698  // piluep scale factors are calculated once per event
699  // therefore vector index is stored in a map for each module that matches to a db_id
700  size_t i = 0;
701  LogDebug("PixelDigitizer ") << " Check PixelEfficiencies -- PUFactors "
702  << "\n";
703  for (const auto& factor : PUFactors) {
704  //
705  LogDebug("PixelDigitizer ") << " factor " << factor.first << " " << factor.second.size() << "\n";
706  for (size_t i = 0, n = factor.second.size(); i < n; i++) {
707  LogDebug("PixelDigitizer ") << " print factor.second for " << i << " " << factor.second[i] << "\n";
708  }
709  //
710  const DetId db_id = DetId(factor.first);
711  for (const auto& it_module : geom->detUnits()) {
712  if (dynamic_cast<PixelGeomDetUnit const*>(it_module) == nullptr)
713  continue;
714  const DetId detid = it_module->geographicalId();
715  if (!matches(detid, db_id, DetIdmasks))
716  continue;
717  if (iPU.count(detid.rawId())) {
718  throw cms::Exception("Database")
719  << "Multiple db_ids match to same module in SiPixelDynamicInefficiency DB Object";
720  } else {
721  iPU[detid.rawId()] = i;
722  }
723  }
724  thePUEfficiency.push_back(factor.second);
725  ++i;
726  }
727  pu_scale.resize(thePUEfficiency.size());
728 }
729 
731  const DetId& db_id,
732  const std::vector<uint32_t>& DetIdmasks) {
733  if (detid.subdetId() != db_id.subdetId())
734  return false;
735  for (size_t i = 0; i < DetIdmasks.size(); ++i) {
736  DetId maskid = DetId(DetIdmasks.at(i));
737  if (maskid.subdetId() != db_id.subdetId())
738  continue;
739  if ((detid.rawId() & maskid.rawId()) != (db_id.rawId() & maskid.rawId()) &&
740  (db_id.rawId() & maskid.rawId()) != DetId(db_id.det(), db_id.subdetId()).rawId())
741  return false;
742  }
743  return true;
744 }
745 
747  bool AddAging,
749  int NumberOfEndcapDisks) {
750  // pixel aging
751  // Don't use Hard coded values, read aging in from python or don't use any
752  if (AddAging) {
753  int NumberOfTotLayers = NumberOfBarrelLayers + NumberOfEndcapDisks;
754  FPixIndex = NumberOfBarrelLayers;
755 
756  int i = 0;
757  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix1");
758  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix2");
759  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix3");
760  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix4");
761 
762  // to be removed when Gaelle will have the phase2 digitizer
763  if (NumberOfBarrelLayers >= 5) {
764  if (NumberOfTotLayers > 20) {
765  throw cms::Exception("Configuration") << "SiPixelDigitizer was given more layers than it can handle";
766  }
767  // For Phase2 tracker layers just set the outermost BPix aging 0.
768  for (int j = 5; j <= NumberOfBarrelLayers; j++) {
769  thePixelPseudoRadDamage[j - 1] = 0.;
770  }
771  }
772  //
773  i = FPixIndex;
774  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix1");
775  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix2");
776  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix3");
777 
778  //To be removed when Phase2 digitizer will be available
779  if (NumberOfEndcapDisks >= 4) {
780  if (NumberOfTotLayers > 20) {
781  throw cms::Exception("Configuration") << "SiPixelDigitizer was given more layers than it can handle";
782  }
783  // For Phase2 tracker layers just set the extra FPix disk aging to 0. BE CAREFUL THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
784  for (int j = 4 + FPixIndex; j <= NumberOfEndcapDisks + NumberOfBarrelLayers; j++) {
785  thePixelPseudoRadDamage[j - 1] = 0.;
786  }
787  }
788  }
789  // the first "NumberOfBarrelLayers" settings [0],[1], ... , [NumberOfBarrelLayers-1] are for the barrel pixels
790  // the next "NumberOfEndcapDisks" settings [NumberOfBarrelLayers],[NumberOfBarrelLayers+1], ... [NumberOfEndcapDisks+NumberOfBarrelLayers-1]
791 }
792 
793 //=========================================================================
794 void SiPixelDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
795  std::vector<PSimHit>::const_iterator inputEnd,
796  const size_t inputBeginGlobalIndex,
797  const unsigned int tofBin,
798  const PixelGeomDetUnit* pixdet,
799  const GlobalVector& bfield,
800  const TrackerTopology* tTopo,
801  CLHEP::HepRandomEngine* engine) {
802  // produce SignalPoint's for all SimHit's in detector
803  // Loop over hits
804 
805  uint32_t detId = pixdet->geographicalId().rawId();
806  size_t simHitGlobalIndex = inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
807  for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; ssbegin != inputEnd; ++ssbegin, ++simHitGlobalIndex) {
808  // skip hits not in this detector.
809  if ((*ssbegin).detUnitId() != detId) {
810  continue;
811  }
812 
813 #ifdef TP_DEBUG
814  LogDebug("Pixel Digitizer") << (*ssbegin).particleType() << " " << (*ssbegin).pabs() << " "
815  << (*ssbegin).energyLoss() << " " << (*ssbegin).tof() << " " << (*ssbegin).trackId()
816  << " " << (*ssbegin).processType() << " " << (*ssbegin).detUnitId()
817  << (*ssbegin).entryPoint() << " " << (*ssbegin).exitPoint();
818 #endif
819 
820  std::vector<EnergyDepositUnit> ionization_points;
821  std::vector<SignalPoint> collection_points;
822 
823  // fill collection_points for this SimHit, indpendent of topology
824  // Check the TOF cut
825  if (((*ssbegin).tof() - pixdet->surface().toGlobal((*ssbegin).localPosition()).mag() / 30.) >= theTofLowerCut &&
826  ((*ssbegin).tof() - pixdet->surface().toGlobal((*ssbegin).localPosition()).mag() / 30.) <= theTofUpperCut) {
827  primary_ionization(*ssbegin, ionization_points, engine); // fills _ionization_points
828  drift(*ssbegin,
829  pixdet,
830  bfield,
831  tTopo,
832  ionization_points,
833  collection_points); // transforms _ionization_points to collection_points
834  // compute induced signal on readout elements and add to _signal
835  induce_signal(inputBegin,
836  inputEnd,
837  *ssbegin,
838  simHitGlobalIndex,
839  inputBeginGlobalIndex,
840  tofBin,
841  pixdet,
842  collection_points); // 1st 3 args needed only for SimHit<-->Digi link
843  } // end if
844  } // end for
845 }
846 
847 //============================================================================
849  //Instlumi scalefactor calculating for dynamic inefficiency
850 
851  if (puInfo) {
852  const std::vector<int>& bunchCrossing = puInfo->getMix_bunchCrossing();
853  const std::vector<float>& TrueInteractionList = puInfo->getMix_TrueInteractions();
854  //const int bunchSpacing = puInfo->getMix_bunchSpacing();
855 
856  int pui = 0, p = 0;
857  std::vector<int>::const_iterator pu;
858  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
859 
860  for (pu = bunchCrossing.begin(); pu != bunchCrossing.end(); ++pu) {
861  if (*pu == 0) {
862  pu0 = pu;
863  p = pui;
864  }
865  pui++;
866  }
867  if (pu0 != bunchCrossing.end()) {
868  for (size_t i = 0, n = pixelEfficiencies_.thePUEfficiency.size(); i < n; i++) {
869  double instlumi = TrueInteractionList.at(p) * pixelEfficiencies_.theInstLumiScaleFactor;
870  double instlumi_pow = 1.;
872  for (size_t j = 0; j < pixelEfficiencies_.thePUEfficiency[i].size(); j++) {
874  instlumi_pow *= instlumi;
875  }
876  }
877  }
878  } else {
879  for (int i = 0, n = pixelEfficiencies_.thePUEfficiency.size(); i < n; i++) {
881  }
882  }
883 }
884 
885 //============================================================================
886 void SiPixelDigitizerAlgorithm::calculateInstlumiFactor(const std::vector<PileupSummaryInfo>& ps, int bunchSpacing) {
887  int p = -1;
888  for (unsigned int i = 0; i < ps.size(); i++)
889  if (ps[i].getBunchCrossing() == 0)
890  p = i;
891 
892  if (p >= 0) {
893  for (size_t i = 0, n = pixelEfficiencies_.thePUEfficiency.size(); i < n; i++) {
894  double instlumi = ps[p].getTrueNumInteractions() * pixelEfficiencies_.theInstLumiScaleFactor;
895  double instlumi_pow = 1.;
897  for (size_t j = 0; j < pixelEfficiencies_.thePUEfficiency[i].size(); j++) {
899  instlumi_pow *= instlumi;
900  }
901  }
902  } else {
903  for (int i = 0, n = pixelEfficiencies_.thePUEfficiency.size(); i < n; i++) {
905  }
906  }
907 }
908 
909 // ========== StuckTBMs
910 
912 
913 std::unique_ptr<PixelFEDChannelCollection> SiPixelDigitizerAlgorithm::chooseScenario(
914  const std::vector<PileupSummaryInfo>& ps, CLHEP::HepRandomEngine* engine) {
915  std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_ = nullptr;
917 
918  std::vector<int> bunchCrossing;
919  std::vector<float> TrueInteractionList;
920 
921  for (unsigned int i = 0; i < ps.size(); i++) {
922  bunchCrossing.push_back(ps[i].getBunchCrossing());
923  TrueInteractionList.push_back(ps[i].getTrueNumInteractions());
924  }
925 
926  int pui = 0, p = 0;
927  std::vector<int>::const_iterator pu;
928  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
929 
930  for (pu = bunchCrossing.begin(); pu != bunchCrossing.end(); ++pu) {
931  if (*pu == 0) {
932  pu0 = pu;
933  p = pui;
934  }
935  pui++;
936  }
937 
938  if (pu0 != bunchCrossing.end()) {
939  unsigned int PUBin = TrueInteractionList.at(p); // case delta PU=1, fix me
940  const auto& theProbabilitiesPerScenario = scenarioProbability_->getProbabilities(PUBin);
941  std::vector<double> probabilities;
942  probabilities.reserve(theProbabilitiesPerScenario.size());
943  for (auto it = theProbabilitiesPerScenario.begin(); it != theProbabilitiesPerScenario.end(); it++) {
944  probabilities.push_back(it->second);
945  }
946 
947  CLHEP::RandGeneral randGeneral(*engine, &(probabilities.front()), probabilities.size());
948  double x = randGeneral.shoot();
949  unsigned int index = x * probabilities.size() - 1;
950  const std::string& scenario = theProbabilitiesPerScenario.at(index).first;
951 
952  PixelFEDChannelCollection_ = std::make_unique<PixelFEDChannelCollection>(quality_map->at(scenario));
954  std::make_unique<PixelFEDChannelCollection>(quality_map->at(scenario));
955  }
956 
957  return PixelFEDChannelCollection_;
958 }
959 
960 std::unique_ptr<PixelFEDChannelCollection> SiPixelDigitizerAlgorithm::chooseScenario(PileupMixingContent* puInfo,
961  CLHEP::HepRandomEngine* engine) {
962  //Determine scenario to use for the current event based on pileup information
963 
964  std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_ = nullptr;
966  if (puInfo) {
967  const std::vector<int>& bunchCrossing = puInfo->getMix_bunchCrossing();
968  const std::vector<float>& TrueInteractionList = puInfo->getMix_TrueInteractions();
969 
970  int pui = 0, p = 0;
971  std::vector<int>::const_iterator pu;
972  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
973 
974  for (pu = bunchCrossing.begin(); pu != bunchCrossing.end(); ++pu) {
975  if (*pu == 0) {
976  pu0 = pu;
977  p = pui;
978  }
979  pui++;
980  }
981 
982  if (pu0 != bunchCrossing.end()) {
983  unsigned int PUBin = TrueInteractionList.at(p); // case delta PU=1, fix me
984  const auto& theProbabilitiesPerScenario = scenarioProbability_->getProbabilities(PUBin);
985  std::vector<double> probabilities;
986  probabilities.reserve(theProbabilitiesPerScenario.size());
987  for (auto it = theProbabilitiesPerScenario.begin(); it != theProbabilitiesPerScenario.end(); it++) {
988  probabilities.push_back(it->second);
989  }
990 
991  CLHEP::RandGeneral randGeneral(*engine, &(probabilities.front()), probabilities.size());
992  double x = randGeneral.shoot();
993  unsigned int index = x * probabilities.size() - 1;
994  const std::string& scenario = theProbabilitiesPerScenario.at(index).first;
995 
996  PixelFEDChannelCollection_ = std::make_unique<PixelFEDChannelCollection>(quality_map->at(scenario));
998  std::make_unique<PixelFEDChannelCollection>(quality_map->at(scenario));
999  }
1000  }
1001  return PixelFEDChannelCollection_;
1002 }
1003 
1004 //============================================================================
1005 void SiPixelDigitizerAlgorithm::setSimAccumulator(const std::map<uint32_t, std::map<int, int> >& signalMap) {
1006  for (const auto& det : signalMap) {
1007  auto& theSignal = _signal[det.first];
1008  for (const auto& chan : det.second) {
1009  theSignal[chan.first].set(chan.second *
1010  theElectronPerADC); // will get divided again by theElectronPerAdc in digitize...
1011  }
1012  }
1013 }
1014 
1015 //============================================================================
1017  std::vector<PixelDigi>& digis,
1018  std::vector<PixelDigiSimLink>& simlinks,
1019  const TrackerTopology* tTopo,
1020  CLHEP::HepRandomEngine* engine) {
1021  // Pixel Efficiency moved from the constructor to this method because
1022  // the information of the det are not available in the constructor
1023  // Efficiency parameters. 0 - no inefficiency, 1-low lumi, 10-high lumi
1024 
1025  uint32_t detID = pixdet->geographicalId().rawId();
1026  const signal_map_type& theSignal = _signal[detID];
1027 
1028  // Noise already defined in electrons
1029  // thePixelThresholdInE = thePixelThreshold * theNoiseInElectrons ;
1030  // Find the threshold in noise units, needed for the noiser.
1031 
1032  float thePixelThresholdInE = 0.;
1033 
1034  if (theNoiseInElectrons > 0.) {
1035  if (pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()) { // Barrel modules
1036  int lay = tTopo->layer(detID);
1037  if (addThresholdSmearing) {
1040  if (lay == 1) {
1041  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(
1042  engine, theThresholdInE_BPix_L1, theThresholdSmearing_BPix_L1); // gaussian smearing
1043  } else if (lay == 2) {
1044  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(
1045  engine, theThresholdInE_BPix_L2, theThresholdSmearing_BPix_L2); // gaussian smearing
1046  } else {
1047  thePixelThresholdInE =
1048  CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix, theThresholdSmearing_BPix); // gaussian smearing
1049  }
1050  }
1051  } else {
1054  if (lay == 1) {
1055  thePixelThresholdInE = theThresholdInE_BPix_L1;
1056  } else if (lay == 2) {
1057  thePixelThresholdInE = theThresholdInE_BPix_L2;
1058  } else {
1059  thePixelThresholdInE = theThresholdInE_BPix; // no smearing
1060  }
1061  }
1062  }
1063  } else if (pixdet->type().isTrackerPixel()) { // Forward disks modules
1064  if (addThresholdSmearing) {
1065  thePixelThresholdInE =
1066  CLHEP::RandGaussQ::shoot(engine, theThresholdInE_FPix, theThresholdSmearing_FPix); // gaussian smearing
1067  } else {
1068  thePixelThresholdInE = theThresholdInE_FPix; // no smearing
1069  }
1070  } else {
1071  throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;
1072  }
1073  }
1074 
1075 #ifdef TP_DEBUG
1076  const PixelTopology* topol = &pixdet->specificTopology();
1077  int numColumns = topol->ncolumns(); // det module number of cols&rows
1078  int numRows = topol->nrows();
1079  // full detector thickness
1080  float moduleThickness = pixdet->specificSurface().bounds().thickness();
1081  LogDebug("PixelDigitizer") << " PixelDigitizer " << numColumns << " " << numRows << " " << moduleThickness;
1082 #endif
1083 
1084  if (addNoise)
1085  add_noise(pixdet, thePixelThresholdInE / theNoiseInElectrons, engine); // generate noise
1086 
1087  // Do only if needed
1088 
1089  if ((AddPixelInefficiency) && (!theSignal.empty()))
1090  pixel_inefficiency(pixelEfficiencies_, pixdet, tTopo, engine); // Kill some pixels
1091 
1092  if (use_ineff_from_db_ && (!theSignal.empty()))
1093  pixel_inefficiency_db(detID);
1094 
1095  if (use_module_killing_) {
1096  if (use_deadmodule_DB_) { // remove dead modules using DB
1097  module_killing_DB(detID);
1098  } else { // remove dead modules using the list in cfg file
1099  module_killing_conf(detID);
1100  }
1101  }
1102 
1103  make_digis(thePixelThresholdInE, detID, pixdet, digis, simlinks, tTopo);
1104 
1105 #ifdef TP_DEBUG
1106  LogDebug("PixelDigitizer") << "[SiPixelDigitizerAlgorithm] converted " << digis.size() << " PixelDigis in DetUnit"
1107  << detID;
1108 #endif
1109 }
1110 
1111 //***********************************************************************/
1112 // Generate primary ionization along the track segment.
1113 // Divide the track into small sub-segments
1115  std::vector<EnergyDepositUnit>& ionization_points,
1116  CLHEP::HepRandomEngine* engine) const {
1117  // Straight line approximation for trajectory inside active media
1118 
1119  const float SegmentLength = 0.0010; //10microns in cm
1120  float energy;
1121 
1122  // Get the 3D segment direction vector
1123  LocalVector direction = hit.exitPoint() - hit.entryPoint();
1124 
1125  float eLoss = hit.energyLoss(); // Eloss in GeV
1126  float length = direction.mag(); // Track length in Silicon
1127 
1128  int NumberOfSegments = int(length / SegmentLength); // Number of segments
1129  if (NumberOfSegments < 1)
1130  NumberOfSegments = 1;
1131 
1132 #ifdef TP_DEBUG
1133  LogDebug("Pixel Digitizer") << " enter primary_ionzation " << NumberOfSegments
1134  << " shift = " << (hit.exitPoint().x() - hit.entryPoint().x()) << " "
1135  << (hit.exitPoint().y() - hit.entryPoint().y()) << " "
1136  << (hit.exitPoint().z() - hit.entryPoint().z()) << " " << hit.particleType() << " "
1137  << hit.pabs();
1138 #endif
1139 
1140  float* elossVector = new float[NumberOfSegments]; // Eloss vector
1141 
1142  if (fluctuateCharge) {
1143  //MP DA RIMUOVERE ASSOLUTAMENTE
1144  int pid = hit.particleType();
1145  //int pid=211; // assume it is a pion
1146 
1147  float momentum = hit.pabs();
1148  // Generate fluctuated charge points
1149  fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments, elossVector, engine);
1150  }
1151 
1152  ionization_points.resize(NumberOfSegments); // set size
1153 
1154  // loop over segments
1155  for (int i = 0; i != NumberOfSegments; i++) {
1156  // Divide the segment into equal length subsegments
1157  Local3DPoint point = hit.entryPoint() + float((i + 0.5) / NumberOfSegments) * direction;
1158 
1159  if (fluctuateCharge)
1160  energy = elossVector[i] / GeVperElectron; // Convert charge to elec.
1161  else
1162  energy = hit.energyLoss() / GeVperElectron / float(NumberOfSegments);
1163 
1164  EnergyDepositUnit edu(energy, point); //define position,energy point
1165  ionization_points[i] = edu; // save
1166 
1167 #ifdef TP_DEBUG
1168  LogDebug("Pixel Digitizer") << i << " " << ionization_points[i].x() << " " << ionization_points[i].y() << " "
1169  << ionization_points[i].z() << " " << ionization_points[i].energy();
1170 #endif
1171 
1172  } // end for loop
1173 
1174  delete[] elossVector;
1175 }
1176 //******************************************************************************
1177 
1178 // Fluctuate the charge comming from a small (10um) track segment.
1179 // Use the G4 routine. For mip pions for the moment.
1181  float particleMomentum,
1182  float eloss,
1183  float length,
1184  int NumberOfSegs,
1185  float elossVector[],
1186  CLHEP::HepRandomEngine* engine) const {
1187  // Get dedx for this track
1188  //float dedx;
1189  //if( length > 0.) dedx = eloss/length;
1190  //else dedx = eloss;
1191 
1192  double particleMass = 139.6; // Mass in MeV, Assume pion
1193  pid = std::abs(pid);
1194  if (pid != 211) { // Mass in MeV
1195  if (pid == 11)
1196  particleMass = 0.511;
1197  else if (pid == 13)
1198  particleMass = 105.7;
1199  else if (pid == 321)
1200  particleMass = 493.7;
1201  else if (pid == 2212)
1202  particleMass = 938.3;
1203  }
1204  // What is the track segment length.
1205  float segmentLength = length / NumberOfSegs;
1206 
1207  // Generate charge fluctuations.
1208  float de = 0.;
1209  float sum = 0.;
1210  double segmentEloss = (1000. * eloss) / NumberOfSegs; //eloss in MeV
1211  for (int i = 0; i < NumberOfSegs; i++) {
1212  // material,*, momentum,energy,*, *, mass
1213  //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
1214  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
1215  // track segment length in mm, segment eloss in MeV
1216  // Returns fluctuated eloss in MeV
1217  double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
1218  de = fluctuate->SampleFluctuations(double(particleMomentum * 1000.),
1219  particleMass,
1220  deltaCutoff,
1221  double(segmentLength * 10.),
1222  segmentEloss,
1223  engine) /
1224  1000.; //convert to GeV
1225  elossVector[i] = de;
1226  sum += de;
1227  }
1228 
1229  if (sum > 0.) { // If fluctuations give eloss>0.
1230  // Rescale to the same total eloss
1231  float ratio = eloss / sum;
1232 
1233  for (int ii = 0; ii < NumberOfSegs; ii++)
1234  elossVector[ii] = ratio * elossVector[ii];
1235  } else { // If fluctuations gives 0 eloss
1236  float averageEloss = eloss / NumberOfSegs;
1237  for (int ii = 0; ii < NumberOfSegs; ii++)
1238  elossVector[ii] = averageEloss;
1239  }
1240  return;
1241 }
1242 
1243 //*******************************************************************************
1244 // Drift the charge segments to the sensor surface (collection plane)
1245 // Include the effect of E-field and B-field
1247  const PixelGeomDetUnit* pixdet,
1248  const GlobalVector& bfield,
1249  const TrackerTopology* tTopo,
1250  const std::vector<EnergyDepositUnit>& ionization_points,
1251  std::vector<SignalPoint>& collection_points) const {
1252 #ifdef TP_DEBUG
1253  LogDebug("Pixel Digitizer") << " enter drift ";
1254 #endif
1255 
1256  collection_points.resize(ionization_points.size()); // set size
1257 
1258  LocalVector driftDir = DriftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
1259  if (driftDir.z() == 0.) {
1260  LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
1261  return;
1262  }
1263 
1264  // tangent of Lorentz angle
1265  //float TanLorenzAngleX = driftDir.x()/driftDir.z();
1266  //float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
1267 
1268  float TanLorenzAngleX, TanLorenzAngleY, dir_z, CosLorenzAngleX, CosLorenzAngleY;
1269  if (alpha2Order) {
1270  TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
1271  TanLorenzAngleY = driftDir.y();
1272  dir_z = driftDir.z(); // The z drift direction
1273  CosLorenzAngleX = 1. / sqrt(1. + TanLorenzAngleX * TanLorenzAngleX); //cosine
1274  CosLorenzAngleY = 1. / sqrt(1. + TanLorenzAngleY * TanLorenzAngleY); //cosine;
1275 
1276  } else {
1277  TanLorenzAngleX = driftDir.x();
1278  TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
1279  dir_z = driftDir.z(); // The z drift direction
1280  CosLorenzAngleX = 1. / sqrt(1. + TanLorenzAngleX * TanLorenzAngleX); //cosine to estimate the path length
1281  CosLorenzAngleY = 1.;
1282  }
1283 
1284  float moduleThickness = pixdet->specificSurface().bounds().thickness();
1285 #ifdef TP_DEBUG
1286  LogDebug("Pixel Digitizer") << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY << " " << CosLorenzAngleX
1287  << " " << CosLorenzAngleY << " " << moduleThickness * TanLorenzAngleX << " " << driftDir;
1288 #endif
1289 
1290  float Sigma_x = 1.; // Charge spread
1291  float Sigma_y = 1.;
1292  float DriftDistance; // Distance between charge generation and collection
1293  float DriftLength; // Actual Drift Lentgh
1294  float Sigma;
1295 
1296  for (unsigned int i = 0; i != ionization_points.size(); i++) {
1297  float SegX, SegY, SegZ; // position
1298  SegX = ionization_points[i].x();
1299  SegY = ionization_points[i].y();
1300  SegZ = ionization_points[i].z();
1301 
1302  // Distance from the collection plane
1303  //DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
1304  // Include explixitely the E drift direction (for CMS dir_z=-1)
1305  DriftDistance = moduleThickness / 2. - (dir_z * SegZ); // Drift to -z
1306 
1307  if (DriftDistance <= 0.)
1308  LogDebug("PixelDigitizer ") << " <=0 " << DriftDistance << " " << i << " " << SegZ << " " << dir_z << " " << SegX
1309  << " " << SegY << " " << (moduleThickness / 2) << " " << ionization_points[i].energy()
1310  << " " << hit.particleType() << " " << hit.pabs() << " " << hit.energyLoss() << " "
1311  << hit.entryPoint() << " " << hit.exitPoint() << "\n";
1312 
1313  if (DriftDistance < 0.) {
1314  DriftDistance = 0.;
1315  } else if (DriftDistance > moduleThickness)
1316  DriftDistance = moduleThickness;
1317 
1318  // Assume full depletion now, partial depletion will come later.
1319  float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
1320  float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
1321 
1322  // Shift cloud center
1323  float CloudCenterX = SegX + XDriftDueToMagField;
1324  float CloudCenterY = SegY + YDriftDueToMagField;
1325 
1326  // Calculate how long is the charge drift path
1327  DriftLength = sqrt(DriftDistance * DriftDistance + XDriftDueToMagField * XDriftDueToMagField +
1328  YDriftDueToMagField * YDriftDueToMagField);
1329 
1330  // What is the charge diffusion after this path
1331  Sigma = sqrt(DriftLength / Dist300) * Sigma0;
1332 
1333  // Project the diffusion sigma on the collection plane
1334  Sigma_x = Sigma / CosLorenzAngleX;
1335  Sigma_y = Sigma / CosLorenzAngleY;
1336 
1337  // Insert a charge loss due to Rad Damage here
1338  float energyOnCollector = ionization_points[i].energy(); // The energy that reaches the collector
1339 
1340  // add pixel aging
1341  if (AddPixelAging) {
1342  float kValue = pixel_aging(pixelAging_, pixdet, tTopo);
1343  energyOnCollector *= exp(-1 * kValue * DriftDistance / moduleThickness);
1344  }
1345 
1346 #ifdef TP_DEBUG
1347  LogDebug("Pixel Digitizer") << " Dift DistanceZ= " << DriftDistance << " module thickness= " << moduleThickness
1348  << " Start Energy= " << ionization_points[i].energy()
1349  << " Energy after loss= " << energyOnCollector;
1350 #endif
1351  SignalPoint sp(CloudCenterX, CloudCenterY, Sigma_x, Sigma_y, hit.tof(), energyOnCollector);
1352 
1353  // Load the Charge distribution parameters
1354  collection_points[i] = (sp);
1355 
1356  } // loop over ionization points, i.
1357 
1358 } // end drift
1359 
1360 //*************************************************************************
1361 // Induce the signal on the collection plane of the active sensor area.
1362 void SiPixelDigitizerAlgorithm::induce_signal(std::vector<PSimHit>::const_iterator inputBegin,
1363  std::vector<PSimHit>::const_iterator inputEnd,
1364  const PSimHit& hit,
1365  const size_t hitIndex,
1366  const size_t FirstHitIndex,
1367  const unsigned int tofBin,
1368  const PixelGeomDetUnit* pixdet,
1369  const std::vector<SignalPoint>& collection_points) {
1370  // X - Rows, Left-Right, 160, (1.6cm) for barrel
1371  // Y - Columns, Down-Up, 416, (6.4cm)
1372 
1373  const PixelTopology* topol = &pixdet->specificTopology();
1374  uint32_t detID = pixdet->geographicalId().rawId();
1375  signal_map_type& theSignal = _signal[detID];
1376 
1377 #ifdef TP_DEBUG
1378  LogDebug("Pixel Digitizer") << " enter induce_signal, " << topol->pitch().first << " " << topol->pitch().second; //OK
1379 #endif
1380 
1381  // local map to store pixels hit by 1 Hit.
1382  typedef std::map<int, float, std::less<int> > hit_map_type;
1383  hit_map_type hit_signal;
1384 
1385  // map to store pixel integrals in the x and in the y directions
1386  std::map<int, float, std::less<int> > x, y;
1387 
1388  // Assign signals to readout channels and store sorted by channel number
1389 
1390  // Iterate over collection points on the collection plane
1391  for (std::vector<SignalPoint>::const_iterator i = collection_points.begin(); i != collection_points.end(); ++i) {
1392  float CloudCenterX = i->position().x(); // Charge position in x
1393  float CloudCenterY = i->position().y(); // in y
1394  float SigmaX = i->sigma_x(); // Charge spread in x
1395  float SigmaY = i->sigma_y(); // in y
1396  float Charge = i->amplitude(); // Charge amplitude
1397 
1398  if (SigmaX == 0 || SigmaY == 0) {
1399  LogDebug("Pixel Digitizer") << SigmaX << " " << SigmaY << " cloud " << i->position().x() << " "
1400  << i->position().y() << " " << i->sigma_x() << " " << i->sigma_y() << " "
1401  << i->amplitude() << "\n";
1402  }
1403 
1404 #ifdef TP_DEBUG
1405  LogDebug("Pixel Digitizer") << " cloud " << i->position().x() << " " << i->position().y() << " " << i->sigma_x()
1406  << " " << i->sigma_y() << " " << i->amplitude();
1407 #endif
1408 
1409  // Find the maximum cloud spread in 2D plane , assume 3*sigma
1410  float CloudRight = CloudCenterX + ClusterWidth * SigmaX;
1411  float CloudLeft = CloudCenterX - ClusterWidth * SigmaX;
1412  float CloudUp = CloudCenterY + ClusterWidth * SigmaY;
1413  float CloudDown = CloudCenterY - ClusterWidth * SigmaY;
1414 
1415  // Define 2D cloud limit points
1416  LocalPoint PointRightUp = LocalPoint(CloudRight, CloudUp);
1417  LocalPoint PointLeftDown = LocalPoint(CloudLeft, CloudDown);
1418 
1419  // This points can be located outside the sensor area.
1420  // The conversion to measurement point does not check for that
1421  // so the returned pixel index might be wrong (outside range).
1422  // We rely on the limits check below to fix this.
1423  // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
1424  // PIXEL AREA IS LOST, it should not be collected.
1425 
1426  // Convert the 2D points to pixel indices
1427  MeasurementPoint mp = topol->measurementPosition(PointRightUp); //OK
1428 
1429  int IPixRightUpX = int(floor(mp.x()));
1430  int IPixRightUpY = int(floor(mp.y()));
1431 
1432 #ifdef TP_DEBUG
1433  LogDebug("Pixel Digitizer") << " right-up " << PointRightUp << " " << mp.x() << " " << mp.y() << " " << IPixRightUpX
1434  << " " << IPixRightUpY;
1435 #endif
1436 
1437  mp = topol->measurementPosition(PointLeftDown); //OK
1438 
1439  int IPixLeftDownX = int(floor(mp.x()));
1440  int IPixLeftDownY = int(floor(mp.y()));
1441 
1442 #ifdef TP_DEBUG
1443  emd::LogDebug("Pixel Digitizer") << " left-down " << PointLeftDown << " " << mp.x() << " " << mp.y() << " "
1444  << IPixLeftDownX << " " << IPixLeftDownY;
1445 #endif
1446 
1447  // Check detector limits to correct for pixels outside range.
1448  int numColumns = topol->ncolumns(); // det module number of cols&rows
1449  int numRows = topol->nrows();
1450 
1451  IPixRightUpX = numRows > IPixRightUpX ? IPixRightUpX : numRows - 1;
1452  IPixRightUpY = numColumns > IPixRightUpY ? IPixRightUpY : numColumns - 1;
1453  IPixLeftDownX = 0 < IPixLeftDownX ? IPixLeftDownX : 0;
1454  IPixLeftDownY = 0 < IPixLeftDownY ? IPixLeftDownY : 0;
1455 
1456  x.clear(); // clear temporary integration array
1457  y.clear();
1458 
1459  // First integrate charge strips in x
1460  int ix; // TT for compatibility
1461  for (ix = IPixLeftDownX; ix <= IPixRightUpX; ix++) { // loop over x index
1462  float xUB, xLB, UpperBound, LowerBound;
1463 
1464  // Why is set to 0 if ix=0, does it meen that we accept charge
1465  // outside the sensor? CHeck How it was done in ORCA?
1466  //if(ix == 0) LowerBound = 0.;
1467  if (ix == 0 || SigmaX == 0.) // skip for surface segemnts
1468  LowerBound = 0.;
1469  else {
1470  mp = MeasurementPoint(float(ix), 0.0);
1471  xLB = topol->localPosition(mp).x();
1472  LowerBound = 1 - calcQ((xLB - CloudCenterX) / SigmaX);
1473  }
1474 
1475  if (ix == numRows - 1 || SigmaX == 0.)
1476  UpperBound = 1.;
1477  else {
1478  mp = MeasurementPoint(float(ix + 1), 0.0);
1479  xUB = topol->localPosition(mp).x();
1480  UpperBound = 1. - calcQ((xUB - CloudCenterX) / SigmaX);
1481  }
1482 
1483  float TotalIntegrationRange = UpperBound - LowerBound; // get strip
1484  x[ix] = TotalIntegrationRange; // save strip integral
1485  if (SigmaX == 0 || SigmaY == 0)
1486  LogDebug("Pixel Digitizer") << TotalIntegrationRange << " " << ix << "\n";
1487  }
1488 
1489  // Now integrate strips in y
1490  int iy; // TT for compatibility
1491  for (iy = IPixLeftDownY; iy <= IPixRightUpY; iy++) { //loope over y ind
1492  float yUB, yLB, UpperBound, LowerBound;
1493 
1494  if (iy == 0 || SigmaY == 0.)
1495  LowerBound = 0.;
1496  else {
1497  mp = MeasurementPoint(0.0, float(iy));
1498  yLB = topol->localPosition(mp).y();
1499  LowerBound = 1. - calcQ((yLB - CloudCenterY) / SigmaY);
1500  }
1501 
1502  if (iy == numColumns - 1 || SigmaY == 0.)
1503  UpperBound = 1.;
1504  else {
1505  mp = MeasurementPoint(0.0, float(iy + 1));
1506  yUB = topol->localPosition(mp).y();
1507  UpperBound = 1. - calcQ((yUB - CloudCenterY) / SigmaY);
1508  }
1509 
1510  float TotalIntegrationRange = UpperBound - LowerBound;
1511  y[iy] = TotalIntegrationRange; // save strip integral
1512  if (SigmaX == 0 || SigmaY == 0)
1513  LogDebug("Pixel Digitizer") << TotalIntegrationRange << " " << iy << "\n";
1514  }
1515 
1516  // Get the 2D charge integrals by folding x and y strips
1517  int chan;
1518  for (ix = IPixLeftDownX; ix <= IPixRightUpX; ix++) { // loop over x index
1519  for (iy = IPixLeftDownY; iy <= IPixRightUpY; iy++) { //loope over y ind
1520 
1521  float ChargeFraction = Charge * x[ix] * y[iy];
1522 
1523  if (ChargeFraction > 0.) {
1524  chan = PixelDigi::pixelToChannel(ix, iy); // Get index
1525  // Load the amplitude
1526  hit_signal[chan] += ChargeFraction;
1527  } // endif
1528 
1529 #ifdef TP_DEBUG
1530  mp = MeasurementPoint(float(ix), float(iy));
1531  LocalPoint lp = topol->localPosition(mp);
1532  chan = topol->channel(lp);
1533  LogDebug("Pixel Digitizer") << " pixel " << ix << " " << iy << " - "
1534  << " " << chan << " " << ChargeFraction << " " << mp.x() << " " << mp.y() << " "
1535  << lp.x() << " " << lp.y() << " " // givex edge position
1536  << chan; // edge belongs to previous ?
1537 #endif
1538 
1539  } // endfor iy
1540  } //endfor ix
1541 
1542  } // loop over charge distributions
1543 
1544  // Fill the global map with all hit pixels from this event
1545 
1546  bool reweighted = false;
1547  // size_t ReferenceIndex4CR=0;
1548  if (UseReweighting) {
1549  if (hit.processType() == 0) {
1550  // ReferenceIndex4CR=hitIndex;
1551  reweighted = TheNewSiPixelChargeReweightingAlgorithmClass->hitSignalReweight(
1552  hit, hit_signal, hitIndex, tofBin, topol, detID, theSignal, hit.processType(), makeDigiSimLinks_);
1553  } else {
1554  std::vector<PSimHit>::const_iterator crSimHit = inputBegin;
1555  // ReferenceIndex4CR=FirstHitIndex;
1556  // if the first hit in the same detId is not associated to the same trackId, try to find a better match
1557  if ((*inputBegin).trackId() != hit.trackId()) {
1558  // loop over all the hit from the 1st in the same detId to the hit itself to find the primary particle of the same trackId
1559  uint32_t detId = pixdet->geographicalId().rawId();
1560  size_t localIndex = FirstHitIndex;
1561  for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; localIndex < hitIndex;
1562  ++ssbegin, ++localIndex) {
1563  if ((*ssbegin).detUnitId() != detId) {
1564  continue;
1565  }
1566  if ((*ssbegin).trackId() == hit.trackId() && (*ssbegin).processType() == 0) {
1567  crSimHit = ssbegin;
1568  // ReferenceIndex4CR=localIndex;
1569  break;
1570  }
1571  }
1572  }
1573 
1574  reweighted = TheNewSiPixelChargeReweightingAlgorithmClass->hitSignalReweight(
1575  (*crSimHit), hit_signal, hitIndex, tofBin, topol, detID, theSignal, hit.processType(), makeDigiSimLinks_);
1576  }
1577  }
1578  if (!reweighted) {
1579  for (hit_map_type::const_iterator im = hit_signal.begin(); im != hit_signal.end(); ++im) {
1580  int chan = (*im).first;
1581  theSignal[chan] += (makeDigiSimLinks_ ? Amplitude((*im).second, &hit, hitIndex, tofBin, (*im).second)
1582  : Amplitude((*im).second, (*im).second));
1583 
1584 #ifdef TP_DEBUG
1585  std::pair<int, int> ip = PixelDigi::channelToPixel(chan);
1586  LogDebug("Pixel Digitizer") << " pixel " << ip.first << " " << ip.second << " " << theSignal[chan];
1587 #endif
1588  }
1589  }
1590 
1591 } // end induce_signal
1592 
1593 /***********************************************************************/
1594 
1595 // Build pixels, check threshold, add misscalibration, ...
1596 void SiPixelDigitizerAlgorithm::make_digis(float thePixelThresholdInE,
1597  uint32_t detID,
1598  const PixelGeomDetUnit* pixdet,
1599  std::vector<PixelDigi>& digis,
1600  std::vector<PixelDigiSimLink>& simlinks,
1601  const TrackerTopology* tTopo) const {
1602 #ifdef TP_DEBUG
1603  LogDebug("Pixel Digitizer") << " make digis "
1604  << " "
1605  << " pixel threshold FPix" << theThresholdInE_FPix << " "
1606  << " pixel threshold BPix" << theThresholdInE_BPix << " "
1607  << " pixel threshold BPix Layer1" << theThresholdInE_BPix_L1 << " "
1608  << " pixel threshold BPix Layer2" << theThresholdInE_BPix_L2 << " "
1609  << " List pixels passing threshold ";
1610 #endif
1611 
1612  // Loop over hit pixels
1613 
1614  signalMaps::const_iterator it = _signal.find(detID);
1615  if (it == _signal.end()) {
1616  return;
1617  }
1618 
1619  const signal_map_type& theSignal = (*it).second;
1620 
1621  // unsigned long is enough to store SimTrack id and EncodedEventId
1622  using TrackEventId = std::pair<decltype(SimTrack().trackId()), decltype(EncodedEventId().rawId())>;
1623  std::map<TrackEventId, float> simi; // re-used
1624 
1625  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1626  float signalInElectrons = (*i).second; // signal in electrons
1627 
1628  // Do the miss calibration for calibration studies only.
1629  //if(doMissCalibrate) signalInElectrons = missCalibrate(signalInElectrons)
1630 
1631  // Do only for pixels above threshold
1632 
1633  if (signalInElectrons >= thePixelThresholdInE &&
1634  signalInElectrons > 0.) { // check threshold, always reject killed (0-charge) digis
1635 
1636  int chan = (*i).first; // channel number
1637  std::pair<int, int> ip = PixelDigi::channelToPixel(chan);
1638  int adc = 0; // ADC count as integer
1639 
1640  // Do the miss calibration for calibration studies only.
1641  if (doMissCalibrate) {
1642  int row = ip.first; // X in row
1643  int col = ip.second; // Y is in col
1644  adc = int(missCalibrate(detID, tTopo, pixdet, col, row, signalInElectrons)); //full misscalib.
1645  } else { // Just do a simple electron->adc conversion
1646  adc = int(signalInElectrons / theElectronPerADC); // calibrate gain
1647  }
1648  adc = std::min(adc, theAdcFullScale); // Check maximum value
1649 #ifdef TP_DEBUG
1650  LogDebug("Pixel Digitizer") << (*i).first << " " << (*i).second << " " << signalInElectrons << " " << adc
1651  << ip.first << " " << ip.second;
1652 #endif
1653 
1654  // Load digis
1655  digis.emplace_back(ip.first, ip.second, adc);
1656 
1657  if (makeDigiSimLinks_ && !(*i).second.hitInfos().empty()) {
1658  //digilink
1659  unsigned int il = 0;
1660  for (const auto& info : (*i).second.hitInfos()) {
1661  // note: according to C++ standard operator[] does
1662  // value-initializiation, which for float means initial value of 0
1663  simi[std::make_pair(info.trackId(), info.eventId().rawId())] += (*i).second.individualampl()[il];
1664  il++;
1665  }
1666 
1667  //sum the contribution of the same trackid
1668  for (const auto& info : (*i).second.hitInfos()) {
1669  // skip if track already processed
1670  auto found = simi.find(std::make_pair(info.trackId(), info.eventId().rawId()));
1671  if (found == simi.end())
1672  continue;
1673 
1674  float sum_samechannel = found->second;
1675  float fraction = sum_samechannel / (*i).second;
1676  if (fraction > 1.f)
1677  fraction = 1.f;
1678 
1679  // Approximation: pick hitIndex and tofBin only from the first SimHit
1680  simlinks.emplace_back((*i).first, info.trackId(), info.hitIndex(), info.tofBin(), info.eventId(), fraction);
1681  simi.erase(found);
1682  }
1683  simi.clear(); // although should be empty already
1684  }
1685  }
1686  }
1687 }
1688 
1689 /***********************************************************************/
1690 
1691 // Add electronic noise to pixel charge
1693  float thePixelThreshold,
1694  CLHEP::HepRandomEngine* engine) {
1695 #ifdef TP_DEBUG
1696  LogDebug("Pixel Digitizer") << " enter add_noise " << theNoiseInElectrons;
1697 #endif
1698 
1699  uint32_t detID = pixdet->geographicalId().rawId();
1700  signal_map_type& theSignal = _signal[detID];
1701 
1702  // First add noise to hit pixels
1703  float theSmearedChargeRMS = 0.0;
1704 
1705  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); i++) {
1706  if (addChargeVCALSmearing) {
1707  if ((*i).second < 3000) {
1708  theSmearedChargeRMS = 543.6 - (*i).second * 0.093;
1709  } else if ((*i).second < 6000) {
1710  theSmearedChargeRMS = 307.6 - (*i).second * 0.01;
1711  } else {
1712  theSmearedChargeRMS = -432.4 + (*i).second * 0.123;
1713  }
1714 
1715  // Noise from Vcal smearing:
1716  float noise_ChargeVCALSmearing = theSmearedChargeRMS * CLHEP::RandGaussQ::shoot(engine, 0., 1.);
1717  // Noise from full readout:
1718  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1719 
1720  if (((*i).second + Amplitude(noise + noise_ChargeVCALSmearing, -1.)) < 0.) {
1721  (*i).second.set(0);
1722  } else {
1723  (*i).second += Amplitude(noise + noise_ChargeVCALSmearing, -1.);
1724  }
1725 
1726  } // End if addChargeVCalSmearing
1727  else {
1728  // Noise: ONLY full READOUT Noise.
1729  // Use here the FULL readout noise, including TBM,ALT,AOH,OPT-REC.
1730  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1731 
1732  if (((*i).second + Amplitude(noise, -1.)) < 0.) {
1733  (*i).second.set(0);
1734  } else {
1735  (*i).second += Amplitude(noise, -1.);
1736  }
1737  } // end if only Noise from full readout
1738  }
1739 
1740  if (!addNoisyPixels) // Option to skip noise in non-hit pixels
1741  return;
1742 
1743  const PixelTopology* topol = &pixdet->specificTopology();
1744  int numColumns = topol->ncolumns(); // det module number of cols&rows
1745  int numRows = topol->nrows();
1746 
1747  // Add noise on non-hit pixels
1748  // Use here the pixel noise
1749  int numberOfPixels = (numRows * numColumns);
1750  std::map<int, float, std::less<int> > otherPixels;
1751  std::map<int, float, std::less<int> >::iterator mapI;
1752 
1753  theNoiser->generate(numberOfPixels,
1754  thePixelThreshold, //thr. in un. of nois
1755  theNoiseInElectrons, // noise in elec.
1756  otherPixels,
1757  engine);
1758 
1759 #ifdef TP_DEBUG
1760  LogDebug("Pixel Digitizer") << " Add noisy pixels " << numRows << " " << numColumns << " " << theNoiseInElectrons
1761  << " " << theThresholdInE_FPix << theThresholdInE_BPix << " " << numberOfPixels << " "
1762  << otherPixels.size();
1763 #endif
1764 
1765  // Add noisy pixels
1766  for (mapI = otherPixels.begin(); mapI != otherPixels.end(); mapI++) {
1767  int iy = ((*mapI).first) / numRows;
1768  int ix = ((*mapI).first) - (iy * numRows);
1769 
1770  // Keep for a while for testing.
1771  if (iy < 0 || iy > (numColumns - 1))
1772  LogWarning("Pixel Geometry") << " error in iy " << iy;
1773  if (ix < 0 || ix > (numRows - 1))
1774  LogWarning("Pixel Geometry") << " error in ix " << ix;
1775 
1776  int chan = PixelDigi::pixelToChannel(ix, iy);
1777 
1778 #ifdef TP_DEBUG
1779  LogDebug("Pixel Digitizer") << " Storing noise = " << (*mapI).first << " " << (*mapI).second << " " << ix << " "
1780  << iy << " " << chan;
1781 #endif
1782 
1783  if (theSignal[chan] == 0) {
1784  // float noise = float( (*mapI).second );
1785  int noise = int((*mapI).second);
1786  theSignal[chan] = Amplitude(noise, -1.);
1787  }
1788  }
1789 }
1790 
1791 /***********************************************************************/
1792 
1793 // Simulate the readout inefficiencies.
1794 // Delete a selected number of single pixels, dcols and rocs.
1796  const PixelGeomDetUnit* pixdet,
1797  const TrackerTopology* tTopo,
1798  CLHEP::HepRandomEngine* engine) {
1799  uint32_t detID = pixdet->geographicalId().rawId();
1800  signal_map_type& theSignal = _signal[detID];
1801  const PixelTopology* topol = &pixdet->specificTopology();
1802  int numColumns = topol->ncolumns(); // det module number of cols&rows
1803  int numRows = topol->nrows();
1804  bool isPhase1 = pixdet->subDetector() == GeomDetEnumerators::SubDetector::P1PXB ||
1806  // Predefined efficiencies
1807  double pixelEfficiency = 1.0;
1808  double columnEfficiency = 1.0;
1809  double chipEfficiency = 1.0;
1810  std::vector<double> pixelEfficiencyROCStdPixels(16, 1);
1811  std::vector<double> pixelEfficiencyROCBigPixels(16, 1);
1812 
1813  auto pIndexConverter = PixelIndices(numColumns, numRows);
1814 
1815  std::vector<int> badRocsFromFEDChannels(16, 0);
1816  if (eff.PixelFEDChannelCollection_ != nullptr) {
1818 
1819  if (it != eff.PixelFEDChannelCollection_->end()) {
1820  const std::vector<CablingPathToDetUnit>& path = map_->pathToDetUnit(detID);
1821  for (const auto& ch : *it) {
1822  for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
1823  for (const auto p : path) {
1824  const PixelROC* myroc = map_->findItem(p);
1825  if (myroc->idInDetUnit() == static_cast<unsigned int>(i_roc)) {
1826  LocalPixel::RocRowCol local = {39, 25}; //corresponding to center of ROC row,col
1827  GlobalPixel global = myroc->toGlobal(LocalPixel(local));
1828  int chipIndex(0), colROC(0), rowROC(0);
1829  pIndexConverter.transformToROC(global.col, global.row, chipIndex, colROC, rowROC);
1830  badRocsFromFEDChannels.at(chipIndex) = 1;
1831  }
1832  }
1833  }
1834  } // loop over channels
1835  } // detID in PixelFEDChannelCollection_
1836  } // has PixelFEDChannelCollection_
1837 
1838  if (eff.FromConfig) {
1839  // setup the chip indices conversion
1841  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P1PXB) { // barrel layers
1842  int layerIndex = tTopo->layer(detID);
1843  pixelEfficiency = eff.thePixelEfficiency[layerIndex - 1];
1844  columnEfficiency = eff.thePixelColEfficiency[layerIndex - 1];
1845  chipEfficiency = eff.thePixelChipEfficiency[layerIndex - 1];
1846  LogDebug("Pixel Digitizer") << "Using BPix columnEfficiency = " << columnEfficiency
1847  << " for layer = " << layerIndex << "\n";
1848  // This should never happen, but only check if it is not an upgrade geometry
1849  if (NumberOfBarrelLayers == 3) {
1850  if (numColumns > 416)
1851  LogWarning("Pixel Geometry") << " wrong columns in barrel " << numColumns;
1852  if (numRows > 160)
1853  LogWarning("Pixel Geometry") << " wrong rows in barrel " << numRows;
1854 
1855  int ladder = tTopo->pxbLadder(detID);
1856  int module = tTopo->pxbModule(detID);
1857  if (module <= 4)
1858  module = 5 - module;
1859  else
1860  module -= 4;
1861 
1862  columnEfficiency *= eff.theLadderEfficiency_BPix[layerIndex - 1][ladder - 1] *
1863  eff.theModuleEfficiency_BPix[layerIndex - 1][module - 1] * eff.pu_scale[layerIndex - 1];
1864  }
1867  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2PXEC) { // forward disks
1868 
1869  unsigned int diskIndex =
1870  tTopo->layer(detID) + eff.FPixIndex; // Use diskIndex-1 later to stay consistent with BPix
1871  unsigned int panelIndex = tTopo->pxfPanel(detID);
1872  unsigned int moduleIndex = tTopo->pxfModule(detID);
1873  //if (eff.FPixIndex>diskIndex-1){throw cms::Exception("Configuration") <<"SiPixelDigitizer is using the wrong efficiency value. index = "
1874  // <<diskIndex-1<<" , MinIndex = "<<eff.FPixIndex<<" ... "<<tTopo->pxfDisk(detID);}
1875  pixelEfficiency = eff.thePixelEfficiency[diskIndex - 1];
1876  columnEfficiency = eff.thePixelColEfficiency[diskIndex - 1];
1877  chipEfficiency = eff.thePixelChipEfficiency[diskIndex - 1];
1878  LogDebug("Pixel Digitizer") << "Using FPix columnEfficiency = " << columnEfficiency
1879  << " for Disk = " << tTopo->pxfDisk(detID) << "\n";
1880  // Sometimes the forward pixels have wrong size,
1881  // this crashes the index conversion, so exit, but only check if it is not an upgrade geometry
1882  if (NumberOfBarrelLayers ==
1883  3) { // whether it is the present or the phase 1 detector can be checked using GeomDetEnumerators::SubDetector
1884  if (numColumns > 260 || numRows > 160) {
1885  if (numColumns > 260)
1886  LogWarning("Pixel Geometry") << " wrong columns in endcaps " << numColumns;
1887  if (numRows > 160)
1888  LogWarning("Pixel Geometry") << " wrong rows in endcaps " << numRows;
1889  return;
1890  }
1891  if ((panelIndex == 1 && (moduleIndex == 1 || moduleIndex == 2)) ||
1892  (panelIndex == 2 && moduleIndex == 1)) { //inner modules
1893  columnEfficiency *= eff.theInnerEfficiency_FPix[diskIndex - 1] * eff.pu_scale[3];
1894  } else { //outer modules
1895  columnEfficiency *= eff.theOuterEfficiency_FPix[diskIndex - 1] * eff.pu_scale[4];
1896  }
1897  } // current detector, forward
1898  } else if (pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2OTB ||
1900  // If phase 2 outer tracker, hardcoded values as they have been so far
1901  pixelEfficiency = 0.999;
1902  columnEfficiency = 0.999;
1903  chipEfficiency = 0.999;
1904  } // if barrel/forward
1905  } else { // Load precomputed factors from Database
1906  pixelEfficiency = eff.PixelGeomFactors.at(detID);
1907  columnEfficiency = eff.ColGeomFactors.at(detID) * eff.pu_scale[eff.iPU.at(detID)];
1908  chipEfficiency = eff.ChipGeomFactors.at(detID);
1909  if (isPhase1) {
1910  for (unsigned int i_roc = 0; i_roc < eff.PixelGeomFactorsROCStdPixels.at(detID).size(); ++i_roc) {
1911  pixelEfficiencyROCStdPixels[i_roc] = eff.PixelGeomFactorsROCStdPixels.at(detID).at(i_roc);
1912  pixelEfficiencyROCBigPixels[i_roc] = eff.PixelGeomFactorsROCBigPixels.at(detID).at(i_roc);
1913  }
1914  } // is Phase 1
1915  }
1916 
1917 #ifdef TP_DEBUG
1918  LogDebug("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " " << columnEfficiency << " "
1919  << chipEfficiency;
1920 #endif
1921 
1922  // Initilize the index converter
1923  //PixelIndices indexConverter(numColumns,numRows);
1924 
1925  int chipIndex = 0;
1926  int rowROC = 0;
1927  int colROC = 0;
1928  std::map<int, int, std::less<int> > chips, columns, pixelStd, pixelBig;
1929  std::map<int, int, std::less<int> >::iterator iter;
1930 
1931  // Find out the number of columns and rocs hits
1932  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1933  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1934  int chan = i->first;
1935  std::pair<int, int> ip = PixelDigi::channelToPixel(chan);
1936  int row = ip.first; // X in row
1937  int col = ip.second; // Y is in col
1938  //transform to ROC index coordinates
1939  pIndexConverter.transformToROC(col, row, chipIndex, colROC, rowROC);
1940  int dColInChip = pIndexConverter.DColumn(colROC); // get ROC dcol from ROC col
1941  //dcol in mod
1942  int dColInDet = pIndexConverter.DColumnInModule(dColInChip, chipIndex);
1943 
1944  chips[chipIndex]++;
1945  columns[dColInDet]++;
1946  if (isPhase1) {
1947  if (topol->isItBigPixelInX(row) || topol->isItBigPixelInY(col))
1948  pixelBig[chipIndex]++;
1949  else
1950  pixelStd[chipIndex]++;
1951  }
1952  }
1953 
1954  // Delete some ROC hits.
1955  for (iter = chips.begin(); iter != chips.end(); iter++) {
1956  //float rand = RandFlat::shoot();
1957  float rand = CLHEP::RandFlat::shoot(engine);
1958  if (rand > chipEfficiency)
1959  chips[iter->first] = 0;
1960  }
1961 
1962  // Delete some Dcol hits.
1963  for (iter = columns.begin(); iter != columns.end(); iter++) {
1964  //float rand = RandFlat::shoot();
1965  float rand = CLHEP::RandFlat::shoot(engine);
1966  if (rand > columnEfficiency)
1967  columns[iter->first] = 0;
1968  }
1969 
1970  // Delete some pixel hits based on DCDC issue damage.
1971  if (isPhase1) {
1972  for (iter = pixelStd.begin(); iter != pixelStd.end(); iter++) {
1973  float rand = CLHEP::RandFlat::shoot(engine);
1974  if (rand > pixelEfficiencyROCStdPixels[iter->first])
1975  pixelStd[iter->first] = 0;
1976  }
1977 
1978  for (iter = pixelBig.begin(); iter != pixelBig.end(); iter++) {
1979  float rand = CLHEP::RandFlat::shoot(engine);
1980  if (rand > pixelEfficiencyROCBigPixels[iter->first])
1981  pixelBig[iter->first] = 0;
1982  }
1983  }
1984 
1985  // Now loop again over pixels to kill some of them.
1986  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1987  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1988  // int chan = i->first;
1989  std::pair<int, int> ip = PixelDigi::channelToPixel(i->first); //get pixel pos
1990  int row = ip.first; // X in row
1991  int col = ip.second; // Y is in col
1992  //transform to ROC index coordinates
1993  pIndexConverter.transformToROC(col, row, chipIndex, colROC, rowROC);
1994  int dColInChip = pIndexConverter.DColumn(colROC); //get ROC dcol from ROC col
1995  //dcol in mod
1996  int dColInDet = pIndexConverter.DColumnInModule(dColInChip, chipIndex);
1997 
1998  //float rand = RandFlat::shoot();
1999  float rand = CLHEP::RandFlat::shoot(engine);
2000  if (chips[chipIndex] == 0 || columns[dColInDet] == 0 || rand > pixelEfficiency ||
2001  (pixelStd.count(chipIndex) && pixelStd[chipIndex] == 0) ||
2002  (pixelBig.count(chipIndex) && pixelBig[chipIndex] == 0)) {
2003  // make pixel amplitude =0, pixel will be lost at clusterization
2004  i->second.set(0.); // reset amplitude,
2005  } // end if
2006  if (isPhase1) {
2007  if ((pixelStd.count(chipIndex) && pixelStd[chipIndex] == 0) ||
2008  (pixelBig.count(chipIndex) && pixelBig[chipIndex] == 0) || (badRocsFromFEDChannels.at(chipIndex) == 1)) {
2009  //============================================================
2010  // make pixel amplitude =0, pixel will be lost at clusterization
2011  i->second.set(0.); // reset amplitude,
2012  } // end if
2013  } // is Phase 1
2014  if (KillBadFEDChannels && badRocsFromFEDChannels.at(chipIndex) == 1) {
2015  i->second.set(0.);
2016  }
2017  } // end pixel loop
2018 } // end pixel_indefficiency
2019 
2020 //***************************************************************************************
2021 // Simulate pixel aging with an exponential function
2022 //**************************************************************************************
2023 
2025  const PixelGeomDetUnit* pixdet,
2026  const TrackerTopology* tTopo) const {
2027  uint32_t detID = pixdet->geographicalId().rawId();
2028 
2029  // Predefined damage parameter (no aging)
2030  float pseudoRadDamage = 0.0f;
2031 
2032  // setup the chip indices conversion
2034  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P1PXB) { // barrel layers
2035  int layerIndex = tTopo->layer(detID);
2036 
2037  pseudoRadDamage = aging.thePixelPseudoRadDamage[layerIndex - 1];
2038 
2039  LogDebug("Pixel Digitizer") << "pixel_aging: "
2040  << "\n";
2041  LogDebug("Pixel Digitizer") << "Subid " << pixdet->subDetector() << " layerIndex " << layerIndex << " ladder "
2042  << tTopo->pxbLadder(detID) << " module " << tTopo->pxbModule(detID) << "\n";
2043 
2046  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2PXEC) { // forward disks
2047  unsigned int diskIndex =
2048  tTopo->layer(detID) + aging.FPixIndex; // Use diskIndex-1 later to stay consistent with BPix
2049 
2050  pseudoRadDamage = aging.thePixelPseudoRadDamage[diskIndex - 1];
2051 
2052  LogDebug("Pixel Digitizer") << "pixel_aging: "
2053  << "\n";
2054  LogDebug("Pixel Digitizer") << "Subid " << pixdet->subDetector() << " diskIndex " << diskIndex << "\n";
2055  } else if (pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2OTB ||
2057  // if phase 2 OT hardcoded value as it has always been
2058  pseudoRadDamage = 0.f;
2059  } // if barrel/forward
2060 
2061  LogDebug("Pixel Digitizer") << " pseudoRadDamage " << pseudoRadDamage << "\n";
2062  LogDebug("Pixel Digitizer") << " end pixel_aging "
2063  << "\n";
2064 
2065  return pseudoRadDamage;
2066 #ifdef TP_DEBUG
2067  LogDebug("Pixel Digitizer") << " enter pixel_aging " << pseudoRadDamage;
2068 #endif
2069 }
2070 
2071 //***********************************************************************
2072 
2073 // Fluctuate the gain and offset for the amplitude calibration
2074 // Use gaussian smearing.
2075 //float SiPixelDigitizerAlgorithm::missCalibrate(const float amp) const {
2076 //float gain = RandGaussQ::shoot(1.,theGainSmearing);
2077 //float offset = RandGaussQ::shoot(0.,theOffsetSmearing);
2078 //float newAmp = amp * gain + offset;
2079 // More complex misscalibration
2081  const TrackerTopology* tTopo,
2082  const PixelGeomDetUnit* pixdet,
2083  int col,
2084  int row,
2085  const float signalInElectrons) const {
2086  // Central values
2087  //const float p0=0.00352, p1=0.868, p2=112., p3=113.; // pix(0,0,0)
2088  // const float p0=0.00382, p1=0.886, p2=112.7, p3=113.0; // average roc=0
2089  //const float p0=0.00492, p1=1.998, p2=90.6, p3=134.1; // average roc=6
2090  // Smeared (rms)
2091  //const float s0=0.00020, s1=0.051, s2=5.4, s3=4.4; // average roc=0
2092  //const float s0=0.00015, s1=0.043, s2=3.2, s3=3.1; // col average roc=0
2093 
2094  // Make 2 sets of parameters for Fpix and BPIx:
2095 
2096  float p0 = 0.0f;
2097  float p1 = 0.0f;
2098  float p2 = 0.0f;
2099  float p3 = 0.0f;
2100 
2101  if (pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()) { // barrel layers
2102  p0 = BPix_p0;
2103  p1 = BPix_p1;
2104  p2 = BPix_p2;
2105  p3 = BPix_p3;
2106  } else if (pixdet->type().isTrackerPixel()) { // forward disks
2107  p0 = FPix_p0;
2108  p1 = FPix_p1;
2109  p2 = FPix_p2;
2110  p3 = FPix_p3;
2111  } else {
2112  throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;
2113  }
2114 
2115  float newAmp = 0.f; //Modified signal
2116 
2117  // Convert electrons to VCAL units
2118  float signal = (signalInElectrons - electronsPerVCAL_Offset) / electronsPerVCAL;
2119 
2120  // New gains/offsets are needed for phase1 L1
2121  int layer = 0;
2122  if (DetId(detID).subdetId() == 1)
2123  layer = tTopo->pxbLayer(detID);
2124  if (layer == 1)
2125  signal = (signalInElectrons - electronsPerVCAL_L1_Offset) / electronsPerVCAL_L1;
2126 
2127  // Simulate the analog response with fixed parametrization
2128  newAmp = p3 + p2 * tanh(p0 * signal - p1);
2129 
2130  // Use the pixel-by-pixel calibrations
2131  //transform to ROC index coordinates
2132  //int chipIndex=0, colROC=0, rowROC=0;
2133  //std::unique_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
2134  //pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
2135 
2136  // Use calibration from a file
2137  //int chanROC = PixelIndices::pixelToChannelROC(rowROC,colROC); // use ROC coordinates
2138  //float pp0=0, pp1=0,pp2=0,pp3=0;
2139  //map<int,CalParameters,std::less<int> >::const_iterator it=calmap.find(chanROC);
2140  //CalParameters y = (*it).second;
2141  //pp0 = y.p0;
2142  //pp1 = y.p1;
2143  //pp2 = y.p2;
2144  //pp3 = y.p3;
2145 
2146  //
2147  // Use random smearing
2148  // Randomize the pixel response
2149  //float pp0 = RandGaussQ::shoot(p0,s0);
2150  //float pp1 = RandGaussQ::shoot(p1,s1);
2151  //float pp2 = RandGaussQ::shoot(p2,s2);
2152  //float pp3 = RandGaussQ::shoot(p3,s3);
2153 
2154  //newAmp = pp3 + pp2 * tanh(pp0*signal - pp1); // Final signal
2155 
2156  LogDebug("Pixel Digitizer") << " misscalibrate " << col << " " << row
2157  << " "
2158  // <<chipIndex<<" " <<colROC<<" " <<rowROC<<" "
2159  << signalInElectrons << " " << signal << " " << newAmp << " "
2160  << (signalInElectrons / theElectronPerADC) << "\n";
2161 
2162  return newAmp;
2163 }
2164 //******************************************************************************
2165 
2166 // Set the drift direction accoring to the Bfield in local det-unit frame
2167 // Works for both barrel and forward pixels.
2168 // Replace the sign convention to fit M.Swartz's formulaes.
2169 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
2170 // parameter value
2171 
2173  const GlobalVector& bfield,
2174  const DetId& detId) const {
2175  Frame detFrame(pixdet->surface().position(), pixdet->surface().rotation());
2176  LocalVector Bfield = detFrame.toLocal(bfield);
2177 
2178  float alpha2_FPix;
2179  float alpha2_BPix;
2180  float alpha2;
2181 
2182  //float dir_x = -tanLorentzAnglePerTesla * Bfield.y();
2183  //float dir_y = +tanLorentzAnglePerTesla * Bfield.x();
2184  //float dir_z = -1.; // E field always in z direction, so electrons go to -z
2185  // The dir_z has to be +/- 1. !
2186  // LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
2187 
2188  float dir_x = 0.0f;
2189  float dir_y = 0.0f;
2190  float dir_z = 0.0f;
2191  float scale = 0.0f;
2192 
2193  uint32_t detID = pixdet->geographicalId().rawId();
2194 
2195  // Read Lorentz angle from cfg file:**************************************************************
2196 
2197  if (!use_LorentzAngle_DB_) {
2198  if (alpha2Order) {
2201  } else {
2202  alpha2_FPix = 0.0f;
2203  alpha2_BPix = 0.0f;
2204  }
2205 
2206  if (pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()) { // barrel layers
2207  dir_x = -(tanLorentzAnglePerTesla_BPix * Bfield.y() + alpha2_BPix * Bfield.z() * Bfield.x());
2208  dir_y = +(tanLorentzAnglePerTesla_BPix * Bfield.x() - alpha2_BPix * Bfield.z() * Bfield.y());
2209  dir_z = -(1 + alpha2_BPix * Bfield.z() * Bfield.z());
2210  scale = -dir_z;
2211  } else if (pixdet->type().isTrackerPixel()) { // forward disks
2212  dir_x = -(tanLorentzAnglePerTesla_FPix * Bfield.y() + alpha2_FPix * Bfield.z() * Bfield.x());
2213  dir_y = +(tanLorentzAnglePerTesla_FPix * Bfield.x() - alpha2_FPix * Bfield.z() * Bfield.y());
2214  dir_z = -(1 + alpha2_FPix * Bfield.z() * Bfield.z());
2215  scale = -dir_z;
2216  } else {
2217  throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;
2218  }
2219  } // end: Read LA from cfg file.
2220 
2221  //Read Lorentz angle from DB:********************************************************************
2222  if (use_LorentzAngle_DB_) {
2223  float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
2224  alpha2 = lorentzAngle * lorentzAngle;
2225  dir_x = -(lorentzAngle * Bfield.y() + alpha2 * Bfield.z() * Bfield.x());
2226  dir_y = +(lorentzAngle * Bfield.x() - alpha2 * Bfield.z() * Bfield.y());
2227  dir_z = -(1 + alpha2 * Bfield.z() * Bfield.z());
2228  scale = -dir_z;
2229  } // end: Read LA from DataBase.
2230 
2231  LocalVector theDriftDirection = LocalVector(dir_x / scale, dir_y / scale, dir_z / scale);
2232 
2233 #ifdef TP_DEBUG
2234  LogDebug("Pixel Digitizer") << " The drift direction in local coordinate is " << theDriftDirection;
2235 #endif
2236 
2237  return theDriftDirection;
2238 }
2239 
2240 //****************************************************************************************************
2241 
2243  signal_map_type& theSignal = _signal[detID];
2244 
2245  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
2246  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
2247  // int chan = i->first;
2248  std::pair<int, int> ip = PixelDigi::channelToPixel(i->first); //get pixel pos
2249  int row = ip.first; // X in row
2250  int col = ip.second; // Y is in col
2251  //transform to ROC index coordinates
2252  if (theSiPixelGainCalibrationService_->isDead(detID, col, row)) {
2253  LogDebug("Pixel Digitizer") << "now in isdead check, row " << detID << " " << col << "," << row << "\n";
2254  // make pixel amplitude =0, pixel will be lost at clusterization
2255  i->second.set(0.); // reset amplitude,
2256  } // end if
2257  } // end pixel loop
2258 } // end pixel_indefficiency
2259 
2260 //****************************************************************************************************
2261 
2263  bool isbad = false;
2264 
2265  Parameters::const_iterator itDeadModules = DeadModules.begin();
2266 
2267  int detid = detID;
2268  for (; itDeadModules != DeadModules.end(); ++itDeadModules) {
2269  int Dead_detID = itDeadModules->getParameter<int>("Dead_detID");
2270  if (detid == Dead_detID) {
2271  isbad = true;
2272  break;
2273  }
2274  }
2275 
2276  if (!isbad)
2277  return;
2278 
2279  signal_map_type& theSignal = _signal[detID];
2280 
2281  std::string Module = itDeadModules->getParameter<std::string>("Module");
2282 
2283  if (Module == "whole") {
2284  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
2285  i->second.set(0.); // reset amplitude
2286  }
2287  }
2288 
2289  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
2290  std::pair<int, int> ip = PixelDigi::channelToPixel(i->first); //get pixel pos
2291 
2292  if (Module == "tbmA" && ip.first >= 80 && ip.first <= 159) {
2293  i->second.set(0.);
2294  }
2295 
2296  if (Module == "tbmB" && ip.first <= 79) {
2297  i->second.set(0.);
2298  }
2299  }
2300 }
2301 //****************************************************************************************************
2303  // Not SLHC safe for now
2304 
2305  bool isbad = false;
2306 
2307  std::vector<SiPixelQuality::disabledModuleType> disabledModules = SiPixelBadModule_->getBadComponentList();
2308 
2310 
2311  for (size_t id = 0; id < disabledModules.size(); id++) {
2312  if (detID == disabledModules[id].DetID) {
2313  isbad = true;
2314  badmodule = disabledModules[id];
2315  break;
2316  }
2317  }
2318 
2319  if (!isbad)
2320  return;
2321 
2322  signal_map_type& theSignal = _signal[detID];
2323 
2324  LogDebug("Pixel Digitizer") << "Hit in: " << detID << " errorType " << badmodule.errorType << " BadRocs=" << std::hex
2325  << SiPixelBadModule_->getBadRocs(detID) << std::dec << " "
2326  << "\n";
2327  if (badmodule.errorType == 0) { // this is a whole dead module.
2328 
2329  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
2330  i->second.set(0.); // reset amplitude
2331  }
2332  } else { // all other module types: half-modules and single ROCs.
2333  // Get Bad ROC position:
2334  //follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
2335  std::vector<GlobalPixel> badrocpositions(0);
2336  for (unsigned int j = 0; j < 16; j++) {
2337  if (SiPixelBadModule_->IsRocBad(detID, j) == true) {
2338  std::vector<CablingPathToDetUnit> path = map_->pathToDetUnit(detID);
2339  typedef std::vector<CablingPathToDetUnit>::const_iterator IT;
2340  for (IT it = path.begin(); it != path.end(); ++it) {
2341  const PixelROC* myroc = map_->findItem(*it);
2342  if (myroc->idInDetUnit() == j) {
2343  LocalPixel::RocRowCol local = {39, 25}; //corresponding to center of ROC row, col
2344  GlobalPixel global = myroc->toGlobal(LocalPixel(local));
2345  badrocpositions.push_back(global);
2346  break;
2347  }
2348  }
2349  }
2350  } // end of getBadRocPositions
2351 
2352  for (signal_map_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
2353  std::pair<int, int> ip = PixelDigi::channelToPixel(i->first); //get pixel pos
2354 
2355  for (std::vector<GlobalPixel>::const_iterator it = badrocpositions.begin(); it != badrocpositions.end(); ++it) {
2356  if (it->row >= 80 && ip.first >= 80) {
2357  if ((std::abs(ip.second - it->col) < 26)) {
2358  i->second.set(0.);
2359  } else if (it->row == 120 && ip.second - it->col == 26) {
2360  i->second.set(0.);
2361  } else if (it->row == 119 && it->col - ip.second == 26) {
2362  i->second.set(0.);
2363  }
2364  } else if (it->row < 80 && ip.first < 80) {
2365  if ((std::abs(ip.second - it->col) < 26)) {
2366  i->second.set(0.);
2367  } else if (it->row == 40 && ip.second - it->col == 26) {
2368  i->second.set(0.);
2369  } else if (it->row == 39 && it->col - ip.second == 26) {
2370  i->second.set(0.);
2371  }
2372  }
2373  }
2374  }
2375  }
2376 }
void init(const edm::EventSetup &es)
static std::pair< int, int > channelToPixelROC(const int chan)
Definition: PixelIndices.h:236
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
void pixel_inefficiency_db(uint32_t detID)
signal_map_type::const_iterator signal_map_const_iterator
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DVector LocalVector
Definition: LocalVector.h:12
float tof() const
deprecated name for timeOfFlight()
Definition: PSimHit.h:76
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
const std::map< unsigned int, double > & getPixelGeomFactors() const
bool isBarrel() const
Definition: GeomDetType.cc:9
std::map< unsigned int, probabilityVec > probabilityMap
T y() const
Definition: PV2DBase.h:44
std::map< int, CalParameters, std::less< int > > initCal() const
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
LocalVector DriftDirection(const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
uint16_t *__restrict__ id
virtual int ncolumns() const =0
const TString p2
Definition: fwPaths.cc:13
const std::vector< float > & getMix_TrueInteractions() const
probabilityVec getProbabilities(const unsigned int puBin) const
const std::unique_ptr< SiG4UniversalFluctuation > fluctuate
PixelEfficiencies(const edm::ParameterSet &conf, bool AddPixelInefficiency, int NumberOfBarrelLayers, int NumberOfEndcapDisks)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
unsigned int pxfDisk(const DetId &id) const
const std::vector< int > & getMix_bunchCrossing() const
tuple chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, &quot;NDC&quot;) lumi.SetBorderSize( 0 ) lumi...
const GeomDetType & type() const override
virtual int rowsperroc() const =0
unsigned int pxbLadder(const DetId &id) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T y() const
Definition: PV3DBase.h:60
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual int nrows() const =0
edm::ESGetToken< SiPixelQualityProbabilities, SiPixelStatusScenarioProbabilityRcd > scenarioProbabilityToken_
std::unique_ptr< PixelFEDChannelCollection > PixelFEDChannelCollection_
const Bounds & bounds() const
Definition: Surface.h:87
std::vector< std::vector< double > > thePUEfficiency
std::vector< sipixelobjects::CablingPathToDetUnit > pathToDetUnit(uint32_t rawDetId) const final
unsigned int pxbModule(const DetId &id) const
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
const std::map< unsigned int, std::vector< double > > & getPUFactors() const
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
Log< level::Error, false > LogError
int ii
Definition: cuy.py:589
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
A arg
Definition: Factorize.h:31
scenario
Definition: constants.h:173
identify pixel inside single ROC
Definition: LocalPixel.h:7
const SiPixelDynamicInefficiency * SiPixelDynamicInefficiency_
std::unique_ptr< SiPixelChargeReweightingAlgorithm > TheNewSiPixelChargeReweightingAlgorithmClass
constexpr std::array< uint8_t, layerIndexSize > layer
void make_digis(float thePixelThresholdInE, uint32_t detID, const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, const TrackerTopology *tTopo) const
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
bool isThere(GeomDetEnumerators::SubDetector subdet) const
virtual float thickness() const =0
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void induce_signal(std::vector< PSimHit >::const_iterator inputBegin, std::vector< PSimHit >::const_iterator inputEnd, const PSimHit &hit, const size_t hitIndex, const size_t FirstHitIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const std::vector< SignalPoint > &collection_points)
const std::vector< uint32_t > getDetIdmasks() const
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > SiPixelBadModuleToken_
T mag() const
Definition: PV3DBase.h:64
SiPixelDigitizerAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
const PixelFEDChannelCollectionMap * quality_map
const std::map< unsigned int, double > & getColGeomFactors() const
float pixel_aging(const PixelAging &aging, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo) const
std::map< uint32_t, std::vector< double > > PixelGeomFactorsROCBigPixels
virtual bool isItBigPixelInX(int ixbin) const =0
virtual int colsperroc() const =0
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
static int pixelToChannelROC(const int rowROC, const int colROC)
Definition: PixelIndices.h:233
virtual bool isItBigPixelInY(int iybin) const =0
constexpr float Bfield
Definition: Config.h:88
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:37
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
void digitize(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &path) const final
const TString p1
Definition: fwPaths.cc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const SiPixelLorentzAngle * SiPixelLorentzAngle_
virtual int channel(const LocalPoint &p) const =0
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
const std::vector< disabledModuleType > getBadComponentList() const
std::vector< LinkConnSpec >::const_iterator IT
void init_from_db(const TrackerGeometry *, const SiPixelDynamicInefficiency *)
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleSimRcd > SiPixelLorentzAngleToken_
bool isTrackerPixel() const
Definition: GeomDetType.cc:15
signal_map_type::iterator signal_map_iterator
void init_DynIneffDB(const edm::EventSetup &)
std::unique_ptr< PixelFEDChannelCollection > chooseScenario(PileupMixingContent *puInfo, CLHEP::HepRandomEngine *)
const TrackerGeomDet * idToDet(DetId) const override
unsigned int pxfModule(const DetId &id) const
void setSimAccumulator(const std::map< uint32_t, std::map< int, int > > &signalMap)
edm::ESGetToken< PixelFEDChannelCollectionMap, SiPixelFEDChannelContainerESProducerRcd > PixelFEDChannelCollectionMapToken_
const std::map< int, CalParameters, std::less< int > > calmap
void fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments, float elossVector[], CLHEP::HepRandomEngine *) const
unsigned int pxbLayer(const DetId &id) const
__shared__ int noise
Log< level::Info, false > LogInfo
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
Definition: DetId.h:17
std::map< int, Amplitude, std::less< int > > signal_map_type
const SiPixelQuality * SiPixelBadModule_
const SiPixelFedCablingMap * map_
std::map< uint32_t, std::vector< double > > PixelGeomFactorsROCStdPixels
void primary_ionization(const PSimHit &hit, std::vector< EnergyDepositUnit > &ionization_points, CLHEP::HepRandomEngine *) const
void calculateInstlumiFactor(PileupMixingContent *puInfo)
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const size_t inputBeginGlobalIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
row and collumn in ROC representation
Definition: LocalPixel.h:13
const std::map< unsigned int, double > & getChipGeomFactors() const
static const GlobalPoint notFound(0, 0, 0)
short getBadRocs(const uint32_t &detid) const
unsigned short processType() const
Definition: PSimHit.h:120
const SiPixelQualityProbabilities * scenarioProbability_
unsigned int layer(const DetId &id) const
float missCalibrate(uint32_t detID, const TrackerTopology *tTopo, const PixelGeomDetUnit *pixdet, int col, int row, float amp) const
bool matches(const DetId &, const DetId &, const std::vector< uint32_t > &)
std::vector< edm::ParameterSet > Parameters
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
int particleType() const
Definition: PSimHit.h:89
tuple columns
Definition: mps_check.py:244
tuple filename
Definition: lut2db_cfg.py:20
unsigned int trackId() const
Definition: PSimHit.h:106
static unsigned int const shift
float getLorentzAngle(const uint32_t &) const
edm::ESGetToken< SiPixelDynamicInefficiency, SiPixelDynamicInefficiencyRcd > SiPixelDynamicInefficiencyToken_
list entry
Definition: mps_splice.py:68
const RotationType & rotation() const
const bool addNoise
virtual std::pair< float, float > pitch() const =0
PixelAging(const edm::ParameterSet &conf, bool AddPixelAging, int NumberOfBarrelLayers, int NumberOfEndcapDisks)
void drift(const PSimHit &hit, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const TrackerTopology *tTopo, const std::vector< EnergyDepositUnit > &ionization_points, std::vector< SignalPoint > &collection_points) const
Log< level::Warning, false > LogWarning
int col
Definition: cuy.py:1009
virtual SubDetector subDetector() const
Which subdetector.
Definition: GeomDet.cc:38
T x() const
Definition: PV2DBase.h:43
T x() const
Definition: PV3DBase.h:59
const PositionType & position() const
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
tuple size
Write out results.
unsigned int pxfPanel(const DetId &id) const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
unsigned int detUnitId() const
Definition: PSimHit.h:97
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
tuple module
Definition: callgraph.py:69
uint16_t *__restrict__ uint16_t const *__restrict__ adc
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:55
void add_noise(const PixelGeomDetUnit *pixdet, float thePixelThreshold, CLHEP::HepRandomEngine *)
#define LogDebug(id)
void pixel_inefficiency(const PixelEfficiencies &eff, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
unsigned transform(const HcalDetId &id, unsigned transformCode)