CMS 3D CMS Logo

HGCalClusteringDummyImpl.cc
Go to the documentation of this file.
2 #include <unordered_map>
3 #include <unordered_set>
6 
7 // class constructor
9  : calibSF_(conf.getParameter<double>("calibSF_cluster")),
10  layerWeights_(conf.getParameter<std::vector<double>>("layerWeights")),
11  applyLayerWeights_(conf.getParameter<bool>("applyLayerCalibration")) {
12  edm::LogInfo("HGCalClusterParameters") << "C2d global calibration factor: " << calibSF_;
13 }
14 
15 // Create one cluster per TC for direct TC->3D clustering
18  std::vector<l1t::HGCalCluster> clustersTmp;
19  for (std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin();
20  tc != triggerCellsPtrs.end();
21  ++tc) {
22  clustersTmp.emplace_back(*tc);
23  }
24 
25  /* store clusters in the persistent collection */
26  clusters.resize(0, clustersTmp.size());
27  for (unsigned i(0); i < clustersTmp.size(); ++i) {
28  calibratePt(clustersTmp.at(i));
29  clusters.set(0, i, clustersTmp.at(i));
30  }
31 }
32 
34  double calibPt = 0.;
35 
36  if (applyLayerWeights_ && !triggerTools_.isNose(cluster.detId())) {
37  unsigned layerN = triggerTools_.layerWithOffset(cluster.detId());
38 
39  if (layerWeights_.at(layerN) == 0.) {
40  throw cms::Exception("BadConfiguration")
41  << "2D cluster energy forced to 0 by calibration coefficients.\n"
42  << "The configuration should be changed. "
43  << "Discarded layers should be defined in hgcalTriggerGeometryESProducer.TriggerGeometry.DisconnectedLayers "
44  "and not with calibration coefficients = 0\n";
45  }
46 
47  calibPt = layerWeights_.at(layerN) * cluster.mipPt();
48 
49  }
50 
51  else {
52  calibPt = cluster.pt() * calibSF_;
53  }
54 
55  cluster.setPt(calibPt);
56 }
double pt() const final
transverse momentum
uint32_t detId() const
Definition: HGCalClusterT.h:92
void setPt(double pt)
Definition: HGCalClusterT.h:94
void clusterizeDummy(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters)
HGCalClusteringDummyImpl(const edm::ParameterSet &conf)
std::vector< double > layerWeights_
void calibratePt(l1t::HGCalCluster &cluster)
bool isNose(const DetId &) const
unsigned layerWithOffset(const DetId &) const
Log< level::Info, false > LogInfo
double mipPt() const
Definition: HGCalClusterT.h:90