CMS 3D CMS Logo

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