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