CMS 3D CMS Logo

Cluster3DPCACalculator.cc
Go to the documentation of this file.
3 
4 #include <cmath>
5 #include <unordered_map>
6 
7 #include "vdt/vdtMath.h"
8 
12 
14 
15 #include "TPrincipal.h"
16 
18  public:
21  pca_(new TPrincipal(3,"D")){
22  }
25 
28 
29  private:
30  std::unique_ptr<TPrincipal> pca_;
31 
33  math::XYZVector& );
34 
36 };
37 
40  "Cluster3DPCACalculator");
41 
44  pca_.reset(new TPrincipal(3,"D"));
46 }
47 
50  for( reco::PFCluster& cluster : clusters ) {
51  pca_.reset(new TPrincipal(3,"D"));
53  }
54 }
55 
58  if( !cluster.seed() ) {
59  throw cms::Exception("ClusterWithNoSeed")
60  << " Found a cluster with no seed: " << cluster;
61  }
62  double cl_energy = 0;
63  double max_e = 0.0;
64  double avg_time = 0.0;
65  double time_norm = 0.0;
66  PFLayer::Layer max_e_layer = PFLayer::NONE;
67  reco::PFRecHitRef refseed;
68  double pcavars[3];
69 
70  for( const reco::PFRecHitFraction& rhf : cluster.recHitFractions() ) {
71  const reco::PFRecHitRef& refhit = rhf.recHitRef();
72  double rh_energy = refhit->energy();
73  double rh_time = refhit->time();
74  cl_energy += rh_energy * rhf.fraction();
75  if( rh_time > 0.0 ) { // time == -1 means no measurement
76  // all times are offset by one nanosecond in digitizer
77  // remove that here so all times of flight
78  // are with respect to (0,0,0)
79  avg_time += (rh_time - 1.0);
80  time_norm += 1.0;
81  }
82  if( rh_energy > max_e ) {
83  max_e = rh_energy;
84  max_e_layer = rhf.recHitRef()->layer();
85  }
86  if( refhit->detId() == cluster.seed() ) refseed = refhit;
87  const double rh_fraction = rhf.fraction();
88  rh_energy = refhit->energy()*rh_fraction;
89  if( edm::isNotFinite(rh_energy) ) {
90  throw cms::Exception("PFClusterAlgo")
91  <<"rechit " << refhit->detId() << " has a NaN energy... "
92  << "The input of the particle flow clustering seems to be corrupted.";
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  }
104  cluster.setEnergy(cl_energy);
105  cluster.setLayer(max_e_layer);
106  // calculate the position
107 
108  pca_->MakePrincipals();
109  const TVectorD& means = *(pca_->GetMeanValues());
110  const TMatrixD& eigens = *(pca_->GetEigenVectors());
111 
112  math::XYZPoint barycenter(means[0],means[1],means[2]);
113  math::XYZVector axis(eigens(0,0),eigens(1,0),eigens(2,0));
114 
115  if( time_norm > 0.0 ) {
116  avg_time = avg_time/time_norm;
117  } else {
119  }
120 
121  if( axis.z()*barycenter.z() < 0.0 ) {
122  axis = math::XYZVector(-eigens(0,0),-eigens(1,0),-eigens(2,0));
123  }
124 
125  cluster.setTime(avg_time);
126  cluster.setPosition(barycenter);
127  cluster.calculatePositionREP();
128 
129 }
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:120
Cluster3DPCACalculator & operator=(const Cluster3DPCACalculator &)=delete
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:115
void setTime(float time, float timeError=0)
Definition: PFCluster.h:92
void showerParameters(const reco::PFCluster &, math::XYZPoint &, math::XYZVector &)
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:111
bool isNotFinite(T x)
Definition: isFinite.h:10
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:100
Cluster3DPCACalculator(const edm::ParameterSet &conf)
T min(T a, T b)
Definition: MathUtil.h:58
void calculateAndSetPositionActual(reco::PFCluster &)
Layer
layer definition
Definition: PFLayer.h:31
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:205
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
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
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
#define DEFINE_EDM_PLUGIN(factory, type, name)