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