CbmRoot
Loading...
Searching...
No Matches
CbmQaCheckerHist1DHandler.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
11
12#include "CbmQaCheckerResult.h"
14#include "Logger.h"
15#include "TCanvas.h"
16#include "TError.h"
17#include "TFolder.h"
18#include "TGraphAsymmErrors.h"
19#include "TH1.h"
20#include "TLegend.h"
21#include "TMath.h"
22#include "TMultiGraph.h"
23#include "TObject.h"
24#include "TStyle.h"
25
26#include <bitset>
27#include <limits>
28
32
33// ---------------------------------------------------------------------------------------------------------------------
34//
35Hist1DHandler::Hist1DHandler(int iObject, int iFile, int iDataset) : ObjectHandler(iObject, iFile, iDataset, "TH1") {}
36
37// ---------------------------------------------------------------------------------------------------------------------
38//
39
40
41// ---------------------------------------------------------------------------------------------------------------------
42// TODO: move somewhere else (probably unite with the CbmQaCompare in a future merger request)
43Result Hist1DHandler::Compare(int iVersion) const
44{
45 auto* pNumerator = static_cast<const TH1*>(fvpObjects[iVersion]);
46 auto* pDenominator = static_cast<const TH1*>(fvpObjects[fpObjDB->GetDefaultID()]);
47 if (pNumerator->GetNbinsX() != pDenominator->GetNbinsX() || pNumerator->GetNbinsY() != pDenominator->GetNbinsY()
48 || pNumerator->GetNbinsZ() != pDenominator->GetNbinsZ()) {
49 return Result(ECmpInference::Different, false, 0., 0., 0.);
50 }
51
52 bool bCompareExact{fCmpBits[static_cast<uint8_t>(ECmpMethod::Exact)]};
53 bool bCompareRatio{fCmpBits[static_cast<uint8_t>(ECmpMethod::Ratio)]};
54 bool bCompareChi2{fCmpBits[static_cast<uint8_t>(ECmpMethod::Chi2)]};
55
56 bool bEqualExactly{true};
57 bool bEqualByRatio{true};
58
59 double ratioLo{1.};
60 double ratioUp{1.};
61
62 if (bCompareExact || bCompareRatio) {
63 bool bRatioMayBeWrong =
64 false; // handling a case, when the histograms contain different empty bins, but the ratio remains 1
65 for (int iBinX{0}; iBinX <= pNumerator->GetNbinsX(); ++iBinX) {
66 for (int iBinY{0}; iBinY <= pNumerator->GetNbinsY(); ++iBinY) {
67 for (int iBinZ{0}; iBinZ <= pNumerator->GetNbinsZ(); ++iBinZ) {
68 auto numBinContent{pNumerator->GetBinContent(iBinX, iBinY, iBinZ)};
69 auto denBinContent{pDenominator->GetBinContent(iBinX, iBinY, iBinZ)};
70 double ratio{static_cast<double>(numBinContent) / denBinContent};
71 // Check bin content
72 if (!TMath::IsNaN(numBinContent) && !TMath::IsNaN(denBinContent)) {
73 if (numBinContent != denBinContent) {
74 bEqualExactly = false;
75 if (bCompareRatio) {
76 bool numEmpty{numBinContent < 1.e-4};
77 bool denEmpty{denBinContent < 1.e-4};
78 if (numEmpty || denEmpty) {
79 bRatioMayBeWrong &= (!numEmpty || !denEmpty);
80 continue; // Ignoring empty bin in ratio estimation (but only, if it is empty in both histograms)
81 }
82 ratioLo = std::min(ratioLo, ratio);
83 ratioUp = std::max(ratioUp, ratio);
84 }
85 }
86 }
87 else {
88 if (TMath::IsNaN(numBinContent) != TMath::IsNaN(denBinContent)) {
89 bEqualExactly = false;
90 }
91 }
92 auto numBinError{pNumerator->GetBinError(iBinX, iBinY, iBinZ)};
93 auto denBinError{pDenominator->GetBinError(iBinX, iBinY, iBinZ)};
94 // Check bin error
95 if (!TMath::IsNaN(numBinError) && !TMath::IsNaN(denBinError)) {
96 if (numBinError != denBinError) {
97 bEqualExactly = false;
98 }
99 }
100 else {
101 if (TMath::IsNaN(numBinError) != TMath::IsNaN(denBinError)) {
102 bEqualExactly = false;
103 }
104 }
105 }
106 }
107 }
108
109 if (bCompareExact && bEqualExactly) {
110
111 return Result(ECmpInference::StronglyEqual, true, ratioLo, ratioUp, 1.);
112 }
113
114 if (bRatioMayBeWrong && (std::fabs(1 - 1.e-4) < ratioLo) && std::fabs(1 + 1.e+4) > ratioUp) {
115 LOG(warn) << "Hist1DHandler::Compare: file " << fpObjDB->GetInputFileName(fDatasetID, fFileID, 0)
116 << ", object: " << pNumerator->GetName()
117 << " has a ratio equal to 1., but some of bins are still different "
118 << "(empty/non-empty bin case)";
119 }
120
121 if (bCompareRatio) {
122 bEqualByRatio = (ratioLo >= fpObjDB->GetRatioRangeMin() && ratioUp <= fpObjDB->GetRatioRangeMax());
123 }
124 }
125
126 bool bEqualByChi2{true};
127 double pVal{0.};
128 if (bCompareChi2) { // perform only, if the histograms were different
129 auto currErrorLevel{gErrorIgnoreLevel};
130 gErrorIgnoreLevel = kFatal;
131 double chi2{0.};
132 int ndf{0};
133 int igood{0};
134 pVal = pDenominator->Chi2TestX(pNumerator, chi2, ndf, igood);
135 bEqualByChi2 = (pVal >= fpObjDB->GetPvalThreshold());
136 gErrorIgnoreLevel = currErrorLevel;
137 }
138
139 // Do not account for method, if it was not required
140 // strongly eq = (!bCmp_0 v bEq_0) ^ (!bCmp_1 v bEq_1) ^ .. ^ (!bCmp_N v bEq_N)
141 // weakly eq = (!bCmp_0 v bEq_0) v (!bCmp_1 v bEq_1) v .. v (!bCmp_N v bEq_N)
142 bEqualExactly = !bCompareExact || bEqualExactly;
143 bEqualByRatio = !bCompareRatio || bEqualByRatio;
144 bEqualByChi2 = !bCompareChi2 || bEqualByChi2;
145 bool bEqualStrongly = bEqualExactly && bEqualByRatio && bEqualByChi2;
146 bool bEqualWeakly = bEqualExactly || bEqualByRatio || bEqualByChi2;
147 ECmpInference res = bEqualStrongly ? ECmpInference::StronglyEqual
149
150 if (ECmpInference::StronglyEqual != res) {
151 auto ResultMsg = [](bool flag) -> std::string {
152 constexpr const char* kEqual = "\e[1;32mconsistent\e[0m";
153 constexpr const char* kDiff = "\e[1;31minconsistent\e[0m";
154 return flag ? kEqual : kDiff;
155 };
156
157 std::stringstream msg;
158 msg << "Found mismatch: file: " << fpObjDB->GetInputFileName(iVersion, fFileID, fDatasetID)
159 << ", object: " << fpObjDB->GetObject(fFileID, fObjectID).first << "\n";
160 if (bCompareExact) {
161 msg << "\texact -> " << ResultMsg(bEqualExactly) << '\n';
162 }
163 if (bCompareRatio) {
164 msg << "\tratio -> " << ResultMsg(bEqualByRatio) << "(lo: " << ratioLo << ", up: " << ratioUp << ")\n";
165 }
166 if (bCompareChi2) {
167 msg << "\tchi2 -> " << ResultMsg(bEqualByChi2) << "(p-val: " << pVal << ")\n";
168 }
169 msg << "Inference: " << ToString(res);
170
171 LOG(info) << msg.str();
172 return Result(res, bEqualExactly, ratioLo, ratioUp, pVal);
173 }
174
175 return Result(res, bEqualStrongly, ratioLo, ratioUp, pVal);
176}
177
178
179// ---------------------------------------------------------------------------------------------------------------------
180//
182{
183 int nVersions = fpObjDB->GetNofVersions();
184 int iDef = fpObjDB->GetDefaultID();
185
186 // ----- Option parsing
187 std::string sOption = opt;
188 for (auto& ch : sOption) {
189 ch = std::tolower(ch);
190 }
191 bool bDrawRatio = sOption.find("r") != std::string::npos;
192 bool bDrawDiff = sOption.find("d") != std::string::npos;
193
194 // Axis titles
195 TString xAxisTitle = static_cast<TH1*>(fvpObjects[0])->GetXaxis()->GetTitle();
196 TString yAxisTitle = static_cast<TH1*>(fvpObjects[0])->GetYaxis()->GetTitle();
197
198 // Original histograms
199 {
200 TString sCanvN = fsBaseName + "_cmp_canvas";
201 TString sCanvT = TString("Comparison result for \"") + fsBaseName + "\"";
202 fpCanvasCmp = std::make_shared<TCanvas>(sCanvN, sCanvT, 500, 500);
203 fpCanvasCmp->SetMargin(0.20, 0.05, 0.20, 0.10);
204 fpCanvasCmp->cd();
205 auto* pMultiGraphOrig = new TMultiGraph(fsBaseName.c_str(), fvpObjects[0]->GetTitle());
206 for (int iV = 0; iV < nVersions; ++iV) {
207 auto* pCopy = new TGraphAsymmErrors((TH1*) fvpObjects[iV]);
208 pCopy->SetMarkerStyle(20);
209 pCopy->SetTitle(fpObjDB->GetVersionLabel(iV).c_str());
210 pMultiGraphOrig->Add(pCopy, "P");
211 }
212 pMultiGraphOrig->GetXaxis()->SetTitle(xAxisTitle);
213 pMultiGraphOrig->GetYaxis()->SetTitle(yAxisTitle);
214 pMultiGraphOrig->Draw("A pmc plc");
215 fpCanvasCmp->BuildLegend();
216 }
217
218 auto* pDefault = static_cast<TH1*>(fvpObjects[iDef]);
219
220 // Histogram ratios to default
221 if (bDrawRatio) {
222 TString sCanvN = fsBaseName + "_ratio_canvas";
223 TString sCanvT = TString("Comparison result for \"") + fsBaseName + "\" (ratio)";
224 fpCanvasRat = std::make_shared<TCanvas>(sCanvN, sCanvT, 500, 500);
225 fpCanvasRat->SetMargin(0.20, 0.05, 0.20, 0.10);
226 fpCanvasRat->cd();
227 TString titleRatio = Form("Ratio to %s", fpObjDB->GetVersionLabel(iDef).c_str());
228 auto* pMultiGraphRatio = new TMultiGraph((fsBaseName + "_ratio").c_str(), titleRatio);
229 for (int iV = 0; iV < nVersions; ++iV) {
230 if (iV == iDef) {
231 continue;
232 }
233 auto* pRatio = static_cast<TH1*>(fvpObjects[iV]->Clone());
234 pRatio->SetDirectory(fpOutDir);
235 auto currErrorLevel = gErrorIgnoreLevel;
236 gErrorIgnoreLevel = kError;
237 auto* pCopy = new TGraphAsymmErrors(pRatio, pDefault, "pois");
238 gErrorIgnoreLevel = currErrorLevel;
239 pCopy->GetYaxis()->SetTitle("ratio");
240 pCopy->SetTitle(fpObjDB->GetVersionLabel(iV).c_str());
241 pCopy->SetMarkerStyle(20);
242 pMultiGraphRatio->Add(pCopy, "P");
243 if (pRatio) {
244 delete pRatio;
245 pRatio = nullptr;
246 }
247 }
248 pMultiGraphRatio->GetYaxis()->SetTitle("ratio");
249 pMultiGraphRatio->GetXaxis()->SetTitle(xAxisTitle);
250 pMultiGraphRatio->Draw("A pmc plc");
251 fpCanvasRat->BuildLegend();
252 }
253
254 // Histogram differences to default
255 if (bDrawDiff) {
256 TString sCanvN = fsBaseName + "_difference_canvas";
257 TString sCanvT = TString("Comparison result for \"") + fsBaseName + "\" (difference)";
258 fpCanvasDiff = std::make_shared<TCanvas>(sCanvN, sCanvT, 500, 500);
259 fpCanvasDiff->SetMargin(0.20, 0.05, 0.20, 0.10);
260 fpCanvasDiff->cd();
261 TString titleDiff = Form("Difference from %s", fpObjDB->GetVersionLabel(iDef).c_str());
262 auto* pMultiGraphDiff = new TMultiGraph((fsBaseName + "_diff").c_str(), titleDiff);
263 for (int iV = 0; iV < nVersions; ++iV) {
264 if (iV == iDef) {
265 continue;
266 }
267 auto* pDiff = static_cast<TH1*>(fvpObjects[iV]->Clone());
268 pDiff->SetDirectory(fpOutDir);
269 pDiff->Add(pDefault, -1.);
270 auto* pCopy = new TGraphAsymmErrors(pDiff);
271 pCopy->GetYaxis()->SetTitle("difference");
272 pCopy->SetTitle(fpObjDB->GetVersionLabel(iV).c_str());
273 pCopy->SetMarkerStyle(20);
274 pMultiGraphDiff->Add(pCopy, "P");
275 if (pDiff) {
276 delete pDiff;
277 pDiff = nullptr;
278 }
279 }
280 pMultiGraphDiff->GetYaxis()->SetTitle("difference");
281 pMultiGraphDiff->GetXaxis()->SetTitle(xAxisTitle);
282 pMultiGraphDiff->Draw("A pmc plc");
283 fpCanvasDiff->BuildLegend();
284 }
285}
Handler class for 1D-histograms (including TProfile objects) (declaration)
Common definitions for QA-Checker framework.
Hist1DHandler(int iObject, int iFile, int iDataset)
Constructor.
void CreateCanvases(Option_t *opt="") override
Creates object comparison canvas.
Result Compare(int iVersion) const override
Compares objects to default.
std::string fsBaseName
Base names of the object.
std::shared_ptr< TCanvas > fpCanvasCmp
Comparison canvas: plots together.
std::shared_ptr< TCanvas > fpCanvasDiff
Comparison canvas: difference with default.
std::shared_ptr< ObjectDB > fpObjDB
Pointer to object database.
std::shared_ptr< TCanvas > fpCanvasRat
Comparison canvas: ratios to default.
std::vector< TNamed * > fvpObjects
Vector of objects.
ObjectHandler(int iObject, int iFile, int iDataset, const char *objType="")
Default constructor.
TDirectory * fpOutDir
Pointer to directory.
std::bitset< static_cast< size_t >(ECmpMethod::END)> fCmpBits
Bitset for comparison methods.
A storable result of the QA-checker comparison routine.
ECmpInference
The object comparison inference.
@ Different
Neither of the comparison methods showed equality.
@ WeaklyEqual
At least one of the comparison methods showed equality.
@ StronglyEqual
All the comparison methods gave equality.
@ Exact
exact equality (point-by-point, error-by-error)
@ Chi2
equality within a chi2 hypothesis test
@ Ratio
ratio equality (max and min ratios do not exceed a required range)
std::string ToString(ECmpInference inference)
String representation of the ECmpInference enum.