Project
Loading...
Searching...
No Matches
ClusterFactory.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
13#include <array>
14#include <gsl/span>
15#include "Rtypes.h"
23#include "EMCALBase/Geometry.h"
24#include "MathUtils/Cartesian.h"
25
27
28using namespace o2::emcal;
29
30template <class InputType>
31ClusterFactory<InputType>::ClusterFactory(gsl::span<const o2::emcal::Cluster> clustersContainer, gsl::span<const InputType> inputsContainer, gsl::span<const int> cellsIndices)
32{
33 setContainer(clustersContainer, inputsContainer, cellsIndices);
34}
35
36template <class InputType>
38{
39 mClustersContainer = gsl::span<const o2::emcal::Cluster>();
40 mInputsContainer = gsl::span<const InputType>();
41 mCellsIndices = gsl::span<int>();
42 mLookUpInit = false;
43 mCellLabelContainer = gsl::span<const o2::emcal::CellLabel>();
44}
45
48//____________________________________________________________________________
49template <class InputType>
51{
52 if (clusterIndex >= mClustersContainer.size()) {
53 throw ClusterRangeException(clusterIndex, mClustersContainer.size());
54 }
55 if (!mGeomPtr) {
57 }
58
59 o2::emcal::AnalysisCluster clusterAnalysis;
60 clusterAnalysis.setID(clusterIndex);
61
62 int firstCellIndex = mClustersContainer[clusterIndex].getCellIndexFirst();
63 int nCells = mClustersContainer[clusterIndex].getNCells();
64
65 gsl::span<const int> inputsIndices = gsl::span<const int>(&mCellsIndices[firstCellIndex], nCells);
66
67 // First calculate the index of input with maximum amplitude and get
68 // the supermodule number where it sits.
69
70 auto [inputIndMax, inputEnergyMax, cellAmp, shared] = getMaximalEnergyIndex(inputsIndices);
71
72 short towerId = mInputsContainer[inputIndMax].getTower();
73
74 float exoticTime = mInputsContainer[inputIndMax].getTimeStamp();
75
76 float fCross = 0.;
77
78 try {
79 clusterAnalysis.setIsExotic(isExoticCell(towerId, inputEnergyMax, exoticTime, fCross));
80 clusterAnalysis.setFCross(fCross);
81 } catch (UninitLookUpTableException& e) {
82 LOG(error) << e.what();
83 }
84
85 clusterAnalysis.setIndMaxInput(inputIndMax);
86
87 clusterAnalysis.setE(cellAmp);
88
89 mSuperModuleNumber = mGeomPtr->GetSuperModuleNumber(towerId);
90
91 clusterAnalysis.setNCells(inputsIndices.size());
92
93 std::vector<unsigned short> cellsIdices;
94
95 bool addClusterLabels = ((clusterLabel != nullptr) && (mCellLabelContainer.size() > 0));
96 for (auto cellIndex : inputsIndices) {
97 cellsIdices.push_back(cellIndex);
98 if (addClusterLabels) {
99 for (size_t iLabel = 0; iLabel < mCellLabelContainer[cellIndex].GetLabelSize(); iLabel++) {
100 if (mCellLabelContainer[cellIndex].GetAmplitudeFraction(iLabel) <= 0.f) {
101 continue; // skip 0 entries
102 }
103 clusterLabel->addValue(mCellLabelContainer[cellIndex].GetLabel(iLabel),
104 mCellLabelContainer[cellIndex].GetAmplitudeFraction(iLabel) * mInputsContainer[cellIndex].getEnergy());
105 }
106 }
107 }
108 if (addClusterLabels) {
109 clusterLabel->orderLabels();
110 clusterLabel->normalize(cellAmp);
111 }
112
113 clusterAnalysis.setCellsIndices(cellsIdices);
114
115 // evaluate global and local position
116 evalGlobalPosition(inputsIndices, clusterAnalysis);
117 evalLocalPosition(inputsIndices, clusterAnalysis);
118
119 // evaluate shower parameters
120 evalElipsAxis(inputsIndices, clusterAnalysis);
121 evalDispersion(inputsIndices, clusterAnalysis);
122
123 // evaluate number of local maxima
124 evalNExMax(inputsIndices, clusterAnalysis);
125
126 evalCoreEnergy(inputsIndices, clusterAnalysis);
127 evalTime(inputsIndices, clusterAnalysis);
128
129 // TODO to be added at a later stage
130 // evalPrimaries(inputsIndices, clusterAnalysis);
131 // evalParents(inputsIndices, clusterAnalysis);
132
133 // TODO to be added at a later stage
134 // Called last because it sets the global position of the cluster?
135 // Do not call it when recalculating clusters out of standard reconstruction
136 // if (!mJustCluster)
137 // evalLocal2TrackingCSTransform();
138
139 return clusterAnalysis;
140}
141
145//____________________________________________________________________________
146template <class InputType>
147void ClusterFactory<InputType>::evalDispersion(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
148{
149 double d = 0., wtot = 0.;
150 int nstat = 0;
151
152 // Calculates the dispersion in cell units
153 double etaMean = 0.0, phiMean = 0.0;
154
155 // Calculate mean values
156 for (auto iInput : inputsIndices) {
157
158 if (clusterAnalysis.E() > 0 && mInputsContainer[iInput].getEnergy() > 0) {
159 auto [nSupMod, nModule, nIphi, nIeta] = mGeomPtr->GetCellIndex(mInputsContainer[iInput].getTower());
160 auto [iphi, ieta] = mGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta);
161
162 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
163 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
164 if (mSharedCluster && nSupMod % 2) {
165 ieta += EMCAL_COLS;
166 }
167
168 double etai = (double)ieta;
169 double phii = (double)iphi;
170 double w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
171
172 if (w > 0.0) {
173 phiMean += phii * w;
174 etaMean += etai * w;
175 wtot += w;
176 }
177 }
178 }
179
180 if (wtot > 0) {
181 phiMean /= wtot;
182 etaMean /= wtot;
183 } else {
184 LOG(error) << Form("Wrong weight %f\n", wtot);
185 }
186
187 // Calculate dispersion
188 for (auto iInput : inputsIndices) {
189
190 if (clusterAnalysis.E() > 0 && mInputsContainer[iInput].getEnergy() > 0) {
191 auto [nSupMod, nModule, nIphi, nIeta] = mGeomPtr->GetCellIndex(mInputsContainer[iInput].getTower());
192 auto [iphi, ieta] = mGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta);
193
194 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
195 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
196 if (mSharedCluster && nSupMod % 2) {
197 ieta += EMCAL_COLS;
198 }
199
200 double etai = (double)ieta;
201 double phii = (double)iphi;
202 double w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
203
204 if (w > 0.0) {
205 nstat++;
206 d += w * ((etai - etaMean) * (etai - etaMean) + (phii - phiMean) * (phii - phiMean));
207 }
208 }
209 }
210
211 if (wtot > 0 && nstat > 1) {
212 d /= wtot;
213 } else {
214 d = 0.;
215 }
216
217 clusterAnalysis.setDispersion(TMath::Sqrt(d));
218}
219
222//____________________________________________________________________________
223template <class InputType>
224void ClusterFactory<InputType>::evalLocalPosition(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
225{
226
227 int nstat = 0;
228
229 double dist = tMaxInCm(double(clusterAnalysis.E()));
230
231 double clXYZ[3] = {0., 0., 0.}, clRmsXYZ[3] = {0., 0., 0.}, xyzi[3], wtot = 0., w = 0.;
232
233 for (auto iInput : inputsIndices) {
234
235 try {
236 mGeomPtr->RelPosCellInSModule(mInputsContainer[iInput].getTower(), dist).GetCoordinates(xyzi[0], xyzi[1], xyzi[2]);
237 } catch (InvalidCellIDException& e) {
238 LOG(error) << e.what();
239 continue;
240 }
241
242 // Temporal patch, due to mapping problem, need to swap "y" in one of the 2 SM, although no effect in position calculation. GCB 05/2010
243 if (mSharedCluster && mSuperModuleNumber != mGeomPtr->GetSuperModuleNumber(mInputsContainer[iInput].getTower())) {
244 xyzi[1] *= -1;
245 }
246
247 if (mLogWeight > 0.0) {
248 w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
249 } else {
250 w = mInputsContainer[iInput].getEnergy(); // just energy
251 }
252
253 if (w > 0.0) {
254 wtot += w;
255 nstat++;
256
257 for (int i = 0; i < 3; i++) {
258 clXYZ[i] += (w * xyzi[i]);
259 clRmsXYZ[i] += (w * xyzi[i] * xyzi[i]);
260 }
261 } // w > 0
262 } // dig loop
263
264 // cout << " wtot " << wtot << endl;
265
266 if (wtot > 0) {
267 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
268 for (int i = 0; i < 3; i++) {
269 clXYZ[i] /= wtot;
270
271 if (nstat > 1) {
272 clRmsXYZ[i] /= (wtot * wtot);
273 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i] * clXYZ[i];
274
275 if (clRmsXYZ[i] > 0.0) {
276 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
277 } else {
278 clRmsXYZ[i] = 0;
279 }
280 } else {
281 clRmsXYZ[i] = 0;
282 }
283 }
284 } else {
285 for (int i = 0; i < 3; i++) {
286 clXYZ[i] = clRmsXYZ[i] = -1.;
287 }
288 }
289
290 clusterAnalysis.setLocalPosition(math_utils::Point3D<float>(clXYZ[0], clXYZ[1], clXYZ[2]));
291}
292
295//____________________________________________________________________________
296template <class InputType>
297void ClusterFactory<InputType>::evalGlobalPosition(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
298{
299
300 int i = 0, nstat = 0;
301
302 double dist = tMaxInCm(double(clusterAnalysis.E()));
303
304 double clXYZ[3] = {0., 0., 0.}, clRmsXYZ[3] = {0., 0., 0.}, lxyzi[3], xyzi[3], wtot = 0., w = 0.;
305
306 for (auto iInput : inputsIndices) {
307
308 // get the local coordinates of the cell
309 try {
310 mGeomPtr->RelPosCellInSModule(mInputsContainer[iInput].getTower(), dist).GetCoordinates(lxyzi[0], lxyzi[1], lxyzi[2]);
311 } catch (InvalidCellIDException& e) {
312 LOG(error) << e.what();
313 continue;
314 }
315
316 // Now get the global coordinate
317 mGeomPtr->GetGlobal(lxyzi, xyzi, mGeomPtr->GetSuperModuleNumber(mInputsContainer[iInput].getTower()));
318
319 if (mLogWeight > 0.0) {
320 w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
321 } else {
322 w = mInputsContainer[iInput].getEnergy(); // just energy
323 }
324
325 if (w > 0.0) {
326 wtot += w;
327 nstat++;
328
329 for (i = 0; i < 3; i++) {
330 clXYZ[i] += (w * xyzi[i]);
331 clRmsXYZ[i] += (w * xyzi[i] * xyzi[i]);
332 }
333 }
334 }
335
336 // cout << " wtot " << wtot << endl;
337
338 if (wtot > 0) {
339 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
340 for (i = 0; i < 3; i++) {
341 clXYZ[i] /= wtot;
342
343 if (nstat > 1) {
344 clRmsXYZ[i] /= (wtot * wtot);
345 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i] * clXYZ[i];
346
347 if (clRmsXYZ[i] > 0.0) {
348 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
349 } else {
350 clRmsXYZ[i] = 0;
351 }
352 } else {
353 clRmsXYZ[i] = 0;
354 }
355 }
356 } else {
357 for (i = 0; i < 3; i++) {
358 clXYZ[i] = clRmsXYZ[i] = -1.;
359 }
360 }
361
362 clusterAnalysis.setGlobalPosition(math_utils::Point3D<float>(clXYZ[0], clXYZ[1], clXYZ[2]));
363}
364
367//____________________________________________________________________________
368template <class InputType>
369void ClusterFactory<InputType>::evalLocalPositionFit(double deff, double mLogWeight,
370 double phiSlope, gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
371{
372 int i = 0, nstat = 0;
373 double clXYZ[3] = {0., 0., 0.}, clRmsXYZ[3] = {0., 0., 0.}, xyzi[3], wtot = 0., w = 0.;
374
375 for (auto iInput : inputsIndices) {
376
377 try {
378 mGeomPtr->RelPosCellInSModule(mInputsContainer[iInput].getTower(), deff).GetCoordinates(xyzi[0], xyzi[1], xyzi[2]);
379 } catch (InvalidCellIDException& e) {
380 LOG(error) << e.what();
381 continue;
382 }
383
384 if (mLogWeight > 0.0) {
385 w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
386 } else {
387 w = mInputsContainer[iInput].getEnergy(); // just energy
388 }
389
390 if (w > 0.0) {
391 wtot += w;
392 nstat++;
393
394 for (i = 0; i < 3; i++) {
395 clXYZ[i] += (w * xyzi[i]);
396 clRmsXYZ[i] += (w * xyzi[i] * xyzi[i]);
397 }
398 }
399 } // loop
400
401 // cout << " wtot " << wtot << endl;
402
403 if (wtot > 0) {
404 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
405 for (i = 0; i < 3; i++) {
406 clXYZ[i] /= wtot;
407
408 if (nstat > 1) {
409 clRmsXYZ[i] /= (wtot * wtot);
410 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i] * clXYZ[i];
411
412 if (clRmsXYZ[i] > 0.0) {
413 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
414 } else {
415 clRmsXYZ[i] = 0;
416 }
417 } else {
418 clRmsXYZ[i] = 0;
419 }
420 }
421 } else {
422 for (i = 0; i < 3; i++) {
423 clXYZ[i] = clRmsXYZ[i] = -1.;
424 }
425 }
426
427 // clRmsXYZ[i] ??
428
429 if (phiSlope != 0.0 && mLogWeight > 0.0 && wtot) {
430 // Correction in phi direction (y - coords here); Aug 16;
431 // May be put to global level or seperate method
432 double ycorr = clXYZ[1] * (1. + phiSlope);
433
434 // printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope);
435 clXYZ[1] = ycorr;
436 }
437
438 clusterAnalysis.setLocalPosition(math_utils::Point3D<float>(clXYZ[0], clXYZ[1], clXYZ[2]));
439}
440
445//_____________________________________________________________________________
446template <class InputType>
447void ClusterFactory<InputType>::getDeffW0(const double esum, double& deff, double& w0)
448{
449 double e = 0.0;
450 const double kdp0 = 9.25147, kdp1 = 1.16700; // Hard coded now
451 const double kwp0 = 4.83713, kwp1 = -2.77970e-01, kwp2 = 4.41116;
452
453 // No extrapolation here
454 e = esum < 0.5 ? 0.5 : esum;
455 e = e > 100. ? 100. : e;
456
457 deff = kdp0 + kdp1 * TMath::Log(e);
458 w0 = kwp0 / (1. + TMath::Exp(kwp1 * (e + kwp2)));
459}
460
468//______________________________________________________________________________
469template <class InputType>
470void ClusterFactory<InputType>::evalCoreEnergy(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
471{
472
473 float coreEnergy = 0.;
474
475 if (!clusterAnalysis.getLocalPosition().Mag2()) {
476 evalLocalPosition(inputsIndices, clusterAnalysis);
477 }
478
479 double phiPoint = clusterAnalysis.getLocalPosition().Phi();
480 double etaPoint = clusterAnalysis.getLocalPosition().Eta();
481 for (auto iInput : inputsIndices) {
482
483 auto [eta, phi] = mGeomPtr->EtaPhiFromIndex(mInputsContainer[iInput].getTower());
484 phi = phi * TMath::DegToRad();
485
486 double distance = TMath::Sqrt((eta - etaPoint) * (eta - etaPoint) + (phi - phiPoint) * (phi - phiPoint));
487
488 if (distance < mCoreRadius) {
489 coreEnergy += mInputsContainer[iInput].getEnergy();
490 }
491 }
492 clusterAnalysis.setCoreEnergy(coreEnergy);
493}
494
497//____________________________________________________________________________
498template <class InputType>
499void ClusterFactory<InputType>::evalNExMax(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
500{
501 // Pre-compute cell indices and energies for all cells in cluster to avoid multiple expensive geometry lookups
502 const size_t n = inputsIndices.size();
503 std::vector<short> rows;
504 std::vector<short> columns;
505 std::vector<double> energies;
506
507 rows.reserve(n);
508 columns.reserve(n);
509 energies.reserve(n);
510
511 for (auto iInput : inputsIndices) {
512 auto [nSupMod, nModule, nIphi, nIeta] = mGeomPtr->GetCellIndex(mInputsContainer[iInput].getTower());
513
514 // get a nice topological indexing that is done in exactly the same way as used by the clusterizer
515 // this way we can handle the shared cluster cases correctly
516 const auto [row, column] = mGeomPtr->GetTopologicalRowColumn(nSupMod, nModule, nIphi, nIeta);
517
518 rows.push_back(row);
519 columns.push_back(column);
520 energies.push_back(mInputsContainer[iInput].getEnergy());
521 }
522
523 // Now find local maxima using pre-computed data
524 int nExMax = 0;
525 for (size_t i = 0; i < n; i++) {
526 // this cell is assumed to be local maximum unless we find a higher energy cell in the neighborhood
527 bool isExMax = true;
528
529 // loop over all other cells in cluster
530 for (size_t j = 0; j < n; j++) {
531 if (i == j) {
532 continue;
533 }
534
535 // adjacent cell is any cell with adjacent phi or eta index
536 if (std::abs(rows[i] - rows[j]) <= 1 &&
537 std::abs(columns[i] - columns[j]) <= 1) {
538
539 // if there is a cell with higher energy than the current cell, it is not a local maximum
540 if (energies[j] > energies[i]) {
541 isExMax = false;
542 break;
543 }
544 }
545 }
546 if (isExMax) {
547 nExMax++;
548 }
549 }
550 clusterAnalysis.setNExMax(nExMax);
551}
552
556//____________________________________________________________________________
557template <class InputType>
558void ClusterFactory<InputType>::evalElipsAxis(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
559{
560 double wtot = 0.;
561 double x = 0.;
562 double z = 0.;
563 double dxx = 0.;
564 double dzz = 0.;
565 double dxz = 0.;
566
567 std::array<float, 2> lambda;
568
569 for (auto iInput : inputsIndices) {
570
571 auto [nSupMod, nModule, nIphi, nIeta] = mGeomPtr->GetCellIndex(mInputsContainer[iInput].getTower());
572 auto [iphi, ieta] = mGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta);
573
574 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
575 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
576 if (mSharedCluster && nSupMod % 2) {
577 ieta += EMCAL_COLS;
578 }
579
580 double etai = (double)ieta;
581 double phii = (double)iphi;
582
583 double w = TMath::Max(0., mLogWeight + TMath::Log(mInputsContainer[iInput].getEnergy() / clusterAnalysis.E()));
584 // clusterAnalysis.E() summed amplitude of inputs, i.e. energy of cluster
585 // Gives smaller value of lambda than log weight
586 // w = mEnergyList[iInput] / clusterAnalysis.E(); // Nov 16, 2006 - try just energy
587
588 dxx += w * etai * etai;
589 x += w * etai;
590 dzz += w * phii * phii;
591 z += w * phii;
592
593 dxz += w * etai * phii;
594
595 wtot += w;
596 }
597
598 if (wtot > 0) {
599 dxx /= wtot;
600 x /= wtot;
601 dxx -= x * x;
602 dzz /= wtot;
603 z /= wtot;
604 dzz -= z * z;
605 dxz /= wtot;
606 dxz -= x * z;
607
608 lambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt(0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz);
609
610 if (lambda[0] > 0) {
611 lambda[0] = TMath::Sqrt(lambda[0]);
612 } else {
613 lambda[0] = 0;
614 }
615
616 lambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt(0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz);
617
618 if (lambda[1] > 0) { // To avoid exception if numerical errors lead to negative lambda.
619 lambda[1] = TMath::Sqrt(lambda[1]);
620 } else {
621 lambda[1] = 0.;
622 }
623 } else {
624 lambda[0] = 0.;
625 lambda[1] = 0.;
626 }
627
628 clusterAnalysis.setM02(lambda[0] * lambda[0]);
629 clusterAnalysis.setM20(lambda[1] * lambda[1]);
630}
631
634//____________________________________________________________________________
635template <class InputType>
636std::tuple<int, float, float, bool> ClusterFactory<InputType>::getMaximalEnergyIndex(gsl::span<const int> inputsIndices) const
637{
638
639 float energy = 0.;
640 int mid = 0;
641 float cellAmp = 0;
642 int iSupMod0 = -1;
643 bool shared = false;
644 for (auto iInput : inputsIndices) {
645 if (iInput >= mInputsContainer.size()) {
646 throw CellIndexRangeException(iInput, mInputsContainer.size());
647 }
648 cellAmp += mInputsContainer[iInput].getEnergy();
649 if (iSupMod0 == -1) {
650 iSupMod0 = mGeomPtr->GetSuperModuleNumber(mInputsContainer[iInput].getTower());
651 } else if (iSupMod0 != mGeomPtr->GetSuperModuleNumber(mInputsContainer[iInput].getTower())) {
652 shared = true;
653 }
654 if (mInputsContainer[iInput].getEnergy() > energy) {
655 energy = mInputsContainer[iInput].getEnergy();
656 mid = iInput;
657 }
658 } // loop on cluster inputs
659
660 return std::make_tuple(mid, energy, cellAmp, shared);
661}
662
665//____________________________________________________________________________
666template <class InputType>
667bool ClusterFactory<InputType>::isExoticCell(short towerId, float ecell, float const exoticTime, float& fCross) const
668{
669 if (ecell < mExoticCellMinAmplitude) {
670 return false; // do not reject low energy cells
671 }
672
673 // if the look up table is not set yet (mostly due to a reset call) then set it up now.
674 if (!getLookUpInit()) {
676 }
677
678 float eCross = getECross(towerId, ecell, exoticTime);
679 fCross = 1.f - eCross / ecell;
680
681 if (fCross > mExoticCellFraction) {
682 LOG(debug) << "EXOTIC CELL id " << towerId << ", eCell " << ecell << ", eCross " << eCross << ", 1-eCross/eCell " << 1 - eCross / ecell;
683 return true;
684 }
685
686 return false;
687}
688
691//____________________________________________________________________________
692template <class InputType>
693float ClusterFactory<InputType>::getECross(short towerId, float energy, float const exoticTime) const
694{
695 auto [iSM, iMod, iIphi, iIeta] = mGeomPtr->GetCellIndex(towerId);
696 auto [iphi, ieta] = mGeomPtr->GetCellPhiEtaIndexInSModule(iSM, iMod, iIphi, iIeta);
697
698 // Get close cells index, energy and time, not in corners
699
700 short towerId1 = -1;
701 short towerId2 = -1;
702
703 if (iphi < o2::emcal::EMCAL_ROWS - 1) {
704 try {
705 towerId1 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi + 1, ieta);
706 } catch (InvalidCellIDException& e) {
707 towerId1 = -1 * e.getCellID();
708 }
709 }
710 if (iphi > 0) {
711 try {
712 towerId2 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi - 1, ieta);
713 } catch (InvalidCellIDException& e) {
714 towerId2 = -1 * e.getCellID();
715 }
716 }
717
718 // In case of cell in eta = 0 border, depending on SM shift the cross cell index
719
720 short towerId3 = -1;
721 short towerId4 = -1;
722
723 if (ieta == o2::emcal::EMCAL_COLS - 1 && !(iSM % 2)) {
724 try {
725 towerId3 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM + 1, iphi, 0);
726 } catch (InvalidCellIDException& e) {
727 towerId3 = -1 * e.getCellID();
728 }
729 try {
730 towerId4 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi, ieta - 1);
731 } catch (InvalidCellIDException& e) {
732 towerId4 = -1 * e.getCellID();
733 }
734 } else if (ieta == 0 && iSM % 2) {
735 try {
736 towerId3 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi, ieta + 1);
737 } catch (InvalidCellIDException& e) {
738 towerId3 = -1 * e.getCellID();
739 }
740 try {
741 towerId4 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM - 1, iphi, o2::emcal::EMCAL_COLS - 1);
742 } catch (InvalidCellIDException& e) {
743 towerId4 = -1 * e.getCellID();
744 }
745 } else {
746 if (ieta < o2::emcal::EMCAL_COLS - 1) {
747 try {
748 towerId3 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi, ieta + 1);
749 } catch (InvalidCellIDException& e) {
750 towerId3 = -1 * e.getCellID();
751 }
752 }
753 if (ieta > 0) {
754 try {
755 towerId4 = mGeomPtr->GetAbsCellIdFromCellIndexes(iSM, iphi, ieta - 1);
756 } catch (InvalidCellIDException& e) {
757 towerId4 = -1 * e.getCellID();
758 }
759 }
760 }
761
762 LOG(debug) << "iSM " << iSM << ", towerId " << towerId << ", a " << towerId1 << ", b " << towerId2 << ", c " << towerId3 << ", e " << towerId3;
763
764 short index1 = (towerId1 > -1) ? mLoolUpTowerToIndex.at(towerId1) : -1;
765 short index2 = (towerId2 > -1) ? mLoolUpTowerToIndex.at(towerId2) : -1;
766 short index3 = (towerId3 > -1) ? mLoolUpTowerToIndex.at(towerId3) : -1;
767 short index4 = (towerId4 > -1) ? mLoolUpTowerToIndex.at(towerId4) : -1;
768
769 std::array<std::pair<float, float>, 4> cellData = {
770 {{(index1 > -1) ? mInputsContainer[index1].getEnergy() : 0., (index1 > -1) ? mInputsContainer[index1].getTimeStamp() : 0.},
771 {(index2 > -1) ? mInputsContainer[index2].getEnergy() : 0., (index2 > -1) ? mInputsContainer[index2].getTimeStamp() : 0.},
772 {(index3 > -1) ? mInputsContainer[index3].getEnergy() : 0., (index3 > -1) ? mInputsContainer[index3].getTimeStamp() : 0.},
773 {(index4 > -1) ? mInputsContainer[index4].getEnergy() : 0., (index4 > -1) ? mInputsContainer[index4].getTimeStamp() : 0.}}};
774
775 for (auto& cell : cellData) {
776 if (std::abs(exoticTime - cell.second) > mExoticCellDiffTime) {
777 cell.first = 0;
778 }
779 }
780
781 float w1 = 1, w2 = 1, w3 = 1, w4 = 1;
782 if (mUseWeightExotic) {
783 w1 = GetCellWeight(cellData[0].first, energy);
784 w2 = GetCellWeight(cellData[1].first, energy);
785 w3 = GetCellWeight(cellData[2].first, energy);
786 w4 = GetCellWeight(cellData[3].first, energy);
787 }
788
789 if (cellData[0].first < mExoticCellInCrossMinAmplitude || w1 <= 0) {
790 cellData[0].first = 0;
791 }
792 if (cellData[1].first < mExoticCellInCrossMinAmplitude || w2 <= 0) {
793 cellData[1].first = 0;
794 }
795 if (cellData[2].first < mExoticCellInCrossMinAmplitude || w3 <= 0) {
796 cellData[2].first = 0;
797 }
798 if (cellData[3].first < mExoticCellInCrossMinAmplitude || w4 <= 0) {
799 cellData[3].first = 0;
800 }
801
802 return cellData[0].first + cellData[1].first + cellData[2].first + cellData[3].first;
803}
804
807//____________________________________________________________________________
808template <class InputType>
809float ClusterFactory<InputType>::GetCellWeight(float eCell, float eCluster) const
810{
811 if (eCell > 0 && eCluster > 0) {
812 if (mLogWeight > 0) {
813 return std::max(0.f, mLogWeight + std::log(eCell / eCluster));
814 } else {
815 return std::log(eCluster / eCell);
816 }
817 } else {
818 return 0.;
819 }
820}
821
824//____________________________________________________________________________
825template <class InputType>
826int ClusterFactory<InputType>::getMultiplicityAtLevel(float H, gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
827{
828 int multipl = 0;
829 for (auto iInput : inputsIndices) {
830 if (mInputsContainer[iInput].getEnergy() > H * clusterAnalysis.E()) {
831 multipl++;
832 }
833 }
834
835 return multipl;
836}
837
840//____________________________________________________________________________
841template <class InputType>
842void ClusterFactory<InputType>::evalTime(gsl::span<const int> inputsIndices, AnalysisCluster& clusterAnalysis) const
843{
844 float maxE = 0;
845 unsigned short maxAt = 0;
846 for (auto iInput : inputsIndices) {
847 if (mInputsContainer[iInput].getEnergy() > maxE) {
848 maxE = mInputsContainer[iInput].getEnergy();
849 maxAt = iInput;
850 }
851 }
852
853 clusterAnalysis.setClusterTime(mInputsContainer[maxAt].getTimeStamp());
854}
855
859//_____________________________________________________________________
860template <class InputType>
861double ClusterFactory<InputType>::tMaxInCm(const double e, const int key) const
862{
863 const double ca = 4.82; // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
864 double tmax = 0.; // position of electromagnetic shower max in cm
865
866 const double x0 = 1.31; // radiation lenght (cm)
867
868 if (e > 0.1) {
869 tmax = TMath::Log(e) + ca;
870 if (key == 0) {
871 tmax += 0.5;
872 } else {
873 tmax -= 0.5;
874 }
875 tmax *= x0; // convert to cm
876 }
877
878 return tmax;
879}
880
883//______________________________________________________________________________
884template <class InputType>
886{
887 return (2. * TMath::ATan(TMath::Exp(-arg)));
888}
889
892//______________________________________________________________________________
893template <class InputType>
895{
896 return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
897}
898
899template <class InputType>
900ClusterFactory<InputType>::ClusterIterator::ClusterIterator(const ClusterFactory& factory, int clusterIndex, bool forward) : mClusterFactory(factory),
901 mCurrentCluster(),
902 mClusterID(clusterIndex),
903 mForward(forward)
904{
905 mCurrentCluster = mClusterFactory.buildCluster(mClusterID);
906}
907
908template <class InputType>
910{
911 return &mClusterFactory == &rhs.mClusterFactory && mClusterID == rhs.mClusterID && mForward == rhs.mForward;
912}
913
914template <class InputType>
916{
917 if (mForward) {
918 mClusterID++;
919 } else {
920 mClusterID--;
921 }
922 mCurrentCluster = mClusterFactory.buildCluster(mClusterID);
923 return *this;
924}
925
926template <class InputType>
928{
929 auto tmp = *this;
930 ++(*this);
931 return tmp;
932}
933
934template <class InputType>
936{
937 if (mForward) {
938 mClusterID--;
939 } else {
940 mClusterID++;
941 }
942 mCurrentCluster = mClusterFactory.buildCluster(mClusterID);
943 return *this;
944}
945
946template <class InputType>
948{
949 auto tmp = *this;
950 --(*this);
951 return tmp;
952}
953
std::ostringstream debug
int32_t i
uint32_t j
Definition RawData.h:0
StringRef key
Cluster class for kinematic cluster parametersported from AliVCluster in AliRoot.
void setNExMax(unsigned char nExMax)
void setGlobalPosition(math_utils::Point3D< float > x)
Set the cluster global position.
void setFCross(float fCross)
void setCoreEnergy(float energy)
math_utils::Point3D< float > getLocalPosition() const
void setCellsIndices(const std::vector< unsigned short > &array)
Set the array of cell indices.
void setIndMaxInput(const int ind)
void setClusterTime(float time)
void setLocalPosition(math_utils::Point3D< float > x)
ClusterIterator & operator--()
Prefix decrementation operator.
ClusterIterator(const ClusterFactory &factory, int clusterIndex, bool forward)
Constructor, initializing the iterator.
bool operator==(const ClusterIterator &rhs) const
Check for equalness.
ClusterIterator & operator++()
Prefix incrementation operator.
Exception handling uninitialized look up table.
const char * what() const noexcept final
Access to error message of the exception.
EMCal clusters factory Ported from class AliEMCALcluster.
void reset()
Reset containers.
std::tuple< int, float, float, bool > getMaximalEnergyIndex(gsl::span< const int > inputsIndices) const
Finds the maximum energy in the cluster and computes the Summed amplitude of digits/cells.
float GetCellWeight(float eCell, float eCluster) const
return weight of cell for shower shape calculation
float thetaToEta(float arg) const
Converts Theta (Radians) to Eta (Radians)
void evalNExMax(gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
Calculate the number of local maxima in the cluster.
float getECross(short absID, float energy, float const exoticTime) const
Calculate the energy in the cross around the energy of a given cell.
void evalGlobalPosition(gsl::span< const int > inputsIndices, AnalysisCluster &cluster) const
Calculates the center of gravity in the global ALICE coordinates.
void evalLocalPosition(gsl::span< const int > inputsIndices, AnalysisCluster &cluster) const
Calculates the center of gravity in the local EMCAL-module coordinates.
Double_t tMaxInCm(const Double_t e=0.0, const int key=0) const
static void getDeffW0(const Double_t esum, Double_t &deff, Double_t &w0)
void evalTime(gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
Time is set to the time of the digit with the maximum energy.
int getMultiplicityAtLevel(float level, gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
Calculates the multiplicity of digits/cells with energy larger than level*energy.
ClusterFactory()=default
Dummy constructor.
void evalLocalPositionFit(Double_t deff, Double_t w0, Double_t phiSlope, gsl::span< const int > inputsIndices, AnalysisCluster &cluster) const
evaluates local position of clusters in SM
void evalElipsAxis(gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
void evalCoreEnergy(gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
void evalDispersion(gsl::span< const int > inputsIndices, AnalysisCluster &clusterAnalysis) const
float etaToTheta(float arg) const
Converts Eta (Radians) to Theta (Radians)
bool isExoticCell(short towerId, float ecell, float const exoticTime, float &fCross) const
Look to cell neighbourhood and reject if it seems exotic.
AnalysisCluster buildCluster(int index, o2::emcal::ClusterLabel *clusterLabel=nullptr) const
evaluates cluster parameters: position, shower shape, primaries ...
cluster class for MC particle IDs and their respective energy fraction
void orderLabels()
Sort the labels and energy fraction in descending order (largest energy fraction to smallest)
void normalize(float factor)
Normalize the energy fraction.
void addValue(int label, float energyFraction)
Add label and energy fraction to the.
Exception handling non-existing cell IDs.
int getCellID() const noexcept
Access to cell ID raising the exception.
const char * what() const noexcept final
Access to error message of the exception.
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLsizei GLsizei GLfloat distance
Definition glcorearb.h:5506
GLuint GLfloat x0
Definition glcorearb.h:5034
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
@ EMCAL_ROWS
Number of rows per module for EMCAL.
Definition Constants.h:25
@ EMCAL_COLS
Number of columns per module for EMCAL.
Definition Constants.h:26
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< int > row
std::vector< ReadoutWindowData > rows