CMS 3D CMS Logo

TracksterLinkingbySkeletons.cc
Go to the documentation of this file.
2 #include <Math/VectorUtil.h>
11 #include "TICLGraph.h"
12 
13 namespace {
14  bool isRoundTrackster(std::array<ticl::Vector, 3> skeleton) { return (skeleton[0].Z() == skeleton[2].Z()); }
15 
16  bool isGoodTrackster(const ticl::Trackster &trackster,
17  const std::array<ticl::Vector, 3> &skeleton,
18  const unsigned int min_num_lcs,
19  const float min_trackster_energy,
20  const float pca_quality_th) {
21  bool isGood = false;
22 
23  if (isRoundTrackster(skeleton) or trackster.vertices().size() < min_num_lcs or
24  trackster.raw_energy() < min_trackster_energy) {
25  isGood = false;
26  } else {
27  auto const &eigenvalues = trackster.eigenvalues();
28  auto const sum = std::accumulate(std::begin(eigenvalues), std::end(eigenvalues), 0.f);
29  float pcaQuality = eigenvalues[0] / sum;
30  if (pcaQuality > pca_quality_th) {
31  isGood = true;
32  }
33  }
34  return isGood;
35  }
36 
37  //distance between skeletons
38  float projective_distance(const ticl::Vector &point1, const ticl::Vector &point2) {
39  // squared projective distance
40  float r1 = std::sqrt(point1.x() * point1.x() + point1.y() * point1.y());
41  float r2_at_z1 =
42  std::sqrt(point2.x() * point2.x() + point2.y() * point2.y()) * std::abs(point1.z()) / std::abs(point2.z());
43  float delta_phi = reco::deltaPhi(point1.Phi(), point2.Phi());
44  float projective_distance = (r1 - r2_at_z1) * (r1 - r2_at_z1) + r2_at_z1 * r2_at_z1 * delta_phi * delta_phi;
45  LogDebug("TracksterLinkingbySkeletons") << "Computing distance between point : " << point1 << " And point "
46  << point2 << " Distance " << projective_distance << std::endl;
47  return projective_distance;
48  }
49 } // namespace
50 
51 using namespace ticl;
52 
54  : TracksterLinkingAlgoBase(conf, iC),
55  timing_quality_threshold_(conf.getParameter<double>("track_time_quality_threshold")),
56  del_(conf.getParameter<double>("wind")),
57  min_num_lcs_(conf.getParameter<unsigned int>("min_num_lcs")),
58  min_trackster_energy_(conf.getParameter<double>("min_trackster_energy")),
59  pca_quality_th_(conf.getParameter<double>("pca_quality_th")),
60  dot_prod_th_(conf.getParameter<double>("dot_prod_th")),
61  max_distance_projective_sqr_(conf.getParameter<std::vector<double>>("max_distance_projective_sqr")),
62  min_distance_z_(conf.getParameter<std::vector<double>>("min_distance_z")),
63  max_distance_projective_sqr_closest_points_(
64  conf.getParameter<std::vector<double>>("max_distance_projective_sqr_closest_points")),
65  max_z_distance_closest_points_(conf.getParameter<std::vector<double>>("max_z_distance_closest_points")),
66  cylinder_radius_sqr_(conf.getParameter<std::vector<double>>("cylinder_radius_sqr"))
67 
68 {}
69 
71  // build disks at HGCal front & EM-Had interface for track propagation
72 
73  float zVal = hgcons_->waferZ(1, true);
74  std::pair<float, float> rMinMax = hgcons_->rangeR(zVal, true);
75 
76  float zVal_interface = rhtools_.getPositionLayer(rhtools_.lastLayerEE()).z();
77  std::pair<float, float> rMinMax_interface = hgcons_->rangeR(zVal_interface, true);
78 
79  for (int iSide = 0; iSide < 2; ++iSide) {
80  float zSide = (iSide == 0) ? (-1. * zVal) : zVal;
81  firstDisk_[iSide] =
82  std::make_unique<GeomDet>(Disk::build(Disk::PositionType(0, 0, zSide),
84  SimpleDiskBounds(rMinMax.first, rMinMax.second, zSide - 0.5, zSide + 0.5))
85  .get());
86 
87  zSide = (iSide == 0) ? (-1. * zVal_interface) : zVal_interface;
88  interfaceDisk_[iSide] = std::make_unique<GeomDet>(
89  Disk::build(Disk::PositionType(0, 0, zSide),
91  SimpleDiskBounds(rMinMax_interface.first, rMinMax_interface.second, zSide - 0.5, zSide + 0.5))
92  .get());
93  }
94 }
95 
97  const hgcal::RecHitTools rhtools,
98  const edm::ESHandle<MagneticField> bfieldH,
99  const edm::ESHandle<Propagator> propH) {
100  hgcons_ = hgcons;
101  rhtools_ = rhtools;
102  buildLayers();
103 
104  bfield_ = bfieldH;
105  propagator_ = propH;
106 }
107 
109  const ticl::Trackster &trackster,
110  float lower_percentage,
111  float upper_percentage,
112  const std::vector<reco::CaloCluster> &layerClusters,
113  const hgcal::RecHitTools &rhtools) {
114  auto const &vertices = trackster.vertices();
115  auto const trackster_raw_energy = trackster.raw_energy();
116  // sort vertices by layerId
117  std::vector<unsigned int> sortedVertices(vertices);
118  std::sort(sortedVertices.begin(), sortedVertices.end(), [&layerClusters](unsigned int i, unsigned int j) {
119  return std::abs(layerClusters[i].z()) < std::abs(layerClusters[j].z());
120  });
121 
122  // now loop over sortedVertices and find the layerId that contains the lower_percentage of the energy
123  // and the layerId that contains the upper_percentage of the energy
124  float cumulativeEnergyFraction = 0.f;
125  int innerLayerId = rhtools.getLayerWithOffset(layerClusters[sortedVertices[0]].hitsAndFractions()[0].first);
126  float innerLayerZ = layerClusters[sortedVertices[0]].z();
127  int outerLayerId = rhtools.getLayerWithOffset(layerClusters[sortedVertices.back()].hitsAndFractions()[0].first);
128  float outerLayerZ = layerClusters[sortedVertices.back()].z();
129  bool foundInnerLayer = false;
130  bool foundOuterLayer = false;
131  for (auto const &v : sortedVertices) {
132  auto const &lc = layerClusters[v];
133  auto const &n_lay = rhtools.getLayerWithOffset(lc.hitsAndFractions()[0].first);
134  cumulativeEnergyFraction += lc.energy() / trackster_raw_energy;
135  if (cumulativeEnergyFraction >= lower_percentage and not foundInnerLayer) {
136  innerLayerId = n_lay;
137  innerLayerZ = lc.z();
138  foundInnerLayer = true;
139  }
140  if (cumulativeEnergyFraction >= upper_percentage and not foundOuterLayer) {
141  outerLayerId = n_lay;
142  outerLayerZ = lc.z();
143  foundOuterLayer = true;
144  }
145  }
146  std::array<ticl::Vector, 3> skeleton;
147  int minimumDistanceInLayers = 4;
148  if (outerLayerId - innerLayerId < minimumDistanceInLayers) {
149  skeleton = {{trackster.barycenter(), trackster.barycenter(), trackster.barycenter()}};
150  } else {
151  auto intersectLineWithSurface = [](float surfaceZ, const Vector &origin, const Vector &direction) -> Vector {
152  auto const t = (surfaceZ - origin.Z()) / direction.Z();
153  auto const iX = t * direction.X() + origin.X();
154  auto const iY = t * direction.Y() + origin.Y();
155  auto const iZ = surfaceZ;
156  const Vector intersection(iX, iY, iZ);
157  return intersection;
158  };
159 
160  auto const &t0_p1 = trackster.barycenter();
161  auto const t0_p0 = intersectLineWithSurface(innerLayerZ, t0_p1, trackster.eigenvectors(0));
162  auto const t0_p2 = intersectLineWithSurface(outerLayerZ, t0_p1, trackster.eigenvectors(0));
163  skeleton = {{t0_p0, t0_p1, t0_p2}};
164  std::sort(skeleton.begin(), skeleton.end(), [](Vector &v1, Vector &v2) { return v1.Z() < v2.Z(); });
165  }
166 
167  return skeleton;
168 }
169 
170 bool isInCylinder(const std::array<ticl::Vector, 3> &mySkeleton,
171  const std::array<ticl::Vector, 3> &otherSkeleton,
172  const float radius_sqr) {
173  const auto &first = mySkeleton[0];
174  const auto &last = mySkeleton[2];
175  const auto &pointToCheck = otherSkeleton[0];
176 
177  const auto &cylAxis = last - first;
178  const auto &vecToPoint = pointToCheck - first;
179 
180  auto axisNorm = cylAxis.Dot(cylAxis);
181  auto projLength = vecToPoint.Dot(cylAxis) / axisNorm;
182  bool isWithinLength = projLength >= 0 && projLength <= 1;
183 
184  if (!isWithinLength)
185  return false;
186 
187  const auto &proj = cylAxis * projLength;
188 
189  const auto &pointOnAxis = first + proj;
190 
191  const auto &distance = pointToCheck - pointOnAxis;
192  auto distance2 = distance.Dot(distance);
193 
194  bool isWithinRadius = distance2 <= radius_sqr;
195 
196  return isWithinRadius;
197 }
198 
200  const ticl::Trackster &otherTrackster,
201  const std::array<ticl::Vector, 3> &mySkeleton,
202  const std::array<ticl::Vector, 3> &otherSkeleton) {
203  //do not start links from small/bad tracksters
204  float zVal_interface = rhtools_.getPositionLayer(rhtools_.lastLayerEE()).z();
205  if (!isGoodTrackster(myTrackster, mySkeleton, min_num_lcs_, min_trackster_energy_, pca_quality_th_)) {
206  LogDebug("TracksterLinkingbySkeletons") << "Inner Trackster with energy " << myTrackster.raw_energy() << " Num LCs "
207  << myTrackster.vertices().size() << " NOT GOOD " << std::endl;
208  return false;
209  } else {
210  LogDebug("TracksterLinkingbySkeletons") << "Inner Trackster wi energy " << myTrackster.raw_energy() << " Num LCs "
211  << myTrackster.vertices().size() << " IS GOOD " << std::endl;
212  float proj_distance = projective_distance(mySkeleton[1], otherSkeleton[1]);
213  auto isEE = mySkeleton[1].z() <= zVal_interface ? 0 : 1;
214  bool areAlignedInProjectiveSpace = proj_distance < max_distance_projective_sqr_[isEE];
215  LogDebug("TracksterLinkingbySkeletons")
216  << "\t Trying to compare with outer Trackster with energy " << otherTrackster.raw_energy() << " Num LCS "
217  << otherTrackster.vertices().size() << " Projective distance " << proj_distance << " areAlignedProjective "
218  << areAlignedInProjectiveSpace << " TH " << max_distance_projective_sqr_[isEE] << std::endl;
219  //check if otherTrackster is good
220  if (isGoodTrackster(otherTrackster, otherSkeleton, min_num_lcs_, min_trackster_energy_, pca_quality_th_)) {
221  // if both tracksters are good, then we can check the projective distance between the barycenters.
222  // if the barycenters are aligned, then we check that the two skeletons are aligned
223  if (areAlignedInProjectiveSpace) {
224  auto dotProdSkeletons =
225  ((mySkeleton[2] - mySkeleton[0]).Unit()).Dot((otherSkeleton[2] - otherSkeleton[0]).Unit());
226  bool alignedSkeletons = dotProdSkeletons > dot_prod_th_;
227  LogDebug("TracksterLinkingbySkeletons")
228  << "\t Outer Trackster is Good, checking for skeleton alignment " << alignedSkeletons << " dotProd "
229  << dotProdSkeletons << " Threshold " << dot_prod_th_ << std::endl;
230  if (alignedSkeletons)
231  LogDebug("TracksterLinkingbySkeletons") << "\t\t Linked! " << std::endl;
232  return alignedSkeletons;
233  } else {
234  // we measure the distance between the two closest nodes in the two skeletons
235  LogDebug("TracksterLinkingbySkeletons")
236  << "\t Outer Trackster is not aligned, check skeletons distances " << std::endl;
237  int myClosestPoint = -1;
238  int otherClosestPoint = -1;
239  float minDistance_z = std::numeric_limits<float>::max();
240  for (int i = 0; i < 3; i++) {
241  for (int j = 0; j < 3; j++) {
242  float dist_z = std::abs(mySkeleton[i].Z() - otherSkeleton[j].Z());
243  if (dist_z < minDistance_z) {
244  myClosestPoint = i;
245  otherClosestPoint = j;
246  minDistance_z = dist_z;
247  }
248  }
249  }
250  if (minDistance_z < min_distance_z_[isEE]) {
251  LogDebug("TracksterLinkingbySkeletons")
252  << "\t Trackster have distance in Z " << minDistance_z
253  << "Checking if they are aligned in projective space "
254  << projective_distance(mySkeleton[myClosestPoint], otherSkeleton[otherClosestPoint]) << " TH "
255  << max_distance_projective_sqr_[isEE] << std::endl;
256  if (projective_distance(mySkeleton[myClosestPoint], otherSkeleton[otherClosestPoint]) <
258  LogDebug("TracksterLinkingbySkeletons") << "\t\t Linked! " << std::endl;
259  }
260  return projective_distance(mySkeleton[myClosestPoint], otherSkeleton[otherClosestPoint]) <
262  } else {
263  LogDebug("TracksterLinkingbySkeletons") << "\t\t Not Linked Distance Z " << minDistance_z << std::endl;
264  return false;
265  }
266  }
267  } else {
268  LogDebug("TracksterLinkingbySkeletons")
269  << "\t Outer Trackster is NOT GOOD, check projective space alignment " << areAlignedInProjectiveSpace
270  << " proj_distance " << max_distance_projective_sqr_[isEE] << std::endl;
271  if (areAlignedInProjectiveSpace) {
272  LogDebug("TracksterLinkingbySkeletons") << "\t\t Linked! " << std::endl;
273  return true;
274  } else {
275  LogDebug("TracksterLinkingbySkeletons")
276  << "\t Not aligned in projective space, check distance between closest points in the two skeletons "
277  << std::endl;
278  // we measure the distance between the two closest nodes in the two skeletons
279  int myClosestPoint = -1;
280  int otherClosestPoint = -1;
281  float minDistance_z = std::numeric_limits<float>::max();
282  // we skip the innermost node of mySkeleton
283  for (int i = 1; i < 3; i++) {
284  for (int j = 0; j < 3; j++) {
285  float dist_z = std::abs(mySkeleton[i].Z() - otherSkeleton[j].Z());
286  if (dist_z < minDistance_z) {
287  myClosestPoint = i;
288  otherClosestPoint = j;
289  minDistance_z = dist_z;
290  }
291  }
292  }
293  float d = projective_distance(mySkeleton[myClosestPoint], otherSkeleton[otherClosestPoint]);
294  LogDebug("TracksterLinkingbySkeletons")
295  << "\t\t Distance between closest points " << d << " TH " << 10.f << " Z Distance " << minDistance_z
296  << " TH " << max_distance_projective_sqr_closest_points_[isEE] << std::endl;
298  minDistance_z < max_z_distance_closest_points_[isEE]) {
299  LogDebug("TracksterLinkingbySkeletons") << "\t\t\t Linked! " << d << std::endl;
300  return true;
301  } else {
302  LogDebug("TracksterLinkingbySkeletons") << "Distance between closest point " << d << " Distance in z "
303  << max_z_distance_closest_points_[isEE] << std::endl;
304  bool isInCyl = isInCylinder(mySkeleton, otherSkeleton, cylinder_radius_sqr_[isEE]);
305  LogDebug("TracksterLinkingbySkeletons") << "Two Points are in Cylinder " << isInCylinder << std::endl;
306  if (isInCyl) {
307  LogDebug("TracksterLinkingbySkeletons") << "\t\t\t Linked! " << d << std::endl;
308  }
309  return isInCyl;
310  }
311  }
312  }
313  }
314 }
315 
317  const Inputs &input,
318  std::vector<Trackster> &resultTracksters,
319  std::vector<std::vector<unsigned int>> &linkedResultTracksters,
320  std::vector<std::vector<unsigned int>> &linkedTracksterIdToInputTracksterId) {
321  const auto &tracksters = input.tracksters;
322  const auto &layerClusters = input.layerClusters;
323 
324  // sort tracksters by energy
325  std::vector<unsigned int> sortedTracksters(tracksters.size());
326  std::iota(sortedTracksters.begin(), sortedTracksters.end(), 0);
327  std::sort(sortedTracksters.begin(), sortedTracksters.end(), [&tracksters](unsigned int i, unsigned int j) {
328  return tracksters[i].raw_energy() > tracksters[j].raw_energy();
329  });
330  // fill tiles for trackster linking
331  // tile 0 for negative eta
332  // tile 1 for positive eta
333  std::array<TICLLayerTile, 2> tracksterTile;
334  // loop over tracksters sorted by energy and calculate skeletons
335  // fill tiles for trackster linking
336  std::vector<std::array<ticl::Vector, 3>> skeletons(tracksters.size());
337  for (auto const t_idx : sortedTracksters) {
338  const auto &trackster = tracksters[t_idx];
339  skeletons[t_idx] = findSkeletonNodes(tracksters[t_idx], 0.1, 0.9, layerClusters, rhtools_);
340  tracksterTile[trackster.barycenter().eta() > 0.f].fill(
341  trackster.barycenter().eta(), trackster.barycenter().phi(), t_idx);
342  }
343  std::vector<int> maskReceivedLink(tracksters.size(), 1);
344  std::vector<int> isRootTracksters(tracksters.size(), 1);
345 
346  std::vector<Node> allNodes;
347  for (size_t it = 0; it < tracksters.size(); ++it) {
348  allNodes.emplace_back(it);
349  }
350 
351  // loop over tracksters sorted by energy and link them
352  for (auto const &t_idx : sortedTracksters) {
353  auto const &trackster = tracksters[t_idx];
354  auto const &skeleton = skeletons[t_idx];
355 
356  auto const bary = trackster.barycenter();
357  float eta_min = std::max(abs(bary.eta()) - del_, TileConstants::minEta);
358  float eta_max = std::min(abs(bary.eta()) + del_, TileConstants::maxEta);
359  int tileIndex = bary.eta() > 0.f;
360  const auto &tiles = tracksterTile[tileIndex];
361  std::array<int, 4> search_box = tiles.searchBoxEtaPhi(eta_min, eta_max, bary.phi() - del_, bary.phi() + del_);
362  if (search_box[2] > search_box[3]) {
363  search_box[3] += TileConstants::nPhiBins;
364  }
365 
366  for (int eta_i = search_box[0]; eta_i <= search_box[1]; ++eta_i) {
367  for (int phi_i = search_box[2]; phi_i <= search_box[3]; ++phi_i) {
368  auto &neighbours = tiles[tiles.globalBin(eta_i, (phi_i % TileConstants::nPhiBins))];
369  for (auto n : neighbours) {
370  if (t_idx == n)
371  continue;
372  if (maskReceivedLink[n] == 0 or allNodes[t_idx].isInnerNeighbour(n))
373  continue;
374  if (isGoodTrackster(trackster, skeleton, min_num_lcs_, min_trackster_energy_, pca_quality_th_)) {
375  LogDebug("TracksterLinkingbySkeletons")
376  << "Trying to Link Trackster " << t_idx << " With Trackster " << n << std::endl;
377  if (areCompatible(trackster, tracksters[n], skeleton, skeletons[n])) {
378  LogDebug("TracksterLinkingbySkeletons")
379  << "\t==== LINK: Trackster " << t_idx << " Linked with Trackster " << n << std::endl;
380  maskReceivedLink[n] = 0;
381  allNodes[t_idx].addOuterNeighbour(n);
382  allNodes[n].addInnerNeighbour(t_idx);
383  isRootTracksters[n] = 0;
384  }
385  }
386  }
387  }
388  }
389  }
390 
391  LogDebug("TracksterLinkingbySkeletons") << "**************** FINAL GRAPH **********************" << std::endl;
392  for (auto const &node : allNodes) {
393  if (isRootTracksters[node.getId()]) {
394  LogDebug("TracksterLinkingbySkeletons")
395  << "ISROOT "
396  << " Node " << node.getId() << " position " << tracksters[node.getId()].barycenter() << " energy "
397  << tracksters[node.getId()].raw_energy() << std::endl;
398  } else {
399  LogDebug("TracksterLinkingbySkeletons")
400  << "Node " << node.getId() << " position " << tracksters[node.getId()].barycenter() << " energy "
401  << tracksters[node.getId()].raw_energy() << std::endl;
402  }
403  }
404  LogDebug("TracksterLinkingbySkeletons") << "********************************************************" << std::endl;
405 
406  TICLGraph graph(allNodes, isRootTracksters);
407 
408  int ic = 0;
409  auto const &components = graph.findSubComponents();
410  linkedTracksterIdToInputTracksterId.resize(components.size());
411  for (auto const &comp : components) {
412  LogDebug("TracksterLinkingbySkeletons") << "Component " << ic << " Node: ";
413  std::vector<unsigned int> linkedTracksters;
414  Trackster outTrackster;
415  if (comp.size() == 1) {
416  if (input.tracksters[comp[0]].vertices().size() <= 3) {
417  continue;
418  }
419  }
420  for (auto const &node : comp) {
421  LogDebug("TracksterLinkingbySkeletons") << node << " ";
422  linkedTracksterIdToInputTracksterId[ic].push_back(node);
423  outTrackster.mergeTracksters(input.tracksters[node]);
424  }
425  linkedTracksters.push_back(resultTracksters.size());
426  resultTracksters.push_back(outTrackster);
427  linkedResultTracksters.push_back(linkedTracksters);
428  LogDebug("TracksterLinkingbySkeletons") << "\n";
429  ++ic;
430  }
431 } // linkTracksters
double waferZ(int layer, bool reco) const
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
static DiskPointer build(Args &&... args)
Definition: BoundDisk.h:38
edm::ESHandle< MagneticField > bfield_
static constexpr int nPhiBins
Definition: Common.h:15
const Vector & barycenter() const
Definition: Trackster.h:159
static constexpr float maxEta
Definition: Common.h:13
void linkTracksters(const Inputs &input, std::vector< Trackster > &resultTracksters, std::vector< std::vector< unsigned int >> &linkedResultTracksters, std::vector< std::vector< unsigned int >> &linkedTracksterIdToInputTracksterId) override
bool areCompatible(const ticl::Trackster &myTrackster, const ticl::Trackster &otherTrackster, const std::array< ticl::Vector, 3 > &mySkeleton, const std::array< ticl::Vector, 3 > &otherSkeleton)
static std::string const input
Definition: EdmProvDump.cc:50
int32_t tileIndex(int32_t layer, int32_t ring, int32_t phi)
const std::array< Vector, 3 > & eigenvectors() const
Definition: Trackster.h:161
TracksterLinkingbySkeletons(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
void mergeTracksters(const Trackster &other)
Definition: Trackster.h:81
std::array< ticl::Vector, 3 > findSkeletonNodes(const ticl::Trackster &trackster, float lower_percentage, float upper_percentage, const std::vector< reco::CaloCluster > &layerClusters, const hgcal::RecHitTools &rhtools)
const float raw_energy() const
Definition: Trackster.h:154
T sqrt(T t)
Definition: SSEVec.h:23
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
std::unique_ptr< GeomDet > interfaceDisk_[2]
static constexpr float minEta
Definition: Common.h:12
d
Definition: ztail.py:151
std::vector< double > max_distance_projective_sqr_closest_points_
std::pair< double, double > rangeR(double z, bool reco) const
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
std::vector< std::vector< unsigned int > > findSubComponents()
Definition: TICLGraph.cc:20
std::unique_ptr< GeomDet > firstDisk_[2]
GlobalPoint getPositionLayer(int layer, bool nose=false) const
Definition: RecHitTools.cc:152
math::XYZVectorF Vector
Definition: Common.h:42
void initialize(const HGCalDDDConstants *hgcons, const hgcal::RecHitTools rhtools, const edm::ESHandle< MagneticField > bfieldH, const edm::ESHandle< Propagator > propH) override
Definition: Common.h:10
const std::array< float, 3 > & eigenvalues() const
Definition: Trackster.h:160
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:76
bool isInCylinder(const std::array< ticl::Vector, 3 > &mySkeleton, const std::array< ticl::Vector, 3 > &otherSkeleton, const float radius_sqr)
#define LogDebug(id)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381