CMS 3D CMS Logo

SiPixelChargeReweightingAlgorithm.cc
Go to the documentation of this file.
1 //class SiPixelChargeReweightingAlgorithm SimTracker/SiPixelDigitizer/src/SiPixelChargeReweightingAlgorithm.cc
2 
3 // Original Author Caroline Collard
4 // September 2020 : Extraction of the code for cluster charge reweighting from SiPixelDigitizerAlgorithm to a new class
5 //
6 #include <iostream>
7 #include <iomanip>
8 
10 
16 
17 //#include "PixelIndices.h"
19 
23 
24 // Accessing dead pixel modules from the DB:
26 
28 
33 
35 
36 // Geometry
40 
41 using namespace edm;
42 using namespace sipixelobjects;
43 
45  // Read template files for charge reweighting
46  if (UseReweighting || applyLateReweighting_) {
47  dbobject_den = &es.getData(SiPixel2DTemp_den_token_);
48  dbobject_num = &es.getData(SiPixel2DTemp_num_token_);
49 
50  int numOfTemplates = dbobject_den->numOfTempl() + dbobject_num->numOfTempl();
51  templateStores_.reserve(numOfTemplates);
52  SiPixelTemplate2D::pushfile(*dbobject_den, templateStores_);
53  SiPixelTemplate2D::pushfile(*dbobject_num, templateStores_);
54 
55  track.resize(6);
56  }
57 }
58 
59 //=========================================================================
60 
63  :
64 
65  templ2D(templateStores_),
66  xdouble(TXSIZE),
67  ydouble(TYSIZE),
68  IDnum(conf.exists("TemplateIDnumerator") ? conf.getParameter<int>("TemplateIDnumerator") : 0),
69  IDden(conf.exists("TemplateIDdenominator") ? conf.getParameter<int>("TemplateIDdenominator") : 0),
70 
71  UseReweighting(conf.getParameter<bool>("UseReweighting")),
72  applyLateReweighting_(conf.exists("applyLateReweighting") ? conf.getParameter<bool>("applyLateReweighting")
73  : false),
74  PrintClusters(conf.getParameter<bool>("PrintClusters")),
75  PrintTemplates(conf.getParameter<bool>("PrintTemplates")) {
77  SiPixel2DTemp_den_token_ = iC.esConsumes(edm::ESInputTag("", "denominator"));
79  }
80  edm::LogVerbatim("PixelDigitizer ") << "SiPixelChargeReweightingAlgorithm constructed"
81  << " with UseReweighting = " << UseReweighting
82  << " and applyLateReweighting_ = " << applyLateReweighting_;
83 }
84 
85 //=========================================================================
87  LogDebug("PixelDigitizer") << "SiPixelChargeReweightingAlgorithm deleted";
88 }
89 
90 //============================================================================
91 
93  std::map<int, float, std::less<int> >& hit_signal,
94  const size_t hitIndex,
95  const size_t hitIndex4CR,
96  const unsigned int tofBin,
97  const PixelTopology* topol,
98  uint32_t detID,
99  signal_map_type& theSignal,
100  unsigned short int processType,
101  const bool& boolmakeDigiSimLinks) {
102  int irow_min = topol->nrows();
103  int irow_max = 0;
104  int icol_min = topol->ncolumns();
105  int icol_max = 0;
106 
107  float chargeBefore = 0;
108  float chargeAfter = 0;
109  signal_map_type hitSignal;
110  LocalVector direction = hit.exitPoint() - hit.entryPoint();
111 
112  for (std::map<int, float, std::less<int> >::const_iterator im = hit_signal.begin(); im != hit_signal.end(); ++im) {
113  int chan = (*im).first;
114  std::pair<int, int> pixelWithCharge = PixelDigi::channelToPixel(chan);
115 
116  hitSignal[chan] += (boolmakeDigiSimLinks ? SiPixelDigitizerAlgorithm::Amplitude(
117  (*im).second, &hit, hitIndex, hitIndex4CR, tofBin, (*im).second)
118  : SiPixelDigitizerAlgorithm::Amplitude((*im).second, (*im).second));
119  chargeBefore += (*im).second;
120 
121  if (pixelWithCharge.first < irow_min)
122  irow_min = pixelWithCharge.first;
123  if (pixelWithCharge.first > irow_max)
124  irow_max = pixelWithCharge.first;
125  if (pixelWithCharge.second < icol_min)
126  icol_min = pixelWithCharge.second;
127  if (pixelWithCharge.second > icol_max)
128  icol_max = pixelWithCharge.second;
129  }
130 
131  LocalPoint hitEntryPoint = hit.entryPoint();
132 
133  float trajectoryScaleToPosition = std::abs(hitEntryPoint.z() / direction.z());
134 
135  LocalPoint hitPosition = hitEntryPoint + trajectoryScaleToPosition * direction;
136 
137  MeasurementPoint hitPositionPixel = topol->measurementPosition(hit.localPosition());
138  std::pair<int, int> hitPixel =
139  std::pair<int, int>(int(floor(hitPositionPixel.x())), int(floor(hitPositionPixel.y())));
140 
141  MeasurementPoint originPixel = MeasurementPoint(hitPixel.first - THX + 0.5, hitPixel.second - THY + 0.5);
142  LocalPoint origin = topol->localPosition(originPixel);
143 
144  MeasurementPoint hitEntryPointPixel = topol->measurementPosition(hit.entryPoint());
145  MeasurementPoint hitExitPointPixel = topol->measurementPosition(hit.exitPoint());
146  std::pair<int, int> entryPixel =
147  std::pair<int, int>(int(floor(hitEntryPointPixel.x())), int(floor(hitEntryPointPixel.y())));
148  std::pair<int, int> exitPixel =
149  std::pair<int, int>(int(floor(hitExitPointPixel.x())), int(floor(hitExitPointPixel.y())));
150 
151  int hitcol_min, hitcol_max, hitrow_min, hitrow_max;
152  if (entryPixel.first > exitPixel.first) {
153  hitrow_min = exitPixel.first;
154  hitrow_max = entryPixel.first;
155  } else {
156  hitrow_min = entryPixel.first;
157  hitrow_max = exitPixel.first;
158  }
159 
160  if (entryPixel.second > exitPixel.second) {
161  hitcol_min = exitPixel.second;
162  hitcol_max = entryPixel.second;
163  } else {
164  hitcol_min = entryPixel.second;
165  hitcol_max = exitPixel.second;
166  }
167 
168 #ifdef TP_DEBUG
169  LocalPoint CMSSWhitPosition = hit.localPosition();
170 
171  LogDebug("Pixel Digitizer") << "\n"
172  << "Particle ID is: " << hit.particleType() << "\n"
173  << "Process type: " << hit.processType() << "\n"
174  << "HitPosition:"
175  << "\n"
176  << "Hit entry x/y/z: " << hit.entryPoint().x() << " " << hit.entryPoint().y() << " "
177  << hit.entryPoint().z() << " "
178  << "Hit exit x/y/z: " << hit.exitPoint().x() << " " << hit.exitPoint().y() << " "
179  << hit.exitPoint().z() << " "
180 
181  << "Pixel Pos - X: " << hitPositionPixel.x() << " Y: " << hitPositionPixel.y() << "\n"
182  << "Cart.Cor. - X: " << CMSSWhitPosition.x() << " Y: " << CMSSWhitPosition.y() << "\n"
183  << "Z=0 Pos - X: " << hitPosition.x() << " Y: " << hitPosition.y() << "\n"
184 
185  << "Origin of the template:"
186  << "\n"
187  << "Pixel Pos - X: " << originPixel.x() << " Y: " << originPixel.y() << "\n"
188  << "Cart.Cor. - X: " << origin.x() << " Y: " << origin.y() << "\n"
189  << "\n"
190  << "Entry/Exit:"
191  << "\n"
192  << "Entry - X: " << hit.entryPoint().x() << " Y: " << hit.entryPoint().y()
193  << " Z: " << hit.entryPoint().z() << "\n"
194  << "Exit - X: " << hit.exitPoint().x() << " Y: " << hit.exitPoint().y()
195  << " Z: " << hit.exitPoint().z() << "\n"
196 
197  << "Entry - X Pixel: " << hitEntryPointPixel.x() << " Y Pixel: " << hitEntryPointPixel.y()
198  << "\n"
199  << "Exit - X Pixel: " << hitExitPointPixel.x() << " Y Pixel: " << hitExitPointPixel.y()
200  << "\n"
201 
202  << "row min: " << irow_min << " col min: " << icol_min << "\n";
203 #endif
204 
205  if (!(irow_min <= hitrow_max && irow_max >= hitrow_min && icol_min <= hitcol_max && icol_max >= hitcol_min)) {
206  // The clusters do not have an overlap, hence the hit is NOT reweighted
207  return false;
208  }
209 
210  track[0] = (hitPosition.x() - origin.x()) * cmToMicrons;
211  track[1] = (hitPosition.y() - origin.y()) * cmToMicrons;
212  track[2] = 0.0f; //Middle of sensor is origin for Z-axis
213  track[3] = direction.x();
214  track[4] = direction.y();
215  track[5] = direction.z();
216 
217  array_2d pixrewgt(boost::extents[TXSIZE][TYSIZE]);
218 
219  for (int row = 0; row < TXSIZE; ++row) {
220  for (int col = 0; col < TYSIZE; ++col) {
221  pixrewgt[row][col] = 0;
222  }
223  }
224 
225  for (int row = 0; row < TXSIZE; ++row) {
226  xdouble[row] = topol->isItBigPixelInX(hitPixel.first + row - THX);
227  }
228 
229  for (int col = 0; col < TYSIZE; ++col) {
230  ydouble[col] = topol->isItBigPixelInY(hitPixel.second + col - THY);
231  }
232 
233  // define loop boundaries that will prevent the row and col loops
234  // from going out of physical bounds of the pixel module
235  int rowmin = std::max(0, THX - hitPixel.first);
236  int rowmax = std::min(TXSIZE, topol->nrows() + THX - hitPixel.first);
237  int colmin = std::max(0, THY - hitPixel.second);
238  int colmax = std::min(TYSIZE, topol->ncolumns() + THY - hitPixel.second);
239 
240  for (int row = rowmin; row < rowmax; ++row) {
241  for (int col = colmin; col < colmax; ++col) {
242  //Fill charges into 21x13 Pixel Array with hitPixel in centre
243  pixrewgt[row][col] =
244  hitSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)];
245  }
246  }
247 
248  if (PrintClusters) {
249  LogDebug("PixelDigitizer ") << "Cluster before reweighting: ";
250  printCluster(pixrewgt);
251  }
252 
253  int ierr;
254  // for unirradiated: 2nd argument is IDden
255  // for irradiated: 2nd argument is IDnum
256  if (UseReweighting == true) {
257  int ID1 = dbobject_num->getTemplateID(detID);
258  int ID0 = dbobject_den->getTemplateID(detID);
259 
260  if (ID0 == ID1) {
261  return false;
262  }
263  ierr = PixelTempRewgt2D(ID0, ID1, pixrewgt);
264  } else {
265  ierr = PixelTempRewgt2D(IDden, IDden, pixrewgt);
266  }
267  if (ierr != 0) {
268 #ifdef TP_DEBUG
269  LogDebug("PixelDigitizer ") << "Cluster Charge Reweighting did not work properly.";
270 #endif
271  return false;
272  }
273 
274  if (PrintClusters) {
275  LogDebug("PixelDigitizer ") << "Cluster after reweighting: ";
276  printCluster(pixrewgt);
277  }
278 
279  for (int row = rowmin; row < rowmax; ++row) {
280  for (int col = colmin; col < colmax; ++col) {
281  float charge = pixrewgt[row][col];
282  if (charge > 0) {
283  chargeAfter += charge;
284  theSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)] +=
285  (boolmakeDigiSimLinks
286  ? SiPixelDigitizerAlgorithm::Amplitude(charge, &hit, hitIndex, hitIndex4CR, tofBin, charge)
288  }
289  }
290  }
291 
292  if (chargeBefore != 0. && chargeAfter == 0.) {
293  return false;
294  }
295 
296  if (PrintClusters) {
297  LogDebug("PixelDigitizer ") << "Charges (before->after): " << chargeBefore << " -> " << chargeAfter;
298  LogDebug("PixelDigitizer ") << "Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 << " % \n";
299  }
300 
301  return true;
302 }
303 
304 // *******************************************************************************************************
312 // *******************************************************************************************************
313 int SiPixelChargeReweightingAlgorithm::PixelTempRewgt2D(int id_in, int id_rewgt, array_2d& cluster) {
314  // Local variables
315  int i, j, k, l, kclose;
316  int nclusx, nclusy, success;
317  float xsize, ysize, q50i, q100i, q50r, q10r, q100r, xhit2D, yhit2D, qclust, dist2, dmin2;
318  float xy_in[BXM2][BYM2], xy_rewgt[BXM2][BYM2], xy_clust[TXSIZE][TYSIZE];
319  int denx_clust[TXSIZE][TYSIZE], deny_clust[TXSIZE][TYSIZE];
320  int goodWeightsUsed, nearbyWeightsUsed, noWeightsUsed;
321  float cotalpha, cotbeta;
322  // success = 0 is returned if everthing is OK
323  success = 0;
324 
325  // Copy the array to remember original charges
326  array_2d clust(cluster);
327 
328  // Take the pixel dimensions from the 2D template
329  templ2D.getid(id_in);
330  xsize = templ2D.xsize();
331  ysize = templ2D.ysize();
332 
333  // Calculate the track angles
334 
335  if (std::abs(track[5]) > 0.f) {
336  cotalpha = track[3] / track[5]; //if track[5] (direction in z) is 0 the hit is not processed by re-weighting
337  cotbeta = track[4] / track[5];
338  } else {
339  LogDebug("Pixel Digitizer") << "Reweighting angle is not good! \n";
340  return 9; //returned value here indicates that no reweighting was done in this case
341  }
342 
343  // The 2-D templates are defined on a shifted coordinate system wrt the 1D templates
344  if (ydouble[0]) {
345  yhit2D = track[1] - cotbeta * track[2] + ysize;
346  } else {
347  yhit2D = track[1] - cotbeta * track[2] + 0.5f * ysize;
348  }
349  if (xdouble[0]) {
350  xhit2D = track[0] - cotalpha * track[2] + xsize;
351  } else {
352  xhit2D = track[0] - cotalpha * track[2] + 0.5f * xsize;
353  }
354 
355  // Zero the input and output templates
356  for (i = 0; i < BYM2; ++i) {
357  for (j = 0; j < BXM2; ++j) {
358  xy_in[j][i] = 0.f;
359  xy_rewgt[j][i] = 0.f;
360  }
361  }
362 
363  // Next, interpolate the CMSSW template needed to analyze this cluster
364 
365  if (!templ2D.xytemp(id_in, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_in)) {
366  success = 1;
367  }
368  if (success != 0) {
369 #ifdef TP_DEBUG
370  LogDebug("Pixel Digitizer") << "No matching template found \n";
371 #endif
372  return 2;
373  }
374 
375  if (PrintTemplates) {
376  LogDebug("Pixel Digitizer") << "Template unirrad: \n";
377  printCluster(xy_in);
378  }
379 
380  q50i = templ2D.s50();
381  //q50i = 0;
382  q100i = 2.f * q50i;
383 
384  // Check that the cluster container is a 13x21 matrix
385 
386  if (cluster.num_dimensions() != 2) {
387  LogWarning("Pixel Digitizer") << "Cluster is not 2-dimensional. Return. \n";
388  return 3;
389  }
390  nclusx = (int)cluster.shape()[0];
391  nclusy = (int)cluster.shape()[1];
392  if (nclusx != TXSIZE || xdouble.size() != TXSIZE) {
393  LogWarning("Pixel Digitizer") << "Sizes in x do not match: nclusx=" << nclusx << " xdoubleSize=" << xdouble.size()
394  << " TXSIZE=" << TXSIZE << ". Return. \n";
395  return 4;
396  }
397  if (nclusy != TYSIZE || ydouble.size() != TYSIZE) {
398  LogWarning("Pixel Digitizer") << "Sizes in y do not match. Return. \n";
399  return 5;
400  }
401 
402  // Sum initial charge in the cluster
403 
404  qclust = 0.f;
405  for (i = 0; i < TYSIZE; ++i) {
406  for (j = 0; j < TXSIZE; ++j) {
407  xy_clust[j][i] = 0.f;
408  denx_clust[j][i] = 0;
409  deny_clust[j][i] = 0;
410  if (cluster[j][i] > q100i) {
411  qclust += cluster[j][i];
412  }
413  }
414  }
415 
416  // Next, interpolate the physical output template needed to reweight
417 
418  if (!templ2D.xytemp(id_rewgt, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_rewgt)) {
419  success = 1;
420  }
421 
422  if (PrintTemplates) {
423  LogDebug("Pixel Digitizer") << "Template irrad: \n";
424  printCluster(xy_rewgt);
425  }
426 
427  q50r = templ2D.s50();
428  q100r = 2.f * q50r;
429  q10r = 0.2f * q50r;
430 
431  // Find all non-zero denominator pixels in the input template and generate "inside" weights
432 
433  int ntpix = 0;
434  int ncpix = 0;
435  std::vector<int> ytclust;
436  std::vector<int> xtclust;
437  std::vector<int> ycclust;
438  std::vector<int> xcclust;
439  qclust = 0.f;
440  for (i = 0; i < TYSIZE; ++i) {
441  for (j = 0; j < TXSIZE; ++j) {
442  if (xy_in[j + 1][i + 1] > q100i) {
443  ++ntpix;
444  ytclust.push_back(i);
445  xtclust.push_back(j);
446  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[j + 1][i + 1];
447  denx_clust[j][i] = j;
448  deny_clust[j][i] = i;
449  }
450  }
451  }
452 
453  // Find all non-zero numerator pixels not matched to denominator in the output template and generate "inside" weights
454 
455  for (i = 0; i < TYSIZE; ++i) {
456  for (j = 0; j < TXSIZE; ++j) {
457  if (xy_rewgt[j + 1][i + 1] > q10r && xy_clust[j][i] == 0.f && ntpix > 0) {
458  // Search for nearest denominator pixel
459  dmin2 = 10000.f;
460  kclose = 0;
461  for (k = 0; k < ntpix; ++k) {
462  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
463  if (dist2 < dmin2) {
464  dmin2 = dist2;
465  kclose = k;
466  }
467  }
468  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[xtclust[kclose] + 1][ytclust[kclose] + 1];
469  denx_clust[j][i] = xtclust[kclose];
470  deny_clust[j][i] = ytclust[kclose];
471  }
472  }
473  }
474 
475  if (PrintTemplates) {
476  LogDebug("Pixel Digitizer") << "Weights: \n";
477  printCluster(xy_clust);
478  }
479 
480  // Do the reweighting
481  goodWeightsUsed = 0;
482  nearbyWeightsUsed = 0;
483  noWeightsUsed = 0;
484 
485  for (i = 0; i < TYSIZE; ++i) {
486  for (j = 0; j < TXSIZE; ++j) {
487  if (xy_clust[j][i] > 0.f) {
488  cluster[j][i] = xy_clust[j][i] * clust[denx_clust[j][i]][deny_clust[j][i]];
489  if (cluster[j][i] > q100r) {
490  qclust += cluster[j][i];
491  }
492  if (cluster[j][i] > 0) {
493  goodWeightsUsed++;
494  }
495  } else {
496  if (clust[j][i] > 0.f) {
497  ++ncpix;
498  ycclust.push_back(i);
499  xcclust.push_back(j);
500  }
501  }
502  }
503  }
504 
505  // Now reweight pixels outside of template footprint using closest weights
506 
507  if (ncpix > 0) {
508  for (l = 0; l < ncpix; ++l) {
509  i = ycclust[l];
510  j = xcclust[l];
511  dmin2 = 10000.f;
512  kclose = 0;
513  for (k = 0; k < ntpix; ++k) {
514  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
515  if (dist2 < dmin2) {
516  dmin2 = dist2;
517  kclose = k;
518  }
519  }
520  if (dmin2 < 5.f) {
521  nearbyWeightsUsed++;
522  cluster[j][i] *= xy_clust[xtclust[kclose]][ytclust[kclose]];
523  if (cluster[j][i] > q100r) {
524  qclust += cluster[j][i];
525  }
526  } else {
527  noWeightsUsed++;
528  cluster[j][i] = 0.f;
529  }
530  }
531  }
532 
533  return success;
534 } // PixelTempRewgt2D
535 
537  for (int col = 0; col < TYSIZE; ++col) {
538  for (int row = 0; row < TXSIZE; ++row) {
539  LogDebug("Pixel Digitizer") << cluster[row][col];
540  }
541  LogDebug("Pixel Digitizer") << "\n";
542  }
543 }
544 
546  for (int col = 0; col < BYM2; ++col) {
547  for (int row = 0; row < BXM2; ++row) {
548  LogDebug("Pixel Digitizer") << arr[row][col];
549  }
550  LogDebug("Pixel Digitizer") << "\n";
551  }
552 }
553 
555  for (int col = 0; col < TYSIZE; ++col) {
556  for (int row = 0; row < TXSIZE; ++row) {
557  LogDebug("Pixel Digitizer") << arr[row][col];
558  }
559  LogDebug("Pixel Digitizer") << "\n";
560  }
561 }
562 
564  std::vector<PixelDigi>& digis,
565  PixelSimHitExtraInfo& loopTempSH,
566  signal_map_type& theNewDigiSignal,
567  const TrackerTopology* tTopo,
568  CLHEP::HepRandomEngine* engine) {
569  uint32_t detID = pixdet->geographicalId().rawId();
570  const PixelTopology* topol = &pixdet->specificTopology();
571 
572  if (UseReweighting) {
573  LogError("Pixel Digitizer") << " ******************************** \n";
574  LogError("Pixel Digitizer") << " ******************************** \n";
575  LogError("Pixel Digitizer") << " ******************************** \n";
576  LogError("Pixel Digitizer") << " ***** INCONSISTENCY !!! ***** \n";
577  LogError("Pixel Digitizer")
578  << " applyLateReweighting_ and UseReweighting can not be true at the same time for PU ! \n";
579  LogError("Pixel Digitizer") << " ---> DO NOT APPLY CHARGE REWEIGHTING TWICE !!! \n";
580  LogError("Pixel Digitizer") << " ******************************** \n";
581  LogError("Pixel Digitizer") << " ******************************** \n";
582  return false;
583  }
584 
585  signal_map_type theDigiSignal;
586 
587  int irow_min = topol->nrows();
588  int irow_max = 0;
589  int icol_min = topol->ncolumns();
590  int icol_max = 0;
591 
592  float chargeBefore = 0;
593  float chargeAfter = 0;
594 
595  //loop on digis
596  std::vector<PixelDigi>::const_iterator loopDigi;
597  for (loopDigi = digis.begin(); loopDigi != digis.end(); ++loopDigi) {
598  unsigned int chan = loopDigi->channel();
599  if (loopTempSH.isInTheList(chan)) {
600  float corresponding_charge = loopDigi->adc();
601  theDigiSignal[chan] += SiPixelDigitizerAlgorithm::Amplitude(corresponding_charge, corresponding_charge);
602  chargeBefore += corresponding_charge;
603  if (loopDigi->row() < irow_min)
604  irow_min = loopDigi->row();
605  if (loopDigi->row() > irow_max)
606  irow_max = loopDigi->row();
607  if (loopDigi->column() < icol_min)
608  icol_min = loopDigi->column();
609  if (loopDigi->column() > icol_max)
610  icol_max = loopDigi->column();
611  }
612  }
613  // end loop on digis
614 
615  LocalVector direction = loopTempSH.exitPoint() - loopTempSH.entryPoint();
616  LocalPoint hitEntryPoint = loopTempSH.entryPoint();
617  float trajectoryScaleToPosition = std::abs(hitEntryPoint.z() / direction.z());
618  LocalPoint hitPosition = hitEntryPoint + trajectoryScaleToPosition * direction;
619 
620  // addition as now the hit himself is not available
621  // https://github.com/cms-sw/cmssw/blob/master/SimDataFormats/TrackingHit/interface/PSimHit.h#L52
622  LocalPoint hitLocalPosition = hitEntryPoint + 0.5 * direction;
623  MeasurementPoint hitPositionPixel = topol->measurementPosition(hitLocalPosition);
624  std::pair<int, int> hitPixel =
625  std::pair<int, int>(int(floor(hitPositionPixel.x())), int(floor(hitPositionPixel.y())));
626 
627  MeasurementPoint originPixel = MeasurementPoint(hitPixel.first - THX + 0.5, hitPixel.second - THY + 0.5);
628  LocalPoint origin = topol->localPosition(originPixel);
629 
630  MeasurementPoint hitEntryPointPixel = topol->measurementPosition(loopTempSH.entryPoint());
631  MeasurementPoint hitExitPointPixel = topol->measurementPosition(loopTempSH.exitPoint());
632  std::pair<int, int> entryPixel =
633  std::pair<int, int>(int(floor(hitEntryPointPixel.x())), int(floor(hitEntryPointPixel.y())));
634  std::pair<int, int> exitPixel =
635  std::pair<int, int>(int(floor(hitExitPointPixel.x())), int(floor(hitExitPointPixel.y())));
636 
637  int hitcol_min, hitcol_max, hitrow_min, hitrow_max;
638  if (entryPixel.first > exitPixel.first) {
639  hitrow_min = exitPixel.first;
640  hitrow_max = entryPixel.first;
641  } else {
642  hitrow_min = entryPixel.first;
643  hitrow_max = exitPixel.first;
644  }
645 
646  if (entryPixel.second > exitPixel.second) {
647  hitcol_min = exitPixel.second;
648  hitcol_max = entryPixel.second;
649  } else {
650  hitcol_min = entryPixel.second;
651  hitcol_max = exitPixel.second;
652  }
653 
654  if (!(irow_min <= hitrow_max && irow_max >= hitrow_min && icol_min <= hitcol_max && icol_max >= hitcol_min)) {
655  // The clusters do not have an overlap, hence the hit is NOT reweighted
656  return false;
657  }
658 
659  track[0] = (hitPosition.x() - origin.x()) * cmToMicrons;
660  track[1] = (hitPosition.y() - origin.y()) * cmToMicrons;
661  track[2] = 0.0f; //Middle of sensor is origin for Z-axis
662  track[3] = direction.x();
663  track[4] = direction.y();
664  track[5] = direction.z();
665 
666  array_2d pixrewgt(boost::extents[TXSIZE][TYSIZE]);
667 
668  /*
669  for (int row = 0; row < TXSIZE; ++row) {
670  for (int col = 0; col < TYSIZE; ++col) {
671  pixrewgt[row][col] = 0;
672  }
673  }
674 */
675 
676  for (int row = 0; row < TXSIZE; ++row) {
677  xdouble[row] = topol->isItBigPixelInX(hitPixel.first + row - THX);
678  }
679 
680  for (int col = 0; col < TYSIZE; ++col) {
681  ydouble[col] = topol->isItBigPixelInY(hitPixel.second + col - THY);
682  }
683 
684  for (int row = 0; row < TXSIZE; ++row) {
685  for (int col = 0; col < TYSIZE; ++col) {
686  //Fill charges into 21x13 Pixel Array with hitPixel in centre
687  pixrewgt[row][col] =
688  theDigiSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)];
689  }
690  }
691 
692  if (PrintClusters) {
693  LogDebug("Pixel Digitizer") << " Cluster before reweighting: ";
694  printCluster(pixrewgt);
695  }
696 
697  int ierr;
698  // for unirradiated: 2nd argument is IDden
699  // for irradiated: 2nd argument is IDnum
700  int ID1 = dbobject_num->getTemplateID(detID);
701  int ID0 = dbobject_den->getTemplateID(detID);
702 
703  if (ID0 == ID1) {
704  LogDebug("Pixel Digitizer") << " same template for num and den ";
705  return false;
706  }
707  ierr = PixelTempRewgt2D(ID0, ID1, pixrewgt);
708  if (ierr != 0) {
709 #ifdef TP_DEBUG
710  LogDebug("PixelDigitizer ") << "Cluster Charge Reweighting did not work properly.";
711 #endif
712  return false;
713  }
714  if (PrintClusters) {
715  LogDebug("Pixel Digitizer") << " Cluster after reweighting: ";
716  printCluster(pixrewgt);
717  }
718 
719  for (int row = 0; row < TXSIZE; ++row) {
720  for (int col = 0; col < TYSIZE; ++col) {
721  float charge = 0;
722  charge = pixrewgt[row][col];
723  if ((hitPixel.first + row - THX) >= 0 && (hitPixel.first + row - THX) < topol->nrows() &&
724  (hitPixel.second + col - THY) >= 0 && (hitPixel.second + col - THY) < topol->ncolumns() && charge > 0) {
725  chargeAfter += charge;
726  theNewDigiSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)] +=
728  }
729  }
730  }
731 
732  if (chargeBefore != 0. && chargeAfter == 0.) {
733  return false;
734  }
735 
736  if (PrintClusters) {
737  LogDebug("Pixel Digitizer") << "Charges (before->after): " << chargeBefore << " -> " << chargeAfter;
738  LogDebug("Pixel Digitizer") << "Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 << " %";
739  }
740 
741  // need to store the digi out of the 21x13 box.
742  for (signal_map_const_iterator i = theDigiSignal.begin(); i != theDigiSignal.end(); ++i) {
743  // check if in the 21x13 box
744  int chanDigi = (*i).first;
745  std::pair<int, int> ip = PixelDigi::channelToPixel(chanDigi);
746  int row_digi = ip.first;
747  int col_digi = ip.second;
748  int i_transformed_row = row_digi - hitPixel.first + THX;
749  int i_transformed_col = col_digi - hitPixel.second + THY;
750  if (i_transformed_row < 0 || i_transformed_row > TXSIZE || i_transformed_col < 0 || i_transformed_col > TYSIZE) {
751  // not in the box
752  if (chanDigi >= 0 && (*i).second > 0) {
753  theNewDigiSignal[chanDigi] += SiPixelDigitizerAlgorithm::Amplitude((*i).second, (*i).second);
754  }
755  }
756  }
757 
758  return true;
759 }
Log< level::Info, true > LogVerbatim
float xsize()
pixel x-size (microns)
bool xytemp(float xhit, float yhit, bool ydouble[21+2], bool xdouble[13+2], float template2d[13+2][21+2], bool dervatives, float dpdx2d[2][13+2][21+2], float &QTemplate)
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_num_token_
virtual int ncolumns() const =0
T z() const
Definition: PV3DBase.h:61
#define TXSIZE
bool lateSignalReweight(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, PixelSimHitExtraInfo &loopTempSH, signal_map_type &theNewDigiSignal, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *engine)
virtual int nrows() const =0
std::map< int, SiPixelDigitizerAlgorithm::Amplitude, std::less< int > > signal_map_type
const SiPixel2DTemplateDBObject * dbobject_den
T x() const
Definition: PV2DBase.h:43
const Int_t ysize
Log< level::Error, false > LogError
float s50()
1/2 of the pixel threshold signal in adc units
#define THX
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
#define BXM2
T y() const
Definition: PV2DBase.h:44
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
bool isInTheList(unsigned int channelToCheck)
virtual bool isItBigPixelInX(int ixbin) const =0
const Local3DPoint & entryPoint() const
virtual bool isItBigPixelInY(int iybin) const =0
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define BYM2
double f[11][100]
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
bool getData(T &iHolder) const
Definition: EventSetup.h:122
#define THY
signal_map_type::const_iterator signal_map_const_iterator
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
const SiPixel2DTemplateDBObject * dbobject_num
int PixelTempRewgt2D(int id_gen, int id_rewgt, array_2d &cluster)
#define TYSIZE
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
SiPixelChargeReweightingAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
float ysize()
pixel y-size (microns)
short getTemplateID(const uint32_t &detid) const
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_den_token_
HLT enums.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
col
Definition: cuy.py:1009
bool hitSignalReweight(const PSimHit &hit, std::map< int, float, std::less< int > > &hit_signal, const size_t hitIndex, const size_t hitIndex4CR, const unsigned int tofBin, const PixelTopology *topol, uint32_t detID, signal_map_type &theSignal, unsigned short int processType, const bool &boolmakeDigiSimLinks)
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore2D > &pixelTemp, std::string dir="CalibTracker/SiPixelESProducers/data/")
Log< level::Warning, false > LogWarning
const Int_t xsize
const Local3DPoint & exitPoint() const
#define LogDebug(id)