CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalClusterLazyTools.h
Go to the documentation of this file.
1 #ifndef RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
2 #define RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
3 
18 
19 
20 class CaloTopology;
21 class CaloGeometry;
22 
24  public:
25  EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, edm::InputTag redEBRecHits, edm::InputTag redEERecHits );
27 
28  // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster
29  //NOTE (29/10/08): we now use an eta/phi coordinate system rather than phi/eta
30  //to minmise possible screwups, for now e5x1 isnt defined all the majority of people who call it actually want e1x5 and
31  //it is thought it is better that their code doesnt compile rather than pick up the wrong function
32  //therefore in this version and later e1x5 = e5x1 in the old version
33  //so 1x5 is 1 crystal in eta and 5 crystals in phi
34  //note e3x2 does not have a definate eta/phi geometry, it takes the maximum 3x2 block containing the
35  //seed regardless of whether that 3 in eta or phi
36  float e1x3( const reco::BasicCluster &cluster );
37  float e3x1( const reco::BasicCluster &cluster );
38  float e1x5( const reco::BasicCluster &cluster );
39  //float e5x1( const reco::BasicCluster &cluster );
40  float e2x2( const reco::BasicCluster &cluster );
41  float e3x2( const reco::BasicCluster &cluster );
42  float e3x3( const reco::BasicCluster &cluster );
43  float e4x4( const reco::BasicCluster &cluster );
44  float e5x5( const reco::BasicCluster &cluster );
45  // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
46  // 2 crystals wide in eta, 5 wide in phi.
47  float e2x5Right( const reco::BasicCluster &cluster );
48  // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
49  float e2x5Left( const reco::BasicCluster &cluster );
50  // energy in the 5x2 strip above the max crystal (does not contain max crystal)
51  // 5 crystals wide in eta, 2 wide in phi.
52  float e2x5Top( const reco::BasicCluster &cluster );
53  // energy in the 5x2 strip below the max crystal (does not contain max crystal)
54  float e2x5Bottom( const reco::BasicCluster &cluster );
55  // energy in a 2x5 strip containing the seed (max) crystal.
56  // 2 crystals wide in eta, 5 wide in phi.
57  // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
58  float e2x5Max( const reco::BasicCluster &cluster );
59  // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
60  float eLeft( const reco::BasicCluster &cluster );
61  float eRight( const reco::BasicCluster &cluster );
62  float eTop( const reco::BasicCluster &cluster );
63  float eBottom( const reco::BasicCluster &cluster );
64  // the energy of the most energetic crystal in the cluster
65  float eMax( const reco::BasicCluster &cluster );
66  // the energy of the second most energetic crystal in the cluster
67  float e2nd( const reco::BasicCluster &cluster );
68  // get the DetId and the energy of the maximum energy crystal of the input cluster
69  std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster );
70  std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster );
71  std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster );
72  // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
73  std::vector<float> lat( const reco::BasicCluster &cluster, bool logW = true, float w0 = 4.7 );
74  // return a vector v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
75  std::vector<float> covariances(const reco::BasicCluster &cluster, float w0 = 4.7 );
76  // return a vector v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
77  //this function calculates differences in eta/phi in units of crystals not global eta/phi
78  //this is gives better performance in the crack regions of the calorimeter but gives otherwise identical results to covariances function
79  //this is only defined for the barrel, it returns covariances when the cluster is in the endcap
80  //Warning: covIEtaIEta has been studied by egamma, but so far covIPhiIPhi hasnt been studied extensively so there could be a bug in
81  // the covIPhiIEta or covIPhiIPhi calculations. I dont think there is but as it hasnt been heavily used, there might be one
82  std::vector<float> localCovariances(const reco::BasicCluster &cluster, float w0 = 4.7);
83  std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, float w0 = 4.7);
84 
85  double zernike20( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
86  double zernike42( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
87 
88  // get the detId's of a matrix centered in the maximum energy crystal = (0,0)
89  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
90  std::vector<DetId> matrixDetId( DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
91  // get the energy deposited in a matrix centered in the maximum energy crystal = (0,0)
92  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
93  float matrixEnergy( const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
94 
95  private:
96  void getGeometry( const edm::EventSetup &es );
97  void getTopology( const edm::EventSetup &es );
98  void getEBRecHits( const edm::Event &ev, edm::InputTag redEBRecHits );
99  void getEERecHits( const edm::Event &ev, edm::InputTag redEERecHits );
101 
106 };
107 
108 #endif
std::vector< float > localCovariances(const reco::BasicCluster &cluster, float w0=4.7)
void getEBRecHits(const edm::Event &ev, edm::InputTag redEBRecHits)
float e2x2(const reco::BasicCluster &cluster)
float e4x4(const reco::BasicCluster &cluster)
std::vector< float > energyBasketFractionPhi(const reco::BasicCluster &cluster)
float eBottom(const reco::BasicCluster &cluster)
EcalClusterLazyTools(const edm::Event &ev, const edm::EventSetup &es, edm::InputTag redEBRecHits, edm::InputTag redEERecHits)
float e2nd(const reco::BasicCluster &cluster)
float eRight(const reco::BasicCluster &cluster)
float e1x5(const reco::BasicCluster &cluster)
float e3x3(const reco::BasicCluster &cluster)
float eMax(const reco::BasicCluster &cluster)
const CaloGeometry * geometry_
void getEERecHits(const edm::Event &ev, edm::InputTag redEERecHits)
const EcalRecHitCollection * ebRecHits_
float matrixEnergy(const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax)
float e2x5Max(const reco::BasicCluster &cluster)
float e2x5Right(const reco::BasicCluster &cluster)
double zernike20(const reco::BasicCluster &cluster, double R0=6.6, bool logW=true, float w0=4.7)
double zernike42(const reco::BasicCluster &cluster, double R0=6.6, bool logW=true, float w0=4.7)
const EcalRecHitCollection * eeRecHits_
std::vector< float > scLocalCovariances(const reco::SuperCluster &cluster, float w0=4.7)
float e3x2(const reco::BasicCluster &cluster)
std::vector< float > covariances(const reco::BasicCluster &cluster, float w0=4.7)
float eTop(const reco::BasicCluster &cluster)
float e2x5Bottom(const reco::BasicCluster &cluster)
std::vector< float > lat(const reco::BasicCluster &cluster, bool logW=true, float w0=4.7)
Definition: DetId.h:20
float e2x5Left(const reco::BasicCluster &cluster)
float e5x5(const reco::BasicCluster &cluster)
float eLeft(const reco::BasicCluster &cluster)
std::pair< DetId, float > getMaximum(const reco::BasicCluster &cluster)
float e2x5Top(const reco::BasicCluster &cluster)
void getGeometry(const edm::EventSetup &es)
float e3x1(const reco::BasicCluster &cluster)
std::vector< DetId > matrixDetId(DetId id, int ixMin, int ixMax, int iyMin, int iyMax)
const CaloTopology * topology_
const EcalRecHitCollection * getEcalRecHitCollection(const reco::BasicCluster &cluster)
void getTopology(const edm::EventSetup &es)
std::vector< float > energyBasketFractionEta(const reco::BasicCluster &cluster)
float e1x3(const reco::BasicCluster &cluster)