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