CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HGCalTriggerClusterInterpretationEM.cc
Go to the documentation of this file.
1 
6 
8 public:
11  void initialize(const edm::ParameterSet& conf) final;
12  void eventSetup(const edm::EventSetup& es) final;
13  void interpret(l1t::HGCalMulticlusterBxCollection& multiclusters) const final;
14 
15 private:
16  std::vector<double> layer_containment_corrs_;
17  std::vector<double> scale_corrections_coeff_;
18  std::vector<double> dr_bylayer_;
19 
21 };
22 
23 DEFINE_HGC_TPG_CLUSTER_INTERPRETER(HGCalTriggerClusterInterpretationEM, "HGCalTriggerClusterInterpretationEM");
24 
26 
28  layer_containment_corrs_ = conf.getParameter<std::vector<double>>("layer_containment_corrs");
29  scale_corrections_coeff_ = conf.getParameter<std::vector<double>>("scale_correction_coeff");
30  dr_bylayer_ = conf.getParameter<std::vector<double>>("dr_bylayer");
31 
32  const unsigned corrections_size = 2;
33  if (scale_corrections_coeff_.size() != corrections_size) {
34  throw cms::Exception("HGCTriggerParameterError")
35  << "HGCalTriggerClusterInterpretationEM::scale_correction_coeff parameter has size: "
36  << scale_corrections_coeff_.size() << " while expected is " << corrections_size;
37  }
38  if (layer_containment_corrs_.size() != dr_bylayer_.size()) {
39  throw cms::Exception("HGCTriggerParameterError")
40  << "HGCalTriggerClusterInterpretationEM::layer_containment_corrs and "
41  "HGCalTriggerClusterInterpretationEM::dr_bylayer have different size!";
42  }
43 }
44 
46 
48  for (unsigned int idx = 0; idx != multiclusters.size(); idx++) {
49  l1t::HGCalMulticluster& cluster3d = multiclusters[idx];
50 
51  const GlobalPoint& cluster3d_position = cluster3d.centreProj();
52  double energy = 0.;
53 
54  for (const auto& cluster2d : cluster3d.constituents()) {
55  const unsigned layer = triggerTools_.triggerLayer(cluster2d.first);
56  if (layer <= layer_containment_corrs_.size() - 1) {
57  double dr = (cluster3d_position - cluster2d.second->centreProj()).mag();
58  if (dr <= dr_bylayer_.at(layer)) {
59  energy += layer_containment_corrs_.at(layer) * cluster2d.second->energy();
60  }
61  }
62  }
63  energy += scale_corrections_coeff_.at(1) * std::abs(cluster3d.eta()) + scale_corrections_coeff_.at(0);
64  cluster3d.saveEnergyInterpretation(l1t::HGCalMulticluster::EnergyInterpretation::EM, max(energy, 0.));
65  }
66 }
T getParameter(std::string const &) const
void interpret(l1t::HGCalMulticlusterBxCollection &multiclusters) const final
unsigned size(int bx) const
double eta() const final
momentum pseudorapidity
void eventSetup(const edm::EventSetup &)
void initialize(const edm::ParameterSet &conf) final
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const GlobalPoint & centreProj() const
#define DEFINE_HGC_TPG_CLUSTER_INTERPRETER(type, name)
unsigned triggerLayer(const unsigned id) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void saveEnergyInterpretation(const HGCalMulticluster::EnergyInterpretation eInt, double energy)
void eventSetup(const edm::EventSetup &es) final
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50