CMS 3D CMS Logo

EcalClusterTools.h
Go to the documentation of this file.
1 #ifndef RecoEcal_EgammaCoreTools_EcalClusterTools_h
2 #define RecoEcal_EgammaCoreTools_EcalClusterTools_h
3 
34 //#include "DataFormats/Math/interface/Point3D.h"
36 //includes for ShowerShape function to work
37 #include <vector>
38 #include <cmath>
39 #include <TMath.h>
40 #include <TMatrixT.h>
41 #include <TMatrixD.h>
42 #include <TVectorT.h>
43 #include <TVectorD.h>
44 
49 
51 
52 #include "CLHEP/Geometry/Transform3D.h"
53 
61 
62 class DetId;
63 class CaloTopology;
64 class CaloGeometry;
65 
67  // major and minor cluster moments wrt principale axes:
68  float sMaj;
69  float sMin;
70  // angle between sMaj and phi:
71  float alpha;
72  Cluster2ndMoments() : sMaj(0.), sMin(0.), alpha(0.) {}
73 };
74 
75 template <bool noZS>
77 public:
78  // various energies in the matrix nxn surrounding the maximum energy crystal of the input cluster
79  //we use an eta/phi coordinate system rather than phi/eta
80  //note e3x2 does not have a definate eta/phi geometry, it takes the maximum 3x2 block containing the
81  //seed regardless of whether that 3 in eta or phi
82  static float e1x3(const reco::BasicCluster &cluster,
84  const CaloTopology *topology);
85 
86  static float e3x1(const reco::BasicCluster &cluster,
88  const CaloTopology *topology);
89 
90  static float e1x5(const reco::BasicCluster &cluster,
92  const CaloTopology *topology);
93 
94  static float e5x1(const reco::BasicCluster &cluster,
96  const CaloTopology *topology);
97 
98  static float e2x2(const reco::BasicCluster &cluster,
100  const CaloTopology *topology);
101 
102  static float e3x2(const reco::BasicCluster &cluster,
104  const CaloTopology *topology);
105 
106  static float e3x3(const reco::BasicCluster &cluster,
108  const CaloTopology *topology);
109 
110  static float e4x4(const reco::BasicCluster &cluster,
112  const CaloTopology *topology);
113 
114  static float e5x5(const reco::BasicCluster &cluster,
116  const CaloTopology *topology);
117  static int n5x5(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology);
118 
119  // energy in the 2x5 strip right of the max crystal (does not contain max crystal)
120  // 2 crystals wide in eta, 5 wide in phi.
121  static float e2x5Right(const reco::BasicCluster &cluster,
123  const CaloTopology *topology);
124  // energy in the 2x5 strip left of the max crystal (does not contain max crystal)
125 
126  static float e2x5Left(const reco::BasicCluster &cluster,
128  const CaloTopology *topology);
129  // energy in the 5x2 strip above the max crystal (does not contain max crystal)
130  // 5 crystals wide in eta, 2 wide in phi.
131 
132  static float e2x5Top(const reco::BasicCluster &cluster,
134  const CaloTopology *topology);
135  // energy in the 5x2 strip below the max crystal (does not contain max crystal)
136 
137  static float e2x5Bottom(const reco::BasicCluster &cluster,
139  const CaloTopology *topology);
140  // energy in a 2x5 strip containing the seed (max) crystal.
141  // 2 crystals wide in eta, 5 wide in phi.
142  // it is the maximum of either (1x5left + 1x5center) or (1x5right + 1x5center)
143  static float e2x5Max(const reco::BasicCluster &cluster,
145  const CaloTopology *topology);
146 
147  // energies in the crystal left, right, top, bottom w.r.t. to the most energetic crystal
148  static float eLeft(const reco::BasicCluster &cluster,
150  const CaloTopology *topology);
151 
152  static float eRight(const reco::BasicCluster &cluster,
154  const CaloTopology *topology);
155 
156  static float eTop(const reco::BasicCluster &cluster,
158  const CaloTopology *topology);
159 
160  static float eBottom(const reco::BasicCluster &cluster,
162  const CaloTopology *topology);
163  // the energy of the most energetic crystal in the cluster
164 
165  static float eMax(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits);
166 
167  // the energy of the second most energetic crystal in the cluster
168  static float e2nd(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits);
169 
170  // get the DetId and the energy of the maximum energy crystal of the input cluster
171  static std::pair<DetId, float> getMaximum(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits);
172 
173  static std::vector<float> energyBasketFractionEta(const reco::BasicCluster &cluster,
175 
176  static std::vector<float> energyBasketFractionPhi(const reco::BasicCluster &cluster,
178 
179  // return a vector v with v[0] = etaLat, v[1] = phiLat, v[2] = lat
180  static std::vector<float> lat(const reco::BasicCluster &cluster,
182  const CaloGeometry *geometry,
183  bool logW = true,
184  float w0 = 4.7);
185 
186  // return an array v with v[0] = covEtaEta, v[1] = covEtaPhi, v[2] = covPhiPhi
187  static std::array<float, 3> covariances(const reco::BasicCluster &cluster,
189  const CaloTopology *topology,
190  const CaloGeometry *geometry,
191  float w0 = 4.7);
192 
193  // return an array v with v[0] = covIEtaIEta, v[1] = covIEtaIPhi, v[2] = covIPhiIPhi
194  //this function calculates differences in eta/phi in units of crystals not global eta/phi
195  //this is gives better performance in the crack regions of the calorimeter but gives otherwise identical results to covariances function
196  // except that it doesnt need an eta based correction funtion in the endcap
197  //it is multipled by an approprate crystal size to ensure it gives similar values to covariances(...)
198  //
199  //Warning: covIEtaIEta has been studied by egamma, but so far covIPhiIPhi hasnt been studied extensively so there could be a bug in
200  // the covIPhiIEta or covIPhiIPhi calculations. I dont think there is but as it hasnt been heavily used, there might be one
201  static std::array<float, 3> localCovariances(const reco::BasicCluster &cluster,
203  const CaloTopology *topology,
205  const EcalPFRecHitThresholds *thresholds = nullptr,
206  float multEB = 0.0,
207  float multEE = 0.0);
208 
209  static std::array<float, 3> scLocalCovariances(const reco::SuperCluster &cluster,
211  const CaloTopology *topology,
212  float w0 = 4.7);
213 
214  // cluster second moments with respect to principal axes:
215  static Cluster2ndMoments cluster2ndMoments(const reco::BasicCluster &basicCluster,
217  double phiCorrectionFactor = 0.8,
218  double w0 = 4.7,
219  bool useLogWeights = true);
220 
221  static Cluster2ndMoments cluster2ndMoments(const reco::SuperCluster &superCluster,
223  double phiCorrectionFactor = 0.8,
224  double w0 = 4.7,
225  bool useLogWeights = true);
226  static Cluster2ndMoments cluster2ndMoments(const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs,
227  double phiCorrectionFactor = 0.8,
228  double w0 = 4.7,
229  bool useLogWeights = true);
230 
231  static double zernike20(const reco::BasicCluster &cluster,
233  const CaloGeometry *geometry,
234  double R0 = 6.6,
235  bool logW = true,
236  float w0 = 4.7);
237  static double zernike42(const reco::BasicCluster &cluster,
239  const CaloGeometry *geometry,
240  double R0 = 6.6,
241  bool logW = true,
242  float w0 = 4.7);
243 
244  // get the detId's of a matrix centered in the maximum energy crystal = (0,0)
245  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
246  static std::vector<DetId> matrixDetId(const CaloTopology *topology, DetId id, CaloRectangle rectangle);
247  static std::vector<DetId> matrixDetId(const CaloTopology *topology, DetId id, int size) {
248  return matrixDetId(topology, id, CaloRectangle{-size, size, -size, size});
249  }
250 
251  // get the energy deposited in a matrix centered in the maximum energy crystal = (0,0)
252  // the size is specified by ixMin, ixMax, iyMin, iyMax in unit of crystals
253  static float matrixEnergy(const reco::BasicCluster &cluster,
255  const CaloTopology *topology,
256  DetId id,
257  CaloRectangle rectangle);
258  static int matrixSize(const reco::BasicCluster &cluster,
260  const CaloTopology *topology,
261  DetId id,
262  CaloRectangle rectangle);
263 
264  static float getFraction(const std::vector<std::pair<DetId, float>> &v_id, DetId id);
265  // get the DetId and the energy of the maximum energy crystal in a vector of DetId
266  static std::pair<DetId, float> getMaximum(const std::vector<std::pair<DetId, float>> &v_id,
268 
269  // get the energy of a DetId, return 0 if the DetId is not in the collection
270  static float recHitEnergy(DetId id, const EcalRecHitCollection *recHits);
271 
272  //Shower shape variables return vector <Roundness, Angle> of a photon
273  static std::vector<float> roundnessBarrelSuperClusters(const reco::SuperCluster &superCluster,
275  int weightedPositionMethod = 0,
276  float energyThreshold = 0.0);
277  static std::vector<float> roundnessBarrelSuperClustersUserExtended(const reco::SuperCluster &superCluster,
279  int ieta_delta = 0,
280  int iphi_delta = 0,
281  float energyRHThresh = 0.00000,
282  int weightedPositionMethod = 0);
283  static std::vector<float> roundnessSelectedBarrelRecHits(
284  const std::vector<std::pair<const EcalRecHit *, float>> &rhVector, int weightedPositionMethod = 0);
285 
286  //works out the number of staturated crystals in 5x5
287  static int nrSaturatedCrysIn5x5(const DetId &id, const EcalRecHitCollection *recHits, const CaloTopology *topology);
288 
289 private:
292  double r;
293  double phi;
294  };
295 
296  static std::vector<EcalClusterEnergyDeposition> getEnergyDepTopology(const reco::BasicCluster &cluster,
298  const CaloGeometry *geometry,
299  bool logW,
300  float w0);
301 
304  const CaloTopology *topology,
305  const CaloGeometry *geometry);
306 
307  //return energy weighted mean distance from the seed crystal in number of crystals
308  //<iEta,iPhi>, iPhi is not defined for endcap and is returned as zero
309  static std::pair<float, float> mean5x5PositionInLocalCrysCoord(const reco::BasicCluster &cluster,
311  const CaloTopology *topology);
312 
313  static std::pair<float, float> mean5x5PositionInXY(const reco::BasicCluster &cluster,
315  const CaloTopology *topology);
316 
317  static double f00(double r) { return 1; }
318  static double f11(double r) { return r; }
319  static double f20(double r) { return 2.0 * r * r - 1.0; }
320  static double f22(double r) { return r * r; }
321  static double f31(double r) { return 3.0 * r * r * r - 2.0 * r; }
322  static double f33(double r) { return r * r * r; }
323  static double f40(double r) { return 6.0 * r * r * r * r - 6.0 * r * r + 1.0; }
324  static double f42(double r) { return 4.0 * r * r * r * r - 3.0 * r * r; }
325  static double f44(double r) { return r * r * r * r; }
326  static double f51(double r) { return 10.0 * pow(r, 5) - 12.0 * pow(r, 3) + 3.0 * r; }
327  static double f53(double r) { return 5.0 * pow(r, 5) - 4.0 * pow(r, 3); }
328  static double f55(double r) { return pow(r, 5); }
329 
330  static double absZernikeMoment(const reco::BasicCluster &cluster,
332  const CaloGeometry *geometry,
333  int n,
334  int m,
335  double R0,
336  bool logW,
337  float w0);
338  static double fast_AbsZernikeMoment(const reco::BasicCluster &cluster,
340  const CaloGeometry *geometry,
341  int n,
342  int m,
343  double R0,
344  bool logW,
345  float w0);
346  static double calc_AbsZernikeMoment(const reco::BasicCluster &cluster,
348  const CaloGeometry *geometry,
349  int n,
350  int m,
351  double R0,
352  bool logW,
353  float w0);
354 
355  static double factorial(int n) {
356  double res = 1.;
357  for (int i = 2; i <= n; ++i)
358  res *= i;
359  return res;
360  }
361 
362  //useful functions for the localCovariances function
363  static float getIEta(const DetId &id);
364  static float getIPhi(const DetId &id);
365  static float getNormedIX(const DetId &id);
366  static float getNormedIY(const DetId &id);
367  static float getDPhiEndcap(const DetId &crysId, float meanX, float meanY);
368  static float getNrCrysDiffInEta(const DetId &crysId, const DetId &orginId);
369  static float getNrCrysDiffInPhi(const DetId &crysId, const DetId &orginId);
370 
371  //useful functions for showerRoundnessBarrel function
372  static int deltaIEta(int seed_ieta, int rh_ieta);
373  static int deltaIPhi(int seed_iphi, int rh_iphi);
374  static std::vector<int> getSeedPosition(const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs);
375  static float getSumEnergy(const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs);
376  static float computeWeight(float eRH, float energyTotal, int weightedPositionMethod);
377 };
378 
379 // implementation
380 template <bool noZS>
381 float EcalClusterToolsT<noZS>::getFraction(const std::vector<std::pair<DetId, float>> &v_id, DetId id) {
382  if (noZS)
383  return 1.0;
384  float frac = 0.0;
385  for (size_t i = 0; i < v_id.size(); ++i) {
386  if (v_id[i].first.rawId() == id.rawId()) {
387  frac = v_id[i].second;
388  break;
389  }
390  }
391  return frac;
392 }
393 
394 template <bool noZS>
395 std::pair<DetId, float> EcalClusterToolsT<noZS>::getMaximum(const std::vector<std::pair<DetId, float>> &v_id,
396  const EcalRecHitCollection *recHits) {
397  float max = 0;
398  DetId id(0);
399  for (size_t i = 0; i < v_id.size(); ++i) {
400  float energy = recHitEnergy(v_id[i].first, recHits) * (noZS ? 1.0 : v_id[i].second);
401  if (energy > max) {
402  max = energy;
403  id = v_id[i].first;
404  }
405  }
406  return std::pair<DetId, float>(id, max);
407 }
408 
409 template <bool noZS>
410 std::pair<DetId, float> EcalClusterToolsT<noZS>::getMaximum(const reco::BasicCluster &cluster,
411  const EcalRecHitCollection *recHits) {
412  return getMaximum(cluster.hitsAndFractions(), recHits);
413 }
414 
415 template <bool noZS>
417  if (id == DetId(0)) {
418  return 0;
419  } else {
421  if (it != recHits->end()) {
422  if (noZS && (it->checkFlag(EcalRecHit::kTowerRecovered) || it->checkFlag(EcalRecHit::kWeird) ||
423  (it->detid().subdetId() == EcalBarrel && it->checkFlag(EcalRecHit::kDiWeird))))
424  return 0.0;
425  else
426  return it->energy();
427  } else {
428  //throw cms::Exception("EcalRecHitNotFound") << "The recHit corresponding to the DetId" << id.rawId() << " not found in the EcalRecHitCollection";
429  // the recHit is not in the collection (hopefully zero suppressed)
430  return 0;
431  }
432  }
433  return 0;
434 }
435 
436 // Returns the energy in a rectangle of crystals
437 // specified in eta by ixMin and ixMax
438 // and in phi by iyMin and iyMax
439 //
440 // Reference picture (X=seed crystal)
441 // iy ___________
442 // 2 |_|_|_|_|_|
443 // 1 |_|_|_|_|_|
444 // 0 |_|_|X|_|_|
445 // -1 |_|_|_|_|_|
446 // -2 |_|_|_|_|_|
447 // -2 -1 0 1 2 ix
448 template <bool noZS>
451  const CaloTopology *topology,
452  DetId id,
453  CaloRectangle rectangle) {
454  float energy = 0;
455  auto const &vId = cluster.hitsAndFractions();
456 
457  for (auto const &detId : rectangle(id, *topology)) {
458  energy += recHitEnergy(detId, recHits) * getFraction(vId, detId);
459  }
460 
461  return energy;
462 }
463 
464 template <bool noZS>
467  const CaloTopology *topology,
468  DetId id,
469  CaloRectangle rectangle) {
470  int result = 0;
471 
472  for (auto const &detId : rectangle(id, *topology)) {
473  const float energy = recHitEnergy(detId, recHits);
474  const float frac = getFraction(cluster.hitsAndFractions(), detId);
475  if (energy * frac > 0)
476  result++;
477  }
478  return result;
479 }
480 
481 template <bool noZS>
482 std::vector<DetId> EcalClusterToolsT<noZS>::matrixDetId(const CaloTopology *topology,
483  DetId id,
484  CaloRectangle rectangle) {
485  std::vector<DetId> v;
486  for (auto const &detId : rectangle(id, *topology)) {
487  if (detId != DetId(0))
488  v.push_back(detId);
489  }
490  return v;
491 }
492 
493 template <bool noZS>
496  const CaloTopology *topology) {
497  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
498  std::list<float> energies;
499  float max_E = matrixEnergy(cluster, recHits, topology, id, {-1, 0, -1, 0});
500  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-1, 0, 0, 1}));
501  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {0, 1, 0, 1}));
502  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {0, 1, -1, 0}));
503  return max_E;
504 }
505 
506 template <bool noZS>
509  const CaloTopology *topology) {
510  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
511  float max_E = matrixEnergy(cluster, recHits, topology, id, {-1, 1, -1, 0});
512  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {0, 1, -1, 1}));
513  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-1, 1, 0, 1}));
514  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-1, 0, -1, 1}));
515  return max_E;
516 }
517 
518 template <bool noZS>
521  const CaloTopology *topology) {
522  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
523  return matrixEnergy(cluster, recHits, topology, id, {-1, 1, -1, 1});
524 }
525 
526 template <bool noZS>
529  const CaloTopology *topology) {
530  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
531  float max_E = matrixEnergy(cluster, recHits, topology, id, {-1, 2, -2, 1});
532  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-2, 1, -2, 1}));
533  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-2, 1, -1, 2}));
534  max_E = std::max(max_E, matrixEnergy(cluster, recHits, topology, id, {-1, 2, -1, 2}));
535  return max_E;
536 }
537 
538 template <bool noZS>
541  const CaloTopology *topology) {
542  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
543  return matrixEnergy(cluster, recHits, topology, id, {-2, 2, -2, 2});
544 }
545 
546 template <bool noZS>
549  const CaloTopology *topology) {
550  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
551  return matrixSize(cluster, recHits, topology, id, {-2, 2, -2, 2});
552 }
553 
554 template <bool noZS>
556  return getMaximum(cluster.hitsAndFractions(), recHits).second;
557 }
558 
559 template <bool noZS>
561  std::vector<float> energies;
562  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
563  energies.reserve(v_id.size());
564  if (v_id.size() < 2)
565  return 0;
566  for (size_t i = 0; i < v_id.size(); ++i) {
567  energies.push_back(recHitEnergy(v_id[i].first, recHits) * (noZS ? 1.0 : v_id[i].second));
568  }
569  std::partial_sort(energies.begin(), energies.begin() + 2, energies.end(), std::greater<float>());
570  return energies[1];
571 }
572 
573 template <bool noZS>
576  const CaloTopology *topology) {
577  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
578  return matrixEnergy(cluster, recHits, topology, id, {1, 2, -2, 2});
579 }
580 
581 template <bool noZS>
584  const CaloTopology *topology) {
585  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
586  return matrixEnergy(cluster, recHits, topology, id, {-2, -1, -2, 2});
587 }
588 
589 template <bool noZS>
592  const CaloTopology *topology) {
593  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
594  return matrixEnergy(cluster, recHits, topology, id, {-2, 2, 1, 2});
595 }
596 
597 template <bool noZS>
600  const CaloTopology *topology) {
601  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
602  return matrixEnergy(cluster, recHits, topology, id, {-2, 2, -2, -1});
603 }
604 
605 // Energy in 2x5 strip containing the max crystal.
606 // Adapted from code by Sam Harper
607 template <bool noZS>
610  const CaloTopology *topology) {
611  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
612 
613  // 1x5 strip left of seed
614  float left = matrixEnergy(cluster, recHits, topology, id, {-1, -1, -2, 2});
615  // 1x5 strip right of seed
616  float right = matrixEnergy(cluster, recHits, topology, id, {1, 1, -2, 2});
617  // 1x5 strip containing seed
618  float centre = matrixEnergy(cluster, recHits, topology, id, {0, 0, -2, 2});
619 
620  // Return the maximum of (left+center) or (right+center) strip
621  return left > right ? left + centre : right + centre;
622 }
623 
624 template <bool noZS>
627  const CaloTopology *topology) {
628  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
629  return matrixEnergy(cluster, recHits, topology, id, {0, 0, -2, 2});
630 }
631 
632 template <bool noZS>
635  const CaloTopology *topology) {
636  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
637  return matrixEnergy(cluster, recHits, topology, id, {-2, 2, 0, 0});
638 }
639 
640 template <bool noZS>
643  const CaloTopology *topology) {
644  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
645  return matrixEnergy(cluster, recHits, topology, id, {0, 0, -1, 1});
646 }
647 
648 template <bool noZS>
651  const CaloTopology *topology) {
652  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
653  return matrixEnergy(cluster, recHits, topology, id, {-1, 1, 0, 0});
654 }
655 
656 template <bool noZS>
659  const CaloTopology *topology) {
660  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
661  return matrixEnergy(cluster, recHits, topology, id, {-1, -1, 0, 0});
662 }
663 
664 template <bool noZS>
667  const CaloTopology *topology) {
668  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
669  return matrixEnergy(cluster, recHits, topology, id, {1, 1, 0, 0});
670 }
671 
672 template <bool noZS>
675  const CaloTopology *topology) {
676  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
677  return matrixEnergy(cluster, recHits, topology, id, {0, 0, 1, 1});
678 }
679 
680 template <bool noZS>
683  const CaloTopology *topology) {
684  DetId id = getMaximum(cluster.hitsAndFractions(), recHits).first;
685  return matrixEnergy(cluster, recHits, topology, id, {0, 0, -1, -1});
686 }
687 
688 template <bool noZS>
690  const EcalRecHitCollection *recHits) {
691  std::vector<float> basketFraction(2 * EBDetId::kModulesPerSM);
692  float clusterEnergy = cluster.energy();
693  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
694  if (v_id[0].first.subdetId() != EcalBarrel) {
695  edm::LogWarning("EcalClusterToolsT<noZS>::energyBasketFractionEta")
696  << "Trying to get basket fraction for endcap basic-clusters. Basket fractions can be obtained ONLY for barrel "
697  "basic-clusters. Returning empty vector.";
698  return basketFraction;
699  }
700  for (size_t i = 0; i < v_id.size(); ++i) {
701  basketFraction[EBDetId(v_id[i].first).im() - 1 + EBDetId(v_id[i].first).positiveZ() * EBDetId::kModulesPerSM] +=
702  recHitEnergy(v_id[i].first, recHits) * v_id[i].second / clusterEnergy;
703  }
704  std::sort(basketFraction.rbegin(), basketFraction.rend());
705  return basketFraction;
706 }
707 
708 template <bool noZS>
710  const EcalRecHitCollection *recHits) {
711  std::vector<float> basketFraction(2 * (EBDetId::MAX_IPHI / EBDetId::kCrystalsInPhi));
712  float clusterEnergy = cluster.energy();
713  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
714  if (v_id[0].first.subdetId() != EcalBarrel) {
715  edm::LogWarning("EcalClusterToolsT<noZS>::energyBasketFractionPhi")
716  << "Trying to get basket fraction for endcap basic-clusters. Basket fractions can be obtained ONLY for barrel "
717  "basic-clusters. Returning empty vector.";
718  return basketFraction;
719  }
720  for (size_t i = 0; i < v_id.size(); ++i) {
721  basketFraction[(EBDetId(v_id[i].first).iphi() - 1) / EBDetId::kCrystalsInPhi +
723  recHitEnergy(v_id[i].first, recHits) * (noZS ? 1.0 : v_id[i].second) / clusterEnergy;
724  }
725  std::sort(basketFraction.rbegin(), basketFraction.rend());
726  return basketFraction;
727 }
728 
729 template <bool noZS>
730 std::vector<typename EcalClusterToolsT<noZS>::EcalClusterEnergyDeposition>
733  const CaloGeometry *geometry,
734  bool logW,
735  float w0) {
736  std::vector<typename EcalClusterToolsT<noZS>::EcalClusterEnergyDeposition> energyDistribution;
737  // init a map of the energy deposition centered on the
738  // cluster centroid. This is for momenta calculation only.
739  CLHEP::Hep3Vector clVect(cluster.position().x(), cluster.position().y(), cluster.position().z());
740  CLHEP::Hep3Vector clDir(clVect);
741  clDir *= 1.0 / clDir.mag();
742  // in the transverse plane, axis perpendicular to clusterDir
743  CLHEP::Hep3Vector theta_axis(clDir.y(), -clDir.x(), 0.0);
744  theta_axis *= 1.0 / theta_axis.mag();
745  CLHEP::Hep3Vector phi_axis = theta_axis.cross(clDir);
746 
747  const std::vector<std::pair<DetId, float>> &clusterDetIds = cluster.hitsAndFractions();
748 
750  EcalRecHit testEcalRecHit;
751  std::vector<std::pair<DetId, float>>::const_iterator posCurrent;
752  // loop over crystals
753  for (posCurrent = clusterDetIds.begin(); posCurrent != clusterDetIds.end(); ++posCurrent) {
754  EcalRecHitCollection::const_iterator itt = recHits->find((*posCurrent).first);
755  testEcalRecHit = *itt;
756 
757  if (((*posCurrent).first != DetId(0)) && (recHits->find((*posCurrent).first) != recHits->end())) {
758  clEdep.deposited_energy = testEcalRecHit.energy() * (noZS ? 1.0 : (*posCurrent).second);
759  // if logarithmic weight is requested, apply cut on minimum energy of the recHit
760  if (logW) {
761  //double w0 = parameterMap_.find("W0")->second;
762 
763  double weight = std::max(0.0, w0 + log(std::abs(clEdep.deposited_energy) / cluster.energy()));
764  if (weight == 0) {
765  LogDebug("ClusterShapeAlgo") << "Crystal has insufficient energy: E = " << clEdep.deposited_energy
766  << " GeV; skipping... ";
767  continue;
768  } else
769  LogDebug("ClusterShapeAlgo") << "===> got crystal. Energy = " << clEdep.deposited_energy << " GeV. ";
770  }
771  DetId id_ = (*posCurrent).first;
772  const CaloSubdetectorGeometry *geo = geometry->getSubdetectorGeometry(id_);
773  auto this_cell = geo->getGeometry(id_);
774  const GlobalPoint &cellPos = this_cell->getPosition();
775  CLHEP::Hep3Vector gblPos(cellPos.x(), cellPos.y(), cellPos.z()); //surface position?
776  // Evaluate the distance from the cluster centroid
777  CLHEP::Hep3Vector diff = gblPos - clVect;
778  // Important: for the moment calculation, only the "lateral distance" is important
779  // "lateral distance" r_i = distance of the digi position from the axis Origin-Cluster Center
780  // ---> subtract the projection on clDir
781  CLHEP::Hep3Vector DigiVect = diff - diff.dot(clDir) * clDir;
782  clEdep.r = DigiVect.mag();
783  LogDebug("ClusterShapeAlgo") << "E = " << clEdep.deposited_energy << "\tdiff = " << diff.mag()
784  << "\tr = " << clEdep.r;
785  clEdep.phi = DigiVect.angle(theta_axis);
786  if (DigiVect.dot(phi_axis) < 0)
787  clEdep.phi = 2 * M_PI - clEdep.phi;
788  energyDistribution.push_back(clEdep);
789  }
790  }
791  return energyDistribution;
792 }
793 
794 template <bool noZS>
795 std::vector<float> EcalClusterToolsT<noZS>::lat(const reco::BasicCluster &cluster,
797  const CaloGeometry *geometry,
798  bool logW,
799  float w0) {
800  std::vector<EcalClusterToolsT::EcalClusterEnergyDeposition> energyDistribution =
801  getEnergyDepTopology(cluster, recHits, geometry, logW, w0);
802 
803  std::vector<float> lat;
804  double r, redmoment = 0;
805  double phiRedmoment = 0;
806  double etaRedmoment = 0;
807  int n, n1, n2, tmp;
808  int clusterSize = energyDistribution.size();
809  float etaLat_, phiLat_, lat_;
810  if (clusterSize < 3) {
811  etaLat_ = 0.0;
812  lat_ = 0.0;
813  lat.push_back(0.);
814  lat.push_back(0.);
815  lat.push_back(0.);
816  return lat;
817  }
818 
819  n1 = 0;
820  n2 = 1;
821  if (energyDistribution[1].deposited_energy > energyDistribution[0].deposited_energy) {
822  tmp = n2;
823  n2 = n1;
824  n1 = tmp;
825  }
826  for (int i = 2; i < clusterSize; i++) {
827  n = i;
828  if (energyDistribution[i].deposited_energy > energyDistribution[n1].deposited_energy) {
829  tmp = n2;
830  n2 = n1;
831  n1 = i;
832  n = tmp;
833  } else {
834  if (energyDistribution[i].deposited_energy > energyDistribution[n2].deposited_energy) {
835  tmp = n2;
836  n2 = i;
837  n = tmp;
838  }
839  }
840 
841  r = energyDistribution[n].r;
842  redmoment += r * r * energyDistribution[n].deposited_energy;
843  double rphi = r * cos(energyDistribution[n].phi);
844  phiRedmoment += rphi * rphi * energyDistribution[n].deposited_energy;
845  double reta = r * sin(energyDistribution[n].phi);
846  etaRedmoment += reta * reta * energyDistribution[n].deposited_energy;
847  }
848  double e1 = energyDistribution[n1].deposited_energy;
849  double e2 = energyDistribution[n2].deposited_energy;
850 
851  lat_ = redmoment / (redmoment + 2.19 * 2.19 * (e1 + e2));
852  phiLat_ = phiRedmoment / (phiRedmoment + 2.19 * 2.19 * (e1 + e2));
853  etaLat_ = etaRedmoment / (etaRedmoment + 2.19 * 2.19 * (e1 + e2));
854 
855  lat.push_back(etaLat_);
856  lat.push_back(phiLat_);
857  lat.push_back(lat_);
858  return lat;
859 }
860 
861 template <bool noZS>
864  const CaloTopology *topology,
865  const CaloGeometry *geometry) {
866  // find mean energy position of a 5x5 cluster around the maximum
867  math::XYZVector meanPosition(0.0, 0.0, 0.0);
868  const std::vector<std::pair<DetId, float>> &hsAndFs = cluster.hitsAndFractions();
869  std::vector<DetId> v_id = matrixDetId(topology, getMaximum(cluster, recHits).first, 2);
870  for (const std::pair<DetId, float> &hitAndFrac : hsAndFs) {
871  for (std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it) {
872  if (hitAndFrac.first != *it && !noZS)
873  continue;
874  const CaloSubdetectorGeometry *geo = geometry->getSubdetectorGeometry(*it);
875  GlobalPoint positionGP = geo->getGeometry(*it)->getPosition();
876  math::XYZVector position(positionGP.x(), positionGP.y(), positionGP.z());
877  meanPosition = meanPosition + recHitEnergy(*it, recHits) * position * hitAndFrac.second;
878  }
879  if (noZS)
880  break;
881  }
882  return meanPosition / e5x5(cluster, recHits, topology);
883 }
884 
885 //returns mean energy weighted eta/phi in crystals from the seed
886 //iPhi is not defined for endcap and is returned as zero
887 //return <eta,phi>
888 //we have an issue in working out what to do for negative energies
889 //I (Sam Harper) think it makes sense to ignore crystals with E<0 in the calculation as they are ignored
890 //in the sigmaIEtaIEta calculation (well they arent fully ignored, they do still contribute to the e5x5 sum
891 //in the sigmaIEtaIEta calculation but not here)
892 template <bool noZS>
895  const CaloTopology *topology) {
896  DetId seedId = getMaximum(cluster, recHits).first;
897  float meanDEta = 0.;
898  float meanDPhi = 0.;
899  float energySum = 0.;
900 
901  const std::vector<std::pair<DetId, float>> &hsAndFs = cluster.hitsAndFractions();
902  std::vector<DetId> v_id = matrixDetId(topology, seedId, 2);
903  for (const std::pair<DetId, float> &hAndF : hsAndFs) {
904  for (std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it) {
905  if (hAndF.first != *it && !noZS)
906  continue;
907  float energy = recHitEnergy(*it, recHits) * hAndF.second;
908  if (energy < 0.)
909  continue; //skipping negative energy crystals
910  meanDEta += energy * getNrCrysDiffInEta(*it, seedId);
911  meanDPhi += energy * getNrCrysDiffInPhi(*it, seedId);
912  energySum += energy;
913  }
914  if (noZS)
915  break;
916  }
917  meanDEta /= energySum;
918  meanDPhi /= energySum;
919  return std::pair<float, float>(meanDEta, meanDPhi);
920 }
921 
922 //returns mean energy weighted x/y in normalised crystal coordinates
923 //only valid for endcap, returns 0,0 for barrel
924 //we have an issue in working out what to do for negative energies
925 //I (Sam Harper) think it makes sense to ignore crystals with E<0 in the calculation as they are ignored
926 //in the sigmaIEtaIEta calculation (well they arent fully ignored, they do still contribute to the e5x5 sum
927 //in the sigmaIEtaIEta calculation but not here)
928 template <bool noZS>
931  const CaloTopology *topology) {
932  DetId seedId = getMaximum(cluster, recHits).first;
933 
934  std::pair<float, float> meanXY(0., 0.);
935  if (seedId.subdetId() == EcalBarrel)
936  return meanXY;
937 
938  float energySum = 0.;
939 
940  const std::vector<std::pair<DetId, float>> &hsAndFs = cluster.hitsAndFractions();
941  std::vector<DetId> v_id = matrixDetId(topology, seedId, 2);
942  for (const std::pair<DetId, float> &hAndF : hsAndFs) {
943  for (std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it) {
944  if (hAndF.first != *it && !noZS)
945  continue;
946  float energy = recHitEnergy(*it, recHits) * hAndF.second;
947  if (energy < 0.)
948  continue; //skipping negative energy crystals
949  meanXY.first += energy * getNormedIX(*it);
950  meanXY.second += energy * getNormedIY(*it);
951  energySum += energy;
952  }
953  if (noZS)
954  break;
955  }
956  meanXY.first /= energySum;
957  meanXY.second /= energySum;
958  return meanXY;
959 }
960 
961 template <bool noZS>
962 std::array<float, 3> EcalClusterToolsT<noZS>::covariances(const reco::BasicCluster &cluster,
964  const CaloTopology *topology,
965  const CaloGeometry *geometry,
966  float w0) {
967  float e_5x5 = e5x5(cluster, recHits, topology);
968  float covEtaEta, covEtaPhi, covPhiPhi;
969  if (e_5x5 >= 0.) {
970  //double w0_ = parameterMap_.find("W0")->second;
971  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
972  math::XYZVector meanPosition = meanClusterPosition(cluster, recHits, topology, geometry);
973 
974  // now we can calculate the covariances
975  double numeratorEtaEta = 0;
976  double numeratorEtaPhi = 0;
977  double numeratorPhiPhi = 0;
978  double denominator = 0;
979 
980  DetId id = getMaximum(v_id, recHits).first;
981  CaloRectangle rectangle{-2, 2, -2, 2};
982  for (auto const &detId : rectangle(id, *topology)) {
983  float frac = getFraction(v_id, detId);
984  float energy = recHitEnergy(detId, recHits) * frac;
985 
986  if (energy <= 0)
987  continue;
988 
989  const CaloSubdetectorGeometry *geo = geometry->getSubdetectorGeometry(detId);
990  GlobalPoint position = geo->getGeometry(detId)->getPosition();
991 
992  double dPhi = position.phi() - meanPosition.phi();
993  if (dPhi > +Geom::pi()) {
994  dPhi = Geom::twoPi() - dPhi;
995  }
996  if (dPhi < -Geom::pi()) {
997  dPhi = Geom::twoPi() + dPhi;
998  }
999 
1000  double dEta = position.eta() - meanPosition.eta();
1001  double w = 0.;
1002  w = std::max(0.0f, w0 + std::log(energy / e_5x5));
1003 
1004  denominator += w;
1005  numeratorEtaEta += w * dEta * dEta;
1006  numeratorEtaPhi += w * dEta * dPhi;
1007  numeratorPhiPhi += w * dPhi * dPhi;
1008  }
1009 
1010  if (denominator != 0.0) {
1011  covEtaEta = numeratorEtaEta / denominator;
1012  covEtaPhi = numeratorEtaPhi / denominator;
1013  covPhiPhi = numeratorPhiPhi / denominator;
1014  } else {
1015  covEtaEta = 999.9;
1016  covEtaPhi = 999.9;
1017  covPhiPhi = 999.9;
1018  }
1019 
1020  } else {
1021  // Warn the user if there was no energy in the cells and return zeroes.
1022  // std::cout << "\ClusterShapeAlgo::Calculate_Covariances: no energy in supplied cells.\n";
1023  covEtaEta = 0;
1024  covEtaPhi = 0;
1025  covPhiPhi = 0;
1026  }
1027  std::array<float, 3> v{{covEtaEta, covEtaPhi, covPhiPhi}};
1028  return v;
1029 }
1030 
1031 //for covIEtaIEta,covIEtaIPhi and covIPhiIPhi are defined but only covIEtaIEta has been actively studied
1032 //instead of using absolute eta/phi it counts crystals normalised so that it gives identical results to normal covariances except near the cracks where of course its better
1033 //it also does not require any eta correction function in the endcap
1034 //it is multipled by an approprate crystal size to ensure it gives similar values to covariances(...)
1035 template <bool noZS>
1038  const CaloTopology *topology,
1039  float w0,
1041  float multEB,
1042  float multEE) {
1043  float e_5x5 = e5x5(cluster, recHits, topology);
1044  float covEtaEta, covEtaPhi, covPhiPhi;
1045 
1046  if (e_5x5 >= 0.) {
1047  //double w0_ = parameterMap_.find("W0")->second;
1048  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
1049  std::pair<float, float> mean5x5PosInNrCrysFromSeed = mean5x5PositionInLocalCrysCoord(cluster, recHits, topology);
1050  std::pair<float, float> mean5x5XYPos = mean5x5PositionInXY(cluster, recHits, topology);
1051 
1052  // now we can calculate the covariances
1053  double numeratorEtaEta = 0;
1054  double numeratorEtaPhi = 0;
1055  double numeratorPhiPhi = 0;
1056  double denominator = 0;
1057 
1058  //these allow us to scale the localCov by the crystal size
1059  //so that the localCovs have the same average value as the normal covs
1060  const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
1061  const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
1062 
1063  DetId seedId = getMaximum(v_id, recHits).first;
1064 
1065  bool isBarrel = seedId.subdetId() == EcalBarrel;
1066  const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
1067  float mult = isBarrel ? multEB : multEE; // we will multiply PF RecHit threshold by mult.
1068  // mult = 1 should work reasonably well for noise cleaning.
1069  // dedicated studies showed mult=1.25 works best for Run3 in endcap, where noise is high.
1070  // If no noise cleaning is intended then put mult=0.
1071 
1072  CaloRectangle rectangle{-2, 2, -2, 2};
1073  for (auto const &detId : rectangle(seedId, *topology)) {
1074  float frac = getFraction(v_id, detId);
1075  float energy = recHitEnergy(detId, recHits) * frac;
1076 
1077  if ((thresholds == nullptr) && (mult != 0.0)) {
1078  throw cms::Exception("EmptyPFRechHitThresColl")
1079  << "In EcalClusterTools::localCovariances, if EcalPFRecHitThresholds==nulptr, then multEB and multEE "
1080  "should be 0 as well.";
1081  }
1082  float rhThres = 0.0;
1083  if (thresholds != nullptr) {
1084  rhThres = (*thresholds)[detId]; // access PFRechit thresholds for noise cleaning
1085  }
1086 
1087  if (energy <= (rhThres * mult))
1088  continue;
1089 
1090  float dEta = getNrCrysDiffInEta(detId, seedId) - mean5x5PosInNrCrysFromSeed.first;
1091  float dPhi = 0;
1092 
1093  if (isBarrel)
1094  dPhi = getNrCrysDiffInPhi(detId, seedId) - mean5x5PosInNrCrysFromSeed.second;
1095  else
1096  dPhi = getDPhiEndcap(detId, mean5x5XYPos.first, mean5x5XYPos.second);
1097 
1098  double w = std::max(0.0f, w0 + std::log(energy / e_5x5));
1099 
1100  denominator += w;
1101  numeratorEtaEta += w * dEta * dEta;
1102  numeratorEtaPhi += w * dEta * dPhi;
1103  numeratorPhiPhi += w * dPhi * dPhi;
1104  }
1105 
1106  //multiplying by crysSize to make the values compariable to normal covariances
1107  if (denominator != 0.0) {
1108  covEtaEta = crysSize * crysSize * numeratorEtaEta / denominator;
1109  covEtaPhi = crysSize * crysSize * numeratorEtaPhi / denominator;
1110  covPhiPhi = crysSize * crysSize * numeratorPhiPhi / denominator;
1111  } else {
1112  covEtaEta = 999.9;
1113  covEtaPhi = 999.9;
1114  covPhiPhi = 999.9;
1115  }
1116 
1117  } else {
1118  // Warn the user if there was no energy in the cells and return zeroes.
1119  // std::cout << "\ClusterShapeAlgo::Calculate_Covariances: no energy in supplied cells.\n";
1120  covEtaEta = 0;
1121  covEtaPhi = 0;
1122  covPhiPhi = 0;
1123  }
1124  std::array<float, 3> v{{covEtaEta, covEtaPhi, covPhiPhi}};
1125  return v;
1126 }
1127 
1128 template <bool noZS>
1131  const CaloGeometry *geometry,
1132  double R0,
1133  bool logW,
1134  float w0) {
1135  return absZernikeMoment(cluster, recHits, geometry, 2, 0, R0, logW, w0);
1136 }
1137 
1138 template <bool noZS>
1141  const CaloGeometry *geometry,
1142  double R0,
1143  bool logW,
1144  float w0) {
1145  return absZernikeMoment(cluster, recHits, geometry, 4, 2, R0, logW, w0);
1146 }
1147 
1148 template <bool noZS>
1151  const CaloGeometry *geometry,
1152  int n,
1153  int m,
1154  double R0,
1155  bool logW,
1156  float w0) {
1157  // 1. Check if n,m are correctly
1158  if ((m > n) || ((n - m) % 2 != 0) || (n < 0) || (m < 0))
1159  return -1;
1160 
1161  // 2. Check if n,R0 are within validity Range :
1162  // n>20 or R0<2.19cm just makes no sense !
1163  if ((n > 20) || (R0 <= 2.19))
1164  return -1;
1165  if (n <= 5)
1166  return fast_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0);
1167  else
1168  return calc_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0);
1169 }
1170 
1171 template <bool noZS>
1174  const CaloGeometry *geometry,
1175  int n,
1176  int m,
1177  double R0,
1178  bool logW,
1179  float w0) {
1180  double r, ph, e, Re = 0, Im = 0;
1181  double TotalEnergy = cluster.energy();
1182  int index = (n / 2) * (n / 2) + (n / 2) + m;
1183  std::vector<EcalClusterEnergyDeposition> energyDistribution =
1184  getEnergyDepTopology(cluster, recHits, geometry, logW, w0);
1185  int clusterSize = energyDistribution.size();
1186  if (clusterSize < 3)
1187  return 0.0;
1188 
1189  for (int i = 0; i < clusterSize; i++) {
1190  r = energyDistribution[i].r / R0;
1191  if (r < 1) {
1192  std::vector<double> pol;
1193  pol.push_back(f00(r));
1194  pol.push_back(f11(r));
1195  pol.push_back(f20(r));
1196  pol.push_back(f22(r));
1197  pol.push_back(f31(r));
1198  pol.push_back(f33(r));
1199  pol.push_back(f40(r));
1200  pol.push_back(f42(r));
1201  pol.push_back(f44(r));
1202  pol.push_back(f51(r));
1203  pol.push_back(f53(r));
1204  pol.push_back(f55(r));
1205  ph = (energyDistribution[i]).phi;
1206  e = energyDistribution[i].deposited_energy;
1207  Re = Re + e / TotalEnergy * pol[index] * cos((double)m * ph);
1208  Im = Im - e / TotalEnergy * pol[index] * sin((double)m * ph);
1209  }
1210  }
1211  return sqrt(Re * Re + Im * Im);
1212 }
1213 
1214 template <bool noZS>
1217  const CaloGeometry *geometry,
1218  int n,
1219  int m,
1220  double R0,
1221  bool logW,
1222  float w0) {
1223  double r, ph, e, Re = 0, Im = 0, f_nm;
1224  double TotalEnergy = cluster.energy();
1225  std::vector<EcalClusterEnergyDeposition> energyDistribution =
1226  getEnergyDepTopology(cluster, recHits, geometry, logW, w0);
1227  int clusterSize = energyDistribution.size();
1228  if (clusterSize < 3)
1229  return 0.0;
1230 
1231  for (int i = 0; i < clusterSize; ++i) {
1232  r = energyDistribution[i].r / R0;
1233  if (r < 1) {
1234  ph = energyDistribution[i].phi;
1235  e = energyDistribution[i].deposited_energy;
1236  f_nm = 0;
1237  for (int s = 0; s <= (n - m) / 2; s++) {
1238  if (s % 2 == 0) {
1239  f_nm = f_nm + factorial(n - s) * pow(r, (double)(n - 2 * s)) /
1240  (factorial(s) * factorial((n + m) / 2 - s) * factorial((n - m) / 2 - s));
1241  } else {
1242  f_nm = f_nm - factorial(n - s) * pow(r, (double)(n - 2 * s)) /
1243  (factorial(s) * factorial((n + m) / 2 - s) * factorial((n - m) / 2 - s));
1244  }
1245  }
1246  Re = Re + e / TotalEnergy * f_nm * cos((double)m * ph);
1247  Im = Im - e / TotalEnergy * f_nm * sin((double)m * ph);
1248  }
1249  }
1250  return sqrt(Re * Re + Im * Im);
1251 }
1252 
1253 //returns the crystal 'eta' from the det id
1254 //it is defined as the number of crystals from the centre in the eta direction
1255 //for the barrel with its eta/phi geometry it is always integer
1256 //for the endcap it is fractional due to the x/y geometry
1257 template <bool noZS>
1259  if (id.det() == DetId::Ecal) {
1260  if (id.subdetId() == EcalBarrel) {
1261  EBDetId ebId(id);
1262  return ebId.ieta();
1263  } else if (id.subdetId() == EcalEndcap) {
1264  float iXNorm = getNormedIX(id);
1265  float iYNorm = getNormedIY(id);
1266 
1267  return std::sqrt(iXNorm * iXNorm + iYNorm * iYNorm);
1268  }
1269  }
1270  return 0.;
1271 }
1272 
1273 //returns the crystal 'phi' from the det id
1274 //it is defined as the number of crystals from the centre in the phi direction
1275 //for the barrel with its eta/phi geometry it is always integer
1276 //for the endcap it is not defined
1277 template <bool noZS>
1279  if (id.det() == DetId::Ecal) {
1280  if (id.subdetId() == EcalBarrel) {
1281  EBDetId ebId(id);
1282  return ebId.iphi();
1283  }
1284  }
1285  return 0.;
1286 }
1287 
1288 //want to map 1=-50,50=-1,51=1 and 100 to 50 so sub off one if zero or neg
1289 template <bool noZS>
1291  if (id.det() == DetId::Ecal && id.subdetId() == EcalEndcap) {
1292  EEDetId eeId(id);
1293  int iXNorm = eeId.ix() - 50;
1294  if (iXNorm <= 0)
1295  iXNorm--;
1296  return iXNorm;
1297  }
1298  return 0;
1299 }
1300 
1301 //want to map 1=-50,50=-1,51=1 and 100 to 50 so sub off one if zero or neg
1302 template <bool noZS>
1304  if (id.det() == DetId::Ecal && id.subdetId() == EcalEndcap) {
1305  EEDetId eeId(id);
1306  int iYNorm = eeId.iy() - 50;
1307  if (iYNorm <= 0)
1308  iYNorm--;
1309  return iYNorm;
1310  }
1311  return 0;
1312 }
1313 
1314 //nr crystals crysId is away from orgin id in eta
1315 template <bool noZS>
1316 float EcalClusterToolsT<noZS>::getNrCrysDiffInEta(const DetId &crysId, const DetId &orginId) {
1317  float crysIEta = getIEta(crysId);
1318  float orginIEta = getIEta(orginId);
1319  bool isBarrel = orginId.subdetId() == EcalBarrel;
1320 
1321  float nrCrysDiff = crysIEta - orginIEta;
1322 
1323  //no iEta=0 in barrel, so if go from positive to negative
1324  //need to reduce abs(detEta) by 1
1325  if (isBarrel) {
1326  if (crysIEta * orginIEta < 0) { // -1 to 1 transition
1327  if (crysIEta > 0)
1328  nrCrysDiff--;
1329  else
1330  nrCrysDiff++;
1331  }
1332  }
1333  return nrCrysDiff;
1334 }
1335 
1336 //nr crystals crysId is away from orgin id in phi
1337 template <bool noZS>
1338 float EcalClusterToolsT<noZS>::getNrCrysDiffInPhi(const DetId &crysId, const DetId &orginId) {
1339  float crysIPhi = getIPhi(crysId);
1340  float orginIPhi = getIPhi(orginId);
1341  bool isBarrel = orginId.subdetId() == EcalBarrel;
1342 
1343  float nrCrysDiff = crysIPhi - orginIPhi;
1344 
1345  if (isBarrel) { //if barrel, need to map into 0-180
1346  if (nrCrysDiff > +180) {
1347  nrCrysDiff = nrCrysDiff - 360;
1348  }
1349  if (nrCrysDiff < -180) {
1350  nrCrysDiff = nrCrysDiff + 360;
1351  }
1352  }
1353  return nrCrysDiff;
1354 }
1355 
1356 //nr crystals crysId is away from mean phi in 5x5 in phi
1357 template <bool noZS>
1358 float EcalClusterToolsT<noZS>::getDPhiEndcap(const DetId &crysId, float meanX, float meanY) {
1359  float iXNorm = getNormedIX(crysId);
1360  float iYNorm = getNormedIY(crysId);
1361 
1362  float hitLocalR2 = (iXNorm - meanX) * (iXNorm - meanX) + (iYNorm - meanY) * (iYNorm - meanY);
1363  float hitR2 = iXNorm * iXNorm + iYNorm * iYNorm;
1364  float meanR2 = meanX * meanX + meanY * meanY;
1365  float hitR = sqrt(hitR2);
1366  float meanR = sqrt(meanR2);
1367 
1368  float tmp = (hitR2 + meanR2 - hitLocalR2) / (2 * hitR * meanR);
1369  if (tmp < -1)
1370  tmp = -1;
1371  if (tmp > 1)
1372  tmp = 1;
1373  float phi = acos(tmp);
1374  float dPhi = hitR * phi;
1375 
1376  return dPhi;
1377 }
1378 
1379 template <bool noZS>
1382  const CaloTopology *topology,
1383  float w0) {
1384  const reco::BasicCluster bcluster = *(cluster.seed());
1385 
1386  float e_5x5 = e5x5(bcluster, recHits, topology);
1387  float covEtaEta, covEtaPhi, covPhiPhi;
1388 
1389  if (e_5x5 >= 0.) {
1390  const std::vector<std::pair<DetId, float>> &v_id = cluster.hitsAndFractions();
1391  std::pair<float, float> mean5x5PosInNrCrysFromSeed = mean5x5PositionInLocalCrysCoord(bcluster, recHits, topology);
1392  std::pair<float, float> mean5x5XYPos = mean5x5PositionInXY(cluster, recHits, topology);
1393  // now we can calculate the covariances
1394  double numeratorEtaEta = 0;
1395  double numeratorEtaPhi = 0;
1396  double numeratorPhiPhi = 0;
1397  double denominator = 0;
1398 
1399  const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
1400  const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
1401 
1402  DetId seedId = getMaximum(v_id, recHits).first;
1403  bool isBarrel = seedId.subdetId() == EcalBarrel;
1404 
1405  const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
1406 
1407  for (size_t i = 0; i < v_id.size(); ++i) {
1408  float frac = getFraction(v_id, v_id[i].first);
1409  float energy = recHitEnergy(v_id[i].first, recHits) * frac;
1410 
1411  if (energy <= 0)
1412  continue;
1413 
1414  float dEta = getNrCrysDiffInEta(v_id[i].first, seedId) - mean5x5PosInNrCrysFromSeed.first;
1415  float dPhi = 0;
1416  if (isBarrel)
1417  dPhi = getNrCrysDiffInPhi(v_id[i].first, seedId) - mean5x5PosInNrCrysFromSeed.second;
1418  else
1419  dPhi = getDPhiEndcap(v_id[i].first, mean5x5XYPos.first, mean5x5XYPos.second);
1420 
1421  double w = 0.;
1422  w = std::max(0.0f, w0 + std::log(energy / e_5x5));
1423 
1424  denominator += w;
1425  numeratorEtaEta += w * dEta * dEta;
1426  numeratorEtaPhi += w * dEta * dPhi;
1427  numeratorPhiPhi += w * dPhi * dPhi;
1428  }
1429 
1430  //multiplying by crysSize to make the values compariable to normal covariances
1431  if (denominator != 0.0) {
1432  covEtaEta = crysSize * crysSize * numeratorEtaEta / denominator;
1433  covEtaPhi = crysSize * crysSize * numeratorEtaPhi / denominator;
1434  covPhiPhi = crysSize * crysSize * numeratorPhiPhi / denominator;
1435  } else {
1436  covEtaEta = 999.9;
1437  covEtaPhi = 999.9;
1438  covPhiPhi = 999.9;
1439  }
1440 
1441  } else {
1442  // Warn the user if there was no energy in the cells and return zeroes.
1443  // std::cout << "\ClusterShapeAlgo::Calculate_Covariances: no energy in supplied cells.\n";
1444  covEtaEta = 0;
1445  covEtaPhi = 0;
1446  covPhiPhi = 0;
1447  }
1448 
1449  std::array<float, 3> v{{covEtaEta, covEtaPhi, covPhiPhi}};
1450  return v;
1451 }
1452 
1453 // compute cluster second moments with respect to principal axes (eigenvectors of sEtaEta, sPhiPhi, sEtaPhi matrix)
1454 // store also angle alpha between major axis and phi.
1455 // takes into account shower elongation in phi direction due to magnetic field effect:
1456 // default value of 0.8 ensures sMaj = sMin for unconverted photons
1457 // (if phiCorrectionFactor=1 sMaj > sMin and alpha=0 also for unconverted photons)
1458 template <bool noZS>
1461  double phiCorrectionFactor,
1462  double w0,
1463  bool useLogWeights) {
1464  if (recHits.empty())
1465  return Cluster2ndMoments();
1466 
1467  std::vector<std::pair<const EcalRecHit *, float>> RH_ptrs_fracs;
1468 
1469  const std::vector<std::pair<DetId, float>> &myHitsPair = basicCluster.hitsAndFractions();
1470 
1471  for (unsigned int i = 0; i < myHitsPair.size(); i++) {
1472  //get pointer to recHit object
1473  EcalRecHitCollection::const_iterator myRH = recHits.find(myHitsPair[i].first);
1474  if (myRH != recHits.end()) {
1475  RH_ptrs_fracs.push_back(std::make_pair(&(*myRH), myHitsPair[i].second));
1476  }
1477  }
1478 
1479  return EcalClusterToolsT<noZS>::cluster2ndMoments(RH_ptrs_fracs, phiCorrectionFactor, w0, useLogWeights);
1480 }
1481 
1482 template <bool noZS>
1485  double phiCorrectionFactor,
1486  double w0,
1487  bool useLogWeights) {
1489  *(superCluster.seed()), recHits, phiCorrectionFactor, w0, useLogWeights);
1490 }
1491 
1492 template <bool noZS>
1494  const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs,
1495  double phiCorrectionFactor,
1496  double w0,
1497  bool useLogWeights) {
1498  if (RH_ptrs_fracs.empty())
1499  return Cluster2ndMoments();
1500 
1501  double mid_eta(0), mid_phi(0), mid_x(0), mid_y(0);
1502 
1503  double Etot = EcalClusterToolsT<noZS>::getSumEnergy(RH_ptrs_fracs);
1504 
1505  double max_phi = -10.;
1506  double min_phi = 100.;
1507 
1508  std::vector<double> etaDetId;
1509  std::vector<double> phiDetId;
1510  std::vector<double> xDetId;
1511  std::vector<double> yDetId;
1512  std::vector<double> wiDetId;
1513 
1514  unsigned int nCry = 0;
1515  double denominator = 0.;
1516  bool isBarrel(true);
1517 
1518  // loop over rechits and compute weights:
1519  for (std::vector<std::pair<const EcalRecHit *, float>>::const_iterator rhf_ptr = RH_ptrs_fracs.begin();
1520  rhf_ptr != RH_ptrs_fracs.end();
1521  rhf_ptr++) {
1522  const EcalRecHit *rh_ptr = rhf_ptr->first;
1523 
1524  //get iEta, iPhi
1525  double temp_eta(0), temp_phi(0), temp_x(0), temp_y(0);
1526  isBarrel = rh_ptr->detid().subdetId() == EcalBarrel;
1527 
1528  if (isBarrel) {
1529  temp_eta = (getIEta(rh_ptr->detid()) > 0. ? getIEta(rh_ptr->detid()) + 84.5 : getIEta(rh_ptr->detid()) + 85.5);
1530  temp_phi = getIPhi(rh_ptr->detid()) - 0.5;
1531  } else {
1532  temp_eta = getIEta(rh_ptr->detid());
1533  temp_x = getNormedIX(rh_ptr->detid());
1534  temp_y = getNormedIY(rh_ptr->detid());
1535  }
1536 
1537  double temp_ene = rh_ptr->energy() * (noZS ? 1.0 : rhf_ptr->second);
1538 
1539  double temp_wi = ((useLogWeights) ? std::max(0.0, w0 + std::log(std::abs(temp_ene) / Etot)) : temp_ene);
1540 
1541  if (temp_phi > max_phi)
1542  max_phi = temp_phi;
1543  if (temp_phi < min_phi)
1544  min_phi = temp_phi;
1545  etaDetId.push_back(temp_eta);
1546  phiDetId.push_back(temp_phi);
1547  xDetId.push_back(temp_x);
1548  yDetId.push_back(temp_y);
1549  wiDetId.push_back(temp_wi);
1550  denominator += temp_wi;
1551  nCry++;
1552  }
1553 
1554  if (isBarrel) {
1555  // correct phi wrap-around:
1556  if (max_phi == 359.5 && min_phi == 0.5) {
1557  for (unsigned int i = 0; i < nCry; i++) {
1558  if (phiDetId[i] - 179. > 0.)
1559  phiDetId[i] -= 360.;
1560  mid_phi += phiDetId[i] * wiDetId[i];
1561  mid_eta += etaDetId[i] * wiDetId[i];
1562  }
1563  } else {
1564  for (unsigned int i = 0; i < nCry; i++) {
1565  mid_phi += phiDetId[i] * wiDetId[i];
1566  mid_eta += etaDetId[i] * wiDetId[i];
1567  }
1568  }
1569  } else {
1570  for (unsigned int i = 0; i < nCry; i++) {
1571  mid_eta += etaDetId[i] * wiDetId[i];
1572  mid_x += xDetId[i] * wiDetId[i];
1573  mid_y += yDetId[i] * wiDetId[i];
1574  }
1575  }
1576 
1577  mid_eta /= denominator;
1578  mid_phi /= denominator;
1579  mid_x /= denominator;
1580  mid_y /= denominator;
1581 
1582  // See = sigma eta eta
1583  // Spp = (B field corrected) sigma phi phi
1584  // See = (B field corrected) sigma eta phi
1585  double See = 0.;
1586  double Spp = 0.;
1587  double Sep = 0.;
1588  double deta(0), dphi(0);
1589  // compute (phi-corrected) covariance matrix:
1590  for (unsigned int i = 0; i < nCry; i++) {
1591  if (isBarrel) {
1592  deta = etaDetId[i] - mid_eta;
1593  dphi = phiDetId[i] - mid_phi;
1594  } else {
1595  deta = etaDetId[i] - mid_eta;
1596  float hitLocalR2 = (xDetId[i] - mid_x) * (xDetId[i] - mid_x) + (yDetId[i] - mid_y) * (yDetId[i] - mid_y);
1597  float hitR2 = xDetId[i] * xDetId[i] + yDetId[i] * yDetId[i];
1598  float meanR2 = mid_x * mid_x + mid_y * mid_y;
1599  float hitR = sqrt(hitR2);
1600  float meanR = sqrt(meanR2);
1601  float phi = acos((hitR2 + meanR2 - hitLocalR2) / (2 * hitR * meanR));
1602  dphi = hitR * phi;
1603  }
1604  See += (wiDetId[i] * deta * deta) / denominator;
1605  Spp += phiCorrectionFactor * (wiDetId[i] * dphi * dphi) / denominator;
1606  Sep += sqrt(phiCorrectionFactor) * (wiDetId[i] * deta * dphi) / denominator;
1607  }
1608 
1609  Cluster2ndMoments returnMoments;
1610 
1611  // compute matrix eigenvalues:
1612  returnMoments.sMaj = ((See + Spp) + sqrt((See - Spp) * (See - Spp) + 4. * Sep * Sep)) / 2.;
1613  returnMoments.sMin = ((See + Spp) - sqrt((See - Spp) * (See - Spp) + 4. * Sep * Sep)) / 2.;
1614 
1615  returnMoments.alpha = atan((See - Spp + sqrt((Spp - See) * (Spp - See) + 4. * Sep * Sep)) / (2. * Sep));
1616 
1617  return returnMoments;
1618 }
1619 
1620 //compute shower shapes: roundness and angle in a vector. Roundness is 0th element, Angle is 1st element.
1621 //description: uses classical mechanics inertia tensor.
1622 // roundness is smaller_eValue/larger_eValue
1623 // angle is the angle from the iEta axis to the smallest eVector (a.k.a. the shower's elongated axis)
1624 // this function uses only recHits belonging to a SC above energyThreshold (default 0)
1625 // you can select linear weighting = energy_recHit/total_energy (weightedPositionMethod=0) default
1626 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
1627 template <bool noZS>
1629  const reco::SuperCluster &superCluster,
1631  int weightedPositionMethod,
1632  float energyThreshold) { //int positionWeightingMethod=0){
1633  std::vector<std::pair<const EcalRecHit *, float>> RH_ptrs_fracs;
1634  const std::vector<std::pair<DetId, float>> &myHitsPair = superCluster.hitsAndFractions();
1635  for (unsigned int i = 0; i < myHitsPair.size(); ++i) {
1636  //get pointer to recHit object
1637  EcalRecHitCollection::const_iterator myRH = recHits.find(myHitsPair[i].first);
1638  if (myRH != recHits.end() && myRH->energy() * (noZS ? 1.0 : myHitsPair[i].second) > energyThreshold) {
1639  //require rec hit to have positive energy
1640  RH_ptrs_fracs.push_back(std::make_pair(&(*myRH), myHitsPair[i].second));
1641  }
1642  }
1643  std::vector<float> temp =
1644  EcalClusterToolsT<noZS>::roundnessSelectedBarrelRecHits(RH_ptrs_fracs, weightedPositionMethod);
1645  return temp;
1646 }
1647 
1648 // this function uses all recHits within specified window ( with default values ieta_delta=2, iphi_delta=5) around SC's highest recHit.
1649 // recHits must pass an energy threshold "energyRHThresh" (default 0)
1650 // you can select linear weighting = energy_recHit/total_energy (weightedPositionMethod=0)
1651 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
1652 template <bool noZS>
1654  const reco::SuperCluster &superCluster,
1656  int ieta_delta,
1657  int iphi_delta,
1658  float energyRHThresh,
1659  int weightedPositionMethod) {
1660  std::vector<std::pair<const EcalRecHit *, float>> RH_ptrs_fracs;
1661  const std::vector<std::pair<DetId, float>> &myHitsPair = superCluster.hitsAndFractions();
1662  for (unsigned int i = 0; i < myHitsPair.size(); ++i) {
1663  //get pointer to recHit object
1664  EcalRecHitCollection::const_iterator myRH = recHits.find(myHitsPair[i].first);
1665  if (myRH != recHits.end() && myRH->energy() * (noZS ? 1.0 : myHitsPair[i].second) > energyRHThresh)
1666  RH_ptrs_fracs.push_back(std::make_pair(&(*myRH), myHitsPair[i].second));
1667  }
1668 
1669  std::vector<int> seedPosition = EcalClusterToolsT<noZS>::getSeedPosition(RH_ptrs_fracs);
1670 
1671  for (EcalRecHitCollection::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
1672  EBDetId EBdetIdi(rh->detid());
1673  float the_fraction = 0;
1674  //if(rh != recHits.end())
1675  bool inEtaWindow = (abs(deltaIEta(seedPosition[0], EBdetIdi.ieta())) <= ieta_delta);
1676  bool inPhiWindow = (abs(deltaIPhi(seedPosition[1], EBdetIdi.iphi())) <= iphi_delta);
1677  bool passEThresh = (rh->energy() > energyRHThresh);
1678  bool alreadyCounted = false;
1679 
1680  // figure out if the rechit considered now is already inside the SC
1681  bool is_SCrh_inside_recHits = false;
1682  for (unsigned int i = 0; i < myHitsPair.size(); i++) {
1683  EcalRecHitCollection::const_iterator SCrh = recHits.find(myHitsPair[i].first);
1684  if (SCrh != recHits.end()) {
1685  the_fraction = myHitsPair[i].second;
1686  is_SCrh_inside_recHits = true;
1687  if (rh->detid() == SCrh->detid())
1688  alreadyCounted = true;
1689  }
1690  } //for loop over SC's recHits
1691 
1692  if (is_SCrh_inside_recHits && !alreadyCounted && passEThresh && inEtaWindow && inPhiWindow) {
1693  RH_ptrs_fracs.push_back(std::make_pair(&(*rh), the_fraction));
1694  }
1695 
1696  } //for loop over rh
1697  return EcalClusterToolsT<noZS>::roundnessSelectedBarrelRecHits(RH_ptrs_fracs, weightedPositionMethod);
1698 }
1699 
1700 // this function computes the roundness and angle variables for vector of pointers to recHits you pass it
1701 // you can select linear weighting = energy_recHit/total_energy (weightedPositionMethod=0)
1702 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
1703 template <bool noZS>
1705  const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs,
1706  int weightedPositionMethod) { //int weightedPositionMethod = 0){
1707  //positionWeightingMethod = 0 linear weighting, 1 log energy weighting
1708 
1709  std::vector<float> shapes; // this is the returning vector
1710 
1711  //make sure photon has more than one crystal; else roundness and angle suck
1712  if (RH_ptrs_fracs.size() < 2) {
1713  shapes.push_back(-3);
1714  shapes.push_back(-3);
1715  return shapes;
1716  }
1717 
1718  //Find highest E RH (Seed) and save info, compute sum total energy used
1719  std::vector<int> seedPosition = EcalClusterToolsT<noZS>::getSeedPosition(RH_ptrs_fracs); // *recHits);
1720  int tempInt = seedPosition[0];
1721  if (tempInt < 0)
1722  tempInt++;
1723  float energyTotal = EcalClusterToolsT<noZS>::getSumEnergy(RH_ptrs_fracs);
1724 
1725  //1st loop over rechits: compute new weighted center position in coordinates centered on seed
1726  float centerIEta = 0.;
1727  float centerIPhi = 0.;
1728  float denominator = 0.;
1729 
1730  for (std::vector<std::pair<const EcalRecHit *, float>>::const_iterator rhf_ptr = RH_ptrs_fracs.begin();
1731  rhf_ptr != RH_ptrs_fracs.end();
1732  rhf_ptr++) {
1733  const EcalRecHit *rh_ptr = rhf_ptr->first;
1734  //get iEta, iPhi
1735  EBDetId EBdetIdi(rh_ptr->detid());
1736  if (fabs(energyTotal) < 0.0001) {
1737  // don't /0, bad!
1738  shapes.push_back(-2);
1739  shapes.push_back(-2);
1740  return shapes;
1741  }
1742  float rh_energy = rh_ptr->energy() * (noZS ? 1.0 : rhf_ptr->second);
1743  float weight = 0;
1744  if (std::abs(weightedPositionMethod) < 0.0001) { //linear
1745  weight = rh_energy / energyTotal;
1746  } else { //logrithmic
1747  weight = std::max(0.0, 4.2 + log(rh_energy / energyTotal));
1748  }
1749  denominator += weight;
1750  centerIEta += weight * deltaIEta(seedPosition[0], EBdetIdi.ieta());
1751  centerIPhi += weight * deltaIPhi(seedPosition[1], EBdetIdi.iphi());
1752  }
1753  if (fabs(denominator) < 0.0001) {
1754  // don't /0, bad!
1755  shapes.push_back(-2);
1756  shapes.push_back(-2);
1757  return shapes;
1758  }
1759  centerIEta = centerIEta / denominator;
1760  centerIPhi = centerIPhi / denominator;
1761 
1762  //2nd loop over rechits: compute inertia tensor
1763  TMatrixDSym inertia(2); //initialize 2d inertia tensor
1764  double inertia00 = 0.;
1765  double inertia01 = 0.; // = inertia10 b/c matrix should be symmetric
1766  double inertia11 = 0.;
1767  int i = 0;
1768  for (std::vector<std::pair<const EcalRecHit *, float>>::const_iterator rhf_ptr = RH_ptrs_fracs.begin();
1769  rhf_ptr != RH_ptrs_fracs.end();
1770  rhf_ptr++) {
1771  const EcalRecHit *rh_ptr = rhf_ptr->first;
1772  //get iEta, iPhi
1773  EBDetId EBdetIdi(rh_ptr->detid());
1774 
1775  if (fabs(energyTotal) < 0.0001) {
1776  // don't /0, bad!
1777  shapes.push_back(-2);
1778  shapes.push_back(-2);
1779  return shapes;
1780  }
1781  float rh_energy = rh_ptr->energy() * (noZS ? 1.0 : rhf_ptr->second);
1782  float weight = 0;
1783  if (std::abs(weightedPositionMethod) < 0.0001) { //linear
1784  weight = rh_energy / energyTotal;
1785  } else { //logrithmic
1786  weight = std::max(0.0, 4.2 + log(rh_energy / energyTotal));
1787  }
1788 
1789  float ieta_rh_to_center = deltaIEta(seedPosition[0], EBdetIdi.ieta()) - centerIEta;
1790  float iphi_rh_to_center = deltaIPhi(seedPosition[1], EBdetIdi.iphi()) - centerIPhi;
1791 
1792  inertia00 += weight * iphi_rh_to_center * iphi_rh_to_center;
1793  inertia01 -= weight * iphi_rh_to_center * ieta_rh_to_center;
1794  inertia11 += weight * ieta_rh_to_center * ieta_rh_to_center;
1795  i++;
1796  }
1797 
1798  inertia[0][0] = inertia00;
1799  inertia[0][1] = inertia01; // use same number here
1800  inertia[1][0] = inertia01; // and here to insure symmetry
1801  inertia[1][1] = inertia11;
1802 
1803  //step 1 find principal axes of inertia
1804  TMatrixD eVectors(2, 2);
1805  TVectorD eValues(2);
1806  //std::cout<<"EcalClusterToolsT<noZS>::showerRoundness- about to compute eVectors"<<std::endl;
1807  eVectors = inertia.EigenVectors(eValues); //ordered highest eV to lowest eV (I checked!)
1808  //and eVectors are in columns of matrix! I checked!
1809  //and they are normalized to 1
1810 
1811  //step 2 select eta component of smaller eVal's eVector
1812  TVectorD smallerAxis(2); //easiest to spin SC on this axis (smallest eVal)
1813  smallerAxis[0] = eVectors[0][1]; //row,col //eta component
1814  smallerAxis[1] = eVectors[1][1]; //phi component
1815 
1816  //step 3 compute interesting quatities
1817  Double_t temp = fabs(smallerAxis[0]); // closer to 1 ->beamhalo, closer to 0 something else
1818  if (fabs(eValues[0]) < 0.0001) {
1819  // don't /0, bad!
1820  shapes.push_back(-2);
1821  shapes.push_back(-2);
1822  return shapes;
1823  }
1824 
1825  float Roundness = eValues[1] / eValues[0];
1826  float Angle = acos(temp);
1827 
1828  if (-0.00001 < Roundness && Roundness < 0)
1829  Roundness = 0.;
1830  if (-0.00001 < Angle && Angle < 0)
1831  Angle = 0.;
1832 
1833  shapes.push_back(Roundness);
1834  shapes.push_back(Angle);
1835  return shapes;
1836 }
1837 
1838 template <bool noZS>
1841  const CaloTopology *topology) {
1842  int nrSat = 0;
1843  CaloRectangle rectangle{-2, 2, -2, 2};
1844  for (auto const &detId : rectangle(id, *topology)) {
1845  auto recHitIt = recHits->find(detId);
1846  if (recHitIt != recHits->end() && recHitIt->checkFlag(EcalRecHit::kSaturated)) {
1847  nrSat++;
1848  }
1849  }
1850  return nrSat;
1851 }
1852 
1853 //private functions useful for roundnessBarrelSuperClusters etc.
1854 //compute delta iphi between a seed and a particular recHit
1855 //iphi [1,360]
1856 //safe gaurds are put in to ensure the difference is between [-180,180]
1857 template <bool noZS>
1858 int EcalClusterToolsT<noZS>::deltaIPhi(int seed_iphi, int rh_iphi) {
1859  int rel_iphi = rh_iphi - seed_iphi;
1860  // take care of cyclic variable iphi [1,360]
1861  if (rel_iphi > 180)
1862  rel_iphi = rel_iphi - 360;
1863  if (rel_iphi < -180)
1864  rel_iphi = rel_iphi + 360;
1865  return rel_iphi;
1866 }
1867 
1868 //compute delta ieta between a seed and a particular recHit
1869 //ieta [-85,-1] and [1,85]
1870 //safe gaurds are put in to shift the negative ieta +1 to make an ieta=0 so differences are computed correctly
1871 template <bool noZS>
1872 int EcalClusterToolsT<noZS>::deltaIEta(int seed_ieta, int rh_ieta) {
1873  // get rid of the fact that there is no ieta=0
1874  if (seed_ieta < 0)
1875  seed_ieta++;
1876  if (rh_ieta < 0)
1877  rh_ieta++;
1878  int rel_ieta = rh_ieta - seed_ieta;
1879  return rel_ieta;
1880 }
1881 
1882 //return (ieta,iphi) of highest energy recHit of the recHits passed to this function
1883 template <bool noZS>
1885  const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs) {
1886  std::vector<int> seedPosition;
1887  float eSeedRH = 0;
1888  int iEtaSeedRH = 0;
1889  int iPhiSeedRH = 0;
1890 
1891  for (auto const &rhf : RH_ptrs_fracs) {
1892  const EcalRecHit *rh_ptr = rhf.first;
1893  //get iEta, iPhi
1894  EBDetId EBdetIdi(rh_ptr->detid());
1895  float rh_energy = rh_ptr->energy() * (noZS ? 1.0 : rhf.second);
1896 
1897  if (eSeedRH < rh_energy) {
1898  eSeedRH = rh_energy;
1899  iEtaSeedRH = EBdetIdi.ieta();
1900  iPhiSeedRH = EBdetIdi.iphi();
1901  }
1902 
1903  } // for loop
1904 
1905  seedPosition.push_back(iEtaSeedRH);
1906  seedPosition.push_back(iPhiSeedRH);
1907  return seedPosition;
1908 }
1909 
1910 // return the total energy of the recHits passed to this function
1911 template <bool noZS>
1912 float EcalClusterToolsT<noZS>::getSumEnergy(const std::vector<std::pair<const EcalRecHit *, float>> &RH_ptrs_fracs) {
1913  float sumE = 0.;
1914  for (const auto &hAndF : RH_ptrs_fracs) {
1915  sumE += hAndF.first->energy() * (noZS ? 1.0 : hAndF.second);
1916  }
1917  return sumE;
1918 }
1919 
1921 
1922 namespace noZS {
1924 }
1925 
1926 #endif
EBDetId::kCrystalsInPhi
static const int kCrystalsInPhi
Definition: EBDetId.h:142
EcalClusterToolsT::EcalClusterEnergyDeposition::phi
double phi
Definition: EcalClusterTools.h:293
change_name.diff
diff
Definition: change_name.py:13
EBDetId::positiveZ
bool positiveZ() const
Definition: EBDetId.h:76
EcalClusterToolsT::n5x5
static int n5x5(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:547
EcalRecHit
Definition: EcalRecHit.h:15
EBDetId::ieta
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
mps_fire.i
i
Definition: mps_fire.py:428
EcalBarrelGeometry.h
edm::SortedCollection< EcalRecHit >::const_iterator
std::vector< EcalRecHit >::const_iterator const_iterator
Definition: SortedCollection.h:80
EcalClusterToolsT::nrSaturatedCrysIn5x5
static int nrSaturatedCrysIn5x5(const DetId &id, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:1839
MessageLogger.h
EcalClusterToolsT::zernike20
static double zernike20(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0=6.6, bool logW=true, float w0=4.7)
Definition: EcalClusterTools.h:1129
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
EcalRecHit::detid
const DetId & detid() const
Definition: EcalRecHit.h:72
particleFlowZeroSuppressionECAL_cff.thresholds
thresholds
Definition: particleFlowZeroSuppressionECAL_cff.py:31
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
EcalClusterTools
EcalClusterToolsT< false > EcalClusterTools
Definition: EcalClusterTools.h:1920
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
reco::SuperCluster
Definition: SuperCluster.h:18
BasicCluster.h
EcalClusterToolsT::calc_AbsZernikeMoment
static double calc_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0)
Definition: EcalClusterTools.h:1215
EcalClusterToolsT::e2x5Bottom
static float e2x5Bottom(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:598
EBDetId
Definition: EBDetId.h:17
mps_merge.weight
weight
Definition: mps_merge.py:88
geometry
Definition: geometry.py:1
EBDetId.h
EEDetId.h
factorial
int factorial(int n)
factorial function
Definition: L1GlobalTriggerFunctions.cc:13
EcalRecHit::kSaturated
Definition: EcalRecHit.h:27
EBDetId::im
int im() const
get the number of module inside the SM (1-4)
Definition: EBDetId.h:64
EcalClusterToolsT::eMax
static float eMax(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:555
EcalClusterToolsT::fast_AbsZernikeMoment
static double fast_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0)
Definition: EcalClusterTools.h:1172
edm::SortedCollection< EcalRecHit >
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
EcalClusterToolsT::roundnessBarrelSuperClustersUserExtended
static std::vector< float > roundnessBarrelSuperClustersUserExtended(const reco::SuperCluster &superCluster, const EcalRecHitCollection &recHits, int ieta_delta=0, int iphi_delta=0, float energyRHThresh=0.00000, int weightedPositionMethod=0)
Definition: EcalClusterTools.h:1653
EcalClusterToolsT::matrixEnergy
static float matrixEnergy(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, DetId id, CaloRectangle rectangle)
Definition: EcalClusterTools.h:449
EcalRecHit::energy
float energy() const
Definition: EcalRecHit.h:68
findQualityFiles.v
v
Definition: findQualityFiles.py:179
EcalClusterToolsT::getEnergyDepTopology
static std::vector< EcalClusterEnergyDeposition > getEnergyDepTopology(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW, float w0)
Definition: EcalClusterTools.h:731
HLT_FULL_cff.dPhi
dPhi
Definition: HLT_FULL_cff.py:13695
EcalClusterToolsT::localCovariances
static std::array< float, 3 > localCovariances(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, float w0=EgammaLocalCovParamDefaults::kRelEnCut, const EcalPFRecHitThresholds *thresholds=nullptr, float multEB=0.0, float multEE=0.0)
Definition: EcalClusterTools.h:1036
EcalRecHit::kTowerRecovered
Definition: EcalRecHit.h:30
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
EcalCondObjectContainer
Definition: EcalCondObjectContainer.h:13
EEDetId::ix
int ix() const
Definition: EEDetId.h:77
groupFilesInBlocks.temp
list temp
Definition: groupFilesInBlocks.py:142
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
CaloTopology
Definition: CaloTopology.h:19
EcalClusterToolsT::EcalClusterEnergyDeposition::deposited_energy
double deposited_energy
Definition: EcalClusterTools.h:291
EcalBarrel
Definition: EcalSubdetector.h:10
EcalRecHitCollections.h
EBDetId::kModulesPerSM
static const int kModulesPerSM
Definition: EBDetId.h:140
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
CastorDataFrameFilter_impl::energySum
double energySum(const DataFrame &df, int fs, int ls)
Definition: CastorDataFrameFilter.cc:22
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
EcalClusterToolsT::getMaximum
static std::pair< DetId, float > getMaximum(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:410
DetId
Definition: DetId.h:17
EcalClusterToolsT::f33
static double f33(double r)
Definition: EcalClusterTools.h:322
alignCSCRings.s
s
Definition: alignCSCRings.py:92
CaloGeometry
Definition: CaloGeometry.h:21
EcalPFRecHitThresholds.h
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
EcalClusterToolsT::e2x5Max
static float e2x5Max(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:608
DivergingColor.frac
float frac
Definition: DivergingColor.py:175
EcalClusterToolsT::f11
static double f11(double r)
Definition: EcalClusterTools.h:318
EcalClusterToolsT::roundnessBarrelSuperClusters
static std::vector< float > roundnessBarrelSuperClusters(const reco::SuperCluster &superCluster, const EcalRecHitCollection &recHits, int weightedPositionMethod=0, float energyThreshold=0.0)
Definition: EcalClusterTools.h:1628
EcalClusterToolsT::e2x5Top
static float e2x5Top(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:590
EcalClusterToolsT::matrixSize
static int matrixSize(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, DetId id, CaloRectangle rectangle)
Definition: EcalClusterTools.h:465
EcalClusterToolsT::getIEta
static float getIEta(const DetId &id)
Definition: EcalClusterTools.h:1258
EcalClusterToolsT::e3x1
static float e3x1(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:649
EcalClusterToolsT::e2x2
static float e2x2(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:494
w
const double w
Definition: UKUtility.cc:23
EcalClusterToolsT::matrixDetId
static std::vector< DetId > matrixDetId(const CaloTopology *topology, DetId id, CaloRectangle rectangle)
Definition: EcalClusterTools.h:482
EcalClusterToolsT::getIPhi
static float getIPhi(const DetId &id)
Definition: EcalClusterTools.h:1278
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:78
EcalClusterToolsT::deltaIEta
static int deltaIEta(int seed_ieta, int rh_ieta)
Definition: EcalClusterTools.h:1872
reco::CaloCluster
Definition: CaloCluster.h:31
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
EcalClusterToolsT::covariances
static std::array< float, 3 > covariances(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, const CaloGeometry *geometry, float w0=4.7)
Definition: EcalClusterTools.h:962
EcalClusterToolsT::matrixDetId
static std::vector< DetId > matrixDetId(const CaloTopology *topology, DetId id, int size)
Definition: EcalClusterTools.h:247
EcalClusterToolsT::EcalClusterEnergyDeposition::r
double r
Definition: EcalClusterTools.h:292
Geom::pi
constexpr double pi()
Definition: Pi.h:31
EcalClusterToolsT::getNrCrysDiffInPhi
static float getNrCrysDiffInPhi(const DetId &crysId, const DetId &orginId)
Definition: EcalClusterTools.h:1338
EcalClusterToolsT::deltaIPhi
static int deltaIPhi(int seed_iphi, int rh_iphi)
Definition: EcalClusterTools.h:1858
Geom::twoPi
constexpr double twoPi()
Definition: Pi.h:32
EBDetId::MAX_IPHI
static const int MAX_IPHI
Definition: EBDetId.h:137
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
Cluster2ndMoments::alpha
float alpha
Definition: EcalClusterTools.h:71
CaloRectangle
Definition: CaloRectangle.h:13
callgraph.n2
n2
Definition: callgraph.py:50
Point3DBase< float, GlobalTag >
EcalClusterToolsT::f40
static double f40(double r)
Definition: EcalClusterTools.h:323
EcalClusterToolsT::mean5x5PositionInXY
static std::pair< float, float > mean5x5PositionInXY(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:929
EcalClusterToolsT::cluster2ndMoments
static Cluster2ndMoments cluster2ndMoments(const reco::BasicCluster &basicCluster, const EcalRecHitCollection &recHits, double phiCorrectionFactor=0.8, double w0=4.7, bool useLogWeights=true)
Definition: EcalClusterTools.h:1459
EcalSeverityLevelAlgo.h
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:125
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
EEDetId
Definition: EEDetId.h:14
EcalClusterToolsT::eTop
static float eTop(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:673
EcalClusterToolsT::e5x1
static float e5x1(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:633
EcalEndcap
Definition: EcalSubdetector.h:10
FastTrackerRecHitMaskProducer_cfi.recHits
recHits
Definition: FastTrackerRecHitMaskProducer_cfi.py:8
EgammaLocalCovParamDefaults::kRelEnCut
static constexpr float kRelEnCut
Definition: EgammaLocalCovParamDefaults.h:5
EcalClusterToolsT::f31
static double f31(double r)
Definition: EcalClusterTools.h:321
PixelPluginsPhase0_cfi.isBarrel
isBarrel
Definition: PixelPluginsPhase0_cfi.py:17
EcalClusterToolsT::getNormedIY
static float getNormedIY(const DetId &id)
Definition: EcalClusterTools.h:1303
EcalClusterToolsT::getSumEnergy
static float getSumEnergy(const std::vector< std::pair< const EcalRecHit *, float >> &RH_ptrs_fracs)
Definition: EcalClusterTools.h:1912
EcalClusterToolsT::f44
static double f44(double r)
Definition: EcalClusterTools.h:325
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
Cluster2ndMoments::sMaj
float sMaj
Definition: EcalClusterTools.h:68
reco::SuperCluster::seed
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:77
PVValHelper::phi
Definition: PVValidationHelpers.h:69
reco::CaloCluster::hitsAndFractions
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
EcalClusterToolsT::eBottom
static float eBottom(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:681
EcalClusterToolsT::e4x4
static float e4x4(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:527
EcalClusterToolsT::getNrCrysDiffInEta
static float getNrCrysDiffInEta(const DetId &crysId, const DetId &orginId)
Definition: EcalClusterTools.h:1316
EcalClusterToolsT::getDPhiEndcap
static float getDPhiEndcap(const DetId &crysId, float meanX, float meanY)
Definition: EcalClusterTools.h:1358
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
EcalClusterToolsT::f55
static double f55(double r)
Definition: EcalClusterTools.h:328
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
EcalPFRecHitThresholdsRcd.h
EcalClusterToolsT::getFraction
static float getFraction(const std::vector< std::pair< DetId, float >> &v_id, DetId id)
Definition: EcalClusterTools.h:381
EcalClusterToolsT::e3x2
static float e3x2(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:507
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
EBDetId::kTowersInPhi
static const int kTowersInPhi
Definition: EBDetId.h:139
EcalClusterToolsT::absZernikeMoment
static double absZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0)
Definition: EcalClusterTools.h:1149
EcalClusterToolsT::mean5x5PositionInLocalCrysCoord
static std::pair< float, float > mean5x5PositionInLocalCrysCoord(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:893
CaloSubdetectorGeometry::getGeometry
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
Definition: CaloSubdetectorGeometry.cc:36
EcalClusterToolsT::e2nd
static float e2nd(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:560
EcalRecHit::kWeird
Definition: EcalRecHit.h:35
DetId::Ecal
Definition: DetId.h:27
EcalClusterToolsT::energyBasketFractionPhi
static std::vector< float > energyBasketFractionPhi(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:709
EcalClusterToolsT::getSeedPosition
static std::vector< int > getSeedPosition(const std::vector< std::pair< const EcalRecHit *, float >> &RH_ptrs)
Definition: EcalClusterTools.h:1884
res
Definition: Electron.h:6
EEDetId::iy
int iy() const
Definition: EEDetId.h:83
EcalClusterToolsT::f00
static double f00(double r)
Definition: EcalClusterTools.h:317
EcalClusterToolsT::eRight
static float eRight(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:665
EcalClusterToolsT::e2x5Right
static float e2x5Right(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:574
EcalClusterToolsT
Definition: EcalClusterTools.h:76
reco::CaloCluster::position
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
alignCSCRings.r
r
Definition: alignCSCRings.py:93
CaloTopology.h
EcalClusterToolsT::e3x3
static float e3x3(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:519
EcalClusterToolsT::f20
static double f20(double r)
Definition: EcalClusterTools.h:319
EcalClusterToolsT::f22
static double f22(double r)
Definition: EcalClusterTools.h:320
noZS
Definition: EcalClusterLazyTools.h:304
EcalClusterToolsT::factorial
static double factorial(int n)
Definition: EcalClusterTools.h:355
HLTTauDQMOffline_cfi.denominator
denominator
Definition: HLTTauDQMOffline_cfi.py:180
CaloCellGeometry.h
Angle
Definition: Angle.h:17
EcalClusterToolsT::getNormedIX
static float getNormedIX(const DetId &id)
Definition: EcalClusterTools.h:1290
EcalClusterToolsT::EcalClusterEnergyDeposition
Definition: EcalClusterTools.h:290
EcalClusterToolsT::recHitEnergy
static float recHitEnergy(DetId id, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:416
DetId.h
HLT_FULL_cff.dEta
dEta
Definition: HLT_FULL_cff.py:13694
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
EcalClusterToolsT::e2x5Left
static float e2x5Left(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:582
SuperCluster.h
Exception
Definition: hltDiff.cc:245
CaloGeometry.h
EcalClusterToolsT::e1x5
static float e1x5(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:625
EcalClusterToolsT::zernike42
static double zernike42(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0=6.6, bool logW=true, float w0=4.7)
Definition: EcalClusterTools.h:1139
CaloRectangle.h
EcalClusterToolsT::f53
static double f53(double r)
Definition: EcalClusterTools.h:327
EcalClusterToolsT::f51
static double f51(double r)
Definition: EcalClusterTools.h:326
EcalClusterToolsT::energyBasketFractionEta
static std::vector< float > energyBasketFractionEta(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
Definition: EcalClusterTools.h:689
EcalClusterToolsT::e5x5
static float e5x5(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:539
CaloSubdetectorGeometry
Definition: CaloSubdetectorGeometry.h:22
VarParsing.mult
mult
Definition: VarParsing.py:658
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
Cluster2ndMoments::sMin
float sMin
Definition: EcalClusterTools.h:69
EcalClusterToolsT::scLocalCovariances
static std::array< float, 3 > scLocalCovariances(const reco::SuperCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, float w0=4.7)
Definition: EcalClusterTools.h:1380
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
EgammaLocalCovParamDefaults.h
EcalClusterToolsT::lat
static std::vector< float > lat(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW=true, float w0=4.7)
Definition: EcalClusterTools.h:795
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
EBDetId::iphi
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
mps_fire.result
result
Definition: mps_fire.py:311
EcalRecHit::kDiWeird
Definition: EcalRecHit.h:36
EcalClusterToolsT::eLeft
static float eLeft(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:657
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
EcalClusterToolsT::f42
static double f42(double r)
Definition: EcalClusterTools.h:324
ClusterTask_cfi.energyThreshold
energyThreshold
Definition: ClusterTask_cfi.py:3
EcalClusterToolsT::computeWeight
static float computeWeight(float eRH, float energyTotal, int weightedPositionMethod)
HLT_FULL_cff.R0
R0
Definition: HLT_FULL_cff.py:8652
noZS::EcalClusterTools
EcalClusterToolsT< true > EcalClusterTools
Definition: EcalClusterTools.h:1923
EcalClusterToolsT::e1x3
static float e1x3(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology)
Definition: EcalClusterTools.h:641
Vector3D.h
Cluster2ndMoments
Definition: EcalClusterTools.h:66
Cluster2ndMoments::Cluster2ndMoments
Cluster2ndMoments()
Definition: EcalClusterTools.h:72
EcalClusterToolsT::meanClusterPosition
static math::XYZVector meanClusterPosition(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, const CaloGeometry *geometry)
Definition: EcalClusterTools.h:862
reco::CaloCluster::energy
double energy() const
cluster energy
Definition: CaloCluster.h:149
weight
Definition: weight.py:1
EcalClusterToolsT::roundnessSelectedBarrelRecHits
static std::vector< float > roundnessSelectedBarrelRecHits(const std::vector< std::pair< const EcalRecHit *, float >> &rhVector, int weightedPositionMethod=0)
Definition: EcalClusterTools.h:1704
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37