CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelGaussianSmearingRecHitConverterAlgorithm.cc
Go to the documentation of this file.
1 
13 // SiPixel Gaussian Smearing
15 
16 // Geometry
17 //#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
19 //#include "Geometry/CommonDetUnit/interface/GeomDetType.h"
23 
24 // Famos
27 
28 // STL
29 
30 // ROOT
31 #include <TFile.h>
32 #include <TH1F.h>
33 //#include <TAxis.h>
34 
35 //#define FAMOS_DEBUG
36 
37 const double microntocm = 0.0001;
38 using namespace std;
39 
41  const edm::ParameterSet& pset,
42  GeomDetType::SubDetector pixelPart)
43 :
44  pset_(pset),
45  thePixelPart(pixelPart)
46 {
47  // Switch between old (ORCA) and new (CMSSW) pixel parameterization
48  useCMSSWPixelParameterization = pset.getParameter<bool>("UseCMSSWPixelParametrization");
49 
52 
54  isForward = false;
55  thePixelResolutionFileName1 = pset_.getParameter<string>( "NewPixelBarrelResolutionFile1" );
56  thePixelResolutionFile1 = new
57  TFile( edm::FileInPath( thePixelResolutionFileName1 ).fullPath().c_str() ,"READ");
58  thePixelResolutionFileName2 = pset_.getParameter<string>( "NewPixelBarrelResolutionFile2" );
60  TFile( edm::FileInPath( thePixelResolutionFileName2 ).fullPath().c_str() ,"READ");
62  tempId = pset_.getParameter<int> ( "templateIdBarrel" );
64  throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm:")
65  <<"SiPixel Barrel Template Not Loaded Correctly!"<<endl;
66 #ifdef FAMOS_DEBUG
67  cout<<"The Barrel map size is "<<theXHistos.size()<<" and "<<theYHistos.size()<<endl;
68 #endif
69  }
70  else
72  isForward = true;
73  thePixelResolutionFileName1 = pset_.getParameter<string>( "NewPixelForwardResolutionFile" );
74  thePixelResolutionFile1 = new
75  TFile( edm::FileInPath( thePixelResolutionFileName1 ).fullPath().c_str() ,"READ");
77  tempId = pset_.getParameter<int> ( "templateIdForward" );
79  throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm:")
80  <<"SiPixel Forward Template Not Loaded Correctly!"<<endl;
81 #ifdef FAMOS_DEBUG
82  cout<<"The Forward map size is "<<theXHistos.size()<<" and "<<theYHistos.size()<<endl;
83 #endif
84  }
85  else
86  throw cms::Exception("SiPixelGaussianSmearingRecHitConverterAlgorithm :")
87  <<"Not a pixel detector"<<endl;
88 
90  thePixelResolutionFile2->Close();
92  }
93  if ( thePixelResolutionFile1) {
94  thePixelResolutionFile1->Close();
96  }
97  thePixelResolutionFile1=0;
99 }
100 
102 {
103 
104  std::map<unsigned,const SimpleHistogramGenerator*>::const_iterator it;
105  for ( it=theXHistos.begin(); it!=theXHistos.end(); ++it )
106  delete it->second;
107  for ( it=theYHistos.begin(); it!=theYHistos.end(); ++it )
108  delete it->second;
109 
110  theXHistos.clear();
111  theYHistos.clear();
112 
113 }
114 
116  const PSimHit& simHit,
117  const PixelGeomDetUnit* detUnit,
118  const double boundX,
119  const double boundY,
121 {
122 
123 #ifdef FAMOS_DEBUG
124  std::cout << " Pixel smearing in " << thePixelPart
125  << std::endl;
126 #endif
127  //
128  // at the beginning the position is the Local Point in the local pixel module reference frame
129  // same code as in PixelCPEBase
130  LocalVector localDir = simHit.momentumAtEntry().unit();
131  float locx = localDir.x();
132  float locy = localDir.y();
133  float locz = localDir.z();
134 
135  // alpha: angle with respect to local x axis in local (x,z) plane
136  float cotalpha = locx/locz;
137  if ( isFlipped( detUnit ) ) { // &&& check for FPIX !!!
138 #ifdef FAMOS_DEBUG
139  std::cout << " isFlipped " << std::endl;
140 #endif
141  }
142  // beta: angle with respect to local y axis in local (y,z) plane
143  float cotbeta = locy/locz;
144  float sign=1.;
145  if( isForward ) {
146  if( cotbeta < 0 ) sign=-1.;
147  cotbeta = sign*cotbeta;
148  }
149 
150 
151  //
152 #ifdef FAMOS_DEBUG
153  std::cout << " Local Direction " << simHit.localDirection()
154  << " cotalpha(x) = " << cotalpha
155  << " cotbeta(y) = " << cotbeta
156  << std::endl;
157 #endif
158 
159  const PixelTopology* theSpecificTopology = &(detUnit->specificType().specificTopology());
160  const RectangularPixelTopology *rectPixelTopology = static_cast<const RectangularPixelTopology*>(theSpecificTopology);
161 
162  const int nrows = theSpecificTopology->nrows();
163  const int ncolumns = theSpecificTopology->ncolumns();
164 
165  const Local3DPoint lp = simHit.localPosition();
166  //Transform local position to measurement position
167  const MeasurementPoint mp = rectPixelTopology->measurementPosition( lp );
168  float mpy = mp.y();
169  float mpx = mp.x();
170  //Get the center of the struck pixel in measurement position
171  float pixelCenterY = 0.5 + (int)mpy;
172  float pixelCenterX = 0.5 + (int)mpx;
173 #ifdef FAMOS_DEBUG
174  cout<<"Struck pixel center at pitch units x: "<<pixelCenterX<<" y: "<<pixelCenterY<<endl;
175 #endif
176 
177  const MeasurementPoint mpCenter(pixelCenterX, pixelCenterY);
178  //Transform the center of the struck pixel back into local position
179  const Local3DPoint lpCenter = rectPixelTopology->localPosition( mpCenter );
180 #ifdef FAMOS_DEBUG
181  cout<<"Struck point at cm x: "<<lp.x()<<" y: "<<lp.y()<<endl;
182  cout<<"Struck pixel center at cm x: "<<lpCenter.x()<<" y: "<<lpCenter.y()<<endl;
183  cout<<"The boundX is "<<boundX<<" boundY is "<<boundY<<endl;
184 #endif
185 
186  //Get the relative position of struck point to the center of the struck pixel
187  float xtrk = lp.x() - lpCenter.x();
188  float ytrk = lp.y() - lpCenter.y();
189  //Pixel Y, X pitch
190  const float ysize={0.015}, xsize={0.01};
191  //Variables for SiPixelTemplate input, see SiPixelTemplate reco
192  float yhit = 20. + 8.*(ytrk/ysize);
193  float xhit = 20. + 8.*(xtrk/xsize);
194  int ybin = (int)yhit;
195  int xbin = (int)xhit;
196  float yfrac= yhit - (float)ybin;
197  float xfrac= xhit - (float)xbin;
198  //Protect againt ybin, xbin being outside of range [0-39]
199  if( ybin < 0 ) ybin = 0;
200  if( ybin > 39 ) ybin = 39;
201  if( xbin < 0 ) xbin = 0;
202  if( xbin > 39 ) xbin = 39;
203 
204  //Variables for SiPixelTemplate output
205  //qBin -- normalized pixel charge deposition
206  float qbin_frac[4];
207  //Single pixel cluster projection possibility
208  float ny1_frac, ny2_frac, nx1_frac, nx2_frac;
209  bool singlex = false, singley = false;
211  templ.interpolate(tempId, cotalpha, cotbeta);
212  templ.qbin_dist(tempId, cotalpha, cotbeta, qbin_frac, ny1_frac, ny2_frac, nx1_frac, nx2_frac );
213  int nqbin;
214 
215  double xsizeProbability = random->flatShoot();
216  double ysizeProbability = random->flatShoot();
217  bool hitbigx = rectPixelTopology->isItBigPixelInX( (int)mpx );
218  bool hitbigy = rectPixelTopology->isItBigPixelInY( (int)mpy );
219 
220  if( hitbigx )
221  if( xsizeProbability < nx2_frac ) singlex = true;
222  else singlex = false;
223  else
224  if( xsizeProbability < nx1_frac ) singlex = true;
225  else singlex = false;
226 
227  if( hitbigy )
228  if( ysizeProbability < ny2_frac ) singley = true;
229  else singley = false;
230  else
231  if( ysizeProbability < ny1_frac ) singley = true;
232  else singley = false;
233 
234 
235 
236  // random multiplicity for alpha and beta
237  double qbinProbability = random->flatShoot();
238  for(int i = 0; i<4; ++i) {
239  nqbin = i;
240  if(qbinProbability < qbin_frac[i]) break;
241  }
242 
243  //Store interpolated pixel cluster profile
244  //BYSIZE, BXSIZE, const definition from SiPixelTemplate
245  float ytempl[41][BYSIZE] = {{0}}, xtempl[41][BXSIZE] = {{0}} ;
246  templ.ytemp(0, 40, ytempl);
247  templ.xtemp(0, 40, xtempl);
248 
249  std::vector<double> ytemp(BYSIZE);
250  for( int i=0; i<BYSIZE; ++i) {
251  ytemp[i]=(1.-yfrac)*ytempl[ybin][i]+yfrac*ytempl[ybin+1][i];
252  }
253 
254  std::vector<double> xtemp(BXSIZE);
255  for(int i=0; i<BXSIZE; ++i) {
256  xtemp[i]=(1.-xfrac)*xtempl[xbin][i]+xfrac*xtempl[xbin+1][i];
257  }
258 
259  //Pixel readout threshold
260  const float qThreshold = templ.s50()*2.0;
261 
262  //Cut away pixels below readout threshold
263  //For cluster lengths calculation
264  int offsetX1=0, offsetX2=0, offsetY1=0, offsetY2=0;
265  int firstY, lastY, firstX, lastX;
266  for( firstY = 0; firstY < BYSIZE; ++firstY ) {
267  bool yCluster = ytemp[firstY] > qThreshold ;
268  if( yCluster )
269  {
270  offsetY1 = BHY -firstY;
271  break;
272  }
273  }
274  for( lastY = firstY; lastY < BYSIZE; ++lastY )
275  {
276  bool yCluster = ytemp[lastY] > qThreshold ;
277  if( !yCluster )
278  {
279  lastY = lastY - 1;
280  offsetY2 = lastY - BHY;
281  break;
282  }
283  }
284 
285  for( firstX = 0; firstX < BXSIZE; ++firstX ) {
286  bool xCluster = xtemp[firstX] > qThreshold ;
287  if( xCluster ) {
288  offsetX1 = BHX - firstX;
289  break;
290  }
291  }
292  for( lastX = firstX; lastX < BXSIZE; ++ lastX ) {
293  bool xCluster = xtemp[lastX] > qThreshold ;
294  if( !xCluster ) {
295  lastX = lastX - 1;
296  offsetX2 = lastX - BHX;
297  break;
298  }
299  }
300 
301  bool edge, edgex, edgey;
302  // bool bigx, bigy;
303  unsigned int clslenx = offsetX1 + offsetX2 + 1;
304  unsigned int clsleny = offsetY1 + offsetY2 + 1;
305 
306  theClslenx = clslenx;
307  theClsleny = clsleny;
308 
309  int firstPixelInX = (int)mpx - offsetX1 ;
310  int firstPixelInY = (int)mpy - offsetY1 ;
311  int lastPixelInX = (int)mpx + offsetX2 ;
312  int lastPixelInY = (int)mpy + offsetY2 ;
313  firstPixelInX = (firstPixelInX >= 0) ? firstPixelInX : 0 ;
314  firstPixelInY = (firstPixelInY >= 0) ? firstPixelInY : 0 ;
315  lastPixelInX = (lastPixelInX < nrows ) ? lastPixelInX : nrows-1 ;
316  lastPixelInY = (lastPixelInY < ncolumns ) ? lastPixelInY : ncolumns-1;
317 
318  edgex = rectPixelTopology->isItEdgePixelInX( firstPixelInX ) || rectPixelTopology->isItEdgePixelInX( lastPixelInX );
319  edgey = rectPixelTopology->isItEdgePixelInY( firstPixelInY ) || rectPixelTopology->isItEdgePixelInY( lastPixelInY );
320  edge = edgex || edgey;
321 
322  // bigx = rectPixelTopology->isItBigPixelInX( firstPixelInX ) || rectPixelTopology->isItBigPixelInX( lastPixelInX );
323  // bigy = rectPixelTopology->isItBigPixelInY( firstPixelInY ) || rectPixelTopology->isItBigPixelInY( lastPixelInY );
324  bool hasBigPixelInX = rectPixelTopology->containsBigPixelInX( firstPixelInX, lastPixelInX );
325  bool hasBigPixelInY = rectPixelTopology->containsBigPixelInY( firstPixelInY, lastPixelInY );
326 
327  //Variables for SiPixelTemplate pixel hit error output
328  float sigmay, sigmax, sy1, sy2, sx1, sx2;
329  templ.temperrors(tempId, cotalpha, cotbeta, nqbin, sigmay, sigmax, sy1, sy2, sx1, sx2 );
330 
331  // define private mebers --> Errors
332  if( edge ) {
333  if( edgex && !edgey ) {
334  theErrorX = 23.0*microntocm;
335  theErrorY = 39.0*microntocm;
336  }
337  else if( !edgex && edgey ) {
338  theErrorX = 24.0*microntocm;
339  theErrorY = 96.0*microntocm;
340  }
341  else
342  {
343  theErrorX = 31.0*microntocm;
344  theErrorY = 90.0*microntocm;
345  }
346 
347  }
348  else {
349  if( singlex )
350  if ( hitbigx )
351  theErrorX = sx2*microntocm;
352  else
353  theErrorX = sx1*microntocm;
354  else theErrorX = sigmax*microntocm;
355  if( singley )
356  if( hitbigy )
357  theErrorY = sy2*microntocm;
358  else
359  theErrorY = sy1*microntocm;
360  else theErrorY = sigmay*microntocm;
361  }
362  theErrorZ = 1e-8; // 1 um means zero
364  // Local Error is 2D: (xx,xy,yy), square of sigma in first an third position
365  // as for resolution matrix
366  //
367 #ifdef FAMOS_DEBUG
368  std::cout << " Pixel Errors "
369  << "\talpha(x) = " << theErrorX
370  << "\tbeta(y) = " << theErrorY
371  << std::endl;
372 #endif
373  // Generate position
374  // get resolution histograms
375  int cotalphaHistBin = (int)( ( cotalpha - rescotAlpha_binMin ) / rescotAlpha_binWidth + 1 );
376  int cotbetaHistBin = (int)( ( cotbeta - rescotBeta_binMin ) / rescotBeta_binWidth + 1 );
377  // protection against out-of-range (undeflows and overflows)
378  if( cotalphaHistBin < 1 ) cotalphaHistBin = 1;
379  if( cotbetaHistBin < 1 ) cotbetaHistBin = 1;
380  if( cotalphaHistBin > (int)rescotAlpha_binN ) cotalphaHistBin = (int)rescotAlpha_binN;
381  if( cotbetaHistBin > (int)rescotBeta_binN ) cotbetaHistBin = (int)rescotBeta_binN;
382  //
383  unsigned int theXHistN;
384  unsigned int theYHistN;
385  if( !isForward ) {
386  if(edge)
387  {
388  theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
389  theYHistN = theXHistN;
390  }
391  else
392  {
393  if(singlex)
394  {
395  if(hitbigx) theXHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
396  else theXHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
397  }
398  else
399  {
400  if(hasBigPixelInX) theXHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
401  else theXHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
402  }
403  if(singley)
404  {
405  if(hitbigy) theYHistN = 1 * 100000 + cotalphaHistBin * 100 + cotbetaHistBin ;
406  else theYHistN = 1 * 100000 + 1 * 1000 + cotalphaHistBin * 100 + cotbetaHistBin ;
407  }
408  else
409  {
410  if(hasBigPixelInY) theYHistN = 1 * 1000000 + 1 * 100000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
411  else theYHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
412  }
413  }
414  }
415  else
416  {
417  if(edge)
418  {
419  theXHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
420  theYHistN = theXHistN;
421  }
422  else
423  {
424  if( singlex )
425  if( hitbigx )
426  theXHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
427  else
428  theXHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
429  else
430  theXHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
431  if( singley )
432  if( hitbigy )
433  theYHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
434  else
435  theYHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
436  else
437  theYHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + (nqbin+1);
438  }
439  }
440  unsigned int counter = 0;
441  do {
442  //
443  // Smear the hit Position
444  thePositionX = theXHistos[theXHistN]->generate(random);
445  thePositionY = theYHistos[theYHistN]->generate(random);
446  if( isForward ) thePositionY *= sign;
447  thePositionZ = 0.0; // set at the centre of the active area
448  //protect from empty resolution histograms
449  //if( thePositionX > 0.0799 ) thePositionX = 0;
450  //if( thePositionY > 0.0799 ) thePositionY = 0;
451  thePosition =
452  Local3DPoint(simHit.localPosition().x() + thePositionX ,
453  simHit.localPosition().y() + thePositionY ,
454  simHit.localPosition().z() + thePositionZ );
455 #ifdef FAMOS_DEBUG
456  std::cout << " Detector bounds: "
457  << "\t\tx = " << boundX
458  << "\ty = " << boundY
459  << std::endl;
460  std::cout << " Generated local position "
461  << "\tx = " << thePosition.x()
462  << "\ty = " << thePosition.y()
463  << std::endl;
464 #endif
465  counter++;
466  if(counter > 20) {
467  thePosition = Local3DPoint(simHit.localPosition().x(),
468  simHit.localPosition().y(),
469  simHit.localPosition().z());
470  break;
471  }
472  } while(fabs(thePosition.x()) > boundX || fabs(thePosition.y()) > boundY);
473 
474 }
475 
476 //-----------------------------------------------------------------------------
477 // I COPIED FROM THE PixelCPEBase BECAUSE IT'S BETTER THAN REINVENT IT
478 // The isFlipped() is a silly way to determine which detectors are inverted.
479 // In the barrel for every 2nd ladder the E field direction is in the
480 // global r direction (points outside from the z axis), every other
481 // ladder has the E field inside. Something similar is in the
482 // forward disks (2 sides of the blade). This has to be recognised
483 // because the charge sharing effect is different.
484 //
485 // The isFliped does it by looking and the relation of the local (z always
486 // in the E direction) to global coordinates. There is probably a much
487 // better way.(PJ: And faster!)
488 //-----------------------------------------------------------------------------
490  // Check the relative position of the local +/- z in global coordinates.
491  float tmp1 = theDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
492  float tmp2 = theDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
493  // std::cout << " 1: " << tmp1 << " 2: " << tmp2 << std::endl;
494  if ( tmp2<tmp1 ) return true;
495  else return false;
496 }
497 
499 {
500  //Hard coded at the moment, can easily be changed to be configurable
501  rescotAlpha_binMin = -0.2;
502  rescotAlpha_binWidth = 0.08 ;
503  rescotAlpha_binN = 5;
504  rescotBeta_binMin = -5.5;
505  rescotBeta_binWidth = 1.0;
506  rescotBeta_binN = 11;
507  resqbin_binMin = 0;
508  resqbin_binWidth = 1;
509  resqbin_binN = 4;
510 
511  // Initialize the barrel histos once and for all, and prepare the random generation
512  for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
513  for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin ) {
514  unsigned int singleBigPixelHistN = 1 * 100000
515  + cotalphaHistBin * 100
516  + cotbetaHistBin ;
517  theXHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
518  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singleBigPixelHistN ) ) );
519  theYHistos[singleBigPixelHistN] = new SimpleHistogramGenerator(
520  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singleBigPixelHistN ) ) );
521  unsigned int singlePixelHistN = 1 * 100000 + 1 * 1000
522  + cotalphaHistBin * 100
523  + cotbetaHistBin ;
524  theXHistos[singlePixelHistN] = new SimpleHistogramGenerator(
525  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , singlePixelHistN ) ) );
526  theYHistos[singlePixelHistN] = new SimpleHistogramGenerator(
527  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , singlePixelHistN ) ) );
528  for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
529  unsigned int edgePixelHistN = cotalphaHistBin * 1000
530  + cotbetaHistBin * 10
531  + qbinBin;
532  theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
533  (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hx0%u" ,edgePixelHistN ) ) );
534  theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
535  (TH1F*) thePixelResolutionFile2->Get( Form( "DQMData/clustBPIX/hy0%u" ,edgePixelHistN ) ) );
536  unsigned int multiPixelBigHistN = 1 * 1000000 + 1 * 100000
537  + cotalphaHistBin * 1000
538  + cotbetaHistBin * 10
539  + qbinBin;
540  theXHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
541  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" ,multiPixelBigHistN ) ) );
542  theYHistos[multiPixelBigHistN] = new SimpleHistogramGenerator(
543  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" ,multiPixelBigHistN ) ) );
544  unsigned int multiPixelHistN = 1 * 1000000 + 1 * 100000 + 1 * 10000
545  + cotalphaHistBin * 1000
546  + cotbetaHistBin * 10
547  + qbinBin;
548  theXHistos[multiPixelHistN] = new SimpleHistogramGenerator(
549  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hx%u" , multiPixelHistN ) ) );
550  theYHistos[multiPixelHistN] = new SimpleHistogramGenerator(
551  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustBPIX/hy%u" , multiPixelHistN ) ) );
552  } //end for qbinBin
553  }//end for cotalphaHistBin, cotbetaHistBin
554 }
555 
557 {
558  //Hard coded at the moment, can easily be changed to be configurable
559  rescotAlpha_binMin = 0.1;
560  rescotAlpha_binWidth = 0.1 ;
561  rescotAlpha_binN = 4;
562  rescotBeta_binMin = 0.;
563  rescotBeta_binWidth = 0.15;
564  rescotBeta_binN = 4;
565  resqbin_binMin = 0;
566  resqbin_binWidth = 1;
567  resqbin_binN = 4;
568 
569  // Initialize the forward histos once and for all, and prepare the random generation
570  for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
571  for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
572  for( unsigned qbinBin=1; qbinBin<=resqbin_binN; ++qbinBin ) {
573  unsigned int edgePixelHistN = cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
574  theXHistos[edgePixelHistN] = new SimpleHistogramGenerator(
575  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx0%u" ,edgePixelHistN ) ) );
576  theYHistos[edgePixelHistN] = new SimpleHistogramGenerator(
577  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy0%u" ,edgePixelHistN ) ) );
578  unsigned int PixelHistN = 100000 + 10000 + cotalphaHistBin * 1000 + cotbetaHistBin * 10 + qbinBin;
579  theXHistos[PixelHistN] = new SimpleHistogramGenerator(
580  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,PixelHistN ) ) );
581  theYHistos[PixelHistN] = new SimpleHistogramGenerator(
582  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,PixelHistN ) ) );
583  }//end cotalphaHistBin, cotbetaHistBin, qbinBin
584 
585  for ( unsigned cotalphaHistBin=1; cotalphaHistBin<=rescotAlpha_binN; ++cotalphaHistBin )
586  for ( unsigned cotbetaHistBin=1; cotbetaHistBin<=rescotBeta_binN; ++cotbetaHistBin )
587  {
588  unsigned int SingleBigPixelHistN = 100000 + cotalphaHistBin * 100 + cotbetaHistBin;
589  theXHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
590  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SingleBigPixelHistN ) ) );
591  theYHistos[SingleBigPixelHistN] = new SimpleHistogramGenerator(
592  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SingleBigPixelHistN ) ) );
593  unsigned int SinglePixelHistN = 100000 + 1000 + cotalphaHistBin * 100 + cotbetaHistBin;
594  theXHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
595  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhx%u" ,SinglePixelHistN ) ) );
596  theYHistos[SinglePixelHistN] = new SimpleHistogramGenerator(
597  (TH1F*) thePixelResolutionFile1->Get( Form( "DQMData/clustFPIX/fhy%u" ,SinglePixelHistN ) ) );
598 
599  }
600 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
T getParameter(std::string const &) const
bool containsBigPixelInY(int iymin, int iymax) const
int i
Definition: DBlmapReader.cc:9
#define BXSIZE
T y() const
Definition: PV2DBase.h:46
T perp() const
Definition: PV3DBase.h:72
bool isItEdgePixelInY(int iybin) const
#define BYSIZE
double flatShoot(double xmin=0.0, double xmax=1.0) const
virtual int ncolumns() const =0
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:47
virtual bool isItBigPixelInY(const int iybin) const
double sign(double x)
bool isItEdgePixelInX(int ixbin) const
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore > &thePixelTemp_)
T y() const
Definition: PV3DBase.h:63
virtual int nrows() const =0
TRandom random
Definition: MVATrainer.cc:138
const Int_t ysize
virtual MeasurementPoint measurementPosition(const LocalPoint &lp) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
std::map< unsigned, const SimpleHistogramGenerator * > theXHistos
bool interpolate(int id, float cotalpha, float cotbeta, float locBz)
void ytemp(int fybin, int lybin, float ytemplate[41][21+4])
Local3DPoint localPosition() const
Definition: PSimHit.h:44
#define BHX
T z() const
Definition: PV3DBase.h:64
void temperrors(int id, float cotalpha, float cotbeta, int qBin, float &sigmay, float &sigmax, float &sy1, float &sy2, float &sx1, float &sx2)
bool containsBigPixelInX(int ixmin, int ixmax) const
float s50()
1/2 of the pixel threshold signal in electrons
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:52
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
Vector3DBase unit() const
Definition: Vector3DBase.h:57
virtual LocalPoint localPosition(const MeasurementPoint &mp) const
void xtemp(int fxbin, int lxbin, float xtemplate[41][13+4])
void qbin_dist(int id, float cotalpha, float cotbeta, float qbin_frac[4], float &ny1_frac, float &ny2_frac, float &nx1_frac, float &nx2_frac)
std::map< unsigned, const SimpleHistogramGenerator * > theYHistos
virtual const TopologyType & specificTopology() const
void smearHit(const PSimHit &simHit, const PixelGeomDetUnit *detUnit, const double boundX, const double boundY, RandomEngineAndDistribution const *)
#define BHY
static std::atomic< unsigned int > counter
virtual bool isItBigPixelInX(const int ixbin) const
tuple cout
Definition: gather_cfg.py:121
SiPixelGaussianSmearingRecHitConverterAlgorithm(const edm::ParameterSet &pset, GeomDetType::SubDetector pixelPart)
for(const auto &isodef:isoDefs)
const Int_t xsize
T x() const
Definition: PV2DBase.h:45
T x() const
Definition: PV3DBase.h:62
virtual const PixelGeomDetType & specificType() const