CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 using namespace edm;
37 using namespace sipixelobjects;
38 
40  // Read template files for charge reweighting
41  if (UseReweighting) {
42  dbobject_den = &es.getData(SiPixel2DTemp_den_token_);
43  dbobject_num = &es.getData(SiPixel2DTemp_num_token_);
44 
45  int numOfTemplates = dbobject_den->numOfTempl() + dbobject_num->numOfTempl();
46  templateStores_.reserve(numOfTemplates);
47  SiPixelTemplate2D::pushfile(*dbobject_den, templateStores_);
48  SiPixelTemplate2D::pushfile(*dbobject_num, templateStores_);
49 
50  track.resize(6);
51  }
52 }
53 
54 //=========================================================================
55 
58  :
59 
60  templ2D(templateStores_),
61  xdouble(TXSIZE),
62  ydouble(TYSIZE),
63  IDnum(conf.exists("TemplateIDnumerator") ? conf.getParameter<int>("TemplateIDnumerator") : 0),
64  IDden(conf.exists("TemplateIDdenominator") ? conf.getParameter<int>("TemplateIDdenominator") : 0),
65 
66  UseReweighting(conf.getParameter<bool>("UseReweighting")),
67  PrintClusters(conf.getParameter<bool>("PrintClusters")),
68  PrintTemplates(conf.getParameter<bool>("PrintTemplates")) {
69  if (UseReweighting) {
70  SiPixel2DTemp_den_token_ = iC.esConsumes(edm::ESInputTag("", "denominator"));
72  }
73  edm::LogVerbatim("PixelDigitizer ") << "SiPixelChargeReweightingAlgorithm constructed"
74  << " with UseReweighting = " << UseReweighting;
75 }
76 
77 //=========================================================================
79  LogDebug("PixelDigitizer") << "SiPixelChargeReweightingAlgorithm deleted";
80 }
81 
82 //============================================================================
83 
85  std::map<int, float, std::less<int> >& hit_signal,
86  const size_t hitIndex,
87  const unsigned int tofBin,
88  const PixelTopology* topol,
89  uint32_t detID,
90  signal_map_type& theSignal,
91  unsigned short int processType,
92  const bool& boolmakeDigiSimLinks) {
93  int irow_min = topol->nrows();
94  int irow_max = 0;
95  int icol_min = topol->ncolumns();
96  int icol_max = 0;
97 
98  float chargeBefore = 0;
99  float chargeAfter = 0;
100  signal_map_type hitSignal;
101  LocalVector direction = hit.exitPoint() - hit.entryPoint();
102 
103  for (std::map<int, float, std::less<int> >::const_iterator im = hit_signal.begin(); im != hit_signal.end(); ++im) {
104  int chan = (*im).first;
105  std::pair<int, int> pixelWithCharge = PixelDigi::channelToPixel(chan);
106 
107  hitSignal[chan] +=
108  (boolmakeDigiSimLinks ? SiPixelDigitizerAlgorithm::Amplitude((*im).second, &hit, hitIndex, tofBin, (*im).second)
109  : SiPixelDigitizerAlgorithm::Amplitude((*im).second, (*im).second));
110  chargeBefore += (*im).second;
111 
112  if (pixelWithCharge.first < irow_min)
113  irow_min = pixelWithCharge.first;
114  if (pixelWithCharge.first > irow_max)
115  irow_max = pixelWithCharge.first;
116  if (pixelWithCharge.second < icol_min)
117  icol_min = pixelWithCharge.second;
118  if (pixelWithCharge.second > icol_max)
119  icol_max = pixelWithCharge.second;
120  }
121 
122  LocalPoint hitEntryPoint = hit.entryPoint();
123 
124  float trajectoryScaleToPosition = hitEntryPoint.z() / direction.z();
125 
126  if ((hitEntryPoint.z() > 0 && direction.z() < 0) || (hitEntryPoint.z() < 0 && direction.z() > 0)) {
127  trajectoryScaleToPosition *= -1;
128  }
129 
130  LocalPoint hitPosition = hitEntryPoint + trajectoryScaleToPosition * direction;
131 
132  MeasurementPoint hitPositionPixel = topol->measurementPosition(hit.localPosition());
133  std::pair<int, int> hitPixel =
134  std::pair<int, int>(int(floor(hitPositionPixel.x())), int(floor(hitPositionPixel.y())));
135 
136  MeasurementPoint originPixel = MeasurementPoint(hitPixel.first - THX + 0.5, hitPixel.second - THY + 0.5);
137  LocalPoint origin = topol->localPosition(originPixel);
138 
139  MeasurementPoint hitEntryPointPixel = topol->measurementPosition(hit.entryPoint());
140  MeasurementPoint hitExitPointPixel = topol->measurementPosition(hit.exitPoint());
141  std::pair<int, int> entryPixel =
142  std::pair<int, int>(int(floor(hitEntryPointPixel.x())), int(floor(hitEntryPointPixel.y())));
143  std::pair<int, int> exitPixel =
144  std::pair<int, int>(int(floor(hitExitPointPixel.x())), int(floor(hitExitPointPixel.y())));
145 
146  int hitcol_min, hitcol_max, hitrow_min, hitrow_max;
147  if (entryPixel.first > exitPixel.first) {
148  hitrow_min = exitPixel.first;
149  hitrow_max = entryPixel.first;
150  } else {
151  hitrow_min = entryPixel.first;
152  hitrow_max = exitPixel.first;
153  }
154 
155  if (entryPixel.second > exitPixel.second) {
156  hitcol_min = exitPixel.second;
157  hitcol_max = entryPixel.second;
158  } else {
159  hitcol_min = entryPixel.second;
160  hitcol_max = exitPixel.second;
161  }
162 
163 #ifdef TP_DEBUG
164  LocalPoint CMSSWhitPosition = hit.localPosition();
165 
166  LogDebug("Pixel Digitizer") << "\n"
167  << "Particle ID is: " << hit.particleType() << "\n"
168  << "Process type: " << hit.processType() << "\n"
169  << "HitPosition:"
170  << "\n"
171  << "Hit entry x/y/z: " << hit.entryPoint().x() << " " << hit.entryPoint().y() << " "
172  << hit.entryPoint().z() << " "
173  << "Hit exit x/y/z: " << hit.exitPoint().x() << " " << hit.exitPoint().y() << " "
174  << hit.exitPoint().z() << " "
175 
176  << "Pixel Pos - X: " << hitPositionPixel.x() << " Y: " << hitPositionPixel.y() << "\n"
177  << "Cart.Cor. - X: " << CMSSWhitPosition.x() << " Y: " << CMSSWhitPosition.y() << "\n"
178  << "Z=0 Pos - X: " << hitPosition.x() << " Y: " << hitPosition.y() << "\n"
179 
180  << "Origin of the template:"
181  << "\n"
182  << "Pixel Pos - X: " << originPixel.x() << " Y: " << originPixel.y() << "\n"
183  << "Cart.Cor. - X: " << origin.x() << " Y: " << origin.y() << "\n"
184  << "\n"
185  << "Entry/Exit:"
186  << "\n"
187  << "Entry - X: " << hit.entryPoint().x() << " Y: " << hit.entryPoint().y()
188  << " Z: " << hit.entryPoint().z() << "\n"
189  << "Exit - X: " << hit.exitPoint().x() << " Y: " << hit.exitPoint().y()
190  << " Z: " << hit.exitPoint().z() << "\n"
191 
192  << "Entry - X Pixel: " << hitEntryPointPixel.x() << " Y Pixel: " << hitEntryPointPixel.y()
193  << "\n"
194  << "Exit - X Pixel: " << hitExitPointPixel.x() << " Y Pixel: " << hitExitPointPixel.y()
195  << "\n"
196 
197  << "row min: " << irow_min << " col min: " << icol_min << "\n";
198 #endif
199 
200  if (!(irow_min <= hitrow_max && irow_max >= hitrow_min && icol_min <= hitcol_max && icol_max >= hitcol_min)) {
201  // The clusters do not have an overlap, hence the hit is NOT reweighted
202  return false;
203  }
204 
205  float cmToMicrons = 10000.f;
206 
207  track[0] = (hitPosition.x() - origin.x()) * cmToMicrons;
208  track[1] = (hitPosition.y() - origin.y()) * cmToMicrons;
209  track[2] = 0.0f; //Middle of sensor is origin for Z-axis
210  track[3] = direction.x();
211  track[4] = direction.y();
212  track[5] = direction.z();
213 
214  array_2d pixrewgt(boost::extents[TXSIZE][TYSIZE]);
215 
216  for (int row = 0; row < TXSIZE; ++row) {
217  for (int col = 0; col < TYSIZE; ++col) {
218  pixrewgt[row][col] = 0;
219  }
220  }
221 
222  for (int row = 0; row < TXSIZE; ++row) {
223  xdouble[row] = topol->isItBigPixelInX(hitPixel.first + row - THX);
224  }
225 
226  for (int col = 0; col < TYSIZE; ++col) {
227  ydouble[col] = topol->isItBigPixelInY(hitPixel.second + col - THY);
228  }
229 
230  // define loop boundaries that will prevent the row and col loops
231  // from going out of physical bounds of the pixel module
232  int rowmin = std::max(0, THX - hitPixel.first);
233  int rowmax = std::min(TXSIZE, topol->nrows() + THX - hitPixel.first);
234  int colmin = std::max(0, THY - hitPixel.second);
235  int colmax = std::min(TYSIZE, topol->ncolumns() + THY - hitPixel.second);
236 
237  for (int row = rowmin; row < rowmax; ++row) {
238  for (int col = colmin; col < colmax; ++col) {
239  //Fill charges into 21x13 Pixel Array with hitPixel in centre
240  pixrewgt[row][col] =
241  hitSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)];
242  }
243  }
244 
245  if (PrintClusters) {
246  std::cout << "Cluster before reweighting: " << std::endl;
247  printCluster(pixrewgt);
248  }
249 
250  int ierr;
251  // for unirradiated: 2nd argument is IDden
252  // for irradiated: 2nd argument is IDnum
253  if (UseReweighting == true) {
254  int ID1 = dbobject_num->getTemplateID(detID);
255  int ID0 = dbobject_den->getTemplateID(detID);
256 
257  if (ID0 == ID1) {
258  return false;
259  }
260  ierr = PixelTempRewgt2D(ID0, ID1, pixrewgt);
261  } else {
262  ierr = PixelTempRewgt2D(IDden, IDden, pixrewgt);
263  }
264  if (ierr != 0) {
265 #ifdef TP_DEBUG
266  LogDebug("PixelDigitizer ") << "Cluster Charge Reweighting did not work properly.";
267 #endif
268  return false;
269  }
270 
271  if (PrintClusters) {
272  std::cout << "Cluster after reweighting: " << std::endl;
273  printCluster(pixrewgt);
274  }
275 
276  for (int row = rowmin; row < rowmax; ++row) {
277  for (int col = colmin; col < colmax; ++col) {
278  float charge = pixrewgt[row][col];
279  if (charge > 0) {
280  chargeAfter += charge;
281  theSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)] +=
282  (boolmakeDigiSimLinks ? SiPixelDigitizerAlgorithm::Amplitude(charge, &hit, hitIndex, tofBin, charge)
283  : SiPixelDigitizerAlgorithm::Amplitude(charge, charge));
284  }
285  }
286  }
287 
288  if (chargeBefore != 0. && chargeAfter == 0.) {
289  return false;
290  }
291 
292  if (PrintClusters) {
293  std::cout << std::endl;
294  std::cout << "Charges (before->after): " << chargeBefore << " -> " << chargeAfter << std::endl;
295  std::cout << "Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 << " %" << std::endl << std::endl;
296  }
297 
298  return true;
299 }
300 
301 // *******************************************************************************************************
309 // *******************************************************************************************************
310 int SiPixelChargeReweightingAlgorithm::PixelTempRewgt2D(int id_in, int id_rewgt, array_2d& cluster) {
311  // Local variables
312  int i, j, k, l, kclose;
313  int nclusx, nclusy, success;
314  float xsize, ysize, q50i, q100i, q50r, q10r, q100r, xhit2D, yhit2D, qclust, dist2, dmin2;
315  float xy_in[BXM2][BYM2], xy_rewgt[BXM2][BYM2], xy_clust[TXSIZE][TYSIZE];
316  int denx_clust[TXSIZE][TYSIZE], deny_clust[TXSIZE][TYSIZE];
317  int goodWeightsUsed, nearbyWeightsUsed, noWeightsUsed;
318  float cotalpha, cotbeta;
319  // success = 0 is returned if everthing is OK
320  success = 0;
321 
322  // Copy the array to remember original charges
323  array_2d clust(cluster);
324 
325  // Take the pixel dimensions from the 2D template
326  templ2D.getid(id_in);
327  xsize = templ2D.xsize();
328  ysize = templ2D.ysize();
329 
330  // Calculate the track angles
331 
332  if (std::abs(track[5]) > 0.f) {
333  cotalpha = track[3] / track[5]; //if track[5] (direction in z) is 0 the hit is not processed by re-weighting
334  cotbeta = track[4] / track[5];
335  } else {
336  LogDebug("Pixel Digitizer") << "Reweighting angle is not good!" << std::endl;
337  return 9; //returned value here indicates that no reweighting was done in this case
338  }
339 
340  // The 2-D templates are defined on a shifted coordinate system wrt the 1D templates
341  if (ydouble[0]) {
342  yhit2D = track[1] - cotbeta * track[2] + ysize;
343  } else {
344  yhit2D = track[1] - cotbeta * track[2] + 0.5f * ysize;
345  }
346  if (xdouble[0]) {
347  xhit2D = track[0] - cotalpha * track[2] + xsize;
348  } else {
349  xhit2D = track[0] - cotalpha * track[2] + 0.5f * xsize;
350  }
351 
352  // Zero the input and output templates
353  for (i = 0; i < BYM2; ++i) {
354  for (j = 0; j < BXM2; ++j) {
355  xy_in[j][i] = 0.f;
356  xy_rewgt[j][i] = 0.f;
357  }
358  }
359 
360  // Next, interpolate the CMSSW template needed to analyze this cluster
361 
362  if (!templ2D.xytemp(id_in, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_in)) {
363  success = 1;
364  }
365  if (success != 0) {
366 #ifdef TP_DEBUG
367  LogDebug("Pixel Digitizer") << "No matching template found" << std::endl;
368 #endif
369  return 2;
370  }
371 
372  if (PrintTemplates) {
373  std::cout << "Template unirrad: " << std::endl;
374  printCluster(xy_in);
375  }
376 
377  q50i = templ2D.s50();
378  //q50i = 0;
379  q100i = 2.f * q50i;
380 
381  // Check that the cluster container is a 13x21 matrix
382 
383  if (cluster.num_dimensions() != 2) {
384  LogWarning("Pixel Digitizer") << "Cluster is not 2-dimensional. Return." << std::endl;
385  return 3;
386  }
387  nclusx = (int)cluster.shape()[0];
388  nclusy = (int)cluster.shape()[1];
389  if (nclusx != TXSIZE || xdouble.size() != TXSIZE) {
390  LogWarning("Pixel Digitizer") << "Sizes in x do not match: nclusx=" << nclusx << " xdoubleSize=" << xdouble.size()
391  << " TXSIZE=" << TXSIZE << ". Return." << std::endl;
392  return 4;
393  }
394  if (nclusy != TYSIZE || ydouble.size() != TYSIZE) {
395  LogWarning("Pixel Digitizer") << "Sizes in y do not match. Return." << std::endl;
396  return 5;
397  }
398 
399  // Sum initial charge in the cluster
400 
401  qclust = 0.f;
402  for (i = 0; i < TYSIZE; ++i) {
403  for (j = 0; j < TXSIZE; ++j) {
404  xy_clust[j][i] = 0.f;
405  denx_clust[j][i] = 0;
406  deny_clust[j][i] = 0;
407  if (cluster[j][i] > q100i) {
408  qclust += cluster[j][i];
409  }
410  }
411  }
412 
413  // Next, interpolate the physical output template needed to reweight
414 
415  if (!templ2D.xytemp(id_rewgt, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_rewgt)) {
416  success = 1;
417  }
418 
419  if (PrintTemplates) {
420  std::cout << "Template irrad: " << std::endl;
421  printCluster(xy_rewgt);
422  }
423 
424  q50r = templ2D.s50();
425  q100r = 2.f * q50r;
426  q10r = 0.2f * q50r;
427 
428  // Find all non-zero denominator pixels in the input template and generate "inside" weights
429 
430  int ntpix = 0;
431  int ncpix = 0;
432  std::vector<int> ytclust;
433  std::vector<int> xtclust;
434  std::vector<int> ycclust;
435  std::vector<int> xcclust;
436  qclust = 0.f;
437  for (i = 0; i < TYSIZE; ++i) {
438  for (j = 0; j < TXSIZE; ++j) {
439  if (xy_in[j + 1][i + 1] > q100i) {
440  ++ntpix;
441  ytclust.push_back(i);
442  xtclust.push_back(j);
443  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[j + 1][i + 1];
444  denx_clust[j][i] = j;
445  deny_clust[j][i] = i;
446  }
447  }
448  }
449 
450  // Find all non-zero numerator pixels not matched to denominator in the output template and generate "inside" weights
451 
452  for (i = 0; i < TYSIZE; ++i) {
453  for (j = 0; j < TXSIZE; ++j) {
454  if (xy_rewgt[j + 1][i + 1] > q10r && xy_clust[j][i] == 0.f && ntpix > 0) {
455  // Search for nearest denominator pixel
456  dmin2 = 10000.f;
457  kclose = 0;
458  for (k = 0; k < ntpix; ++k) {
459  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
460  if (dist2 < dmin2) {
461  dmin2 = dist2;
462  kclose = k;
463  }
464  }
465  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[xtclust[kclose] + 1][ytclust[kclose] + 1];
466  denx_clust[j][i] = xtclust[kclose];
467  deny_clust[j][i] = ytclust[kclose];
468  }
469  }
470  }
471 
472  if (PrintTemplates) {
473  std::cout << "Weights:" << std::endl;
474  printCluster(xy_clust);
475  }
476 
477  // Do the reweighting
478  goodWeightsUsed = 0;
479  nearbyWeightsUsed = 0;
480  noWeightsUsed = 0;
481 
482  for (i = 0; i < TYSIZE; ++i) {
483  for (j = 0; j < TXSIZE; ++j) {
484  if (xy_clust[j][i] > 0.f) {
485  cluster[j][i] = xy_clust[j][i] * clust[denx_clust[j][i]][deny_clust[j][i]];
486  if (cluster[j][i] > q100r) {
487  qclust += cluster[j][i];
488  }
489  if (cluster[j][i] > 0) {
490  goodWeightsUsed++;
491  }
492  } else {
493  if (clust[j][i] > 0.f) {
494  ++ncpix;
495  ycclust.push_back(i);
496  xcclust.push_back(j);
497  }
498  }
499  }
500  }
501 
502  // Now reweight pixels outside of template footprint using closest weights
503 
504  if (ncpix > 0) {
505  for (l = 0; l < ncpix; ++l) {
506  i = ycclust[l];
507  j = xcclust[l];
508  dmin2 = 10000.f;
509  kclose = 0;
510  for (k = 0; k < ntpix; ++k) {
511  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
512  if (dist2 < dmin2) {
513  dmin2 = dist2;
514  kclose = k;
515  }
516  }
517  if (dmin2 < 5.f) {
518  nearbyWeightsUsed++;
519  cluster[j][i] *= xy_clust[xtclust[kclose]][ytclust[kclose]];
520  if (cluster[j][i] > q100r) {
521  qclust += cluster[j][i];
522  }
523  } else {
524  noWeightsUsed++;
525  cluster[j][i] = 0.f;
526  }
527  }
528  }
529 
530  return success;
531 } // PixelTempRewgt2D
532 
534  for (int col = 0; col < TYSIZE; ++col) {
535  for (int row = 0; row < TXSIZE; ++row) {
536  std::cout << std::setw(10) << std::setprecision(0) << std::fixed;
537  std::cout << cluster[row][col];
538  }
539  std::cout << std::endl;
540  }
541  std::cout.copyfmt(std::ios(nullptr));
542 }
543 
545  for (int col = 0; col < BYM2; ++col) {
546  for (int row = 0; row < BXM2; ++row) {
547  std::cout << std::setw(10) << std::setprecision(0) << std::fixed;
548  std::cout << arr[row][col];
549  }
550  std::cout << std::endl;
551  }
552  std::cout.copyfmt(std::ios(nullptr));
553 }
554 
556  for (int col = 0; col < TYSIZE; ++col) {
557  for (int row = 0; row < TXSIZE; ++row) {
558  std::cout << std::setw(10) << std::fixed;
559  std::cout << arr[row][col];
560  }
561  std::cout << std::endl;
562  }
563  std::cout.copyfmt(std::ios(nullptr));
564 }
Log< level::Info, true > LogVerbatim
float xsize()
pixel x-size (microns)
T y() const
Definition: PV2DBase.h:44
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
#define TXSIZE
tuple chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, &quot;NDC&quot;) lumi.SetBorderSize( 0 ) lumi...
T y() const
Definition: PV3DBase.h:60
virtual int nrows() const =0
std::map< int, SiPixelDigitizerAlgorithm::Amplitude, std::less< int > > signal_map_type
bool hitSignalReweight(const PSimHit &hit, std::map< int, float, std::less< int > > &hit_signal, const size_t hitIndex, const unsigned int tofBin, const PixelTopology *topol, uint32_t detID, signal_map_type &theSignal, unsigned short int processType, const bool &boolmakeDigiSimLinks)
const SiPixel2DTemplateDBObject * dbobject_den
const Int_t ysize
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
bool getData(T &iHolder) const
Definition: EventSetup.h:128
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
Local3DPoint localPosition() const
Definition: PSimHit.h:52
virtual bool isItBigPixelInX(int ixbin) const =0
T z() const
Definition: PV3DBase.h:61
virtual bool isItBigPixelInY(int iybin) const =0
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define BYM2
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
T min(T a, T b)
Definition: MathUtil.h:58
#define THY
const SiPixel2DTemplateDBObject * dbobject_num
int PixelTempRewgt2D(int id_gen, int id_rewgt, array_2d &cluster)
#define TYSIZE
short getTemplateID(const uint32_t &detid) const
SiPixelChargeReweightingAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
unsigned short processType() const
Definition: PSimHit.h:120
float ysize()
pixel y-size (microns)
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_den_token_
int particleType() const
Definition: PSimHit.h:89
tuple cout
Definition: gather_cfg.py:144
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore2D > &pixelTemp, std::string dir="CalibTracker/SiPixelESProducers/data/")
Log< level::Warning, false > LogWarning
int col
Definition: cuy.py:1009
const Int_t xsize
T x() const
Definition: PV2DBase.h:43
T x() const
Definition: PV3DBase.h:59
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
#define LogDebug(id)