CMS 3D CMS Logo

ClusterShapeAlgo.h
Go to the documentation of this file.
1 #ifndef RecoEcal_EgammaCoreTools_ClusterShapeAlgo_h
2 #define RecoEcal_EgammaCoreTools_ClusterShapeAlgo_h
3 
13 #include <map>
14 
16 
26 
28 
30 {
32  double r;
33  double phi;
34 };
35 
37 {
38 
39  public:
42  reco::ClusterShape Calculate(const reco::BasicCluster &passedCluster,
46 
47  private:
48  void Calculate_TopEnergy(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits);
49  void Calculate_2ndEnergy(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits);
50  void Create_Map(const EcalRecHitCollection *hits, const CaloSubdetectorTopology* topology);
51  void Calculate_e2x2();
52  void Calculate_e3x2();
53  void Calculate_e3x3();
54  void Calculate_e4x4();
55  void Calculate_e5x5();
56  void Calculate_e2x5Right();
57  void Calculate_e2x5Left();
58  void Calculate_e2x5Top();
59  void Calculate_e2x5Bottom();
60  void Calculate_Covariances(const reco::BasicCluster &passedCluster,
61  const EcalRecHitCollection* hits,
62  const CaloSubdetectorGeometry* geometry);
63  void Calculate_BarrelBasketEnergyFraction(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits,
64  const int EtaPhi,const CaloSubdetectorGeometry * geometry);
65  // defines a energy deposition topology in a reference system centered on the cluster
66  void Calculate_EnergyDepTopology(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits, const CaloSubdetectorGeometry * geometry, bool logW=true);
67  void Calculate_Polynomials(double rho);
68  double factorial(int n) const;
69  void Calculate_lat(const reco::BasicCluster &passedCluster);
70  void Calculate_ComplexZernikeMoments(const reco::BasicCluster &passedCluster);
71  // explicit implementation of polynomial part of
72  // Zernike-Functions for n<=5;
73  double f00(double r);
74  double f11(double r);
75  double f20(double r);
76  double f22(double r);
77  double f31(double r);
78  double f33(double r);
79  double f40(double r);
80  double f42(double r);
81  double f44(double r);
82  double f51(double r);
83  double f53(double r);
84  double f55(double r);
85  double absZernikeMoment(const reco::BasicCluster &passedCluster, int n, int m, double R0=6.6);
86  double fast_AbsZernikeMoment(const reco::BasicCluster &passedCluster, int n, int m, double R0);
87  // Calculation of Zernike-Moments for general values of (n,m)
88  double calc_AbsZernikeMoment(const reco::BasicCluster &passedCluster, int n, int m, double R0);
89 
91 
92  std::pair<DetId, double> energyMap_[5][5];
93  int e2x2_Diagonal_X_, e2x2_Diagonal_Y_;
94 
95  double covEtaEta_, covEtaPhi_, covPhiPhi_;
96  double eMax_, e2nd_, e2x2_, e3x2_, e3x3_, e4x4_, e5x5_;
97  double e2x5Right_, e2x5Left_, e2x5Top_, e2x5Bottom_;
98  double e3x2Ratio_;
99  double lat_;
100  double etaLat_ ;
101  double phiLat_ ;
102  double A20_, A42_;
103  std::vector<double> energyBasketFractionEta_;
104  std::vector<double> energyBasketFractionPhi_;
105  DetId eMaxId_, e2ndId_;
106  std::vector<EcalClusterEnergyDeposition> energyDistribution_;
107  std::vector<double> fcn_;
108 
109  enum { Eta, Phi };
110 
111 };
112 
113 #endif
std::vector< double > fcn_
CaloTopology const * topology(0)
std::vector< EcalClusterEnergyDeposition > energyDistribution_
edm::ParameterSet parameterSet_
Helper class for the calculation of a top and a W boson mass estime.
Definition: DetId.h:18
std::vector< double > energyBasketFractionPhi_
int factorial(int n)
factorial function
std::vector< double > energyBasketFractionEta_