CbmRoot
Loading...
Searching...
No Matches
CbmQaCompare.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
10#include "CbmQaCompare.h"
11
12#include "CbmQaCanvas.h"
13#include "Logger.h"
14#include "TAxis.h"
15#include "TCanvas.h"
16#include "TError.h"
17#include "TGraphAsymmErrors.h"
18#include "TH1.h"
19#include "TH2.h"
20#include "TLegend.h"
21#include "TMath.h"
22#include "TMultiGraph.h"
23#include "TObject.h"
24#include "TProfile.h"
25#include "TProfile2D.h"
26#include "TStyle.h"
27
28#include <sstream>
29#include <type_traits>
30
31
32// ---------------------------------------------------------------------------------------------------------------------
33//
34template<class Obj>
35CbmQaCompare<Obj>::CbmQaCompare(const Obj* pHistL, const Obj* pHistR, int verbose)
36 : fpObjL(pHistL)
37 , fpObjR(pHistR)
38 , fVerbose(verbose)
39{
40}
41
42// ---------------------------------------------------------------------------------------------------------------------
43//
44template<class Obj>
48
49// ---------------------------------------------------------------------------------------------------------------------
50//
51template<class Obj>
53 const std::string& optStat) const
54{
55 // ---- Parse options
56 bool bCmpPointByPoint = opt.find("p") != std::string::npos;
57 bool bCmpStatTest = opt.find("s") != std::string::npos;
58 bool bCmpRatio = opt.find("r") != std::string::npos;
59
60 Result res;
61
62
63 if constexpr (std::is_base_of_v<TH1, Obj>) {
64 auto CheckAxes = [&](const TAxis* pL, const TAxis* pR, const char* name) {
65 bool ret = true;
66 if (pL->GetNbins() != pR->GetNbins()) {
67 LOG(warn) << "histogram " << name << " has inconsistent bin number with the default histogram";
68 ret = false;
69 }
70 if (pL->GetXmin() != pR->GetXmin()) {
71 LOG(warn) << "histogram " << name << " has inconsistent min of an axis x, y or z with the default histogram";
72 ret = false;
73 }
74 if (pL->GetXmax() != pR->GetXmax()) {
75 LOG(warn) << "histogram " << name << " has inconsistent max of an axis x, y or z with the default histogram";
76 ret = false;
77 }
78 return ret;
79 };
80
81 // ---- Consistency check ------------------------------------------------------------------------------------------
82 if (!(CheckAxes(fpObjL->GetXaxis(), fpObjR->GetXaxis(), fpObjL->GetName())
83 && CheckAxes(fpObjL->GetYaxis(), fpObjR->GetYaxis(), fpObjL->GetName())
84 && CheckAxes(fpObjL->GetZaxis(), fpObjR->GetZaxis(), fpObjL->GetName()))) {
85 res.fConsistent = false;
86 return res;
87 }
88 else {
89 res.fConsistent = true;
90 }
91
92 // ---- Point-by-point comparison ----------------------------------------------------------------------------------
93 if (bCmpPointByPoint || bCmpRatio) {
94 res.fPointByPoint = true;
95 res.fRatioLo = 1.;
96 res.fRatioUp = 1.;
97 for (int iBinX = 0; iBinX <= fpObjL->GetNbinsX(); ++iBinX) {
98 for (int iBinY = 0; iBinY <= fpObjL->GetNbinsY(); ++iBinY) {
99 for (int iBinZ = 0; iBinZ <= fpObjL->GetNbinsZ(); ++iBinZ) {
100 auto numBinContent = fpObjL->GetBinContent(iBinX, iBinY, iBinZ);
101 auto denBinContent = fpObjR->GetBinContent(iBinX, iBinY, iBinZ);
102 double ratio = static_cast<double>(numBinContent) / denBinContent;
103 // Check bin content
104 if (!TMath::IsNaN(numBinContent) && !TMath::IsNaN(denBinContent)) {
105 if (numBinContent != denBinContent) {
106 res.fPointByPoint = false;
107 if (bCmpRatio) {
108 if (res.fRatioLo > ratio) {
109 res.fRatioLo = ratio;
110 }
111 if (res.fRatioUp < ratio) {
112 res.fRatioUp = ratio;
113 }
114 // NOTE: case numBinContent == denBinContent == 0 is accounted as true
115 }
116 }
117 }
118 else {
119 if (TMath::IsNaN(numBinContent) != TMath::IsNaN(denBinContent)) {
120 res.fPointByPoint = false;
121 }
122 }
123 auto numBinError = fpObjL->GetBinError(iBinX, iBinY, iBinZ);
124 auto denBinError = fpObjR->GetBinError(iBinX, iBinY, iBinZ);
125 // Check bin error
126 if (!TMath::IsNaN(numBinError) && !TMath::IsNaN(denBinError)) {
127 if (numBinError != denBinError) {
128 res.fPointByPoint = false;
129 }
130 }
131 else {
132 if (TMath::IsNaN(numBinError) != TMath::IsNaN(denBinError)) {
133 res.fPointByPoint = false;
134 }
135 }
136 }
137 }
138 }
139 }
140 // ---- Stat test comparison ---------------------------------------------------------------------------------------
141 if (bCmpStatTest) {
142 TString optParam = optStat + "CHI2/NDF"; // Set CHI2/NDF as an output
143 res.fChi2NDF = fpObjL->Chi2Test(fpObjR, optParam);
144 }
145 }
146 return res;
147}
148
149// ---------------------------------------------------------------------------------------------------------------------
150//
151template<class Obj>
152TCanvas* CbmQaCompare<Obj>::GetComparisonCanvas(const std::string& opt)
153{
154 // ---- Constant parameters
155 constexpr int kPadPX = 500; // Pad x-size in pixels
156 constexpr int kPadPY = 500; // Pad y-size in pixels
157 constexpr double kPadMarginL = 0.15; // Left margin of a pad
158 constexpr double kPadMarginR = 0.05; // Right margin of a pad
159 constexpr double kPadMarginB = 0.12; // Bottom margin of a pad
160 constexpr double kPadMarginT = 0.10; // Top margin of a pad
161 constexpr Style_t kMarkerL = 24;
162 constexpr Style_t kMarkerR = 25;
163 constexpr Style_t kMarkerDif = 20;
164 constexpr Style_t kMarkerRat = 20;
165
166 // ---- Parse options
167 bool bDrawCmp = opt.find("c") != std::string::npos;
168 bool bDrawRat = opt.find("r") != std::string::npos;
169 bool bDrawDif = opt.find("d") != std::string::npos;
170 int nPads = bDrawCmp + bDrawDif + bDrawRat;
171
172 if (!nPads) {
173 return nullptr;
174 }
175
176 auto* pCanv = new TCanvas("comparison", "comparison", kPadPX * nPads, kPadPY);
177 pCanv->Divide(nPads, 1);
178
179 if constexpr (std::is_base_of_v<TH1, Obj>) {
180 if (bDrawCmp) {
181 pCanv->cd(1);
182 gPad->SetMargin(kPadMarginL, kPadMarginR, kPadMarginB, kPadMarginT);
183 auto* pMGraph = new TMultiGraph(Form("%s_cmp", fpObjL->GetName()), fpObjL->GetTitle());
184 auto* pGrL = new TGraphAsymmErrors(fpObjL);
185 pGrL->SetLineColor(kOrange + 5);
186 pGrL->SetMarkerColor(kOrange + 5);
187 pGrL->SetMarkerStyle(kMarkerL);
188 auto* pGrR = new TGraphAsymmErrors(fpObjR);
189 pGrR->SetLineColor(kCyan - 2);
190 pGrR->SetMarkerColor(kCyan - 2);
191 pGrR->SetMarkerStyle(kMarkerR);
192 pMGraph->Add(pGrL, "P");
193 pMGraph->Add(pGrR, "P");
194 pMGraph->GetXaxis()->SetTitle(fpObjL->GetXaxis()->GetTitle());
195 pMGraph->GetYaxis()->SetTitle(fpObjL->GetYaxis()->GetTitle());
196 pMGraph->Draw("AP");
197 auto* pLegend = new TLegend(0.55, 0.80, (1. - kPadMarginR), (1. - kPadMarginT));
198 pLegend->SetTextSize(0.035);
199 //pLegend->SetBorderSize(0);
200 pLegend->SetFillStyle(0);
201 pLegend->SetMargin(0.2);
202 pLegend->AddEntry(pGrL, fsLabelL.c_str(), "P");
203 pLegend->AddEntry(pGrR, fsLabelR.c_str(), "P");
204 pLegend->Draw();
205 }
206
207 TH1* pHistL = nullptr; // Temporary histogram
208 TH1* pHistR = nullptr; // Temporary histogram
209 if constexpr (std::is_base_of_v<TProfile2D, Obj>) {
210 pHistL = fpObjL->ProjectionXY(Form("tmp_%s", fsLabelL.c_str()));
211 pHistR = fpObjR->ProjectionXY(Form("tmp_%s", fsLabelR.c_str()));
212 }
213 else if constexpr (std::is_base_of_v<TProfile, Obj>) {
214 pHistL = fpObjL->ProjectionX(Form("tmp_%s", fsLabelL.c_str()));
215 pHistR = fpObjR->ProjectionX(Form("tmp_%s", fsLabelR.c_str()));
216 }
217 else {
218 pHistL = static_cast<Obj*>(fpObjL->Clone(Form("tmp_%s", fsLabelL.c_str())));
219 pHistR = static_cast<Obj*>(fpObjR->Clone(Form("tmp_%s", fsLabelR.c_str())));
220 }
221
222 // Ratio plot
223 if (bDrawRat) {
224 pCanv->cd(static_cast<int>(bDrawCmp) + static_cast<int>(bDrawDif) + 1);
225 gPad->SetMargin(kPadMarginL, kPadMarginR, kPadMarginB, kPadMarginT);
226 auto currErrorLevel = gErrorIgnoreLevel; // Suppress log error about skipping null points
227 gErrorIgnoreLevel = kError;
228 auto* pGrRat = new TGraphAsymmErrors(pHistL, pHistR, "pois");
229 gErrorIgnoreLevel = currErrorLevel;
230 pGrRat->GetXaxis()->SetTitle(pHistL->GetXaxis()->GetTitle());
231 pGrRat->GetYaxis()->SetTitle("ratio");
232 pGrRat->SetTitle(Form("Ratio of \"%s\" and \"%s\"", fsLabelL.c_str(), fsLabelR.c_str()));
233 pGrRat->SetMarkerStyle(kMarkerRat);
234 pGrRat->Draw("AP");
235 }
236
237 // Difference plot
238 if (bDrawDif) {
239 pCanv->cd(static_cast<int>(bDrawCmp) + 1);
240 gPad->SetMargin(kPadMarginL, kPadMarginR, kPadMarginB, kPadMarginT);
241 auto* pDif = static_cast<Obj*>(pHistL->Clone());
242 pDif->Add(pHistR, -1.);
243 auto* pGrDif = new TGraphAsymmErrors(pDif);
244 pGrDif->GetXaxis()->SetTitle(pHistL->GetXaxis()->GetTitle());
245 pGrDif->GetYaxis()->SetTitle("difference");
246 pGrDif->SetTitle(Form("Difference of \"%s\" and \"%s\"", fsLabelL.c_str(), fsLabelR.c_str()));
247 pGrDif->SetMarkerStyle(kMarkerDif);
248 pGrDif->Draw("AP");
249 delete pDif;
250 }
251
252 if (pHistL) {
253 delete pHistL;
254 }
255 if (pHistR) {
256 delete pHistR;
257 }
258 }
259
260 return pCanv;
261}
262
263// ---------------------------------------------------------------------------------------------------------------------
264//
265template<class Obj>
266void CbmQaCompare<Obj>::SetObjectLabels(const std::string& labelL, const std::string& labelR)
267{
268 fsLabelL = labelL;
269 fsLabelR = labelR;
270}
271
272template class CbmQaCompare<TH1>;
273template class CbmQaCompare<TH2>;
274template class CbmQaCompare<TProfile>;
Definition of the CbmQaCanvas class.
A histogram comparison module for the QA task (declaration)
Class to compare histograms of the QA task with default ones.
Result operator()(const std::string &opt, const std::string &optStat="UU") const
Compares two histograms, returns a comparison status.
~CbmQaCompare()
Destructor.
void SetObjectLabels(const std::string &labelL, const std::string &labelR)
Set version lables.
TCanvas * GetComparisonCanvas(const std::string &opt)
Creates a comparison canvas.
CbmQaCompare(const Obj *pHistL, const Obj *pHistR, int verbose)
Constructor.
double fRatioUp
Upper bound of the ratio.
bool fConsistent
Consistency flag.
double fRatioLo
Lower bound of the ratio.
double fChi2NDF
Chi2/NDF value.