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 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  for (int row = 0; row < TXSIZE; ++row) {
231  for (int col = 0; col < TYSIZE; ++col) {
232  //Fill charges into 21x13 Pixel Array with hitPixel in centre
233  pixrewgt[row][col] =
234  hitSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)];
235  }
236  }
237 
238  if (PrintClusters) {
239  std::cout << "Cluster before reweighting: " << std::endl;
240  printCluster(pixrewgt);
241  }
242 
243  int ierr;
244  // for unirradiated: 2nd argument is IDden
245  // for irradiated: 2nd argument is IDnum
246  if (UseReweighting == true) {
247  int ID1 = dbobject_num->getTemplateID(detID);
248  int ID0 = dbobject_den->getTemplateID(detID);
249 
250  if (ID0 == ID1) {
251  return false;
252  }
253  ierr = PixelTempRewgt2D(ID0, ID1, pixrewgt);
254  } else {
255  ierr = PixelTempRewgt2D(IDden, IDden, pixrewgt);
256  }
257  if (ierr != 0) {
258 #ifdef TP_DEBUG
259  LogDebug("PixelDigitizer ") << "Cluster Charge Reweighting did not work properly.";
260 #endif
261  return false;
262  }
263 
264  if (PrintClusters) {
265  std::cout << "Cluster after reweighting: " << std::endl;
266  printCluster(pixrewgt);
267  }
268 
269  for (int row = 0; row < TXSIZE; ++row) {
270  for (int col = 0; col < TYSIZE; ++col) {
271  float charge = 0;
272  charge = pixrewgt[row][col];
273  if ((hitPixel.first + row - THX) >= 0 && (hitPixel.first + row - THX) < topol->nrows() &&
274  (hitPixel.second + col - THY) >= 0 && (hitPixel.second + col - THY) < topol->ncolumns() && charge > 0) {
275  chargeAfter += charge;
276  theSignal[PixelDigi::pixelToChannel(hitPixel.first + row - THX, hitPixel.second + col - THY)] +=
277  (boolmakeDigiSimLinks ? SiPixelDigitizerAlgorithm::Amplitude(charge, &hit, hitIndex, tofBin, charge)
279  }
280  }
281  }
282 
283  if (chargeBefore != 0. && chargeAfter == 0.) {
284  return false;
285  }
286 
287  if (PrintClusters) {
288  std::cout << std::endl;
289  std::cout << "Charges (before->after): " << chargeBefore << " -> " << chargeAfter << std::endl;
290  std::cout << "Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 << " %" << std::endl << std::endl;
291  }
292 
293  return true;
294 }
295 
296 // *******************************************************************************************************
304 // *******************************************************************************************************
305 int SiPixelChargeReweightingAlgorithm::PixelTempRewgt2D(int id_in, int id_rewgt, array_2d& cluster) {
306  // Local variables
307  int i, j, k, l, kclose;
308  int nclusx, nclusy, success;
309  float xsize, ysize, q50i, q100i, q50r, q10r, q100r, xhit2D, yhit2D, qclust, dist2, dmin2;
310  float xy_in[BXM2][BYM2], xy_rewgt[BXM2][BYM2], xy_clust[TXSIZE][TYSIZE];
311  int denx_clust[TXSIZE][TYSIZE], deny_clust[TXSIZE][TYSIZE];
312  int goodWeightsUsed, nearbyWeightsUsed, noWeightsUsed;
313  float cotalpha, cotbeta;
314  // success = 0 is returned if everthing is OK
315  success = 0;
316 
317  // Copy the array to remember original charges
318  array_2d clust(cluster);
319 
320  // Take the pixel dimensions from the 2D template
321  templ2D.getid(id_in);
322  xsize = templ2D.xsize();
323  ysize = templ2D.ysize();
324 
325  // Calculate the track angles
326 
327  if (std::abs(track[5]) > 0.f) {
328  cotalpha = track[3] / track[5]; //if track[5] (direction in z) is 0 the hit is not processed by re-weighting
329  cotbeta = track[4] / track[5];
330  } else {
331  LogDebug("Pixel Digitizer") << "Reweighting angle is not good!" << std::endl;
332  return 9; //returned value here indicates that no reweighting was done in this case
333  }
334 
335  // The 2-D templates are defined on a shifted coordinate system wrt the 1D templates
336  if (ydouble[0]) {
337  yhit2D = track[1] - cotbeta * track[2] + ysize;
338  } else {
339  yhit2D = track[1] - cotbeta * track[2] + 0.5f * ysize;
340  }
341  if (xdouble[0]) {
342  xhit2D = track[0] - cotalpha * track[2] + xsize;
343  } else {
344  xhit2D = track[0] - cotalpha * track[2] + 0.5f * xsize;
345  }
346 
347  // Zero the input and output templates
348  for (i = 0; i < BYM2; ++i) {
349  for (j = 0; j < BXM2; ++j) {
350  xy_in[j][i] = 0.f;
351  xy_rewgt[j][i] = 0.f;
352  }
353  }
354 
355  // Next, interpolate the CMSSW template needed to analyze this cluster
356 
357  if (!templ2D.xytemp(id_in, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_in)) {
358  success = 1;
359  }
360  if (success != 0) {
361 #ifdef TP_DEBUG
362  LogDebug("Pixel Digitizer") << "No matching template found" << std::endl;
363 #endif
364  return 2;
365  }
366 
367  if (PrintTemplates) {
368  std::cout << "Template unirrad: " << std::endl;
369  printCluster(xy_in);
370  }
371 
372  q50i = templ2D.s50();
373  //q50i = 0;
374  q100i = 2.f * q50i;
375 
376  // Check that the cluster container is a 13x21 matrix
377 
378  if (cluster.num_dimensions() != 2) {
379  LogWarning("Pixel Digitizer") << "Cluster is not 2-dimensional. Return." << std::endl;
380  return 3;
381  }
382  nclusx = (int)cluster.shape()[0];
383  nclusy = (int)cluster.shape()[1];
384  if (nclusx != TXSIZE || xdouble.size() != TXSIZE) {
385  LogWarning("Pixel Digitizer") << "Sizes in x do not match: nclusx=" << nclusx << " xdoubleSize=" << xdouble.size()
386  << " TXSIZE=" << TXSIZE << ". Return." << std::endl;
387  return 4;
388  }
389  if (nclusy != TYSIZE || ydouble.size() != TYSIZE) {
390  LogWarning("Pixel Digitizer") << "Sizes in y do not match. Return." << std::endl;
391  return 5;
392  }
393 
394  // Sum initial charge in the cluster
395 
396  qclust = 0.f;
397  for (i = 0; i < TYSIZE; ++i) {
398  for (j = 0; j < TXSIZE; ++j) {
399  xy_clust[j][i] = 0.f;
400  denx_clust[j][i] = 0;
401  deny_clust[j][i] = 0;
402  if (cluster[j][i] > q100i) {
403  qclust += cluster[j][i];
404  }
405  }
406  }
407 
408  // Next, interpolate the physical output template needed to reweight
409 
410  if (!templ2D.xytemp(id_rewgt, cotalpha, cotbeta, xhit2D, yhit2D, ydouble, xdouble, xy_rewgt)) {
411  success = 1;
412  }
413 
414  if (PrintTemplates) {
415  std::cout << "Template irrad: " << std::endl;
416  printCluster(xy_rewgt);
417  }
418 
419  q50r = templ2D.s50();
420  q100r = 2.f * q50r;
421  q10r = 0.2f * q50r;
422 
423  // Find all non-zero denominator pixels in the input template and generate "inside" weights
424 
425  int ntpix = 0;
426  int ncpix = 0;
427  std::vector<int> ytclust;
428  std::vector<int> xtclust;
429  std::vector<int> ycclust;
430  std::vector<int> xcclust;
431  qclust = 0.f;
432  for (i = 0; i < TYSIZE; ++i) {
433  for (j = 0; j < TXSIZE; ++j) {
434  if (xy_in[j + 1][i + 1] > q100i) {
435  ++ntpix;
436  ytclust.push_back(i);
437  xtclust.push_back(j);
438  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[j + 1][i + 1];
439  denx_clust[j][i] = j;
440  deny_clust[j][i] = i;
441  }
442  }
443  }
444 
445  // Find all non-zero numerator pixels not matched to denominator in the output template and generate "inside" weights
446 
447  for (i = 0; i < TYSIZE; ++i) {
448  for (j = 0; j < TXSIZE; ++j) {
449  if (xy_rewgt[j + 1][i + 1] > q10r && xy_clust[j][i] == 0.f && ntpix > 0) {
450  // Search for nearest denominator pixel
451  dmin2 = 10000.f;
452  kclose = 0;
453  for (k = 0; k < ntpix; ++k) {
454  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
455  if (dist2 < dmin2) {
456  dmin2 = dist2;
457  kclose = k;
458  }
459  }
460  xy_clust[j][i] = xy_rewgt[j + 1][i + 1] / xy_in[xtclust[kclose] + 1][ytclust[kclose] + 1];
461  denx_clust[j][i] = xtclust[kclose];
462  deny_clust[j][i] = ytclust[kclose];
463  }
464  }
465  }
466 
467  if (PrintTemplates) {
468  std::cout << "Weights:" << std::endl;
469  printCluster(xy_clust);
470  }
471 
472  // Do the reweighting
473  goodWeightsUsed = 0;
474  nearbyWeightsUsed = 0;
475  noWeightsUsed = 0;
476 
477  for (i = 0; i < TYSIZE; ++i) {
478  for (j = 0; j < TXSIZE; ++j) {
479  if (xy_clust[j][i] > 0.f) {
480  cluster[j][i] = xy_clust[j][i] * clust[denx_clust[j][i]][deny_clust[j][i]];
481  if (cluster[j][i] > q100r) {
482  qclust += cluster[j][i];
483  }
484  if (cluster[j][i] > 0) {
485  goodWeightsUsed++;
486  }
487  } else {
488  if (clust[j][i] > 0.f) {
489  ++ncpix;
490  ycclust.push_back(i);
491  xcclust.push_back(j);
492  }
493  }
494  }
495  }
496 
497  // Now reweight pixels outside of template footprint using closest weights
498 
499  if (ncpix > 0) {
500  for (l = 0; l < ncpix; ++l) {
501  i = ycclust[l];
502  j = xcclust[l];
503  dmin2 = 10000.f;
504  kclose = 0;
505  for (k = 0; k < ntpix; ++k) {
506  dist2 = (i - ytclust[k]) * (i - ytclust[k]) + 0.44444f * (j - xtclust[k]) * (j - xtclust[k]);
507  if (dist2 < dmin2) {
508  dmin2 = dist2;
509  kclose = k;
510  }
511  }
512  if (dmin2 < 5.f) {
513  nearbyWeightsUsed++;
514  cluster[j][i] *= xy_clust[xtclust[kclose]][ytclust[kclose]];
515  if (cluster[j][i] > q100r) {
516  qclust += cluster[j][i];
517  }
518  } else {
519  noWeightsUsed++;
520  cluster[j][i] = 0.f;
521  }
522  }
523  }
524 
525  return success;
526 } // PixelTempRewgt2D
527 
529  for (int col = 0; col < TYSIZE; ++col) {
530  for (int row = 0; row < TXSIZE; ++row) {
531  std::cout << std::setw(10) << std::setprecision(0) << std::fixed;
532  std::cout << cluster[row][col];
533  }
534  std::cout << std::endl;
535  }
536  std::cout.copyfmt(std::ios(nullptr));
537 }
538 
540  for (int col = 0; col < BYM2; ++col) {
541  for (int row = 0; row < BXM2; ++row) {
542  std::cout << std::setw(10) << std::setprecision(0) << std::fixed;
543  std::cout << arr[row][col];
544  }
545  std::cout << std::endl;
546  }
547  std::cout.copyfmt(std::ios(nullptr));
548 }
549 
551  for (int col = 0; col < TYSIZE; ++col) {
552  for (int row = 0; row < TXSIZE; ++row) {
553  std::cout << std::setw(10) << std::fixed;
554  std::cout << arr[row][col];
555  }
556  std::cout << std::endl;
557  }
558  std::cout.copyfmt(std::ios(nullptr));
559 }
Vector3DBase< float, LocalTag >
SiPixelTemplate2D::xsize
float xsize()
pixel x-size (microns)
Definition: SiPixelTemplate2D.h:262
Point2DBase
Definition: Point2DBase.h:9
alignBH_cfg.fixed
fixed
Definition: alignBH_cfg.py:54
electrons_cff.bool
bool
Definition: electrons_cff.py:366
SiPixelTemplate2D::xytemp
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)
Definition: SiPixelTemplate2D.cc:989
mps_fire.i
i
Definition: mps_fire.py:428
edm::ESInputTag
Definition: ESInputTag.h:87
PixelSubdetector.h
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11724
SiPixel2DTemplateDBObject.h
MessageLogger.h
SiPixelChargeReweightingAlgorithm::dbobject_num
const SiPixel2DTemplateDBObject * dbobject_num
Definition: SiPixelChargeReweightingAlgorithm.h:62
hit::y
double y
Definition: SiStripHitEffFromCalibTree.cc:90
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
SiPixelSimParameters_cfi.PrintTemplates
PrintTemplates
Definition: SiPixelSimParameters_cfi.py:53
SiPixelChargeReweightingAlgorithm::signal_map_type
std::map< int, SiPixelDigitizerAlgorithm::Amplitude, std::less< int > > signal_map_type
Definition: SiPixelChargeReweightingAlgorithm.h:25
edm
HLT enums.
Definition: AlignableModifier.h:19
SiPixelChargeReweightingAlgorithm::templ2D
SiPixelTemplate2D templ2D
Definition: SiPixelChargeReweightingAlgorithm.h:47
cuy.col
col
Definition: cuy.py:1009
gather_cfg.cout
cout
Definition: gather_cfg.py:144
PSimHitContainer.h
THX
#define THX
Definition: SiPixelTemplateDefs.h:31
SiPixelTemplate2D::getid
bool getid(int id)
Definition: SiPixelTemplate2D.cc:536
SiPixelChargeReweightingAlgorithm::PixelTempRewgt2D
int PixelTempRewgt2D(int id_gen, int id_rewgt, array_2d &cluster)
Definition: SiPixelChargeReweightingAlgorithm.cc:305
GlobalPixel.h
SiPixelSimParameters_cfi.PrintClusters
PrintClusters
Definition: SiPixelSimParameters_cfi.py:52
SiPixelChargeReweightingAlgorithm::init
void init(const edm::EventSetup &es)
Definition: SiPixelChargeReweightingAlgorithm.cc:39
ysize
const Int_t ysize
Definition: trackSplitPlot.h:43
Topology::localPosition
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
SiPixelChargeReweightingAlgorithm::UseReweighting
const bool UseReweighting
Definition: SiPixelChargeReweightingAlgorithm.h:53
PixelDigi.h
xsize
const Int_t xsize
Definition: trackSplitPlot.h:42
edm::ConsumesCollector::esConsumes
auto esConsumes()
Definition: ConsumesCollector.h:97
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
GaussianTailNoiseGenerator.h
TXSIZE
#define TXSIZE
Definition: SiPixelTemplateDefs.h:30
hit::x
double x
Definition: SiStripHitEffFromCalibTree.cc:89
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
SiG4UniversalFluctuation.h
PixelTopology::isItBigPixelInX
virtual bool isItBigPixelInX(int ixbin) const =0
SiPixelChargeReweightingAlgorithm::hitSignalReweight
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)
Definition: SiPixelChargeReweightingAlgorithm.cc:84
SiPixelTemplate2D::s50
float s50()
1/2 of the pixel threshold signal in adc units
Definition: SiPixelTemplate2D.h:206
BYM2
#define BYM2
Definition: SiPixelTemplateDefs.h:28
SiPixelChargeReweightingAlgorithm::xdouble
std::vector< bool > xdouble
Definition: SiPixelChargeReweightingAlgorithm.h:48
SiPixelChargeReweightingAlgorithm::PrintTemplates
const bool PrintTemplates
Definition: SiPixelChargeReweightingAlgorithm.h:55
sipixelobjects
Definition: CablingPathToDetUnit.h:4
BXM2
#define BXM2
Definition: SiPixelTemplateDefs.h:36
SiPixelChargeReweightingAlgorithm::ydouble
std::vector< bool > ydouble
Definition: SiPixelChargeReweightingAlgorithm.h:49
summarizeEdmComparisonLogfiles.success
success
Definition: summarizeEdmComparisonLogfiles.py:114
PixelTopology::ncolumns
virtual int ncolumns() const =0
dqmdumpme.k
k
Definition: dqmdumpme.py:60
PixelDigi::pixelToChannel
static int pixelToChannel(int row, int col)
Definition: PixelDigi.h:75
Point3DBase< float, LocalTag >
PixelTopology
Definition: PixelTopology.h:10
SiPixelChargeReweightingAlgorithm::SiPixel2DTemp_num_token_
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_num_token_
Definition: SiPixelChargeReweightingAlgorithm.h:60
SiPixelChargeReweightingAlgorithm::track
std::vector< float > track
Definition: SiPixelChargeReweightingAlgorithm.h:50
hit::z
double z
Definition: SiStripHitEffFromCalibTree.cc:91
SiPixelChargeReweightingAlgorithm::SiPixel2DTemp_den_token_
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_den_token_
Definition: SiPixelChargeReweightingAlgorithm.h:59
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
MeasurementPoint
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
Definition: MeasurementPoint.h:12
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
PixelROC.h
SiPixel2DTemplateDBObject::getTemplateID
short getTemplateID(const uint32_t &detid) const
Definition: SiPixel2DTemplateDBObject.h:90
PV2DBase::y
T y() const
Definition: PV2DBase.h:44
PV2DBase::x
T x() const
Definition: PV2DBase.h:43
PixelTopology::isItBigPixelInY
virtual bool isItBigPixelInY(int iybin) const =0
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
createfilelist.int
int
Definition: createfilelist.py:10
TYSIZE
#define TYSIZE
Definition: SiPixelTemplateDefs.h:21
LocalPixel.h
SiPixelChargeReweightingAlgorithm::printCluster
void printCluster(array_2d &cluster)
Definition: SiPixelChargeReweightingAlgorithm.cc:528
edm::EventSetup
Definition: EventSetup.h:58
SiPixelChargeReweightingAlgorithm::array_2d
boost::multi_array< float, 2 > array_2d
Definition: SiPixelChargeReweightingAlgorithm.h:44
Topology::measurementPosition
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
SiPixelChargeReweightingAlgorithm::IDden
int IDden
Definition: SiPixelChargeReweightingAlgorithm.h:51
PixelIndices.h
SiPixelChargeReweightingAlgorithm::~SiPixelChargeReweightingAlgorithm
~SiPixelChargeReweightingAlgorithm()
Definition: SiPixelChargeReweightingAlgorithm.cc:78
SiPixelTemplate2D::ysize
float ysize()
pixel y-size (microns)
Definition: SiPixelTemplate2D.h:263
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
THY
#define THY
Definition: SiPixelTemplateDefs.h:22
SiPixelChargeReweightingAlgorithm.h
LaserClient_cfi.Amplitude
Amplitude
Definition: LaserClient_cfi.py:39
DetId.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
officialStyle.chan
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi....
Definition: officialStyle.py:106
ConsumesCollector.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
genParticles_cff.map
map
Definition: genParticles_cff.py:11
ParameterSet.h
SiPixelSimParameters_cfi.UseReweighting
UseReweighting
Definition: SiPixelSimParameters_cfi.py:51
PSimHit
Definition: PSimHit.h:15
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
SiPixelChargeReweightingAlgorithm::PrintClusters
const bool PrintClusters
Definition: SiPixelChargeReweightingAlgorithm.h:54
PixelTopology::nrows
virtual int nrows() const =0
SiPixelTemplate2D::pushfile
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore2D > &pixelTemp, std::string dir="CalibTracker/SiPixelESProducers/data/")
Definition: SiPixelTemplate2D.cc:62
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
hit
Definition: SiStripHitEffFromCalibTree.cc:88
SiPixelChargeReweightingAlgorithm::dbobject_den
const SiPixel2DTemplateDBObject * dbobject_den
Definition: SiPixelChargeReweightingAlgorithm.h:61
PixelDigi::channelToPixel
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
SiPixelChargeReweightingAlgorithm::SiPixelChargeReweightingAlgorithm
SiPixelChargeReweightingAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
Definition: SiPixelChargeReweightingAlgorithm.cc:56