CbmRoot
Loading...
Searching...
No Matches
CbmQaOnlineInterface.h
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
10// TODO: move to more suitable place (histoserv, for example)
11
12#pragma once
13
14#include "Histogram.h"
15#include "Logger.h"
16#include "TH1D.h"
17#include "TH2D.h"
18#include "TProfile.h"
19#include "TProfile2D.h"
20
21//#include <log.hpp>
22#include <type_traits>
23
24namespace
25{
30} // namespace
31
32namespace cbm::qa
33{
40 template<class RootHistogram>
41 class RootHistogramAccessor : public RootHistogram {
42
43 public:
44 template<class... Args>
45 RootHistogramAccessor(Args... args) : RootHistogram(args...)
46 {
47 }
48
53 template<class SourceQaHistogram>
54 void AddSliceFromQaHistogram(const SourceQaHistogram& histo, double val);
55
58 template<class QaHistogram>
59 void SetFromQaHistogram(const QaHistogram& histo);
60 };
61
65 public:
66 // TODO: generalize
71 static void AddSlice(const H1D& src, double value, TH2D* dst);
72
77 static void AddSlice(const Prof1D& src, double value, TProfile2D* dst);
78
82 static TH1D* ROOTHistogram(const H1D& hist);
83
87 static TH2D* ROOTHistogram(const H2D& hist);
88
92 static TProfile* ROOTHistogram(const Prof1D& prof);
93
97 static TProfile2D* ROOTHistogram(const Prof2D& prof);
98 };
99} // namespace cbm::qa
100
101
102namespace cbm::qa
103{
104 // -------------------------------------------------------------------------------------------------------------------
105 //
106 template<class RootHistogram>
107 template<class QaHistogram>
109 {
110 constexpr bool IsDstH2D = std::is_same_v<RootHistogram, TH2D>;
111 constexpr bool IsDstProf2D = std::is_same_v<RootHistogram, TProfile2D>;
112 constexpr bool IsSrcH1D = std::is_same_v<QaHistogram, H1D>;
113 constexpr bool IsSrcProf1D = std::is_same_v<QaHistogram, Prof1D>;
114
115 static_assert((IsDstH2D && IsSrcH1D) || (IsDstProf2D && IsSrcProf1D),
116 "function not implemented for a given (RootHistogram, QaHistogram) pair");
117
118 auto iBinX = this->GetXaxis()->FindBin(val);
119 for (int iBinY = 0; iBinY <= this->GetNbinsY() + 1; ++iBinY) {
120 int iGlobBin = this->GetBin(iBinX, iBinY);
121 auto binSrc = src.GetBinAccumulator(iBinY);
122 if constexpr (IsSrcH1D) {
123 this->fArray[iGlobBin] += binSrc.value();
124 if (this->fSumw2.fN) {
125 this->fSumw2.fArray[iGlobBin] += binSrc.variance();
126 }
127 }
128 else if constexpr (IsSrcProf1D) {
129 this->fArray[iGlobBin] += binSrc.GetSumWV();
130 this->fBinEntries.fArray[iGlobBin] += binSrc.GetSumW();
131 this->fSumw2.fArray[iGlobBin] += binSrc.GetSumWV2();
132 this->fBinSumw2.fArray[iGlobBin] += binSrc.GetSumW2();
133 this->fTsumwz += binSrc.GetSumWV();
134 this->fTsumwz2 += binSrc.GetSumWV2();
135 }
136 }
137 this->fTsumw += src.GetTotSumW();
138 this->fTsumw2 += src.GetTotSumW2();
139 this->fTsumwx += src.GetTotSumWX();
140 this->fTsumwx2 += src.GetTotSumWX2();
141 this->SetEntries(this->GetEntries() + src.GetEntries());
142 }
143
144 // -------------------------------------------------------------------------------------------------------------------
145 //
146 template<class RootHistogram>
147 template<class QaHistogram>
149 {
150 // Set individual bin properties
151 if constexpr (std::is_same_v<RootHistogram, TProfile> || std::is_same_v<RootHistogram, TProfile2D>) {
152 this->Sumw2();
153 auto SetBin = [&](int iBinRoot, const auto& binQa) {
154 this->fArray[iBinRoot] = binQa.GetSumWV();
155 this->fBinEntries.fArray[iBinRoot] = binQa.GetSumW();
156 this->fSumw2.fArray[iBinRoot] = binQa.GetSumWV2();
157 this->fBinSumw2.fArray[iBinRoot] = binQa.GetSumW2();
158 };
159 if constexpr (std::is_same_v<RootHistogram, TProfile>) {
160 for (int iBinX = 0; iBinX <= this->GetNbinsX() + 1; ++iBinX) {
161 int iBin = this->GetBin(iBinX);
162 auto bin = hist.GetBinAccumulator(iBinX);
163 SetBin(iBin, bin);
164 this->fTsumwy += bin.GetSumWV();
165 this->fTsumwy2 += bin.GetSumWV2();
166 }
167 }
168 else if constexpr (std::is_same_v<RootHistogram, TProfile2D>) {
169 for (int iBinX = 0; iBinX <= this->GetNbinsX() + 1; ++iBinX) {
170 for (int iBinY = 0; iBinY <= this->GetNbinsY() + 1; ++iBinY) {
171 int iBin = this->GetBin(iBinX, iBinY);
172 auto bin = hist.GetBinAccumulator(iBinX, iBinY);
173 SetBin(iBin, bin);
174 this->fTsumwz += bin.GetSumWV();
175 this->fTsumwz2 += bin.GetSumWV2();
176 }
177 }
178 }
179 }
180 else if constexpr (std::is_same_v<RootHistogram, TH1D> || std::is_same_v<RootHistogram, TH2D>) {
181 if (hist.GetTotSumW() != hist.GetTotSumW2()) { // some of weights were not equal to 1.
182 if (!this->fSumw2.fN) {
183 this->Sumw2();
184 }
185 }
186 auto SetBin = [&](int iBinRoot, const auto& binQa) {
187 this->fArray[iBinRoot] = binQa.value();
188 if (this->fSumw2.fN) {
189 this->fSumw2.fArray[iBinRoot] = binQa.variance();
190 }
191 };
192 if constexpr (std::is_same_v<RootHistogram, TH1D>) {
193 for (int iBinX = 0; iBinX <= this->GetNbinsX() + 1; ++iBinX) {
194 SetBin(this->GetBin(iBinX), hist.GetBinAccumulator(iBinX));
195 }
196 }
197 else if constexpr (std::is_same_v<RootHistogram, TH2D>) {
198 for (int iBinX = 0; iBinX <= this->GetNbinsX() + 1; ++iBinX) {
199 for (int iBinY = 0; iBinY <= this->GetNbinsY() + 1; ++iBinY) {
200 SetBin(this->GetBin(iBinX, iBinY), hist.GetBinAccumulator(iBinX, iBinY));
201 }
202 }
203 }
204 }
205
206 // Set total sums
207 if constexpr (std::is_base_of_v<TH2D, RootHistogram>) {
208 this->fTsumwy = hist.GetTotSumWY();
209 this->fTsumwy2 = hist.GetTotSumWY2();
210 this->fTsumwxy = hist.GetTotSumWXY();
211 }
212 this->fTsumw = hist.GetTotSumW();
213 this->fTsumw2 = hist.GetTotSumW2();
214 this->fTsumwx = hist.GetTotSumWX();
215 this->fTsumwx2 = hist.GetTotSumWX2();
216 this->SetEntries(hist.GetEntries());
217 }
218} // namespace cbm::qa
ROOT-free implementation of a histogram.
1D-histogram
2D-histogram
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 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 SetFromQaHistogram(const QaHistogram &histo)
Sets fields from qa histogram.
void AddSliceFromQaHistogram(const SourceQaHistogram &histo, double val)
Sets slice from the lower-dimension histogram.