CMS 3D CMS Logo

EgammaPCAHelper.h
Go to the documentation of this file.
1 //--------------------------------------------------------------------------------------------------
2 //
3 // EGammaPCAHelper
4 //
5 // Helper Class to compute PCA
6 //
7 //
8 //--------------------------------------------------------------------------------------------------
9 #ifndef RecoEgamma_EgammaTools_EGammaPCAHelper_h
10 #define RecoEgamma_EgammaTools_EGammaPCAHelper_h
11 
13 
19 
24 #include "Math/Transform3D.h"
25 #include <unordered_map>
26 
27 #include "TPrincipal.h"
28 
29 class HGCalRecHit;
30 
31 namespace hgcal {
32 
34  public:
35  typedef ROOT::Math::Transform3D Transform3D;
37 
39 
40  // for the GsfElectrons
41  void storeRecHits(const reco::CaloCluster &theCluster);
42  void storeRecHits(const reco::HGCalMultiCluster &cluster);
43 
44  const TPrincipal &pcaResult();
46  void setHitMap(const std::unordered_map<DetId, const unsigned int> *hitMap);
47  const std::unordered_map<DetId, const unsigned int> *getHitMap() { return hitMap_; }
48 
49  void setRecHitTools(const hgcal::RecHitTools *recHitTools);
50  void setRecHits(edm::Handle<HGCRecHitCollection> recHitHandleEE,
51  edm::Handle<HGCRecHitCollection> recHitHandleFH,
52  edm::Handle<HGCRecHitCollection> recHitHandleBH);
53 
54  inline void setdEdXWeights(const std::vector<double> &dEdX) { dEdXWeights_ = dEdX; }
55 
56  void pcaInitialComputation() { computePCA(-1., false); }
57 
58  void computePCA(float radius, bool withHalo = true);
59  const math::XYZPoint &barycenter() const { return barycenter_; }
60  const math::XYZVector &axis() const { return axis_; }
61 
62  void computeShowerWidth(float radius, bool withHalo = true);
63 
64  inline double sigmaUU() const { return checkIteration() ? sigu_ : -1.; }
65  inline double sigmaVV() const { return checkIteration() ? sigv_ : -1.; }
66  inline double sigmaEE() const { return checkIteration() ? sige_ : -1.; }
67  inline double sigmaPP() const { return checkIteration() ? sigp_ : -1.; }
68 
69  inline const TVectorD &eigenValues() const { return *pca_->GetEigenValues(); }
70  inline const TVectorD &sigmas() const { return *pca_->GetSigmas(); }
71  // contains maxlayer+1 values, first layer is [1]
72  LongDeps energyPerLayer(float radius, bool withHalo = true);
73 
74  float clusterDepthCompatibility(const LongDeps &, float &measuredDepth, float &expectedDepth, float &expectedSigma);
75  void printHits(float radius) const;
76  void clear();
77 
78  private:
79  bool checkIteration() const;
80  void storeRecHits(const std::vector<std::pair<DetId, float>> &hf);
81  float findZFirstLayer(const LongDeps &) const;
82 
84  bool debug_;
85 
86  //parameters
87  std::vector<double> dEdXWeights_;
88  std::vector<double> invThicknessCorrection_;
89 
91  const std::unordered_map<DetId, const unsigned int> *hitMap_;
92  std::vector<Spot> theSpots_;
94  unsigned int maxlayer_;
95 
96  // output quantities
99 
101  double sigu_, sigv_, sige_, sigp_;
102 
103  // helper
104  std::unique_ptr<TPrincipal> pca_;
107 
108  std::vector<const HGCRecHit *> hits_;
109  };
110 
111 } // namespace hgcal
112 
113 #endif
void computeShowerWidth(float radius, bool withHalo=true)
const TVectorD & sigmas() const
std::vector< const HGCRecHit * > hits_
void printHits(float radius) const
math::XYZPoint barycenter_
void computePCA(float radius, bool withHalo=true)
LongDeps energyPerLayer(float radius, bool withHalo=true)
double sigmaEE() const
std::unique_ptr< TPrincipal > pca_
double sigmaVV() const
void setdEdXWeights(const std::vector< double > &dEdX)
const math::XYZPoint & barycenter() const
const TVectorD & eigenValues() const
void storeRecHits(const reco::CaloCluster &theCluster)
math::XYZPoint Point
std::vector< Spot > theSpots_
const std::unordered_map< DetId, const unsigned int > * hitMap_
ROOT::Math::Transform3D Transform3D
const std::unordered_map< DetId, const unsigned int > * getHitMap()
double sigmaPP() const
const math::XYZVector & axis() const
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
float clusterDepthCompatibility(const LongDeps &, float &measuredDepth, float &expectedDepth, float &expectedSigma)
float findZFirstLayer(const LongDeps &) const
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 setHitMap(const std::unordered_map< DetId, const unsigned int > *hitMap)
to set once per event
math::XYZVector axis_
std::vector< double > invThicknessCorrection_
std::vector< double > dEdXWeights_
const reco::CaloCluster * theCluster_
const TPrincipal & pcaResult()
ROOT::Math::Transform3D::Point Point
double sigmaUU() const
void setRecHits(edm::Handle< HGCRecHitCollection > recHitHandleEE, edm::Handle< HGCRecHitCollection > recHitHandleFH, edm::Handle< HGCRecHitCollection > recHitHandleBH)
const hgcal::RecHitTools * recHitTools_