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 
20 
27 
28 
29 class CaloTopology;
30 class CaloGeometry;
32 
34  public:
35  EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, const edm::InputTag& redEBRecHits, const edm::InputTag& redEERecHits,const edm::ParameterSet& config );
36  EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, const edm::InputTag& redEBRecHits, const edm::InputTag& redEERecHits);
37  EcalClusterLazyTools( const edm::Event &ev, const edm::EventSetup &es, const edm::InputTag& redEBRecHits, const edm::InputTag& redEERecHits, const edm::InputTag& redESRecHits);
38 
40 
41  // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster
42  //NOTE (29/10/08): we now use an eta/phi coordinate system rather than phi/eta
43  //to minmise possible screwups, for now e5x1 isnt defined all the majority of people who call it actually want e1x5 and
44  //it is thought it is better that their code doesnt compile rather than pick up the wrong function
45  //therefore in this version and later e1x5 = e5x1 in the old version
46  //so 1x5 is 1 crystal in eta and 5 crystals in phi
47  //note e3x2 does not have a definate eta/phi geometry, it takes the maximum 3x2 block containing the
48  //seed regardless of whether that 3 in eta or phi
49  float e1x3( const reco::BasicCluster &cluster );
50  float e1x3( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
51 
52  float e3x1( const reco::BasicCluster &cluster );
53  float e3x1( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
54 
55  float e1x5( const reco::BasicCluster &cluster );
56  float e1x5( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
57 
58  float e5x1( const reco::BasicCluster &cluster );
59  float e5x1( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
60 
61  float e2x2( const reco::BasicCluster &cluster );
62  float e2x2( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
63 
64  float e3x2( const reco::BasicCluster &cluster );
65  float e3x2( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
66 
67  float e3x3( const reco::BasicCluster &cluster );
68  float e3x3( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
69 
70  float e4x4( const reco::BasicCluster &cluster );
71  float e4x4( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
72 
73  float e5x5( const reco::BasicCluster &cluster );
74  float e5x5( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
75  // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
76  // 2 crystals wide in eta, 5 wide in phi.
77  float e2x5Right( const reco::BasicCluster &cluster );
78  float e2x5Right( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
79  // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
80  float e2x5Left( const reco::BasicCluster &cluster );
81  float e2x5Left( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
82  // energy in the 5x2 strip above the max crystal (does not contain max crystal)
83  // 5 crystals wide in eta, 2 wide in phi.
84  float e2x5Top( const reco::BasicCluster &cluster );
85  float e2x5Top( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
86 
87  // energy in the 5x2 strip below the max crystal (does not contain max crystal)
88  float e2x5Bottom( const reco::BasicCluster &cluster );
89  float e2x5Bottom( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
90  // energy in a 2x5 strip containing the seed (max) crystal.
91  // 2 crystals wide in eta, 5 wide in phi.
92  // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
93  float e2x5Max( const reco::BasicCluster &cluster );
94  float e2x5Max( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
95  // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
96  float eLeft( const reco::BasicCluster &cluster );
97  float eLeft( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
98 
99  float eRight( const reco::BasicCluster &cluster );
100  float eRight( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
101 
102  float eTop( const reco::BasicCluster &cluster );
103  float eTop( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
104 
105  float eBottom( const reco::BasicCluster &cluster );
106  float eBottom( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
107  // the energy of the most energetic crystal in the cluster
108  float eMax( const reco::BasicCluster &cluster );
109  float eMax( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
110  // the energy of the second most energetic crystal in the cluster
111  float e2nd( const reco::BasicCluster &cluster );
112  float e2nd( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
113 
114  // get the DetId and the energy of the maximum energy crystal of the input cluster
115  std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster );
116  std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
117 
118  std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster );
119  std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster,const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
120 
121  std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster );
122  std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
123 
124  // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
125  std::vector<float> lat( const reco::BasicCluster &cluster, bool logW = true, float w0 = 4.7 );
126  std::vector<float> lat( const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv, bool logW = true, float w0 = 4.7 );
127 
128  // return a vector v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
129  std::vector<float> covariances(const reco::BasicCluster &cluster, float w0 = 4.7 );
130  std::vector<float> covariances(const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7 );
131 
132  // return a vector v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
133  //this function calculates differences in eta/phi in units of crystals not global eta/phi
134  //this is gives better performance in the crack regions of the calorimeter but gives otherwise identical results to covariances function
135  //this is only defined for the barrel, it returns covariances when the cluster is in the endcap
136  //Warning: covIEtaIEta has been studied by egamma, but so far covIPhiIPhi hasnt been studied extensively so there could be a bug in
137  // the covIPhiIEta or covIPhiIPhi calculations. I dont think there is but as it hasnt been heavily used, there might be one
138  std::vector<float> localCovariances(const reco::BasicCluster &cluster, float w0 = 4.7);
139  std::vector<float> localCovariances(const reco::BasicCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7);
140 
141  std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, float w0 = 4.7);
142  std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0 = 4.7);
143 
144  double zernike20( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
145  double zernike42( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 );
146 
147  // get the detId's of a matrix centered in the maximum energy crystal = (0,0)
148  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
149  std::vector<DetId> matrixDetId( DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
150  // get the energy deposited in a matrix centered in the maximum energy crystal = (0,0)
151  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
152  float matrixEnergy( const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax );
153  float matrixEnergy( const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax, const std::vector<int>& flagsexcl, const std::vector<int>& severitiesexcl, const EcalSeverityLevelAlgo *sevLv );
154 
155  // get time of basic cluster seed crystal
156  float BasicClusterSeedTime(const reco::BasicCluster &cluster);
157  // error-weighted average of time from constituents of basic cluster
158  float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev);
159  // get BasicClusterSeedTime of the seed basic cluser of the supercluster
160  float SuperClusterSeedTime(const reco::SuperCluster &cluster);
161  // get BasicClusterTime of the seed basic cluser of the supercluster
162  float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev);
163 
164  // mapping for preshower rechits
165  std::map<DetId, EcalRecHit> rechits_map_;
166  // get Preshower hit array
167  std::vector<float> getESHits(double X, double Y, double Z, std::map<DetId, EcalRecHit> rechits_map, const CaloGeometry* geometry, CaloSubdetectorTopology *topology_p, int row=0, int plane=1);
168  // get Preshower hit shape
169  float getESShape(const std::vector<float>& ESHits0);
170  // get Preshower effective sigmaRR
171  float eseffsirir( const reco::SuperCluster &cluster );
172  float eseffsixix( const reco::SuperCluster &cluster );
173  float eseffsiyiy( const reco::SuperCluster &cluster );
174 
175 // std::vector<int> flagsexcl_;
176  //std::vector<int> severitiesexcl_;
177  // const EcalSeverityLevelAlgo *sevLv;
178 
179  private:
180  void getGeometry( const edm::EventSetup &es );
181  void getTopology( const edm::EventSetup &es );
182  void getEBRecHits( const edm::Event &ev, const edm::InputTag& redEBRecHits );
183  void getEERecHits( const edm::Event &ev, const edm::InputTag& redEERecHits );
184  void getESRecHits( const edm::Event &ev, const edm::InputTag& redESRecHits );
186 
192 
193  //const EcalIntercalibConstantMap& icalMap;
198  void getIntercalibConstants( const edm::EventSetup &es );
199  void getADCToGeV ( const edm::EventSetup &es );
200  void getLaserDbService ( const edm::EventSetup &es );
201 
202 // std::vector<int> flagsexcl_;
203 // std::vector<int> severitiesexcl_;
204 
205 };
206 
207 #endif
const double Z[kNumberCalorimeter]
std::vector< float > localCovariances(const reco::BasicCluster &cluster, float w0=4.7)
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)
float e2nd(const reco::BasicCluster &cluster)
const EcalRecHitCollection * esRecHits_
std::map< DetId, EcalRecHit > rechits_map_
float eseffsiyiy(const reco::SuperCluster &cluster)
void getEBRecHits(const edm::Event &ev, const edm::InputTag &redEBRecHits)
edm::ESHandle< EcalADCToGeVConstant > agc
float eRight(const reco::BasicCluster &cluster)
float e1x5(const reco::BasicCluster &cluster)
void getEERecHits(const edm::Event &ev, const edm::InputTag &redEERecHits)
float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev)
#define X(str)
Definition: MuonsGrabber.cc:49
float e3x3(const reco::BasicCluster &cluster)
float eMax(const reco::BasicCluster &cluster)
const CaloGeometry * geometry_
float getESShape(const std::vector< float > &ESHits0)
float SuperClusterSeedTime(const reco::SuperCluster &cluster)
void getLaserDbService(const edm::EventSetup &es)
edm::ESHandle< EcalLaserDbService > laser
const EcalRecHitCollection * ebRecHits_
float eseffsirir(const reco::SuperCluster &cluster)
float matrixEnergy(const reco::BasicCluster &cluster, DetId id, int ixMin, int ixMax, int iyMin, int iyMax)
void getIntercalibConstants(const edm::EventSetup &es)
float e2x5Max(const reco::BasicCluster &cluster)
float e5x1(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)
void getADCToGeV(const edm::EventSetup &es)
double zernike42(const reco::BasicCluster &cluster, double R0=6.6, bool logW=true, float w0=4.7)
const EcalRecHitCollection * eeRecHits_
float BasicClusterSeedTime(const reco::BasicCluster &cluster)
std::vector< float > getESHits(double X, double Y, double Z, std::map< DetId, EcalRecHit > rechits_map, const CaloGeometry *geometry, CaloSubdetectorTopology *topology_p, int row=0, int plane=1)
std::vector< float > scLocalCovariances(const reco::SuperCluster &cluster, float w0=4.7)
float eseffsixix(const reco::SuperCluster &cluster)
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)
edm::ESHandle< EcalIntercalibConstants > ical
std::pair< DetId, float > getMaximum(const reco::BasicCluster &cluster)
EcalClusterLazyTools(const edm::Event &ev, const edm::EventSetup &es, const edm::InputTag &redEBRecHits, const edm::InputTag &redEERecHits, const edm::ParameterSet &config)
ESHandle< TrackerGeometry > geometry
EcalIntercalibConstantMap icalMap
float e2x5Top(const reco::BasicCluster &cluster)
void getGeometry(const edm::EventSetup &es)
float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev)
float e3x1(const reco::BasicCluster &cluster)
std::vector< DetId > matrixDetId(DetId id, int ixMin, int ixMax, int iyMin, int iyMax)
const CaloTopology * topology_
void getESRecHits(const edm::Event &ev, const edm::InputTag &redESRecHits)
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)