CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalMulticluster.h
Go to the documentation of this file.
1 #ifndef DataFormats_L1Trigger_HGCalMulticluster_h
2 #define DataFormats_L1Trigger_HGCalMulticluster_h
3 
8 #include <boost/iterator/transform_iterator.hpp>
9 #include <functional>
10 
11 namespace l1t {
12 
13  class HGCalMulticluster : public HGCalClusterT<l1t::HGCalCluster> {
14  public:
16  HGCalMulticluster(const LorentzVector p4, int pt = 0, int eta = 0, int phi = 0);
17 
19 
20  ~HGCalMulticluster() override;
21 
22  float hOverE() const {
23  // --- this below would be faster when reading old objects, as HoE will only be computed once,
24  // --- but it may not be allowed by CMS rules because of the const_cast
25  // --- and could potentially cause a data race
26  // if (!hOverEValid_) (const_cast<HGCalMulticluster*>(this))->saveHOverE();
27  // --- this below is safe in any case
29  }
30 
31  void saveHOverE() {
33  hOverEValid_ = true;
34  }
35 
36  enum EnergyInterpretation { EM = 0 };
37 
39 
41  return energy() * interpretationFraction(eInt);
42  }
43 
44  double iPt(const HGCalMulticluster::EnergyInterpretation eInt) const { return pt() * interpretationFraction(eInt); }
45 
47  return p4() * interpretationFraction(eInt);
48  }
49 
51  return math::PtEtaPhiMLorentzVector(pt() * interpretationFraction(eInt), eta(), phi(), 0.);
52  }
53 
54  private:
55  template <typename Iter>
56  struct KeyGetter : std::unary_function<typename Iter::value_type, typename Iter::value_type::first_type> {
57  const typename Iter::value_type::first_type& operator()(const typename Iter::value_type& p) const {
58  return p.first;
59  }
60  };
61 
62  template <typename Iter>
63  boost::transform_iterator<KeyGetter<Iter>, Iter> key_iterator(Iter itr) const {
64  return boost::make_transform_iterator<KeyGetter<Iter>, Iter>(itr, KeyGetter<Iter>());
65  }
66 
67  public:
68  typedef boost::transform_iterator<KeyGetter<std::map<EnergyInterpretation, double>::const_iterator>,
69  std::map<EnergyInterpretation, double>::const_iterator>
71 
72  std::pair<EnergyInterpretation_const_iterator, EnergyInterpretation_const_iterator> energyInterpretations() const {
73  return std::make_pair(key_iterator(energyInterpretationFractions_.cbegin()),
75  }
76 
79  }
80 
83  }
84 
86 
87  private:
89 
90  float hOverE_;
92  std::map<EnergyInterpretation, double> energyInterpretationFractions_;
93  };
94 
96 
97 } // namespace l1t
98 
99 #endif
math::XYZTLorentzVector iP4(const HGCalMulticluster::EnergyInterpretation eInt) const
EnergyInterpretation_const_iterator interpretations_end() const
boost::transform_iterator< KeyGetter< std::map< EnergyInterpretation, double >::const_iterator >, std::map< EnergyInterpretation, double >::const_iterator > EnergyInterpretation_const_iterator
const Iter::value_type::first_type & operator()(const typename Iter::value_type &p) const
double pt() const final
transverse momentum
double iEnergy(const HGCalMulticluster::EnergyInterpretation eInt) const
size_t size_type
Definition: Candidate.h:29
double hOverE() const
double iPt(const HGCalMulticluster::EnergyInterpretation eInt) const
std::map< EnergyInterpretation, double > energyInterpretationFractions_
const LorentzVector & p4() const final
four-momentum Lorentz vector
EnergyInterpretation_const_iterator interpretations_begin() const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double p() const final
magnitude of momentum vector
std::pair< EnergyInterpretation_const_iterator, EnergyInterpretation_const_iterator > energyInterpretations() const
BXVector< HGCalMulticluster > HGCalMulticlusterBxCollection
math::PtEtaPhiMLorentzVector iPolarP4(const HGCalMulticluster::EnergyInterpretation eInt) const
math::XYZTLorentzVector LorentzVector
boost::transform_iterator< KeyGetter< Iter >, Iter > key_iterator(Iter itr) const
double interpretationFraction(const HGCalMulticluster::EnergyInterpretation eInt) const
void saveEnergyInterpretation(const HGCalMulticluster::EnergyInterpretation eInt, double energy)
size_type interpretations_size() const
double phi() const final
momentum azimuthal angle
double energy() const final
energy
double eta() const final
momentum pseudorapidity