Project
Loading...
Searching...
No Matches
DataContainer3D.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
16#include "TPCBase/Mapper.h"
17#include "Framework/Logger.h"
18#include "TFile.h"
19#include "ROOT/RDataFrame.hxx"
20
21#include <memory>
22#include <iomanip>
23#include <algorithm>
24
25using namespace o2::tpc;
26
27template <typename DataT>
28template <typename DataTOut>
29int DataContainer3D<DataT>::writeToFile(TFile& outf, const char* name) const
30{
31 if (outf.IsZombie()) {
32 LOGP(error, "Failed to write to file: {}", outf.GetName());
33 return -1;
34 }
35
36 DataContainer3D<DataTOut> containerTmp(mZVertices, mRVertices, mPhiVertices);
37 containerTmp.getData() = std::vector<DataTOut>(mData.begin(), mData.end());
38
39 outf.WriteObjectAny(&containerTmp, DataContainer3D<DataTOut>::Class(), name);
40 return 0;
41}
42
43template <typename DataT>
44int DataContainer3D<DataT>::writeToFile(std::string_view file, std::string_view option, std::string_view name, const int nthreads) const
45{
46 // max number of floats per Entry
47 const size_t maxvalues = sizeof(float) * 1024 * 1024;
48
49 // total number of values to be stored
50 const size_t nsize = getNDataPoints();
51
52 // calculate number of entries in the tree and restrict if the number of values per threads exceeds max size
53 size_t entries = ((nsize / nthreads) > maxvalues) ? (nsize / maxvalues) : nthreads;
54
55 if (entries > nsize) {
56 entries = nsize;
57 }
58
59 // calculate numbers to store per entry
60 const size_t values_per_entry = nsize / entries;
61
62 // in case of remainder add additonal entry
63 const size_t values_lastEntry = nsize % entries;
64 if (values_lastEntry) {
65 entries += 1;
66 }
67
68 // in case EnableImplicitMT was already called with different number of threads, perform reset
69 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
70 ROOT::DisableImplicitMT();
71 }
72 ROOT::EnableImplicitMT(nthreads);
73
74 // define dataframe which will be stored in the TTree
75 ROOT::RDataFrame dFrame(entries);
76
77 // define function which is used to fill the data frame
78 auto dfStore = dFrame.DefineSlotEntry(name, [&data = std::as_const(mData), entries, values_per_entry](unsigned int, ULong64_t entry) { return DataContainer3D<DataT>::getDataSlice(data, entries, values_per_entry, entry); });
79 dfStore = dfStore.Define("nz", [mZVertices = mZVertices]() { return mZVertices; });
80 dfStore = dfStore.Define("nr", [mRVertices = mRVertices]() { return mRVertices; });
81 dfStore = dfStore.Define("nphi", [mPhiVertices = mPhiVertices]() { return mPhiVertices; });
82
83 // define options of TFile
84 ROOT::RDF::RSnapshotOptions opt;
85 opt.fMode = option;
86 opt.fOverwriteIfExists = true; // overwrite if already exists
87
88 TStopwatch timer;
89 // note: first call has some overhead (~2s)
90 dfStore.Snapshot(name, file, {name.data(), "nz", "nr", "nphi"}, opt);
91 timer.Print("u");
92 return 0;
94
95template <typename DataT>
96bool DataContainer3D<DataT>::initFromFile(std::string_view file, std::string_view name, const int nthreads)
98 // in case EnableImplicitMT was already called with different number of threads, perform reset
99 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
100 ROOT::DisableImplicitMT();
101 }
102 ROOT::EnableImplicitMT(nthreads);
103
104 // compare first the meta data (is the number of vertices the same)
105 // define data frame from imput file
106 ROOT::RDataFrame dFrame(name, file);
107
108 // compare vertices
109 auto comp = [mZVertices = mZVertices, mRVertices = mRVertices, mPhiVertices = mPhiVertices](const unsigned short nz, const unsigned short nr, const unsigned short nphi) {
110 if ((nz == mZVertices) && (nr == mRVertices) && (nphi == mPhiVertices)) {
111 return false;
113 return true;
114 };
115
116 auto count = dFrame.Filter(comp, {"nz", "nr", "nphi"}).Count();
117 if (*count != 0) {
118 LOGP(error, "Data from input file has different number of vertices! Found {} same vertices", *count);
119 return false;
121
122 // define lambda function which is used to copy the data
123 auto readData = [&mData = mData](const std::pair<long, std::vector<float>>& data) {
124 std::copy(data.second.begin(), data.second.end(), mData.begin() + data.first);
125 };
126
127 LOGP(info, "Reading {} from file {}", name, file);
128
129 // fill data from RDataFrame
130 TStopwatch timer;
131 dFrame.Foreach(readData, {name.data()});
132 timer.Print("u");
133 return true;
134}
135
137template <typename DataT>
138template <typename DataTIn>
139bool DataContainer3D<DataT>::initFromFile(TFile& inpf, const char* name)
140{
141 if (inpf.IsZombie()) {
142 LOGP(error, "Failed to read from file: {}", inpf.GetName());
143 return false;
144 }
145 DataContainer3D<DataTIn>* dataCont{nullptr};
146 dataCont = reinterpret_cast<DataContainer3D<DataTIn>*>(inpf.GetObjectChecked(name, DataContainer3D<DataTIn>::Class()));
147
148 if (!dataCont) {
149 LOGP(error, "Failed to load {} from {}", name, inpf.GetName());
150 return false;
151 }
153 if (mZVertices != dataCont->getNZ() || mRVertices != dataCont->getNR() || mPhiVertices != dataCont->getNPhi()) {
154 LOGP(error, "Data from input file has different definition of vertices!");
155 LOGP(error, "set vertices before creating the sc object to: SpaceCharge<>::setGrid({}, {}, {})", dataCont->getNZ(), dataCont->getNR(), dataCont->getNPhi());
156 delete dataCont;
157 return false;
158 }
159
160 mData = std::vector<DataT>(dataCont->getData().begin(), dataCont->getData().end());
161 delete dataCont;
162 return true;
163}
164
165template <typename DataT>
167{
168 if (inpf.IsZombie()) {
169 LOGP(error, "Failed to read from file {}", inpf.GetName());
170 return nullptr;
171 }
172 DataContainer3D<DataT>* dataCont{nullptr};
173
174 dataCont = reinterpret_cast<DataContainer3D<DataT>*>(inpf.GetObjectChecked(name, DataContainer3D<DataT>::Class()));
175 if (!dataCont) {
176 LOGP(error, "Failed to load {} from {}", name, inpf.GetName());
177 return nullptr;
178 }
179 return dataCont;
180}
181
182template <typename DataT>
185 std::stringstream stream;
186 stream.precision(3);
187 auto&& w = std::setw(9);
188 stream << std::endl;
189
190 for (unsigned int iz = 0; iz < mPhiVertices; ++iz) {
191 stream << "z layer: " << iz << "\n";
192 // print top x row
193 stream << "⎡" << w << (*this)(0, 0, iz);
194 for (unsigned int ix = 1; ix < mZVertices; ++ix) {
195 stream << ", " << w << (*this)(ix, 0, iz);
196 }
197 stream << " ⎤ \n";
198
199 for (unsigned int iy = 1; iy < mRVertices - 1; ++iy) {
200 stream << "⎢" << w << (*this)(0, iy, iz);
201 for (unsigned int ix = 1; ix < mZVertices; ++ix) {
202 stream << ", " << w << (*this)(ix, iy, iz);
203 }
204 stream << " ⎥ \n";
205 }
206
207 stream << "⎣" << w << (*this)(0, mRVertices - 1, iz);
208 for (unsigned int ix = 1; ix < mZVertices; ++ix) {
209 stream << ", " << w << (*this)(ix, mRVertices - 1, iz);
210 }
211 stream << " ⎦ \n \n";
212 }
213 LOGP(info, "{} \n \n", stream.str());
214}
215
216template <typename DataT>
217auto DataContainer3D<DataT>::getDataSlice(const std::vector<DataT>& data, size_t entries, const size_t values_per_entry, ULong64_t entry)
218{
219 const long indStart = entry * values_per_entry;
220 if (entry < (entries - 1)) {
221 return std::pair(indStart, std::vector<float>(data.begin() + indStart, data.begin() + indStart + values_per_entry));
222 } else if (entry == (entries - 1)) {
223 // last entry might have different number of values. just copy the rest...
224 return std::pair(indStart, std::vector<float>(data.begin() + indStart, data.end()));
225 }
226 return std::pair(indStart, std::vector<float>());
227};
228
229template <typename DataT>
231{
232 std::transform(mData.begin(), mData.end(), mData.begin(), [value = value](auto& val) { return val * value; });
233 return *this;
234}
235
236template <typename DataT>
238{
239 std::transform(mData.begin(), mData.end(), other.mData.begin(), mData.begin(), std::plus<>());
240 return *this;
241}
242
243template <typename DataT>
245{
246 std::transform(mData.begin(), mData.end(), other.mData.begin(), mData.begin(), std::minus<>());
247 return *this;
248}
249
250template <typename DataT>
252{
253 std::transform(mData.begin(), mData.end(), other.mData.begin(), mData.begin(), std::multiplies<>());
254 return *this;
255}
256
257template <typename DataT>
258size_t DataContainer3D<DataT>::getIndexZ(size_t index, const int nz, const int nr, const int nphi)
259{
260 const size_t iphi = index / (nz * nr);
261 index -= (iphi * nz * nr);
262 const size_t iz = index % nz;
263 return iz;
264}
265
266template <typename DataT>
267size_t DataContainer3D<DataT>::getIndexR(size_t index, const int nz, const int nr, const int nphi)
268{
269 const size_t iphi = index / (nz * nr);
270 index -= (iphi * nz * nr);
271 const size_t ir = index / nz;
272 return ir;
273}
274
275template <typename DataT>
276size_t DataContainer3D<DataT>::getIndexPhi(size_t index, const int nz, const int nr, const int nphi)
277{
278 return index / (nz * nr);
279}
280
281template <typename DataT>
283{
284 tree->SetAlias("ir", "o2::tpc::DataContainer3D<float>::getIndexR(first + Iteration$, nz, nr, nphi)");
285 tree->SetAlias("iz", "o2::tpc::DataContainer3D<float>::getIndexZ(first + Iteration$, nz, nr, nphi)");
286 tree->SetAlias("iphi", "o2::tpc::DataContainer3D<float>::getIndexPhi(first + Iteration$, nz, nr, nphi)");
287 tree->SetAlias("r", "o2::tpc::GridProperties<float>::getRMin() + o2::tpc::GridProperties<float>::getGridSpacingR(nr) * ir");
288 tree->SetAlias("z", "o2::tpc::GridProperties<float>::getZMin() + o2::tpc::GridProperties<float>::getGridSpacingZ(nz) * iz");
289 tree->SetAlias("phi", "o2::tpc::GridProperties<float>::getPhiMin() + o2::tpc::GridProperties<float>::getGridSpacingPhi(nphi) * iphi");
290}
291
292template <typename DataT>
294{
295 // actuall stored value
296 tree->SetAlias("val", "_0");
297
298 // some meta data
299 tree->SetAlias("iz", "_1");
300 tree->SetAlias("ir", "_2");
301 tree->SetAlias("iphi", "_3");
302 tree->SetAlias("z", "_4");
303 tree->SetAlias("r", "_5");
304 tree->SetAlias("phi", "_6");
305 tree->SetAlias("lpos", "_7");
306 tree->SetAlias("lx", "lpos.fCoordinates.fX");
307 tree->SetAlias("ly", "lpos.fCoordinates.fY");
308 tree->SetAlias("index", "_8");
309}
310
311template <typename DataT>
312void DataContainer3D<DataT>::setGrid(unsigned short nZ, unsigned short nR, unsigned short nPhi, const bool resize)
313{
314 mZVertices = nZ;
315 mRVertices = nR;
316 mPhiVertices = nPhi;
317 if (resize) {
318 mData.resize(nZ * nR * static_cast<size_t>(nPhi));
319 }
320}
321
322template <typename DataT>
323void DataContainer3D<DataT>::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)
324{
325 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
326 ROOT::DisableImplicitMT();
327 }
328 ROOT::EnableImplicitMT(nthreads);
329 ROOT::RDataFrame dFrame(treename, fileIn);
330
331 auto df = dFrame.Define("slice", [rangeiZ, rangeiR, rangeiPhi](const std::pair<long, std::vector<float>>& values, unsigned short nz, unsigned short nr, unsigned short nphi) {
332 const bool simOneSectorOnly = MGParameters::normalizeGridToOneSector;
333 std::vector<size_t> ir;
334 std::vector<size_t> iphi;
335 std::vector<size_t> iz;
336 std::vector<float> r;
337 std::vector<float> phi;
338 std::vector<float> z;
339 std::vector<float> vals;
340 std::vector<size_t> globalIdx;
341 std::vector<LocalPosition3D> lPos;
342 const auto nvalues = values.second.size();
343 ir.reserve(nvalues);
344 iphi.reserve(nvalues);
345 iz.reserve(nvalues);
346 r.reserve(nvalues);
347 phi.reserve(nvalues);
348 z.reserve(nvalues);
349 vals.reserve(nvalues);
350 lPos.reserve(nvalues);
351 globalIdx.reserve(nvalues);
352 for (size_t i = 0; i < nvalues; ++i) {
353 const size_t idx = values.first + i;
354 const auto iZTmp = o2::tpc::DataContainer3D<float>::getIndexZ(idx, nz, nr, nphi);
355 if ((rangeiZ.first < rangeiZ.second) && ((iZTmp < rangeiZ.first) || (iZTmp > rangeiZ.second))) {
356 continue;
357 }
358
359 const auto iRTmp = o2::tpc::DataContainer3D<float>::getIndexR(idx, nz, nr, nphi);
360 if ((rangeiR.first < rangeiR.second) && ((iRTmp < rangeiR.first) || (iRTmp > rangeiR.second))) {
361 continue;
362 }
363
364 const auto iPhiTmp = o2::tpc::DataContainer3D<float>::getIndexPhi(idx, nz, nr, nphi);
365 if ((rangeiPhi.first < rangeiPhi.second) && ((iPhiTmp < rangeiPhi.first) || (iPhiTmp > rangeiPhi.second))) {
366 continue;
367 }
368
371 const float phiTmp = o2::tpc::GridProperties<float>::getPhiMin() + o2::tpc::GridProperties<float>::getGridSpacingPhi(nphi) / (simOneSectorOnly ? SECTORSPERSIDE : 1) * iPhiTmp;
372
373 const float x = rTmp * std::cos(phiTmp);
374 const float y = rTmp * std::sin(phiTmp);
375 const LocalPosition3D pos(x, y, zTmp);
376 unsigned char secNum = simOneSectorOnly ? 0 : std::floor(phiTmp / SECPHIWIDTH);
377 Sector sector(secNum + (pos.Z() < 0) * SECTORSPERSIDE);
378 LocalPosition3D lPosTmp = Mapper::GlobalToLocal(pos, sector);
379
380 lPos.emplace_back(lPosTmp);
381 ir.emplace_back(iRTmp);
382 iphi.emplace_back(iPhiTmp);
383 iz.emplace_back(iZTmp);
384 r.emplace_back(rTmp);
385 phi.emplace_back(phiTmp);
386 z.emplace_back(zTmp);
387 vals.emplace_back(values.second[i]);
388 globalIdx.emplace_back(idx);
389 }
390 return std::make_tuple(vals, iz, ir, iphi, z, r, phi, lPos, globalIdx);
391 },
392 {treename.data(), "nz", "nr", "nphi"});
393
394 // define options of TFile
395 ROOT::RDF::RSnapshotOptions opt;
396 opt.fMode = option;
397 df.Snapshot(treename, fileOut, {"slice"}, opt);
398}
399
400template <typename DataT>
402{
403 TriCubicInterpolator<DataT> interpolator(*this, grid);
404 return interpolator(z, r, phi);
405}
406
407template <typename DataT>
408void DataContainer3D<DataT>::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)
409{
410 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
411 ROOT::DisableImplicitMT();
412 }
413 ROOT::EnableImplicitMT(nthreads);
414 ROOT::RDataFrame dFrame(nPhi);
415
416 // get vertices for input TTree which is needed to define the grid for interpolation
417 unsigned short nr, nz, nphi;
418 if (!getVertices(treename, fileIn, nr, nz, nphi)) {
419 return;
420 }
421
422 // load data from input TTree
424 data.setGrid(nz, nr, nphi, true);
425 data.initFromFile(fileIn, treename, nthreads);
426
427 // define grid for interpolation
428 using GridProp = GridProperties<DataT>;
429 const RegularGrid3D<DataT> mGrid3D(GridProp::ZMIN, GridProp::RMIN, GridProp::PHIMIN, GridProp::getGridSpacingZ(nz), GridProp::getGridSpacingR(nr), o2::tpc::GridProperties<float>::getGridSpacingPhi(nphi) / (MGParameters::normalizeGridToOneSector ? SECTORSPERSIDE : 1), ParamSpaceCharge{nr, nz, nphi});
430
431 auto interpolate = [&mGrid3D = std::as_const(mGrid3D), &data = std::as_const(data), rangeR, rangeZ, rangePhi, nR, nZ, nPhi](unsigned int, ULong64_t iPhi) {
432 const bool simOneSectorOnly = MGParameters::normalizeGridToOneSector;
433 std::vector<size_t> ir;
434 std::vector<size_t> iphi;
435 std::vector<size_t> iz;
436 std::vector<float> r;
437 std::vector<float> phi;
438 std::vector<float> z;
439 std::vector<float> vals;
440 std::vector<size_t> globalIdx;
441 std::vector<LocalPosition3D> lPos;
442 const auto nvalues = nR * nZ;
443 ir.reserve(nvalues);
444 iphi.reserve(nvalues);
445 iz.reserve(nvalues);
446 r.reserve(nvalues);
447 phi.reserve(nvalues);
448 z.reserve(nvalues);
449 vals.reserve(nvalues);
450 lPos.reserve(nvalues);
451 globalIdx.reserve(nvalues);
452
453 const float rSpacing = (rangeR.second - rangeR.first) / (nR - 1);
454 const float zSpacing = (rangeZ.second - rangeZ.first) / (nZ - 1);
455 const float phiSpacing = (rangePhi.second - rangePhi.first) / (nPhi - 1);
456 const DataT phiPos = rangePhi.first + iPhi * phiSpacing;
457 // loop over grid and interpolate values
458 for (int iR = 0; iR < nR; ++iR) {
459 const DataT rPos = rangeR.first + iR * rSpacing;
460 for (int iZ = 0; iZ < nZ; ++iZ) {
461 const size_t idx = (iZ + nZ * (iR + iPhi * nR)); // unique index to Build index with other friend TTrees
462 const DataT zPos = rangeZ.first + iZ * zSpacing;
463 ir.emplace_back(iR);
464 iphi.emplace_back(iPhi);
465 iz.emplace_back(iZ);
466 r.emplace_back(rPos);
467 phi.emplace_back(phiPos);
468 z.emplace_back(zPos);
469 vals.emplace_back(data.interpolate(zPos, rPos, phiPos, mGrid3D)); // interpolated values
470 globalIdx.emplace_back(idx);
471 const float x = rPos * std::cos(phiPos);
472 const float y = rPos * std::sin(phiPos);
473 const LocalPosition3D pos(x, y, zPos);
474 unsigned char secNum = simOneSectorOnly ? 0 : std::floor(phiPos / SECPHIWIDTH);
475 Sector sector(secNum + (pos.Z() < 0) * SECTORSPERSIDE);
476 LocalPosition3D lPosTmp = Mapper::GlobalToLocal(pos, sector);
477 lPos.emplace_back(lPosTmp);
478 }
479 }
480 return std::make_tuple(vals, iz, ir, iphi, z, r, phi, lPos, globalIdx);
481 };
482
483 // define RDataFrame entry
484 auto dfStore = dFrame.DefineSlotEntry(treename, interpolate);
485
486 // define options of TFile
487 ROOT::RDF::RSnapshotOptions opt;
488 opt.fMode = option;
489
490 TStopwatch timer;
491 // note: first call has some overhead (~2s)
492 dfStore.Snapshot(treename, fileOut, {treename.data()}, opt);
493 timer.Print("u");
494}
495
496template <typename DataT>
497bool DataContainer3D<DataT>::getVertices(std::string_view treename, std::string_view fileIn, unsigned short& nR, unsigned short& nZ, unsigned short& nPhi)
498{
499 TFile fTmp(fileIn.data(), "READ");
500 TTree* tree = (TTree*)fTmp.Get(treename.data());
501 if (!tree) {
502 LOGP(warning, "Tree {} not found in input file {}", treename, fileIn);
503 return false;
504 }
505 tree->SetBranchAddress("nz", &nZ);
506 tree->SetBranchAddress("nr", &nR);
507 tree->SetBranchAddress("nphi", &nPhi);
508 tree->GetEntry(0);
509 delete tree;
510 return true;
511}
512
515
516// deprecated functions (to be removed...)
517template int o2::tpc::DataContainer3D<float>::writeToFile<float>(TFile&, const char*) const;
518template int o2::tpc::DataContainer3D<float>::writeToFile<double>(TFile&, const char*) const;
519template int o2::tpc::DataContainer3D<double>::writeToFile<float>(TFile&, const char*) const;
520template int o2::tpc::DataContainer3D<double>::writeToFile<double>(TFile&, const char*) const;
521template bool o2::tpc::DataContainer3D<float>::initFromFile<float>(TFile&, const char*);
522template bool o2::tpc::DataContainer3D<float>::initFromFile<double>(TFile&, const char*);
523template bool o2::tpc::DataContainer3D<double>::initFromFile<float>(TFile&, const char*);
524template bool o2::tpc::DataContainer3D<double>::initFromFile<double>(TFile&, const char*);
This class provides a simple method to store values on a large 3-Dim grid with ROOT io functionality.
int32_t i
uint16_t pos
Definition RawData.h:3
Definition of RegularGrid3D class.
Definition of TriCubic class.
static LocalPosition3D GlobalToLocal(const GlobalPosition3D &pos, const double alpha)
Definition Mapper.h:468
GLint GLenum GLint x
Definition glcorearb.h:403
GLint GLsizei count
Definition glcorearb.h:399
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
GLenum GLsizei GLsizei GLint * values
Definition glcorearb.h:1576
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat * val
Definition glcorearb.h:1582
GLboolean r
Definition glcorearb.h:1233
GLuint GLuint stream
Definition glcorearb.h:1806
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
Global TPC definitions and constants.
Definition SimTraits.h:167
constexpr double SECPHIWIDTH
Definition Defs.h:45
constexpr unsigned char SECTORSPERSIDE
Definition Defs.h:40
void readData(o2::tpc::GBTFrameContainer &container, std::vector< std::ofstream * > &outfiles, int &run, int &done)
const auto & getData() const
static size_t getIndexR(size_t index, const int nz, const int nr, const int nphi)
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)
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)
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)
static bool getVertices(std::string_view treename, std::string_view fileIn, unsigned short &nR, unsigned short &nZ, unsigned short &nPhi)
static DataContainer3D< DataT > * loadFromFile(TFile &inpf, const char *name="data")
get pointer to object from file (deprecated!)
DataContainer3D< DataT > & operator-=(const DataContainer3D< DataT > &other)
DataContainer3D< DataT > & operator+=(const DataContainer3D< DataT > &other)
static constexpr DataT getRMin()
static constexpr DataT getZMin()
static constexpr DataT getGridSpacingR(const unsigned int nR)
static constexpr DataT getPhiMin()
static constexpr DataT getGridSpacingZ(const unsigned int nZ)
static constexpr DataT getGridSpacingPhi(const unsigned int nPhi)
static bool normalizeGridToOneSector
the grid in phi direction is squashed from 2 Pi to (2 Pi / SECTORSPERSIDE). This can used to get the ...
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)