Project
Loading...
Searching...
No Matches
GPUTPCClusterStatistics.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
16#include "GPULogging.h"
17#include "GPUO2DataTypes.h"
18#include <algorithm>
19#include <cstring>
20#include <map>
21#include <queue>
22
23using namespace o2::gpu;
24
25// Small helper to compute Huffman probabilities
26namespace o2::gpu
27{
28namespace // anonymous
29{
30typedef std::vector<bool> HuffCode;
31typedef std::map<uint32_t, HuffCode> HuffCodeMap;
32
33class INode
34{
35 public:
36 const double f;
37
38 virtual ~INode() = default;
39
40 protected:
41 INode(double v) : f(v) {}
42};
43
44class InternalNode : public INode
45{
46 public:
47 INode* const left;
48 INode* const right;
49
50 InternalNode(INode* c0, INode* c1) : INode(c0->f + c1->f), left(c0), right(c1) {}
51 ~InternalNode() override
52 {
53 delete left;
54 delete right;
55 }
56};
57
58class LeafNode : public INode
59{
60 public:
61 const uint32_t c;
62
63 LeafNode(double v, uint32_t w) : INode(v), c(w) {}
64};
65
66struct NodeCmp {
67 bool operator()(const INode* lhs, const INode* rhs) const { return lhs->f > rhs->f; }
68};
69
70INode* BuildTree(const double* frequencies, uint32_t UniqueSymbols)
71{
72 std::priority_queue<INode*, std::vector<INode*>, NodeCmp> trees;
73
74 for (uint32_t i = 0; i < UniqueSymbols; i++) {
75 if (frequencies[i] != 0) {
76 trees.push(new LeafNode(frequencies[i], i));
77 }
78 }
79 while (trees.size() > 1) {
80 INode* childR = trees.top();
81 trees.pop();
82
83 INode* childL = trees.top();
84 trees.pop();
85
86 INode* parent = new InternalNode(childR, childL);
87 trees.push(parent);
88 }
89 return trees.top();
90}
91
92void GenerateCodes(const INode* node, const HuffCode& prefix, HuffCodeMap& outCodes)
93{
94 if (const LeafNode* lf = dynamic_cast<const LeafNode*>(node)) {
95 outCodes[lf->c] = prefix;
96 } else if (const InternalNode* in = dynamic_cast<const InternalNode*>(node)) {
97 HuffCode leftPrefix = prefix;
98 leftPrefix.push_back(false);
99 GenerateCodes(in->left, leftPrefix, outCodes);
100
101 HuffCode rightPrefix = prefix;
102 rightPrefix.push_back(true);
103 GenerateCodes(in->right, rightPrefix, outCodes);
104 }
105}
106} // anonymous namespace
107} // namespace o2::gpu
108
110{
111 uint32_t decodingErrors = 0;
112 o2::tpc::ClusterNativeAccess clustersNativeDecoded;
113 std::vector<o2::tpc::ClusterNative> clusterBuffer;
114 GPUInfo("Compression statistics, decoding: %d attached (%d tracks), %d unattached", clustersCompressed->nAttachedClusters, clustersCompressed->nTracks, clustersCompressed->nUnattachedClusters);
115 auto allocator = [&clusterBuffer](size_t size) {clusterBuffer.resize(size); return clusterBuffer.data(); };
116 mDecoder.decompress(clustersCompressed, clustersNativeDecoded, allocator, param, true);
117 std::vector<o2::tpc::ClusterNative> tmpClusters;
118 if (param.rec.tpc.rejectionStrategy == GPUSettings::RejectionNone) { // verification does not make sense if we reject clusters during compression
119 for (uint32_t i = 0; i < NSECTORS; i++) {
120 for (uint32_t j = 0; j < GPUTPCGeometry::NROWS; j++) {
121 if (clustersNative->nClusters[i][j] != clustersNativeDecoded.nClusters[i][j]) {
122 GPUError("Number of clusters mismatch sector %u row %u: expected %d v.s. decoded %d", i, j, clustersNative->nClusters[i][j], clustersNativeDecoded.nClusters[i][j]);
123 decodingErrors++;
124 continue;
125 }
126 tmpClusters.resize(clustersNative->nClusters[i][j]);
127 for (uint32_t k = 0; k < clustersNative->nClusters[i][j]; k++) {
128 tmpClusters[k] = clustersNative->clusters[i][j][k];
129 if (param.rec.tpc.compressionTypeMask & GPUSettings::CompressionTruncate) {
130 GPUTPCCompression::truncateSignificantBitsChargeMax(tmpClusters[k].qMax, param);
131 GPUTPCCompression::truncateSignificantBitsWidth(tmpClusters[k].sigmaPadPacked, param);
132 if (!tmpClusters[k].isSaturated()) [[likely]] {
133 GPUTPCCompression::truncateSignificantBitsCharge(tmpClusters[k].qTot, param);
134 GPUTPCCompression::truncateSignificantBitsWidth(tmpClusters[k].sigmaTimePacked, param);
135 }
136 }
137 }
138 std::sort(tmpClusters.begin(), tmpClusters.end());
139 for (uint32_t k = 0; k < clustersNative->nClusters[i][j]; k++) {
140 const o2::tpc::ClusterNative& c1 = tmpClusters[k];
141 const o2::tpc::ClusterNative& c2 = clustersNativeDecoded.clusters[i][j][k];
142 if (c1.timeFlagsPacked != c2.timeFlagsPacked || c1.padPacked != c2.padPacked || c1.sigmaTimePacked != c2.sigmaTimePacked || c1.sigmaPadPacked != c2.sigmaPadPacked || c1.qMax != c2.qMax || c1.qTot != c2.qTot) {
143 if (decodingErrors++ < 100) {
144 GPUWarning("Cluster mismatch: sector %2u row %3u hit %5u: %6d %3d %4d %3d %3d %4d %4d", i, j, k, (int32_t)c1.getTimePacked(), (int32_t)c1.getFlags(), (int32_t)c1.padPacked, (int32_t)c1.sigmaTimePacked, (int32_t)c1.sigmaPadPacked, (int32_t)c1.qMax, (int32_t)c1.qTot);
145 GPUWarning("%45s %6d %3d %4d %3d %3d %4d %4d", "", (int32_t)c2.getTimePacked(), (int32_t)c2.getFlags(), (int32_t)c2.padPacked, (int32_t)c2.sigmaTimePacked, (int32_t)c2.sigmaPadPacked, (int32_t)c2.qMax, (int32_t)c2.qTot);
146 }
147 }
148 }
149 }
150 }
151 if (decodingErrors) {
152 mDecodingError = true;
153 GPUWarning("Errors during cluster decoding %u\n", decodingErrors);
154 } else {
155 GPUInfo("Cluster decoding verification: PASSED");
156 }
157 }
158
159 FillStatistic(mPqTotA, clustersCompressed->qTotA, clustersCompressed->nAttachedClusters);
160 FillStatistic(mPqMaxA, clustersCompressed->qMaxA, clustersCompressed->nAttachedClusters);
161 FillStatistic(mPflagsA, clustersCompressed->flagsA, clustersCompressed->nAttachedClusters);
162 FillStatistic(mProwDiffA, clustersCompressed->rowDiffA, clustersCompressed->nAttachedClustersReduced);
163 FillStatistic(mPsectorLegDiffA, clustersCompressed->sliceLegDiffA, clustersCompressed->nAttachedClustersReduced);
164 FillStatistic(mPpadResA, clustersCompressed->padResA, clustersCompressed->nAttachedClustersReduced);
165 FillStatistic(mPtimeResA, clustersCompressed->timeResA, clustersCompressed->nAttachedClustersReduced);
166 FillStatistic(mPsigmaPadA, clustersCompressed->sigmaPadA, clustersCompressed->nAttachedClusters);
167 FillStatistic(mPsigmaTimeA, clustersCompressed->sigmaTimeA, clustersCompressed->nAttachedClusters);
168 FillStatistic(mPqPtA, clustersCompressed->qPtA, clustersCompressed->nTracks);
169 FillStatistic(mProwA, clustersCompressed->rowA, clustersCompressed->nTracks);
170 FillStatistic(mPsectorA, clustersCompressed->sliceA, clustersCompressed->nTracks);
171 FillStatistic(mPtimeA, clustersCompressed->timeA, clustersCompressed->nTracks);
172 FillStatistic(mPpadA, clustersCompressed->padA, clustersCompressed->nTracks);
173 FillStatistic(mPqTotU, clustersCompressed->qTotU, clustersCompressed->nUnattachedClusters);
174 FillStatistic(mPqMaxU, clustersCompressed->qMaxU, clustersCompressed->nUnattachedClusters);
175 FillStatistic(mPflagsU, clustersCompressed->flagsU, clustersCompressed->nUnattachedClusters);
176 FillStatistic(mPpadDiffU, clustersCompressed->padDiffU, clustersCompressed->nUnattachedClusters);
177 FillStatistic(mPtimeDiffU, clustersCompressed->timeDiffU, clustersCompressed->nUnattachedClusters);
178 FillStatistic(mPsigmaPadU, clustersCompressed->sigmaPadU, clustersCompressed->nUnattachedClusters);
179 FillStatistic(mPsigmaTimeU, clustersCompressed->sigmaTimeU, clustersCompressed->nUnattachedClusters);
180 FillStatistic<uint16_t, 1>(mPnTrackClusters, clustersCompressed->nTrackClusters, clustersCompressed->nTracks);
181 FillStatistic<uint32_t, 1>(mPnSectorRowClusters, clustersCompressed->nSliceRowClusters, clustersCompressed->nSliceRows);
182 FillStatisticCombined(mPsigmaA, clustersCompressed->sigmaPadA, clustersCompressed->sigmaTimeA, clustersCompressed->nAttachedClusters, P_MAX_SIGMA);
183 FillStatisticCombined(mPsigmaU, clustersCompressed->sigmaPadU, clustersCompressed->sigmaTimeU, clustersCompressed->nUnattachedClusters, P_MAX_SIGMA);
184 FillStatisticCombined(mPQA, clustersCompressed->qMaxA, clustersCompressed->qTotA, clustersCompressed->nAttachedClusters, P_MAX_QMAX);
185 FillStatisticCombined(mPQU, clustersCompressed->qMaxU, clustersCompressed->qTotU, clustersCompressed->nUnattachedClusters, P_MAX_QMAX);
186 FillStatisticCombined(mProwSectorA, clustersCompressed->rowDiffA, clustersCompressed->sliceLegDiffA, clustersCompressed->nAttachedClustersReduced, GPUTPCGeometry::NROWS);
187 mNTotalClusters += clustersCompressed->nAttachedClusters + clustersCompressed->nUnattachedClusters;
188}
189
191{
192 if (mDecodingError) {
193 GPUError("-----------------------------------------\nERROR - INCORRECT CLUSTER DECODING!\n-----------------------------------------");
194 }
195 if (mNTotalClusters == 0) {
196 return;
197 }
198
199 GPUInfo("\nRunning cluster compression entropy statistics");
200 double eQ = Analyze(mPqTotA, "qTot Attached", false);
201 eQ += Analyze(mPqMaxA, "qMax Attached", false);
202 Analyze(mPflagsA, "flags Attached");
203 double eRowSector = Analyze(mProwDiffA, "rowDiff Attached", false);
204 eRowSector += Analyze(mPsectorLegDiffA, "sectorDiff Attached", false);
205 Analyze(mPpadResA, "padRes Attached");
206 Analyze(mPtimeResA, "timeRes Attached");
207 double eSigma = Analyze(mPsigmaPadA, "sigmaPad Attached", false);
208 eSigma += Analyze(mPsigmaTimeA, "sigmaTime Attached", false);
209 Analyze(mPqPtA, "qPt Attached");
210 Analyze(mProwA, "row Attached");
211 Analyze(mPsectorA, "sector Attached");
212 Analyze(mPtimeA, "time Attached");
213 Analyze(mPpadA, "pad Attached");
214 eQ += Analyze(mPqTotU, "qTot Unattached", false);
215 eQ += Analyze(mPqMaxU, "qMax Unattached", false);
216 Analyze(mPflagsU, "flags Unattached");
217 Analyze(mPpadDiffU, "padDiff Unattached");
218 Analyze(mPtimeDiffU, "timeDiff Unattached");
219 eSigma += Analyze(mPsigmaPadU, "sigmaPad Unattached", false);
220 eSigma += Analyze(mPsigmaTimeU, "sigmaTime Unattached", false);
221 Analyze(mPnTrackClusters, "nClusters in Track");
222 Analyze(mPnSectorRowClusters, "nClusters in Row");
223 double eSigmaCombined = Analyze(mPsigmaA, "combined sigma Attached");
224 eSigmaCombined += Analyze(mPsigmaU, "combined sigma Unattached");
225 double eQCombined = Analyze(mPQA, "combined Q Attached");
226 eQCombined += Analyze(mPQU, "combined Q Unattached");
227 double eRowSectorCombined = Analyze(mProwSectorA, "combined row/sector Attached");
228
229 GPUInfo("Combined Row/Sector: %6.4f --> %6.4f (%6.4f%%)", eRowSector, eRowSectorCombined, eRowSector > 1e-1 ? (100. * (eRowSector - eRowSectorCombined) / eRowSector) : 0.f);
230 GPUInfo("Combined Sigma: %6.4f --> %6.4f (%6.4f%%)", eSigma, eSigmaCombined, eSigma > 1e-3 ? (100. * (eSigma - eSigmaCombined) / eSigma) : 0.f);
231 GPUInfo("Combined Q: %6.4f --> %6.4f (%6.4f%%)", eQ, eQCombined, eQ > 1e-3 ? (100. * (eQ - eQCombined) / eQ) : 0.f);
232
233 printf("\nCombined Entropy: %7.4f (Size %'13.0f, %'zu clusters)\nCombined Huffman: %7.4f (Size %'13.0f, %f%%)\n\n", mEntropy / mNTotalClusters, mEntropy / 8., mNTotalClusters, mHuffman / mNTotalClusters, mHuffman / 8., 100. * (mHuffman - mEntropy) / mHuffman);
234}
235
236float GPUTPCClusterStatistics::Analyze(std::vector<int32_t>& p, const char* name, bool count)
237{
238 double entropy = 0.;
239 double huffmanSize = 0;
240
241 std::vector<double> prob(p.size());
242 double log2 = log(2.);
243 size_t total = 0;
244 for (uint32_t i = 0; i < p.size(); i++) {
245 total += p[i];
246 }
247 if (total) {
248 for (uint32_t i = 0; i < prob.size(); i++) {
249 if (p[i]) {
250 prob[i] = (double)p[i] / total;
251 double I = -log(prob[i]) / log2;
252 double H = I * prob[i];
253
254 entropy += H;
255 }
256 }
257
258 INode* root = BuildTree(prob.data(), prob.size());
259
260 HuffCodeMap codes;
261 GenerateCodes(root, HuffCode(), codes);
262 delete root;
263
264 for (HuffCodeMap::const_iterator it = codes.begin(); it != codes.end(); it++) {
265 huffmanSize += it->second.size() * prob[it->first];
266 }
267
268 if (count) {
269 mEntropy += entropy * total;
270 mHuffman += huffmanSize * total;
271 }
272 }
273 GPUInfo("Size: %30s: Entropy %7.4f Huffman %7.4f (Count) %9ld", name, entropy, huffmanSize, (int64_t)total);
274 return entropy;
275}
276
277template <class T, int32_t I>
278void GPUTPCClusterStatistics::FillStatistic(std::vector<int32_t>& p, const T* ptr, size_t n)
279{
280 for (size_t i = 0; i < n; i++) {
281 uint32_t val = ptr[i];
282 if (val >= p.size()) {
283 if (I) {
284 p.resize(val + 1);
285 } else {
286 GPUError("Invalid Value: %d >= %d", val, (int32_t)p.size());
287 continue;
288 }
289 }
290 p[val]++;
291 }
292}
293
294template <class T, class S, int32_t I>
295void GPUTPCClusterStatistics::FillStatisticCombined(std::vector<int32_t>& p, const T* ptr1, const S* ptr2, size_t n, int32_t max1)
296{
297 for (size_t i = 0; i < n; i++) {
298 uint32_t val = ptr1[i] + ptr2[i] * max1;
299 if (val >= p.size()) {
300 if (I) {
301 p.resize(val + 1);
302 } else {
303 GPUError("Invalid Value: %d >= %d", val, (int32_t)p.size());
304 continue;
305 }
306 }
307 p[val]++;
308 }
309}
std::unique_ptr< expressions::Node > node
int32_t i
uint32_t j
Definition RawData.h:0
TBranch * ptr
float Analyze(std::vector< int32_t > &p, const char *name, bool count=true)
void FillStatisticCombined(std::vector< int32_t > &p, const T *ptr1, const S *ptr2, size_t n, int32_t max1)
void FillStatistic(std::vector< int32_t > &p, const T *ptr, size_t n)
void RunStatistics(const o2::tpc::ClusterNativeAccess *clustersNative, const o2::tpc::CompressedClusters *clustersCompressed, const GPUParam &param)
static constexpr uint32_t NROWS
static int32_t decompress(const o2::tpc::CompressedClustersFlat *clustersCompressed, o2::tpc::ClusterNativeAccess &clustersNative, std::function< o2::tpc::ClusterNative *(size_t)> allocator, const GPUParam &param, bool deterministicRec)
GLdouble n
Definition glcorearb.h:1982
GLint GLsizei count
Definition glcorearb.h:399
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLuint const GLchar * name
Definition glcorearb.h:781
GLdouble GLdouble right
Definition glcorearb.h:4077
GLdouble f
Definition glcorearb.h:310
GLuint GLfloat * val
Definition glcorearb.h:1582
GLenum GLfloat param
Definition glcorearb.h:271
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
D const SVectorGPU< T, D > & rhs
Definition SMatrixGPU.h:193
unsigned int nClusters[constants::MAXSECTOR][constants::MAXGLOBALPADROW]
const ClusterNative * clusters[constants::MAXSECTOR][constants::MAXGLOBALPADROW]