12 #include "oneapi/tbb/task_arena.h" 13 #include "oneapi/tbb.h" 19 template <
typename T,
typename STRATEGY>
22 numberOfClustersPerLayer_.clear();
23 cells_.resize(2 * (maxlayer_ + 1));
24 numberOfClustersPerLayer_.resize(2 * (maxlayer_ + 1), 0);
27 template <
typename T,
typename STRATEGY>
36 for (
unsigned int i = 0;
i <
hits.size(); ++
i) {
39 unsigned int layerOnSide = (rhtools_.getLayerWithOffset(detid) - 1);
43 float sigmaNoise = 1.f;
45 int thickness_index = rhtools_.getSiThickIndex(detid);
46 if (thickness_index == -1)
47 thickness_index = maxNumberOfThickIndices_;
49 double storedThreshold = thresholds_[layerOnSide][thickness_index];
51 storedThreshold = thresholds_[layerOnSide][thickness_index + deltasi_index_regemfac_];
53 sigmaNoise = v_sigmaNoise_[layerOnSide][thickness_index];
55 if (hgrh.
energy() < storedThreshold)
59 if (!dependSensor_ && hgrh.
energy() < ecut_)
62 int offset = ((rhtools_.zside(detid) + 1) >> 1) * maxlayer_;
63 int layer = layerOnSide +
offset;
65 cells_[layer].detid.emplace_back(detid);
66 if constexpr (std::is_same_v<STRATEGY, HGCalScintillatorStrategy>) {
67 cells_[layer].dim1.emplace_back(
position.eta());
68 cells_[layer].dim2.emplace_back(
position.phi());
71 cells_[layer].dim1.emplace_back(
position.x());
72 cells_[layer].dim2.emplace_back(
position.y());
74 cells_[layer].weight.emplace_back(hgrh.
energy());
75 cells_[layer].sigmaNoise.emplace_back(sigmaNoise);
79 template <
typename T,
typename STRATEGY>
81 auto cellsSize = cells_[
l].detid.size();
82 cells_[
l].rho.resize(cellsSize, 0.
f);
83 cells_[
l].delta.resize(cellsSize, 9999999);
84 cells_[
l].nearestHigher.resize(cellsSize, -1);
85 cells_[
l].clusterIndex.resize(cellsSize, -1);
86 cells_[
l].followers.resize(cellsSize);
87 cells_[
l].isSeed.resize(cellsSize,
false);
94 template <
typename T,
typename STRATEGY>
97 tbb::this_task_arena::isolate([&] {
98 tbb::parallel_for(
size_t(0),
size_t(2 * maxlayer_ + 2), [&](
size_t i) {
99 prepareDataStructures(
i);
102 lt.fill(cells_[
i].dim1, cells_[
i].dim2);
105 if constexpr (std::is_same_v<STRATEGY, HGCalSiliconStrategy>) {
108 if (
i % maxlayer_ < lastLayerEE_)
109 delta_c = vecDeltas_[0];
110 else if (
i % maxlayer_ < (firstLayerBH_ - 1))
111 delta_c = vecDeltas_[1];
113 delta_c = vecDeltas_[2];
119 LogDebug(
"HGCalCLUEAlgo") <<
"maxlayer: " << maxlayer_ <<
" lastLayerEE: " << lastLayerEE_
120 <<
" firstLayerBH: " << firstLayerBH_ <<
"\n";
122 calculateLocalDensity(lt,
i,
delta);
123 calculateDistanceToHigher(lt,
i,
delta);
124 numberOfClustersPerLayer_[
i] = findAndAssignClusters(
i,
delta);
129 template <
typename T,
typename STRATEGY>
131 std::vector<int>
offsets(numberOfClustersPerLayer_.size(), 0);
133 int maxClustersOnLayer = numberOfClustersPerLayer_[0];
135 for (
unsigned layerId = 1; layerId <
offsets.size(); ++layerId) {
136 offsets[layerId] =
offsets[layerId - 1] + numberOfClustersPerLayer_[layerId - 1];
138 maxClustersOnLayer =
std::max(maxClustersOnLayer, numberOfClustersPerLayer_[layerId]);
141 auto totalNumberOfClusters =
offsets.back() + numberOfClustersPerLayer_.back();
142 clusters_v_.resize(totalNumberOfClusters);
143 std::vector<std::vector<int>> cellsIdInCluster;
144 cellsIdInCluster.reserve(maxClustersOnLayer);
146 for (
unsigned int layerId = 0; layerId < 2 * maxlayer_ + 2; ++layerId) {
147 cellsIdInCluster.resize(numberOfClustersPerLayer_[layerId]);
148 auto& cellsOnLayer = cells_[layerId];
149 unsigned int numberOfCells = cellsOnLayer.detid.size();
150 auto firstClusterIdx =
offsets[layerId];
152 for (
unsigned int i = 0;
i < numberOfCells; ++
i) {
153 auto clusterIndex = cellsOnLayer.clusterIndex[
i];
154 if (clusterIndex != -1)
155 cellsIdInCluster[clusterIndex].push_back(
i);
158 std::vector<std::pair<DetId, float>> thisCluster;
160 for (
auto&
cl : cellsIdInCluster) {
165 for (
auto cellIdx :
cl) {
166 energy += cellsOnLayer.weight[cellIdx];
167 thisCluster.emplace_back(cellsOnLayer.detid[cellIdx], 1.f);
168 if (cellsOnLayer.isSeed[cellIdx]) {
169 seedDetId = cellsOnLayer.detid[cellIdx];
172 auto globalClusterIndex = cellsOnLayer.clusterIndex[
cl[0]] + firstClusterIdx;
174 clusters_v_[globalClusterIndex] =
176 clusters_v_[globalClusterIndex].setSeed(seedDetId);
180 cellsIdInCluster.clear();
184 template <
typename T,
typename STRATEGY>
186 const unsigned int layerId,
189 auto& cellsOnLayer = cells_[layerId];
190 unsigned int numberOfCells = cellsOnLayer.detid.size();
191 for (
unsigned int i = 0;
i < numberOfCells;
i++) {
192 std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[
i] -
delta,
193 cellsOnLayer.dim1[
i] +
delta,
194 cellsOnLayer.dim2[
i] -
delta,
195 cellsOnLayer.dim2[
i] +
delta);
197 for (
int xBin = search_box[0];
xBin < search_box[1] + 1; ++
xBin) {
198 for (
int yBin = search_box[2];
yBin < search_box[3] + 1; ++
yBin) {
199 int binId = lt.getGlobalBinByBin(
xBin,
yBin);
200 size_t binSize = lt[binId].size();
202 for (
unsigned int j = 0;
j < binSize;
j++) {
203 unsigned int otherId = lt[binId][
j];
205 cellsOnLayer.rho[
i] += (
i == otherId ? 1.f : 0.5f) * cellsOnLayer.weight[otherId];
210 LogDebug(
"HGCalCLUEAlgo") <<
"Debugging calculateLocalDensity: \n" 211 <<
" cell: " <<
i <<
" eta: " << cellsOnLayer.dim1[
i] <<
" phi: " << cellsOnLayer.dim2[
i]
212 <<
" energy: " << cellsOnLayer.weight[
i] <<
" density: " << cellsOnLayer.rho[
i] <<
"\n";
215 template <
typename T,
typename STRATEGY>
217 const unsigned int layerId,
220 auto& cellsOnLayer = cells_[layerId];
221 unsigned int numberOfCells = cellsOnLayer.detid.size();
222 for (
unsigned int i = 0;
i < numberOfCells;
i++) {
223 std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[
i] -
delta,
224 cellsOnLayer.dim1[
i] +
delta,
225 cellsOnLayer.dim2[
i] -
delta,
226 cellsOnLayer.dim2[
i] +
delta);
227 cellsOnLayer.rho[
i] += cellsOnLayer.weight[
i];
228 float northeast(0), northwest(0), southeast(0), southwest(0),
all(0);
232 int binId = lt.getGlobalBinByBin(
etaBin, phi);
233 size_t binSize = lt[binId].size();
234 for (
unsigned int j = 0;
j < binSize;
j++) {
235 unsigned int otherId = lt[binId][
j];
241 int dIPhi = otherIPhi - iPhi;
242 dIPhi +=
abs(dIPhi) < 2 ? 0
243 : dIPhi < 0 ? scintMaxIphi_
245 int dIEta = otherIEta -
iEta;
246 LogDebug(
"HGCalCLUEAlgo") <<
" Debugging calculateLocalDensity for Scintillator: \n" 247 <<
" cell: " << otherId <<
" energy: " << cellsOnLayer.weight[otherId]
248 <<
" otherIPhi: " << otherIPhi <<
" iPhi: " << iPhi <<
" otherIEta: " << otherIEta
249 <<
" iEta: " <<
iEta <<
"\n";
252 auto neighborCellContribution = 0.5f * cellsOnLayer.weight[otherId];
253 all += neighborCellContribution;
254 if (dIPhi >= 0 && dIEta >= 0)
255 northeast += neighborCellContribution;
256 if (dIPhi <= 0 && dIEta >= 0)
257 southeast += neighborCellContribution;
258 if (dIPhi >= 0 && dIEta <= 0)
259 northwest += neighborCellContribution;
260 if (dIPhi <= 0 && dIEta <= 0)
261 southwest += neighborCellContribution;
263 LogDebug(
"HGCalCLUEAlgo") <<
" Debugging calculateLocalDensity for Scintillator: \n" 264 <<
" northeast: " << northeast <<
" southeast: " << southeast
265 <<
" northwest: " << northwest <<
" southwest: " << southwest <<
"\n";
270 float neighborsval = (
std::max(northeast, northwest) >
std::max(southeast, southwest))
274 cellsOnLayer.rho[
i] += neighborsval;
276 cellsOnLayer.rho[
i] +=
all;
277 LogDebug(
"HGCalCLUEAlgo") <<
"Debugging calculateLocalDensity: \n" 278 <<
" cell: " <<
i <<
" eta: " << cellsOnLayer.dim1[
i] <<
" phi: " << cellsOnLayer.dim2[
i]
279 <<
" energy: " << cellsOnLayer.weight[
i] <<
" density: " << cellsOnLayer.rho[
i] <<
"\n";
282 template <
typename T,
typename STRATEGY>
284 if constexpr (std::is_same_v<STRATEGY, HGCalSiliconStrategy>) {
291 template <
typename T,
typename STRATEGY>
293 auto& cellsOnLayer = cells_[layerId];
294 unsigned int numberOfCells = cellsOnLayer.detid.size();
296 for (
unsigned int i = 0;
i < numberOfCells;
i++) {
300 int i_nearestHigher = -1;
302 std::array<int, 4> search_box = lt.searchBox(cellsOnLayer.dim1[
i] -
range,
303 cellsOnLayer.dim1[
i] +
range,
304 cellsOnLayer.dim2[
i] -
range,
305 cellsOnLayer.dim2[
i] +
range);
307 for (
int dim1Bin = search_box[0]; dim1Bin < search_box[1] + 1; ++dim1Bin) {
308 for (
int dim2Bin = search_box[2]; dim2Bin < search_box[3] + 1; ++dim2Bin) {
310 size_t binId = lt.getGlobalBinByBin(dim1Bin, dim2Bin);
311 if constexpr (std::is_same_v<STRATEGY, HGCalScintillatorStrategy>)
312 binId = lt.getGlobalBinByBin(dim1Bin, (dim2Bin % T::type::nRows));
314 size_t binSize = lt[binId].size();
317 for (
unsigned int j = 0;
j < binSize;
j++) {
318 unsigned int otherId = lt[binId][
j];
319 float dist =
distance(lt,
i, otherId, layerId);
321 (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[
i]) ||
322 (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[
i] && cellsOnLayer.detid[otherId] > cellsOnLayer.detid[
i]);
323 if (foundHigher && dist <= i_delta) {
327 i_nearestHigher = otherId;
332 bool foundNearestHigherInSearchBox = (i_delta !=
maxDelta);
333 if (foundNearestHigherInSearchBox) {
334 cellsOnLayer.delta[
i] = i_delta;
335 cellsOnLayer.nearestHigher[
i] = i_nearestHigher;
340 cellsOnLayer.nearestHigher[
i] = -1;
343 LogDebug(
"HGCalCLUEAlgo") <<
"Debugging calculateDistanceToHigher: \n" 344 <<
" cell: " <<
i <<
" eta: " << cellsOnLayer.dim1[
i] <<
" phi: " << cellsOnLayer.dim2[
i]
345 <<
" energy: " << cellsOnLayer.weight[
i] <<
" density: " << cellsOnLayer.rho[
i]
346 <<
" nearest higher: " << cellsOnLayer.nearestHigher[
i]
347 <<
" distance: " << cellsOnLayer.delta[
i] <<
"\n";
351 template <
typename T,
typename STRATEGY>
357 unsigned int nClustersOnLayer = 0;
358 auto& cellsOnLayer = cells_[layerId];
359 unsigned int numberOfCells = cellsOnLayer.detid.size();
360 std::vector<int> localStack;
362 for (
unsigned int i = 0;
i < numberOfCells;
i++) {
363 float rho_c = kappa_ * cellsOnLayer.sigmaNoise[
i];
365 cellsOnLayer.clusterIndex[
i] = -1;
366 bool isSeed = (cellsOnLayer.delta[
i] >
delta) && (cellsOnLayer.rho[
i] >= rho_c);
367 bool isOutlier = (cellsOnLayer.delta[
i] > outlierDeltaFactor_ *
delta) && (cellsOnLayer.rho[
i] < rho_c);
369 cellsOnLayer.clusterIndex[
i] = nClustersOnLayer;
370 cellsOnLayer.isSeed[
i] =
true;
372 localStack.push_back(
i);
374 }
else if (!isOutlier) {
375 cellsOnLayer.followers[cellsOnLayer.nearestHigher[
i]].push_back(
i);
380 while (!localStack.empty()) {
381 int endStack = localStack.back();
382 auto& thisSeed = cellsOnLayer.followers[endStack];
383 localStack.pop_back();
386 for (
int j : thisSeed) {
388 cellsOnLayer.clusterIndex[
j] = cellsOnLayer.clusterIndex[endStack];
390 localStack.push_back(
j);
393 return nClustersOnLayer;
396 template <
typename T,
typename STRATEGY>
411 std::vector<double>
dummy;
413 dummy.resize(maxNumberOfThickIndices_ + !isNose_, 0);
414 thresholds_.resize(maxlayer_,
dummy);
415 v_sigmaNoise_.resize(maxlayer_,
dummy);
417 for (
unsigned ilayer = 1; ilayer <= maxlayer_; ++ilayer) {
418 for (
unsigned ithick = 0; ithick < maxNumberOfThickIndices_; ++ithick) {
419 float sigmaNoise = 0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
420 (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
421 thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
422 v_sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
423 LogDebug(
"HGCalCLUEAlgo") <<
"ilayer: " << ilayer <<
" nonAgedNoises: " << nonAgedNoises_[ithick]
424 <<
" fcPerEle: " << fcPerEle_ <<
" fcPerMip: " << fcPerMip_[ithick]
425 <<
" noiseMip: " << fcPerEle_ * nonAgedNoises_[ithick] / fcPerMip_[ithick]
426 <<
" sigmaNoise: " << sigmaNoise <<
"\n";
430 float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer] / sciThicknessCorrection_;
431 thresholds_[ilayer - 1][maxNumberOfThickIndices_] = ecut_ * scintillators_sigmaNoise;
432 v_sigmaNoise_[ilayer - 1][maxNumberOfThickIndices_] = scintillators_sigmaNoise;
433 LogDebug(
"HGCalCLUEAlgo") <<
"ilayer: " << ilayer <<
" noiseMip: " << noiseMip_
434 <<
" scintillators_sigmaNoise: " << scintillators_sigmaNoise <<
"\n";
constexpr const DetId & detid() const
void calculateDistanceToHigher(const TILE <, const unsigned int layerId, float delta)
constexpr Detector det() const
get the detector field from this detid
void populate(const HGCRecHitCollection &hits) override
constexpr float energy() const
int iphi() const
get the phi index
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Abs< T >::type abs(const T &t)
void makeClusters() override
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum) ...
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
void calculateLocalDensity(const TILE <, const unsigned int layerId, float delta)
XYZPointD XYZPoint
point in space with cartesian internal representation
void prepareDataStructures(const unsigned int layerId)
static int position[264][3]
int findAndAssignClusters(const unsigned int layerId, float delta)
std::vector< reco::BasicCluster > getClusters(bool) override