CMS 3D CMS Logo

TrackstersPCA.cc
Go to the documentation of this file.
4 #include "TrackstersPCA.h"
5 
6 #include <iostream>
7 #include <set>
8 
9 #include <Eigen/Core>
10 #include <Eigen/Dense>
11 #include <vector>
12 #include <functional>
13 
14 void ticl::assignPCAtoTracksters(std::vector<Trackster> &tracksters,
15  const std::vector<reco::CaloCluster> &layerClusters,
16  const edm::ValueMap<std::pair<float, float>> &layerClustersTime,
17  double z_limit_em,
18  const hgcal::RecHitTools &rhtools,
19  bool computeLocalTime,
20  bool energyWeight,
21  bool clean,
22  int minLayer,
23  int maxLayer) {
24  LogDebug("TrackstersPCA_Eigen") << "------- Eigen -------" << std::endl;
25 
26  for (auto &trackster : tracksters) {
27  LogDebug("TrackstersPCA_Eigen") << "start testing teackster with size:" << trackster.vertices().size() << std::endl;
28 
30  point << 0., 0., 0.;
31  Eigen::Vector3f barycenter;
32  barycenter << 0., 0., 0.;
33  Eigen::Vector3f filtered_barycenter;
34  filtered_barycenter << 0., 0., 0.;
35 
36  auto fillPoint = [&](const reco::CaloCluster &c, const float weight = 1.f) {
37  point[0] = weight * c.x();
38  point[1] = weight * c.y();
39  point[2] = weight * c.z();
40  };
41 
42  // Initialize this trackster with default, dummy values
43  trackster.setRawEnergy(0.f);
44  trackster.setRawEmEnergy(0.f);
45  trackster.setRawPt(0.f);
46  trackster.setRawEmPt(0.f);
47 
48  size_t N = trackster.vertices().size();
49  if (N == 0)
50  continue;
51  float weight = 1.f / N;
52  float weights2_sum = 0.f;
53 
54  std::vector<float> layerClusterEnergies;
55 
56  for (size_t i = 0; i < N; ++i) {
57  auto fraction = 1.f / trackster.vertex_multiplicity(i);
58  trackster.addToRawEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
59  if (std::abs(layerClusters[trackster.vertices(i)].z()) <= z_limit_em)
60  trackster.addToRawEmEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
61 
62  // Compute the weighted barycenter.
63  if (energyWeight)
64  weight = layerClusters[trackster.vertices(i)].energy() * fraction;
65  fillPoint(layerClusters[trackster.vertices(i)], weight);
66  for (size_t j = 0; j < 3; ++j)
67  barycenter[j] += point[j];
68 
69  layerClusterEnergies.push_back(layerClusters[trackster.vertices(i)].energy());
70  }
71  float raw_energy = trackster.raw_energy();
72  float inv_raw_energy = 1.f / raw_energy;
73  if (energyWeight)
74  barycenter *= inv_raw_energy;
75  trackster.setBarycenter(ticl::Trackster::Vector(barycenter));
76 
77  LogDebug("TrackstersPCA_Eigen") << "cleaning is :" << clean << std::endl;
78 
79  std::vector<unsigned> filtered_idx; // indices of layer clusters to consider for cleaned PCA
80  float filtered_energy = 0;
81  float inv_filtered_energy = 0;
82  if (clean) {
83  // Filter layerclusters for the cleaned PCA
84  auto maxE_vertex = std::distance(layerClusterEnergies.begin(),
85  std::max_element(layerClusterEnergies.begin(), layerClusterEnergies.end()));
86  auto maxE_layer = getLayerFromLC(layerClusters[trackster.vertices(maxE_vertex)], rhtools);
87 
88  auto vertices_by_layer = sortByLayer(trackster, layerClusters, rhtools);
89 
90  for (unsigned i = 1; i <= rhtools.lastLayer(); ++i) {
91  auto vertices_in_layer = vertices_by_layer[i];
92  if (vertices_in_layer.empty())
93  continue;
94 
95  std::vector<float> energies_in_layer;
96  for (auto vrt : vertices_in_layer)
97  energies_in_layer.push_back(layerClusters[trackster.vertices(vrt)].energy());
98 
99  unsigned maxEid_inLayer = std::distance(energies_in_layer.begin(),
100  std::max_element(energies_in_layer.begin(), energies_in_layer.end()));
101 
102  // layer based filtering of what goes into the PCA
103  if ((int)i >= (int)maxE_layer - minLayer && (int)i <= (int)maxE_layer + maxLayer) {
104  auto filtered_vert = vertices_in_layer[maxEid_inLayer];
105  filtered_idx.push_back(filtered_vert);
106 
107  const auto &maxE_LC = layerClusters[trackster.vertices(filtered_vert)];
108  fillPoint(maxE_LC, maxE_LC.energy() * (1.f / trackster.vertex_multiplicity(filtered_vert)));
109  for (size_t j = 0; j < 3; ++j)
110  filtered_barycenter[j] += point[j];
111  filtered_energy += maxE_LC.energy();
112  }
113  }
114  inv_filtered_energy = 1. / filtered_energy;
115  filtered_barycenter *= inv_filtered_energy;
116  }
117  LogDebug("TrackstersPCA_Eigen") << "min, max " << minLayer << " " << maxLayer << std::endl;
118 
119  std::pair<float, float> timeTrackster;
120  if (computeLocalTime)
121  timeTrackster = ticl::computeLocalTracksterTime(trackster, layerClusters, layerClustersTime, barycenter, N);
122  else
123  timeTrackster = ticl::computeTracksterTime(trackster, layerClustersTime, N);
124 
125  trackster.setTimeAndError(timeTrackster.first, timeTrackster.second);
126  LogDebug("TrackstersPCA") << "Use energy weighting: " << energyWeight << std::endl;
127  LogDebug("TrackstersPCA") << "\nTrackster characteristics: " << std::endl;
128  LogDebug("TrackstersPCA") << "Size: " << N << std::endl;
129  LogDebug("TrackstersPCA") << "Energy: " << trackster.raw_energy() << std::endl;
130  LogDebug("TrackstersPCA") << "raw_pt: " << trackster.raw_pt() << std::endl;
131  LogDebug("TrackstersPCA") << "Means: " << barycenter[0] << ", " << barycenter[1] << ", " << barycenter[2]
132  << std::endl;
133  LogDebug("TrackstersPCA") << "Time: " << trackster.time() << " +/- " << trackster.timeError() << std::endl;
134 
135  if (N > 2) {
136  Eigen::Vector3f sigmas;
137  sigmas << 0., 0., 0.;
138  Eigen::Vector3f sigmasEigen;
139  sigmasEigen << 0., 0., 0.;
140  Eigen::Matrix3f covM = Eigen::Matrix3f::Zero();
141  // Compute the Covariance Matrix and the sum of the squared weights, used
142  // to compute the correct normalization.
143  // The barycenter has to be known.
144 
145  auto calc_covM = [&](size_t i) {
146  fillPoint(layerClusters[trackster.vertices(i)]);
147  if (energyWeight && trackster.raw_energy()) {
148  weight = (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) *
149  (clean ? inv_filtered_energy : inv_raw_energy);
150  if (trackster.vertex_multiplicity(i) > 1)
151  LogDebug("TrackstersPCA_Eigen")
152  << "trackster.vertex_multiplicity(i) :" << trackster.vertex_multiplicity(i);
153  }
154  weights2_sum += weight * weight;
155  for (size_t x = 0; x < 3; ++x) {
156  for (size_t y = 0; y <= x; ++y) {
157  covM(x, y) += weight * (point[x] - (clean ? filtered_barycenter[x] : barycenter[x])) *
158  (point[y] - (clean ? filtered_barycenter[y] : barycenter[y]));
159  covM(y, x) = covM(x, y);
160  }
161  }
162  };
163  if (clean) {
164  for (size_t i : filtered_idx) {
165  calc_covM(i);
166  }
167  } else {
168  for (size_t i = 0; i < N; ++i) {
169  calc_covM(i);
170  }
171  }
172 
173  covM *= 1.f / (1.f - weights2_sum);
174 
175  // Perform the actual decomposition
176  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f>::RealVectorType eigenvalues_fromEigen;
177  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f>::EigenvectorsType eigenvectors_fromEigen;
178  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigensolver(covM);
179  if (eigensolver.info() != Eigen::Success) {
180  eigenvalues_fromEigen = eigenvalues_fromEigen.Zero();
181  eigenvectors_fromEigen = eigenvectors_fromEigen.Zero();
182  } else {
183  eigenvalues_fromEigen = eigensolver.eigenvalues();
184  eigenvectors_fromEigen = eigensolver.eigenvectors();
185  }
186 
187  // Compute the spread in the both spaces.
188  auto calc_spread = [&](size_t i) {
189  fillPoint(layerClusters[trackster.vertices(i)]);
190  sigmas += weight * (point - (clean ? filtered_barycenter : barycenter)).cwiseAbs2();
191  Eigen::Vector3f point_transformed =
192  eigenvectors_fromEigen * (point - (clean ? filtered_barycenter : barycenter));
193  if (energyWeight && raw_energy)
194  weight = (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) *
195  (clean ? inv_filtered_energy : inv_raw_energy);
196  sigmasEigen += weight * (point_transformed.cwiseAbs2());
197  };
198 
199  if (clean) {
200  for (size_t i : filtered_idx) {
201  calc_spread(i);
202  }
203  } else {
204  for (size_t i = 0; i < N; ++i) {
205  calc_spread(i);
206  }
207  }
208 
209  sigmas /= (1.f - weights2_sum);
210  sigmasEigen /= (1.f - weights2_sum);
211 
212  trackster.fillPCAVariables(eigenvalues_fromEigen,
213  eigenvectors_fromEigen,
214  sigmas,
215  sigmasEigen,
216  3,
218 
219  LogDebug("TrackstersPCA") << "EigenValues from Eigen/Tr(cov): " << eigenvalues_fromEigen[2] / covM.trace() << ", "
220  << eigenvalues_fromEigen[1] / covM.trace() << ", "
221  << eigenvalues_fromEigen[0] / covM.trace() << std::endl;
222  LogDebug("TrackstersPCA") << "EigenValues from Eigen: " << eigenvalues_fromEigen[2] << ", "
223  << eigenvalues_fromEigen[1] << ", " << eigenvalues_fromEigen[0] << std::endl;
224  LogDebug("TrackstersPCA") << "EigenVector 3 from Eigen: " << eigenvectors_fromEigen(0, 2) << ", "
225  << eigenvectors_fromEigen(1, 2) << ", " << eigenvectors_fromEigen(2, 2) << std::endl;
226  LogDebug("TrackstersPCA") << "EigenVector 2 from Eigen: " << eigenvectors_fromEigen(0, 1) << ", "
227  << eigenvectors_fromEigen(1, 1) << ", " << eigenvectors_fromEigen(2, 1) << std::endl;
228  LogDebug("TrackstersPCA") << "EigenVector 1 from Eigen: " << eigenvectors_fromEigen(0, 0) << ", "
229  << eigenvectors_fromEigen(1, 0) << ", " << eigenvectors_fromEigen(2, 0) << std::endl;
230  LogDebug("TrackstersPCA") << "Original sigmas: " << sigmas[0] << ", " << sigmas[1] << ", " << sigmas[2]
231  << std::endl;
232  LogDebug("TrackstersPCA") << "SigmasEigen in PCA space: " << sigmasEigen[2] << ", " << sigmasEigen[1] << ", "
233  << sigmasEigen[0] << std::endl;
234  LogDebug("TrackstersPCA") << "covM: \n" << covM << std::endl;
235  }
236  }
237 }
238 
239 std::pair<float, float> ticl::computeLocalTracksterTime(const Trackster &trackster,
240  const std::vector<reco::CaloCluster> &layerClusters,
241  const edm::ValueMap<std::pair<float, float>> &layerClustersTime,
242  const Eigen::Vector3f &barycenter,
243  size_t N) {
244  float tracksterTime = 0.;
245  float tracksterTimeErr = 0.;
246 
247  auto project_lc_to_pca = [](const std::array<float, 3> &point, const std::array<float, 3> &segment_end) {
248  float dot_product = 0.0;
249  float segment_dot = 0.0;
250 
251  for (int i = 0; i < 3; ++i) {
252  dot_product += point[i] * segment_end[i];
253  segment_dot += segment_end[i] * segment_end[i];
254  }
255 
256  float projection = 0.0;
257  if (segment_dot != 0.0) {
258  projection = dot_product / segment_dot;
259  }
260 
261  std::array<float, 3> closest_point;
262  for (int i = 0; i < 3; ++i) {
263  closest_point[i] = projection * segment_end[i];
264  }
265 
266  float distanceSquared = 0.f;
267  for (int i = 0; i < 3; ++i) {
268  distanceSquared += std::pow(point[i] - closest_point[i], 2);
269  }
270  return distanceSquared;
271  };
272 
273  constexpr float c = 29.9792458; // cm/ns
274  for (size_t i = 0; i < N; ++i) {
275  // Add timing from layerClusters not already used
276  float timeE = layerClustersTime.get(trackster.vertices(i)).second;
277  if (timeE > 0.f) {
278  float time = layerClustersTime.get(trackster.vertices(i)).first;
279  timeE = 1.f / pow(timeE, 2);
280  float x = layerClusters[trackster.vertices(i)].x();
281  float y = layerClusters[trackster.vertices(i)].y();
282  float z = layerClusters[trackster.vertices(i)].z();
283 
284  if (project_lc_to_pca({{x, y, z}}, {{barycenter[0], barycenter[1], barycenter[2]}}) < 9.f) { // set MR to 3
285  float invz = 1.f / z;
286  float deltaT = 1.f / c *
287  std::sqrt(((barycenter[2] * invz - 1.f) * x) * ((barycenter[2] * invz - 1.f) * x) +
288  ((barycenter[2] * invz - 1.f) * y) * ((barycenter[2] * invz - 1.f) * y) +
289  (barycenter[2] - z) * (barycenter[2] - z));
290  time = std::abs(barycenter[2]) < std::abs(z) ? time - deltaT : time + deltaT;
291 
292  tracksterTime += time * timeE;
293  tracksterTimeErr += timeE;
294  }
295  }
296  }
297  if (tracksterTimeErr > 0.f)
298  return {tracksterTime / tracksterTimeErr, 1.f / std::sqrt(tracksterTimeErr)};
299  else
300  return {-99.f, -1.f};
301 }
302 
303 std::pair<float, float> ticl::computeTracksterTime(const Trackster &trackster,
304  const edm::ValueMap<std::pair<float, float>> &layerClustersTime,
305  size_t N) {
306  std::vector<float> times;
307  std::vector<float> timeErrors;
308 
309  for (size_t i = 0; i < N; ++i) {
310  float timeE = layerClustersTime.get(trackster.vertices(i)).second;
311  if (timeE > 0.f) {
312  times.push_back(layerClustersTime.get(trackster.vertices(i)).first);
313  timeErrors.push_back(1.f / pow(timeE, 2));
314  }
315  }
316 
318  return timeEstimator.fixSizeHighestDensity(times, timeErrors);
319 }
std::vector< std::vector< unsigned > > sortByLayer(const Trackster &ts, const std::vector< reco::CaloCluster > &layerClusters, const hgcal::RecHitTools &rhtools)
Definition: TrackstersPCA.h:47
unsigned getLayerFromLC(const reco::CaloCluster &LC, const hgcal::RecHitTools &rhtools)
Definition: TrackstersPCA.h:40
math::XYZVectorF Vector
Definition: Trackster.h:21
void assignPCAtoTracksters(std::vector< Trackster > &tracksters, const std::vector< reco::CaloCluster > &layerClusters, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, double z_limit_em, hgcal::RecHitTools const &rhTools, bool computeLocalTime=false, bool energyWeight=true, bool clean=false, int minLayer=10, int maxLayer=10)
Definition: weight.py:1
static void clean(char *s)
U second(std::pair< T, U > const &p)
std::pair< float, float > fixSizeHighestDensity(std::vector< float > &time, std::vector< float > weight=std::vector< float >(), unsigned int minNhits=3, float deltaT=0.210, float timeWidthBy=0.5)
T sqrt(T t)
Definition: SSEVec.h:23
std::pair< float, float > computeLocalTracksterTime(const Trackster &trackster, const std::vector< reco::CaloCluster > &layerClusters, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, const Eigen::Vector3f &barycenter, size_t N)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
Eigen::Matrix3f Matrix3f
Definition: FitUtils.h:51
std::pair< float, float > computeTracksterTime(const Trackster &trackster, const edm::ValueMap< std::pair< float, float >> &layerClustersTime, size_t N)
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
#define N
Definition: blowfish.cc:9
float x
Eigen::Vector3f Vector3f
Definition: FitUtils.h:52
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
*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
#define LogDebug(id)
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80