CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Classes | Typedefs | Functions
ticl Namespace Reference

Classes

class  CacheBase
 
class  ClusterFilterBase
 
class  ClusterFilterByAlgo
 
class  ClusterFilterByAlgoAndSize
 
class  ClusterFilterByAlgoAndSizeAndLayerRange
 
class  ClusterFilterBySize
 
class  PatternRecognitionAlgoBaseT
 
class  PatternRecognitionbyCA
 
class  PatternRecognitionbyCLUE3D
 
class  PatternRecognitionbyMultiClusters
 
class  SeedingRegionAlgoBase
 
class  SeedingRegionByHF
 
class  SeedingRegionByL1
 
class  SeedingRegionByTracks
 
class  SeedingRegionGlobal
 
struct  TileConstants
 
struct  TileConstantsHFNose
 
class  Trackster
 
class  TracksterMomentumPluginBase
 
class  TracksterP4FromEnergySum
 
class  TracksterP4FromTrackAndPCA
 
class  TracksterRecoTrackPlugin
 
class  TrackstersCache
 
class  TracksterTrackPluginBase
 

Typedefs

typedef std::vector< std::pair
< unsigned int, float > > 
HgcalClusterFilterMask
 
using TICLLayerTile = TICLLayerTileT< TileConstants >
 
using TICLLayerTileHFNose = TICLLayerTileT< TileConstantsHFNose >
 
using Tiles = std::array< TICLLayerTile, TileConstants::nLayers >
 
using TilesHFNose = std::array< TICLLayerTileHFNose, TileConstantsHFNose::nLayers >
 
typedef std::vector< TracksterTracksterCollection
 
using TracksterTiles = std::array< TICLLayerTile, TileConstants::iterations >
 
using TracksterTilesHFNose = std::array< TICLLayerTileHFNose, TileConstantsHFNose::iterations >
 

Functions

static void addTrackster (const int &index, const std::vector< std::pair< edm::Ref< reco::CaloClusterCollection >, std::pair< float, float >>> &lcVec, const std::vector< float > &inputClusterMask, const float &fractionCut_, const float &energy, const int &pdgId, const int &charge, const edm::ProductID &seed, std::vector< float > &output_mask, std::vector< Trackster > &result)
 
void assignPCAtoTracksters (std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool energyWeight=true)
 
Trackster::ParticleType tracksterParticleTypeFromPdgId (int pdgId, int charge)
 

Typedef Documentation

typedef std::vector<std::pair<unsigned int, float> > ticl::HgcalClusterFilterMask

Definition at line 32 of file Common.h.

Definition at line 61 of file TICLLayerTile.h.

Definition at line 65 of file TICLLayerTile.h.

using ticl::Tiles = typedef std::array<TICLLayerTile, TileConstants::nLayers>

Definition at line 62 of file TICLLayerTile.h.

Definition at line 66 of file TICLLayerTile.h.

typedef std::vector<Trackster> ticl::TracksterCollection

Definition at line 203 of file Trackster.h.

Definition at line 63 of file TICLLayerTile.h.

Definition at line 67 of file TICLLayerTile.h.

Function Documentation

static void ticl::addTrackster ( const int &  index,
const std::vector< std::pair< edm::Ref< reco::CaloClusterCollection >, std::pair< float, float >>> &  lcVec,
const std::vector< float > &  inputClusterMask,
const float &  fractionCut_,
const float &  energy,
const int &  pdgId,
const int &  charge,
const edm::ProductID seed,
std::vector< float > &  output_mask,
std::vector< Trackster > &  result 
)
static

Definition at line 36 of file commons.h.

References validate-o2o-wbm::f, HLT_FULL_cff::fraction, and tracksterParticleTypeFromPdgId().

Referenced by SimTrackstersProducer::produce().

46  {
47  if (lcVec.empty())
48  return;
49 
50  Trackster tmpTrackster;
51  tmpTrackster.zeroProbabilities();
52  tmpTrackster.vertices().reserve(lcVec.size());
53  tmpTrackster.vertex_multiplicity().reserve(lcVec.size());
54  for (auto const& [lc, energyScorePair] : lcVec) {
55  if (inputClusterMask[lc.index()] > 0) {
56  double fraction = energyScorePair.first / lc->energy();
57  if (fraction < fractionCut_)
58  continue;
59  tmpTrackster.vertices().push_back(lc.index());
60  output_mask[lc.index()] -= fraction;
61  tmpTrackster.vertex_multiplicity().push_back(1. / fraction);
62  }
63  }
64 
65  tmpTrackster.setIdProbability(tracksterParticleTypeFromPdgId(pdgId, charge), 1.f);
66  tmpTrackster.setRegressedEnergy(energy);
67  tmpTrackster.setSeed(seed, index);
68  result.emplace_back(tmpTrackster);
69  }
tuple result
Definition: mps_fire.py:311
Trackster::ParticleType tracksterParticleTypeFromPdgId(int pdgId, int charge)
Definition: commons.h:10
void ticl::assignPCAtoTracksters ( std::vector< Trackster > &  tracksters,
const std::vector< reco::CaloCluster > &  layerClusters,
const edm::ValueMap< std::pair< float, float >> &  layerClustersTime,
double  z_limit_em,
bool  energyWeight = true 
)

