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 
30 class CaloTopology;
31 class CaloGeometry;
33 
35  public:
38 
39  // get time of basic cluster seed crystal
40  float BasicClusterSeedTime(const reco::BasicCluster &cluster);
41  // error-weighted average of time from constituents of basic cluster
42  float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev);
43  // get BasicClusterSeedTime of the seed basic cluser of the supercluster
44  float SuperClusterSeedTime(const reco::SuperCluster &cluster);
45  // get BasicClusterTime of the seed basic cluser of the supercluster
46  float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev);
47 
48  // mapping for preshower rechits
49  std::map<DetId, EcalRecHit> rechits_map_;
50  // get Preshower hit array
51  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);
52  // get Preshower hit shape
53  float getESShape(const std::vector<float>& ESHits0);
54  // get Preshower effective sigmaRR
55  float eseffsirir( const reco::SuperCluster &cluster );
56  float eseffsixix( const reco::SuperCluster &cluster );
57  float eseffsiyiy( const reco::SuperCluster &cluster );
58 
64 
65  protected:
67 
73 
75 
76  std::shared_ptr<CaloSubdetectorTopology const> ecalPS_topology_;
77 
82  void getIntercalibConstants( const edm::EventSetup &es );
83  void getADCToGeV ( const edm::EventSetup &es );
84  void getLaserDbService ( const edm::EventSetup &es );
85 
86  private:
87  void getESRecHits( const edm::Event &ev );
88 
89 }; // class EcalClusterLazyToolsBase
90 
91 template<class ClusterTools>
93  public:
94 
96  const edm::EventSetup &es,
99  : EcalClusterLazyToolsBase(ev,es,token1,token2) {}
100 
102  const edm::EventSetup &es,
106  : EcalClusterLazyToolsBase(ev,es,token1,token2,token3) {}
107 
108  // Get the rec hit energies in a rectangle matrix around the seed.
109  std::vector<float> energyMatrix(reco::BasicCluster const& cluster, int size) const {
110 
111  auto recHits = getEcalRecHitCollection(cluster);
112  DetId maxId = ClusterTools::getMaximum(cluster, recHits).first;
113 
114  std::vector<float> energies;
115  for (auto const& detId : CaloRectangleRange(size, maxId, *topology_)) {
116  energies.push_back(ClusterTools::recHitEnergy( detId, recHits));
117  }
118 
119  return energies;
120  }
121 
122  // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster
123  //
124  // NOTE (29/10/08): we now use an eta/phi coordinate system rather than
125  // phi/eta to minmise possible screwups, for now e5x1 isnt
126  // defined all the majority of people who call it actually
127  // want e1x5 and it is thought it is better that their code
128  // doesnt compile rather than pick up the wrong function
129  // therefore in this version and later e1x5 = e5x1 in the old
130  // version so 1x5 is 1 crystal in eta and 5 crystals in phi
131  // note e3x2 does not have a definate eta/phi geometry, it
132  // takes the maximum 3x2 block containing the seed regardless
133  // of whether that 3 in eta or phi
134  float e1x3( const reco::BasicCluster &cluster )
135  { return ClusterTools::e1x3( cluster, getEcalRecHitCollection(cluster), topology_ ); }
136  float e3x1( const reco::BasicCluster &cluster )
137  { return ClusterTools::e3x1( cluster, getEcalRecHitCollection(cluster), topology_ ); }
138  float e1x5( const reco::BasicCluster &cluster )
139  { return ClusterTools::e1x5( cluster, getEcalRecHitCollection(cluster), topology_ ); }
140  float e5x1( const reco::BasicCluster &cluster )
141  { return ClusterTools::e5x1( cluster, getEcalRecHitCollection(cluster), topology_ ); }
142  float e2x2( const reco::BasicCluster &cluster )
143  { return ClusterTools::e2x2( cluster, getEcalRecHitCollection(cluster), topology_ ); }
144  float e3x2( const reco::BasicCluster &cluster )
145  { return ClusterTools::e3x2( cluster, getEcalRecHitCollection(cluster), topology_ ); }
146  float e3x3( const reco::BasicCluster &cluster )
147  { return ClusterTools::e3x3( cluster, getEcalRecHitCollection(cluster), topology_ ); }
148  float e4x4( const reco::BasicCluster &cluster )
149  { return ClusterTools::e4x4( cluster, getEcalRecHitCollection(cluster), topology_ ); }
150  float e5x5( const reco::BasicCluster &cluster )
151  { return ClusterTools::e5x5( cluster, getEcalRecHitCollection(cluster), topology_ ); }
152  int n5x5( const reco::BasicCluster &cluster )
153  { return ClusterTools::n5x5( cluster, getEcalRecHitCollection(cluster), topology_ ); }
154  // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
155  // 2 crystals wide in eta, 5 wide in phi.
156  float e2x5Right( const reco::BasicCluster &cluster )
157  { return ClusterTools::e2x5Right( cluster, getEcalRecHitCollection(cluster), topology_ ); }
158  // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
159  float e2x5Left( const reco::BasicCluster &cluster )
160  { return ClusterTools::e2x5Left( cluster, getEcalRecHitCollection(cluster), topology_ ); }
161  // energy in the 5x2 strip above the max crystal (does not contain max crystal)
162  // 5 crystals wide in eta, 2 wide in phi.
163  float e2x5Top( const reco::BasicCluster &cluster )
164  { return ClusterTools::e2x5Top( cluster, getEcalRecHitCollection(cluster), topology_ ); }
165 
166  // energy in the 5x2 strip below the max crystal (does not contain max crystal)
167  float e2x5Bottom( const reco::BasicCluster &cluster )
168  { return ClusterTools::e2x5Bottom( cluster, getEcalRecHitCollection(cluster), topology_ ); }
169  // energy in a 2x5 strip containing the seed (max) crystal.
170  // 2 crystals wide in eta, 5 wide in phi.
171  // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
172  float e2x5Max( const reco::BasicCluster &cluster )
173  { return ClusterTools::e2x5Max( cluster, getEcalRecHitCollection(cluster), topology_ ); }
174  // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
175  float eLeft( const reco::BasicCluster &cluster )
176  { return ClusterTools::eLeft( cluster, getEcalRecHitCollection(cluster), topology_ ); }
177  float eRight( const reco::BasicCluster &cluster )
178  { return ClusterTools::eRight( cluster, getEcalRecHitCollection(cluster), topology_ ); }
179  float eTop( const reco::BasicCluster &cluster )
180  { return ClusterTools::eTop( cluster, getEcalRecHitCollection(cluster), topology_ ); }
181  float eBottom( const reco::BasicCluster &cluster )
182  { return ClusterTools::eBottom( cluster, getEcalRecHitCollection(cluster), topology_ ); }
183  // the energy of the most energetic crystal in the cluster
184  float eMax( const reco::BasicCluster &cluster )
185  { return ClusterTools::eMax( cluster, getEcalRecHitCollection(cluster) ); }
186  // the energy of the second most energetic crystal in the cluster
187  float e2nd( const reco::BasicCluster &cluster )
188  { return ClusterTools::e2nd( cluster, getEcalRecHitCollection(cluster) ); }
189 
190  // get the DetId and the energy of the maximum energy crystal of the input cluster
191  std::pair<DetId, float> getMaximum( const reco::BasicCluster &cluster )
192  { return ClusterTools::getMaximum( cluster, getEcalRecHitCollection(cluster) ); }
193  std::vector<float> energyBasketFractionEta( const reco::BasicCluster &cluster )
194  { return ClusterTools::energyBasketFractionEta( cluster, getEcalRecHitCollection(cluster) ); }
195  std::vector<float> energyBasketFractionPhi( const reco::BasicCluster &cluster )
196  { return ClusterTools::energyBasketFractionPhi( cluster, getEcalRecHitCollection(cluster) ); }
197 
198  // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
199  std::vector<float> lat( const reco::BasicCluster &cluster, bool logW = true, float w0 = 4.7 )
200  { return ClusterTools::lat( cluster, getEcalRecHitCollection(cluster), geometry_, logW, w0 ); }
201 
202  // return a vector v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
203  std::vector<float> covariances(const reco::BasicCluster &cluster, float w0 = 4.7 )
204  { return ClusterTools::covariances( cluster, getEcalRecHitCollection(cluster), topology_, geometry_, w0 ); }
205 
206  // return a vector v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
207  // this function calculates differences in eta/phi in units of crystals not
208  // global eta/phi this is gives better performance in the crack regions of
209  // the calorimeter but gives otherwise identical results to covariances
210  // function this is only defined for the barrel, it returns covariances when
211  // the cluster is in the endcap Warning: covIEtaIEta has been studied by
212  // egamma, but so far covIPhiIPhi hasnt been studied extensively so there
213  // could be a bug in the covIPhiIEta or covIPhiIPhi calculations. I dont
214  // think there is but as it hasnt been heavily used, there might be one
215  std::vector<float> localCovariances(const reco::BasicCluster &cluster, float w0 = 4.7)
216  { return ClusterTools::localCovariances( cluster, getEcalRecHitCollection(cluster), topology_, w0 ); }
217  std::vector<float> scLocalCovariances(const reco::SuperCluster &cluster, float w0 = 4.7)
218  { return ClusterTools::scLocalCovariances( cluster, getEcalRecHitCollection(cluster), topology_, w0 ); }
219  double zernike20( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 )
220  { return ClusterTools::zernike20( cluster, getEcalRecHitCollection(cluster), geometry_, R0, logW, w0 ); }
221  double zernike42( const reco::BasicCluster &cluster, double R0 = 6.6, bool logW = true, float w0 = 4.7 )
222  { return ClusterTools::zernike42( cluster, getEcalRecHitCollection(cluster), geometry_, R0, logW, w0 ); }
223 
224 }; // class EcalClusterLazyToolsT
225 
226 
229 
230 
231 #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:48
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:18
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)