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 
10 #include "Math/Vector3D.h"
11 
12 
13 namespace l1t
14 {
15  template <class C> class HGCalClusterT : public L1Candidate
16  {
17 
18  public:
20 
21  public:
24  int pt=0,
25  int eta=0,
26  int phi=0
27  )
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 
37  valid_(true),
38  detId_( c->detId() ),
39  centre_(0., 0., 0.),
40  centreProj_(0., 0., 0.),
41  mipPt_(0.),
42  seedMipPt_(0.)
43  {
44  addConstituent(c);
45  }
46 
47  ~HGCalClusterT() override {};
48 
49  const edm::PtrVector<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 )
55  {
56  if( constituents_.empty() )
57  {
58  detId_ = HGCalDetId(c->detId());
59  seedMipPt_ = c->mipPt();
60  }
61 
62  /* update cluster positions */
63  Basic3DVector<float> constituentCentre( c->position() );
64  Basic3DVector<float> clusterCentre( centre_ );
65 
66  clusterCentre = clusterCentre*mipPt_ + constituentCentre*c->mipPt();
67  if( mipPt_ + c->mipPt()!=0 )
68  {
69  clusterCentre /= ( mipPt_ + c->mipPt() ) ;
70  }
71  centre_ = GlobalPoint( clusterCentre );
72 
73  if( clusterCentre.z()!=0 )
74  {
75  centreProj_= GlobalPoint( clusterCentre / clusterCentre.z() );
76  }
77  /* update cluster energies */
78  mipPt_ += c->mipPt();
79 
80  int updatedPt = hwPt() + c->hwPt();
81  setHwPt(updatedPt);
82 
83  math::PtEtaPhiMLorentzVector updatedP4 ( p4() );
84  updatedP4 += c->p4();
85  setP4( updatedP4 );
86 
88 
89  }
90 
91  bool valid() const { return valid_;}
92  void setValid(bool valid) { valid_ = valid;}
93 
94  double mipPt() const { return mipPt_; }
95  double seedMipPt() const { return seedMipPt_; }
96  uint32_t detId() const { return detId_.rawId(); }
97 
98  /* distance in 'cm' */
99  double distance( const l1t::HGCalTriggerCell &tc ) const
100  {
101  return ( tc.position() - centre_ ).mag();
102  }
103 
104  const GlobalPoint& position() const { return centre_; }
105  const GlobalPoint& centre() const { return centre_; }
106  const GlobalPoint& centreProj() const { return centreProj_; }
107 
108  // FIXME: will need to fix places where the shapes are directly accessed
109  // Right now keep shapes() getter as non-const
111  double hOverE() const
112  {
113  double pt_em = 0.;
114  double pt_had = 0.;
115  double hOe = 0.;
116 
117  for(const auto& constituent : constituents())
118  {
119  switch( constituent->subdetId() )
120  {
121  case HGCEE:
122  pt_em += constituent->pt();
123  break;
124  case HGCHEF:
125  pt_had += constituent->pt();
126  break;
127  case HGCHEB:
128  pt_had += constituent->pt();
129  break;
130  default:
131  break;
132  }
133  }
134  if(pt_em>0) hOe = pt_had / pt_em ;
135  else hOe = -1.;
136  return hOe;
137  }
138 
139  uint32_t subdetId() const {return detId_.subdetId();}
140  uint32_t layer() const {return detId_.layer();}
141  int32_t zside() const {return detId_.zside();}
142 
143 
144  /* operators */
145  bool operator<(const HGCalClusterT<C>& cl) const {return mipPt() < cl.mipPt();}
146  bool operator>(const HGCalClusterT<C>& cl) const { return cl<*this; }
147  bool operator<=(const HGCalClusterT<C>& cl) const { return !(cl>*this); }
148  bool operator>=(const HGCalClusterT<C>& cl) const { return !(cl<*this); }
149 
150  private:
151 
152  bool valid_;
156  GlobalPoint centreProj_; // centre projected onto the first HGCal layer
157 
158  double mipPt_;
159  double seedMipPt_;
160 
162 
163  };
164 
165 }
166 
167 #endif
double distance(const l1t::HGCalTriggerCell &tc) const
Definition: HGCalClusterT.h:99
const_iterator constituents_end() const
Definition: HGCalClusterT.h:51
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:74
double eta() const final
momentum pseudorapidity
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
unsigned size() const
Definition: HGCalClusterT.h:52
HGCalClusterT(const LorentzVector p4, int pt=0, int eta=0, int phi=0)
Definition: HGCalClusterT.h:23
double hOverE() const
bool empty() const
Is the RefVector empty.
Definition: PtrVectorBase.h:71
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
uint32_t subdetId() const
ClusterShapes shapes_
double pt() const final
transverse momentum
ClusterShapes & shapes()
delete x;
Definition: CaloConfig.h:22
const GlobalPoint & centreProj() const
const GlobalPoint & position() const
const_iterator begin() const
Definition: PtrVector.h:129
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
const_iterator constituents_begin() const
Definition: HGCalClusterT.h:50
int zside() const
get the z-side of the cell (1/-1)
Definition: HGCalDetId.h:51
T z() const
Cartesian z coordinate.
~HGCalClusterT() override
Definition: HGCalClusterT.h:47
bool operator>(const HGCalClusterT< C > &cl) const
GlobalPoint centreProj_
double seedMipPt() const
Definition: HGCalClusterT.h:95
bool valid() const
Definition: HGCalClusterT.h:91
const GlobalPoint & centre() const
const_iterator end() const
Definition: PtrVector.h:134
double mipPt() const
Definition: HGCalClusterT.h:94
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
uint32_t detId() const
Definition: HGCalClusterT.h:96
HGCalClusterT(const edm::Ptr< C > &c)
Definition: HGCalClusterT.h:36
edm::PtrVector< C >::const_iterator const_iterator
Definition: HGCalClusterT.h:19
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
GlobalPoint centre_
uint32_t layer() const
void addConstituent(const edm::Ptr< C > &c)
Definition: HGCalClusterT.h:54
int hwPt() const
Definition: L1Candidate.h:48
const edm::PtrVector< C > & constituents() const
Definition: HGCalClusterT.h:49
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void setHwPt(int pt)
Definition: L1Candidate.h:41
const GlobalPoint & position() const
int32_t zside() const
edm::PtrVector< C > constituents_
bool operator>=(const HGCalClusterT< C > &cl) const
double phi() const final
momentum azimuthal angle
void setValid(bool valid)
Definition: HGCalClusterT.h:92
void setP4(const LorentzVector &p4) final
set 4-momentum
int layer() const
get the layer #
Definition: HGCalDetId.h:48