Definition at line 12 of file TrackstersPCA.cc.

References funct::abs(), ticl::Trackster::ascending, c, validate-o2o-wbm::f, first, hgcalsimclustertime::ComputeClusterTime::fixSizeHighestDensity(), HLT_FULL_cff::fraction, mps_fire::i, dqmiolumiharvest::j, LogDebug, N, point, funct::pow(), edm::second(), histoStyle::weight, gpuClustering::x, reco::CaloCluster::x(), detailsBasic3DVector::y, reco::CaloCluster::y(), and reco::CaloCluster::z().

Referenced by ticl::PatternRecognitionbyCLUE3D< TILES >::makeTracksters(), ticl::PatternRecognitionbyCA< TILES >::makeTracksters(), TrackstersMergeProducer::produce(), and SimTrackstersProducer::produce().

16  {
17  LogDebug("TrackstersPCA_Eigen") << "------- Eigen -------" << std::endl;
18 
19  for (auto &trackster : tracksters) {
21  point << 0., 0., 0.;
22  Eigen::Vector3d barycenter;
23  barycenter << 0., 0., 0.;
24 
25  auto fillPoint = [&](const reco::CaloCluster &c, const float weight = 1.f) {
26  point[0] = weight * c.x();
27  point[1] = weight * c.y();
28  point[2] = weight * c.z();
29  };
30 
31  // Initialize this trackster with default, dummy values
32  trackster.setRawEnergy(0.f);
33  trackster.setRawEmEnergy(0.f);
34  trackster.setRawPt(0.f);
35  trackster.setRawEmPt(0.f);
36 
37  size_t N = trackster.vertices().size();
38  float weight = 1.f / N;
39  float weights2_sum = 0.f;
40  Eigen::Vector3d sigmas;
41  sigmas << 0., 0., 0.;
42  Eigen::Vector3d sigmasEigen;
43  sigmasEigen << 0., 0., 0.;
44  Eigen::Matrix3d covM = Eigen::Matrix3d::Zero();
45 
46  std::vector<float> times;
47  std::vector<float> timeErrors;
48  std::set<uint32_t> usedLC;
49 
50  for (size_t i = 0; i < N; ++i) {
51  auto fraction = 1.f / trackster.vertex_multiplicity(i);
52  trackster.addToRawEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
53  if (std::abs(layerClusters[trackster.vertices(i)].z()) <= z_limit_em)
54  trackster.addToRawEmEnergy(layerClusters[trackster.vertices(i)].energy() * fraction);
55 
56  // Compute the weighted barycenter.
57  if (energyWeight)
58  weight = layerClusters[trackster.vertices(i)].energy() * fraction;
59  fillPoint(layerClusters[trackster.vertices(i)], weight);
60  for (size_t j = 0; j < 3; ++j)
61  barycenter[j] += point[j];
62 
63  // Add timing from layerClusters not already used
64  if ((usedLC.insert(trackster.vertices(i))).second) {
65  float timeE = layerClustersTime.get(trackster.vertices(i)).second;
66  if (timeE > 0.f) {
67  times.push_back(layerClustersTime.get(trackster.vertices(i)).first);
68  timeErrors.push_back(1. / pow(timeE, 2));
69  }
70  }
71  }
72  if (energyWeight && trackster.raw_energy())
73  barycenter /= trackster.raw_energy();
74 
76  std::pair<float, float> timeTrackster = timeEstimator.fixSizeHighestDensity(times, timeErrors);
77 
78  // Compute the Covariance Matrix and the sum of the squared weights, used
79  // to compute the correct normalization.
80  // The barycenter has to be known.
81  for (size_t i = 0; i < N; ++i) {
82  fillPoint(layerClusters[trackster.vertices(i)]);
83  if (energyWeight && trackster.raw_energy())
84  weight =
85  (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) / trackster.raw_energy();
86  weights2_sum += weight * weight;
87  for (size_t x = 0; x < 3; ++x)
88  for (size_t y = 0; y <= x; ++y) {
89  covM(x, y) += weight * (point[x] - barycenter[x]) * (point[y] - barycenter[y]);
90  covM(y, x) = covM(x, y);
91  }
92  }
93  covM *= 1. / (1. - weights2_sum);
94 
95  // Perform the actual decomposition
96  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d>::RealVectorType eigenvalues_fromEigen;
97  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d>::EigenvectorsType eigenvectors_fromEigen;
98  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eigensolver(covM);
99  if (eigensolver.info() != Eigen::Success) {
100  eigenvalues_fromEigen = eigenvalues_fromEigen.Zero();
101  eigenvectors_fromEigen = eigenvectors_fromEigen.Zero();
102  } else {
103  eigenvalues_fromEigen = eigensolver.eigenvalues();
104  eigenvectors_fromEigen = eigensolver.eigenvectors();
105  }
106 
107  // Compute the spread in the both spaces.
108  for (size_t i = 0; i < N; ++i) {
109  fillPoint(layerClusters[trackster.vertices(i)]);
110  sigmas += weight * (point - barycenter).cwiseAbs2();
111  Eigen::Vector3d point_transformed = eigenvectors_fromEigen * (point - barycenter);
112  if (energyWeight && trackster.raw_energy())
113  weight =
114  (layerClusters[trackster.vertices(i)].energy() / trackster.vertex_multiplicity(i)) / trackster.raw_energy();
115  sigmasEigen += weight * (point_transformed.cwiseAbs2());
116  }
117  sigmas /= (1. - weights2_sum);
118  sigmasEigen /= (1. - weights2_sum);
119 
120  // Add trackster attributes
121  trackster.setBarycenter(ticl::Trackster::Vector(barycenter));
122  trackster.setTimeAndError(timeTrackster.first, timeTrackster.second);
123  trackster.fillPCAVariables(
124  eigenvalues_fromEigen, eigenvectors_fromEigen, sigmas, sigmasEigen, 3, ticl::Trackster::PCAOrdering::ascending);
125 
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  LogDebug("TrackstersPCA") << "EigenValues from Eigen/Tr(cov): " << eigenvalues_fromEigen[2] / covM.trace() << ", "
135  << eigenvalues_fromEigen[1] / covM.trace() << ", "
136  << eigenvalues_fromEigen[0] / covM.trace() << std::endl;
137  LogDebug("TrackstersPCA") << "EigenValues from Eigen: " << eigenvalues_fromEigen[2] << ", "
138  << eigenvalues_fromEigen[1] << ", " << eigenvalues_fromEigen[0] << std::endl;
139  LogDebug("TrackstersPCA") << "EigenVector 3 from Eigen: " << eigenvectors_fromEigen(0, 2) << ", "
140  << eigenvectors_fromEigen(1, 2) << ", " << eigenvectors_fromEigen(2, 2) << std::endl;
141  LogDebug("TrackstersPCA") << "EigenVector 2 from Eigen: " << eigenvectors_fromEigen(0, 1) << ", "
142  << eigenvectors_fromEigen(1, 1) << ", " << eigenvectors_fromEigen(2, 1) << std::endl;
143  LogDebug("TrackstersPCA") << "EigenVector 1 from Eigen: " << eigenvectors_fromEigen(0, 0) << ", "
144  << eigenvectors_fromEigen(1, 0) << ", " << eigenvectors_fromEigen(2, 0) << std::endl;
145  LogDebug("TrackstersPCA") << "Original sigmas: " << sigmas[0] << ", " << sigmas[1] << ", " << sigmas[2]
146  << std::endl;
147  LogDebug("TrackstersPCA") << "SigmasEigen in PCA space: " << sigmasEigen[2] << ", " << sigmasEigen[1] << ", "
148  << sigmasEigen[0] << std::endl;
149  LogDebug("TrackstersPCA") << "covM: \n" << covM << std::endl;
150  }
151 }
const edm::EventSetup & c
math::XYZVector Vector
Definition: Trackster.h:21
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:178
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)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint16_t const *__restrict__ x
Definition: gpuClustering.h:39
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:172
const_reference_type get(ProductID id, size_t idx) const
Definition: ValueMap.h:144
#define N
Definition: blowfish.cc:9
double y() const
y coordinate of cluster centroid
Definition: CaloCluster.h:175
int weight
Definition: histoStyle.py:51
Eigen::Matrix3d Matrix3d
Definition: FitResult.h:18
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
Eigen::Vector3d Vector3d
Definition: FitResult.h:14
#define LogDebug(id)
Trackster::ParticleType ticl::tracksterParticleTypeFromPdgId ( int  pdgId,
int  charge 
)
inline

Definition at line 10 of file commons.h.

References funct::abs(), ticl::Trackster::charged_hadron, ticl::Trackster::electron, ticl::Trackster::muon, ticl::Trackster::neutral_hadron, ticl::Trackster::neutral_pion, ticl::Trackster::photon, and ticl::Trackster::unknown.

Referenced by addTrackster().

10  {
11  if (pdgId == 111) {
12  return Trackster::ParticleType::neutral_pion;
13  } else {
14  pdgId = std::abs(pdgId);
15  if (pdgId == 22) {
16  return Trackster::ParticleType::photon;
17  } else if (pdgId == 11) {
19  } else if (pdgId == 13) {
21  } else {
22  bool isHadron = (pdgId > 100 and pdgId < 900) or (pdgId > 1000 and pdgId < 9000);
23  if (isHadron) {
24  if (charge != 0) {
25  return Trackster::ParticleType::charged_hadron;
26  } else {
27  return Trackster::ParticleType::neutral_hadron;
28  }
29  } else {
31  }
32  }
33  }
34  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22