CMS 3D CMS Logo

HGCalTriggerClusterInterpretationEM.cc
Go to the documentation of this file.
1 
6 
8 public:
11  void initialize(const edm::ParameterSet& conf) final;
12  void setGeometry(const HGCalTriggerGeometryBase* const) 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 
47 }
48 
50  for (unsigned int idx = 0; idx != multiclusters.size(); idx++) {
51  l1t::HGCalMulticluster& cluster3d = multiclusters[idx];
52 
53  const GlobalPoint& cluster3d_position = cluster3d.centreProj();
54  double energy = 0.;
55 
56  for (const auto& cluster2d : cluster3d.constituents()) {
57  const unsigned layer = triggerTools_.triggerLayer(cluster2d.first);
58  if (layer <= layer_containment_corrs_.size() - 1) {
59  double dr = (cluster3d_position - cluster2d.second->centreProj()).mag();
60  if (dr <= dr_bylayer_.at(layer)) {
61  energy += layer_containment_corrs_.at(layer) * cluster2d.second->energy();
62  }
63  }
64  }
66  cluster3d.saveEnergyInterpretation(l1t::HGCalMulticluster::EnergyInterpretation::EM, max(energy, 0.));
67  }
68 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void initialize(const edm::ParameterSet &conf) final
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:49
#define DEFINE_HGC_TPG_CLUSTER_INTERPRETER(type, name)
void setGeometry(const HGCalTriggerGeometryBase *const)
unsigned size(int bx) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setGeometry(const HGCalTriggerGeometryBase *const) final
unsigned triggerLayer(const unsigned id) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void saveEnergyInterpretation(const HGCalMulticluster::EnergyInterpretation eInt, double energy)
const GlobalPoint & centreProj() const
void interpret(l1t::HGCalMulticlusterBxCollection &multiclusters) const final
double eta() const final
momentum pseudorapidity