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"
28#include <boost/property_tree/ptree.hpp>
29
30namespace o2::tpc
31{
32
33template <class T>
34class CalDet;
35
37 std::array<std::vector<float>, CRU::MaxCRU> idcs{};
38 int relTF{};
39 ClassDefNV(IDCFactorizeSplit, 1)
40};
41
43{
44 public:
52 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);
53
57 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){};
58
60 IDCFactorization() = default;
61
64
67
70 static std::vector<o2::tpc::Side> getSides(const std::vector<uint32_t>& crus);
71
77 void factorizeIDCs(const bool norm, const bool calcDeltas);
78
81 void calcIDCZero(const bool norm);
82
85
88 void createStatusMapOutlier(const bool debug = false);
89
91 void createStatusMap();
92
95 void checkFECs(const float maxOutliersPerFEC = 0.7);
96
101 void checkNeighbourOutliers(const int maxIter = 10, const int nOutliersNeighbours = 8);
102
104 void calcIDCOne();
105
107 void calcIDCDelta();
108
110 template <typename DataVec>
111 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);
112
119 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]; }
120
127 float getIDCValUngrouped(const unsigned int sector, const unsigned int region, unsigned int urow, unsigned int upad, unsigned int integrationInterval) const;
128
134 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)); }
135
143 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)); }
144
149 void getLocalIntegrationInterval(const unsigned int integrationInterval, unsigned int& chunk, unsigned int& localintegrationInterval) const;
150
152 unsigned int getNTimeframes() const { return mTimeFrames; }
153
156 unsigned long getNIntegrationIntervalsInChunk(const unsigned int chunk) const;
157
160 unsigned long getNIntegrationIntervalsToChunk(const unsigned int chunk) const;
161
163 unsigned long getNIntegrationIntervals(const int cru) const;
164
166 unsigned long getNIntegrationIntervals() const;
167
170 const std::vector<float>& getIDCZeroVec(const o2::tpc::Side side) const { return mIDCZero[mSideIndex[side]].mIDCZero; }
171
174 const IDCZero& getIDCZero(const o2::tpc::Side side) const& { return mIDCZero[mSideIndex[side]]; }
175
178 auto getIDCZero(const o2::tpc::Side side) && { return mIDCZero[mSideIndex[side]]; }
179
182 const std::vector<float>& getIDCOneVec(const o2::tpc::Side side) const { return mIDCOne[mSideIndex[side]].mIDCOne; }
183
186 const IDCOne& getIDCOne(const o2::tpc::Side side) const { return mIDCOne[mSideIndex[side]]; }
187
191 const std::vector<float>& getIDCDeltaValuesUncompressed(const unsigned int chunk, const o2::tpc::Side side) const { return mIDCDelta[mSideIndex[side]][chunk].getIDCDelta(); }
192
195 const auto& getIDCDeltaUncompressed(const unsigned int chunk, const Side side) const& { return mIDCDelta[mSideIndex[side]][chunk]; }
196
198 auto getIDCDeltaUncompressed(const unsigned int chunk, const Side side) && { return std::move(mIDCDelta[mSideIndex[side]][chunk]); }
199
202 auto getIDCDeltaMediumCompressed(const unsigned int chunk, const Side side) const { return IDCDeltaCompressionHelper<unsigned short>::getCompressedIDCs(mIDCDelta[mSideIndex[side]][chunk]); }
203
206 auto getIDCDeltaHighCompressed(const unsigned int chunk, const Side side) const { return IDCDeltaCompressionHelper<unsigned char>::getCompressedIDCs(mIDCDelta[mSideIndex[side]][chunk]); }
207
209 unsigned int getNChunks(const Side side) const { return mIDCDelta[mSideIndex[side]].size(); }
210
212 const auto& getIDCZero() const { return mIDCZero; }
213
215 const auto& getIDCOne() const& { return mIDCOne; }
216
218 auto getIDCOne() && { return std::move(mIDCOne); }
219
221 const auto& getIDCs() const { return mIDCs; }
222
224 auto getCRUs() const { return mCRUs; }
225
227 void setTimeStamp(const long timeStamp) { mTimeStamp = timeStamp; }
228
230 long getTimeStamp() const { return mTimeStamp; }
231
233 void setRun(const int run) { mRun = run; }
234
236 int getRun() const { return mRun; }
237
238 // get number of TFs in which the DeltaIDCs are split/stored
239 unsigned int getTimeFramesDeltaIDC() const { return mTimeFramesDeltaIDC; }
240
242 static int getNThreads() { return sNThreads; }
243
248 void setIDCs(std::vector<float>&& idcs, const unsigned int cru, const unsigned int timeframe) { mIDCs[cru][timeframe] = std::move(idcs); }
249
252 static void setNThreads(const int nThreads) { sNThreads = nThreads; }
253
255 static void setMinCompressedIDCDelta(const float minIDCDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCIDCCompressionParam", "minIDCDeltaValue", minIDCDeltaValue); }
256
258 static void setMaxCompressedIDCDelta(const float maxIDCDeltaValue) { o2::conf::ConfigurableParam::setValue<float>("TPCIDCCompressionParam", "maxIDCDeltaValue", maxIDCDeltaValue); }
259
264 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); }
265
269 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); }
270
276 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); }
277
282 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); }
283
286 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); }
287
291 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); }
292
298 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); }
299
304 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); }
305
310 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); }
311
315 void dumpToFile(const char* outFileName = "IDCFactorized.root", const char* outName = "IDCFactorized") const;
316
320 void dumpLargeObjectToFile(const char* outFileName = "IDCFactorized.root", const char* outName = "IDCFactorized") const;
321
326 static std::unique_ptr<IDCFactorization> getLargeObjectFromFile(const char* inpFileName = "IDCFactorized.root", const char* inName = "IDCFactorized");
327
329 void dumpIDCZeroToFile(const Side side, const char* outFileName = "IDCZero.root", const char* outName = "IDC0") const;
330
332 void dumpIDCOneToFile(const Side side, const char* outFileName = "IDCOne.root", const char* outName = "IDC1") const;
333
336 void dumpToTree(const Side side, const char* outFileName = "IDCTree.root");
337
341 void dumpToTreeIDC1(const float integrationTimeOrbits = 12, const char* outFileName = "IDC1Tree.root") const;
342
344 void dumpIDCDeltaToTree(const Side side, const int chunk = 0, const char* outFileName = "IDCDeltaTree.root");
345
348 std::vector<unsigned int> getIntegrationIntervalsPerTF(const int cru = -1) const;
349
351 std::vector<unsigned int> getAllIntegrationIntervalsPerTF() const;
352
354 void reset();
355
359 void setGainMap(const char* inpFile, const char* mapName);
360
362 void setGainMap(const CalDet<float>& gainmap);
363
367 void setPadFlagMap(const char* inpFile, const char* mapName);
368
370 void setPadFlagMap(const CalDet<PadFlags>& flagmap);
371
373 void setUsePadStatusMap(const bool usePadStatusMap) { mUsePadStatusMap = usePadStatusMap; }
374
376 bool getUsePadStatusMap() const { return mUsePadStatusMap; }
377
381 void dumpPadFlagMap(const char* outFile, const char* mapName);
382
384 CalDet<PadFlags>* getPadStatusMapPtr() const { return mPadFlagsMap.get(); }
385
387 std::unique_ptr<CalDet<PadFlags>> getPadStatusMap() { return std::move(mPadFlagsMap); }
388
390 const std::vector<Side>& getSides() const { return mSides; }
391
393 void setIDCZero(const Side side, const IDCZero& idcZero) { mIDCZero[mSideIndex[side]] = idcZero; }
394
396 void setIDCOne(const Side side, const IDCOne& idcOne) { mIDCOne[mSideIndex[side]] = idcOne; }
397
399 void setIDCDelta(const Side side, const IDCDelta<float>& idcDelta, const int index = 0) { mIDCDelta[mSideIndex[side]][index] = idcDelta; }
400
402 float getMeanZ(const Side side) const;
403
405 void normIDCZeroGain() { normIDCZero(0); }
406
408 void normIDCZeroStackMedian() { normIDCZero(1); }
409
411 std::array<float, o2::tpc::GEMSTACKS> getStackMedian() const;
412
413 private:
414 const unsigned int mTimeFrames{};
415 const unsigned int mTimeFramesDeltaIDC{};
416 std::array<std::vector<std::vector<float>>, CRU::MaxCRU> mIDCs{};
417 std::vector<IDCZero> mIDCZero{};
418 std::vector<IDCOne> mIDCOne{};
419 std::vector<std::vector<IDCDelta<float>>> mIDCDelta{};
420 inline static int sNThreads{1};
421 std::unique_ptr<CalDet<float>> mGainMap;
422 std::unique_ptr<CalDet<PadFlags>> mPadFlagsMap;
423 bool mInputGrouped{false};
424 bool mUsePadStatusMap{true};
425 const std::vector<uint32_t> mCRUs{};
426 std::array<unsigned int, SIDES> mSideIndex{0, 1};
427 std::vector<Side> mSides{};
428 std::vector<unsigned int> mIntegrationIntervalsPerTF{};
429 long mTimeStamp{0};
430 int mRun{0};
431
433 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;
434
436 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;
437
439 void drawIDCZeroHelper(const bool type, const Sector sector, const std::string filename, const float minZ, const float maxZ) const;
440
442 void getTF(unsigned int integrationInterval, unsigned int& timeFrame, unsigned int& interval) const;
443
445 unsigned int getChunk(const unsigned int timeframe) const;
446
448 unsigned int getNTFsPerChunk(const unsigned int chunk) const;
449
451 void drawPadFlagMap(const bool type, const Sector sector, const std::string filename, const PadFlags flag) const;
452
454 void normIDCZero(const int type);
455
458 bool checkReceivedIDCs();
459
460 ClassDefNV(IDCFactorization, 2)
461};
462
463} // namespace o2::tpc
464
465#endif
std::ostringstream debug
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
@ 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:168
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
IDCDeltaCompression
IDC Delta IDC Compression types.
@ NO
no compression using floats
@ 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)