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 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  }
64  cluster3d.saveEnergyInterpretation(l1t::HGCalMulticluster::EnergyInterpretation::EM, max(energy, 0.));
65  }
66 }
HGCalTriggerTools.h
HGCalTriggerClusterInterpretationEM::interpret
void interpret(l1t::HGCalMulticlusterBxCollection &multiclusters) const final
Definition: HGCalTriggerClusterInterpretationEM.cc:47
HGCalTriggerTools::eventSetup
void eventSetup(const edm::EventSetup &)
Definition: HGCalTriggerTools.cc:35
DEFINE_HGC_TPG_CLUSTER_INTERPRETER
#define DEFINE_HGC_TPG_CLUSTER_INTERPRETER(type, name)
Definition: HGCalTriggerClusterInterpreterBase.h:19
HGCalTriggerClusterInterpretationEM::eventSetup
void eventSetup(const edm::EventSetup &es) final
Definition: HGCalTriggerClusterInterpretationEM.cc:45
HGCalTriggerClusterInterpretationEM::initialize
void initialize(const edm::ParameterSet &conf) final
Definition: HGCalTriggerClusterInterpretationEM.cc:27
HGCalTriggerClusterInterpretationEM::HGCalTriggerClusterInterpretationEM
HGCalTriggerClusterInterpretationEM()
Definition: HGCalTriggerClusterInterpretationEM.cc:25
HGCalTriggerClusterInterpreterBase.h
HGCalTriggerClusterInterpretationEM::~HGCalTriggerClusterInterpretationEM
~HGCalTriggerClusterInterpretationEM() override
Definition: HGCalTriggerClusterInterpretationEM.cc:10
HGCalTriggerClusterInterpretationEM::scale_corrections_coeff_
std::vector< double > scale_corrections_coeff_
Definition: HGCalTriggerClusterInterpretationEM.cc:17
HGCalTriggerClusterInterpretationEM::dr_bylayer_
std::vector< double > dr_bylayer_
Definition: HGCalTriggerClusterInterpretationEM.cc:18
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
BXVector
Definition: BXVector.h:15
l1t::HGCalMulticluster
Definition: HGCalMulticluster.h:13
HGCalMulticluster.h
l1t::HGCalMulticluster::saveEnergyInterpretation
void saveEnergyInterpretation(const HGCalMulticluster::EnergyInterpretation eInt, double energy)
Definition: HGCalMulticluster.cc:13
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
l1t::HGCalClusterT::constituents
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:50
Point3DBase< float, GlobalTag >
edm::ParameterSet
Definition: ParameterSet.h:47
reco::LeafCandidate::eta
double eta() const final
momentum pseudorapidity
Definition: LeafCandidate.h:152
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
HGCalTriggerClusterInterpretationEM::triggerTools_
HGCalTriggerTools triggerTools_
Definition: HGCalTriggerClusterInterpretationEM.cc:20
edm::EventSetup
Definition: EventSetup.h:57
HGCalTriggerTools::triggerLayer
unsigned triggerLayer(const unsigned id) const
Definition: HGCalTriggerTools.h:84
HGCalTriggerClusterInterpreterBase
Definition: HGCalTriggerClusterInterpreterBase.h:7
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
HGCalTriggerClusterInterpretationEM
Definition: HGCalTriggerClusterInterpretationEM.cc:7
Frameworkfwd.h
flavorHistoryFilter_cfi.dr
dr
Definition: flavorHistoryFilter_cfi.py:37
Exception
Definition: hltDiff.cc:246
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HGCalTriggerClusterInterpretationEM::layer_containment_corrs_
std::vector< double > layer_containment_corrs_
Definition: HGCalTriggerClusterInterpretationEM.cc:16
HGCalTriggerTools
Definition: HGCalTriggerTools.h:32
BXVector::size
unsigned size(int bx) const
l1t::HGCalClusterT::centreProj
const GlobalPoint & centreProj() const
Definition: HGCalClusterT.h:102
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22