Project
Loading...
Searching...
No Matches
IDCFactorization.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
17
18#ifndef ALICEO2_IDCFACTORIZATION_H_
19#define ALICEO2_IDCFACTORIZATION_H_
20
21#include <vector>
22#include "Rtypes.h"
23#include "TPCBase/Mapper.h"
26#include "DataFormatsTPC/Defs.h"
27#include <boost/property_tree/ptree.hpp>
28
29namespace o2::tpc
30{
31
32template <class T>
33class CalDet;
34
36 std::array<std::vector<float>, CRU::MaxCRU> idcs{};
37 int relTF{};
38 ClassDefNV(IDCFactorizeSplit, 1)
39};
40
42{
43 public:
51 IDCFactorization(const std::array<unsigned char, Mapper::NREGIONS>& groupPads, const std::array<unsigned char, Mapper::NREGIONS>& groupRows, const std::array<unsigned char, Mapper::NREGIONS>& groupLastRowsThreshold, const std::array<unsigned char, Mapper::NREGIONS>& groupLastPadsThreshold, const unsigned int groupNotnPadsSectorEdges, const unsigned int timeFrames, const unsigned int timeframesDeltaIDC, const std::vector<uint32_t>& crus);
52
56 IDCFactorization(const unsigned int timeFrames, const unsigned int timeframesDeltaIDC, const std::vector<uint32_t>& crus) : IDCFactorization(std::array<unsigned char, Mapper::NREGIONS>{1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, std::array<unsigned char, Mapper::NREGIONS>{1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, std::array<unsigned char, Mapper::NREGIONS>{1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, std::array<unsigned char, Mapper::NREGIONS>{1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, 0, timeFrames, timeframesDeltaIDC, crus){};
57
59 IDCFactorization() = default;
60
63
66
69 static std::vector<o2::tpc::Side> getSides(const std::vector<uint32_t>& crus);
70
76 void factorizeIDCs(const bool norm, const bool calcDeltas);
77
80 void calcIDCZero(const bool norm);
81
84
87 void createStatusMapOutlier(const bool debug = false);
88
90 void createStatusMap();
91
94 void checkFECs(const float maxOutliersPerFEC = 0.7);
95
100 void checkNeighbourOutliers(const int maxIter = 10, const int nOutliersNeighbours = 8);
101
103 void calcIDCOne();
104
106 void calcIDCDelta();
107
109 template <typename DataVec>
110 static void calcIDCOne(const DataVec& idcsData, const int idcsPerCRU, const int integrationIntervalOffset, const unsigned int indexOffset, const CRU cru, std::vector<std::vector<float>>& idcOneTmp, const IDCZero* idcZero, const CalDet<PadFlags>* flagMap = nullptr, const bool usePadStatusMap = false);
111
118 float getIDCValGrouped(const unsigned int sector, const unsigned int region, const unsigned int grow, unsigned int gpad, unsigned int integrationInterval) const { return mIDCs[sector * Mapper::NREGIONS + region][integrationInterval][mOffsRow[region][grow] + gpad]; }
119
126 float getIDCValUngrouped(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const;
127
133 float getIDCZeroVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad) const { return mIDCZero[mSideIndex[Sector(sector).side()]].getValueIDCZero(getIndexUngrouped(sector, region, urow, upad, 0)); }
134
142 float getIDCDeltaVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int chunk, unsigned int localintegrationInterval) const { return mIDCDelta[mSideIndex[Sector(sector).side()]][chunk].getValue(getIndexUngrouped(sector, region, urow, upad, localintegrationInterval)); }
143
148 void getLocalIntegrationInterval(const unsigned int integrationInterval, unsigned int& chunk, unsigned int& localintegrationInterval) const;
149
151 unsigned int getNTimeframes() const { return mTimeFrames; }
152
155 unsigned long getNIntegrationIntervalsInChunk(const unsigned int chunk) const;
156
159 unsigned long getNIntegrationIntervalsToChunk(const unsigned int chunk) const;
160
162 unsigned long getNIntegrationIntervals(const int cru) const;
163
165 unsigned long getNIntegrationIntervals() const;
166
169 const std::vector<float>& getIDCZeroVec(const o2::tpc::Side side) const { return mIDCZero[mSideIndex[side]].mIDCZero; }
170
173 const IDCZero& getIDCZero(const o2::tpc::Side side) const& { return mIDCZero[mSideIndex[side]]; }
174
177 auto getIDCZero(const o2::tpc::Side side) && { return mIDCZero[mSideIndex[side]]; }
178
181 const std::vector<float>& getIDCOneVec(const o2::tpc::Side side) const { return mIDCOne[mSideIndex[side]].mIDCOne; }
182
185 const IDCOne& getIDCOne(const o2::tpc::Side side) const { return mIDCOne[mSideIndex[side]]; }
186
190 const std::vector<float>& getIDCDeltaValuesUncompressed(const unsigned int chunk, const o2::tpc::Side side) const { return mIDCDelta[mSideIndex[side]][chunk].getIDCDelta(); }
191
194 const auto& getIDCDeltaUncompressed(const unsigned int chunk, const Side side) const& { return mIDCDelta[mSideIndex[side]][chunk]; }
195
197 auto getIDCDeltaUncompressed(const unsigned int chunk, const Side side) && { return std::move(mIDCDelta[mSideIndex[side]][chunk]); }
198
201 auto getIDCDeltaMediumCompressed(const unsigned int chunk, const Side side) const { return IDCDeltaCompressionHelper<unsigned short>::getCompressedIDCs(mIDCDelta[mSideIndex[side]][chunk]); }
202
205 auto getIDCDeltaHighCompressed(const unsigned int chunk, const Side side) const { return IDCDeltaCompressionHelper<unsigned char>::getCompressedIDCs(mIDCDelta[mSideIndex[side]][chunk]); }
206
208 unsigned int getNChunks(const Side side) const { return mIDCDelta[mSideIndex[side]].size(); }
209
211 const auto& getIDCZero() const { return mIDCZero; }
212
214 const auto& getIDCOne() const& { return mIDCOne; }
215
217 auto getIDCOne() && { return std::move(mIDCOne); }
218
220 const auto& getIDCs() const { return mIDCs; }
221
223 auto getCRUs() const { return mCRUs; }
224
226 void setTimeStamp(const long timeStamp) { mTimeStamp = timeStamp; }
227
229 long getTimeStamp() const { return mTimeStamp; }
230
232 void setRun(const int run) { mRun = run; }
233
235 int getRun() const { return mRun; }
236
237 // get number of TFs in which the DeltaIDCs are split/stored
238 unsigned int getTimeFramesDeltaIDC() const { return mTimeFramesDeltaIDC; }
239
241 static int getNThreads() { return sNThreads; }
242
247 void setIDCs(std::vector<float>&& idcs, const unsigned int cru, const unsigned int timeframe) { mIDCs[cru][timeframe] = std::move(idcs); }
248
251 static void setNThreads(const int nThreads) { sNThreads = nThreads; }
252
254 static void setMinCompressedIDCDelta(const float minIDCDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCIDCCompressionParam", "minIDCDeltaValue", minIDCDeltaValue); }
255
257 static void setMaxCompressedIDCDelta(const float maxIDCDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCIDCCompressionParam", "maxIDCDeltaValue", maxIDCDeltaValue); }
258
263 void drawIDCsSector(const unsigned int sector, const unsigned int integrationInterval, const float minZ = 0, const float maxZ = -1, const std::string filename = "IDCsSector.pdf") const { drawIDCHelper(false, Sector(sector), integrationInterval, filename, minZ, maxZ); }
264
268 void drawIDCZeroSector(const unsigned int sector, const float minZ = 0, const float maxZ = -1, const std::string filename = "IDCZeroSector.pdf") const { drawIDCZeroHelper(false, Sector(sector), filename, minZ, maxZ); }
269
275 void drawIDCDeltaSector(const unsigned int sector, const unsigned int integrationInterval, const float minZ = 0, const float maxZ = -1, const IDCDeltaCompression compression = IDCDeltaCompression::NO, const std::string filename = "IDCDeltaSector.pdf") const { drawIDCDeltaHelper(false, Sector(sector), integrationInterval, compression, filename, minZ, maxZ); }
276
281 void drawIDCsSide(const o2::tpc::Side side, const unsigned int integrationInterval, const float minZ = 0, const float maxZ = -1, const std::string filename = "IDCsSide.pdf") const { drawIDCHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), integrationInterval, filename, minZ, maxZ); }
282
285 void drawIDCsSideGIF(const unsigned int integrationIntervals = 0, const float minZ = 0, const float maxZ = -1, const int run = -1, const std::string filename = "IDCsSideGIF") const { drawIDCHelper(true, Sector(0), integrationIntervals, filename, minZ, maxZ, true, run); }
286
290 void drawIDCZeroSide(const o2::tpc::Side side, const float minZ = 0, const float maxZ = -1, const std::string filename = "IDCZeroSide.pdf") const { drawIDCZeroHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), filename, minZ, maxZ); }
291
297 void drawIDCDeltaSide(const o2::tpc::Side side, const unsigned int integrationInterval, const float minZ = 0, const float maxZ = -1, const IDCDeltaCompression compression = IDCDeltaCompression::NO, const std::string filename = "IDCDeltaSide.pdf") const { drawIDCDeltaHelper(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), integrationInterval, compression, filename, minZ, maxZ); }
298
303 void drawPadStatusFlagsMapSector(const unsigned int sector, const PadFlags flag = PadFlags::flagSkip, const std::string filename = "PadStatusFlags_Sector.pdf") const { drawPadFlagMap(false, Sector(sector), filename, flag); }
304
309 void drawPadStatusFlagsMapSide(const o2::tpc::Side side, const PadFlags flag = PadFlags::flagSkip, const std::string filename = "PadStatusFlags_Side.pdf") const { drawPadFlagMap(true, side == Side::A ? Sector(0) : Sector(Sector::MAXSECTOR - 1), filename, flag); }
310
314 void dumpToFile(const char* outFileName = "IDCFactorized.root", const char* outName = "IDCFactorized") const;
315
319 void dumpLargeObjectToFile(const char* outFileName = "IDCFactorized.root", const char* outName = "IDCFactorized") const;
320
325 static std::unique_ptr<IDCFactorization> getLargeObjectFromFile(const char* inpFileName = "IDCFactorized.root", const char* inName = "IDCFactorized");
326
328 void dumpIDCZeroToFile(const Side side, const char* outFileName = "IDCZero.root", const char* outName = "IDC0") const;
329
331 void dumpIDCOneToFile(const Side side, const char* outFileName = "IDCOne.root", const char* outName = "IDC1") const;
332
335 void dumpToTree(const Side side, const char* outFileName = "IDCTree.root");
336
340 void dumpToTreeIDC1(const float integrationTimeOrbits = 12, const char* outFileName = "IDC1Tree.root") const;
341
343 void dumpIDCDeltaToTree(const Side side, const int chunk = 0, const char* outFileName = "IDCDeltaTree.root");
344
347 std::vector<unsigned int> getIntegrationIntervalsPerTF(const int cru = -1) const;
348
350 std::vector<unsigned int> getAllIntegrationIntervalsPerTF() const;
351
353 void reset();
354
358 void setGainMap(const char* inpFile, const char* mapName);
359
361 void setGainMap(const CalDet<float>& gainmap);
362
366 void setPadFlagMap(const char* inpFile, const char* mapName);
367
369 void setPadFlagMap(const CalDet<PadFlags>& flagmap);
370
372 void setUsePadStatusMap(const bool usePadStatusMap) { mUsePadStatusMap = usePadStatusMap; }
373
375 bool getUsePadStatusMap() const { return mUsePadStatusMap; }
376
380 void dumpPadFlagMap(const char* outFile, const char* mapName);
381
383 CalDet<PadFlags>* getPadStatusMapPtr() const { return mPadFlagsMap.get(); }
384
386 std::unique_ptr<CalDet<PadFlags>> getPadStatusMap() { return std::move(mPadFlagsMap); }
387
389 const std::vector<Side>& getSides() const { return mSides; }
390
392 void setIDCZero(const Side side, const IDCZero& idcZero) { mIDCZero[mSideIndex[side]] = idcZero; }
393
395 void setIDCOne(const Side side, const IDCOne& idcOne) { mIDCOne[mSideIndex[side]] = idcOne; }
396
398 void setIDCDelta(const Side side, const IDCDelta<float>& idcDelta, const int index = 0) { mIDCDelta[mSideIndex[side]][index] = idcDelta; }
399
401 float getMeanZ(const Side side) const;
402
404 void normIDCZeroGain() { normIDCZero(0); }
405
407 void normIDCZeroStackMedian() { normIDCZero(1); }
408
410 std::array<float, o2::tpc::GEMSTACKS> getStackMedian() const;
411
412 private:
413 const unsigned int mTimeFrames{};
414 const unsigned int mTimeFramesDeltaIDC{};
415 std::array<std::vector<std::vector<float>>, CRU::MaxCRU> mIDCs{};
416 std::vector<IDCZero> mIDCZero{};
417 std::vector<IDCOne> mIDCOne{};
418 std::vector<std::vector<IDCDelta<float>>> mIDCDelta{};
419 inline static int sNThreads{1};
420 std::unique_ptr<CalDet<float>> mGainMap;
421 std::unique_ptr<CalDet<PadFlags>> mPadFlagsMap;
422 bool mInputGrouped{false};
423 bool mUsePadStatusMap{true};
424 const std::vector<uint32_t> mCRUs{};
425 std::array<unsigned int, SIDES> mSideIndex{0, 1};
426 std::vector<Side> mSides{};
427 std::vector<unsigned int> mIntegrationIntervalsPerTF{};
428 long mTimeStamp{0};
429 int mRun{0};
430
432 void drawIDCDeltaHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const IDCDeltaCompression compression, const std::string filename, const float minZ, const float maxZ) const;
433
435 void drawIDCHelper(const bool type, const Sector sector, const unsigned int integrationInterval, const std::string filename, const float minZ, const float maxZ, const bool drawGIF = false, const int run = 0) const;
436
438 void drawIDCZeroHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const;
439
441 void getTF(unsigned int integrationInterval, unsigned int& timeFrame, unsigned int& interval) const;
442
444 unsigned int getChunk(const unsigned int timeframe) const;
445
447 unsigned int getNTFsPerChunk(const unsigned int chunk) const;
448
450 void drawPadFlagMap(const bool type, const Sector sector, const std::string filename, const PadFlags flag) const;
451
453 void normIDCZero(const int type);
454
457 bool checkReceivedIDCs();
458
459 ClassDefNV(IDCFactorization, 2)
460};
461
462} // namespace o2::tpc
463
464#endif
This file provides the structs for storing the factorized IDC values and fourier coefficients to be s...
helper class for grouping of pads and rows for one sector
uint32_t side
Definition RawData.h:0
std::ostringstream debug
@ MaxCRU
Definition CRU.h:31
static IDCDelta< DataT > getCompressedIDCs(const IDCDelta< float > &idcDeltaUncompressed)
std::array< float, o2::tpc::GEMSTACKS > getStackMedian() const
void setGainMap(const char *inpFile, const char *mapName)
const std::vector< Side > & getSides() const
void calcIDCOne()
calculate I_1(t) = <I(r,\phi,t) / I_0(r,\phi)>_{r,\phi}
void dumpIDCDeltaToTree(const Side side, const int chunk=0, const char *outFileName="IDCDeltaTree.root")
dump IDCDelta to a TTree
unsigned int getNTimeframes() const
void setIDCDelta(const Side side, const IDCDelta< float > &idcDelta, const int index=0)
set the IDCDelta
void setRun(const int run)
void drawIDCZeroSide(const o2::tpc::Side side, const float minZ=0, const float maxZ=-1, const std::string filename="IDCZeroSide.pdf") const
void reset()
resetting aggregated IDCs
void fillIDCZeroDeadPads()
fill I_0 values in case of dead pads,FECs etc.
void drawPadStatusFlagsMapSide(const o2::tpc::Side side, const PadFlags flag=PadFlags::flagSkip, const std::string filename="PadStatusFlags_Side.pdf") const
const auto & getIDCOne() const &
const auto & getIDCs() const
unsigned long getNIntegrationIntervalsInChunk(const unsigned int chunk) const
void dumpToTree(const Side side, const char *outFileName="IDCTree.root")
unsigned int getTimeFramesDeltaIDC() const
IDCFactorization()=default
default constructor for ROOT I/O
float getIDCDeltaVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int chunk, unsigned int localintegrationInterval) const
void drawIDCsSideGIF(const unsigned int integrationIntervals=0, const float minZ=0, const float maxZ=-1, const int run=-1, const std::string filename="IDCsSideGIF") const
IDCFactorization(IDCFactorization &&)=default
default move constructor
static void setMinCompressedIDCDelta(const float minIDCDeltaValue)
static void setNThreads(const int nThreads)
float getIDCZeroVal(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad) const
void drawIDCDeltaSector(const unsigned int sector, const unsigned int integrationInterval, const float minZ=0, const float maxZ=-1, const IDCDeltaCompression compression=IDCDeltaCompression::NO, const std::string filename="IDCDeltaSector.pdf") const
const auto & getIDCZero() const
unsigned long getNIntegrationIntervalsToChunk(const unsigned int chunk) const
const std::vector< float > & getIDCOneVec(const o2::tpc::Side side) const
void factorizeIDCs(const bool norm, const bool calcDeltas)
void dumpPadFlagMap(const char *outFile, const char *mapName)
void dumpIDCZeroToFile(const Side side, const char *outFileName="IDCZero.root", const char *outName="IDC0") const
dump the IDC0 to file
float getIDCValUngrouped(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const
void checkNeighbourOutliers(const int maxIter=10, const int nOutliersNeighbours=8)
void setIDCZero(const Side side, const IDCZero &idcZero)
set the IDCZero object
auto getIDCDeltaMediumCompressed(const unsigned int chunk, const Side side) const
std::vector< unsigned int > getIntegrationIntervalsPerTF(const int cru=-1) const
unsigned int getNChunks(const Side side) const
auto getIDCDeltaUncompressed(const unsigned int chunk, const Side side) &&
std::vector< unsigned int > getAllIntegrationIntervalsPerTF() const
void drawIDCZeroSector(const unsigned int sector, const float minZ=0, const float maxZ=-1, const std::string filename="IDCZeroSector.pdf") const
float getIDCValGrouped(const unsigned int sector, const unsigned int region, const unsigned int grow, unsigned int gpad, unsigned int integrationInterval) const
unsigned long getNIntegrationIntervals() const
void dumpToFile(const char *outFileName="IDCFactorized.root", const char *outName="IDCFactorized") const
void calcIDCZero(const bool norm)
void setIDCs(std::vector< float > &&idcs, const unsigned int cru, const unsigned int timeframe)
void normIDCZeroGain()
normalize IDC0 with the set gain map
static std::unique_ptr< IDCFactorization > getLargeObjectFromFile(const char *inpFileName="IDCFactorized.root", const char *inName="IDCFactorized")
static void setMaxCompressedIDCDelta(const float maxIDCDeltaValue)
void dumpLargeObjectToFile(const char *outFileName="IDCFactorized.root", const char *outName="IDCFactorized") const
void drawIDCsSector(const unsigned int sector, const unsigned int integrationInterval, const float minZ=0, const float maxZ=-1, const std::string filename="IDCsSector.pdf") const
void setUsePadStatusMap(const bool usePadStatusMap)
setting the usage of the pad-by-pad status map during the factorization of the IDCs
const auto & getIDCDeltaUncompressed(const unsigned int chunk, const Side side) const &
void calcIDCDelta()
calculate \Delta I(r,\phi,t) = I(r,\phi,t) / ( I_0(r,\phi) * I_1(t) )
void checkFECs(const float maxOutliersPerFEC=0.7)
void getLocalIntegrationInterval(const unsigned int integrationInterval, unsigned int &chunk, unsigned int &localintegrationInterval) const
void drawIDCDeltaSide(const o2::tpc::Side side, const unsigned int integrationInterval, const float minZ=0, const float maxZ=-1, const IDCDeltaCompression compression=IDCDeltaCompression::NO, const std::string filename="IDCDeltaSide.pdf") const
void normIDCZeroStackMedian()
normalize IDC0 with the median per stack
const std::vector< float > & getIDCDeltaValuesUncompressed(const unsigned int chunk, const o2::tpc::Side side) const
auto getIDCZero(const o2::tpc::Side side) &&
CalDet< PadFlags > * getPadStatusMapPtr() const
void setIDCOne(const Side side, const IDCOne &idcOne)
set the IDCOne
float getMeanZ(const Side side) const
void drawIDCsSide(const o2::tpc::Side side, const unsigned int integrationInterval, const float minZ=0, const float maxZ=-1, const std::string filename="IDCsSide.pdf") const
const IDCZero & getIDCZero(const o2::tpc::Side side) const &
void dumpToTreeIDC1(const float integrationTimeOrbits=12, const char *outFileName="IDC1Tree.root") const
void dumpIDCOneToFile(const Side side, const char *outFileName="IDCOne.root", const char *outName="IDC1") const
dump the IDC1 to file
IDCFactorization(const unsigned int timeFrames, const unsigned int timeframesDeltaIDC, const std::vector< uint32_t > &crus)
void createStatusMapOutlier(const bool debug=false)
const IDCOne & getIDCOne(const o2::tpc::Side side) const
std::unique_ptr< CalDet< PadFlags > > getPadStatusMap()
void setTimeStamp(const long timeStamp)
void setPadFlagMap(const char *inpFile, const char *mapName)
const std::vector< float > & getIDCZeroVec(const o2::tpc::Side side) const
void drawPadStatusFlagsMapSector(const unsigned int sector, const PadFlags flag=PadFlags::flagSkip, const std::string filename="PadStatusFlags_Sector.pdf") const
auto getIDCDeltaHighCompressed(const unsigned int chunk, const Side side) const
Helper class for accessing grouped pads for one sector.
std::array< std::vector< unsigned int >, Mapper::NREGIONS > mOffsRow
offset to calculate the index in the data from row and pad per region
unsigned int getIndexUngrouped(const unsigned int sector, const unsigned int region, unsigned int ulrow, unsigned int upad, unsigned int integrationInterval) const
static constexpr unsigned int NREGIONS
total number of regions in one sector
Definition Mapper.h:527
Side side() const
Definition Sector.h:96
static constexpr int MAXSECTOR
Definition Sector.h:44
GLenum array
Definition glcorearb.h:4274
GLuint index
Definition glcorearb.h:781
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
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
IDCDeltaCompression
IDC Delta IDC Compression types.
@ NO
no compression using floats
PadFlags
Definition Defs.h:100
@ flagSkip
flag for defining a pad which is just ignored during the calculation of I1 and IDCDelta
Defining DataPointCompositeObject explicitly as copiable.
std::string filename()
struct to access and set Delta IDCs
std::array< std::vector< float >, CRU::MaxCRU > idcs
struct containing the ITPC0 values (integrated TPC clusters)