CMS 3D CMS Logo

HGCalLayerClustersFromSoAProducer.cc
Go to the documentation of this file.
2 
4 
12 
19 
21 
22 #define DEBUG_CLUSTERS_ALPAKA 0
23 
24 #if DEBUG_CLUSTERS_ALPAKA
26 #endif
27 
29 public:
31  : getTokenSoAClusters_(consumes(config.getParameter<edm::InputTag>("src"))),
32  getTokenSoACells_(consumes(config.getParameter<edm::InputTag>("hgcalRecHitsSoA"))),
33  getTokenSoARecHitsExtra_(consumes(config.getParameter<edm::InputTag>("hgcalRecHitsLayerClustersSoA"))),
34  detector_(config.getParameter<std::string>("detector")),
35  hitsTime_(config.getParameter<unsigned int>("nHitsTime")),
36  timeClname_(config.getParameter<std::string>("timeClname")) {
37 #if DEBUG_CLUSTERS_ALPAKA
38  moduleLabel_ = config.getParameter<std::string>("@module_label");
39 #endif
40  if (detector_ == "HFNose") {
42  } else if (detector_ == "EE") {
44  } else { //for FH or BH
46  }
47 
48  produces<std::vector<float>>("InitialLayerClustersMask");
49  produces<std::vector<reco::BasicCluster>>();
50  produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
51  }
52 
53  ~HGCalLayerClustersFromSoAProducer() override = default;
54 
55  void produce(edm::Event& iEvent, edm::EventSetup const& iSetup) override {
56  auto const& deviceData = iEvent.get(getTokenSoAClusters_);
57 
58  auto const& deviceSoARecHitsExtra = iEvent.get(getTokenSoARecHitsExtra_);
59  auto const soaRecHitsExtra_v = deviceSoARecHitsExtra.view();
60 
61  auto const& deviceSoACells = iEvent.get(getTokenSoACells_);
62  auto const soaCells_v = deviceSoACells.view();
63 
64  auto const deviceView = deviceData.view();
65 
66  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
67  clusters->reserve(deviceData->metadata().size());
68 
69  // Create a vector of <clusters> locations, where each location holds a
70  // vector of <nCells> floats. These vectors are used to compute the time for
71  // each cluster.
72  std::vector<std::vector<float>> times(deviceData->metadata().size());
73  std::vector<std::vector<float>> timeErrors(deviceData->metadata().size());
74 
75  for (int i = 0; i < deviceData->metadata().size(); ++i) {
76  std::vector<std::pair<DetId, float>> thisCluster;
77  thisCluster.reserve(deviceView.cells(i));
78  clusters->emplace_back(deviceView.energy(i),
79  math::XYZPoint(deviceView.x(i), deviceView.y(i), deviceView.z(i)),
81  std::move(thisCluster),
82  algoId_);
83  clusters->back().setSeed(deviceView.seed(i));
84  times[i].reserve(deviceView.cells(i));
85  timeErrors[i].reserve(deviceView.cells(i));
86  }
87 
88  // Populate hits and fractions required to compute the cluster's time.
89  // This procedure is complex and involves two SoAs: the original RecHits
90  // SoA and the clustering algorithm's output SoA. Both SoAs have the same
91  // cardinality, and crucially, the output SoA includes the cluster index.
92  for (int32_t i = 0; i < soaRecHitsExtra_v.metadata().size(); ++i) {
93  if (soaRecHitsExtra_v[i].clusterIndex() == -1) {
94  continue;
95  }
96  assert(soaRecHitsExtra_v[i].clusterIndex() < (int)clusters->size());
97  (*clusters)[soaRecHitsExtra_v[i].clusterIndex()].addHitAndFraction(soaCells_v[i].detid(), 1.f);
98  if (soaCells_v[i].timeError() < 0.f) {
99  continue;
100  }
101  times[soaRecHitsExtra_v[i].clusterIndex()].push_back(soaCells_v[i].time());
102  timeErrors[soaRecHitsExtra_v[i].clusterIndex()].push_back(
103  1.f / (soaCells_v[i].timeError() * soaCells_v[i].timeError()));
104  }
105 
106  // Finally, compute and assign the time to each cluster.
107  std::vector<std::pair<float, float>> cluster_times;
108  cluster_times.reserve(clusters->size());
110  for (unsigned i = 0; i < clusters->size(); ++i) {
111  if (detector_ != "BH") {
112  cluster_times.push_back(timeEstimator.fixSizeHighestDensity(times[i], timeErrors[i], hitsTime_));
113  } else {
114  cluster_times.push_back(std::pair<float, float>(-99.f, -1.f));
115  }
116  }
117 
118 #if DEBUG_CLUSTERS_ALPAKA
119  auto runNumber = iEvent.eventAuxiliary().run();
120  auto lumiNumber = iEvent.eventAuxiliary().luminosityBlock();
121  auto evtNumber = iEvent.eventAuxiliary().id().event();
122 
123  hgcalUtils::DumpCellsSoA dumperCellsSoA;
124  dumperCellsSoA.dumpInfos(deviceSoACells, moduleLabel_, runNumber, lumiNumber, evtNumber);
125 
127  dumper.dumpInfos(*clusters, moduleLabel_, runNumber, lumiNumber, evtNumber, true);
128 
129  hgcalUtils::DumpClustersSoA dumperSoA;
130  dumperSoA.dumpInfos(deviceSoARecHitsExtra, moduleLabel_, runNumber, lumiNumber, evtNumber);
131 #endif
132 
133  auto clusterHandle = iEvent.put(std::move(clusters));
134 
135  auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
137  filler.insert(clusterHandle, cluster_times.begin(), cluster_times.end());
138  filler.fill();
139  iEvent.put(std::move(timeCl), timeClname_);
140 
141  // The layerClusterMask for the HGCAL detector is created at a later
142  // stage, when the layer clusters from the different components of HGCAL
143  // are merged together into a unique collection. For the case of HFNose,
144  // since there is no further merging step needed, we create the
145  // layerClustersMask directly here.
146  if (detector_ == "HFNose") {
147  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
148  layerClustersMask->resize(clusterHandle->size(), 1.0);
149  iEvent.put(std::move(layerClustersMask), "InitialLayerClustersMask");
150  }
151  }
152 
155  desc.add<edm::InputTag>("src", edm::InputTag("hltHgcalSoALayerClustersProducer"));
156  desc.add<edm::InputTag>("hgcalRecHitsLayerClustersSoA", edm::InputTag("hltHgcalSoARecHitsLayerClustersProducer"));
157  desc.add<edm::InputTag>("hgcalRecHitsSoA", edm::InputTag("hltHgcalSoARecHitsProducer"));
158  desc.add<unsigned int>("nHitsTime", 3);
159  desc.add<std::string>("timeClname", "timeLayerCluster");
161  "detector", "EE", true, edm::Comment("the HGCAL component used to create clusters.")),
162  edm::allowedValues<std::string>("EE", "FH"));
163  descriptions.addWithDefaultLabel(desc);
164  }
165 
166 private:
171  unsigned int hitsTime_;
174 #if DEBUG_CLUSTERS_ALPAKA
175  std::string moduleLabel_;
176 #endif
177 };
edm::EDGetTokenT< HGCalSoARecHitsExtraHostCollection > const getTokenSoARecHitsExtra_
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
HGCalLayerClustersFromSoAProducer(edm::ParameterSet const &config)
void produce(edm::Event &iEvent, edm::EventSetup const &iSetup) override
void dumpInfos(const T &clusters, const std::string &moduleLabel, edm::RunNumber_t run, edm::LuminosityBlockNumber_t lumi, edm::EventNumber_t event, bool dumpCellsDetId=false) const
Definition: config.py:1
edm::EDGetTokenT< HGCalSoAClustersHostCollection > const getTokenSoAClusters_
assert(be >=bs)
std::pair< float, float > fixSizeHighestDensity(std::vector< float > &time, std::vector< float > weight=std::vector< float >(), unsigned int minNhits=3, float deltaT=0.210, float timeWidthBy=0.5)
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< HGCalSoARecHitsHostCollection > const getTokenSoACells_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~HGCalLayerClustersFromSoAProducer() override=default
void dumpInfos(const T &cells, const std::string &moduleLabel, edm::RunNumber_t run, edm::LuminosityBlockNumber_t lumi, edm::EventNumber_t event) const
void dumpInfos(const T &clustersSoA, const std::string &moduleLabel, edm::RunNumber_t run, edm::LuminosityBlockNumber_t lumi, edm::EventNumber_t event) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
HLT enums.
def move(src, dest)
Definition: eostools.py:511