CMS 3D CMS Logo

Cluster3DPCACalculator.cc
Go to the documentation of this file.
3 
4 #include <cmath>
5 #include <memory>
6 
7 #include <unordered_map>
8 
9 #include "vdt/vdtMath.h"
10 
15 
17 
18 #include "TPrincipal.h"
19 
21 public:
24  updateTiming_(conf.getParameter<bool>("updateTiming")),
25  pca_(new TPrincipal(3, "D")) {}
28 
31 
32 private:
33  const bool updateTiming_;
34  std::unique_ptr<TPrincipal> pca_;
35 
37 
39 };
40 
42 
44  pca_ = std::make_unique<TPrincipal>(3, "D");
46 }
47 
49  for (reco::PFCluster& cluster : clusters) {
50  pca_ = std::make_unique<TPrincipal>(3, "D");
52  }
53 }
54 
56  if (!cluster.seed()) {
57  throw cms::Exception("ClusterWithNoSeed") << " Found a cluster with no seed: " << cluster;
58  }
59  double cl_energy = 0;
60  double max_e = 0.0;
61  double avg_time = 0.0;
62  double time_norm = 0.0;
63  PFLayer::Layer max_e_layer = PFLayer::NONE;
64  reco::PFRecHitRef refseed;
65  double pcavars[3];
66 
67  for (const reco::PFRecHitFraction& rhf : cluster.recHitFractions()) {
68  const reco::PFRecHitRef& refhit = rhf.recHitRef();
69  double rh_energy = refhit->energy();
70  double rh_time = refhit->time();
71  cl_energy += rh_energy * rhf.fraction();
72  if (rh_time > 0.0) { // time == -1 means no measurement
73  // all times are offset by one nanosecond in digitizer
74  // remove that here so all times of flight
75  // are with respect to (0,0,0)
76  avg_time += (rh_time - 1.0);
77  time_norm += 1.0;
78  }
79  if (rh_energy > max_e) {
80  max_e = rh_energy;
81  max_e_layer = rhf.recHitRef()->layer();
82  }
83  if (refhit->detId() == cluster.seed())
84  refseed = refhit;
85  const double rh_fraction = rhf.fraction();
86  rh_energy = refhit->energy() * rh_fraction;
87  if (edm::isNotFinite(rh_energy)) {
88  //temporarily changed exception to warning
89  // throw cms::Exception("PFClusterAlgo")
90  edm::LogWarning("PFClusterAlgo") << "rechit " << refhit->detId() << " has a NaN energy... "
91  << "The input of the particle flow clustering seems to be corrupted.";
92  continue;
93  }
94  pcavars[0] = refhit->position().x();
95  pcavars[1] = refhit->position().y();
96  pcavars[2] = refhit->position().z();
97  int nhit = int(rh_energy * 100); // put rec_hit energy in units of 10 MeV
98 
99  for (int i = 0; i < nhit; ++i) {
100  pca_->AddRow(pcavars);
101  }
102  }
103  cluster.setEnergy(cl_energy);
104  cluster.setLayer(max_e_layer);
105  // calculate the position
106 
107  pca_->MakePrincipals();
108  const TVectorD& means = *(pca_->GetMeanValues());
109  const TMatrixD& eigens = *(pca_->GetEigenVectors());
110 
111  math::XYZPoint barycenter(means[0], means[1], means[2]);
112  math::XYZVector axis(eigens(0, 0), eigens(1, 0), eigens(2, 0));
113 
114  if (time_norm > 0.0) {
115  avg_time = avg_time / time_norm;
116  } else {
118  }
119 
120  if (axis.z() * barycenter.z() < 0.0) {
121  axis = math::XYZVector(-eigens(0, 0), -eigens(1, 0), -eigens(2, 0));
122  }
123 
124  if (updateTiming_) {
125  cluster.setTime(avg_time);
126  }
127  cluster.setPosition(barycenter);
128  cluster.calculatePositionREP();
129 }
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:49
Cluster3DPCACalculator & operator=(const Cluster3DPCACalculator &)=delete
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:139
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:218
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
void showerParameters(const reco::PFCluster &, math::XYZPoint &, math::XYZVector &)
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
void setEnergy(double energy)
Definition: CaloCluster.h:135
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
void calculateAndSetPositionActual(reco::PFCluster &)
Layer
layer definition
Definition: PFLayer.h:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void calculateAndSetPositions(reco::PFClusterCollection &) override
std::unique_ptr< TPrincipal > pca_
void calculateAndSetPosition(reco::PFCluster &) override
Cluster3DPCACalculator(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
Log< level::Warning, false > LogWarning