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