Project
Loading...
Searching...
No Matches
SACFactorization.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
16
17#ifndef ALICEO2_SACFACTORIZATION_H_
18#define ALICEO2_SACFACTORIZATION_H_
19
20#include <vector>
21#include "Rtypes.h"
23#include "DataFormatsTPC/Defs.h"
24#include <boost/property_tree/ptree.hpp>
25
26namespace o2::tpc
27{
28
30{
31 public:
33
36 SACFactorization(const unsigned int timeFrames) : mTimeFrames{timeFrames} {};
37
39 SACFactorization() = default;
40
44 void factorizeSACs();
45
47 void calcSACZero();
48
50 void calcSACOne();
51
53 void calcSACDelta();
54
56 void createStatusMap();
57
61 auto getSACValue(const unsigned int stack, const unsigned int interval) const { return mSACs[stack][interval]; }
62
65 float getSACZeroVal(const unsigned int stack) const { return mSACZero.getValueIDCZero(getSide(stack), stack % GEMSTACKSPERSIDE); }
66
70 float getSACDeltaVal(const unsigned int stack, unsigned int interval) const { return mSACDelta.getValue(getSide(stack), getSACDeltaIndex(stack, interval)); }
71
75 float getSACOneVal(const Side side, unsigned int interval) const { return mSACOne.getValue(side, interval); }
76
80 static unsigned int getSACDeltaIndex(const unsigned int stack, unsigned int interval) { return stack % GEMSTACKSPERSIDE + GEMSTACKSPERSIDE * interval; }
81
83 unsigned int getNTimeframes() const { return mTimeFrames; }
84
86 unsigned long getNintervals(const int stack = 0) const { return mSACs[stack].size(); }
87
89 unsigned long getNIntegrationIntervals(const int stack = 0) const { return mSACs[stack].size(); }
90
93 const std::vector<float>& getSACZero(const o2::tpc::Side side) const { return mSACZero.mSACZero[side].mIDCZero; }
94
97 const std::vector<float>& getSACOne(const o2::tpc::Side side) const { return mSACOne.mSACOne[side].mIDCOne; }
98
101 const std::vector<float>& getSACDeltaUncompressed(const o2::tpc::Side side) const { return mSACDelta.mSACDelta[side].getIDCDelta(); }
102
104 const auto& getSACDeltaUncompressed() const& { return mSACDelta; }
105
107 auto getSACDeltaUncompressed() && { return std::move(mSACDelta); }
108
111
114
116 const auto& getSACZero() const { return mSACZero; }
117
119 const auto& getSACOne() const& { return mSACOne; }
120
122 auto getSACOne() && { return std::move(mSACOne); }
123
125 const auto& getSACs() const { return mSACs; }
126
128 const auto& getOutlierMap() const { return mOutlierMap; }
129
131 static int getNThreads() { return sNThreads; }
132
137 void setSACs(std::vector<int32_t>&& SACs, const unsigned int stack) { mSACs[stack] = std::move(SACs); }
138
141 static void setNThreads(const int nThreads) { sNThreads = nThreads; }
142
144 static void setMinCompressedSACDelta(const float minSACDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCSACCompressionParam", "minSACDeltaValue", minSACDeltaValue); }
145
147 static void setMaxCompressedSACDelta(const float maxSACDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCSACCompressionParam", "maxSACDeltaValue", maxSACDeltaValue); }
148
153 void drawSACsSector(const unsigned int sector, const unsigned int interval, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACsSector.pdf") const { drawSACHelper(false, Sector(sector), interval, filename, minZ, maxZ); }
154
158 void drawSACZeroSector(const unsigned int sector, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACZeroSector.pdf") const { drawSACZeroHelper(false, Sector(sector), filename, minZ, maxZ); }
159
163 void drawSACZeroOutlierSector(const unsigned int sector, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACZeroOutlierSector.pdf") const { drawSACZeroOutlierHelper(false, Sector(sector), filename, minZ, maxZ); }
164
170 void drawSACDeltaSector(const unsigned int sector, const unsigned int interval, const float minZ = 0, const float maxZ = -1, const SACDeltaCompression compression = SACDeltaCompression::NO, const std::string filename = "SACDeltaSector.pdf") const { drawSACDeltaHelper(false, Sector(sector), interval, compression, filename, minZ, maxZ); }
171
176 void drawSACsSide(const o2::tpc::Side side, const unsigned int interval, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACsSide.pdf") const { drawSACHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), interval, filename, minZ, maxZ); }
177
181 void drawSACZeroSide(const o2::tpc::Side side, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACZeroSide.pdf") const { drawSACZeroHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), filename, minZ, maxZ); }
182
186 void drawSACZeroOutlierSide(const o2::tpc::Side side, const float minZ = 0, const float maxZ = -1, const std::string filename = "SACZeroOutlierSide.pdf") const { drawSACZeroOutlierHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), filename, minZ, maxZ); }
187
193 void drawSACDeltaSide(const o2::tpc::Side side, const unsigned int interval, const float minZ = 0, const float maxZ = -1, const SACDeltaCompression compression = SACDeltaCompression::NO, const std::string filename = "SACDeltaSide.pdf") const { drawSACDeltaHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), interval, compression, filename, minZ, maxZ); }
194
198 void dumpToFile(const char* outFileName = "SACFactorized.root", const char* outName = "SACFactorized") const;
199
202 void dumpToTree(int intervals = -1, const char* outFileName = "SACTree.root") const;
203
205 void reset();
206
208 static Side getSide(const unsigned int gemStack) { return (gemStack < GEMSTACKSPERSIDE) ? Side::A : Side::C; }
209
211 static unsigned int getStack(const unsigned int sector, const unsigned int stack) { return static_cast<unsigned int>(stack + sector * GEMSTACKSPERSECTOR); }
212
214 static unsigned int getStackInSide(const unsigned int sector, const unsigned int stack) { return getStack(sector, stack) % (GEMSTACKS / 2); }
215
216 void setSACZero(const SACZero& sacZero) { mSACZero = sacZero; }
217
218 private:
219 const unsigned int mTimeFrames{};
220 std::array<std::vector<int32_t>, o2::tpc::GEMSTACKS> mSACs{};
221 SACZero mSACZero{};
222 SACOne mSACOne{};
223 std::array<int, GEMSTACKS> mOutlierMap{};
224 SACDelta<float> mSACDelta{};
225 inline static int sNThreads{1};
226
228 void drawSACDeltaHelper(const bool type, const Sector sector, const unsigned int interval, const SACDeltaCompression compression, const std::string filename, const float minZ, const float maxZ) const;
229
231 void drawSACHelper(const bool type, const Sector sector, const unsigned int interval, const std::string filename, const float minZ, const float maxZ) const;
232
234 void drawSACZeroHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const;
235
237 void drawSACZeroOutlierHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const;
238
239 ClassDefNV(SACFactorization, 1)
240};
241
242} // namespace o2::tpc
243
244#endif
This file provides the structs for storing the factorized IDC values and fourier coefficients to be s...
uint32_t side
Definition RawData.h:0
uint32_t stack
Definition RawData.h:1
static SACDelta< DataT > getCompressedSACs(const SACDelta< float > &sacDeltaUncompressed)
const std::vector< float > & getSACZero(const o2::tpc::Side side) const
const auto & getSACOne() const &
void reset()
resetting aggregated SACs
void drawSACDeltaSector(const unsigned int sector, const unsigned int interval, const float minZ=0, const float maxZ=-1, const SACDeltaCompression compression=SACDeltaCompression::NO, const std::string filename="SACDeltaSector.pdf") const
const auto & getSACZero() const
static unsigned int getStackInSide(const unsigned int sector, const unsigned int stack)
float getSACDeltaVal(const unsigned int stack, unsigned int interval) const
void setSACZero(const SACZero &sacZero)
void drawSACZeroOutlierSector(const unsigned int sector, const float minZ=0, const float maxZ=-1, const std::string filename="SACZeroOutlierSector.pdf") const
static unsigned int getStack(const unsigned int sector, const unsigned int stack)
unsigned long getNintervals(const int stack=0) const
void setSACs(std::vector< int32_t > &&SACs, const unsigned int stack)
void drawSACDeltaSide(const o2::tpc::Side side, const unsigned int interval, const float minZ=0, const float maxZ=-1, const SACDeltaCompression compression=SACDeltaCompression::NO, const std::string filename="SACDeltaSide.pdf") const
void dumpToTree(int intervals=-1, const char *outFileName="SACTree.root") const
void drawSACZeroSector(const unsigned int sector, const float minZ=0, const float maxZ=-1, const std::string filename="SACZeroSector.pdf") const
void drawSACsSector(const unsigned int sector, const unsigned int interval, const float minZ=0, const float maxZ=-1, const std::string filename="SACsSector.pdf") const
void calcSACDelta()
calculate \Delta I(r,\phi,t) = I(r,\phi,t) / ( I_0(r,\phi) * I_1(t) )
void createStatusMap()
create outlier map using the SAC0
void calcSACOne()
calculate I_1(t) = <I(r,\phi,t) / I_0(r,\phi)>_{r,\phi}
unsigned long getNIntegrationIntervals(const int stack=0) const
SACFactorization()=default
default constructor
void drawSACsSide(const o2::tpc::Side side, const unsigned int interval, const float minZ=0, const float maxZ=-1, const std::string filename="SACsSide.pdf") const
static void setMinCompressedSACDelta(const float minSACDeltaValue)
float getSACZeroVal(const unsigned int stack) const
const std::vector< float > & getSACDeltaUncompressed(const o2::tpc::Side side) const
const auto & getSACs() const
static Side getSide(const unsigned int gemStack)
const auto & getSACDeltaUncompressed() const &
void drawSACZeroOutlierSide(const o2::tpc::Side side, const float minZ=0, const float maxZ=-1, const std::string filename="SACZeroOutlierSide.pdf") const
auto getSACDeltaMediumCompressed() const
const std::vector< float > & getSACOne(const o2::tpc::Side side) const
auto getSACValue(const unsigned int stack, const unsigned int interval) const
unsigned int getNTimeframes() const
void calcSACZero()
calculate I_0(r,\phi) = <I(r,\phi,t)>_t
static unsigned int getSACDeltaIndex(const unsigned int stack, unsigned int interval)
const auto & getOutlierMap() const
float getSACOneVal(const Side side, unsigned int interval) const
SACFactorization(const unsigned int timeFrames)
IDCDeltaCompression SACDeltaCompression
void dumpToFile(const char *outFileName="SACFactorized.root", const char *outName="SACFactorized") const
static void setMaxCompressedSACDelta(const float maxSACDeltaValue)
static void setNThreads(const int nThreads)
void drawSACZeroSide(const o2::tpc::Side side, const float minZ=0, const float maxZ=-1, const std::string filename="SACZeroSide.pdf") const
static constexpr int MAXSECTOR
Definition Sector.h:44
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat maxZ
Definition glcorearb.h:2910
GLfloat GLfloat minZ
Definition glcorearb.h:2910
Global TPC definitions and constants.
Definition SimTraits.h:167
constexpr unsigned short GEMSTACKS
Definition Defs.h:60
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
@ C
Definition Defs.h:36
IDCDeltaCompression
IDC Delta IDC Compression types.
constexpr unsigned short GEMSTACKSPERSECTOR
Definition Defs.h:57
constexpr unsigned short GEMSTACKSPERSIDE
Definition Defs.h:59
std::string filename()
float getValue(const Side side, const int index) const
std::array< IDCDelta< DataT >, SIDES > mSACDelta
std::array< IDCOne, SIDES > mSACOne
float getValue(const Side side, const int interval) const
float getValueIDCZero(const Side side, const int stackInSector) const
std::array< IDCZero, SIDES > mSACZero