Project
Loading...
Searching...
No Matches
Clusterizer.h
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#ifndef ALICEO2_EMCAL_CLUSTERIZER_H
15#define ALICEO2_EMCAL_CLUSTERIZER_H
16
17#include <array>
18#include <gsl/span>
19#include "Rtypes.h"
23#include "EMCALBase/Geometry.h"
24
25namespace o2
26{
27
28namespace emcal
29{
30
31// Define numbers rows/columns for topological representation of cells
32constexpr unsigned int NROWS = (24 + 1) * (6 + 4); // 10x supermodule rows (6 for EMCAL, 4 for DCAL). +1 accounts for topological gap between two supermodules
33constexpr unsigned int NCOLS = 48 * 2 + 1; // 2x supermodule columns + 1 empty space in between for DCAL (not used for EMCAL)
34
36
44
45template <class InputType>
47{
50 struct cellWithE {
51
53 cellWithE() : energy(0.), row(0), column(0) {}
54
59 cellWithE(float e, int r, int c) : energy(e), row(r), column(c) {}
60
67 bool operator<(cellWithE const& rhs) const
68 {
69 return energy < rhs.energy;
70 }
71
72 float energy;
73 int row;
74 int column;
75 };
76
79 struct InputwithIndex {
80 const InputType* mInput;
81 ClusterIndex mIndex;
82 };
83
84 public:
93 Clusterizer(double timeCut, double timeMin, double timeMax, double gradientCut, bool doEnergyGradientCut, double thresholdSeedE, double thresholdCellE);
94
97
99 ~Clusterizer() = default;
100
102 void clear()
103 {
104 mFoundClusters.clear();
105 mInputIndices.clear();
106 }
107
116 void initialize(double timeCut, double timeMin, double timeMax, double gradientCut, bool doEnergyGradientCut, double thresholdSeedE, double thresholdCellE);
117
123 void findClusters(const gsl::span<InputType const>& inputArray);
124
127 const std::vector<Cluster>* getFoundClusters() const { return &mFoundClusters; }
128
131 const std::vector<ClusterIndex>* getFoundClustersInputIndices() const { return &mInputIndices; }
132
135 void setGeometry(Geometry* geometry) { mEMCALGeometry = geometry; }
136
139 Geometry* getGeometry() { return mEMCALGeometry; }
140
141 private:
146 void getClusterFromNeighbours(std::vector<InputwithIndex>& clusterInputs, int row, int column);
147
152 void getTopologicalRowColumn(const InputType& input, int& row, int& column);
153
154 Geometry* mEMCALGeometry = nullptr;
155 std::array<cellWithE, NROWS * NCOLS> mSeedList;
156 std::array<std::array<InputwithIndex, NCOLS>, NROWS> mInputMap;
157 std::array<std::array<bool, NCOLS>, NROWS> mCellMask;
158
159 std::vector<Cluster> mFoundClusters;
160 std::vector<ClusterIndex> mInputIndices;
161
162 double mTimeCut;
163 double mTimeMin;
164 double mTimeMax;
165 double mGradientCut;
166 bool mDoEnergyGradientCut;
167 double mThresholdSeedEnergy;
168 double mThresholdCellEnergy;
169 ClassDefNV(Clusterizer, 1);
170};
171
174
175} // namespace emcal
176} // namespace o2
177#endif /* ALICEO2_EMCAL_CLUSTERIZER_H */
uint32_t c
Definition RawData.h:2
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.
void clear()
Clear internal buffers of found clusters and cell indices.
Clusterizer()
Default constructor.
const std::vector< Cluster > * getFoundClusters() const
Get list of found clusters.
Geometry * getGeometry()
Get pointer to geometry.
~Clusterizer()=default
Destructor.
const std::vector< ClusterIndex > * getFoundClustersInputIndices() const
Get list of found cell indices.
void setGeometry(Geometry *geometry)
Set EMCAL geometry.
EMCAL geometry definition.
Definition Geometry.h:40
GLboolean r
Definition glcorearb.h:1233
constexpr unsigned int NCOLS
Definition Clusterizer.h:33
constexpr unsigned int NROWS
Definition Clusterizer.h:32
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::vector< int > row