1 #ifndef SimTracker_Common_SiPixelChargeReweightingAlgorithm_h 2 #define SimTracker_Common_SiPixelChargeReweightingAlgorithm_h 13 #include <type_traits> 14 #include <gsl/gsl_sf_erf.h> 16 #include "boost/multi_array.hpp" 60 template <
class AmplitudeType,
typename SignalType>
62 std::map<
int,
float, std::less<int> >& hit_signal,
63 const size_t hitIndex,
64 const size_t hitIndex4CR,
65 const unsigned int tofBin,
68 SignalType& theSignal,
69 unsigned short int processType,
70 const bool& boolmakeDigiSimLinks);
72 template <
class AmplitudeType,
typename SignalType>
74 std::vector<PixelDigi>& digis,
76 SignalType& theNewDigiSignal,
78 CLHEP::HepRandomEngine* engine);
126 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Init SiPixelChargeReweightingAlgorithm";
136 templ2D(templateStores_),
139 IDnum(conf.exists(
"TemplateIDnumerator") ? conf.getParameter<
int>(
"TemplateIDnumerator") : 0),
140 IDden(conf.exists(
"TemplateIDdenominator") ? conf.getParameter<
int>(
"TemplateIDdenominator") : 0),
143 applyLateReweighting_(conf.exists(
"applyLateReweighting") ? conf.getParameter<
bool>(
"applyLateReweighting")
151 LogDebug(
"SiPixelChargeReweightingAlgorithm")
152 <<
"SiPixelChargeReweightingAlgorithm constructed" 158 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"SiPixelChargeReweightingAlgorithm deleted";
162 template <
class AmplitudeType,
typename SignalType>
164 std::map<
int,
float, std::less<int> >& hit_signal,
165 const size_t hitIndex,
166 const size_t hitIndex4CR,
167 const unsigned int tofBin,
170 SignalType& theSignal,
171 unsigned short int processType,
172 const bool& boolmakeDigiSimLinks) {
173 int irow_min = topol->
nrows();
178 float chargeBefore = 0;
179 float chargeAfter = 0;
180 SignalType hitSignal;
183 for (
std::map<
int,
float, std::less<int> >::const_iterator im = hit_signal.begin(); im != hit_signal.end(); ++im) {
184 int chan = (*im).first;
186 if constexpr (std::is_same_v<AmplitudeType, digitizerUtility::Ph2Amplitude>) {
187 hitSignal[
chan] += AmplitudeType((*im).second, &
hit, (*im).second, 0., hitIndex, tofBin);
190 (boolmakeDigiSimLinks ? AmplitudeType((*im).second, &
hit, hitIndex, hitIndex4CR, tofBin, (*im).second)
191 : AmplitudeType((*im).second, (*im).second));
193 chargeBefore += (*im).second;
195 if (pixelWithCharge.first < irow_min)
196 irow_min = pixelWithCharge.first;
197 if (pixelWithCharge.first > irow_max)
198 irow_max = pixelWithCharge.first;
199 if (pixelWithCharge.second < icol_min)
200 icol_min = pixelWithCharge.second;
201 if (pixelWithCharge.second > icol_max)
202 icol_max = pixelWithCharge.second;
207 float trajectoryScaleToPosition =
std::abs(hitEntryPoint.
z() / direction.
z());
209 LocalPoint hitPosition = hitEntryPoint + trajectoryScaleToPosition * direction;
212 std::pair<int, int> hitPixel =
213 std::pair<int, int>(
int(floor(hitPositionPixel.
x())),
int(floor(hitPositionPixel.
y())));
220 std::pair<int, int> entryPixel =
221 std::pair<int, int>(
int(floor(hitEntryPointPixel.
x())),
int(floor(hitEntryPointPixel.
y())));
222 std::pair<int, int> exitPixel =
223 std::pair<int, int>(
int(floor(hitExitPointPixel.
x())),
int(floor(hitExitPointPixel.
y())));
225 int hitcol_min, hitcol_max, hitrow_min, hitrow_max;
226 if (entryPixel.first > exitPixel.first) {
227 hitrow_min = exitPixel.first;
228 hitrow_max = entryPixel.first;
230 hitrow_min = entryPixel.first;
231 hitrow_max = exitPixel.first;
234 if (entryPixel.second > exitPixel.second) {
235 hitcol_min = exitPixel.second;
236 hitcol_max = entryPixel.second;
238 hitcol_min = entryPixel.second;
239 hitcol_max = exitPixel.second;
243 LogDebug(
"SiPixelChargeReweightingAlgorithm")
245 <<
"Particle ID is: " <<
hit.particleType() <<
"\n" 246 <<
"Process type: " <<
hit.processType() <<
"\n" 249 <<
"Hit entry x/y/z: " <<
hit.entryPoint().
x() <<
" " <<
hit.entryPoint().
y() <<
" " <<
hit.entryPoint().
z()
251 <<
"Hit exit x/y/z: " <<
hit.exitPoint().
x() <<
" " <<
hit.exitPoint().
y() <<
" " <<
hit.exitPoint().
z() <<
" " 253 <<
"Pixel Pos - X: " << hitPositionPixel.
x() <<
" Y: " << hitPositionPixel.
y() <<
"\n" 254 <<
"Cart.Cor. - X: " << CMSSWhitPosition.
x() <<
" Y: " << CMSSWhitPosition.
y() <<
"\n" 255 <<
"Z=0 Pos - X: " << hitPosition.x() <<
" Y: " << hitPosition.y() <<
"\n" 257 <<
"Origin of the template:" 259 <<
"Pixel Pos - X: " << originPixel.x() <<
" Y: " << originPixel.y() <<
"\n" 260 <<
"Cart.Cor. - X: " << origin.
x() <<
" Y: " << origin.
y() <<
"\n" 264 <<
"Entry - X: " <<
hit.entryPoint().
x() <<
" Y: " <<
hit.entryPoint().
y() <<
" Z: " <<
hit.entryPoint().
z()
266 <<
"Exit - X: " <<
hit.exitPoint().
x() <<
" Y: " <<
hit.exitPoint().
y() <<
" Z: " <<
hit.exitPoint().
z() <<
"\n" 268 <<
"Entry - X Pixel: " << hitEntryPointPixel.
x() <<
" Y Pixel: " << hitEntryPointPixel.
y() <<
"\n" 269 <<
"Exit - X Pixel: " << hitExitPointPixel.
x() <<
" Y Pixel: " << hitExitPointPixel.
y() <<
"\n" 271 <<
"row min: " << irow_min <<
" col min: " << icol_min <<
"\n";
273 if (!(irow_min <= hitrow_max && irow_max >= hitrow_min && icol_min <= hitcol_max && icol_max >= hitcol_min)) {
275 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Outside the row/col boundary, exit charge reweighting";
282 track[3] = direction.x();
283 track[4] = direction.y();
284 track[5] = direction.z();
286 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Init array of size x = " <<
TXSIZE <<
" and y = " <<
TYSIZE;
289 for (
int row = 0; row <
TXSIZE; ++row) {
291 pixrewgt[row][
col] = 0;
295 for (
int row = 0; row <
TXSIZE; ++row) {
310 for (
int row = rowmin; row < rowmax; ++row) {
311 for (
int col = colmin;
col < colmax; ++
col) {
319 LogDebug(
"SiPixelChargeReweightingAlgorithm ") <<
"Cluster before reweighting: ";
338 LogDebug(
"SiPixelChargeReweightingAlgorithm ") <<
"Cluster Charge Reweighting did not work properly.";
343 LogDebug(
"SiPixelChargeReweightingAlgorithm ") <<
"Cluster after reweighting: ";
347 for (
int row = rowmin; row < rowmax; ++row) {
348 for (
int col = colmin;
col < colmax; ++
col) {
352 if constexpr (std::is_same_v<AmplitudeType, digitizerUtility::Ph2Amplitude>) {
357 (boolmakeDigiSimLinks ? AmplitudeType(
charge, &
hit, hitIndex, hitIndex4CR, tofBin,
charge)
364 if (chargeBefore != 0. && chargeAfter == 0.) {
369 LogDebug(
"SiPixelChargeReweightingAlgorithm ")
370 <<
"Charges (before->after): " << chargeBefore <<
" -> " << chargeAfter;
371 LogDebug(
"SiPixelChargeReweightingAlgorithm ")
372 <<
"Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 <<
" % \n";
389 int i,
j,
k,
l, kclose;
391 float xsize,
ysize, q50i, q100i, q50r, q10r, q100r, xhit2D, yhit2D, qclust, dist2, dmin2;
394 int goodWeightsUsed, nearbyWeightsUsed, noWeightsUsed;
395 float cotalpha, cotbeta;
413 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Reweighting angle is not good! \n";
433 xy_rewgt[
j][
i] = 0.f;
443 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"No matching template found \n";
448 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Template unirrad: \n";
458 if (cluster.num_dimensions() != 2) {
459 edm::LogWarning(
"SiPixelChargeReweightingAlgorithm") <<
"Cluster is not 2-dimensional. Return. \n";
462 nclusx = (
int)cluster.shape()[0];
463 nclusy = (
int)cluster.shape()[1];
466 <<
"Sizes in x do not match: nclusx=" << nclusx <<
" xdoubleSize=" <<
xdouble.size() <<
" TXSIZE=" <<
TXSIZE 471 edm::LogWarning(
"SiPixelChargeReweightingAlgorithm") <<
"Sizes in y do not match. Return. \n";
480 xy_clust[
j][
i] = 0.f;
481 denx_clust[
j][
i] = 0;
482 deny_clust[
j][
i] = 0;
483 if (cluster[
j][
i] > q100i) {
484 qclust += cluster[
j][
i];
496 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Template irrad: \n";
508 std::vector<int> ytclust;
509 std::vector<int> xtclust;
510 std::vector<int> ycclust;
511 std::vector<int> xcclust;
515 if (xy_in[
j + 1][
i + 1] > q100i) {
517 ytclust.push_back(
i);
518 xtclust.push_back(
j);
519 xy_clust[
j][
i] = xy_rewgt[
j + 1][
i + 1] / xy_in[
j + 1][
i + 1];
520 denx_clust[
j][
i] =
j;
521 deny_clust[
j][
i] =
i;
530 if (xy_rewgt[
j + 1][
i + 1] > q10r && xy_clust[
j][
i] == 0.
f && ntpix > 0) {
534 for (
k = 0;
k < ntpix; ++
k) {
535 dist2 = (
i - ytclust[
k]) * (
i - ytclust[
k]) + 0.44444f * (
j - xtclust[
k]) * (
j - xtclust[
k]);
541 xy_clust[
j][
i] = xy_rewgt[
j + 1][
i + 1] / xy_in[xtclust[kclose] + 1][ytclust[kclose] + 1];
542 denx_clust[
j][
i] = xtclust[kclose];
543 deny_clust[
j][
i] = ytclust[kclose];
549 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Weights: \n";
555 nearbyWeightsUsed = 0;
560 if (xy_clust[
j][
i] > 0.
f) {
561 cluster[
j][
i] = xy_clust[
j][
i] * clust[denx_clust[
j][
i]][deny_clust[
j][
i]];
562 if (cluster[
j][
i] > q100r) {
563 qclust += cluster[
j][
i];
565 if (cluster[
j][
i] > 0) {
569 if (clust[
j][
i] > 0.
f) {
571 ycclust.push_back(
i);
572 xcclust.push_back(
j);
581 for (
l = 0;
l < ncpix; ++
l) {
586 for (
k = 0;
k < ntpix; ++
k) {
587 dist2 = (
i - ytclust[
k]) * (
i - ytclust[
k]) + 0.44444f * (
j - xtclust[
k]) * (
j - xtclust[
k]);
595 cluster[
j][
i] *= xy_clust[xtclust[kclose]][ytclust[kclose]];
596 if (cluster[
j][
i] > q100r) {
597 qclust += cluster[
j][
i];
611 for (
int row = 0; row <
TXSIZE; ++row) {
612 LogDebug(
"SiPixelChargeReweightingAlgorithm") << cluster[row][
col];
614 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"\n";
620 for (
int row = 0; row <
BXM2; ++row) {
621 LogDebug(
"SiPixelChargeReweightingAlgorithm") << arr[row][
col];
623 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"\n";
629 for (
int row = 0; row <
TXSIZE; ++row) {
630 LogDebug(
"SiPixelChargeReweightingAlgorithm") << arr[row][
col];
632 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"\n";
636 template <
class AmplitudeType,
typename SignalType>
638 std::vector<PixelDigi>& digis,
640 SignalType& theNewDigiSignal,
642 CLHEP::HepRandomEngine* engine) {
647 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ******************************** \n";
648 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ******************************** \n";
649 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ******************************** \n";
650 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ***** INCONSISTENCY !!! ***** \n";
652 <<
" applyLateReweighting_ and UseReweighting can not be true at the same time for PU ! \n";
653 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ---> DO NOT APPLY CHARGE REWEIGHTING TWICE !!! \n";
654 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ******************************** \n";
655 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
" ******************************** \n";
659 SignalType theDigiSignal;
661 int irow_min = topol->
nrows();
666 float chargeBefore = 0;
667 float chargeAfter = 0;
670 std::vector<PixelDigi>::const_iterator loopDigi;
671 for (loopDigi = digis.begin(); loopDigi != digis.end(); ++loopDigi) {
672 unsigned int chan = loopDigi->channel();
674 float corresponding_charge = loopDigi->adc();
675 if constexpr (std::is_same_v<AmplitudeType, digitizerUtility::Ph2Amplitude>) {
677 <<
"Phase-2 late charge reweighting not supported (not sure we need it at all)";
680 theDigiSignal[
chan] += AmplitudeType(corresponding_charge, corresponding_charge);
682 chargeBefore += corresponding_charge;
683 if (loopDigi->row() < irow_min)
684 irow_min = loopDigi->row();
685 if (loopDigi->row() > irow_max)
686 irow_max = loopDigi->row();
687 if (loopDigi->column() < icol_min)
688 icol_min = loopDigi->column();
689 if (loopDigi->column() > icol_max)
690 icol_max = loopDigi->column();
697 float trajectoryScaleToPosition =
std::abs(hitEntryPoint.
z() / direction.
z());
698 LocalPoint hitPosition = hitEntryPoint + trajectoryScaleToPosition * direction;
702 LocalPoint hitLocalPosition = hitEntryPoint + 0.5 * direction;
704 std::pair<int, int> hitPixel =
705 std::pair<int, int>(
int(floor(hitPositionPixel.
x())),
int(floor(hitPositionPixel.
y())));
712 std::pair<int, int> entryPixel =
713 std::pair<int, int>(
int(floor(hitEntryPointPixel.
x())),
int(floor(hitEntryPointPixel.
y())));
714 std::pair<int, int> exitPixel =
715 std::pair<int, int>(
int(floor(hitExitPointPixel.
x())),
int(floor(hitExitPointPixel.
y())));
717 int hitcol_min, hitcol_max, hitrow_min, hitrow_max;
718 if (entryPixel.first > exitPixel.first) {
719 hitrow_min = exitPixel.first;
720 hitrow_max = entryPixel.first;
722 hitrow_min = entryPixel.first;
723 hitrow_max = exitPixel.first;
726 if (entryPixel.second > exitPixel.second) {
727 hitcol_min = exitPixel.second;
728 hitcol_max = entryPixel.second;
730 hitcol_min = entryPixel.second;
731 hitcol_max = exitPixel.second;
734 if (!(irow_min <= hitrow_max && irow_max >= hitrow_min && icol_min <= hitcol_max && icol_max >= hitcol_min)) {
742 track[3] = direction.x();
743 track[4] = direction.y();
744 track[5] = direction.z();
756 for (
int row = 0; row <
TXSIZE; ++row) {
764 for (
int row = 0; row <
TXSIZE; ++row) {
773 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
" Cluster before reweighting: ";
784 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
" same template for num and den ";
789 edm::LogError(
"SiPixelChargeReweightingAlgorithm") <<
"Cluster Charge Reweighting did not work properly.";
793 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
" Cluster after reweighting: ";
797 for (
int row = 0; row <
TXSIZE; ++row) {
801 if ((hitPixel.first + row -
THX) >= 0 && (hitPixel.first + row -
THX) < topol->
nrows() &&
804 if constexpr (std::is_same_v<AmplitudeType, digitizerUtility::Ph2Amplitude>) {
806 <<
"Phase-2 late charge reweighting not supported (not sure we need it at all)";
816 if (chargeBefore != 0. && chargeAfter == 0.) {
821 LogDebug(
"SiPixelChargeReweightingAlgorithm")
822 <<
"Charges (before->after): " << chargeBefore <<
" -> " << chargeAfter;
823 LogDebug(
"SiPixelChargeReweightingAlgorithm") <<
"Charge loss: " << (1 - chargeAfter / chargeBefore) * 100 <<
" %";
827 for (
auto&
i : theDigiSignal) {
829 int chanDigi =
i.first;
831 int row_digi = ip.first;
832 int col_digi = ip.second;
833 int i_transformed_row = row_digi - hitPixel.first +
THX;
834 int i_transformed_col = col_digi - hitPixel.second +
THY;
835 if (i_transformed_row < 0 || i_transformed_row >
TXSIZE || i_transformed_col < 0 || i_transformed_col >
TYSIZE) {
837 if (chanDigi >= 0 &&
i.second > 0) {
838 if constexpr (std::is_same_v<AmplitudeType, digitizerUtility::Ph2Amplitude>) {
840 <<
"Phase-2 late charge reweighting not supported (not sure we need it at all)";
843 theNewDigiSignal[chanDigi] += AmplitudeType(
i.second,
i.second);
boost::multi_array< float, 2 > array_2d
const bool PrintTemplates
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, SignalType &theSignal, unsigned short int processType, const bool &boolmakeDigiSimLinks)
float xsize()
pixel x-size (microns)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
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)
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
edm::ESGetToken< SiPixel2DTemplateDBObject, SiPixel2DTemplateDBObjectRcd > SiPixel2DTemp_num_token_
virtual int ncolumns() const =0
std::vector< bool > xdouble
virtual int nrows() const =0
const bool UseReweighting
const SiPixel2DTemplateDBObject * dbobject_den
Log< level::Error, false > LogError
float s50()
1/2 of the pixel threshold signal in adc units
bool lateSignalReweight(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, PixelSimHitExtraInfo &loopTempSH, SignalType &theNewDigiSignal, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *engine)
static int pixelToChannel(int row, int col)
bool applyLateReweighting_
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
SiPixelTemplate2D templ2D
virtual bool isItBigPixelInX(int ixbin) const =0
void init(const edm::EventSetup &es)
GloballyPositioned< double > Frame
std::vector< edm::ParameterSet > Parameters
virtual bool isItBigPixelInY(int iybin) const =0
Abs< T >::type abs(const T &t)
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
DetId geographicalId() const
The label of this GeomDet.
const SiPixel2DTemplateDBObject * dbobject_num
void printCluster(array_2d &cluster)
int PixelTempRewgt2D(int id_gen, int id_rewgt, array_2d &cluster)
constexpr uint32_t rawId() const
get the raw id
SiPixelChargeReweightingAlgorithm(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
~SiPixelChargeReweightingAlgorithm()
std::vector< bool > ydouble
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_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
std::vector< SiPixelTemplateStore2D > templateStores_
std::vector< float > track
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore2D > &pixelTemp, std::string dir="CalibTracker/SiPixelESProducers/data/")
Log< level::Warning, false > LogWarning
static constexpr float cmToMicrons