Project
Loading...
Searching...
No Matches
KrBoxClusterFinder.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
15
74
75#ifndef ALICEO2_TPC_KrBoxClusterFinder_H_
76#define ALICEO2_TPC_KrBoxClusterFinder_H_
77
81
82#include "TPCBase/Mapper.h"
83
84#include "TPCBase/CalDet.h"
85
86#include <tuple>
87#include <vector>
88#include <array>
89#include <deque>
90#include <gsl/span>
91
92namespace o2
93{
94namespace tpc
95{
96
104
106{
107 private:
110 static constexpr size_t MaxPads = 138;
111 static constexpr size_t MaxRows = 152;
112 size_t mMaxTimes = 114048;
113
114 public:
118 explicit KrBoxClusterFinder() = default;
119
120 explicit KrBoxClusterFinder(int sector) : mSector(sector){};
121
124 void loadGainMapFromFile(const std::string_view calDetFileName, const std::string_view gainMapName = "GainMap");
125
127 std::vector<std::tuple<int, int, int>> findLocalMaxima(bool directFilling = true, const int timeOffset = 0);
128
131 KrCluster buildCluster(int clusterCenterPad, int clusterCenterRow, int clusterCenterTime, bool directFilling = false, const int timeOffset = 0);
132
134 void resetClusters() { mClusters.clear(); }
135
136 std::vector<KrCluster>& getClusters() { return mClusters; }
137 const std::vector<KrCluster>& getClusters() const { return mClusters; }
138
140 void setSector(int sector) { mSector = sector; }
141
143 int getSector() const { return mSector; }
144
146 void init();
147
149 void setMinNumberOfNeighbours(int minNumberOfNeighbours) { mMinNumberOfNeighbours = minNumberOfNeighbours; }
150
152 void setMinQTreshold(int minQThreshold) { mQThresholdMax = minQThreshold; }
153
155 void setMaxTimes(int maxTimes) { mMaxTimes = maxTimes; }
156
158 void setMaxClusterSize(int maxClusterSizeRowIROC, int maxClusterSizeRowOROC1, int maxClusterSizeRowOROC2, int maxClusterSizeRowOROC3,
159 int maxClusterSizePadIROC, int maxClusterSizePadOROC1, int maxClusterSizePadOROC2, int maxClusterSizePadOROC3,
160 int maxClusterSizeTime)
161 {
162 mMaxClusterSizeRowIROC = maxClusterSizeRowIROC;
163 mMaxClusterSizeRowOROC1 = maxClusterSizeRowOROC1;
164 mMaxClusterSizeRowOROC2 = maxClusterSizeRowOROC2;
165 mMaxClusterSizeRowOROC3 = maxClusterSizeRowOROC3;
166
167 mMaxClusterSizePadIROC = maxClusterSizePadIROC;
168 mMaxClusterSizePadOROC1 = maxClusterSizePadOROC1;
169 mMaxClusterSizePadOROC2 = maxClusterSizePadOROC2;
170 mMaxClusterSizePadOROC3 = maxClusterSizePadOROC3;
171
172 mMaxClusterSizeTime = maxClusterSizeTime;
173 }
174
175 void loopOverSector(const gsl::span<const Digit> eventSector, const int sector);
176
177 void loopOverSector(const std::vector<Digit>& eventSector, const int sector)
178 {
179 loopOverSector(gsl::span(eventSector.data(), eventSector.size()), sector);
180 }
181
182 private:
183 // These variables can be varied
184 // They were choses such that the box in each readout chamber is approx. the same size
185 // NOTE: They will be overwritten by the values in KrBoxClusterFinderParam in case the init() function is called
186 int mMaxClusterSizeTime = 3;
187 int mMaxClusterSizeRow;
188 int mMaxClusterSizePad;
189
190 int mMaxClusterSizeRowIROC = 3;
191 int mMaxClusterSizeRowOROC1 = 2;
192 int mMaxClusterSizeRowOROC2 = 2;
193 int mMaxClusterSizeRowOROC3 = 1;
194
195 int mMaxClusterSizePadIROC = 5;
196 int mMaxClusterSizePadOROC1 = 3;
197 int mMaxClusterSizePadOROC2 = 3;
198 int mMaxClusterSizePadOROC3 = 3;
199
200 float mQThresholdMax = 30.0;
201 float mQThreshold = 1.0;
202 int mMinNumberOfNeighbours = 2;
203
204 float mCutMinSigmaTime{0};
205 float mCutMaxSigmaTime{1000};
206 float mCutMinSigmaPad{0};
207 float mCutMaxSigmaPad{1000};
208 float mCutMinSigmaRow{0};
209 float mCutMaxSigmaRow{1000};
210 float mCutMaxQtot{1e10};
211 float mCutQtot0{1e10};
212 float mCutQtotSizeSlope{0};
213 unsigned char mCutMaxSize{255};
214 bool mApplyCuts{false};
215
216 int mSector = -1;
217 std::unique_ptr<CalPad> mGainMap;
218
220 static constexpr size_t MaxRowsIROC = 63;
221 static constexpr size_t MaxRowsOROC1 = 34;
222 static constexpr size_t MaxRowsOROC2 = 30;
223 static constexpr size_t MaxRowsOROC3 = 25;
224
227 std::vector<o2::tpc::KrCluster> mClusters;
228
230 const Mapper& mMapperInstance = o2::tpc::Mapper::instance();
231
232 KrCluster mTempCluster;
233
234 using TimeSliceSector = std::array<std::array<float, MaxPads>, MaxRows>;
235
237 // TimeSliceSector mTempTimeSlice;
238
240 size_t mFirstDigit = 0;
241
260 std::deque<TimeSliceSector> mSetOfTimeSlices{};
261
263 struct ThresholdInfo {
264 bool digitAboveThreshold{};
265 std::array<bool, MaxRows> rowAboveThreshold{};
266 };
267
268 std::deque<ThresholdInfo> mThresholdInfo{};
269
270 void createInitialMap(const gsl::span<const Digit> eventSector);
271 void popFirstTimeSliceFromMap()
272 {
273 mSetOfTimeSlices.pop_front();
274 mThresholdInfo.pop_front();
275 }
276 void fillADCValueInLastSlice(int cru, int rowInSector, int padInRow, float adcValue);
277 void addTimeSlice(const gsl::span<const Digit> eventSector, const int timeSlice);
278
280 void setMaxClusterSize(int row);
281
283 void updateTempCluster(float tempCharge, int tempPad, int tempRow, int tempTime);
285 void updateTempClusterFinal(const int timeOffset = 0);
286
288 int signnum(int val) { return (0 < val) - (val < 0); }
289
291 bool acceptCluster(const KrCluster& cl);
292
293 ClassDefNV(KrBoxClusterFinder, 0);
294};
295
296} // namespace tpc
297} // namespace o2
298
299#endif
Definition of the TPC Digit.
Parameter class for the Kr box cluster finder.
Struct for Krypton and X-ray clusters.
void loopOverSector(const std::vector< Digit > &eventSector, const int sector)
void resetClusters()
reset cluster vector
void loadGainMapFromFile(const std::string_view calDetFileName, const std::string_view gainMapName="GainMap")
void init()
initialize the parameters from KrBoxClusterFinderParam
void setMaxTimes(int maxTimes)
Set Function for minimal charge required for maxCharge of a cluster.
void setMinNumberOfNeighbours(int minNumberOfNeighbours)
Set Function for minimum number of direct neighbours required.
void setSector(int sector)
set sector of this instance
std::vector< KrCluster > & getClusters()
void loopOverSector(const gsl::span< const Digit > eventSector, const int sector)
std::vector< std::tuple< int, int, int > > findLocalMaxima(bool directFilling=true, const int timeOffset=0)
After the map is created, we look for local maxima with this function:
void setMinQTreshold(int minQThreshold)
Set Function for minimal charge required for maxCharge of a cluster.
int getSector() const
ger sector of this instance
KrCluster buildCluster(int clusterCenterPad, int clusterCenterRow, int clusterCenterTime, bool directFilling=false, const int timeOffset=0)
void setMaxClusterSize(int maxClusterSizeRowIROC, int maxClusterSizeRowOROC1, int maxClusterSizeRowOROC2, int maxClusterSizeRowOROC3, int maxClusterSizePadIROC, int maxClusterSizePadOROC1, int maxClusterSizePadOROC2, int maxClusterSizePadOROC3, int maxClusterSizeTime)
Set Function for maximal cluster sizes in different ROCs.
const std::vector< KrCluster > & getClusters() const
static Mapper & instance(const std::string mappingDir="")
Definition Mapper.h:44
GLuint GLfloat * val
Definition glcorearb.h:1582
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::vector< int > row