CbmRoot
Loading...
Searching...
No Matches
CbmQaOnlineInterface.cxx
Go to the documentation of this file.
1/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
11//#include "Histogram.h"
12
13#include "TFile.h"
14#include "TH1D.h"
15#include "TH2D.h"
16#include "TProfile.h"
17#include "TProfile2D.h"
19
20#include <boost/filesystem.hpp>
21
28
29// ---------------------------------------------------------------------------------------------------------------------
30//
31void OnlineInterface::AddSlice(const H1D& src, double value, TH2D* dst)
32{
33 auto* pDst = static_cast<RootHistogramAccessor<TH2D>*>(dst);
34 pDst->AddSliceFromQaHistogram<H1D>(src, value);
35}
36
37// ---------------------------------------------------------------------------------------------------------------------
38//
39void OnlineInterface::AddSlice(const Prof1D& src, double value, TProfile2D* dst)
40{
41 auto* pDst = static_cast<RootHistogramAccessor<TProfile2D>*>(dst);
42 pDst->AddSliceFromQaHistogram<Prof1D>(src, value);
43}
44
45// ---------------------------------------------------------------------------------------------------------------------
46//
48{
49 int nBins = hist.GetNbinsX();
50 double xMin = hist.GetMinX();
51 double xMax = hist.GetMaxX();
52 bool add = TH1::AddDirectoryStatus();
53 TH1::AddDirectory(false);
54 auto* pRes = new RootHistogramAccessor<TH1D>(hist.GetName().c_str(), hist.GetTitle().c_str(), nBins, xMin, xMax);
55 TH1::AddDirectory(add);
56 pRes->SetFromQaHistogram<H1D>(hist);
57 return pRes;
58}
59
60// ---------------------------------------------------------------------------------------------------------------------
61//
63{
64 int nBinsX = hist.GetNbinsX();
65 double xMin = hist.GetMinX();
66 double xMax = hist.GetMaxX();
67 int nBinsY = hist.GetNbinsY();
68 double yMin = hist.GetMinY();
69 double yMax = hist.GetMaxY();
70 bool add = TH1::AddDirectoryStatus();
71 TH1::AddDirectory(false);
72 auto* pRes = new RootHistogramAccessor<TH2D>(hist.GetName().c_str(), hist.GetTitle().c_str(), nBinsX, xMin, xMax,
73 nBinsY, yMin, yMax);
74 TH1::AddDirectory(add);
75 pRes->SetFromQaHistogram<H2D>(hist);
76 return pRes;
77}
78
79// ---------------------------------------------------------------------------------------------------------------------
80//
82{
83 const char* name = prof.GetName().c_str();
84 const char* titl = prof.GetTitle().c_str();
85 int nBinsX = prof.GetNbinsX();
86 double xMin = prof.GetMinX();
87 double xMax = prof.GetMaxX();
88 double yMin = prof.GetMinY();
89 double yMax = prof.GetMaxY();
90 bool add = TH1::AddDirectoryStatus();
91 TH1::AddDirectory(false);
92 auto* pRes = new RootHistogramAccessor<TProfile>(name, titl, nBinsX, xMin, xMax, yMin, yMax);
93 TH1::AddDirectory(add);
94 pRes->SetFromQaHistogram<Prof1D>(prof);
95 return pRes;
96}
97
98// ---------------------------------------------------------------------------------------------------------------------
99//
101{
102 const char* name = prof.GetName().c_str();
103 const char* titl = prof.GetTitle().c_str();
104 int nBinsX = prof.GetNbinsX();
105 double xMin = prof.GetMinX();
106 double xMax = prof.GetMaxX();
107 int nBinsY = prof.GetNbinsY();
108 double yMin = prof.GetMinY();
109 double yMax = prof.GetMaxY();
110 double zMin = prof.GetMinZ();
111 double zMax = prof.GetMaxZ();
112 bool add = TH1::AddDirectoryStatus();
113 TH1::AddDirectory(false);
114 auto* pRes = new RootHistogramAccessor<TProfile2D>(name, titl, nBinsX, xMin, xMax, nBinsY, yMin, yMax, zMin, zMax);
115 TH1::AddDirectory(add);
116 pRes->SetFromQaHistogram<Prof2D>(prof);
117 return pRes;
118}
119
120// ---------------------------------------------------------------------------------------------------------------------
121//
122bool OnlineInterface::ConvertOutput(const std::string& inFilename, std::string outFilename)
123{
124 //* Read instance of the qa::StorableData
125 auto archive = std::make_unique<cbm::algo::qa::DataInputArchive>(inFilename);
126 LOG(info) << "ConvertOutput: reading from input archive " << inFilename;
127 auto desc = archive->descriptor();
128 LOG(info) << " - Time created: " << desc.time_created();
129 LOG(info) << " - Host name : " << desc.hostname();
130 LOG(info) << " - User name : " << desc.username();
131
132 auto data = archive->get();
133 if (!data) {
134 LOG(error) << "ConvertOutput: failed to read qa::Data from archive";
135 return false;
136 }
137
138 if (outFilename.empty()) {
139 auto inPath = boost::filesystem::path(inFilename);
140 outFilename = (inPath.parent_path() / (inPath.stem().string() + ".root")).string();
141 }
142 LOG(info) << "ConvertOutput: storing converted histograms in " << outFilename;
143
144 TFile fileOut(outFilename.c_str(), "RECREATE");
145
146 //* Loop over data contents and store the histograms
147 auto ProcessHistograms = [&fileOut](const auto& container) {
148 for (const auto& histogram : container) {
149 std::string dirName = "";
150 std::string histName = histogram.GetName();
151 LOG(info) << " - processing histogram " << histName;
152
153 // FIXME: optimize directory access in outFile
154 auto slashPos = histogram.GetName().rfind('/');
155 if (slashPos != std::string::npos) { // the histogram is in a subdirectory
156 dirName = histName.substr(0, slashPos);
157 histName = histName.substr(slashPos + 1, std::string::npos);
158 }
159
160 auto* pDir = fileOut.GetDirectory(dirName.c_str());
161 if (!pDir) {
162 pDir = fileOut.mkdir(dirName.c_str());
163 }
164
165 auto* pRootHist = ROOTHistogram(histogram);
166 pRootHist->SetName(histName.c_str());
167 pRootHist->SetDirectory(pDir);
168 }
169 };
170
171 ProcessHistograms(data->H1());
172 ProcessHistograms(data->H2());
173 ProcessHistograms(data->P1());
174 ProcessHistograms(data->P2());
175
176 fileOut.Write();
177 fileOut.Close();
178
179 return true;
180}
Set of tools for online->ROOT QA-objects conversions (header)
Archive reading handler for QaData objects.
static void AddSlice(const H1D &src, double value, TH2D *dst)
Fills a slice of a histogram of a higher dimension for a given value (....)
1D-histogram
Definition Histogram.h:442
2D-histogram
Definition Histogram.h:511
const std::string & GetTitle() const
Gets title.
Definition Histogram.h:349
uint32_t GetNbinsY() const
Gets number of bins for y axis.
Definition Histogram.h:317
double GetMaxY() const
Gets y-axis lower bound.
Definition Histogram.h:331
const std::string & GetName() const
Gets name.
Definition Histogram.h:298
double GetMinY() const
Gets y-axis lower bound.
Definition Histogram.h:324
double GetMaxX() const
Gets x-axis lower bound.
Definition Histogram.h:313
double GetMinX() const
Gets x-axis lower bound.
Definition Histogram.h:310
uint32_t GetNbinsX() const
Gets number of bins for x axis.
Definition Histogram.h:307
double GetMaxY() const
Gets y-axis lower bound.
Definition Histogram.h:683
double GetMinY() const
Gets y-axis lower bound.
Definition Histogram.h:680
double GetMaxZ() const
Gets z-axis lower bound.
Definition Histogram.h:806
double GetMinZ() const
Gets z-axis lower bound.
Definition Histogram.h:803
A collection of tools for online QA object conversions.
static void AddSlice(const H1D &src, double value, TH2D *dst)
Fills a slice of a histogram of a higher dimension for a given value (....)
static bool ConvertOutput(const std::string &inFilename, std::string outFilename="")
Converts histograms from a .qa.out file from the online binary to the ROOT format and stores them int...
static TH1D * ROOTHistogram(const H1D &hist)
Converts histogram H1D to ROOT histogram TH1D.
Helper class to access protected fields of TH1, TH2, TProfile and TProfile2D.
void AddSliceFromQaHistogram(const SourceQaHistogram &histo, double val)
Sets slice from the lower-dimension histogram.