CMS 3D CMS Logo

HGCalShowerShape.h
Go to the documentation of this file.
1 #ifndef __L1Trigger_L1THGCal_HGCalShowerShape_h__
2 #define __L1Trigger_L1THGCal_HGCalShowerShape_h__
3 #include <vector>
4 #include <functional>
5 #include <cmath>
11 
13 public:
15 
18 
20 
22 
23  int firstLayer(const l1t::HGCalMulticluster& c3d) const;
24  int lastLayer(const l1t::HGCalMulticluster& c3d) const;
25  int maxLayer(const l1t::HGCalMulticluster& c3d) const;
26  int showerLength(const l1t::HGCalMulticluster& c3d) const {
27  return lastLayer(c3d) - firstLayer(c3d) + 1;
28  } //in number of layers
29  // Maximum number of consecutive layers in the cluster
30  int coreShowerLength(const l1t::HGCalMulticluster& c3d, const HGCalTriggerGeometryBase& triggerGeometry) const;
31  float percentileLayer(const l1t::HGCalMulticluster& c3d,
32  const HGCalTriggerGeometryBase& triggerGeometry,
33  float quantile = 0.5) const;
34 
35  float percentileTriggerCells(const l1t::HGCalMulticluster& c3d, float quantile = 0.5) const;
36 
37  float eMax(const l1t::HGCalMulticluster& c3d) const;
38 
39  float meanZ(const l1t::HGCalMulticluster& c3d) const;
40  float sigmaZZ(const l1t::HGCalMulticluster& c3d) const;
41 
42  float sigmaEtaEtaTot(const l1t::HGCalMulticluster& c3d) const;
43  float sigmaEtaEtaTot(const l1t::HGCalCluster& c2d) const;
44  float sigmaEtaEtaMax(const l1t::HGCalMulticluster& c3d) const;
45 
46  float sigmaPhiPhiTot(const l1t::HGCalMulticluster& c3d) const;
47  float sigmaPhiPhiTot(const l1t::HGCalCluster& c2d) const;
48  float sigmaPhiPhiMax(const l1t::HGCalMulticluster& c3d) const;
49 
50  float sigmaRRTot(const l1t::HGCalMulticluster& c3d) const;
51  float sigmaRRTot(const l1t::HGCalCluster& c2d) const;
52  float sigmaRRMax(const l1t::HGCalMulticluster& c3d) const;
53  float sigmaRRMean(const l1t::HGCalMulticluster& c3d, float radius = 5.) const;
54 
56 
57 private:
58  float meanX(const std::vector<pair<float, float>>& energy_X_tc) const;
59  // Compute energy-weighted RMS of any variable X in the cluster
60  // Delta(a,b) functor as template argument. Default is (a-b)
61  template <typename Delta = std::minus<float>>
62  float sigmaXX(const std::vector<pair<float, float>>& energy_X_tc, const float X_cluster) const {
63  Delta delta;
64  float Etot = 0;
65  float deltaX2_sum = 0;
66  for (const auto& energy_X : energy_X_tc) {
67  deltaX2_sum += energy_X.first * pow(delta(energy_X.second, X_cluster), 2);
68  Etot += energy_X.first;
69  }
70  float X_MSE = 0;
71  if (Etot > 0)
72  X_MSE = deltaX2_sum / Etot;
73  float X_RMS = sqrt(X_MSE);
74  return X_RMS;
75  }
76  // Special case of delta for phi
77  template <class T>
78  struct DeltaPhi {
79  T operator()(const T& x, const T& y) const { return deltaPhi(x, y); }
80  };
81  float sigmaPhiPhi(const std::vector<pair<float, float>>& energy_phi_tc, const float phi_cluster) const {
82  return sigmaXX<DeltaPhi<float>>(energy_phi_tc, phi_cluster);
83  }
84  template <typename T, typename Tref>
85  bool pass(const T& obj, const Tref& ref) const {
86  bool pass_threshold = (obj.mipPt() > threshold_);
87  GlobalPoint proj(Basic3DVector<float>(obj.position()) / std::abs(obj.position().z()));
88  bool pass_distance = ((proj - ref.centreProj()).mag() < distance_);
89  return pass_threshold && pass_distance;
90  }
91 
93  double threshold_ = 0.;
94  double distance_ = 1.;
95 };
96 
97 #endif
float meanX(const std::vector< pair< float, float >> &energy_X_tc) const
constexpr int pow(int x)
Definition: conifer.h:24
HGCalTriggerTools triggerTools_
float percentileTriggerCells(const l1t::HGCalMulticluster &c3d, float quantile=0.5) const
int lastLayer(const l1t::HGCalMulticluster &c3d) const
float percentileLayer(const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry, float quantile=0.5) const
int showerLength(const l1t::HGCalMulticluster &c3d) const
float sigmaEtaEtaTot(const l1t::HGCalMulticluster &c3d) const
float sigmaZZ(const l1t::HGCalMulticluster &c3d) const
void setGeometry(const HGCalTriggerGeometryBase *const)
float eMax(const l1t::HGCalMulticluster &c3d) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T operator()(const T &x, const T &y) const
float sigmaPhiPhiMax(const l1t::HGCalMulticluster &c3d) const
T sqrt(T t)
Definition: SSEVec.h:19
void setGeometry(const HGCalTriggerGeometryBase *const geom)
void fillShapes(l1t::HGCalMulticluster &, const HGCalTriggerGeometryBase &) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int coreShowerLength(const l1t::HGCalMulticluster &c3d, const HGCalTriggerGeometryBase &triggerGeometry) const
math::XYZTLorentzVector LorentzVector
bool pass(const T &obj, const Tref &ref) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
int maxLayer(const l1t::HGCalMulticluster &c3d) const
float sigmaRRMax(const l1t::HGCalMulticluster &c3d) const
float sigmaRRTot(const l1t::HGCalMulticluster &c3d) const
float sigmaEtaEtaMax(const l1t::HGCalMulticluster &c3d) const
float sigmaPhiPhiTot(const l1t::HGCalMulticluster &c3d) const
long double T
int firstLayer(const l1t::HGCalMulticluster &c3d) const
float sigmaXX(const std::vector< pair< float, float >> &energy_X_tc, const float X_cluster) const
float meanZ(const l1t::HGCalMulticluster &c3d) const
float sigmaRRMean(const l1t::HGCalMulticluster &c3d, float radius=5.) const
float sigmaPhiPhi(const std::vector< pair< float, float >> &energy_phi_tc, const float phi_cluster) const