Project
Loading...
Searching...
No Matches
DataContainer3D.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_TPC_DATACONTAINER3D_H_
19#define ALICEO2_TPC_DATACONTAINER3D_H_
20
21#include "Rtypes.h"
22#include <vector>
23
24class TFile;
25class TTree;
26
27namespace o2
28{
29namespace tpc
30{
31
32template <class T>
33class RegularGrid3D;
34
37
39template <typename DataT = double>
41
46 DataContainer3D(unsigned short nZ, unsigned short nR, unsigned short nPhi) : mZVertices{nZ}, mRVertices{nR}, mPhiVertices{nPhi}, mData(nZ * nR * nPhi){};
47
49 DataContainer3D() = default;
50
52 const DataT& operator[](size_t i) const { return mData[i]; }
53 DataT& operator[](size_t i) { return mData[i]; }
54
55 const auto& getData() const { return mData; }
56 auto& getData() { return mData; }
57
62 const DataT& operator()(size_t iz, size_t ir, size_t iphi) const { return mData[getDataIndex(iz, ir, iphi)]; }
63
68 DataT& operator()(size_t iz, size_t ir, size_t iphi) { return mData[getDataIndex(iz, ir, iphi)]; }
69
71 DataT interpolate(const DataT z, const DataT r, const DataT phi, const o2::tpc::RegularGrid3D<DataT>& grid) const;
72
77 size_t getDataIndex(const size_t iz, const size_t ir, const size_t iphi) const { return (iz + mZVertices * (ir + iphi * mRVertices)); }
78
80 size_t getNDataPoints() const { return mData.size(); }
81
83 unsigned short getNZ() const { return mZVertices; }
84
86 unsigned short getNR() const { return mRVertices; }
87
89 unsigned short getNPhi() const { return mPhiVertices; }
90
93 static void setAliases(TTree* tree);
94
97 static void setAliasesForDump(TTree* tree);
98
104 static size_t getIndexZ(size_t index, const int nz, const int nr, const int nphi);
105 size_t getIndexZ(size_t index) const { return getIndexZ(index, getNZ(), getNR(), getNPhi()); }
106
112 static size_t getIndexR(size_t index, const int nz, const int nr, const int nphi);
113 size_t getIndexR(size_t index) const { return getIndexR(index, getNZ(), getNR(), getNPhi()); }
114
120 static size_t getIndexPhi(size_t index, const int nz, const int nr, const int nphi);
121 size_t getIndexPhi(size_t index) const { return getIndexPhi(index, getNZ(), getNR(), getNPhi()); }
122
124 void setGrid(unsigned short nZ, unsigned short nR, unsigned short nPhi, const bool resize);
125
130 template <typename DataTOut = DataT>
131 int writeToFile(TFile& outf, const char* name = "data") const;
132
138 int writeToFile(std::string_view file, std::string_view option, std::string_view name = "data", const int nthreads = 1) const;
139
142 template <typename DataTIn = DataT>
143 bool initFromFile(TFile& inpf, const char* name = "data");
144
149 bool initFromFile(std::string_view file, std::string_view name = "data", const int nthreads = 1);
150
152 inline static DataContainer3D<DataT>* loadFromFile(TFile& inpf, const char* name = "data");
153
162 static void dumpSlice(std::string_view treename, std::string_view fileIn, std::string_view fileOut, std::string_view option, std::pair<unsigned short, unsigned short> rangeiR, std::pair<unsigned short, unsigned short> rangeiZ, std::pair<unsigned short, unsigned short> rangeiPhi, const int nthreads = 1);
163
175 static void dumpInterpolation(std::string_view treename, std::string_view fileIn, std::string_view fileOut, std::string_view option, std::pair<float, float> rangeR, std::pair<float, float> rangeZ, std::pair<float, float> rangePhi, const int nR, const int nZ, const int nPhi, const int nthreads = 1);
176
184 static bool getVertices(std::string_view treename, std::string_view fileIn, unsigned short& nR, unsigned short& nZ, unsigned short& nPhi);
185
187 void print() const;
188
194
195 private:
196 unsigned short mZVertices{};
197 unsigned short mRVertices{};
198 unsigned short mPhiVertices{};
199 std::vector<DataT> mData{};
200
201 static auto getDataSlice(const std::vector<DataT>& data, size_t entries, const size_t values_per_entry, ULong64_t entry);
202
203 ClassDefNV(DataContainer3D, 1)
204};
205
206} // namespace tpc
207} // namespace o2
208
209#endif
int32_t i
GLuint entry
Definition glcorearb.h:5735
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLboolean * data
Definition glcorearb.h:298
GLboolean r
Definition glcorearb.h:1233
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
const auto & getData() const
static size_t getIndexR(size_t index, const int nz, const int nr, const int nphi)
unsigned short getNZ() const
unsigned short getNR() const
static void dumpSlice(std::string_view treename, std::string_view fileIn, std::string_view fileOut, std::string_view option, std::pair< unsigned short, unsigned short > rangeiR, std::pair< unsigned short, unsigned short > rangeiZ, std::pair< unsigned short, unsigned short > rangeiPhi, const int nthreads=1)
void print() const
print the matrix
static void dumpInterpolation(std::string_view treename, std::string_view fileIn, std::string_view fileOut, std::string_view option, std::pair< float, float > rangeR, std::pair< float, float > rangeZ, std::pair< float, float > rangePhi, const int nR, const int nZ, const int nPhi, const int nthreads=1)
unsigned short getNPhi() const
bool initFromFile(TFile &inpf, const char *name="data")
set values from file
int writeToFile(TFile &outf, const char *name="data") const
void setGrid(unsigned short nZ, unsigned short nR, unsigned short nPhi, const bool resize)
set the grid points
static void setAliasesForDump(TTree *tree)
static void setAliases(TTree *tree)
size_t getNDataPoints() const
DataT interpolate(const DataT z, const DataT r, const DataT phi, const o2::tpc::RegularGrid3D< DataT > &grid) const
DataContainer3D< DataT > & operator*=(const DataT value)
operator overload
static size_t getIndexZ(size_t index, const int nz, const int nr, const int nphi)
static size_t getIndexPhi(size_t index, const int nz, const int nr, const int nphi)
size_t getIndexPhi(size_t index) const
static bool getVertices(std::string_view treename, std::string_view fileIn, unsigned short &nR, unsigned short &nZ, unsigned short &nPhi)
DataContainer3D()=default
default constructor for Root I/O
static DataContainer3D< DataT > * loadFromFile(TFile &inpf, const char *name="data")
get pointer to object from file (deprecated!)
DataContainer3D< DataT > & operator-=(const DataContainer3D< DataT > &other)
size_t getIndexZ(size_t index) const
DataContainer3D(unsigned short nZ, unsigned short nR, unsigned short nPhi)
DataContainer3D< DataT > & operator+=(const DataContainer3D< DataT > &other)
DataT & operator()(size_t iz, size_t ir, size_t iphi)
size_t getDataIndex(const size_t iz, const size_t ir, const size_t iphi) const
size_t getIndexR(size_t index) const
const DataT & operator()(size_t iz, size_t ir, size_t iphi) const
DataT & operator[](size_t i)
const DataT & operator[](size_t i) const
operator to directly access the values
VectorOfTObjectPtrs other
o2::InteractionRecord ir(0, 0)
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))
pedsdata resize(norbits)