File indexing completed on 2024-04-06 12:04:32
0001 #ifndef DataFormats_L1Trigger_HGCalMulticluster_h
0002 #define DataFormats_L1Trigger_HGCalMulticluster_h
0003
0004 #include "DataFormats/Common/interface/Ptr.h"
0005 #include "DataFormats/L1Trigger/interface/BXVector.h"
0006 #include "DataFormats/L1THGCal/interface/HGCalClusterT.h"
0007 #include "DataFormats/L1THGCal/interface/HGCalCluster.h"
0008 #include <boost/iterator/transform_iterator.hpp>
0009 #include <functional>
0010
0011 namespace l1t {
0012
0013 class HGCalMulticluster : public HGCalClusterT<l1t::HGCalCluster> {
0014 public:
0015 HGCalMulticluster() : hOverEValid_(false) {}
0016 HGCalMulticluster(const LorentzVector p4, int pt = 0, int eta = 0, int phi = 0);
0017
0018 HGCalMulticluster(const edm::Ptr<l1t::HGCalCluster>& tc, float fraction = 1);
0019
0020 ~HGCalMulticluster() override;
0021
0022 float hOverE() const {
0023
0024
0025
0026
0027
0028 return hOverEValid_ ? hOverE_ : l1t::HGCalClusterT<l1t::HGCalCluster>::hOverE();
0029 }
0030
0031 void saveHOverE() {
0032 hOverE_ = l1t::HGCalClusterT<l1t::HGCalCluster>::hOverE();
0033 hOverEValid_ = true;
0034 }
0035
0036 enum EnergyInterpretation { EM = 0 };
0037
0038 void saveEnergyInterpretation(const HGCalMulticluster::EnergyInterpretation eInt, double energy);
0039
0040 double iEnergy(const HGCalMulticluster::EnergyInterpretation eInt) const {
0041 return energy() * interpretationFraction(eInt);
0042 }
0043
0044 double iPt(const HGCalMulticluster::EnergyInterpretation eInt) const { return pt() * interpretationFraction(eInt); }
0045
0046 math::XYZTLorentzVector iP4(const HGCalMulticluster::EnergyInterpretation eInt) const {
0047 return p4() * interpretationFraction(eInt);
0048 }
0049
0050 math::PtEtaPhiMLorentzVector iPolarP4(const HGCalMulticluster::EnergyInterpretation eInt) const {
0051 return math::PtEtaPhiMLorentzVector(pt() * interpretationFraction(eInt), eta(), phi(), 0.);
0052 }
0053
0054 private:
0055 template <typename Iter>
0056 struct KeyGetter {
0057 const typename Iter::value_type::first_type& operator()(const typename Iter::value_type& p) const {
0058 return p.first;
0059 }
0060 };
0061
0062 template <typename Iter>
0063 boost::transform_iterator<KeyGetter<Iter>, Iter> key_iterator(Iter itr) const {
0064 return boost::make_transform_iterator<KeyGetter<Iter>, Iter>(itr, KeyGetter<Iter>());
0065 }
0066
0067 public:
0068 typedef boost::transform_iterator<KeyGetter<std::map<EnergyInterpretation, double>::const_iterator>,
0069 std::map<EnergyInterpretation, double>::const_iterator>
0070 EnergyInterpretation_const_iterator;
0071
0072 std::pair<EnergyInterpretation_const_iterator, EnergyInterpretation_const_iterator> energyInterpretations() const {
0073 return std::make_pair(key_iterator(energyInterpretationFractions_.cbegin()),
0074 key_iterator(energyInterpretationFractions_.cend()));
0075 }
0076
0077 EnergyInterpretation_const_iterator interpretations_begin() const {
0078 return key_iterator(energyInterpretationFractions_.cbegin());
0079 }
0080
0081 EnergyInterpretation_const_iterator interpretations_end() const {
0082 return key_iterator(energyInterpretationFractions_.cend());
0083 }
0084
0085 size_type interpretations_size() const { return energyInterpretationFractions_.size(); }
0086
0087 private:
0088 double interpretationFraction(const HGCalMulticluster::EnergyInterpretation eInt) const;
0089
0090 float hOverE_;
0091 bool hOverEValid_;
0092 std::map<EnergyInterpretation, double> energyInterpretationFractions_;
0093 };
0094
0095 typedef BXVector<HGCalMulticluster> HGCalMulticlusterBxCollection;
0096
0097 }
0098
0099 #endif