Project
Loading...
Searching...
No Matches
Clusterizer.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
14#include <cstring>
15#include <gsl/span>
16#include <fairlogger/Logger.h> // for LOG
18
19using namespace o2::emcal;
20
21//____________________________________________________________________________
22template <class InputType>
23Clusterizer<InputType>::Clusterizer(double timeCut, double timeMin, double timeMax, double gradientCut, bool doEnergyGradientCut, double thresholdSeedE, double thresholdCellE) : mSeedList(), mInputMap(), mCellMask(), mTimeCut(timeCut), mTimeMin(timeMin), mTimeMax(timeMax), mGradientCut(gradientCut), mDoEnergyGradientCut(doEnergyGradientCut), mThresholdSeedEnergy(thresholdSeedE), mThresholdCellEnergy(thresholdCellE)
24{
25}
26
27//____________________________________________________________________________
28template <class InputType>
29Clusterizer<InputType>::Clusterizer() : mSeedList(), mInputMap(), mCellMask(), mTimeCut(0), mTimeMin(0), mTimeMax(0), mGradientCut(0), mDoEnergyGradientCut(false), mThresholdSeedEnergy(0), mThresholdCellEnergy(0)
30{
31}
32
33//____________________________________________________________________________
34template <class InputType>
35void Clusterizer<InputType>::initialize(double timeCut, double timeMin, double timeMax, double gradientCut, bool doEnergyGradientCut, double thresholdSeedE, double thresholdCellE)
36{
37 mTimeCut = timeCut;
38 mTimeMin = timeMin;
39 mTimeMax = timeMax;
40 mGradientCut = gradientCut;
41 mDoEnergyGradientCut = doEnergyGradientCut;
42 mThresholdSeedEnergy = thresholdSeedE;
43 mThresholdCellEnergy = thresholdCellE;
44}
45
46//____________________________________________________________________________
47template <class InputType>
48void Clusterizer<InputType>::getClusterFromNeighbours(std::vector<InputwithIndex>& clusterInputs, int row, int column)
49{
50 // Recursion 0, add seed cell/digit to cluster
51 if (!clusterInputs.size()) {
52 clusterInputs.emplace_back(mInputMap[row][column]);
53 }
54
55 // Mark the current cell as clustered
56 mCellMask[row][column] = kTRUE;
57
58 // Now go recursively to the next 4 neighbours and add them to the cluster if they fulfill the conditions
59 constexpr int rowDiffs[4] = {-1, 0, 0, 1};
60 constexpr int colDiffs[4] = {0, -1, 1, 0};
61 for (int dir = 0; dir < 4; dir++) {
62 if ((row + rowDiffs[dir] < 0) || (row + rowDiffs[dir] >= NROWS)) {
63 continue;
64 }
65 if ((column + colDiffs[dir] < 0) || (column + colDiffs[dir] >= NCOLS)) {
66 continue;
67 }
68
69 if (mInputMap[row + rowDiffs[dir]][column + colDiffs[dir]].mInput) {
70 if (!mCellMask[row + rowDiffs[dir]][column + colDiffs[dir]]) {
71 if (mDoEnergyGradientCut && (mInputMap[row + rowDiffs[dir]][column + colDiffs[dir]].mInput->getEnergy() > mInputMap[row][column].mInput->getEnergy() + mGradientCut)) {
72 continue;
73 }
74 if (not(TMath::Abs(mInputMap[row + rowDiffs[dir]][column + colDiffs[dir]].mInput->getTimeStamp() - mInputMap[row][column].mInput->getTimeStamp()) > mTimeCut)) {
75 getClusterFromNeighbours(clusterInputs, row + rowDiffs[dir], column + colDiffs[dir]);
76 // Add the cell/digit to the current cluster -- if we end up here, the selected cluster fulfills the condition
77 clusterInputs.emplace_back(mInputMap[row + rowDiffs[dir]][column + colDiffs[dir]]);
78 }
79 }
80 }
81 }
82}
83
84//____________________________________________________________________________
85template <class InputType>
86void Clusterizer<InputType>::getTopologicalRowColumn(const InputType& input, int& row, int& column)
87{
88 // Get SM number and relative row/column for SM
89 auto cellIndex = mEMCALGeometry->GetCellIndex(input.getTower());
90 int nSupMod = std::get<0>(cellIndex);
91
92 auto phiEtaIndex = mEMCALGeometry->GetCellPhiEtaIndexInSModule(nSupMod, std::get<1>(cellIndex), std::get<2>(cellIndex), std::get<3>(cellIndex));
93 row = std::get<0>(phiEtaIndex);
94 column = std::get<1>(phiEtaIndex);
95
96 // Add shifts wrt. supermodule and type of calorimeter
97 // NOTE:
98 // * Rows (phi) are arranged that one space is left empty between supermodules in phi
99 // This is due to the physical gap that forbids clustering
100 // * For DCAL, there is an additional empty column between two supermodules in eta
101 // Again, this is to account for the gap in DCAL
102
103 row += nSupMod / 2 * (24 + 1);
104 // In DCAL, leave a gap between two SMs with same phi
105 if (!mEMCALGeometry->IsDCALSM(nSupMod)) { // EMCAL
106 column += nSupMod % 2 * 48;
107 } else {
108 column += nSupMod % 2 * (48 + 1);
109 }
110}
111
112//____________________________________________________________________________
113template <class InputType>
114void Clusterizer<InputType>::findClusters(const gsl::span<InputType const>& inputArray)
115{
116 clear();
117
118 // Algorithm
119 // - Fill cells/digits in 2D topological map
120 // - Fill struct arrays (energy,x,y) (to get mapping energy -> (x,y))
121 // - Create 2D bitmap (cell/digit is already clustered or not)
122 // - Sort struct arrays with descending energy
123 //
124 // - Loop over arrays:
125 // --> Check 2D bitmap (don't use cell/digit which are already clustered)
126 // --> Take valid cell/digit with highest energy as seed (they are already sorted)
127 // --> Recursive to neighboughs and create cluster
128 // --> Seed cell and all neighbours belonging to cluster will be put in 2D bitmap
129
130 // Reset cell/digit maps and cell masks
131 // Loop over one array dim, then reset each array
132 for (auto iArr = 0; iArr < NROWS; iArr++) {
133 mCellMask[iArr].fill(kFALSE);
134 mInputMap[iArr].fill({nullptr, -1});
135 }
136
137 // Calibrate cells/digits and fill the maps/arrays
138 int nCells = 0;
139 double ehs = 0.0;
140 //for (auto dig : inputArray) {
141 for (int iIndex = 0; iIndex < inputArray.size(); iIndex++) {
142
143 auto& dig = inputArray[iIndex];
144
145 Float_t inputEnergy = dig.getEnergy();
146 Float_t time = dig.getTimeStamp();
147
148 if (inputEnergy < mThresholdCellEnergy || time > mTimeMax || time < mTimeMin) {
149 continue;
150 }
151 if (!mEMCALGeometry->CheckAbsCellId(dig.getTower())) {
152 continue;
153 }
154 ehs += inputEnergy;
155
156 // Put cell/digit to 2D map
157 int row = 0, column = 0;
158 getTopologicalRowColumn(dig, row, column);
159 // not referencing dig here to get proper reference and not local copy
160 mInputMap[row][column].mInput = inputArray.data() + iIndex; //
161 mInputMap[row][column].mIndex = iIndex; // mInputMap saves the position of cells/digits in the input array
162 mSeedList[nCells].energy = inputEnergy;
163 mSeedList[nCells].row = row;
164 mSeedList[nCells].column = column;
165 nCells++;
166 }
167
168 // Sort struct arrays with ascending energy
169 std::sort(mSeedList.begin(), std::next(std::begin(mSeedList), nCells));
170
171 // Take next valid cell/digit in calorimeter as seed (in descending energy order)
172 for (int i = nCells - 1; i >= 0; i--) {
173 int row = mSeedList[i].row, column = mSeedList[i].column;
174 // Continue if the cell is already masked (i.e. was already clustered)
175 if (mCellMask[row][column]) {
176 continue;
177 }
178 // Continue if energy constraints are not fulfilled
179 if (mSeedList[i].energy <= mThresholdSeedEnergy) {
180 continue;
181 }
182
183 // Seed is found, form cluster recursively
184 std::vector<InputwithIndex> clusterInputs;
185 getClusterFromNeighbours(clusterInputs, row, column);
186
187 // Add cells/digits for current cluster to cell/digit index vector
188 int inputIndexStart = mInputIndices.size();
189 for (auto dig : clusterInputs) {
190 mInputIndices.emplace_back(dig.mIndex);
191 }
192 int inputIndexSize = mInputIndices.size() - inputIndexStart;
193
194 // Now form cluster object from cells/digits
195 mFoundClusters.emplace_back(mInputMap[row][column].mInput->getTimeStamp(), inputIndexStart, inputIndexSize); // Cluster object initialized w/ time of seed cell, start + size of associated cells
196 }
197 LOG(debug) << mFoundClusters.size() << "clusters found from " << nCells << " cells/digits (total=" << inputArray.size() << ")-> ehs " << ehs << " (minE " << mThresholdCellEnergy << ")";
198}
199
Definition of the EMCAL clusterizer.
int16_t time
Definition RawEventData.h:4
int32_t i
std::ostringstream debug
Meta class for recursive clusterizer.
Definition Clusterizer.h:47
void initialize(double timeCut, double timeMin, double timeMax, double gradientCut, bool doEnergyGradientCut, double thresholdSeedE, double thresholdCellE)
Initialize class member vars if not done in constructor.
void findClusters(const gsl::span< InputType const > &inputArray)
Find clusters based on a give input collection.
Clusterizer()
Default constructor.
constexpr unsigned int NCOLS
Definition Clusterizer.h:33
constexpr unsigned int NROWS
Definition Clusterizer.h:32
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
vec clear()
std::vector< int > row