CMS 3D CMS Logo

HGCalClusterT.h
Go to the documentation of this file.
1 #ifndef DataFormats_L1Trigger_HGCalClusterT_h
2 #define DataFormats_L1Trigger_HGCalClusterT_h
3 
4 /* CMSSW */
13 
14 /* ROOT */
15 #include "Math/Vector3D.h"
16 
17 #include <unordered_map>
18 
19 namespace l1t {
20  template <class C>
21  class HGCalClusterT : public L1Candidate {
22  public:
23  typedef typename std::unordered_map<uint32_t, edm::Ptr<C>>::const_iterator const_iterator;
24 
25  public:
27  HGCalClusterT(const LorentzVector p4, int pt = 0, int eta = 0, int phi = 0)
28  : L1Candidate(p4, pt, eta, phi),
29  valid_(true),
30  detId_(0),
31  centre_(0, 0, 0),
32  centreProj_(0., 0., 0.),
33  mipPt_(0),
34  seedMipPt_(0),
35  sumPt_(0) {}
36 
37  HGCalClusterT(const edm::Ptr<C>& c, float fraction = 1.)
38  : valid_(true),
39  detId_(c->detId()),
40  centre_(0., 0., 0.),
41  centreProj_(0., 0., 0.),
42  mipPt_(0.),
43  seedMipPt_(0.),
44  sumPt_(0.) {
45  addConstituent(c, true, fraction);
46  }
47 
48  ~HGCalClusterT() override{};
49 
50  const std::unordered_map<uint32_t, edm::Ptr<C>>& constituents() const { return constituents_; }
51  const_iterator constituents_begin() const { return constituents_.begin(); }
52  const_iterator constituents_end() const { return constituents_.end(); }
53  unsigned size() const { return constituents_.size(); }
54 
55  void addConstituent(const edm::Ptr<C>& c, bool updateCentre = true, float fraction = 1.) {
56  double cMipt = c->mipPt() * fraction;
57 
58  if (constituents_.empty()) {
59  detId_ = DetId(c->detId());
60  seedMipPt_ = cMipt;
61  /* if the centre will not be dynamically calculated
62  the seed centre is considere as cluster centre */
63  if (!updateCentre) {
64  centre_ = c->position();
65  }
66  }
67  updateP4AndPosition(c, updateCentre, fraction);
68 
69  constituents_.emplace(c->detId(), c);
70  constituentsFraction_.emplace(c->detId(), fraction);
71  }
72 
73  void removeConstituent(const edm::Ptr<C>& c, bool updateCentre = true) {
74  /* remove the pointer to c from the std::vector */
75  double fraction = 0;
76  const auto& constituent_itr = constituents_.find(c->detId());
77  const auto& fraction_itr = constituentsFraction_.find(c->detId());
78  if (constituent_itr != constituents_.end()) {
79  // remove constituent and get its fraction in the cluster
80  fraction = fraction_itr->second;
81  constituents_.erase(constituent_itr);
82  constituentsFraction_.erase(fraction_itr);
83 
84  updateP4AndPosition(c, updateCentre, -fraction);
85  }
86  }
87 
88  bool valid() const { return valid_; }
89  void setValid(bool valid) { valid_ = valid; }
90 
91  double mipPt() const { return mipPt_; }
92  double seedMipPt() const { return seedMipPt_; }
93  uint32_t detId() const { return detId_.rawId(); }
94  void setDetId(uint32_t id) { detId_ = id; }
95  void setPt(double pt) { setP4(math::PtEtaPhiMLorentzVector(pt, eta(), phi(), mass())); }
96  double sumPt() const { return sumPt_; }
97  /* distance in 'cm' */
98  double distance(const l1t::HGCalTriggerCell& tc) const { return (tc.position() - centre_).mag(); }
99 
100  const GlobalPoint& position() const { return centre_; }
101  const GlobalPoint& centre() const { return centre_; }
102  const GlobalPoint& centreProj() const { return centreProj_; }
103 
104  // FIXME: will need to fix places where the shapes are directly accessed
105  // Right now keep shapes() getter as non-const
106  ClusterShapes& shapes() { return shapes_; }
107  double hOverE() const {
108  double pt_em = 0.;
109  double pt_had = 0.;
110  double hOe = 0.;
111 
112  for (const auto& id_constituent : constituents()) {
113  DetId id(id_constituent.first);
114  auto id_fraction = constituentsFraction_.find(id_constituent.first);
115  double fraction = (id_fraction != constituentsFraction_.end() ? id_fraction->second : 1.);
116  if ((id.det() == DetId::Forward && id.subdetId() == HGCEE) || (id.det() == DetId::HGCalEE) ||
117  (id.det() == DetId::HGCalTrigger &&
119  pt_em += id_constituent.second->pt() * fraction;
120  } else if ((id.det() == DetId::Forward && id.subdetId() == HGCHEF) ||
121  (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) || (id.det() == DetId::HGCalHSi) ||
122  (id.det() == DetId::HGCalHSc) ||
123  (id.det() == DetId::HGCalTrigger &&
125  pt_had += id_constituent.second->pt() * fraction;
126  }
127  }
128  if (pt_em > 0)
129  hOe = pt_had / pt_em;
130  else
131  hOe = -1.;
132  return hOe;
133  }
134 
135  uint32_t subdetId() const { return detId_.subdetId(); }
136 
137  //shower shape
138 
139  int showerLength() const { return showerLength_; }
140  int coreShowerLength() const { return coreShowerLength_; }
141  int firstLayer() const { return firstLayer_; }
142  int maxLayer() const { return maxLayer_; }
143  float eMax() const { return eMax_; }
144  float sigmaEtaEtaMax() const { return sigmaEtaEtaMax_; }
145  float sigmaPhiPhiMax() const { return sigmaPhiPhiMax_; }
146  float sigmaEtaEtaTot() const { return sigmaEtaEtaTot_; }
147  float sigmaPhiPhiTot() const { return sigmaPhiPhiTot_; }
148  float sigmaZZ() const { return sigmaZZ_; }
149  float sigmaRRTot() const { return sigmaRRTot_; }
150  float sigmaRRMax() const { return sigmaRRMax_; }
151  float sigmaRRMean() const { return sigmaRRMean_; }
152  float zBarycenter() const { return zBarycenter_; }
153  float layer10percent() const { return layer10percent_; }
154  float layer50percent() const { return layer50percent_; }
155  float layer90percent() const { return layer90percent_; }
158 
163  void eMax(float eMax) { eMax_ = eMax; }
171  void sigmaZZ(float sigmaZZ) { sigmaZZ_ = sigmaZZ; }
178 
179  /* operators */
180  bool operator<(const HGCalClusterT<C>& cl) const { return mipPt() < cl.mipPt(); }
181  bool operator>(const HGCalClusterT<C>& cl) const { return cl < *this; }
182  bool operator<=(const HGCalClusterT<C>& cl) const { return !(cl > *this); }
183  bool operator>=(const HGCalClusterT<C>& cl) const { return !(cl < *this); }
184 
185  private:
186  bool valid_;
188 
189  std::unordered_map<uint32_t, edm::Ptr<C>> constituents_;
190  std::unordered_map<uint32_t, double> constituentsFraction_;
191 
193  GlobalPoint centreProj_; // centre projected onto the first HGCal layer
194 
195  double mipPt_;
196  double seedMipPt_;
197  double sumPt_;
198 
199  //shower shape
200 
201  int showerLength_ = 0;
203  int firstLayer_ = 0;
204  int maxLayer_ = 0;
205  float eMax_ = 0.;
206  float sigmaEtaEtaMax_ = 0.;
207  float sigmaPhiPhiMax_ = 0.;
208  float sigmaRRMax_ = 0.;
209  float sigmaEtaEtaTot_ = 0.;
210  float sigmaPhiPhiTot_ = 0.;
211  float sigmaRRTot_ = 0.;
212  float sigmaRRMean_ = 0.;
213  float sigmaZZ_ = 0.;
214  float zBarycenter_ = 0.;
215  float layer10percent_ = 0.;
216  float layer50percent_ = 0.;
217  float layer90percent_ = 0.;
220 
222 
223  void updateP4AndPosition(const edm::Ptr<C>& c, bool updateCentre = true, float fraction = 1.) {
224  double cMipt = c->mipPt() * fraction;
225  double cPt = c->pt() * fraction;
226  /* update cluster positions (IF requested) */
227  if (updateCentre) {
228  Basic3DVector<float> constituentCentre(c->position());
229  Basic3DVector<float> clusterCentre(centre_);
230 
231  clusterCentre = clusterCentre * mipPt_ + constituentCentre * cMipt;
232  if ((mipPt_ + cMipt) > 0) {
233  clusterCentre /= (mipPt_ + cMipt);
234  }
235  centre_ = GlobalPoint(clusterCentre);
236 
237  if (clusterCentre.z() != 0) {
238  centreProj_ = GlobalPoint(clusterCentre / std::abs(clusterCentre.z()));
239  }
240  }
241 
242  /* update cluster energies */
243  mipPt_ += cMipt;
244  sumPt_ += cPt;
245  int updatedPt = hwPt() + (int)(c->hwPt() * fraction);
246  setHwPt(updatedPt);
247 
248  math::PtEtaPhiMLorentzVector updatedP4(p4());
249  updatedP4 += (c->p4() * fraction);
250  setP4(updatedP4);
251  }
252  };
253 
254 } // namespace l1t
255 
256 #endif
float layer10percent() const
void addConstituent(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
Definition: HGCalClusterT.h:55
double distance(const l1t::HGCalTriggerCell &tc) const
Definition: HGCalClusterT.h:98
void layer90percent(float layer90percent)
int maxLayer() const
void sigmaRRTot(float sigmaRRTot)
const_iterator constituents_end() const
Definition: HGCalClusterT.h:52
void sigmaEtaEtaMax(float sigmaEtaEtaMax)
double eta() const final
momentum pseudorapidity
void showerLength(int showerLength)
HGCalTriggerSubdetector subdet() const
get the subdetector
void firstLayer(int firstLayer)
void setDetId(uint32_t id)
Definition: HGCalClusterT.h:94
float sigmaPhiPhiTot() const
float triggerCells90percent() const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void zBarycenter(float zBarycenter)
unsigned size() const
Definition: HGCalClusterT.h:53
HGCalClusterT(const LorentzVector p4, int pt=0, int eta=0, int phi=0)
Definition: HGCalClusterT.h:27
double hOverE() const
float sigmaEtaEtaMax() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
uint32_t subdetId() const
float sigmaRRTot() const
ClusterShapes shapes_
double pt() const final
transverse momentum
double sumPt() const
Definition: HGCalClusterT.h:96
std::unordered_map< uint32_t, edm::Ptr< C > > constituents_
ClusterShapes & shapes()
delete x;
Definition: CaloConfig.h:22
const GlobalPoint & centreProj() const
float sigmaZZ() const
void setPt(double pt)
Definition: HGCalClusterT.h:95
const GlobalPoint & position() const
std::unordered_map< uint32_t, edm::Ptr< C > >::const_iterator const_iterator
Definition: HGCalClusterT.h:23
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void triggerCells90percent(float triggerCells90percent)
void triggerCells67percent(float triggerCells67percent)
void sigmaRRMax(float sigmaRRMax)
const_iterator constituents_begin() const
Definition: HGCalClusterT.h:51
void removeConstituent(const edm::Ptr< C > &c, bool updateCentre=true)
Definition: HGCalClusterT.h:73
int showerLength() const
T z() const
Cartesian z coordinate.
~HGCalClusterT() override
Definition: HGCalClusterT.h:48
void updateP4AndPosition(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
bool operator>(const HGCalClusterT< C > &cl) const
void sigmaPhiPhiMax(float sigmaPhiPhiMax)
GlobalPoint centreProj_
double seedMipPt() const
Definition: HGCalClusterT.h:92
bool valid() const
Definition: HGCalClusterT.h:88
std::unordered_map< uint32_t, double > constituentsFraction_
const GlobalPoint & centre() const
float layer90percent() const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double mipPt() const
Definition: HGCalClusterT.h:91
void sigmaZZ(float sigmaZZ)
const LorentzVector & p4() const final
four-momentum Lorentz vector
void sigmaPhiPhiTot(float sigmaPhiPhiTot)
uint32_t detId() const
Definition: HGCalClusterT.h:93
void eMax(float eMax)
GlobalPoint centre_
float sigmaEtaEtaTot() const
float eMax() const
Definition: DetId.h:17
void layer10percent(float layer10percent)
HGCalClusterT(const edm::Ptr< C > &c, float fraction=1.)
Definition: HGCalClusterT.h:37
int coreShowerLength() const
float layer50percent() const
int hwPt() const
Definition: L1Candidate.h:35
int firstLayer() const
float sigmaRRMean() const
float triggerCells67percent() const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void layer50percent(float layer50percent)
void setHwPt(int pt)
Definition: L1Candidate.h:28
const GlobalPoint & position() const
float sigmaRRMax() const
void maxLayer(int maxLayer)
bool operator>=(const HGCalClusterT< C > &cl) const
float zBarycenter() const
float sigmaPhiPhiMax() const
void sigmaEtaEtaTot(float sigmaEtaEtaTot)
void sigmaRRMean(float sigmaRRMean)
double phi() const final
momentum azimuthal angle
void setValid(bool valid)
Definition: HGCalClusterT.h:89
void setP4(const LorentzVector &p4) final
set 4-momentum
void coreShowerLength(int coreShowerLength)
double mass() const final
mass
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50