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 size_t inputBeginGlobalIndex,
455  const unsigned int tofBin,
456  const PixelGeomDetUnit* pixdet,
457  const GlobalVector& bfield,
458  CLHEP::HepRandomEngine* engine) {
459  // produce SignalPoint's for all SimHit's in detector
460  // Loop over hits
461 
462  uint32_t detId = pixdet->geographicalId().rawId();
463  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
464  for (std::vector<PSimHit>::const_iterator ssbegin = inputBegin; ssbegin != inputEnd; ++ssbegin, ++simHitGlobalIndex) {
465  // skip hits not in this detector.
466  if((*ssbegin).detUnitId() != detId) {
467  continue;
468  }
469 
470 #ifdef TP_DEBUG
471  LogDebug ("Pixel Digitizer")
472  << (*ssbegin).particleType() << " " << (*ssbegin).pabs() << " "
473  << (*ssbegin).energyLoss() << " " << (*ssbegin).tof() << " "
474  << (*ssbegin).trackId() << " " << (*ssbegin).processType() << " "
475  << (*ssbegin).detUnitId()
476  << (*ssbegin).entryPoint() << " " << (*ssbegin).exitPoint() ;
477 #endif
478 
479 
480  std::vector<EnergyDepositUnit> ionization_points;
481  std::vector<SignalPoint> collection_points;
482 
483  // fill collection_points for this SimHit, indpendent of topology
484  // Check the TOF cut
485  if ( ((*ssbegin).tof() - pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.)>= theTofLowerCut &&
486  ((*ssbegin).tof()- pixdet->surface().toGlobal((*ssbegin).localPosition()).mag()/30.) <= theTofUpperCut ) {
487  primary_ionization(*ssbegin, ionization_points, engine); // fills _ionization_points
488  drift(*ssbegin, pixdet, bfield, ionization_points, collection_points); // transforms _ionization_points to collection_points
489  // compute induced signal on readout elements and add to _signal
490  induce_signal(*ssbegin, simHitGlobalIndex, tofBin, pixdet, collection_points); // 1st 3 args needed only for SimHit<-->Digi link
491  } // end if
492  } // end for
493 
494 }
495 
496 //============================================================================
498  //Instlumi scalefactor calculating for dynamic inefficiency
499 
500  if (puInfo) {
501 
502  const std::vector<int> bunchCrossing = puInfo->getMix_bunchCrossing();
503  const std::vector<float> TrueInteractionList = puInfo->getMix_TrueInteractions();
504 
505  int pui = 0, p = 0;
506  std::vector<int>::const_iterator pu;
507  std::vector<int>::const_iterator pu0 = bunchCrossing.end();
508 
509  for (pu=bunchCrossing.begin(); pu!=bunchCrossing.end(); ++pu) {
510  if (*pu==0) {
511  pu0 = pu;
512  p = pui;
513  }
514  pui++;
515  }
516 
517  if (pu0!=bunchCrossing.end()) {
518  for (size_t i=0; i<3; i++) {
519  double instlumi = TrueInteractionList.at(p)*221.95;
520  double instlumi_pow=1.;
521  _pu_scale[i] = 0;
522  for (size_t j=0; j<pixelEfficiencies_.thePUEfficiency_BPix[i].size(); j++){
524  instlumi_pow*=instlumi;
525  }
526  }
527  }
528  }
529  else {
530  for (int i=0; i<3;i++) {
531  _pu_scale[i] = 1.;
532  }
533  }
534 }
535 
536 //============================================================================
538  std::vector<PixelDigi>& digis,
539  std::vector<PixelDigiSimLink>& simlinks,
540  const TrackerTopology *tTopo,
541  CLHEP::HepRandomEngine* engine) {
542 
543  // Pixel Efficiency moved from the constructor to this method because
544  // the information of the det are not available in the constructor
545  // Effciency parameters. 0 - no inefficiency, 1-low lumi, 10-high lumi
546 
547  uint32_t detID = pixdet->geographicalId().rawId();
548  const signal_map_type& theSignal = _signal[detID];
549 
550  const PixelTopology* topol=&pixdet->specificTopology();
551  int numColumns = topol->ncolumns(); // det module number of cols&rows
552  int numRows = topol->nrows();
553 
554  // Noise already defined in electrons
555  //thePixelThresholdInE = thePixelThreshold * theNoiseInElectrons ;
556  // Find the threshold in noise units, needed for the noiser.
557 
558  unsigned int Sub_detid=DetId(detID).subdetId();
559 
560  float thePixelThresholdInE = 0.;
561 
562  if(theNoiseInElectrons>0.){
563  if(Sub_detid == PixelSubdetector::PixelBarrel){ // Barrel modules
564  int lay = tTopo->pxbLayer(detID);
566  if(lay==1) {
567  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix_L1, theThresholdSmearing_BPix_L1); // gaussian smearing
568  } else {
569  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_BPix , theThresholdSmearing_BPix); // gaussian smearing
570  }
571  } else {
572  if(lay==1) {
573  thePixelThresholdInE = theThresholdInE_BPix_L1;
574  } else {
575  thePixelThresholdInE = theThresholdInE_BPix; // no smearing
576  }
577  }
578  } else { // Forward disks modules
580  thePixelThresholdInE = CLHEP::RandGaussQ::shoot(engine, theThresholdInE_FPix, theThresholdSmearing_FPix); // gaussian smearing
581  } else {
582  thePixelThresholdInE = theThresholdInE_FPix; // no smearing
583  }
584  }
585  }
586 
587 
588 #ifdef TP_DEBUG
589  // full detector thickness
590  float moduleThickness = pixdet->specificSurface().bounds().thickness();
591  LogDebug ("PixelDigitizer")
592  << " PixelDigitizer "
593  << numColumns << " " << numRows << " " << moduleThickness;
594 #endif
595 
596  if(addNoise) add_noise(pixdet, thePixelThresholdInE/theNoiseInElectrons, engine); // generate noise
597 
598  // Do only if needed
599 
600  if((AddPixelInefficiency) && (theSignal.size()>0))
601  pixel_inefficiency(pixelEfficiencies_, pixdet, tTopo, engine); // Kill some pixels
602 
603  if(use_ineff_from_db_ && (theSignal.size()>0))
604  pixel_inefficiency_db(detID);
605 
606  if(use_module_killing_) {
607  if (use_deadmodule_DB_) { // remove dead modules using DB
608  module_killing_DB(detID);
609  } else { // remove dead modules using the list in cfg file
610  module_killing_conf(detID);
611  }
612  }
613 
614  make_digis(thePixelThresholdInE, detID, digis, simlinks, tTopo);
615 
616 #ifdef TP_DEBUG
617  LogDebug ("PixelDigitizer") << "[SiPixelDigitizerAlgorithm] converted " << digis.size() << " PixelDigis in DetUnit" << detID;
618 #endif
619 }
620 
621 //***********************************************************************/
622 // Generate primary ionization along the track segment.
623 // Divide the track into small sub-segments
624 void SiPixelDigitizerAlgorithm::primary_ionization(const PSimHit& hit, std::vector<EnergyDepositUnit>& ionization_points, CLHEP::HepRandomEngine* engine) const {
625 
626 // Straight line approximation for trajectory inside active media
627 
628  const float SegmentLength = 0.0010; //10microns in cm
629  float energy;
630 
631  // Get the 3D segment direction vector
632  LocalVector direction = hit.exitPoint() - hit.entryPoint();
633 
634  float eLoss = hit.energyLoss(); // Eloss in GeV
635  float length = direction.mag(); // Track length in Silicon
636 
637  int NumberOfSegments = int ( length / SegmentLength); // Number of segments
638  if(NumberOfSegments < 1) NumberOfSegments = 1;
639 
640 #ifdef TP_DEBUG
641  LogDebug ("Pixel Digitizer")
642  << " enter primary_ionzation " << NumberOfSegments
643  << " shift = "
644  << (hit.exitPoint().x()-hit.entryPoint().x()) << " "
645  << (hit.exitPoint().y()-hit.entryPoint().y()) << " "
646  << (hit.exitPoint().z()-hit.entryPoint().z()) << " "
647  << hit.particleType() <<" "<< hit.pabs() ;
648 #endif
649 
650  float* elossVector = new float[NumberOfSegments]; // Eloss vector
651 
652  if( fluctuateCharge ) {
653  //MP DA RIMUOVERE ASSOLUTAMENTE
654  int pid = hit.particleType();
655  //int pid=211; // assume it is a pion
656 
657  float momentum = hit.pabs();
658  // Generate fluctuated charge points
659  fluctuateEloss(pid, momentum, eLoss, length, NumberOfSegments,
660  elossVector, engine);
661  }
662 
663  ionization_points.resize( NumberOfSegments); // set size
664 
665  // loop over segments
666  for ( int i = 0; i != NumberOfSegments; i++) {
667  // Divide the segment into equal length subsegments
668  Local3DPoint point = hit.entryPoint() +
669  float((i+0.5)/NumberOfSegments) * direction;
670 
671  if( fluctuateCharge )
672  energy = elossVector[i]/GeVperElectron; // Convert charge to elec.
673  else
674  energy = hit.energyLoss()/GeVperElectron/float(NumberOfSegments);
675 
676  EnergyDepositUnit edu( energy, point); //define position,energy point
677  ionization_points[i] = edu; // save
678 
679 #ifdef TP_DEBUG
680  LogDebug ("Pixel Digitizer")
681  << i << " " << ionization_points[i].x() << " "
682  << ionization_points[i].y() << " "
683  << ionization_points[i].z() << " "
684  << ionization_points[i].energy();
685 #endif
686 
687  } // end for loop
688 
689  delete[] elossVector;
690 
691 }
692 //******************************************************************************
693 
694 // Fluctuate the charge comming from a small (10um) track segment.
695 // Use the G4 routine. For mip pions for the moment.
696 void SiPixelDigitizerAlgorithm::fluctuateEloss(int pid, float particleMomentum,
697  float eloss, float length,
698  int NumberOfSegs,float elossVector[],
699  CLHEP::HepRandomEngine* engine) const {
700 
701  // Get dedx for this track
702  //float dedx;
703  //if( length > 0.) dedx = eloss/length;
704  //else dedx = eloss;
705 
706  double particleMass = 139.6; // Mass in MeV, Assume pion
707  pid = std::abs(pid);
708  if(pid!=211) { // Mass in MeV
709  if(pid==11) particleMass = 0.511;
710  else if(pid==13) particleMass = 105.7;
711  else if(pid==321) particleMass = 493.7;
712  else if(pid==2212) particleMass = 938.3;
713  }
714  // What is the track segment length.
715  float segmentLength = length/NumberOfSegs;
716 
717  // Generate charge fluctuations.
718  float de=0.;
719  float sum=0.;
720  double segmentEloss = (1000.*eloss)/NumberOfSegs; //eloss in MeV
721  for (int i=0;i<NumberOfSegs;i++) {
722  // material,*, momentum,energy,*, *, mass
723  //myglandz_(14.,segmentLength,2.,2.,dedx,de,0.14);
724  // The G4 routine needs momentum in MeV, mass in Mev, delta-cut in MeV,
725  // track segment length in mm, segment eloss in MeV
726  // Returns fluctuated eloss in MeV
727  double deltaCutoff = tMax; // the cutoff is sometimes redefined inside, so fix it.
728  de = fluctuate->SampleFluctuations(double(particleMomentum*1000.),
729  particleMass, deltaCutoff,
730  double(segmentLength*10.),
731  segmentEloss, engine )/1000.; //convert to GeV
732 
733  elossVector[i]=de;
734  sum +=de;
735  }
736 
737  if(sum>0.) { // If fluctuations give eloss>0.
738  // Rescale to the same total eloss
739  float ratio = eloss/sum;
740 
741  for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= ratio*elossVector[ii];
742  } else { // If fluctuations gives 0 eloss
743  float averageEloss = eloss/NumberOfSegs;
744  for (int ii=0;ii<NumberOfSegs;ii++) elossVector[ii]= averageEloss;
745  }
746  return;
747 }
748 
749 //*******************************************************************************
750 // Drift the charge segments to the sensor surface (collection plane)
751 // Include the effect of E-field and B-field
753  const PixelGeomDetUnit* pixdet,
754  const GlobalVector& bfield,
755  const std::vector<EnergyDepositUnit>& ionization_points,
756  std::vector<SignalPoint>& collection_points) const {
757 
758 #ifdef TP_DEBUG
759  LogDebug ("Pixel Digitizer") << " enter drift " ;
760 #endif
761 
762  collection_points.resize(ionization_points.size()); // set size
763 
764  LocalVector driftDir=DriftDirection(pixdet, bfield, hit.detUnitId()); // get the charge drift direction
765  if(driftDir.z() ==0.) {
766  LogWarning("Magnetic field") << " pxlx: drift in z is zero ";
767  return;
768  }
769 
770  // tangent of Lorentz angle
771  //float TanLorenzAngleX = driftDir.x()/driftDir.z();
772  //float TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
773 
774  float TanLorenzAngleX, TanLorenzAngleY,dir_z, CosLorenzAngleX,
775  CosLorenzAngleY;
776  if( alpha2Order) {
777 
778  TanLorenzAngleX = driftDir.x(); // tangen of Lorentz angle
779  TanLorenzAngleY = driftDir.y();
780  dir_z = driftDir.z(); // The z drift direction
781  CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine
782  CosLorenzAngleY = 1./sqrt(1.+TanLorenzAngleY*TanLorenzAngleY); //cosine;
783 
784  } else{
785 
786  TanLorenzAngleX = driftDir.x();
787  TanLorenzAngleY = 0.; // force to 0, driftDir.y()/driftDir.z();
788  dir_z = driftDir.z(); // The z drift direction
789  CosLorenzAngleX = 1./sqrt(1.+TanLorenzAngleX*TanLorenzAngleX); //cosine to estimate the path length
790  CosLorenzAngleY = 1.;
791  }
792 
793  float moduleThickness = pixdet->specificSurface().bounds().thickness();
794 #ifdef TP_DEBUG
795  LogDebug ("Pixel Digitizer")
796  << " Lorentz Tan " << TanLorenzAngleX << " " << TanLorenzAngleY <<" "
797  << CosLorenzAngleX << " " << CosLorenzAngleY << " "
798  << moduleThickness*TanLorenzAngleX << " " << driftDir;
799 #endif
800 
801  float Sigma_x = 1.; // Charge spread
802  float Sigma_y = 1.;
803  float DriftDistance; // Distance between charge generation and collection
804  float DriftLength; // Actual Drift Lentgh
805  float Sigma;
806 
807 
808  for (unsigned int i = 0; i != ionization_points.size(); i++) {
809 
810  float SegX, SegY, SegZ; // position
811  SegX = ionization_points[i].x();
812  SegY = ionization_points[i].y();
813  SegZ = ionization_points[i].z();
814 
815  // Distance from the collection plane
816  //DriftDistance = (moduleThickness/2. + SegZ); // Drift to -z
817  // Include explixitely the E drift direction (for CMS dir_z=-1)
818  DriftDistance = moduleThickness/2. - (dir_z * SegZ); // Drift to -z
819 
820  //if( DriftDistance <= 0.)
821  //cout<<" <=0 "<<DriftDistance<<" "<<i<<" "<<SegZ<<" "<<dir_z<<" "
822  // <<SegX<<" "<<SegY<<" "<<(moduleThickness/2)<<" "
823  // <<ionization_points[i].energy()<<" "
824  // <<hit.particleType()<<" "<<hit.pabs()<<" "<<hit.energyLoss()<<" "
825  // <<hit.entryPoint()<<" "<<hit.exitPoint()
826  // <<std::endl;
827 
828  if( DriftDistance < 0.) {
829  DriftDistance = 0.;
830  } else if( DriftDistance > moduleThickness )
831  DriftDistance = moduleThickness;
832 
833  // Assume full depletion now, partial depletion will come later.
834  float XDriftDueToMagField = DriftDistance * TanLorenzAngleX;
835  float YDriftDueToMagField = DriftDistance * TanLorenzAngleY;
836 
837  // Shift cloud center
838  float CloudCenterX = SegX + XDriftDueToMagField;
839  float CloudCenterY = SegY + YDriftDueToMagField;
840 
841  // Calculate how long is the charge drift path
842  DriftLength = sqrt( DriftDistance*DriftDistance +
843  XDriftDueToMagField*XDriftDueToMagField +
844  YDriftDueToMagField*YDriftDueToMagField );
845 
846  // What is the charge diffusion after this path
847  Sigma = sqrt(DriftLength/Dist300) * Sigma0;
848 
849  // Project the diffusion sigma on the collection plane
850  Sigma_x = Sigma / CosLorenzAngleX ;
851  Sigma_y = Sigma / CosLorenzAngleY ;
852 
853  // Insert a charge loss due to Rad Damage here
854  float energyOnCollector = ionization_points[i].energy(); // The energy that reaches the collector
855  // pseudoRadDamage
856  if (pseudoRadDamage>=0.001){
857  float moduleRadius = pixdet->surface().position().perp();
858  if (moduleRadius<=pseudoRadDamageRadius){
859  float kValue = pseudoRadDamage /(moduleRadius*moduleRadius);
860  //std::cout <<"\nmoduleRadius= "<<moduleRadius<<" WITH K="<<kValue <<" wich scales the energy by "<<
861  // exp( -1*kValue*DriftDistance/moduleThickness )<<"\n" ;
862  energyOnCollector = energyOnCollector * exp( -1*kValue*DriftDistance/moduleThickness );
863  }
864  }
865 #ifdef TP_DEBUG
866  LogDebug ("Pixel Digitizer")
867  <<" Dift DistanceZ= "<<DriftDistance<<" module thickness= "<<moduleThickness
868  <<" Start Energy= "<<ionization_points[i].energy()<<" Energy after loss= "<<energyOnCollector;
869 #endif
870  SignalPoint sp( CloudCenterX, CloudCenterY,
871  Sigma_x, Sigma_y, hit.tof(), energyOnCollector );
872 
873  // Load the Charge distribution parameters
874  collection_points[i] = (sp);
875 
876  } // loop over ionization points, i.
877 
878 } // end drift
879 
880 //*************************************************************************
881 // Induce the signal on the collection plane of the active sensor area.
883  const size_t hitIndex,
884  const unsigned int tofBin,
885  const PixelGeomDetUnit* pixdet,
886  const std::vector<SignalPoint>& collection_points) {
887 
888  // X - Rows, Left-Right, 160, (1.6cm) for barrel
889  // Y - Columns, Down-Up, 416, (6.4cm)
890 
891  const PixelTopology* topol=&pixdet->specificTopology();
892  uint32_t detID= pixdet->geographicalId().rawId();
893  signal_map_type& theSignal = _signal[detID];
894 
895 #ifdef TP_DEBUG
896  LogDebug ("Pixel Digitizer")
897  << " enter induce_signal, "
898  << topol->pitch().first << " " << topol->pitch().second; //OK
899 #endif
900 
901  // local map to store pixels hit by 1 Hit.
902  typedef std::map< int, float, std::less<int> > hit_map_type;
903  hit_map_type hit_signal;
904 
905  // map to store pixel integrals in the x and in the y directions
906  std::map<int, float, std::less<int> > x,y;
907 
908  // Assign signals to readout channels and store sorted by channel number
909 
910  // Iterate over collection points on the collection plane
911  for ( std::vector<SignalPoint>::const_iterator i=collection_points.begin();
912  i != collection_points.end(); ++i) {
913 
914  float CloudCenterX = i->position().x(); // Charge position in x
915  float CloudCenterY = i->position().y(); // in y
916  float SigmaX = i->sigma_x(); // Charge spread in x
917  float SigmaY = i->sigma_y(); // in y
918  float Charge = i->amplitude(); // Charge amplitude
919 
920 
921  //if(SigmaX==0 || SigmaY==0) {
922  //cout<<SigmaX<<" "<<SigmaY
923  // << " cloud " << i->position().x() << " " << i->position().y() << " "
924  // << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude()<<std::endl;
925  //}
926 
927 #ifdef TP_DEBUG
928  LogDebug ("Pixel Digitizer")
929  << " cloud " << i->position().x() << " " << i->position().y() << " "
930  << i->sigma_x() << " " << i->sigma_y() << " " << i->amplitude();
931 #endif
932 
933  // Find the maximum cloud spread in 2D plane , assume 3*sigma
934  float CloudRight = CloudCenterX + ClusterWidth*SigmaX;
935  float CloudLeft = CloudCenterX - ClusterWidth*SigmaX;
936  float CloudUp = CloudCenterY + ClusterWidth*SigmaY;
937  float CloudDown = CloudCenterY - ClusterWidth*SigmaY;
938 
939  // Define 2D cloud limit points
940  LocalPoint PointRightUp = LocalPoint(CloudRight,CloudUp);
941  LocalPoint PointLeftDown = LocalPoint(CloudLeft,CloudDown);
942 
943  // This points can be located outside the sensor area.
944  // The conversion to measurement point does not check for that
945  // so the returned pixel index might be wrong (outside range).
946  // We rely on the limits check below to fix this.
947  // But remember whatever we do here THE CHARGE OUTSIDE THE ACTIVE
948  // PIXEL AREA IS LOST, it should not be collected.
949 
950  // Convert the 2D points to pixel indices
951  MeasurementPoint mp = topol->measurementPosition(PointRightUp ); //OK
952 
953  int IPixRightUpX = int( floor( mp.x()));
954  int IPixRightUpY = int( floor( mp.y()));
955 
956 #ifdef TP_DEBUG
957  LogDebug ("Pixel Digitizer") << " right-up " << PointRightUp << " "
958  << mp.x() << " " << mp.y() << " "
959  << IPixRightUpX << " " << IPixRightUpY ;
960 #endif
961 
962  mp = topol->measurementPosition(PointLeftDown ); //OK
963 
964  int IPixLeftDownX = int( floor( mp.x()));
965  int IPixLeftDownY = int( floor( mp.y()));
966 
967 #ifdef TP_DEBUG
968  LogDebug ("Pixel Digitizer") << " left-down " << PointLeftDown << " "
969  << mp.x() << " " << mp.y() << " "
970  << IPixLeftDownX << " " << IPixLeftDownY ;
971 #endif
972 
973  // Check detector limits to correct for pixels outside range.
974  int numColumns = topol->ncolumns(); // det module number of cols&rows
975  int numRows = topol->nrows();
976 
977  IPixRightUpX = numRows>IPixRightUpX ? IPixRightUpX : numRows-1 ;
978  IPixRightUpY = numColumns>IPixRightUpY ? IPixRightUpY : numColumns-1 ;
979  IPixLeftDownX = 0<IPixLeftDownX ? IPixLeftDownX : 0 ;
980  IPixLeftDownY = 0<IPixLeftDownY ? IPixLeftDownY : 0 ;
981 
982  x.clear(); // clear temporary integration array
983  y.clear();
984 
985  // First integrate charge strips in x
986  int ix; // TT for compatibility
987  for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) { // loop over x index
988  float xUB, xLB, UpperBound, LowerBound;
989 
990  // Why is set to 0 if ix=0, does it meen that we accept charge
991  // outside the sensor? CHeck How it was done in ORCA?
992  //if(ix == 0) LowerBound = 0.;
993  if(ix == 0 || SigmaX==0. ) // skip for surface segemnts
994  LowerBound = 0.;
995  else {
996  mp = MeasurementPoint( float(ix), 0.0);
997  xLB = topol->localPosition(mp).x();
998  LowerBound = 1-calcQ((xLB-CloudCenterX)/SigmaX);
999  }
1000 
1001  if(ix == numRows-1 || SigmaX==0. )
1002  UpperBound = 1.;
1003  else {
1004  mp = MeasurementPoint( float(ix+1), 0.0);
1005  xUB = topol->localPosition(mp).x();
1006  UpperBound = 1. - calcQ((xUB-CloudCenterX)/SigmaX);
1007  }
1008 
1009  float TotalIntegrationRange = UpperBound - LowerBound; // get strip
1010  x[ix] = TotalIntegrationRange; // save strip integral
1011  //if(SigmaX==0 || SigmaY==0)
1012  //cout<<TotalIntegrationRange<<" "<<ix<<std::endl;
1013 
1014  }
1015 
1016  // Now integrate strips in y
1017  int iy; // TT for compatibility
1018  for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
1019  float yUB, yLB, UpperBound, LowerBound;
1020 
1021  if(iy == 0 || SigmaY==0.)
1022  LowerBound = 0.;
1023  else {
1024  mp = MeasurementPoint( 0.0, float(iy) );
1025  yLB = topol->localPosition(mp).y();
1026  LowerBound = 1. - calcQ((yLB-CloudCenterY)/SigmaY);
1027  }
1028 
1029  if(iy == numColumns-1 || SigmaY==0. )
1030  UpperBound = 1.;
1031  else {
1032  mp = MeasurementPoint( 0.0, float(iy+1) );
1033  yUB = topol->localPosition(mp).y();
1034  UpperBound = 1. - calcQ((yUB-CloudCenterY)/SigmaY);
1035  }
1036 
1037  float TotalIntegrationRange = UpperBound - LowerBound;
1038  y[iy] = TotalIntegrationRange; // save strip integral
1039  //if(SigmaX==0 || SigmaY==0)
1040  //cout<<TotalIntegrationRange<<" "<<iy<<std::endl;
1041  }
1042 
1043  // Get the 2D charge integrals by folding x and y strips
1044  int chan;
1045  for (ix=IPixLeftDownX; ix<=IPixRightUpX; ix++) { // loop over x index
1046  for (iy=IPixLeftDownY; iy<=IPixRightUpY; iy++) { //loope over y ind
1047 
1048  float ChargeFraction = Charge*x[ix]*y[iy];
1049 
1050  if( ChargeFraction > 0. ) {
1051  chan = PixelDigi::pixelToChannel( ix, iy); // Get index
1052  // Load the amplitude
1053  hit_signal[chan] += ChargeFraction;
1054  } // endif
1055 
1056 
1057  mp = MeasurementPoint( float(ix), float(iy) );
1058  LocalPoint lp = topol->localPosition(mp);
1059  chan = topol->channel(lp);
1060 
1061 #ifdef TP_DEBUG
1062  LogDebug ("Pixel Digitizer")
1063  << " pixel " << ix << " " << iy << " - "<<" "
1064  << chan << " " << ChargeFraction<<" "
1065  << mp.x() << " " << mp.y() <<" "
1066  << lp.x() << " " << lp.y() << " " // givex edge position
1067  << chan; // edge belongs to previous ?
1068 #endif
1069 
1070  } // endfor iy
1071  } //endfor ix
1072 
1073 
1074  // Test conversions (THIS IS FOR TESTING ONLY) comment-out.
1075  // mp = topol->measurementPosition( i->position() ); //OK
1076  // LocalPoint lp = topol->localPosition(mp); //OK
1077  // std::pair<float,float> p = topol->pixel( i->position() ); //OK
1078  // chan = PixelDigi::pixelToChannel( int(p.first), int(p.second));
1079  // std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1080  // MeasurementPoint mp1 = MeasurementPoint( float(ip.first),
1081  // float(ip.second) );
1082  // LogDebug ("Pixel Digitizer") << " Test "<< mp.x() << " " << mp.y()
1083  // << " "<< lp.x() << " " << lp.y() << " "<<" "
1084  // <<p.first <<" "<<p.second<<" "<<chan<< " "
1085  // <<" " << ip.first << " " << ip.second << " "
1086  // << mp1.x() << " " << mp1.y() << " " //OK
1087  // << topol->localPosition(mp1).x() << " " //OK
1088  // << topol->localPosition(mp1).y() << " "
1089  // << topol->channel( i->position() ); //OK
1090 
1091 
1092  } // loop over charge distributions
1093 
1094  // Fill the global map with all hit pixels from this event
1095 
1096  for ( hit_map_type::const_iterator im = hit_signal.begin();
1097  im != hit_signal.end(); ++im) {
1098  int chan = (*im).first;
1099  theSignal[chan] += (makeDigiSimLinks_ ? Amplitude( (*im).second, &hit, hitIndex, tofBin, (*im).second) : Amplitude( (*im).second, (*im).second) ) ;
1100 
1101 #ifdef TP_DEBUG
1102  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1103  LogDebug ("Pixel Digitizer")
1104  << " pixel " << ip.first << " " << ip.second << " "
1105  << theSignal[chan];
1106 #endif
1107  }
1108 
1109 } // end induce_signal
1110 
1111 /***********************************************************************/
1112 
1113 // Build pixels, check threshold, add misscalibration, ...
1114 void SiPixelDigitizerAlgorithm::make_digis(float thePixelThresholdInE,
1115  uint32_t detID,
1116  std::vector<PixelDigi>& digis,
1117  std::vector<PixelDigiSimLink>& simlinks,
1118  const TrackerTopology *tTopo) const {
1119 
1120 #ifdef TP_DEBUG
1121  LogDebug ("Pixel Digitizer") << " make digis "<<" "
1122  << " pixel threshold FPix" << theThresholdInE_FPix << " "
1123  << " pixel threshold BPix" << theThresholdInE_BPix << " "
1124  << " pixel threshold BPix Layer1" << theThresholdInE_BPix_L1 << " "
1125  << " List pixels passing threshold ";
1126 #endif
1127 
1128  // Loop over hit pixels
1129 
1130  signalMaps::const_iterator it = _signal.find(detID);
1131  if (it == _signal.end()) {
1132  return;
1133  }
1134 
1135  const signal_map_type& theSignal = (*it).second;
1136 
1137  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1138 
1139  float signalInElectrons = (*i).second ; // signal in electrons
1140 
1141  // Do the miss calibration for calibration studies only.
1142  //if(doMissCalibrate) signalInElectrons = missCalibrate(signalInElectrons)
1143 
1144  // Do only for pixels above threshold
1145 
1146  if( signalInElectrons >= thePixelThresholdInE) { // check threshold
1147 
1148  int chan = (*i).first; // channel number
1149  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1150  int adc=0; // ADC count as integer
1151 
1152  // Do the miss calibration for calibration studies only.
1153  if(doMissCalibrate) {
1154  int row = ip.first; // X in row
1155  int col = ip.second; // Y is in col
1156  adc = int(missCalibrate(detID, col, row, signalInElectrons)); //full misscalib.
1157  } else { // Just do a simple electron->adc conversion
1158  adc = int( signalInElectrons / theElectronPerADC ); // calibrate gain
1159  }
1160  adc = std::min(adc, theAdcFullScale); // Check maximum value
1161 // Calculate layerIndex
1163  unsigned int Sub_detid=DetId(detID).subdetId();
1164  if(Sub_detid == PixelSubdetector::PixelBarrel){ // Barrel modules
1165  int lay = tTopo->pxbLayer(detID);
1166  if (lay>=theFirstStackLayer) {
1167  // Set to 1 if over the threshold
1168  if (theAdcFullScaleStack==1) {adc=1;}
1169  // Make it a linear fit to the full scale of the normal adc count. Start new adc from 1 not zero.
1171  }
1172  }
1173  } // Only enter this if the Adc changes for the outer layers
1174 #ifdef TP_DEBUG
1175  LogDebug ("Pixel Digitizer")
1176  << (*i).first << " " << (*i).second << " " << signalInElectrons
1177  << " " << adc << ip.first << " " << ip.second ;
1178 #endif
1179 
1180  // Load digis
1181  digis.emplace_back(ip.first, ip.second, adc);
1182 
1183  if (makeDigiSimLinks_ && (*i).second.hitInfo()!=0) {
1184  //digilink
1185  if((*i).second.trackIds().size()>0){
1186  simlink_map simi;
1187  unsigned int il=0;
1188  for( std::vector<unsigned int>::const_iterator itid = (*i).second.trackIds().begin();
1189  itid != (*i).second.trackIds().end(); ++itid) {
1190  simi[*itid].push_back((*i).second.individualampl()[il]);
1191  il++;
1192  }
1193 
1194  //sum the contribution of the same trackid
1195  for( simlink_map::iterator simiiter=simi.begin();
1196  simiiter!=simi.end();
1197  simiiter++){
1198 
1199  float sum_samechannel=0;
1200  for (unsigned int iii=0;iii<(*simiiter).second.size();iii++){
1201  sum_samechannel+=(*simiiter).second[iii];
1202  }
1203  float fraction=sum_samechannel/(*i).second;
1204  if(fraction>1.) fraction=1.;
1205  simlinks.emplace_back((*i).first, (*simiiter).first, (*i).second.hitIndex(), (*i).second.tofBin(), (*i).second.eventId(), fraction);
1206  }
1207  }
1208  }
1209  }
1210  }
1211 }
1212 
1213 /***********************************************************************/
1214 
1215 // Add electronic noise to pixel charge
1217  float thePixelThreshold,
1218  CLHEP::HepRandomEngine* engine) {
1219 
1220 #ifdef TP_DEBUG
1221  LogDebug ("Pixel Digitizer") << " enter add_noise " << theNoiseInElectrons;
1222 #endif
1223 
1224  uint32_t detID= pixdet->geographicalId().rawId();
1225  signal_map_type& theSignal = _signal[detID];
1226 
1227 
1228  // First add noise to hit pixels
1229  float theSmearedChargeRMS = 0.0;
1230 
1231  for ( signal_map_iterator i = theSignal.begin(); i != theSignal.end(); i++) {
1232 
1234  {
1235  if((*i).second < 3000)
1236  {
1237  theSmearedChargeRMS = 543.6 - (*i).second * 0.093;
1238  } else if((*i).second < 6000){
1239  theSmearedChargeRMS = 307.6 - (*i).second * 0.01;
1240  } else{
1241  theSmearedChargeRMS = -432.4 +(*i).second * 0.123;
1242  }
1243 
1244  // Noise from Vcal smearing:
1245  float noise_ChargeVCALSmearing = theSmearedChargeRMS * CLHEP::RandGaussQ::shoot(engine, 0., 1.);
1246  // Noise from full readout:
1247  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1248 
1249  if(((*i).second + Amplitude(noise+noise_ChargeVCALSmearing, -1.)) < 0. ) {
1250  (*i).second.set(0);}
1251  else{
1252  (*i).second +=Amplitude(noise+noise_ChargeVCALSmearing, -1.);
1253  }
1254 
1255  } // End if addChargeVCalSmearing
1256  else
1257  {
1258  // Noise: ONLY full READOUT Noise.
1259  // Use here the FULL readout noise, including TBM,ALT,AOH,OPT-REC.
1260  float noise = CLHEP::RandGaussQ::shoot(engine, 0., theReadoutNoise);
1261 
1262  if(((*i).second + Amplitude(noise, -1.)) < 0. ) {
1263  (*i).second.set(0);}
1264  else{
1265  (*i).second +=Amplitude(noise, -1.);
1266  }
1267  } // end if only Noise from full readout
1268 
1269  }
1270 
1271  if(!addNoisyPixels) // Option to skip noise in non-hit pixels
1272  return;
1273 
1274  const PixelTopology* topol=&pixdet->specificTopology();
1275  int numColumns = topol->ncolumns(); // det module number of cols&rows
1276  int numRows = topol->nrows();
1277 
1278  // Add noise on non-hit pixels
1279  // Use here the pixel noise
1280  int numberOfPixels = (numRows * numColumns);
1281  std::map<int,float, std::less<int> > otherPixels;
1282  std::map<int,float, std::less<int> >::iterator mapI;
1283 
1284  theNoiser->generate(numberOfPixels,
1285  thePixelThreshold, //thr. in un. of nois
1286  theNoiseInElectrons, // noise in elec.
1287  otherPixels,
1288  engine );
1289 
1290 #ifdef TP_DEBUG
1291  LogDebug ("Pixel Digitizer")
1292  << " Add noisy pixels " << numRows << " "
1293  << numColumns << " " << theNoiseInElectrons << " "
1294  << theThresholdInE_FPix << theThresholdInE_BPix <<" "<< numberOfPixels<<" "
1295  << otherPixels.size() ;
1296 #endif
1297 
1298  // Add noisy pixels
1299  for (mapI = otherPixels.begin(); mapI!= otherPixels.end(); mapI++) {
1300  int iy = ((*mapI).first) / numRows;
1301  int ix = ((*mapI).first) - (iy*numRows);
1302 
1303  // Keep for a while for testing.
1304  if( iy < 0 || iy > (numColumns-1) )
1305  LogWarning ("Pixel Geometry") << " error in iy " << iy ;
1306  if( ix < 0 || ix > (numRows-1) )
1307  LogWarning ("Pixel Geometry") << " error in ix " << ix ;
1308 
1309  int chan = PixelDigi::pixelToChannel(ix, iy);
1310 
1311 #ifdef TP_DEBUG
1312  LogDebug ("Pixel Digitizer")
1313  <<" Storing noise = " << (*mapI).first << " " << (*mapI).second
1314  << " " << ix << " " << iy << " " << chan ;
1315 #endif
1316 
1317  if(theSignal[chan] == 0){
1318  // float noise = float( (*mapI).second );
1319  int noise=int( (*mapI).second );
1320  theSignal[chan] = Amplitude (noise, -1.);
1321  }
1322  }
1323 }
1324 
1325 /***********************************************************************/
1326 
1327 // Simulate the readout inefficiencies.
1328 // Delete a selected number of single pixels, dcols and rocs.
1330  const PixelGeomDetUnit* pixdet,
1331  const TrackerTopology *tTopo,
1332  CLHEP::HepRandomEngine* engine) {
1333 
1334  uint32_t detID= pixdet->geographicalId().rawId();
1335 
1336  signal_map_type& theSignal = _signal[detID];
1337 
1338  const PixelTopology* topol=&pixdet->specificTopology();
1339  int numColumns = topol->ncolumns(); // det module number of cols&rows
1340  int numRows = topol->nrows();
1341 
1342  // Predefined efficiencies
1343  float pixelEfficiency = 1.0;
1344  float columnEfficiency = 1.0;
1345  float chipEfficiency = 1.0;
1346 
1347  // setup the chip indices conversion
1348  unsigned int Subid=DetId(detID).subdetId();
1349  if (Subid== PixelSubdetector::PixelBarrel){// barrel layers
1350  int layerIndex=tTopo->pxbLayer(detID);
1351  pixelEfficiency = eff.thePixelEfficiency[layerIndex-1];
1352  columnEfficiency = eff.thePixelColEfficiency[layerIndex-1];
1353  chipEfficiency = eff.thePixelChipEfficiency[layerIndex-1];
1354  //std::cout <<"Using BPix columnEfficiency = "<<columnEfficiency<< " for layer = "<<layerIndex <<"\n";
1355  // This should never happen, but only check if it is not an upgrade geometry
1356  if (NumberOfBarrelLayers==3){
1357  if(numColumns>416) LogWarning ("Pixel Geometry") <<" wrong columns in barrel "<<numColumns;
1358  if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in barrel "<<numRows;
1359 
1360  int ladder=tTopo->pxbLadder(detID);
1361  int module=tTopo->pxbModule(detID);
1362  if (module<=4) module=5-module;
1363  else module-=4;
1364 
1365  columnEfficiency *= eff.theLadderEfficiency_BPix[layerIndex-1][ladder-1]*eff.theModuleEfficiency_BPix[layerIndex-1][module-1]*_pu_scale[layerIndex-1];
1366  }
1367  } else { // forward disks
1368  unsigned int diskIndex=tTopo->pxfDisk(detID)+eff.FPixIndex; // Use diskIndex-1 later to stay consistent with BPix
1369  //if (eff.FPixIndex>diskIndex-1){throw cms::Exception("Configuration") <<"SiPixelDigitizer is using the wrong efficiency value. index = "
1370  // <<diskIndex-1<<" , MinIndex = "<<eff.FPixIndex<<" ... "<<tTopo->pxfDisk(detID);}
1371  pixelEfficiency = eff.thePixelEfficiency[diskIndex-1];
1372  columnEfficiency = eff.thePixelColEfficiency[diskIndex-1];
1373  chipEfficiency = eff.thePixelChipEfficiency[diskIndex-1];
1374  //std::cout <<"Using FPix columnEfficiency = "<<columnEfficiency<<" for Disk = "<< tTopo->pxfDisk(detID)<<"\n";
1375  // Sometimes the forward pixels have wrong size,
1376  // this crashes the index conversion, so exit, but only check if it is not an upgrade geometry
1377  if (NumberOfBarrelLayers==3){
1378  if(numColumns>260 || numRows>160) {
1379  if(numColumns>260) LogWarning ("Pixel Geometry") <<" wrong columns in endcaps "<<numColumns;
1380  if(numRows>160) LogWarning ("Pixel Geometry") <<" wrong rows in endcaps "<<numRows;
1381  return;
1382  }
1383  }
1384  } // if barrel/forward
1385 
1386 #ifdef TP_DEBUG
1387  LogDebug ("Pixel Digitizer") << " enter pixel_inefficiency " << pixelEfficiency << " "
1388  << columnEfficiency << " " << chipEfficiency;
1389 #endif
1390 
1391  // Initilize the index converter
1392  //PixelIndices indexConverter(numColumns,numRows);
1393  std::auto_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
1394 
1395  int chipIndex = 0;
1396  int rowROC = 0;
1397  int colROC = 0;
1398  std::map<int, int, std::less<int> >chips, columns;
1399  std::map<int, int, std::less<int> >::iterator iter;
1400 
1401  // Find out the number of columns and rocs hits
1402  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1403  for (signal_map_const_iterator i = theSignal.begin(); i != theSignal.end(); ++i) {
1404 
1405  int chan = i->first;
1406  std::pair<int,int> ip = PixelDigi::channelToPixel(chan);
1407  int row = ip.first; // X in row
1408  int col = ip.second; // Y is in col
1409  //transform to ROC index coordinates
1410  pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1411  int dColInChip = pIndexConverter->DColumn(colROC); // get ROC dcol from ROC col
1412  //dcol in mod
1413  int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
1414 
1415  chips[chipIndex]++;
1416  columns[dColInDet]++;
1417  }
1418 
1419  // Delete some ROC hits.
1420  for ( iter = chips.begin(); iter != chips.end() ; iter++ ) {
1421  //float rand = RandFlat::shoot();
1422  float rand = CLHEP::RandFlat::shoot(engine);
1423  if( rand > chipEfficiency ) chips[iter->first]=0;
1424  }
1425 
1426  // Delete some Dcol hits.
1427  for ( iter = columns.begin(); iter != columns.end() ; iter++ ) {
1428  //float rand = RandFlat::shoot();
1429  float rand = CLHEP::RandFlat::shoot(engine);
1430  if( rand > columnEfficiency ) columns[iter->first]=0;
1431  }
1432 
1433  // Now loop again over pixels to kill some of them.
1434  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1435  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1436 
1437  // int chan = i->first;
1438  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1439  int row = ip.first; // X in row
1440  int col = ip.second; // Y is in col
1441  //transform to ROC index coordinates
1442  pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1443  int dColInChip = pIndexConverter->DColumn(colROC); //get ROC dcol from ROC col
1444  //dcol in mod
1445  int dColInDet = pIndexConverter->DColumnInModule(dColInChip,chipIndex);
1446 
1447 
1448  //float rand = RandFlat::shoot();
1449  float rand = CLHEP::RandFlat::shoot(engine);
1450  if( chips[chipIndex]==0 || columns[dColInDet]==0
1451  || rand>pixelEfficiency ) {
1452  // make pixel amplitude =0, pixel will be lost at clusterization
1453  i->second.set(0.); // reset amplitude,
1454  } // end if
1455 
1456  } // end pixel loop
1457 } // end pixel_indefficiency
1458 
1459 //***********************************************************************
1460 
1461 // Fluctuate the gain and offset for the amplitude calibration
1462 // Use gaussian smearing.
1463 //float SiPixelDigitizerAlgorithm::missCalibrate(const float amp) const {
1464  //float gain = RandGaussQ::shoot(1.,theGainSmearing);
1465  //float offset = RandGaussQ::shoot(0.,theOffsetSmearing);
1466  //float newAmp = amp * gain + offset;
1467  // More complex misscalibration
1468 float SiPixelDigitizerAlgorithm::missCalibrate(uint32_t detID, int col,int row,
1469  const float signalInElectrons) const {
1470  // Central values
1471  //const float p0=0.00352, p1=0.868, p2=112., p3=113.; // pix(0,0,0)
1472  // const float p0=0.00382, p1=0.886, p2=112.7, p3=113.0; // average roc=0
1473  //const float p0=0.00492, p1=1.998, p2=90.6, p3=134.1; // average roc=6
1474  // Smeared (rms)
1475  //const float s0=0.00020, s1=0.051, s2=5.4, s3=4.4; // average roc=0
1476  //const float s0=0.00015, s1=0.043, s2=3.2, s3=3.1; // col average roc=0
1477 
1478  // Make 2 sets of parameters for Fpix and BPIx:
1479 
1480  float p0=0.0;
1481  float p1=0.0;
1482  float p2=0.0;
1483  float p3=0.0;
1484 
1485  unsigned int Sub_detid=DetId(detID).subdetId();
1486 
1487  if(Sub_detid == PixelSubdetector::PixelBarrel){// barrel layers
1488  p0 = BPix_p0;
1489  p1 = BPix_p1;
1490  p2 = BPix_p2;
1491  p3 = BPix_p3;
1492  } else {// forward disks
1493  p0 = FPix_p0;
1494  p1 = FPix_p1;
1495  p2 = FPix_p2;
1496  p3 = FPix_p3;
1497  }
1498 
1499  // const float electronsPerVCAL = 65.5; // our present VCAL calibration (feb 2009)
1500  // const float electronsPerVCAL_Offset = -414.0; // our present VCAL calibration (feb 2009)
1501  float newAmp = 0.; //Modified signal
1502 
1503  // Convert electrons to VCAL units
1504  float signal = (signalInElectrons-electronsPerVCAL_Offset)/electronsPerVCAL;
1505 
1506  // Simulate the analog response with fixed parametrization
1507  newAmp = p3 + p2 * tanh(p0*signal - p1);
1508 
1509 
1510  // Use the pixel-by-pixel calibrations
1511  //transform to ROC index coordinates
1512  //int chipIndex=0, colROC=0, rowROC=0;
1513  //std::auto_ptr<PixelIndices> pIndexConverter(new PixelIndices(numColumns,numRows));
1514  //pIndexConverter->transformToROC(col,row,chipIndex,colROC,rowROC);
1515 
1516  // Use calibration from a file
1517  //int chanROC = PixelIndices::pixelToChannelROC(rowROC,colROC); // use ROC coordinates
1518  //float pp0=0, pp1=0,pp2=0,pp3=0;
1519  //map<int,CalParameters,std::less<int> >::const_iterator it=calmap.find(chanROC);
1520  //CalParameters y = (*it).second;
1521  //pp0 = y.p0;
1522  //pp1 = y.p1;
1523  //pp2 = y.p2;
1524  //pp3 = y.p3;
1525 
1526  //
1527  // Use random smearing
1528  // Randomize the pixel response
1529  //float pp0 = RandGaussQ::shoot(p0,s0);
1530  //float pp1 = RandGaussQ::shoot(p1,s1);
1531  //float pp2 = RandGaussQ::shoot(p2,s2);
1532  //float pp3 = RandGaussQ::shoot(p3,s3);
1533 
1534  //newAmp = pp3 + pp2 * tanh(pp0*signal - pp1); // Final signal
1535 
1536  //cout<<" misscalibrate "<<col<<" "<<row<<" "<<chipIndex<<" "<<colROC<<" "
1537  // <<rowROC<<" "<<signalInElectrons<<" "<<signal<<" "<<newAmp<<" "
1538  // <<(signalInElectrons/theElectronPerADC)<<std::endl;
1539 
1540  return newAmp;
1541 }
1542 //******************************************************************************
1543 
1544 // Set the drift direction accoring to the Bfield in local det-unit frame
1545 // Works for both barrel and forward pixels.
1546 // Replace the sign convention to fit M.Swartz's formulaes.
1547 // Configurations for barrel and foward pixels possess different tanLorentzAngleperTesla
1548 // parameter value
1549 
1551  const GlobalVector& bfield,
1552  const DetId& detId) const {
1553  Frame detFrame(pixdet->surface().position(),pixdet->surface().rotation());
1554  LocalVector Bfield=detFrame.toLocal(bfield);
1555 
1556  float alpha2_FPix;
1557  float alpha2_BPix;
1558  float alpha2;
1559 
1560  //float dir_x = -tanLorentzAnglePerTesla * Bfield.y();
1561  //float dir_y = +tanLorentzAnglePerTesla * Bfield.x();
1562  //float dir_z = -1.; // E field always in z direction, so electrons go to -z
1563  // The dir_z has to be +/- 1. !
1564  // LocalVector theDriftDirection = LocalVector(dir_x,dir_y,dir_z);
1565 
1566  float dir_x = 0.0;
1567  float dir_y = 0.0;
1568  float dir_z = 0.0;
1569  float scale = 0.0;
1570 
1571  uint32_t detID= pixdet->geographicalId().rawId();
1572 
1573  unsigned int Sub_detid=DetId(detID).subdetId();
1574 
1575  // Read Lorentz angle from cfg file:**************************************************************
1576 
1577  if(!use_LorentzAngle_DB_){
1578 
1579  if( alpha2Order) {
1582  }else {
1583  alpha2_FPix = 0.0;
1584  alpha2_BPix = 0.0;
1585  }
1586 
1587  if(Sub_detid == PixelSubdetector::PixelBarrel){// barrel layers
1588  dir_x = -( tanLorentzAnglePerTesla_BPix * Bfield.y() + alpha2_BPix* Bfield.z()* Bfield.x() );
1589  dir_y = +( tanLorentzAnglePerTesla_BPix * Bfield.x() - alpha2_BPix* Bfield.z()* Bfield.y() );
1590  dir_z = -(1 + alpha2_BPix* Bfield.z()*Bfield.z() );
1591  scale = (1 + alpha2_BPix* Bfield.z()*Bfield.z() );
1592 
1593  } else {// forward disks
1594  dir_x = -( tanLorentzAnglePerTesla_FPix * Bfield.y() + alpha2_FPix* Bfield.z()* Bfield.x() );
1595  dir_y = +( tanLorentzAnglePerTesla_FPix * Bfield.x() - alpha2_FPix* Bfield.z()* Bfield.y() );
1596  dir_z = -(1 + alpha2_FPix* Bfield.z()*Bfield.z() );
1597  scale = (1 + alpha2_FPix* Bfield.z()*Bfield.z() );
1598  }
1599  } // end: Read LA from cfg file.
1600 
1601  //Read Lorentz angle from DB:********************************************************************
1603  float lorentzAngle = SiPixelLorentzAngle_->getLorentzAngle(detId);
1604  alpha2 = lorentzAngle * lorentzAngle;
1605  //std::cout << "detID is: "<< it->first <<"The LA per tesla is: "<< it->second << std::std::endl;
1606  dir_x = -( lorentzAngle * Bfield.y() + alpha2 * Bfield.z()* Bfield.x() );
1607  dir_y = +( lorentzAngle * Bfield.x() - alpha2 * Bfield.z()* Bfield.y() );
1608  dir_z = -(1 + alpha2 * Bfield.z()*Bfield.z() );
1609  scale = (1 + alpha2 * Bfield.z()*Bfield.z() );
1610  }// end: Read LA from DataBase.
1611 
1612  LocalVector theDriftDirection = LocalVector(dir_x/scale, dir_y/scale, dir_z/scale );
1613 
1614 #ifdef TP_DEBUG
1615  LogDebug ("Pixel Digitizer") << " The drift direction in local coordinate is "
1616  << theDriftDirection ;
1617 #endif
1618 
1619  return theDriftDirection;
1620 }
1621 
1622 //****************************************************************************************************
1623 
1625  if(!use_ineff_from_db_)
1626  return;
1627 
1628  signal_map_type& theSignal = _signal[detID];
1629 
1630  // Loop over hit pixels, amplitude in electrons, channel = coded row,col
1631  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1632 
1633  // int chan = i->first;
1634  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1635  int row = ip.first; // X in row
1636  int col = ip.second; // Y is in col
1637  //transform to ROC index coordinates
1638  if(theSiPixelGainCalibrationService_->isDead(detID, col, row)){
1639  // std::cout << "now in isdead check, row " << detID << " " << col << "," << row << std::std::endl;
1640  // make pixel amplitude =0, pixel will be lost at clusterization
1641  i->second.set(0.); // reset amplitude,
1642  } // end if
1643  } // end pixel loop
1644 } // end pixel_indefficiency
1645 
1646 
1647 //****************************************************************************************************
1648 
1650  if(!use_module_killing_)
1651  return;
1652 
1653  bool isbad=false;
1654 
1655  Parameters::const_iterator itDeadModules=DeadModules.begin();
1656 
1657  int detid = detID;
1658  for(; itDeadModules != DeadModules.end(); ++itDeadModules){
1659  int Dead_detID = itDeadModules->getParameter<int>("Dead_detID");
1660  if(detid == Dead_detID){
1661  isbad=true;
1662  break;
1663  }
1664  }
1665 
1666  if(!isbad)
1667  return;
1668 
1669  signal_map_type& theSignal = _signal[detID];
1670 
1671  std::string Module = itDeadModules->getParameter<std::string>("Module");
1672 
1673  if(Module=="whole"){
1674  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1675  i->second.set(0.); // reset amplitude
1676  }
1677  }
1678 
1679  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1680  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1681 
1682  if(Module=="tbmA" && ip.first>=80 && ip.first<=159){
1683  i->second.set(0.);
1684  }
1685 
1686  if( Module=="tbmB" && ip.first<=79){
1687  i->second.set(0.);
1688  }
1689  }
1690 }
1691 //****************************************************************************************************
1693 // Not SLHC safe for now
1694  if(!use_module_killing_)
1695  return;
1696 
1697  bool isbad=false;
1698 
1699  std::vector<SiPixelQuality::disabledModuleType>disabledModules = SiPixelBadModule_->getBadComponentList();
1700 
1702 
1703  for (size_t id=0;id<disabledModules.size();id++)
1704  {
1705  if(detID==disabledModules[id].DetID){
1706  isbad=true;
1707  badmodule = disabledModules[id];
1708  break;
1709  }
1710  }
1711 
1712  if(!isbad)
1713  return;
1714 
1715  signal_map_type& theSignal = _signal[detID];
1716 
1717  //std::cout<<"Hit in: "<< detID <<" errorType "<< badmodule.errorType<<" BadRocs="<<std::hex<<SiPixelBadModule_->getBadRocs(detID)<<dec<<" "<<std::endl;
1718  if(badmodule.errorType == 0){ // this is a whole dead module.
1719 
1720  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1721  i->second.set(0.); // reset amplitude
1722  }
1723  }
1724  else { // all other module types: half-modules and single ROCs.
1725  // Get Bad ROC position:
1726  //follow the example of getBadRocPositions in CondFormats/SiPixelObjects/src/SiPixelQuality.cc
1727  std::vector<GlobalPixel> badrocpositions (0);
1728  for(unsigned int j = 0; j < 16; j++){
1729  if(SiPixelBadModule_->IsRocBad(detID, j) == true){
1730 
1731  std::vector<CablingPathToDetUnit> path = map_.product()->pathToDetUnit(detID);
1732  typedef std::vector<CablingPathToDetUnit>::const_iterator IT;
1733  for (IT it = path.begin(); it != path.end(); ++it) {
1734  const PixelROC* myroc = map_.product()->findItem(*it);
1735  if( myroc->idInDetUnit() == j) {
1736  LocalPixel::RocRowCol local = { 39, 25}; //corresponding to center of ROC row, col
1737  GlobalPixel global = myroc->toGlobal( LocalPixel(local) );
1738  badrocpositions.push_back(global);
1739  break;
1740  }
1741  }
1742  }
1743  }// end of getBadRocPositions
1744 
1745 
1746  for(signal_map_iterator i = theSignal.begin();i != theSignal.end(); ++i) {
1747  std::pair<int,int> ip = PixelDigi::channelToPixel(i->first);//get pixel pos
1748 
1749  for(std::vector<GlobalPixel>::const_iterator it = badrocpositions.begin(); it != badrocpositions.end(); ++it){
1750  if(it->row >= 80 && ip.first >= 80 ){
1751  if((fabs(ip.second - it->col) < 26) ) {i->second.set(0.);}
1752  else if(it->row==120 && ip.second-it->col==26){i->second.set(0.);}
1753  else if(it->row==119 && it->col-ip.second==26){i->second.set(0.);}
1754  }
1755  else if(it->row < 80 && ip.first < 80 ){
1756  if((fabs(ip.second - it->col) < 26) ){i->second.set(0.);}
1757  else if(it->row==40 && ip.second-it->col==26){i->second.set(0.);}
1758  else if(it->row==39 && it->col-ip.second==26){i->second.set(0.);}
1759  }
1760  }
1761  }
1762  }
1763 }
1764 
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 accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, const size_t inputBeginGlobalIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const GlobalVector &bfield, CLHEP::HepRandomEngine *)
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:40
identify pixel inside single ROC
Definition: LocalPixel.h:7
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:68
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
virtual float thickness() const =0
void induce_signal(const PSimHit &hit, const size_t hitIndex, const unsigned int tofBin, const PixelGeomDetUnit *pixdet, const std::vector< SignalPoint > &collection_points)
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
T mag() const
Definition: PV3DBase.h:67
tuple path
else: Piece not in the list, fine.
const std::map< int, CalParameters, std::less< int > > calmap
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:77
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
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:86
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:43
GlobalPixel toGlobal(const LocalPixel &loc) const
Definition: PixelROC.h:58
void add_noise(const PixelGeomDetUnit *pixdet, float thePixelThreshold, CLHEP::HepRandomEngine *)
double p3[4]
Definition: TauolaWrapper.h:91
void pixel_inefficiency(const PixelEfficiencies &eff, const PixelGeomDetUnit *pixdet, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)