CMS 3D CMS Logo

TrackstersPCA.cc
Go to the documentation of this file.
2 #include "TrackstersPCA.h"
3 #include "TPrincipal.h"
4 
5 #include <iostream>
6 
7 #include <Eigen/Core>
8 #include <Eigen/Dense>
9 
10 void ticl::assignPCAtoTracksters(std::vector<Trackster> &tracksters,
11  const std::vector<reco::CaloCluster> &layerClusters,
12  double z_limit_em,
13  bool energyWeight) {
14  LogDebug("TrackstersPCA_Eigen") << "------- Eigen -------" << std::endl;
15 
16  for (auto &trackster : tracksters) {
17  Eigen::Vector3d point;
18  point << 0., 0., 0.;
19  Eigen::Vector3d barycenter;
20  barycenter << 0., 0., 0.;
21 
22  auto fillPoint = [&](const reco::CaloCluster &c, const float weight = 1.f) {
23  point[0] = weight * c.x();
24  point[1] = weight * c.y();
25  point[2] = weight * c.z();
26  };
27 
28  // Initialize this trackster with default, dummy values
29  trackster.setRawEnergy(0.f);
30  trackster.setRawEmEnergy(0.f);
31  trackster.setRawPt(0.f);
32  trackster.setRawEmPt(0.f);
33 
34  size_t N = trackster.vertices().size();
35  float weight = 1.f / N;
36  float weights2_sum = 0.f;
37  Eigen::Vector3d sigmas;
38  sigmas << 0., 0., 0.;
39  Eigen::Vector3d sigmasEigen;
40  sigmasEigen << 0., 0., 0.;
41  Eigen::Matrix3d covM = Eigen::Matrix3d::Zero();
42 
43  for (size_t i = 0; i < N; ++i) {
44  auto fraction = 1.f / trackster.vertex_multiplicity(i);
45  trackster.addToRawEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
46  if (std::abs(layerClusters[trackster.vertices(i)].z()) <= z_limit_em)
47  trackster.addToRawEmEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
48 
49  // Compute the weighted barycenter.
50  if (energyWeight)
51  weight = layerClusters[trackster.vertices(i)].energy() * fraction;
52  fillPoint(layerClusters[trackster.vertices(i)], weight);
53  for (size_t j = 0; j < 3; ++j)
54  barycenter[j] += point[j];
55  }
56  if (energyWeight && trackster.raw_energy())
57  barycenter /= trackster.raw_energy();
58 
59  // Compute the Covariance Matrix and the sum of the squared weights, used
60  // to compute the correct normalization.
61  // The barycenter has to be known.
62  for (size_t i = 0; i < N; ++i) {
63  fillPoint(layerClusters[trackster.vertices(i)]);
64  if (energyWeight && trackster.raw_energy())
65  weight =
66  (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) / trackster.raw_energy();
67  weights2_sum += weight * weight;
68  for (size_t x = 0; x < 3; ++x)
69  for (size_t y = 0; y <= x; ++y) {
70  covM(x, y) += weight * (point[x] - barycenter[x]) * (point[y] - barycenter[y]);
71  covM(y, x) = covM(x, y);
72  }
73  }
74  covM *= 1. / (1. - weights2_sum);
75 
76  // Perform the actual decomposition
77  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d>::RealVectorType eigenvalues_fromEigen;
78  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d>::EigenvectorsType eigenvectors_fromEigen;
79  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eigensolver(covM);
80  if (eigensolver.info() != Eigen::Success) {
81  eigenvalues_fromEigen = eigenvalues_fromEigen.Zero();
82  eigenvectors_fromEigen = eigenvectors_fromEigen.Zero();
83  } else {
84  eigenvalues_fromEigen = eigensolver.eigenvalues();
85  eigenvectors_fromEigen = eigensolver.eigenvectors();
86  }
87 
88  // Compute the spread in the both spaces.
89  for (size_t i = 0; i < N; ++i) {
90  fillPoint(layerClusters[trackster.vertices(i)]);
91  sigmas += weight * (point - barycenter).cwiseAbs2();
92  Eigen::Vector3d point_transformed = eigenvectors_fromEigen * (point - barycenter);
93  if (energyWeight && trackster.raw_energy())
94  weight =
95  (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) / trackster.raw_energy();
96  sigmasEigen += weight * (point_transformed.cwiseAbs2());
97  }
98  sigmas /= (1. - weights2_sum);
99  sigmasEigen /= (1. - weights2_sum);
100 
101  // Add trackster attributes
102  trackster.setBarycenter(ticl::Trackster::Vector(barycenter));
103  trackster.fillPCAVariables(
104  eigenvalues_fromEigen, eigenvectors_fromEigen, sigmas, sigmasEigen, 3, ticl::Trackster::PCAOrdering::ascending);
105 
106  LogDebug("TrackstersPCA") << "Use energy weighting: " << energyWeight << std::endl;
107  LogDebug("TrackstersPCA") << "\nTrackster characteristics: " << std::endl;
108  LogDebug("TrackstersPCA") << "Size: " << N << std::endl;
109  LogDebug("TrackstersPCA") << "Energy: " << trackster.raw_energy() << std::endl;
110  LogDebug("TrackstersPCA") << "raw_pt: " << trackster.raw_pt() << std::endl;
111  LogDebug("TrackstersPCA") << "Means: " << barycenter[0] << ", " << barycenter[1] << ", " << barycenter[2]
112  << std::endl;
113  LogDebug("TrackstersPCA") << "EigenValues from Eigen/Tr(cov): " << eigenvalues_fromEigen[2] / covM.trace() << ", "
114  << eigenvalues_fromEigen[1] / covM.trace() << ", "
115  << eigenvalues_fromEigen[0] / covM.trace() << std::endl;
116  LogDebug("TrackstersPCA") << "EigenValues from Eigen: " << eigenvalues_fromEigen[2] << ", "
117  << eigenvalues_fromEigen[1] << ", " << eigenvalues_fromEigen[0] << std::endl;
118  LogDebug("TrackstersPCA") << "EigenVector 3 from Eigen: " << eigenvectors_fromEigen(0, 2) << ", "
119  << eigenvectors_fromEigen(1, 2) << ", " << eigenvectors_fromEigen(2, 2) << std::endl;
120  LogDebug("TrackstersPCA") << "EigenVector 2 from Eigen: " << eigenvectors_fromEigen(0, 1) << ", "
121  << eigenvectors_fromEigen(1, 1) << ", " << eigenvectors_fromEigen(2, 1) << std::endl;
122  LogDebug("TrackstersPCA") << "EigenVector 1 from Eigen: " << eigenvectors_fromEigen(0, 0) << ", "
123  << eigenvectors_fromEigen(1, 0) << ", " << eigenvectors_fromEigen(2, 0) << std::endl;
124  LogDebug("TrackstersPCA") << "Original sigmas: " << sigmas[0] << ", " << sigmas[1] << ", " << sigmas[2]
125  << std::endl;
126  LogDebug("TrackstersPCA") << "SigmasEigen in PCA space: " << sigmasEigen[2] << ", " << sigmasEigen[1] << ", "
127  << sigmasEigen[0] << std::endl;
128  LogDebug("TrackstersPCA") << "covM: \n" << covM << std::endl;
129  }
130 }
ticl::Trackster::Vector
math::XYZVector Vector
Definition: Trackster.h:21
TrackstersPCA.h
ticl::assignPCAtoTracksters
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, double, bool energyWeight=true)
Definition: TrackstersPCA.cc:10
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
mps_merge.weight
weight
Definition: mps_merge.py:88
ticl::Trackster::PCAOrdering::ascending
reco::CaloCluster
Definition: CaloCluster.h:31
vertices_cff.x
x
Definition: vertices_cff.py:29
N
#define N
Definition: blowfish.cc:9
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
point
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
HLT_2018_cff.fraction
fraction
Definition: HLT_2018_cff.py:51317
weight
Definition: weight.py:1