CMS 3D CMS Logo

EcalClusterLazyTools.h
Go to the documentation of this file.
1 #ifndef RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
2 #define RecoEcal_EgammaCoreTools_EcalClusterLazyTools_h
3 
18 
20 
28 
29 class CaloTopology;
30 class CaloGeometry;
32 
34 public:
36  const edm::EventSetup &es,
40  const edm::EventSetup &es,
44 
45  // get time of basic cluster seed crystal
46  float BasicClusterSeedTime(const reco::BasicCluster &cluster);
47  // error-weighted average of time from constituents of basic cluster
48  float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev);
49  // get BasicClusterSeedTime of the seed basic cluser of the supercluster
50  float SuperClusterSeedTime(const reco::SuperCluster &cluster);
51  // get BasicClusterTime of the seed basic cluser of the supercluster
52  float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev);
53 
54  // mapping for preshower rechits
55  std::map<DetId, EcalRecHit> rechits_map_;
56  // get Preshower hit array
57  std::vector<float> getESHits(double X,
58  double Y,
59  double Z,
60  const std::map<DetId, EcalRecHit> &rechits_map,
61  const CaloGeometry *geometry,
62  CaloSubdetectorTopology const *topology_p,
63  int row = 0,
64  int plane = 1);
65  // get Preshower hit shape
66  float getESShape(const std::vector<float> &ESHits0);
67  // get Preshower effective sigmaRR
68  float eseffsirir(const reco::SuperCluster &cluster);
69  float eseffsixix(const reco::SuperCluster &cluster);
70  float eseffsiyiy(const reco::SuperCluster &cluster);
71 
77 
78 protected:
80 
86 
88 
89  std::shared_ptr<CaloSubdetectorTopology const> ecalPS_topology_;
90 
96  void getADCToGeV(const edm::EventSetup &es);
97  void getLaserDbService(const edm::EventSetup &es);
98 
99 private:
100  void getESRecHits(const edm::Event &ev);
101 
102 }; // class EcalClusterLazyToolsBase
103 
104 template <class ClusterTools>
106 public:
108  const edm::EventSetup &es,
111  : EcalClusterLazyToolsBase(ev, es, token1, token2) {}
112 
114  const edm::EventSetup &es,
118  : EcalClusterLazyToolsBase(ev, es, token1, token2, token3) {}
119 
120  // Get the rec hit energies in a rectangle matrix around the seed.
121  std::vector<float> energyMatrix(reco::BasicCluster const &cluster, int size) const {
122  auto recHits = getEcalRecHitCollection(cluster);
123  DetId maxId = ClusterTools::getMaximum(cluster, recHits).first;
124 
125  std::vector<float> energies;
126  for (auto const &detId : CaloRectangleRange(size, maxId, *topology_)) {
127  energies.push_back(ClusterTools::recHitEnergy(detId, recHits));
128  }
129 
130  return energies;
131  }
132 
133  // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster
134  //
135  // NOTE (29/10/08): we now use an eta/phi coordinate system rather than
136  // phi/eta to minmise possible screwups, for now e5x1 isnt
137  // defined all the majority of people who call it actually
138  // want e1x5 and it is thought it is better that their code
139  // doesnt compile rather than pick up the wrong function
140  // therefore in this version and later e1x5 = e5x1 in the old
141  // version so 1x5 is 1 crystal in eta and 5 crystals in phi
142  // note e3x2 does not have a definate eta/phi geometry, it
143  // takes the maximum 3x2 block containing the seed regardless
144  // of whether that 3 in eta or phi
145  float e1x3(const reco::BasicCluster &cluster) {
146  return ClusterTools::e1x3(cluster, getEcalRecHitCollection(cluster), topology_);
147  }
148  float e3x1(const reco::BasicCluster &cluster) {
149  return ClusterTools::e3x1(cluster, getEcalRecHitCollection(cluster), topology_);
150  }
151  float e1x5(const reco::BasicCluster &cluster) {
152  return ClusterTools::e1x5(cluster, getEcalRecHitCollection(cluster), topology_);
153  }
154  float e5x1(const reco::BasicCluster &cluster) {
155  return ClusterTools::e5x1(cluster, getEcalRecHitCollection(cluster), topology_);
156  }
157  float e2x2(const reco::BasicCluster &cluster) {
158  return ClusterTools::e2x2(cluster, getEcalRecHitCollection(cluster), topology_);
159  }
160  float e3x2(const reco::BasicCluster &cluster) {
161  return ClusterTools::e3x2(cluster, getEcalRecHitCollection(cluster), topology_);
162  }
163  float e3x3(const reco::BasicCluster &cluster) {
164  return ClusterTools::e3x3(cluster, getEcalRecHitCollection(cluster), topology_);
165  }
166  float e4x4(const reco::BasicCluster &cluster) {
167  return ClusterTools::e4x4(cluster, getEcalRecHitCollection(cluster), topology_);
168  }
169  float e5x5(const reco::BasicCluster &cluster) {
170  return ClusterTools::e5x5(cluster, getEcalRecHitCollection(cluster), topology_);
171  }
172  int n5x5(const reco::BasicCluster &cluster) {
173  return ClusterTools::n5x5(cluster, getEcalRecHitCollection(cluster), topology_);
174  }
175  // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
176  // 2 crystals wide in eta, 5 wide in phi.
177  float e2x5Right(const reco::BasicCluster &cluster) {
178  return ClusterTools::e2x5Right(cluster, getEcalRecHitCollection(cluster), topology_);
179  }
180  // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
181  float e2x5Left(const reco::BasicCluster &cluster) {
182  return ClusterTools::e2x5Left(cluster, getEcalRecHitCollection(cluster), topology_);
183  }
184  // energy in the 5x2 strip above the max crystal (does not contain max crystal)
185  // 5 crystals wide in eta, 2 wide in phi.
186  float e2x5Top(const reco::BasicCluster &cluster) {
187  return ClusterTools::e2x5Top(cluster, getEcalRecHitCollection(cluster), topology_);
188  }
189 
190  // energy in the 5x2 strip below the max crystal (does not contain max crystal)
191  float e2x5Bottom(const reco::BasicCluster &cluster) {
192  return ClusterTools::e2x5Bottom(cluster, getEcalRecHitCollection(cluster), topology_);
193  }
194  // energy in a 2x5 strip containing the seed (max) crystal.
195  // 2 crystals wide in eta, 5 wide in phi.
196  // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
197  float e2x5Max(const reco::BasicCluster &cluster) {
198  return ClusterTools::e2x5Max(cluster, getEcalRecHitCollection(cluster), topology_);
199  }
200  // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
201  float eLeft(const reco::BasicCluster &cluster) {
202  return ClusterTools::eLeft(cluster, getEcalRecHitCollection(cluster), topology_);
203  }
204  float eRight(const reco::BasicCluster &cluster) {
205  return ClusterTools::eRight(cluster, getEcalRecHitCollection(cluster), topology_);
206  }
207  float eTop(const reco::BasicCluster &cluster) {
208  return ClusterTools::eTop(cluster, getEcalRecHitCollection(cluster), topology_);
209  }
210  float eBottom(const reco::BasicCluster &cluster) {
211  return ClusterTools::eBottom(cluster, getEcalRecHitCollection(cluster), topology_);
212  }
213  // the energy of the most energetic crystal in the cluster
214  float eMax(const reco::BasicCluster &cluster) {
215  return ClusterTools::eMax(cluster, getEcalRecHitCollection(cluster));
216  }
217  // the energy of the second most energetic crystal in the cluster
218  float e2nd(const reco::BasicCluster &cluster) {
219  return ClusterTools::e2nd(cluster, getEcalRecHitCollection(cluster));
220  }
221 
222  // get the DetId and the energy of the maximum energy crystal of the input cluster
223  std::pair<DetId, float> getMaximum(const reco::BasicCluster &cluster) {
224  return ClusterTools::getMaximum(cluster, getEcalRecHitCollection(cluster));
225  }
226  std::vector<float> energyBasketFractionEta(const reco::BasicCluster &cluster) {
227  return ClusterTools::energyBasketFractionEta(cluster, getEcalRecHitCollection(cluster));
228  }
229  std::vector<float> energyBasketFractionPhi(const reco::BasicCluster &cluster) {
230  return ClusterTools::energyBasketFractionPhi(cluster, getEcalRecHitCollection(cluster));
231  }
232 
233  // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
234  std::vector<float> lat(const reco::BasicCluster &cluster, bool logW = true, float w0 = 4.7) {
235  return ClusterTools::lat(cluster, getEcalRecHitCollection(cluster), geometry_, logW, w0);
236  }
237 
238  // return a vector v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
239  std::vector<float> covariances(const reco::BasicCluster &cluster, float w0 = 4.7) {
240  return ClusterTools::covariances(cluster, getEcalRecHitCollection(cluster), topology_, geometry_, w0);
241  }
242 
243  // return a vector v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
244  // this function calculates differences in eta/phi in units of crystals not
245  // global eta/phi this is gives better performance in the crack regions of
246  // the calorimeter but gives otherwise identical results to covariances
247  // function this is only defined for the barrel, it returns covariances when
248  // the cluster is in the endcap Warning: covIEtaIEta has been studied by
249  // egamma, but so far covIPhiIPhi hasnt been studied extensively so there
250  // could be a bug in the covIPhiIEta or covIPhiIPhi calculations. I dont
251  // think there is but as it hasnt been heavily used, there might be one
252  std::vector<float> localCovariances(const reco::BasicCluster &cluster, float w0 = 4.7) {
253  return ClusterTools::localCovariances(cluster, getEcalRecHitCollection(cluster), topology_, w0);
254  }
255  std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, float w0 = 4.7) {
256  return ClusterTools::scLocalCovariances(cluster, getEcalRecHitCollection(cluster), topology_, w0);
257  }
258  double zernike20(const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7) {
259  return ClusterTools::zernike20(cluster, getEcalRecHitCollection(cluster), geometry_, R0, logW, w0);
260  }
261  double zernike42(const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7) {
262  return ClusterTools::zernike42(cluster, getEcalRecHitCollection(cluster), geometry_, R0, logW, w0);
263  }
264 
265 }; // class EcalClusterLazyToolsT
266 
267 namespace noZS {
269 }
271 
272 #endif
float e2x5Top(const reco::BasicCluster &cluster)
size
Write out results.
edm::ESHandle< EcalIntercalibConstants > ical
edm::ESHandle< EcalLaserDbService > const & getLaserHandle() const
const CaloGeometry * geometry_
std::vector< float > energyBasketFractionEta(const reco::BasicCluster &cluster)
const EcalRecHitCollection * esRecHits_
int n5x5(const reco::BasicCluster &cluster)
std::vector< float > getESHits(double X, double Y, double Z, const std::map< DetId, EcalRecHit > &rechits_map, const CaloGeometry *geometry, CaloSubdetectorTopology const *topology_p, int row=0, int plane=1)
EcalClusterLazyToolsT<::EcalClusterTools > EcalClusterLazyTools
float SuperClusterSeedTime(const reco::SuperCluster &cluster)
EcalRecHitCollection const * getEcalESRecHitCollection() const
float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev)
float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev)
float e3x1(const reco::BasicCluster &cluster)
double zernike20(const reco::BasicCluster &cluster, double R0=6.6, bool logW=true, float w0=4.7)
std::vector< float > covariances(const reco::BasicCluster &cluster, float w0=4.7)
float getESShape(const std::vector< float > &ESHits0)
EcalClusterLazyToolsT(const edm::Event &ev, const edm::EventSetup &es, edm::EDGetTokenT< EcalRecHitCollection > token1, edm::EDGetTokenT< EcalRecHitCollection > token2, edm::EDGetTokenT< EcalRecHitCollection > token3)
EcalClusterLazyToolsT(const edm::Event &ev, const edm::EventSetup &es, edm::EDGetTokenT< EcalRecHitCollection > token1, edm::EDGetTokenT< EcalRecHitCollection > token2)
const EcalIntercalibConstantMap * icalMap
#define X(str)
Definition: MuonsGrabber.cc:38
void getADCToGeV(const edm::EventSetup &es)
std::pair< DetId, float > getMaximum(const reco::BasicCluster &cluster)
bool ev
const EcalRecHitCollection * ebRecHits_
float eseffsixix(const reco::SuperCluster &cluster)
void getIntercalibConstants(const edm::EventSetup &es)
float eLeft(const reco::BasicCluster &cluster)
float eRight(const reco::BasicCluster &cluster)
edm::ESHandle< EcalADCToGeVConstant > agc
float e3x2(const reco::BasicCluster &cluster)
float e1x5(const reco::BasicCluster &cluster)
std::vector< float > energyMatrix(reco::BasicCluster const &cluster, int size) const
float e5x5(const reco::BasicCluster &cluster)
std::vector< float > energyBasketFractionPhi(const reco::BasicCluster &cluster)
std::vector< float > scLocalCovariances(const reco::SuperCluster &cluster, float w0=4.7)
EcalClusterLazyToolsBase(const edm::Event &ev, const edm::EventSetup &es, edm::EDGetTokenT< EcalRecHitCollection > token1, edm::EDGetTokenT< EcalRecHitCollection > token2)
float eTop(const reco::BasicCluster &cluster)
float e2x5Max(const reco::BasicCluster &cluster)
EcalRecHitCollection const * getEcalEERecHitCollection() const
void getESRecHits(const edm::Event &ev)
edm::EDGetTokenT< EcalRecHitCollection > esRHToken_
Definition: DetId.h:17
float eseffsiyiy(const reco::SuperCluster &cluster)
std::shared_ptr< CaloSubdetectorTopology const > ecalPS_topology_
std::vector< float > lat(const reco::BasicCluster &cluster, bool logW=true, float w0=4.7)
EcalRecHitCollection const * getEcalRecHitCollection(const reco::BasicCluster &cluster) const
float e3x3(const reco::BasicCluster &cluster)
double zernike42(const reco::BasicCluster &cluster, double R0=6.6, bool logW=true, float w0=4.7)
float e4x4(const reco::BasicCluster &cluster)
float e2nd(const reco::BasicCluster &cluster)
float e5x1(const reco::BasicCluster &cluster)
float e1x3(const reco::BasicCluster &cluster)
void getLaserDbService(const edm::EventSetup &es)
const EcalRecHitCollection * eeRecHits_
float eBottom(const reco::BasicCluster &cluster)
float e2x5Left(const reco::BasicCluster &cluster)
float eMax(const reco::BasicCluster &cluster)
EcalIntercalibConstants const & getEcalIntercalibConstants() const
float e2x5Right(const reco::BasicCluster &cluster)
std::map< DetId, EcalRecHit > rechits_map_
std::vector< float > localCovariances(const reco::BasicCluster &cluster, float w0=4.7)
EcalRecHitCollection const * getEcalEBRecHitCollection() const
edm::ESHandle< EcalLaserDbService > laser
float e2x2(const reco::BasicCluster &cluster)
float BasicClusterSeedTime(const reco::BasicCluster &cluster)
const CaloTopology * topology_
float e2x5Bottom(const reco::BasicCluster &cluster)
float eseffsirir(const reco::SuperCluster &cluster)