CMS 3D CMS Logo

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 #include <iostream>
43 
45 
51 
52 #include <gsl/gsl_sf_erf.h>
54 #include "CLHEP/Random/RandGaussQ.h"
55 #include "CLHEP/Random/RandFlat.h"
56 
57 //#include "PixelIndices.h"
61 
67 
68 // Accessing dead pixel modules from the DB:
70 
72 
87 
94 
95 // Geometry
100 
102 
103 using namespace edm;
104 using namespace sipixelobjects;
105 
106 #define TP_DEBUG // protect all LogDebug with ifdef. Takes too much CPU
107 
108 
110  if(use_ineff_from_db_){// load gain calibration service fromdb...
111  theSiPixelGainCalibrationService_->setESObjects( es );
112  }
113  if(use_deadmodule_DB_) {
114  es.get<SiPixelQualityRcd>().get(SiPixelBadModule_);
115  }
116  if(use_LorentzAngle_DB_) {
117  // Get Lorentz angle from DB record
118  es.get<SiPixelLorentzAngleSimRcd>().get(SiPixelLorentzAngle_);
119  }
120  //gets the map and geometry from the DB (to kill ROCs)
121  es.get<SiPixelFedCablingMapRcd>().get(map_);
122  es.get<TrackerDigiGeometryRecord>().get(geom_);
123 }
124 
125 //=========================================================================
126 
128 
129  _signal(),
130  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", true)),
131  use_ineff_from_db_(conf.getParameter<bool>("useDB")),
132  use_module_killing_(conf.getParameter<bool>("killModules")), // boolean to kill or not modules
133  use_deadmodule_DB_(conf.getParameter<bool>("DeadModules_DB")), // boolean to access dead modules from DB
134  use_LorentzAngle_DB_(conf.getParameter<bool>("LorentzAngle_DB")), // boolean to access Lorentz angle from DB
135 
136  DeadModules(use_deadmodule_DB_ ? Parameters() : conf.getParameter<Parameters>("DeadModules")), // get dead module from cfg file
137 
138  // Common pixel parameters
139  // These are parameters which are not likely to be changed
140  GeVperElectron(3.61E-09), // 1 electron(3.61eV, 1keV(277e, mod 9/06 d.k.
141  Sigma0(0.00037), // Charge diffusion constant 7->3.7
142  Dist300(0.0300), // normalized to 300micron Silicon
143  alpha2Order(conf.getParameter<bool>("Alpha2Order")), // switch on/off of E.B effect
144  ClusterWidth(3.), // Charge integration spread on the collection plane
145 
146  // get external parameters:
147  // To account for upgrade geometries do not assume the number
148  // of layers or disks.
149  NumberOfBarrelLayers(conf.exists("NumPixelBarrel")?conf.getParameter<int>("NumPixelBarrel"):3),
150  NumberOfEndcapDisks(conf.exists("NumPixelEndcap")?conf.getParameter<int>("NumPixelEndcap"):2),
151 
152  // ADC calibration 1adc count(135e.
153  // Corresponds to 2adc/kev, 270[e/kev]/135[e/adc](2[adc/kev]
154  // Be carefull, this parameter is also used in SiPixelDet.cc to
155  // calculate the noise in adc counts from noise in electrons.
156  // Both defaults should be the same.
157  theElectronPerADC(conf.getParameter<double>("ElectronPerAdc")),
158 
159  // ADC saturation value, 255(8bit adc.
160  //theAdcFullScale(conf.getUntrackedParameter<int>("AdcFullScale",255)),
161  theAdcFullScale(conf.getParameter<int>("AdcFullScale")),
162  theAdcFullScaleStack(conf.exists("AdcFullScaleStack")?conf.getParameter<int>("AdcFullScaleStack"):255),
163 
164  // Noise in electrons:
165  // Pixel cell noise, relevant for generating noisy pixels
166  theNoiseInElectrons(conf.getParameter<double>("NoiseInElectrons")),
167 
168  // Fill readout noise, including all readout chain, relevant for smearing
169  //theReadoutNoise(conf.getUntrackedParameter<double>("ReadoutNoiseInElec",500.)),
170  theReadoutNoise(conf.getParameter<double>("ReadoutNoiseInElec")),
171 
172  // Pixel threshold in units of noise:
173  // thePixelThreshold(conf.getParameter<double>("ThresholdInNoiseUnits")),
174  // Pixel threshold in electron units.
175  theThresholdInE_FPix(conf.getParameter<double>("ThresholdInElectrons_FPix")),
176  theThresholdInE_BPix(conf.getParameter<double>("ThresholdInElectrons_BPix")),
177  theThresholdInE_BPix_L1(conf.exists("ThresholdInElectrons_BPix_L1")?conf.getParameter<double>("ThresholdInElectrons_BPix_L1"):theThresholdInE_BPix),
178  theThresholdInE_BPix_L2(conf.exists("ThresholdInElectrons_BPix_L2")?conf.getParameter<double>("ThresholdInElectrons_BPix_L2"):theThresholdInE_BPix),
179 
180  // Add threshold gaussian smearing:
181  theThresholdSmearing_FPix(conf.getParameter<double>("ThresholdSmearing_FPix")),
182  theThresholdSmearing_BPix(conf.getParameter<double>("ThresholdSmearing_BPix")),
183  theThresholdSmearing_BPix_L1(conf.exists("ThresholdSmearing_BPix_L1")?conf.getParameter<double>("ThresholdSmearing_BPix_L1"):theThresholdSmearing_BPix),
184  theThresholdSmearing_BPix_L2(conf.exists("ThresholdSmearing_BPix_L2")?conf.getParameter<double>("ThresholdSmearing_BPix_L2"):theThresholdSmearing_BPix),
185 
186  // electrons to VCAL conversion needed in misscalibrate()
187  electronsPerVCAL(conf.getParameter<double>("ElectronsPerVcal")),
188  electronsPerVCAL_Offset(conf.getParameter<double>("ElectronsPerVcal_Offset")),
189  electronsPerVCAL_L1(conf.exists("ElectronsPerVcal_L1")?conf.getParameter<double>("ElectronsPerVcal_L1"):electronsPerVCAL),
190  electronsPerVCAL_L1_Offset(conf.exists("ElectronsPerVcal_L1_Offset")?conf.getParameter<double>("ElectronsPerVcal_L1_Offset"):electronsPerVCAL_Offset),
191 
192  //theTofCut 12.5, cut in particle TOD +/- 12.5ns
193  //theTofCut(conf.getUntrackedParameter<double>("TofCut",12.5)),
194  theTofLowerCut(conf.getParameter<double>("TofLowerCut")),
195  theTofUpperCut(conf.getParameter<double>("TofUpperCut")),
196 
197  // Get the Lorentz angle from the cfg file:
198  tanLorentzAnglePerTesla_FPix(use_LorentzAngle_DB_ ? 0.0 : conf.getParameter<double>("TanLorentzAnglePerTesla_FPix")),
199  tanLorentzAnglePerTesla_BPix(use_LorentzAngle_DB_ ? 0.0 : conf.getParameter<double>("TanLorentzAnglePerTesla_BPix")),
200 
201  // signal response new parameterization: split Fpix and BPix
202  FPix_p0(conf.getParameter<double>("FPix_SignalResponse_p0")),
203  FPix_p1(conf.getParameter<double>("FPix_SignalResponse_p1")),
204  FPix_p2(conf.getParameter<double>("FPix_SignalResponse_p2")),
205  FPix_p3(conf.getParameter<double>("FPix_SignalResponse_p3")),
206 
207  BPix_p0(conf.getParameter<double>("BPix_SignalResponse_p0")),
208  BPix_p1(conf.getParameter<double>("BPix_SignalResponse_p1")),
209  BPix_p2(conf.getParameter<double>("BPix_SignalResponse_p2")),
210  BPix_p3(conf.getParameter<double>("BPix_SignalResponse_p3")),
211 
212  // Add noise
213  addNoise(conf.getParameter<bool>("AddNoise")),
214 
215  // Smear the pixel charge with a gaussian which RMS is a function of the
216  // pixel charge (Danek's study)
217  addChargeVCALSmearing(conf.getParameter<bool>("ChargeVCALSmearing")),
218 
219  // Add noisy pixels
220  addNoisyPixels(conf.getParameter<bool>("AddNoisyPixels")),
221 
222  // Fluctuate charge in track subsegments
223  fluctuateCharge(conf.getUntrackedParameter<bool>("FluctuateCharge",true)),
224 
225  // Control the pixel inefficiency
226  AddPixelInefficiency(conf.getParameter<bool>("AddPixelInefficiency")),
227 
228  // Add threshold gaussian smearing:
229  addThresholdSmearing(conf.getParameter<bool>("AddThresholdSmearing")),
230 
231  // Get the constants for the miss-calibration studies
232  doMissCalibrate(conf.getParameter<bool>("MissCalibrate")), // Enable miss-calibration
233  theGainSmearing(conf.getParameter<double>("GainSmearing")), // sigma of the gain smearing
234  theOffsetSmearing(conf.getParameter<double>("OffsetSmearing")), //sigma of the offset smearing
235 
236  // Add pixel radiation damage for upgrade studies
237  AddPixelAging(conf.getParameter<bool>("DoPixelAging")),
238 
239  // delta cutoff in MeV, has to be same as in OSCAR(0.030/cmsim=1.0 MeV
240  //tMax(0.030), // In MeV.
241  //tMax(conf.getUntrackedParameter<double>("deltaProductionCut",0.030)),
242  tMax(conf.getParameter<double>("deltaProductionCut")),
243 
244  fluctuate(fluctuateCharge ? new SiG4UniversalFluctuation() : nullptr),
245  theNoiser(addNoise ? new GaussianTailNoiseGenerator() : nullptr),
246  calmap(doMissCalibrate ? initCal() : std::map<int,CalParameters,std::less<int> >()),
247  theSiPixelGainCalibrationService_(use_ineff_from_db_ ? new SiPixelGainCalibrationOfflineSimService(conf) : nullptr),
248  pixelEfficiencies_(conf, AddPixelInefficiency,NumberOfBarrelLayers,NumberOfEndcapDisks),
249  pixelAging_(conf,AddPixelAging,NumberOfBarrelLayers,NumberOfEndcapDisks)
250 {
251  LogInfo ("PixelDigitizer ") <<"SiPixelDigitizerAlgorithm constructed"
252  <<"Configuration parameters:"
253  << "Threshold/Gain = "
254  << "threshold in electron FPix = "
256  << "threshold in electron BPix = "
258  << "threshold in electron BPix Layer1 = "
260  << "threshold in electron BPix Layer2 = "
262  <<" " << theElectronPerADC << " " << theAdcFullScale
263  << " The delta cut-off is set to " << tMax
264  << " pix-inefficiency "<<AddPixelInefficiency;
265 
266 }
267 
268 std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> >
270 
271  using std::cerr;
272  using std::cout;
273  using std::endl;
274 
275  std::map<int, SiPixelDigitizerAlgorithm::CalParameters, std::less<int> > calmap;
276  // Prepare for the analog amplitude miss-calibration
277  LogDebug ("PixelDigitizer ")
278  << " miss-calibrate the pixel amplitude ";
279 
280  const bool ReadCalParameters = false;
281  if(ReadCalParameters) { // Read the calibration files from file
282  // read the calibration constants from a file (testing only)
283  std::ifstream in_file; // data file pointer
284  char filename[80] = "phCalibrationFit_C0.dat";
285 
286  in_file.open(filename, std::ios::in ); // in C++
287  if(in_file.bad()) {
288  cout << " File not found " << endl;
289  return calmap; // signal error
290  }
291  cout << " file opened : " << filename << endl;
292 
293  char line[500];
294  for (int i = 0; i < 3; i++) {
295  in_file.getline(line, 500,'\n');
296  cout<<line<<endl;
297  }
298 
299  cout << " test map" << endl;
300 
301  float par0,par1,par2,par3;
302  int colid,rowid;
304  // Read MC tracks
305  for(int i=0;i<(52*80);i++) { // loop over tracks
306  in_file >> par0 >> par1 >> par2 >> par3 >> name >> colid >> rowid;
307  if(in_file.bad()) { // check for errors
308  cerr << "Cannot read data file" << endl;
309  return calmap;
310  }
311  if( in_file.eof() != 0 ) {
312  cerr << in_file.eof() << " " << in_file.gcount() << " "
313  << in_file.fail() << " " << in_file.good() << " end of file "
314  << endl;
315  return calmap;
316  }
317 
318  //cout << " line " << i << " " <<par0<<" "<<par1<<" "<<par2<<" "<<par3<<" "
319  // <<colid<<" "<<rowid<<endl;
320 
321  CalParameters onePix;
322  onePix.p0=par0;
323  onePix.p1=par1;
324  onePix.p2=par2;
325  onePix.p3=par3;
326 
327  // Convert ROC pixel index to channel
328  int chan = PixelIndices::pixelToChannelROC(rowid,colid);
329  calmap.insert(std::pair<int,CalParameters>(chan,onePix));
330 
331  // Testing the index conversion, can be skipped
332  std::pair<int,int> p = PixelIndices::channelToPixelROC(chan);
333  if(rowid!=p.first) cout<<" wrong channel row "<<rowid<<" "<<p.first<<endl;
334  if(colid!=p.second) cout<<" wrong channel col "<<colid<<" "<<p.second<<endl;
335 
336  } // pixel loop in a ROC
337 
338  cout << " map size " << calmap.size() <<" max "<<calmap.max_size() << " "
339  <<calmap.empty()<< endl;
340 
341 // cout << " map size " << calmap.size() << endl;
342 // map<int,CalParameters,std::less<int> >::iterator ix,it;
343 // map<int,CalParameters,std::less<int> >::const_iterator ip;
344 // for (ix = calmap.begin(); ix != calmap.end(); ++ix) {
345 // int i = (*ix).first;
346 // std::pair<int,int> p = channelToPixelROC(i);
347 // it = calmap.find(i);
348 // CalParameters y = (*it).second;
349 // CalParameters z = (*ix).second;
350 // cout << i <<" "<<p.first<<" "<<p.second<<" "<<y.p0<<" "<<z.p0<<" "<<calmap[i].p0<<endl;
351 
352 // //int dummy=0;
353 // //cin>>dummy;
354 // }
355 
356  } // end if readparameters
357  return calmap;
358 } // end initCal()
359 
360 //=========================================================================
362  LogDebug ("PixelDigitizer")<<"SiPixelDigitizerAlgorithm deleted";
363 }
364 
365 // Read DynIneff Scale factors from Configuration
367  // pixel inefficiency
368  // Don't use Hard coded values, read inefficiencies in from DB/python config or don't use any
369  int NumberOfTotLayers = NumberOfBarrelLayers + NumberOfEndcapDisks;
370  FPixIndex=NumberOfBarrelLayers;
371  if (AddPixelInefficiency){
372  FromConfig =
373  conf.exists("thePixelColEfficiency_BPix1") && conf.exists("thePixelColEfficiency_BPix2") && conf.exists("thePixelColEfficiency_BPix3") &&
374  conf.exists("thePixelColEfficiency_FPix1") && conf.exists("thePixelColEfficiency_FPix2") &&
375  conf.exists("thePixelEfficiency_BPix1") && conf.exists("thePixelEfficiency_BPix2") && conf.exists("thePixelEfficiency_BPix3") &&
376  conf.exists("thePixelEfficiency_FPix1") && conf.exists("thePixelEfficiency_FPix2") &&
377  conf.exists("thePixelChipEfficiency_BPix1") && conf.exists("thePixelChipEfficiency_BPix2") && conf.exists("thePixelChipEfficiency_BPix3") &&
378  conf.exists("thePixelChipEfficiency_FPix1") && conf.exists("thePixelChipEfficiency_FPix2");
379  if (NumberOfBarrelLayers==3) FromConfig = FromConfig && conf.exists("theLadderEfficiency_BPix1") && conf.exists("theLadderEfficiency_BPix2") && conf.exists("theLadderEfficiency_BPix3") &&
380  conf.exists("theModuleEfficiency_BPix1") && conf.exists("theModuleEfficiency_BPix2") && conf.exists("theModuleEfficiency_BPix3") &&
381  conf.exists("thePUEfficiency_BPix1") && conf.exists("thePUEfficiency_BPix2") && conf.exists("thePUEfficiency_BPix3") &&
382  conf.exists("theInnerEfficiency_FPix1") && conf.exists("theInnerEfficiency_FPix2") &&
383  conf.exists("theOuterEfficiency_FPix1") && conf.exists("theOuterEfficiency_FPix2") &&
384  conf.exists("thePUEfficiency_FPix_Inner") && conf.exists("thePUEfficiency_FPix_Outer") &&
385  conf.exists("theInstLumiScaleFactor");
386  if (NumberOfBarrelLayers>=4) FromConfig = FromConfig && conf.exists("thePixelColEfficiency_BPix4") &&
387  conf.exists("thePixelEfficiency_BPix4") && conf.exists("thePixelChipEfficiency_BPix4");
388  if (NumberOfEndcapDisks>=3) FromConfig = FromConfig && conf.exists("thePixelColEfficiency_FPix4") &&
389  conf.exists("thePixelEfficiency_FPix3") && conf.exists("thePixelChipEfficiency_FPix3");
390  if (FromConfig) {
391  LogInfo ("PixelDigitizer ") <<"The PixelDigitizer inefficiency configuration is read from the config file.\n";
392  theInstLumiScaleFactor = conf.getParameter<double>("theInstLumiScaleFactor");
393  int i=0;
394  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix1");
395  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix2");
396  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix3");
397  if (NumberOfBarrelLayers>=4){thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_BPix4");}
398  //
399  i=0;
400  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix1");
401  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix2");
402  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix3");
403  if (NumberOfBarrelLayers>=4){thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_BPix4");}
404  //
405  i=0;
406  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix1");
407  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix2");
408  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix3");
409  if (NumberOfBarrelLayers>=4){thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_BPix4");}
410  //
411  if (NumberOfBarrelLayers==3){
412  i=0;
413  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix1");
414  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix2");
415  theLadderEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theLadderEfficiency_BPix3");
416  if ( ((theLadderEfficiency_BPix[0].size()!=20) || (theLadderEfficiency_BPix[1].size()!=32) ||
417  (theLadderEfficiency_BPix[2].size()!=44)) && (NumberOfBarrelLayers==3) )
418  throw cms::Exception("Configuration") << "Wrong ladder number in efficiency config!";
419  //
420  i=0;
421  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix1");
422  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix2");
423  theModuleEfficiency_BPix[i++] = conf.getParameter<std::vector<double> >("theModuleEfficiency_BPix3");
424  if ( ((theModuleEfficiency_BPix[0].size()!=4) || (theModuleEfficiency_BPix[1].size()!=4) ||
425  (theModuleEfficiency_BPix[2].size()!=4)) && (NumberOfBarrelLayers==3) )
426  throw cms::Exception("Configuration") << "Wrong module number in efficiency config!";
427  //
428  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix1"));
429  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix2"));
430  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_BPix3"));
431  if ( ((thePUEfficiency[0].empty()) || (thePUEfficiency[1].empty()) ||
432  (thePUEfficiency[2].empty())) && (NumberOfBarrelLayers==3) )
433  throw cms::Exception("Configuration") << "At least one PU efficiency (BPix) number is needed in efficiency config!";
434  }
435  // The next is needed for Phase2 Tracker studies
436  if (NumberOfBarrelLayers>=5){
437  if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
438  // For Phase2 tracker layers just set the outermost BPix inefficiency to 99.9% THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
439  for (int j=5 ; j<=NumberOfBarrelLayers ; j++){
440  thePixelColEfficiency[j-1]=0.999;
441  thePixelEfficiency[j-1]=0.999;
442  thePixelChipEfficiency[j-1]=0.999;
443  }
444  }
445  //
446  i=FPixIndex;
447  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix1");
448  thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix2");
449  if (NumberOfEndcapDisks>=3){thePixelColEfficiency[i++] = conf.getParameter<double>("thePixelColEfficiency_FPix3");}
450  i=FPixIndex;
451  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix1");
452  thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix2");
453  if (NumberOfEndcapDisks>=3){thePixelEfficiency[i++] = conf.getParameter<double>("thePixelEfficiency_FPix3");}
454  i=FPixIndex;
455  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix1");
456  thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix2");
457  if (NumberOfEndcapDisks>=3){thePixelChipEfficiency[i++] = conf.getParameter<double>("thePixelChipEfficiency_FPix3");}
458  // The next is needed for Phase2 Tracker studies
459  if (NumberOfEndcapDisks>=4){
460  if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
461  // For Phase2 tracker layers just set the extra FPix disk inefficiency to 99.9% THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
462  for (int j=4+FPixIndex ; j<=NumberOfEndcapDisks+NumberOfBarrelLayers ; j++){
463  thePixelColEfficiency[j-1]=0.999;
464  thePixelEfficiency[j-1]=0.999;
465  thePixelChipEfficiency[j-1]=0.999;
466  }
467  }
468  //FPix Dynamic Inefficiency
469  if (NumberOfBarrelLayers==3){
470  i=FPixIndex;
471  theInnerEfficiency_FPix[i++] = conf.getParameter<double>("theInnerEfficiency_FPix1");
472  theInnerEfficiency_FPix[i++] = conf.getParameter<double>("theInnerEfficiency_FPix2");
473  i=FPixIndex;
474  theOuterEfficiency_FPix[i++] = conf.getParameter<double>("theOuterEfficiency_FPix1");
475  theOuterEfficiency_FPix[i++] = conf.getParameter<double>("theOuterEfficiency_FPix2");
476  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_FPix_Inner"));
477  thePUEfficiency.push_back(conf.getParameter<std::vector<double> >("thePUEfficiency_FPix_Outer"));
478  if ( ((thePUEfficiency[3].empty()) || (thePUEfficiency[4].empty())) && (NumberOfEndcapDisks==2) )
479  throw cms::Exception("Configuration") << "At least one (FPix) PU efficiency number is needed in efficiency config!";
480  pu_scale.resize(thePUEfficiency.size());
481  }
482  }
483  else LogInfo ("PixelDigitizer ") <<"The PixelDigitizer inefficiency configuration is read from the database.\n";
484  }
485  // the first "NumberOfBarrelLayers" settings [0],[1], ... , [NumberOfBarrelLayers-1] are for the barrel pixels
486  // the next "NumberOfEndcapDisks" settings [NumberOfBarrelLayers],[NumberOfBarrelLayers+1], ... [NumberOfEndcapDisks+NumberOfBarrelLayers-1]
487 }
488 
489 // Read DynIneff Scale factors from DB
492  if (bunchspace == 50) es.get<SiPixelDynamicInefficiencyRcd>().get("50ns",SiPixelDynamicInefficiency_);
495  }
496 }
497 
499 
500  theInstLumiScaleFactor = SiPixelDynamicInefficiency->gettheInstLumiScaleFactor();
501  const std::map<uint32_t, double>& PixelGeomFactorsDB = SiPixelDynamicInefficiency->getPixelGeomFactors();
502  const std::map<uint32_t, double>& ColGeomFactorsDB = SiPixelDynamicInefficiency->getColGeomFactors();
503  const std::map<uint32_t, double>& ChipGeomFactorsDB = SiPixelDynamicInefficiency->getChipGeomFactors();
504  const std::map<uint32_t, std::vector<double> >& PUFactors = SiPixelDynamicInefficiency->getPUFactors();
505  std::vector<uint32_t > DetIdmasks = SiPixelDynamicInefficiency->getDetIdmasks();
506 
507  // Loop on all modules, calculate geometrical scale factors and store in map for easy access
508  for( const auto& it_module : geom->detUnits()) {
509  if( dynamic_cast<PixelGeomDetUnit const*>(it_module)==nullptr) continue;
510  const DetId detid = it_module->geographicalId();
511  uint32_t rawid = detid.rawId();
512  PixelGeomFactors[rawid] = 1;
513  ColGeomFactors[rawid] = 1;
514  ChipGeomFactors[rawid] = 1;
515  for (auto db_factor : PixelGeomFactorsDB) if (matches(detid, DetId(db_factor.first), DetIdmasks)) PixelGeomFactors[rawid] *= db_factor.second;
516  for (auto db_factor : ColGeomFactorsDB) if (matches(detid, DetId(db_factor.first), DetIdmasks)) ColGeomFactors[rawid] *= db_factor.second;
517  for (auto db_factor : ChipGeomFactorsDB) if (matches(detid, DetId(db_factor.first), DetIdmasks)) ChipGeomFactors[rawid] *= db_factor.second;
518  }
519 
520  // piluep scale factors are calculated once per event
521  // therefore vector index is stored in a map for each module that matches to a db_id
522  size_t i=0;
523  for (auto factor : PUFactors) {
524  const DetId db_id = DetId(factor.first);
525  for( const auto& it_module : geom->detUnits()) {
526  if( dynamic_cast<PixelGeomDetUnit const*>(it_module)==nullptr) continue;
527  const DetId detid = it_module->geographicalId();
528  if (!matches(detid, db_id, DetIdmasks)) continue;
529  if (iPU.count(detid.rawId())) {
530  throw cms::Exception("Database")<<"Multiple db_ids match to same module in SiPixelDynamicInefficiency DB Object";
531  } else {
532  iPU[detid.rawId()] = i;
533  }
534  }
535  thePUEfficiency.push_back(factor.second);
536  ++i;
537  }
538  pu_scale.resize(thePUEfficiency.size());
539 }
540 
541 bool SiPixelDigitizerAlgorithm::PixelEfficiencies::matches(const DetId& detid, const DetId& db_id, const std::vector<uint32_t >& DetIdmasks) {
542  if (detid.subdetId() != db_id.subdetId()) return false;
543  for (size_t i=0; i<DetIdmasks.size(); ++i) {
544  DetId maskid = DetId(DetIdmasks.at(i));
545  if (maskid.subdetId() != db_id.subdetId()) continue;
546  if ((detid.rawId()&maskid.rawId()) != (db_id.rawId()&maskid.rawId()) &&
547  (db_id.rawId()&maskid.rawId()) != DetId(db_id.det(), db_id.subdetId()).rawId()) return false;
548  }
549  return true;
550 }
551 
553  // pixel aging
554  // Don't use Hard coded values, read aging in from python or don't use any
555  if(AddAging) {
556  int NumberOfTotLayers = NumberOfBarrelLayers + NumberOfEndcapDisks;
557  FPixIndex=NumberOfBarrelLayers;
558 
559  int i=0;
560  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix1");
561  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix2");
562  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix3");
563  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_BPix4");
564 
565  // to be removed when Gaelle will have the phase2 digitizer
566  if (NumberOfBarrelLayers>=5){
567  if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
568  // For Phase2 tracker layers just set the outermost BPix aging 0.
569  for (int j=5 ; j<=NumberOfBarrelLayers ; j++){
570  thePixelPseudoRadDamage[j-1]=0.;
571  }
572  }
573  //
574  i=FPixIndex;
575  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix1");
576  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix2");
577  thePixelPseudoRadDamage[i++] = conf.getParameter<double>("thePixelPseudoRadDamage_FPix3");
578 
579  //To be removed when Phase2 digitizer will be available
580  if (NumberOfEndcapDisks>=4){
581  if (NumberOfTotLayers>20){throw cms::Exception("Configuration") <<"SiPixelDigitizer was given more layers than it can handle";}
582  // For Phase2 tracker layers just set the extra FPix disk aging to 0. BE CAREFUL THESE VALUES ARE HARDCODED ALSO ELSEWHERE IN THIS FILE
583  for (int j=4+FPixIndex ; j<=NumberOfEndcapDisks+NumberOfBarrelLayers ; j++){
584  thePixelPseudoRadDamage[j-1]=0.;
585  }
586  }
587  }
588  // the first "NumberOfBarrelLayers" settings [0],[1], ... , [NumberOfBarrelLayers-1] are for the barrel pixels
589  // the next "NumberOfEndcapDisks" settings [NumberOfBarrelLayers],[NumberOfBarrelLayers+1], ... [NumberOfEndcapDisks+NumberOfBarrelLayers-1]
590 }
591 
592 //=========================================================================
593 void SiPixelDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
594  std::vector<PSimHit>::const_iterator inputEnd,
595  const size_t inputBeginGlobalIndex,
596  const unsigned int tofBin,
597  const PixelGeomDetUnit* pixdet,
598  const GlobalVector& bfield,
599  const TrackerTopology *tTopo,
600  CLHEP::HepRandomEngine* engine) {
601  // produce SignalPoint's for all SimHit's in detector
602  // Loop over hits
603 
604  uint32_t detId = pixdet->geographicalId().rawId();
605  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
606  for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; ssbegin != inputEnd; ++ssbegin, ++simHitGlobalIndex) {
607  // skip hits not in this detector.
608  if((*ssbegin).detUnitId() != detId) {
609  continue;
610  }
611 
612 #ifdef TP_DEBUG
613  LogDebug ("Pixel Digitizer")
614  << (*ssbegin).particleType() << " " << (*ssbegin).pabs() << " "
615  << (*ssbegin).energyLoss() << " " << (*ssbegin).tof() << " "
616  << (*ssbegin).trackId() << " " << (*ssbegin).processType() << " "
617  << (*ssbegin).detUnitId()
618  << (*ssbegin).entryPoint() << " " << (*ssbegin).exitPoint() ;
619 #endif
620 
621 
622  std::vector<EnergyDepositUnit> ionization_points;
623  std::vector<SignalPoint> collection_points;
624 
625  // fill collection_points for this SimHit, indpendent of topology
626  // Check the TOF cut
627  if ( ((*ssbegin).tof() - pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.)>= theTofLowerCut &&
628  ((*ssbegin).tof()- pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.) <= theTofUpperCut ) {
629  primary_ionization(*ssbegin, ionization_points, engine); // fills _ionization_points
630  drift(*ssbegin, pixdet, bfield, tTopo, ionization_points, collection_points); // transforms _ionization_points to collection_points
631  // compute induced signal on readout elements and add to _signal
632  induce_signal(*ssbegin, simHitGlobalIndex, tofBin, pixdet, collection_points); // 1st 3 args needed only for SimHit<-->Digi link
633  } // end if
634  } // end for
635 
636 }
637 
638 //============================================================================
640  //Instlumi scalefactor calculating for dynamic inefficiency
641 
642  if (puInfo) {
643  const std::vector<int>& bunchCrossing = puInfo->getMix_bunchCrossing();
644  const std::vector<float>& TrueInteractionList = puInfo->getMix_TrueInteractions();
645  //const int bunchSpacing = puInfo->getMix_bunchSpacing();
646 
647  int pui = 0, p = 0;
648  std::vector<int>::const_iterator pu;
649  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
650 
651  for (pu=bunchCrossing.begin(); pu!=bunchCrossing.end(); ++pu) {
652  if (*pu==0) {
653  pu0 = pu;
654  p = pui;
655  }
656  pui++;
657  }
658  if (pu0!=bunchCrossing.end()) {
659  for (size_t i=0, n = pixelEfficiencies_.thePUEfficiency.size(); i<n; i++) {
660  double instlumi = TrueInteractionList.at(p)*pixelEfficiencies_.theInstLumiScaleFactor;
661  double instlumi_pow=1.;
663  for (size_t j=0; j<pixelEfficiencies_.thePUEfficiency[i].size(); j++){
665  instlumi_pow*=instlumi;
666  }
667  }
668  }
669  }
670  else {
671  for (int i=0, n = pixelEfficiencies_.thePUEfficiency.size(); i<n; i++) {
673  }
674  }
675 }
676 
677 //============================================================================
679  std::vector<PixelDigi>& digis,
680  std::vector<PixelDigiSimLink>& simlinks,
681  const TrackerTopology *tTopo,
682  CLHEP::HepRandomEngine* engine) {
683 
684  // Pixel Efficiency moved from the constructor to this method because
685  // the information of the det are not available in the constructor
686  // Efficiency parameters. 0 - no inefficiency, 1-low lumi, 10-high lumi
687 
688  uint32_t detID = pixdet->geographicalId().rawId();
689  const signal_map_type& theSignal = _signal[detID];
690 
691  const PixelTopology* topol=&pixdet->specificTopology();
692  int numColumns = topol->ncolumns(); // det module number of cols&rows
693  int numRows = topol->nrows();
694 
695  // Noise already defined in electrons
696  // thePixelThresholdInE = thePixelThreshold * theNoiseInElectrons ;
697  // Find the threshold in noise units, needed for the noiser.
698 
699 
700  float thePixelThresholdInE = 0.;
701 
702  if(theNoiseInElectrons>0.){
703  if(pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()){ // Barrel modules
704  int lay = tTopo->layer(detID);
707  if (lay==1) {
708  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix_L1, theThresholdSmearing_BPix_L1); // gaussian smearing
709  } else if (lay==2) {
710  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix_L2, theThresholdSmearing_BPix_L2); // gaussian smearing
711  } else {
712  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix , theThresholdSmearing_BPix); // gaussian smearing
713  }
714  }
715  } else {
717  if (lay==1) {
718  thePixelThresholdInE = theThresholdInE_BPix_L1;
719  } else if (lay==2) {
720  thePixelThresholdInE = theThresholdInE_BPix_L2;
721  } else {
722  thePixelThresholdInE = theThresholdInE_BPix; // no smearing
723  }
724  }
725  }
726  } else if(pixdet->type().isTrackerPixel()) { // Forward disks modules
728  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_FPix, theThresholdSmearing_FPix); // gaussian smearing
729  } else {
730  thePixelThresholdInE = theThresholdInE_FPix; // no smearing
731  }
732  }
733  else {throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;}
734  }
735 
736 
737 #ifdef TP_DEBUG
738  // full detector thickness
739  float moduleThickness = pixdet->specificSurface().bounds().thickness();
740  LogDebug ("PixelDigitizer")
741  << " PixelDigitizer "
742  << numColumns << " " << numRows << " " << moduleThickness;
743 #endif
744 
745  if(addNoise) add_noise(pixdet, thePixelThresholdInE/theNoiseInElectrons, engine); // generate noise
746 
747  // Do only if needed
748 
749  if((AddPixelInefficiency) && (!theSignal.empty()))
750  pixel_inefficiency(pixelEfficiencies_, pixdet, tTopo, engine); // Kill some pixels
751 
752  if(use_ineff_from_db_ && (!theSignal.empty()))
753  pixel_inefficiency_db(detID);
754 
755  if(use_module_killing_) {
756  if (use_deadmodule_DB_) { // remove dead modules using DB
757  module_killing_DB(detID);
758  } else { // remove dead modules using the list in cfg file
759  module_killing_conf(detID);
760  }
761  }
762 
763  make_digis(thePixelThresholdInE, detID, pixdet, digis, simlinks, tTopo);
764 
765 #ifdef TP_DEBUG
766  LogDebug ("PixelDigitizer") << "[SiPixelDigitizerAlgorithm] converted " << digis.size() << " PixelDigis in DetUnit" << detID;
767 #endif
768 }
769 
770 //***********************************************************************/
771 // Generate primary ionization along the track segment.
772 // Divide the track into small sub-segments
773 void SiPixelDigitizerAlgorithm::primary_ionization(const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points, CLHEP::HepRandomEngine* engine) const {
774 
775  // Straight line approximation for trajectory inside active media
776 
777  const float SegmentLength = 0.0010; //10microns in cm
778  float energy;
779 
780  // Get the 3D segment direction vector
781  LocalVector direction = hit.exitPoint() - hit.entryPoint();
782 
783  float eLoss = hit.energyLoss(); // Eloss in GeV
784  float length = direction.mag(); // Track length in Silicon
785 
786  int NumberOfSegments = int ( length / SegmentLength); // Number of segments
787  if(NumberOfSegments < 1) NumberOfSegments = 1;
788 
789 #ifdef TP_DEBUG
790  LogDebug ("Pixel Digitizer")
791  << " enter primary_ionzation " << NumberOfSegments
792  << " shift = "
793  << (hit.exitPoint().x()-hit.entryPoint().x()) << " "
794  << (hit.exitPoint().y()-hit.entryPoint().y()) << " "
795  << (hit.exitPoint().z()-hit.entryPoint().z()) << " "
796  << hit.particleType() <<" "<< hit.pabs() ;
797 #endif
798 
799  float* elossVector = new float[NumberOfSegments]; // Eloss vector
800 
801  if( fluctuateCharge ) {
802  //MP DA RIMUOVERE ASSOLUTAMENTE
803  int pid = hit.particleType();
804  //int pid=211; // assume it is a pion
805 
806  float momentum = hit.pabs();
807  // Generate fluctuated charge points
808  fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments,
809  elossVector, engine);
810  }
811 
812  ionization_points.resize( NumberOfSegments); // set size
813 
814  // loop over segments
815  for ( int i = 0; i != NumberOfSegments; i++) {
816  // Divide the segment into equal length subsegments
817  Local3DPoint point = hit.entryPoint() +
818  float((i+0.5)/NumberOfSegments) * direction;
819 
820  if( fluctuateCharge )
821  energy = elossVector[i]/GeVperElectron; // Convert charge to elec.
822  else
823  energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
824 
825  EnergyDepositUnit edu( energy, point); //define position,energy point
826  ionization_points[i] = edu; // save
827 
828 #ifdef TP_DEBUG
829  LogDebug ("Pixel Digitizer")
830  << i << " " << ionization_points[i].x() << " "
831  << ionization_points[i].y() << " "
832  << ionization_points[i].z() << " "
833  << ionization_points[i].energy();
834 #endif
835 
836  } // end for loop
837 
838  delete[] elossVector;
839 
840 }
841 //******************************************************************************
842 
843 // Fluctuate the charge comming from a small (10um) track segment.
844 // Use the G4 routine. For mip pions for the moment.
845 void SiPixelDigitizerAlgorithm::fluctuateEloss(int pid, float particleMomentum,
846  float eloss, float length,
847  int NumberOfSegs,float elossVector[],
848  CLHEP::HepRandomEngine* engine) const {
849 
850  // Get dedx for this track
851  //float dedx;
852  //if( length > 0.) dedx = eloss/length;
853  //else dedx = eloss;
854 
855  double particleMass = 139.6; // Mass in MeV, Assume pion
856  pid = std::abs(pid);
857  if(pid!=211) { // Mass in MeV
858  if(pid==11) particleMass = 0.511;
859  else if(pid==13) particleMass = 105.7;
860  else if(pid==321) particleMass = 493.7;
861  else if(pid==2212) particleMass = 938.3;
862  }
863  // What is the track segment length.
864  float segmentLength = length/NumberOfSegs;
865 
866  // Generate charge fluctuations.
867  float de=0.;
868  float sum=0.;
869  double segmentEloss = (1000.*eloss)/NumberOfSegs; //eloss in MeV
870  for (int i=0;i<NumberOfSegs;i++) {
871  // material,*, momentum,energy,*, *, mass
872  //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
873  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
874  // track segment length in mm, segment eloss in MeV
875  // Returns fluctuated eloss in MeV
876  double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
877  de = fluctuate->SampleFluctuations(double(particleMomentum*1000.),
878  particleMass, deltaCutoff,
879  double(segmentLength*10.),
880  segmentEloss, engine )/1000.; //convert to GeV
881  elossVector[i]=de;
882  sum +=de;
883  }
884 
885  if(sum>0.) { // If fluctuations give eloss>0.
886  // Rescale to the same total eloss
887  float ratio = eloss/sum;
888 
889  for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= ratio*elossVector[ii];
890  } else { // If fluctuations gives 0 eloss
891  float averageEloss = eloss/NumberOfSegs;
892  for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= averageEloss;
893  }
894  return;
895 }
896 
897 //*******************************************************************************
898 // Drift the charge segments to the sensor surface (collection plane)
899 // Include the effect of E-field and B-field
901  const PixelGeomDetUnit* pixdet,
902  const GlobalVector& bfield,
903  const TrackerTopology *tTopo,
904  const std::vector<EnergyDepositUnit>& ionization_points,
905  std::vector<SignalPoint>& collection_points) const {
906 
907 #ifdef TP_DEBUG
908  LogDebug ("Pixel Digitizer") << " enter drift " ;
909 #endif
910 
911  collection_points.resize(ionization_points.size()); // set size
912 
913  LocalVector driftDir=DriftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
914  if(driftDir.z() ==0.) {
915  LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
916  return;
917  }
918 
919  // tangent of Lorentz angle
920  //float TanLorenzAngleX = driftDir.x()/driftDir.z();
921  //float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
922 
923  float TanLorenzAngleX, TanLorenzAngleY,dir_z, CosLorenzAngleX,
924  CosLorenzAngleY;
925  if( alpha2Order) {
926  TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
927  TanLorenzAngleY = driftDir.y();
928  dir_z = driftDir.z(); // The z drift direction
929  CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine
930  CosLorenzAngleY = 1./sqrt(1.+TanLorenzAngleY*TanLorenzAngleY); //cosine;
931 
932  } else{
933  TanLorenzAngleX = driftDir.x();
934  TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
935  dir_z = driftDir.z(); // The z drift direction
936  CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine to estimate the path length
937  CosLorenzAngleY = 1.;
938  }
939 
940  float moduleThickness = pixdet->specificSurface().bounds().thickness();
941 #ifdef TP_DEBUG
942  LogDebug ("Pixel Digitizer")
943  << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY <<" "
944  << CosLorenzAngleX << " " << CosLorenzAngleY << " "
945  << moduleThickness*TanLorenzAngleX << " " << driftDir;
946 #endif
947 
948  float Sigma_x = 1.; // Charge spread
949  float Sigma_y = 1.;
950  float DriftDistance; // Distance between charge generation and collection
951  float DriftLength; // Actual Drift Lentgh
952  float Sigma;
953 
954  for (unsigned int i = 0; i != ionization_points.size(); i++) {
955 
956  float SegX, SegY, SegZ; // position
957  SegX = ionization_points[i].x();
958  SegY = ionization_points[i].y();
959  SegZ = ionization_points[i].z();
960 
961  // Distance from the collection plane
962  //DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
963  // Include explixitely the E drift direction (for CMS dir_z=-1)
964  DriftDistance = moduleThickness/2. - (dir_z * SegZ); // Drift to -z
965 
966  //if( DriftDistance <= 0.)
967  //cout<<" <=0 "<<DriftDistance<<" "<<i<<" "<<SegZ<<" "<<dir_z<<" "
968  // <<SegX<<" "<<SegY<<" "<<(moduleThickness/2)<<" "
969  // <<ionization_points[i].energy()<<" "
970  // <<hit.particleType()<<" "<<hit.pabs()<<" "<<hit.energyLoss()<<" "
971  // <<hit.entryPoint()<<" "<<hit.exitPoint()
972  // <<std::endl;
973 
974  if( DriftDistance < 0.) {
975  DriftDistance = 0.;
976  } else if( DriftDistance > moduleThickness )
977  DriftDistance = moduleThickness;
978 
979  // Assume full depletion now, partial depletion will come later.
980  float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
981  float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
982 
983  // Shift cloud center
984  float CloudCenterX = SegX + XDriftDueToMagField;
985  float CloudCenterY = SegY + YDriftDueToMagField;
986 
987  // Calculate how long is the charge drift path
988  DriftLength = sqrt( DriftDistance*DriftDistance +
989  XDriftDueToMagField*XDriftDueToMagField +
990  YDriftDueToMagField*YDriftDueToMagField );
991 
992  // What is the charge diffusion after this path
993  Sigma = sqrt(DriftLength/Dist300) * Sigma0;
994 
995  // Project the diffusion sigma on the collection plane
996  Sigma_x = Sigma / CosLorenzAngleX ;
997  Sigma_y = Sigma / CosLorenzAngleY ;
998 
999  // Insert a charge loss due to Rad Damage here
1000  float energyOnCollector = ionization_points[i].energy(); // The energy that reaches the collector
1001 
1002  // add pixel aging
1003  if (AddPixelAging) {
1004  float kValue = pixel_aging(pixelAging_,pixdet,tTopo);
1005  energyOnCollector *= exp( -1*kValue*DriftDistance/moduleThickness );
1006  }
1007 
1008 #ifdef TP_DEBUG
1009  LogDebug ("Pixel Digitizer")
1010  <<" Dift DistanceZ= "<<DriftDistance<<" module thickness= "<<moduleThickness
1011  <<" Start Energy= "<<ionization_points[i].energy()<<" Energy after loss= "<<energyOnCollector;
1012 #endif
1013  SignalPoint sp( CloudCenterX, CloudCenterY,
1014  Sigma_x, Sigma_y, hit.tof(), energyOnCollector );
1015 
1016  // Load the Charge distribution parameters
1017  collection_points[i] = (sp);
1018 
1019  } // loop over ionization points, i.
1020 
1021 } // end drift
1022 
1023 //*************************************************************************
1024 // Induce the signal on the collection plane of the active sensor area.
1026  const size_t hitIndex,
1027  const unsigned int tofBin,
1028  const PixelGeomDetUnit* pixdet,
1029  const std::vector<SignalPoint>& collection_points) {
1030 
1031  // X - Rows, Left-Right, 160, (1.6cm) for barrel
1032  // Y - Columns, Down-Up, 416, (6.4cm)
1033 
1034  const PixelTopology* topol=&pixdet->specificTopology();
1035  uint32_t detID= pixdet->geographicalId().rawId();
1036  signal_map_type& theSignal = _signal[detID];
1037 
1038 #ifdef TP_DEBUG
1039  LogDebug ("Pixel Digitizer")
1040  << " enter induce_signal, "
1041  << topol->pitch().first << " " << topol->pitch().second; //OK
1042 #endif
1043 
1044  // local map to store pixels hit by 1 Hit.
1045  typedef std::map< int, float, std::less<int> > hit_map_type;
1046  hit_map_type hit_signal;
1047 
1048  // map to store pixel integrals in the x and in the y directions
1049  std::map<int, float, std::less<int> > x,y;
1050 
1051  // Assign signals to readout channels and store sorted by channel number
1052 
1053  // Iterate over collection points on the collection plane
1054  for ( std::vector<SignalPoint>::const_iterator i=collection_points.begin();
1055  i != collection_points.end(); ++i) {
1056 
1057  float CloudCenterX = i->position().x(); // Charge position in x
1058  float CloudCenterY = i->position().y(); // in y
1059  float SigmaX = i->sigma_x(); // Charge spread in x
1060  float SigmaY = i->sigma_y(); // in y
1061  float Charge = i->amplitude(); // Charge amplitude
1062 
1063 
1064  //if(SigmaX==0 || SigmaY==0) {
1065  //cout<<SigmaX<<" "<<SigmaY
1066  // << " cloud " << i->position().x() << " " << i->position().y() << " "
1067  // << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude()<<std::endl;
1068  //}
1069 
1070 #ifdef TP_DEBUG
1071  LogDebug ("Pixel Digitizer")
1072  << " cloud " << i->position().x() << " " << i->position().y() << " "
1073  << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude();
1074 #endif
1075 
1076  // Find the maximum cloud spread in 2D plane , assume 3*sigma
1077  float CloudRight = CloudCenterX + ClusterWidth*SigmaX;
1078  float CloudLeft = CloudCenterX - ClusterWidth*SigmaX;
1079  float CloudUp = CloudCenterY + ClusterWidth*SigmaY;
1080  float CloudDown = CloudCenterY - ClusterWidth*SigmaY;
1081 
1082  // Define 2D cloud limit points
1083  LocalPoint PointRightUp = LocalPoint(CloudRight,CloudUp);
1084  LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
1085 
1086  // This points can be located outside the sensor area.
1087  // The conversion to measurement point does not check for that
1088  // so the returned pixel index might be wrong (outside range).
1089  // We rely on the limits check below to fix this.
1090  // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
1091  // PIXEL AREA IS LOST, it should not be collected.
1092 
1093  // Convert the 2D points to pixel indices
1094  MeasurementPoint mp = topol->measurementPosition(PointRightUp ); //OK
1095 
1096  int IPixRightUpX = int( floor( mp.x()));
1097  int IPixRightUpY = int( floor( mp.y()));
1098 
1099 #ifdef TP_DEBUG
1100  LogDebug ("Pixel Digitizer") << " right-up " << PointRightUp << " "
1101  << mp.x() << " " << mp.y() << " "
1102  << IPixRightUpX << " " << IPixRightUpY ;
1103 #endif
1104 
1105  mp = topol->measurementPosition(PointLeftDown ); //OK
1106 
1107  int IPixLeftDownX = int( floor( mp.x()));
1108  int IPixLeftDownY = int( floor( mp.y()));
1109 
1110 #ifdef TP_DEBUG
1111  LogDebug ("Pixel Digitizer") << " left-down " << PointLeftDown << " "
1112  << mp.x() << " " << mp.y() << " "
1113  << IPixLeftDownX << " " << IPixLeftDownY ;
1114 #endif
1115 
1116  // Check detector limits to correct for pixels outside range.
1117  int numColumns = topol->ncolumns(); // det module number of cols&rows
1118  int numRows = topol->nrows();
1119 
1120  IPixRightUpX = numRows>IPixRightUpX ? IPixRightUpX : numRows-1 ;
1121  IPixRightUpY = numColumns>IPixRightUpY ? IPixRightUpY : numColumns-1 ;
1122  IPixLeftDownX = 0<IPixLeftDownX ? IPixLeftDownX : 0 ;
1123  IPixLeftDownY = 0<IPixLeftDownY ? IPixLeftDownY : 0 ;
1124 
1125  x.clear(); // clear temporary integration array
1126  y.clear();
1127 
1128  // First integrate charge strips in x
1129  int ix; // TT for compatibility
1130  for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) { // loop over x index
1131  float xUB, xLB, UpperBound, LowerBound;
1132 
1133  // Why is set to 0 if ix=0, does it meen that we accept charge
1134  // outside the sensor? CHeck How it was done in ORCA?
1135  //if(ix == 0) LowerBound = 0.;
1136  if(ix == 0 || SigmaX==0. ) // skip for surface segemnts
1137  LowerBound = 0.;
1138  else {
1139  mp = MeasurementPoint( float(ix), 0.0);
1140  xLB = topol->localPosition(mp).x();
1141  LowerBound = 1-calcQ((xLB-CloudCenterX)/SigmaX);
1142  }
1143 
1144  if(ix == numRows-1 || SigmaX==0. )
1145  UpperBound = 1.;
1146  else {
1147  mp = MeasurementPoint( float(ix+1), 0.0);
1148  xUB = topol->localPosition(mp).x();
1149  UpperBound = 1. - calcQ((xUB-CloudCenterX)/SigmaX);
1150  }
1151 
1152  float TotalIntegrationRange = UpperBound - LowerBound; // get strip
1153  x[ix] = TotalIntegrationRange; // save strip integral
1154  //if(SigmaX==0 || SigmaY==0)
1155  //cout<<TotalIntegrationRange<<" "<<ix<<std::endl;
1156 
1157  }
1158 
1159  // Now integrate strips in y
1160  int iy; // TT for compatibility
1161  for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
1162  float yUB, yLB, UpperBound, LowerBound;
1163 
1164  if(iy == 0 || SigmaY==0.)
1165  LowerBound = 0.;
1166  else {
1167  mp = MeasurementPoint( 0.0, float(iy) );
1168  yLB = topol->localPosition(mp).y();
1169  LowerBound = 1. - calcQ((yLB-CloudCenterY)/SigmaY);
1170  }
1171 
1172  if(iy == numColumns-1 || SigmaY==0. )
1173  UpperBound = 1.;
1174  else {
1175  mp = MeasurementPoint( 0.0, float(iy+1) );
1176  yUB = topol->localPosition(mp).y();
1177  UpperBound = 1. - calcQ((yUB-CloudCenterY)/SigmaY);
1178  }
1179 
1180  float TotalIntegrationRange = UpperBound - LowerBound;
1181  y[iy] = TotalIntegrationRange; // save strip integral
1182  //if(SigmaX==0 || SigmaY==0)
1183  //cout<<TotalIntegrationRange<<" "<<iy<<std::endl;
1184  }
1185 
1186  // Get the 2D charge integrals by folding x and y strips
1187  int chan;
1188  for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) { // loop over x index
1189  for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
1190 
1191  float ChargeFraction = Charge*x[ix]*y[iy];
1192 
1193  if( ChargeFraction > 0. ) {
1194  chan = PixelDigi::pixelToChannel( ix, iy); // Get index
1195  // Load the amplitude
1196  hit_signal[chan] += ChargeFraction;
1197  } // endif
1198 
1199 
1200  mp = MeasurementPoint( float(ix), float(iy) );
1201  LocalPoint lp = topol->localPosition(mp);
1202  chan = topol->channel(lp);
1203 
1204 #ifdef TP_DEBUG
1205  LogDebug ("Pixel Digitizer")
1206  << " pixel " << ix << " " << iy << " - "<<" "
1207  << chan << " " << ChargeFraction<<" "
1208  << mp.x() << " " << mp.y() <<" "
1209  << lp.x() << " " << lp.y() << " " // givex edge position
1210  << chan; // edge belongs to previous ?
1211 #endif
1212 
1213  } // endfor iy
1214  } //endfor ix
1215 
1216 
1217  // Test conversions (THIS IS FOR TESTING ONLY) comment-out.
1218  // mp = topol->measurementPosition( i->position() ); //OK
1219  // LocalPoint lp = topol->localPosition(mp); //OK
1220  // std::pair<float,float> p = topol->pixel( i->position() ); //OK
1221  // chan = PixelDigi::pixelToChannel( int(p.first), int(p.second));
1222  // std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1223  // MeasurementPoint mp1 = MeasurementPoint( float(ip.first),
1224  // float(ip.second) );
1225  // LogDebug ("Pixel Digitizer") << " Test "<< mp.x() << " " << mp.y()
1226  // << " "<< lp.x() << " " << lp.y() << " "<<" "
1227  // <<p.first <<" "<<p.second<<" "<<chan<< " "
1228  // <<" " << ip.first << " " << ip.second << " "
1229  // << mp1.x() << " " << mp1.y() << " " //OK
1230  // << topol->localPosition(mp1).x() << " " //OK
1231  // << topol->localPosition(mp1).y() << " "
1232  // << topol->channel( i->position() ); //OK
1233 
1234 
1235  } // loop over charge distributions
1236 
1237  // Fill the global map with all hit pixels from this event
1238 
1239  for ( hit_map_type::const_iterator im = hit_signal.begin();
1240  im != hit_signal.end(); ++im) {
1241  int chan = (*im).first;
1242  theSignal[chan] += (makeDigiSimLinks_ ? Amplitude( (*im).second, &hit, hitIndex, tofBin, (*im).second) : Amplitude( (*im).second, (*im).second) ) ;
1243 
1244 #ifdef TP_DEBUG
1245  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1246  LogDebug ("Pixel Digitizer")
1247  << " pixel " << ip.first << " " << ip.second << " "
1248  << theSignal[chan];
1249 #endif
1250  }
1251 
1252 } // end induce_signal
1253 
1254 /***********************************************************************/
1255 
1256 // Build pixels, check threshold, add misscalibration, ...
1257 void SiPixelDigitizerAlgorithm::make_digis(float thePixelThresholdInE,
1258  uint32_t detID,
1259  const PixelGeomDetUnit* pixdet,
1260  std::vector<PixelDigi>& digis,
1261  std::vector<PixelDigiSimLink>& simlinks,
1262  const TrackerTopology *tTopo) const {
1263 
1264 #ifdef TP_DEBUG
1265  LogDebug ("Pixel Digitizer") << " make digis "<<" "
1266  << " pixel threshold FPix" << theThresholdInE_FPix << " "
1267  << " pixel threshold BPix" << theThresholdInE_BPix << " "
1268  << " pixel threshold BPix Layer1" << theThresholdInE_BPix_L1 << " "
1269  << " pixel threshold BPix Layer2" << theThresholdInE_BPix_L2 << " "
1270  << " List pixels passing threshold ";
1271 #endif
1272 
1273  // Loop over hit pixels
1274 
1275  signalMaps::const_iterator it = _signal.find(detID);
1276  if (it == _signal.end()) {
1277  return;
1278  }
1279 
1280  const signal_map_type& theSignal = (*it).second;
1281 
1282  // unsigned long is enough to store SimTrack id and EncodedEventId
1283  using TrackEventId = std::pair<decltype(SimTrack().trackId()), decltype(EncodedEventId().rawId())>;
1284  std::map<TrackEventId, float> simi; // re-used
1285 
1286  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1287 
1288  float signalInElectrons = (*i).second ; // signal in electrons
1289 
1290  // Do the miss calibration for calibration studies only.
1291  //if(doMissCalibrate) signalInElectrons = missCalibrate(signalInElectrons)
1292 
1293  // Do only for pixels above threshold
1294 
1295  if( signalInElectrons >= thePixelThresholdInE) { // check threshold
1296 
1297  int chan = (*i).first; // channel number
1298  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1299  int adc=0; // ADC count as integer
1300 
1301  // Do the miss calibration for calibration studies only.
1302  if(doMissCalibrate) {
1303  int row = ip.first; // X in row
1304  int col = ip.second; // Y is in col
1305  adc = int(missCalibrate(detID, tTopo, pixdet, col, row, signalInElectrons)); //full misscalib.
1306  } else { // Just do a simple electron->adc conversion
1307  adc = int( signalInElectrons / theElectronPerADC ); // calibrate gain
1308  }
1309  adc = std::min(adc, theAdcFullScale); // Check maximum value
1310 // Calculate layerIndex
1312  if(pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2OTB) { // Phase 2 OT Barrel only
1313  // Set to 1 if over the threshold
1314  if (theAdcFullScaleStack==1) {adc=1;}
1315  // Make it a linear fit to the full scale of the normal adc count. Start new adc from 1 not zero.
1317  }
1318  } // Only enter this if the Adc changes for the outer layers
1319 #ifdef TP_DEBUG
1320  LogDebug ("Pixel Digitizer")
1321  << (*i).first << " " << (*i).second << " " << signalInElectrons
1322  << " " << adc << ip.first << " " << ip.second ;
1323 #endif
1324 
1325  // Load digis
1326  digis.emplace_back(ip.first, ip.second, adc);
1327 
1328  if (makeDigiSimLinks_ && !(*i).second.hitInfos().empty()) {
1329  //digilink
1330  unsigned int il=0;
1331  for(const auto& info: (*i).second.hitInfos()) {
1332  // note: according to C++ standard operator[] does
1333  // value-initializiation, which for float means initial value of 0
1334  simi[std::make_pair(info.trackId(), info.eventId().rawId())] += (*i).second.individualampl()[il];
1335  il++;
1336  }
1337 
1338  //sum the contribution of the same trackid
1339  for(const auto& info: (*i).second.hitInfos()) {
1340  // skip if track already processed
1341  auto found = simi.find(std::make_pair(info.trackId(), info.eventId().rawId()));
1342  if(found == simi.end())
1343  continue;
1344 
1345  float sum_samechannel = found->second;
1346  float fraction=sum_samechannel/(*i).second;
1347  if(fraction>1.f) fraction=1.f;
1348 
1349  // Approximation: pick hitIndex and tofBin only from the first SimHit
1350  simlinks.emplace_back((*i).first, info.trackId(), info.hitIndex(), info.tofBin(), info.eventId(), fraction);
1351  simi.erase(found);
1352  }
1353  simi.clear(); // although should be empty already
1354  }
1355  }
1356  }
1357 }
1358 
1359 /***********************************************************************/
1360 
1361 // Add electronic noise to pixel charge
1363  float thePixelThreshold,
1364  CLHEP::HepRandomEngine* engine) {
1365 
1366 #ifdef TP_DEBUG
1367  LogDebug ("Pixel Digitizer") << " enter add_noise " << theNoiseInElectrons;
1368 #endif
1369 
1370  uint32_t detID= pixdet->geographicalId().rawId();
1371  signal_map_type& theSignal = _signal[detID];
1372 
1373 
1374  // First add noise to hit pixels
1375  float theSmearedChargeRMS = 0.0;
1376 
1377  for ( signal_map_iterator i = theSignal.begin(); i != theSignal.end(); i++) {
1378 
1380  {
1381  if((*i).second < 3000)
1382  {
1383  theSmearedChargeRMS = 543.6 - (*i).second * 0.093;
1384  } else if((*i).second < 6000){
1385  theSmearedChargeRMS = 307.6 - (*i).second * 0.01;
1386  } else{
1387  theSmearedChargeRMS = -432.4 +(*i).second * 0.123;
1388  }
1389 
1390  // Noise from Vcal smearing:
1391  float noise_ChargeVCALSmearing = theSmearedChargeRMS * CLHEP::RandGaussQ::shoot(engine, 0., 1.);
1392  // Noise from full readout:
1393  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1394 
1395  if(((*i).second + Amplitude(noise+noise_ChargeVCALSmearing, -1.)) < 0. ) {
1396  (*i).second.set(0);}
1397  else{
1398  (*i).second +=Amplitude(noise+noise_ChargeVCALSmearing, -1.);
1399  }
1400 
1401  } // End if addChargeVCalSmearing
1402  else
1403  {
1404  // Noise: ONLY full READOUT Noise.
1405  // Use here the FULL readout noise, including TBM,ALT,AOH,OPT-REC.
1406  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1407 
1408  if(((*i).second + Amplitude(noise, -1.)) < 0. ) {
1409  (*i).second.set(0);}
1410  else{
1411  (*i).second +=Amplitude(noise, -1.);
1412  }
1413  } // end if only Noise from full readout
1414 
1415  }
1416 
1417  if(!addNoisyPixels) // Option to skip noise in non-hit pixels
1418  return;
1419 
1420  const PixelTopology* topol=&pixdet->specificTopology();
1421  int numColumns = topol->ncolumns(); // det module number of cols&rows
1422  int numRows = topol->nrows();
1423 
1424  // Add noise on non-hit pixels
1425  // Use here the pixel noise
1426  int numberOfPixels = (numRows * numColumns);
1427  std::map<int,float, std::less<int> > otherPixels;
1428  std::map<int,float, std::less<int> >::iterator mapI;
1429 
1430  theNoiser->generate(numberOfPixels,
1431  thePixelThreshold, //thr. in un. of nois
1432  theNoiseInElectrons, // noise in elec.
1433  otherPixels,
1434  engine );
1435 
1436 #ifdef TP_DEBUG
1437  LogDebug ("Pixel Digitizer")
1438  << " Add noisy pixels " << numRows << " "
1439  << numColumns << " " << theNoiseInElectrons << " "
1440  << theThresholdInE_FPix << theThresholdInE_BPix <<" "<< numberOfPixels<<" "
1441  << otherPixels.size() ;
1442 #endif
1443 
1444  // Add noisy pixels
1445  for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
1446  int iy = ((*mapI).first) / numRows;
1447  int ix = ((*mapI).first) - (iy*numRows);
1448 
1449  // Keep for a while for testing.
1450  if( iy < 0 || iy > (numColumns-1) )
1451  LogWarning ("Pixel Geometry") << " error in iy " << iy ;
1452  if( ix < 0 || ix > (numRows-1) )
1453  LogWarning ("Pixel Geometry") << " error in ix " << ix ;
1454 
1455  int chan = PixelDigi::pixelToChannel(ix, iy);
1456 
1457 #ifdef TP_DEBUG
1458  LogDebug ("Pixel Digitizer")
1459  <<" Storing noise = " << (*mapI).first << " " << (*mapI).second
1460  << " " << ix << " " << iy << " " << chan ;
1461 #endif
1462 
1463  if(theSignal[chan] == 0){
1464  // float noise = float( (*mapI).second );
1465  int noise=int( (*mapI).second );
1466  theSignal[chan] = Amplitude (noise, -1.);
1467  }
1468  }
1469 }
1470 
1471 /***********************************************************************/
1472 
1473 // Simulate the readout inefficiencies.
1474 // Delete a selected number of single pixels, dcols and rocs.
1476  const PixelGeomDetUnit* pixdet,
1477  const TrackerTopology *tTopo,
1478  CLHEP::HepRandomEngine* engine) {
1479 
1480  uint32_t detID= pixdet->geographicalId().rawId();
1481  signal_map_type& theSignal = _signal[detID];
1482  const PixelTopology* topol=&pixdet->specificTopology();
1483  int numColumns = topol->ncolumns(); // det module number of cols&rows
1484  int numRows = topol->nrows();
1485 
1486  // Predefined efficiencies
1487  double pixelEfficiency = 1.0;
1488  double columnEfficiency = 1.0;
1489  double chipEfficiency = 1.0;
1490 
1491  if (eff.FromConfig) {
1492  // setup the chip indices conversion
1494  pixdet->subDetector()==GeomDetEnumerators::SubDetector::P1PXB){// barrel layers
1495  int layerIndex=tTopo->layer(detID);
1496  pixelEfficiency = eff.thePixelEfficiency[layerIndex-1];
1497  columnEfficiency = eff.thePixelColEfficiency[layerIndex-1];
1498  chipEfficiency = eff.thePixelChipEfficiency[layerIndex-1];
1499  //std::cout <<"Using BPix columnEfficiency = "<<columnEfficiency<< " for layer = "<<layerIndex <<"\n";
1500  // This should never happen, but only check if it is not an upgrade geometry
1501  if (NumberOfBarrelLayers==3){
1502  if(numColumns>416) LogWarning ("Pixel Geometry") <<" wrong columns in barrel "<<numColumns;
1503  if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in barrel "<<numRows;
1504 
1505  int ladder=tTopo->pxbLadder(detID);
1506  int module=tTopo->pxbModule(detID);
1507  if (module<=4) module=5-module;
1508  else module-=4;
1509 
1510  columnEfficiency *= eff.theLadderEfficiency_BPix[layerIndex-1][ladder-1]*eff.theModuleEfficiency_BPix[layerIndex-1][module-1]*eff.pu_scale[layerIndex-1];
1511  }
1514  pixdet->subDetector()==GeomDetEnumerators::SubDetector::P2PXEC){ // forward disks
1515 
1516  unsigned int diskIndex=tTopo->layer(detID)+eff.FPixIndex; // Use diskIndex-1 later to stay consistent with BPix
1517  unsigned int panelIndex=tTopo->pxfPanel(detID);
1518  unsigned int moduleIndex=tTopo->pxfModule(detID);
1519  //if (eff.FPixIndex>diskIndex-1){throw cms::Exception("Configuration") <<"SiPixelDigitizer is using the wrong efficiency value. index = "
1520  // <<diskIndex-1<<" , MinIndex = "<<eff.FPixIndex<<" ... "<<tTopo->pxfDisk(detID);}
1521  pixelEfficiency = eff.thePixelEfficiency[diskIndex-1];
1522  columnEfficiency = eff.thePixelColEfficiency[diskIndex-1];
1523  chipEfficiency = eff.thePixelChipEfficiency[diskIndex-1];
1524  //std::cout <<"Using FPix columnEfficiency = "<<columnEfficiency<<" for Disk = "<< tTopo->pxfDisk(detID)<<"\n";
1525  // Sometimes the forward pixels have wrong size,
1526  // this crashes the index conversion, so exit, but only check if it is not an upgrade geometry
1527  if (NumberOfBarrelLayers==3){ // whether it is the present or the phase 1 detector can be checked using GeomDetEnumerators::SubDetector
1528  if(numColumns>260 || numRows>160) {
1529  if(numColumns>260) LogWarning ("Pixel Geometry") <<" wrong columns in endcaps "<<numColumns;
1530  if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in endcaps "<<numRows;
1531  return;
1532  }
1533  if ((panelIndex==1 && (moduleIndex==1 || moduleIndex==2)) || (panelIndex==2 && moduleIndex==1)) { //inner modules
1534  columnEfficiency*=eff.theInnerEfficiency_FPix[diskIndex-1]*eff.pu_scale[3];
1535  } else { //outer modules
1536  columnEfficiency*=eff.theOuterEfficiency_FPix[diskIndex-1]*eff.pu_scale[4];
1537  }
1538  } // current detector, forward
1540  // If phase 2 outer tracker, hardcoded values as they have been so far
1541  pixelEfficiency = 0.999;
1542  columnEfficiency = 0.999;
1543  chipEfficiency = 0.999;
1544  } // if barrel/forward
1545  } else { // Load precomputed factors from Database
1546  pixelEfficiency = eff.PixelGeomFactors.at(detID);
1547  columnEfficiency = eff.ColGeomFactors.at(detID)*eff.pu_scale[eff.iPU.at(detID)];
1548  chipEfficiency = eff.ChipGeomFactors.at(detID);
1549  }
1550 
1551 #ifdef TP_DEBUG
1552  LogDebug ("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " "
1553  << columnEfficiency << " " << chipEfficiency;
1554 #endif
1555 
1556  // Initilize the index converter
1557  //PixelIndices indexConverter(numColumns,numRows);
1558  std::unique_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
1559 
1560  int chipIndex = 0;
1561  int rowROC = 0;
1562  int colROC = 0;
1563  std::map<int, int, std::less<int> >chips, columns;
1564  std::map<int, int, std::less<int> >::iterator iter;
1565 
1566  // Find out the number of columns and rocs hits
1567  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1568  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1569 
1570  int chan = i->first;
1571  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1572  int row = ip.first; // X in row
1573  int col = ip.second; // Y is in col
1574  //transform to ROC index coordinates
1575  pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1576  int dColInChip = pIndexConverter->DColumn(colROC); // get ROC dcol from ROC col
1577  //dcol in mod
1578  int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
1579 
1580  chips[chipIndex]++;
1581  columns[dColInDet]++;
1582  }
1583 
1584  // Delete some ROC hits.
1585  for ( iter = chips.begin(); iter != chips.end() ; iter++ ) {
1586  //float rand = RandFlat::shoot();
1587  float rand = CLHEP::RandFlat::shoot(engine);
1588  if( rand > chipEfficiency ) chips[iter->first]=0;
1589  }
1590 
1591  // Delete some Dcol hits.
1592  for ( iter = columns.begin(); iter != columns.end() ; iter++ ) {
1593  //float rand = RandFlat::shoot();
1594  float rand = CLHEP::RandFlat::shoot(engine);
1595  if( rand > columnEfficiency ) columns[iter->first]=0;
1596  }
1597 
1598  // Now loop again over pixels to kill some of them.
1599  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1600  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1601 
1602  // int chan = i->first;
1603  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1604  int row = ip.first; // X in row
1605  int col = ip.second; // Y is in col
1606  //transform to ROC index coordinates
1607  pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1608  int dColInChip = pIndexConverter->DColumn(colROC); //get ROC dcol from ROC col
1609  //dcol in mod
1610  int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
1611 
1612  //float rand = RandFlat::shoot();
1613  float rand = CLHEP::RandFlat::shoot(engine);
1614  if( chips[chipIndex]==0 || columns[dColInDet]==0
1615  || rand>pixelEfficiency ) {
1616  // make pixel amplitude =0, pixel will be lost at clusterization
1617  i->second.set(0.); // reset amplitude,
1618  } // end if
1619 
1620  } // end pixel loop
1621 } // end pixel_indefficiency
1622 
1623 //***************************************************************************************
1624 // Simulate pixel aging with an exponential function
1625 //**************************************************************************************
1626 
1628  const PixelGeomDetUnit *pixdet,
1629  const TrackerTopology *tTopo) const {
1630 
1631  uint32_t detID= pixdet->geographicalId().rawId();
1632 
1633 
1634  // Predefined damage parameter (no aging)
1635  float pseudoRadDamage = 0.0;
1636 
1637  // setup the chip indices conversion
1639  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P1PXB){// barrel layers
1640  int layerIndex=tTopo->layer(detID);
1641 
1642  pseudoRadDamage = aging.thePixelPseudoRadDamage[layerIndex-1];
1643 
1644  // std::cout << "pixel_aging: " << std::endl;
1645  // std::cout << "Subid " << Subid << " layerIndex " << layerIndex << " ladder " << tTopo->pxbLadder(detID) << " module " << tTopo->pxbModule(detID) << std::endl;
1646 
1649  pixdet->subDetector() == GeomDetEnumerators::SubDetector::P2PXEC) { // forward disks
1650  unsigned int diskIndex=tTopo->layer(detID)+aging.FPixIndex; // Use diskIndex-1 later to stay consistent with BPix
1651 
1652  pseudoRadDamage = aging.thePixelPseudoRadDamage[diskIndex-1];
1653 
1654  // std::cout << "pixel_aging: " << std::endl;
1655  // std::cout << "Subid " << Subid << " diskIndex " << diskIndex << std::endl;
1657  // if phase 2 OT hardcoded value as it has always been
1658  pseudoRadDamage = 0.;
1659  } // if barrel/forward
1660 
1661  // std::cout << " pseudoRadDamage " << pseudoRadDamage << std::endl;
1662  // std::cout << " end pixel_aging " << std::endl;
1663 
1664  return pseudoRadDamage;
1665 #ifdef TP_DEBUG
1666  LogDebug ("Pixel Digitizer") << " enter pixel_aging " << pseudoRadDamage;
1667 #endif
1668 
1669 }
1670 
1671 //***********************************************************************
1672 
1673 // Fluctuate the gain and offset for the amplitude calibration
1674 // Use gaussian smearing.
1675 //float SiPixelDigitizerAlgorithm::missCalibrate(const float amp) const {
1676  //float gain = RandGaussQ::shoot(1.,theGainSmearing);
1677  //float offset = RandGaussQ::shoot(0.,theOffsetSmearing);
1678  //float newAmp = amp * gain + offset;
1679  // More complex misscalibration
1680 float SiPixelDigitizerAlgorithm::missCalibrate(uint32_t detID, const TrackerTopology *tTopo, const PixelGeomDetUnit* pixdet, int col,int row,
1681  const float signalInElectrons) const {
1682  // Central values
1683  //const float p0=0.00352, p1=0.868, p2=112., p3=113.; // pix(0,0,0)
1684  // const float p0=0.00382, p1=0.886, p2=112.7, p3=113.0; // average roc=0
1685  //const float p0=0.00492, p1=1.998, p2=90.6, p3=134.1; // average roc=6
1686  // Smeared (rms)
1687  //const float s0=0.00020, s1=0.051, s2=5.4, s3=4.4; // average roc=0
1688  //const float s0=0.00015, s1=0.043, s2=3.2, s3=3.1; // col average roc=0
1689 
1690  // Make 2 sets of parameters for Fpix and BPIx:
1691 
1692  float p0=0.0;
1693  float p1=0.0;
1694  float p2=0.0;
1695  float p3=0.0;
1696 
1697  if(pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()){// barrel layers
1698  p0 = BPix_p0;
1699  p1 = BPix_p1;
1700  p2 = BPix_p2;
1701  p3 = BPix_p3;
1702  } else if(pixdet->type().isTrackerPixel()) {// forward disks
1703  p0 = FPix_p0;
1704  p1 = FPix_p1;
1705  p2 = FPix_p2;
1706  p3 = FPix_p3;
1707  } else {
1708  throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;
1709  }
1710 
1711  float newAmp = 0.; //Modified signal
1712 
1713  // Convert electrons to VCAL units
1714  float signal = (signalInElectrons-electronsPerVCAL_Offset)/electronsPerVCAL;
1715 
1716  // New gains/offsets are needed for phase1 L1
1717  int layer = 0;
1718  if (DetId(detID).subdetId()==1) layer = tTopo->pxbLayer(detID);
1719  if (layer==1) signal = (signalInElectrons-electronsPerVCAL_L1_Offset)/electronsPerVCAL_L1;
1720 
1721  // Simulate the analog response with fixed parametrization
1722  newAmp = p3 + p2 * tanh(p0*signal - p1);
1723 
1724 
1725  // Use the pixel-by-pixel calibrations
1726  //transform to ROC index coordinates
1727  //int chipIndex=0, colROC=0, rowROC=0;
1728  //std::unique_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
1729  //pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1730 
1731  // Use calibration from a file
1732  //int chanROC = PixelIndices::pixelToChannelROC(rowROC,colROC); // use ROC coordinates
1733  //float pp0=0, pp1=0,pp2=0,pp3=0;
1734  //map<int,CalParameters,std::less<int> >::const_iterator it=calmap.find(chanROC);
1735  //CalParameters y = (*it).second;
1736  //pp0 = y.p0;
1737  //pp1 = y.p1;
1738  //pp2 = y.p2;
1739  //pp3 = y.p3;
1740 
1741  //
1742  // Use random smearing
1743  // Randomize the pixel response
1744  //float pp0 = RandGaussQ::shoot(p0,s0);
1745  //float pp1 = RandGaussQ::shoot(p1,s1);
1746  //float pp2 = RandGaussQ::shoot(p2,s2);
1747  //float pp3 = RandGaussQ::shoot(p3,s3);
1748 
1749  //newAmp = pp3 + pp2 * tanh(pp0*signal - pp1); // Final signal
1750 
1751  //cout<<" misscalibrate "<<col<<" "<<row<<" "<<chipIndex<<" "<<colROC<<" "
1752  // <<rowROC<<" "<<signalInElectrons<<" "<<signal<<" "<<newAmp<<" "
1753  // <<(signalInElectrons/theElectronPerADC)<<std::endl;
1754 
1755  return newAmp;
1756 }
1757 //******************************************************************************
1758 
1759 // Set the drift direction accoring to the Bfield in local det-unit frame
1760 // Works for both barrel and forward pixels.
1761 // Replace the sign convention to fit M.Swartz's formulaes.
1762 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
1763 // parameter value
1764 
1766  const GlobalVector& bfield,
1767  const DetId& detId) const {
1768  Frame detFrame(pixdet->surface().position(),pixdet->surface().rotation());
1769  LocalVector Bfield=detFrame.toLocal(bfield);
1770 
1771  float alpha2_FPix;
1772  float alpha2_BPix;
1773  float alpha2;
1774 
1775  //float dir_x = -tanLorentzAnglePerTesla * Bfield.y();
1776  //float dir_y = +tanLorentzAnglePerTesla * Bfield.x();
1777  //float dir_z = -1.; // E field always in z direction, so electrons go to -z
1778  // The dir_z has to be +/- 1. !
1779  // LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
1780 
1781  float dir_x = 0.0;
1782  float dir_y = 0.0;
1783  float dir_z = 0.0;
1784  float scale = 0.0;
1785 
1786  uint32_t detID= pixdet->geographicalId().rawId();
1787 
1788 
1789  // Read Lorentz angle from cfg file:**************************************************************
1790 
1791  if(!use_LorentzAngle_DB_){
1792 
1793  if( alpha2Order) {
1796  }else {
1797  alpha2_FPix = 0.0;
1798  alpha2_BPix = 0.0;
1799  }
1800 
1801  if(pixdet->type().isTrackerPixel() && pixdet->type().isBarrel()){// barrel layers
1802  dir_x = -( tanLorentzAnglePerTesla_BPix * Bfield.y() + alpha2_BPix* Bfield.z()* Bfield.x() );
1803  dir_y = +( tanLorentzAnglePerTesla_BPix * Bfield.x() - alpha2_BPix* Bfield.z()* Bfield.y() );
1804  dir_z = -(1 + alpha2_BPix* Bfield.z()*Bfield.z() );
1805  scale = -dir_z;
1806  } else if (pixdet->type().isTrackerPixel()) {// forward disks
1807  dir_x = -( tanLorentzAnglePerTesla_FPix * Bfield.y() + alpha2_FPix* Bfield.z()* Bfield.x() );
1808  dir_y = +( tanLorentzAnglePerTesla_FPix * Bfield.x() - alpha2_FPix* Bfield.z()* Bfield.y() );
1809  dir_z = -(1 + alpha2_FPix* Bfield.z()*Bfield.z() );
1810  scale = -dir_z;
1811  } else {
1812  throw cms::Exception("NotAPixelGeomDetUnit") << "Not a pixel geomdet unit" << detID;
1813  }
1814  } // end: Read LA from cfg file.
1815 
1816  //Read Lorentz angle from DB:********************************************************************
1818  float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
1819  alpha2 = lorentzAngle * lorentzAngle;
1820  //std::cout << "detID is: "<< it->first <<"The LA per tesla is: "<< it->second << std::std::endl;
1821  dir_x = -( lorentzAngle * Bfield.y() + alpha2 * Bfield.z()* Bfield.x() );
1822  dir_y = +( lorentzAngle * Bfield.x() - alpha2 * Bfield.z()* Bfield.y() );
1823  dir_z = -(1 + alpha2 * Bfield.z()*Bfield.z() );
1824  scale = -dir_z;
1825  }// end: Read LA from DataBase.
1826 
1827  LocalVector theDriftDirection = LocalVector(dir_x/scale, dir_y/scale, dir_z/scale );
1828 
1829 #ifdef TP_DEBUG
1830  LogDebug ("Pixel Digitizer") << " The drift direction in local coordinate is "
1831  << theDriftDirection ;
1832 #endif
1833 
1834  return theDriftDirection;
1835 }
1836 
1837 //****************************************************************************************************
1838 
1840 
1841  signal_map_type& theSignal = _signal[detID];
1842 
1843  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1844  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1845 
1846  // int chan = i->first;
1847  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1848  int row = ip.first; // X in row
1849  int col = ip.second; // Y is in col
1850  //transform to ROC index coordinates
1851  if(theSiPixelGainCalibrationService_->isDead(detID, col, row)){
1852  // std::cout << "now in isdead check, row " << detID << " " << col << "," << row << std::std::endl;
1853  // make pixel amplitude =0, pixel will be lost at clusterization
1854  i->second.set(0.); // reset amplitude,
1855  } // end if
1856  } // end pixel loop
1857 } // end pixel_indefficiency
1858 
1859 
1860 //****************************************************************************************************
1861 
1863 
1864  bool isbad=false;
1865 
1866  Parameters::const_iterator itDeadModules=DeadModules.begin();
1867 
1868  int detid = detID;
1869  for(; itDeadModules != DeadModules.end(); ++itDeadModules){
1870  int Dead_detID = itDeadModules->getParameter<int>("Dead_detID");
1871  if(detid == Dead_detID){
1872  isbad=true;
1873  break;
1874  }
1875  }
1876 
1877  if(!isbad)
1878  return;
1879 
1880  signal_map_type& theSignal = _signal[detID];
1881 
1882  std::string Module = itDeadModules->getParameter<std::string>("Module");
1883 
1884  if(Module=="whole"){
1885  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1886  i->second.set(0.); // reset amplitude
1887  }
1888  }
1889 
1890  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1891  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1892 
1893  if(Module=="tbmA" && ip.first>=80 && ip.first<=159){
1894  i->second.set(0.);
1895  }
1896 
1897  if( Module=="tbmB" && ip.first<=79){
1898  i->second.set(0.);
1899  }
1900  }
1901 }
1902 //****************************************************************************************************
1904 // Not SLHC safe for now
1905 
1906  bool isbad=false;
1907 
1908  std::vector<SiPixelQuality::disabledModuleType>disabledModules = SiPixelBadModule_->getBadComponentList();
1909 
1911 
1912  for (size_t id=0;id<disabledModules.size();id++)
1913  {
1914  if(detID==disabledModules[id].DetID){
1915  isbad=true;
1916  badmodule = disabledModules[id];
1917  break;
1918  }
1919  }
1920 
1921  if(!isbad)
1922  return;
1923 
1924  signal_map_type& theSignal = _signal[detID];
1925 
1926  //std::cout<<"Hit in: "<< detID <<" errorType "<< badmodule.errorType<<" BadRocs="<<std::hex<<SiPixelBadModule_->getBadRocs(detID)<<dec<<" "<<std::endl;
1927  if(badmodule.errorType == 0){ // this is a whole dead module.
1928 
1929  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1930  i->second.set(0.); // reset amplitude
1931  }
1932  }
1933  else { // all other module types: half-modules and single ROCs.
1934  // Get Bad ROC position:
1935  //follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
1936  std::vector<GlobalPixel> badrocpositions (0);
1937  for(unsigned int j = 0; j < 16; j++){
1938  if(SiPixelBadModule_->IsRocBad(detID, j) == true){
1939 
1940  std::vector<CablingPathToDetUnit> path = map_.product()->pathToDetUnit(detID);
1941  typedef std::vector<CablingPathToDetUnit>::const_iterator IT;
1942  for (IT it = path.begin(); it != path.end(); ++it) {
1943  const PixelROC* myroc = map_.product()->findItem(*it);
1944  if( myroc->idInDetUnit() == j) {
1945  LocalPixel::RocRowCol local = { 39, 25}; //corresponding to center of ROC row, col
1946  GlobalPixel global = myroc->toGlobal( LocalPixel(local) );
1947  badrocpositions.push_back(global);
1948  break;
1949  }
1950  }
1951  }
1952  }// end of getBadRocPositions
1953 
1954 
1955  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1956  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1957 
1958  for(std::vector<GlobalPixel>::const_iterator it = badrocpositions.begin(); it != badrocpositions.end(); ++it){
1959  if(it->row >= 80 && ip.first >= 80 ){
1960  if((std::abs(ip.second - it->col) < 26) ) {i->second.set(0.);}
1961  else if(it->row==120 && ip.second-it->col==26){i->second.set(0.);}
1962  else if(it->row==119 && it->col-ip.second==26){i->second.set(0.);}
1963  }
1964  else if(it->row < 80 && ip.first < 80 ){
1965  if((std::abs(ip.second - it->col) < 26) ){i->second.set(0.);}
1966  else if(it->row==40 && ip.second-it->col==26){i->second.set(0.);}
1967  else if(it->row==39 && it->col-ip.second==26){i->second.set(0.);}
1968  }
1969  }
1970  }
1971  }
1972 }
1973 
int adc(sample_type sample)
get the ADC sample (12 bits)
#define LogDebug(id)
void init(const edm::EventSetup &es)
size
Write out results.
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
Map map_
virtual int nrows() const =0
const GeomDetType & type() const override
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:72
edm::ESHandle< SiPixelFedCablingMap > map_
void init_DynIneffDB(const edm::EventSetup &, const unsigned int &)
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
bool isBarrel() const
Definition: GeomDetType.cc:13
T y() const
Definition: PV2DBase.h:46
std::map< int, CalParameters, std::less< int > > initCal() const
const std::unique_ptr< SiPixelGainCalibrationOfflineSimService > theSiPixelGainCalibrationService_
LocalVector DriftDirection(const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, const DetId &detId) const
SiPixelDigitizerAlgorithm(const edm::ParameterSet &conf)
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
const std::vector< float > & getMix_TrueInteractions() 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())
const std::vector< int > & getMix_bunchCrossing() const
unsigned int pxbLadder(const DetId &id) const
T y() const
Definition: PV3DBase.h:63
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual std::pair< float, float > pitch() const =0
const Bounds & bounds() const
Definition: Surface.h:120
std::vector< std::vector< double > > thePUEfficiency
unsigned int pxbModule(const DetId &id) const
bool IsRocBad(const uint32_t &detid, const short &rocNb) const
#define nullptr
edm::ESHandle< TrackerGeometry > geom_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
SigmaX
Definition: gun_cff.py:25
identify pixel inside single ROC
Definition: LocalPixel.h:7
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:68
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void induce_signal(const PSimHit &hit, const size_t hitIndex, 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:38
T mag() const
Definition: PV3DBase.h:67
bunchspace
in terms of 25 ns
const std::map< int, CalParameters, std::less< int > > calmap
float pixel_aging(const PixelAging &aging, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo) const
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
static int pixelToChannelROC(const int rowROC, const int colROC)
Definition: PixelIndices.h:250
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
static std::pair< int, int > channelToPixelROC(const int chan)
Definition: PixelIndices.h:253
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:40
void digitize(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
double f[11][100]
virtual int channel(const LocalPoint &p) const =0
const std::vector< disabledModuleType > getBadComponentList() const
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< LinkConnSpec >::const_iterator IT
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
bool isTrackerPixel() const
Definition: GeomDetType.cc:29
signal_map_type::iterator signal_map_iterator
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
double p2[4]
Definition: TauolaWrapper.h:90
edm::ESHandle< SiPixelLorentzAngle > SiPixelLorentzAngle_
unsigned int pxfModule(const DetId &id) const
ii
Definition: cuy.py:588
void fluctuateEloss(int particleId, float momentum, float eloss, float length, int NumberOfSegments, float elossVector[], CLHEP::HepRandomEngine *) const
unsigned int pxbLayer(const DetId &id) const
const std::unique_ptr< GaussianTailNoiseGenerator > theNoiser
Definition: DetId.h:18
std::map< int, Amplitude, std::less< int > > signal_map_type
edm::ESHandle< SiPixelQuality > SiPixelBadModule_
virtual float thickness() const =0
SigmaY
Definition: gun_cff.py:26
void primary_ionization(const PSimHit &hit, std::vector< EnergyDepositUnit > &ionization_points, CLHEP::HepRandomEngine *) const
void init_from_db(const edm::ESHandle< TrackerGeometry > &, const edm::ESHandle< SiPixelDynamicInefficiency > &)
void calculateInstlumiFactor(PileupMixingContent *puInfo)
const T & get() const
Definition: EventSetup.h:58
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:62
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
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:15
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
edm::ESHandle< SiPixelDynamicInefficiency > SiPixelDynamicInefficiency_
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
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:75
HLT enums.
virtual int ncolumns() const =0
int particleType() const
Definition: PSimHit.h:85
double p1[4]
Definition: TauolaWrapper.h:89
col
Definition: cuy.py:1008
Signal rand(Signal arg)
Definition: vlib.cc:442
float getLorentzAngle(const uint32_t &) const
const std::map< unsigned int, std::vector< double > > & getPUFactors() const
const RotationType & rotation() const
std::vector< sipixelobjects::CablingPathToDetUnit > pathToDetUnit(uint32_t rawDetId) const override
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
const sipixelobjects::PixelROC * findItem(const sipixelobjects::CablingPathToDetUnit &path) const override
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
virtual SubDetector subDetector() const
Which subdetector.
Definition: GeomDet.cc:44
T x() const
Definition: PV2DBase.h:45
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
T const * product() const
Definition: ESHandle.h:86
Definition: vlib.h:208
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
unsigned int pxfPanel(const DetId &id) const
const std::map< unsigned int, double > & getColGeomFactors() 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:93
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:45
const std::map< unsigned int, double > & getPixelGeomFactors() const
Definition: aging.py:1
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:59
void add_noise(const PixelGeomDetUnit *pixdet, float thePixelThreshold, CLHEP::HepRandomEngine *)
double p3[4]
Definition: TauolaWrapper.h:91
void pixel_inefficiency(const PixelEfficiencies &eff, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
const std::map< unsigned int, double > & getChipGeomFactors() const