19#include "ROOT/RDataFrame.hxx"
27template <
typename DataT>
28template <
typename DataTOut>
31 if (outf.IsZombie()) {
32 LOGP(error,
"Failed to write to file: {}", outf.GetName());
37 containerTmp.
getData() = std::vector<DataTOut>(mData.begin(), mData.end());
43template <
typename DataT>
47 const size_t maxvalues =
sizeof(float) * 1024 * 1024;
50 const size_t nsize = getNDataPoints();
53 size_t entries = ((nsize / nthreads) > maxvalues) ? (nsize / maxvalues) : nthreads;
55 if (entries > nsize) {
60 const size_t values_per_entry = nsize / entries;
63 const size_t values_lastEntry = nsize % entries;
64 if (values_lastEntry) {
69 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
70 ROOT::DisableImplicitMT();
72 ROOT::EnableImplicitMT(nthreads);
75 ROOT::RDataFrame dFrame(entries);
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; });
84 ROOT::RDF::RSnapshotOptions opt;
86 opt.fOverwriteIfExists =
true;
90 dfStore.Snapshot(
name,
file, {
name.data(),
"nz",
"nr",
"nphi"}, opt);
95template <
typename DataT>
99 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
100 ROOT::DisableImplicitMT();
102 ROOT::EnableImplicitMT(nthreads);
106 ROOT::RDataFrame dFrame(
name,
file);
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)) {
116 auto count = dFrame.Filter(comp, {
"nz",
"nr",
"nphi"}).Count();
118 LOGP(error,
"Data from input file has different number of vertices! Found {} same vertices", *
count);
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);
127 LOGP(info,
"Reading {} from file {}",
name,
file);
137template <
typename DataT>
138template <
typename DataTIn>
141 if (inpf.IsZombie()) {
142 LOGP(error,
"Failed to read from file: {}", inpf.GetName());
149 LOGP(error,
"Failed to load {} from {}",
name, inpf.GetName());
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());
160 mData = std::vector<DataT>(dataCont->getData().begin(), dataCont->getData().end());
165template <
typename DataT>
168 if (inpf.IsZombie()) {
169 LOGP(error,
"Failed to read from file {}", inpf.GetName());
176 LOGP(error,
"Failed to load {} from {}",
name, inpf.GetName());
182template <
typename DataT>
187 auto&&
w = std::setw(9);
190 for (
unsigned int iz = 0; iz < mPhiVertices; ++iz) {
194 for (
unsigned int ix = 1; ix < mZVertices; ++ix) {
195 stream <<
", " <<
w << (*this)(ix, 0, iz);
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);
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);
213 LOGP(info,
"{} \n \n",
stream.str());
216template <
typename DataT>
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)) {
224 return std::pair(indStart, std::vector<float>(
data.begin() + indStart,
data.end()));
226 return std::pair(indStart, std::vector<float>());
229template <
typename DataT>
232 std::transform(mData.begin(), mData.end(), mData.begin(), [
value =
value](
auto&
val) { return val * value; });
236template <
typename DataT>
239 std::transform(mData.begin(), mData.end(),
other.mData.begin(), mData.begin(), std::plus<>());
243template <
typename DataT>
246 std::transform(mData.begin(), mData.end(),
other.mData.begin(), mData.begin(), std::minus<>());
250template <
typename DataT>
253 std::transform(mData.begin(), mData.end(),
other.mData.begin(), mData.begin(), std::multiplies<>());
257template <
typename DataT>
260 const size_t iphi =
index / (nz * nr);
261 index -= (iphi * nz * nr);
262 const size_t iz =
index % nz;
266template <
typename DataT>
269 const size_t iphi =
index / (nz * nr);
270 index -= (iphi * nz * nr);
275template <
typename DataT>
278 return index / (nz * nr);
281template <
typename DataT>
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");
292template <
typename DataT>
296 tree->SetAlias(
"val",
"_0");
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");
311template <
typename DataT>
318 mData.resize(nZ * nR *
static_cast<size_t>(nPhi));
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)
325 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
326 ROOT::DisableImplicitMT();
328 ROOT::EnableImplicitMT(nthreads);
329 ROOT::RDataFrame dFrame(treename, fileIn);
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) {
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();
344 iphi.reserve(nvalues);
347 phi.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;
355 if ((rangeiZ.first < rangeiZ.second) && ((iZTmp < rangeiZ.first) || (iZTmp > rangeiZ.second))) {
360 if ((rangeiR.first < rangeiR.second) && ((iRTmp < rangeiR.first) || (iRTmp > rangeiR.second))) {
365 if ((rangeiPhi.first < rangeiPhi.second) && ((iPhiTmp < rangeiPhi.first) || (iPhiTmp > rangeiPhi.second))) {
373 const float x = rTmp * std::cos(phiTmp);
374 const float y = rTmp * std::sin(phiTmp);
376 unsigned char secNum = simOneSectorOnly ? 0 : std::floor(phiTmp /
SECPHIWIDTH);
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);
390 return std::make_tuple(vals, iz,
ir, iphi,
z,
r, phi, lPos, globalIdx);
392 {treename.data(),
"nz",
"nr",
"nphi"});
395 ROOT::RDF::RSnapshotOptions opt;
397 df.Snapshot(treename, fileOut, {
"slice"}, opt);
400template <
typename DataT>
404 return interpolator(
z,
r, phi);
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)
410 if (ROOT::IsImplicitMTEnabled() && (ROOT::GetThreadPoolSize() != nthreads)) {
411 ROOT::DisableImplicitMT();
413 ROOT::EnableImplicitMT(nthreads);
414 ROOT::RDataFrame dFrame(nPhi);
417 unsigned short nr, nz, nphi;
418 if (!getVertices(treename, fileIn, nr, nz, nphi)) {
425 data.initFromFile(fileIn, treename, nthreads);
431 auto interpolate = [&mGrid3D = std::as_const(mGrid3D), &
data = std::as_const(
data), rangeR, rangeZ, rangePhi, nR, nZ, nPhi](
unsigned int, ULong64_t iPhi) {
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;
444 iphi.reserve(nvalues);
447 phi.reserve(nvalues);
449 vals.reserve(nvalues);
450 lPos.reserve(nvalues);
451 globalIdx.reserve(nvalues);
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;
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));
462 const DataT zPos = rangeZ.first + iZ * zSpacing;
464 iphi.emplace_back(iPhi);
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));
470 globalIdx.emplace_back(idx);
471 const float x = rPos * std::cos(phiPos);
472 const float y = rPos * std::sin(phiPos);
474 unsigned char secNum = simOneSectorOnly ? 0 : std::floor(phiPos /
SECPHIWIDTH);
477 lPos.emplace_back(lPosTmp);
480 return std::make_tuple(vals, iz,
ir, iphi,
z,
r, phi, lPos, globalIdx);
484 auto dfStore = dFrame.DefineSlotEntry(treename, interpolate);
487 ROOT::RDF::RSnapshotOptions opt;
492 dfStore.Snapshot(treename, fileOut, {treename.data()}, opt);
496template <
typename DataT>
499 TFile fTmp(fileIn.data(),
"READ");
500 TTree*
tree = (TTree*)fTmp.Get(treename.data());
502 LOGP(warning,
"Tree {} not found in input file {}", treename, fileIn);
505 tree->SetBranchAddress(
"nz", &nZ);
506 tree->SetBranchAddress(
"nr", &nR);
507 tree->SetBranchAddress(
"nphi", &nPhi);
This class provides a simple method to store values on a large 3-Dim grid with ROOT io functionality.
Definition of RegularGrid3D class.
Definition of TriCubic class.
static LocalPosition3D GlobalToLocal(const GlobalPosition3D &pos, const double alpha)
GLuint const GLchar * name
GLsizei const GLfloat * value
GLenum GLsizei GLsizei GLint * values
GLubyte GLubyte GLubyte GLubyte w
GLdouble GLdouble GLdouble z
Global TPC definitions and constants.
constexpr double SECPHIWIDTH
constexpr unsigned char SECTORSPERSIDE
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()))