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 */
12 
13 /* ROOT */
14 #include "Math/Vector3D.h"
15 
16 #include <unordered_map>
17 
18 namespace l1t {
19  template <class C>
20  class HGCalClusterT : public L1Candidate {
21  public:
22  typedef typename std::unordered_map<uint32_t, edm::Ptr<C>>::const_iterator const_iterator;
23 
24  public:
26  HGCalClusterT(const LorentzVector p4, int pt = 0, int eta = 0, int phi = 0)
27  : L1Candidate(p4, pt, eta, phi),
28  valid_(true),
29  detId_(0),
30  centre_(0, 0, 0),
31  centreProj_(0., 0., 0.),
32  mipPt_(0),
33  seedMipPt_(0),
34  sumPt_(0) {}
35 
36  HGCalClusterT(const edm::Ptr<C>& c, float fraction = 1.)
37  : valid_(true),
38  detId_(c->detId()),
39  centre_(0., 0., 0.),
40  centreProj_(0., 0., 0.),
41  mipPt_(0.),
42  seedMipPt_(0.),
43  sumPt_(0.) {
44  addConstituent(c, true, fraction);
45  }
46 
47  ~HGCalClusterT() override {}
48 
49  const std::unordered_map<uint32_t, edm::Ptr<C>>& constituents() const { return constituents_; }
50  const_iterator constituents_begin() const { return constituents_.begin(); }
51  const_iterator constituents_end() const { return constituents_.end(); }
52  unsigned size() const { return constituents_.size(); }
53 
54  void addConstituent(const edm::Ptr<C>& c, bool updateCentre = true, float fraction = 1.) {
55  double cMipt = c->mipPt() * fraction;
56 
57  if (constituents_.empty()) {
58  detId_ = DetId(c->detId());
59  seedMipPt_ = cMipt;
60  /* if the centre will not be dynamically calculated
61  the seed centre is considere as cluster centre */
62  if (!updateCentre) {
63  centre_ = c->position();
64  }
65  }
66  updateP4AndPosition(c, updateCentre, fraction);
67 
68  constituents_.emplace(c->detId(), c);
69  constituentsFraction_.emplace(c->detId(), fraction);
70  }
71 
72  void removeConstituent(const edm::Ptr<C>& c, bool updateCentre = true) {
73  /* remove the pointer to c from the std::vector */
74  double fraction = 0;
75  const auto& constituent_itr = constituents_.find(c->detId());
76  const auto& fraction_itr = constituentsFraction_.find(c->detId());
77  if (constituent_itr != constituents_.end()) {
78  // remove constituent and get its fraction in the cluster
79  fraction = fraction_itr->second;
80  constituents_.erase(constituent_itr);
81  constituentsFraction_.erase(fraction_itr);
82 
83  updateP4AndPosition(c, updateCentre, -fraction);
84  }
85  }
86 
87  bool valid() const { return valid_; }
88  void setValid(bool valid) { valid_ = valid; }
89 
90  double mipPt() const { return mipPt_; }
91  double seedMipPt() const { return seedMipPt_; }
92  uint32_t detId() const { return detId_.rawId(); }
93  void setDetId(uint32_t id) { detId_ = id; }
94  void setPt(double pt) { setP4(math::PtEtaPhiMLorentzVector(pt, eta(), phi(), mass())); }
95  double sumPt() const { return sumPt_; }
96  /* distance in 'cm' */
97  double distance(const l1t::HGCalTriggerCell& tc) const { return (tc.position() - centre_).mag(); }
98 
99  const GlobalPoint& position() const { return centre_; }
100  const GlobalPoint& centre() const { return centre_; }
101  const GlobalPoint& centreProj() const { return centreProj_; }
102 
103  double hOverE() const {
104  double pt_em = 0.;
105  double pt_had = 0.;
106  double hOe = 0.;
107 
108  for (const auto& id_constituent : constituents()) {
109  DetId id(id_constituent.first);
110  auto id_fraction = constituentsFraction_.find(id_constituent.first);
111  double fraction = (id_fraction != constituentsFraction_.end() ? id_fraction->second : 1.);
112  if ((id.det() == DetId::HGCalEE) ||
113  (id.det() == DetId::HGCalTrigger &&
115  (id.det() == DetId::Forward && id.subdetId() == ForwardSubdetector::HFNose && HFNoseDetId(id).isEE()) ||
116  (id.det() == DetId::HGCalTrigger &&
118  HFNoseTriggerDetId(id).isEE())) {
119  pt_em += id_constituent.second->pt() * fraction;
120  } else if ((id.det() == DetId::HGCalHSi) || (id.det() == DetId::HGCalHSc) ||
121  (id.det() == DetId::HGCalTrigger &&
123  (id.det() == DetId::Forward && id.subdetId() == ForwardSubdetector::HFNose &&
124  HFNoseDetId(id).isHE()) ||
125  (id.det() == DetId::HGCalTrigger &&
127  HFNoseTriggerDetId(id).isHSilicon())) {
128  pt_had += id_constituent.second->pt() * fraction;
129  }
130  }
131  if (pt_em > 0)
132  hOe = pt_had / pt_em;
133  else
134  hOe = -1.;
135  return hOe;
136  }
137 
138  uint32_t subdetId() const { return detId_.subdetId(); }
139 
140  //shower shape
141 
142  int showerLength() const { return showerLength_; }
143  int coreShowerLength() const { return coreShowerLength_; }
144  int firstLayer() const { return firstLayer_; }
145  int maxLayer() const { return maxLayer_; }
146  float eMax() const { return eMax_; }
147  float sigmaEtaEtaMax() const { return sigmaEtaEtaMax_; }
148  float sigmaPhiPhiMax() const { return sigmaPhiPhiMax_; }
149  float sigmaEtaEtaTot() const { return sigmaEtaEtaTot_; }
150  float sigmaPhiPhiTot() const { return sigmaPhiPhiTot_; }
151  float sigmaZZ() const { return sigmaZZ_; }
152  float sigmaRRTot() const { return sigmaRRTot_; }
153  float varRR() const { return varRR_; }
154  float varZZ() const { return varZZ_; }
155  float varEtaEta() const { return varEtaEta_; }
156  float varPhiPhi() const { return varPhiPhi_; }
157  float sigmaRRMax() const { return sigmaRRMax_; }
158  float sigmaRRMean() const { return sigmaRRMean_; }
159  float zBarycenter() const { return zBarycenter_; }
160  float layer10percent() const { return layer10percent_; }
161  float layer50percent() const { return layer50percent_; }
162  float layer90percent() const { return layer90percent_; }
165  float first1layers() const { return first1layers_; }
166  float first3layers() const { return first3layers_; }
167  float first5layers() const { return first5layers_; }
168  float firstHcal1layers() const { return firstHcal1layers_; }
169  float firstHcal3layers() const { return firstHcal3layers_; }
170  float firstHcal5layers() const { return firstHcal5layers_; }
171  float last1layers() const { return last1layers_; }
172  float last3layers() const { return last3layers_; }
173  float last5layers() const { return last5layers_; }
174  float emax1layers() const { return emax1layers_; }
175  float emax3layers() const { return emax3layers_; }
176  float emax5layers() const { return emax5layers_; }
177  float eot() const { return eot_; }
178  int ebm0() const { return ebm0_; }
179  int ebm1() const { return ebm1_; }
180  int hbm() const { return hbm_; }
181 
186  void setEMax(float eMax) { eMax_ = eMax; }
193  void setVarRR(float varRR) { varRR_ = varRR; }
194  void setVarZZ(float varZZ) { varZZ_ = varZZ; }
198  void setSigmaZZ(float sigmaZZ) { sigmaZZ_ = sigmaZZ; }
217  void setEot(float eot) { eot_ = eot; }
218  void setEbm0(int ebm0) { ebm0_ = ebm0; }
219  void setEbm1(int ebm1) { ebm1_ = ebm1; }
220  void setHbm(int hbm) { hbm_ = hbm; }
221 
222  /* operators */
223  bool operator<(const HGCalClusterT<C>& cl) const { return mipPt() < cl.mipPt(); }
224  bool operator>(const HGCalClusterT<C>& cl) const { return cl < *this; }
225  bool operator<=(const HGCalClusterT<C>& cl) const { return !(cl > *this); }
226  bool operator>=(const HGCalClusterT<C>& cl) const { return !(cl < *this); }
227 
228  private:
229  bool valid_ = false;
231 
232  std::unordered_map<uint32_t, edm::Ptr<C>> constituents_;
233  std::unordered_map<uint32_t, double> constituentsFraction_;
234 
236  GlobalPoint centreProj_; // centre projected onto the first HGCal layer
237 
238  double mipPt_ = 0.;
239  double seedMipPt_ = 0.;
240  double sumPt_ = 0.;
241 
242  //shower shape
243 
244  int showerLength_ = 0;
246  int firstLayer_ = 0;
247  int maxLayer_ = 0;
248  float eMax_ = 0.;
249  float sigmaEtaEtaMax_ = 0.;
250  float sigmaPhiPhiMax_ = 0.;
251  float sigmaRRMax_ = 0.;
252  float sigmaEtaEtaTot_ = 0.;
253  float sigmaPhiPhiTot_ = 0.;
254  float sigmaRRTot_ = 0.;
255  float varRR_ = 0.;
256  float varZZ_ = 0.;
257  float varEtaEta_ = 0.;
258  float varPhiPhi_ = 0.;
259  float sigmaRRMean_ = 0.;
260  float sigmaZZ_ = 0.;
261  float zBarycenter_ = 0.;
262  float layer10percent_ = 0.;
263  float layer50percent_ = 0.;
264  float layer90percent_ = 0.;
267  float first1layers_ = 0.;
268  float first3layers_ = 0.;
269  float first5layers_ = 0.;
270  float firstHcal1layers_ = 0.;
271  float firstHcal3layers_ = 0.;
272  float firstHcal5layers_ = 0.;
273  float last1layers_ = 0.;
274  float last3layers_ = 0.;
275  float last5layers_ = 0.;
276  float emax1layers_ = 0.;
277  float emax3layers_ = 0.;
278  float emax5layers_ = 0.;
279  float eot_ = 0.;
280  int ebm0_ = 0;
281  int ebm1_ = 0;
282  int hbm_ = 0;
283 
284  void updateP4AndPosition(const edm::Ptr<C>& c, bool updateCentre = true, float fraction = 1.) {
285  double cMipt = c->mipPt() * fraction;
286  double cPt = c->pt() * fraction;
287  /* update cluster positions (IF requested) */
288  if (updateCentre) {
289  Basic3DVector<float> constituentCentre(c->position());
290  Basic3DVector<float> clusterCentre(centre_);
291 
292  clusterCentre = clusterCentre * mipPt_ + constituentCentre * cMipt;
293  if ((mipPt_ + cMipt) > 0) {
294  clusterCentre /= (mipPt_ + cMipt);
295  }
296  centre_ = GlobalPoint(clusterCentre);
297 
298  if (clusterCentre.z() != 0) {
299  centreProj_ = GlobalPoint(clusterCentre / std::abs(clusterCentre.z()));
300  }
301  }
302 
303  /* update cluster energies */
304  mipPt_ += cMipt;
305  sumPt_ += cPt;
306  int updatedPt = hwPt() + (int)(c->hwPt() * fraction);
307  setHwPt(updatedPt);
308 
309  math::PtEtaPhiMLorentzVector updatedP4(p4());
310  updatedP4 += (c->p4() * fraction);
311  setP4(updatedP4);
312  }
313  };
314 
315 } // namespace l1t
316 
317 #endif
float firstHcal3layers() const
float sigmaEtaEtaTot() const
void setFirst1layers(float first1layers)
void addConstituent(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
Definition: HGCalClusterT.h:54
float eMax() const
void setMaxLayer(int maxLayer)
float firstHcal5layers() const
float last5layers() const
float eot() const
float layer50percent() const
double pt() const final
transverse momentum
bool operator>=(const HGCalClusterT< C > &cl) const
void setDetId(uint32_t id)
Definition: HGCalClusterT.h:93
void setSigmaPhiPhiTot(float sigmaPhiPhiTot)
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:49
void setHbm(int hbm)
HGCalClusterT(const LorentzVector p4, int pt=0, int eta=0, int phi=0)
Definition: HGCalClusterT.h:26
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
int coreShowerLength() const
void setSigmaZZ(float sigmaZZ)
void setVarEtaEta(float varEtaEta)
bool isHE(int etabin, int depth)
uint32_t detId() const
Definition: HGCalClusterT.h:92
void setSigmaRRMean(float sigmaRRMean)
uint32_t subdetId() const
std::unordered_map< uint32_t, edm::Ptr< C > > constituents_
delete x;
Definition: CaloConfig.h:22
void setFirst5layers(float first5layers)
void setSigmaRRTot(float sigmaRRTot)
const_iterator constituents_begin() const
Definition: HGCalClusterT.h:50
void setPt(double pt)
Definition: HGCalClusterT.h:94
float emax5layers() const
float triggerCells90percent() const
void setFirstHcal5layers(float firstHcal5layers)
const LorentzVector & p4() const final
four-momentum Lorentz vector
void setCoreShowerLength(int coreShowerLength)
void setEmax5layers(float emax5layers)
float sigmaPhiPhiTot() const
void setLayer10percent(float layer10percent)
std::unordered_map< uint32_t, edm::Ptr< C > >::const_iterator const_iterator
Definition: HGCalClusterT.h:22
const GlobalPoint & centre() const
float varEtaEta() const
float varRR() const
bool valid() const
Definition: HGCalClusterT.h:87
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void setEbm1(int ebm1)
bool operator>(const HGCalClusterT< C > &cl) const
bool isEE() const
consistency check : no bits left => no overhead
const GlobalPoint & position() const
Definition: HGCalClusterT.h:99
void setEbm0(int ebm0)
void removeConstituent(const edm::Ptr< C > &c, bool updateCentre=true)
Definition: HGCalClusterT.h:72
float sigmaPhiPhiMax() const
int firstLayer() const
float zBarycenter() const
~HGCalClusterT() override
Definition: HGCalClusterT.h:47
void updateP4AndPosition(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
float varPhiPhi() const
float first5layers() const
GlobalPoint centreProj_
float last3layers() const
void setLast1layers(float last1layers)
float emax1layers() const
std::unordered_map< uint32_t, double > constituentsFraction_
HGCalTriggerSubdetector subdet() const
get the subdetector
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned size() const
Definition: HGCalClusterT.h:52
void setVarPhiPhi(float varPhiPhi)
double seedMipPt() const
Definition: HGCalClusterT.h:91
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
void setFirst3layers(float first3layers)
float sigmaRRMax() const
float first3layers() const
void setEMax(float eMax)
GlobalPoint centre_
void setLast5layers(float last5layers)
void setEot(float eot)
int showerLength() const
const GlobalPoint & position() const
void setEmax1layers(float emax1layers)
double distance(const l1t::HGCalTriggerCell &tc) const
Definition: HGCalClusterT.h:97
double hOverE() const
Definition: DetId.h:17
void setLast3layers(float last3layers)
int maxLayer() const
void setLayer50percent(float layer50percent)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
HGCalClusterT(const edm::Ptr< C > &c, float fraction=1.)
Definition: HGCalClusterT.h:36
T z() const
Cartesian z coordinate.
int hwPt() const
Definition: L1Candidate.h:35
float last1layers() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
float emax3layers() const
void setTriggerCells67percent(float triggerCells67percent)
float firstHcal1layers() const
void setFirstLayer(int firstLayer)
float sigmaRRTot() const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
void setTriggerCells90percent(float triggerCells90percent)
void setLayer90percent(float layer90percent)
bool isEE() const
consistency check : no bits left => no overhead
Definition: HFNoseDetId.h:105
float layer90percent() const
void setFirstHcal3layers(float firstHcal3layers)
double mass() const final
mass
float sigmaRRMean() const
void setSigmaEtaEtaTot(float sigmaEtaEtaTot)
void setHwPt(int pt)
Definition: L1Candidate.h:28
void setSigmaEtaEtaMax(float sigmaEtaEtaMax)
float sigmaZZ() const
void setSigmaRRMax(float sigmaRRMax)
float layer10percent() const
void setEmax3layers(float emax3layers)
float triggerCells67percent() const
void setFirstHcal1layers(float firstHcal1layers)
float sigmaEtaEtaMax() const
void setSigmaPhiPhiMax(float sigmaPhiPhiMax)
float first1layers() const
double phi() const final
momentum azimuthal angle
float varZZ() const
void setVarZZ(float varZZ)
void setVarRR(float varRR)
double mipPt() const
Definition: HGCalClusterT.h:90
void setValid(bool valid)
Definition: HGCalClusterT.h:88
const_iterator constituents_end() const
Definition: HGCalClusterT.h:51
const GlobalPoint & centreProj() const
void setP4(const LorentzVector &p4) final
set 4-momentum
double sumPt() const
Definition: HGCalClusterT.h:95
void setZBarycenter(float zBarycenter)
void setShowerLength(int showerLength)
double eta() const final
momentum pseudorapidity