CbmRoot
Loading...
Searching...
No Matches
CbmQaCheckerProfile1DHandler.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
14#include "Logger.h"
15#include "TCanvas.h"
16#include "TError.h"
17#include "TFolder.h"
18#include "TGraphAsymmErrors.h"
19#include "TLegend.h"
20#include "TMath.h"
21#include "TMultiGraph.h"
22#include "TObject.h"
23#include "TProfile.h"
24#include "TRatioPlot.h"
25#include "TStyle.h"
26
27#include <limits>
28
30
31// ---------------------------------------------------------------------------------------------------------------------
32//
33Profile1DHandler::Profile1DHandler(int iObject, int iFile, int iDataset) : Hist1DHandler(iObject, iFile, iDataset) {}
34
35// ---------------------------------------------------------------------------------------------------------------------
36//
38{
39 int nVersions = fpObjDB->GetNofVersions();
40 int iDef = fpObjDB->GetDefaultID();
41
42 // ----- Option parsing
43 std::string sOption = opt;
44 for (auto& ch : sOption) {
45 ch = std::tolower(ch);
46 }
47 bool bDrawRatio = sOption.find("r") != std::string::npos;
48 bool bDrawDiff = sOption.find("d") != std::string::npos;
49 int nPads = static_cast<int>(bDrawRatio) + static_cast<int>(bDrawDiff) + 1;
50
51 // ----- Canvas: comparison of original histograms, ratios and subtractions
52 std::string sCanvName = fsBaseName + "_cmp_canvas";
53 std::string sCanvTitle = "Comparison result for \"" + fsBaseName + "\"";
54 fpCanvas = std::make_shared<TCanvas>(sCanvName.data(), sCanvTitle.data(), 1500, 500);
55 fpCanvas->cd();
56
57 TPad* pPadOrig = nullptr;
58 TPad* pPadRatio = nullptr;
59 TPad* pPadDiff = nullptr;
60
61 double xPadSize = 1.0 / nPads;
62 double xNext = 0.;
63
64 pPadOrig = new TPad("p1", "p1", 0.0, 0.0, xPadSize, 1.0);
65 pPadOrig->SetMargin(0.20, 0.05, 0.20, 0.10);
66 fpCanvas->cd();
67 xNext += xPadSize;
68
69 if (bDrawRatio) {
70 pPadRatio = new TPad("p2", "p2", xNext, 0.0, xNext + xPadSize, 1.0);
71 pPadRatio->SetMargin(0.20, 0.05, 0.20, 0.10);
72 fpCanvas->cd();
73 xNext += xPadSize;
74 }
75
76 if (bDrawDiff) {
77 pPadDiff = new TPad("p3", "p3", xNext, 0.0, xNext + xPadSize, 1.0);
78 pPadDiff->SetMargin(0.20, 0.05, 0.20, 0.10);
79 fpCanvas->cd();
80 xNext += xPadSize;
81 }
82
83 const char* title = fvpObjects[0]->GetTitle();
84 const char* titleRatio = Form("Ratio to %s", fpObjDB->GetVersionLabel(iDef).data());
85 const char* titleDiff = Form("Difference from %s", fpObjDB->GetVersionLabel(iDef).data());
86 const char* xAxisTitle = static_cast<TProfile*>(fvpObjects[0])->GetXaxis()->GetTitle();
87 const char* yAxisTitle = static_cast<TProfile*>(fvpObjects[0])->GetYaxis()->GetTitle();
88
89 // Original histograms
90 pPadOrig->cd();
91 auto* pMultiGraphOrig = new TMultiGraph(fsBaseName.data(), title);
92 //auto* pLegend = new TLegend(kLegendSize[0], kLegendSize[1]);
93 for (int iV = 0; iV < nVersions; ++iV) {
94 auto* pCopy = new TGraphAsymmErrors((TH1*) fvpObjects[iV]);
95 pCopy->SetMarkerStyle(20);
96 pCopy->SetTitle(fpObjDB->GetVersionLabel(iV).data());
97 pMultiGraphOrig->Add(pCopy, "P");
98 }
99 pMultiGraphOrig->GetXaxis()->SetTitle(xAxisTitle);
100 pMultiGraphOrig->GetYaxis()->SetTitle(yAxisTitle);
101 pMultiGraphOrig->Draw("A pmc plc");
102 pPadOrig->BuildLegend();
103 //pLegend->Draw();
104
105 auto* pDefault = static_cast<TProfile*>(fvpObjects[iDef])->ProjectionX();
106
107 // Histogram ratios to default
108 if (bDrawRatio) {
109 pPadRatio->cd();
110 auto* pMultiGraphRatio = new TMultiGraph((fsBaseName + "_ratio").data(), titleRatio);
111 for (int iV = 0; iV < nVersions; ++iV) {
112 if (iV == iDef) {
113 continue;
114 }
115 auto* pRatio = static_cast<TProfile*>(fvpObjects[iV])->ProjectionX();
116 auto currErrorLevel = gErrorIgnoreLevel;
117 gErrorIgnoreLevel = kError;
118 auto* pCopy = new TGraphAsymmErrors(pRatio, pDefault, "pois");
119 gErrorIgnoreLevel = currErrorLevel;
120 pCopy->SetMarkerStyle(20);
121 pMultiGraphRatio->Add(pCopy, "P");
122
123 if (pRatio) {
124 delete pRatio;
125 pRatio = nullptr;
126 }
127 }
128 pMultiGraphRatio->GetXaxis()->SetTitle(xAxisTitle);
129 pMultiGraphRatio->GetYaxis()->SetTitle("ratio");
130 pMultiGraphRatio->Draw("A pmc plc");
131 }
132
133 // Histogram ratios to default
134 if (bDrawDiff) {
135 pPadDiff->cd();
136 auto* pMultiGraphDiff = new TMultiGraph((fsBaseName + "_diff").data(), titleDiff);
137 for (int iV = 0; iV < nVersions; ++iV) {
138 if (iV == iDef) {
139 continue;
140 }
141 auto* pDiff = static_cast<TProfile*>(fvpObjects[iV])->ProjectionX();
142 pDiff->Add(pDefault, -1.);
143 auto* pCopy = new TGraphAsymmErrors(pDiff);
144 pCopy->GetYaxis()->SetTitle("difference");
145 pCopy->SetMarkerStyle(20);
146 pMultiGraphDiff->Add(pCopy, "P");
147
148 if (pDiff) {
149 delete pDiff;
150 pDiff = nullptr;
151 }
152 }
153 pMultiGraphDiff->GetXaxis()->SetTitle(xAxisTitle);
154 pMultiGraphDiff->GetYaxis()->SetTitle("difference");
155 pMultiGraphDiff->Draw("A pmc plc");
156 }
157
158 fpCanvas->cd();
159 pPadOrig->Draw();
160 if (pPadRatio) pPadRatio->Draw();
161 if (pPadDiff) pPadDiff->Draw();
162
163 if (pDefault) {
164 delete pDefault;
165 pDefault = nullptr;
166 }
167}
Handler class for 1D-histograms (including TProfile objects) (declaration)
Handler class for 1D-profiles (implementation)
Common definitions for QA-Checker framework.
std::string fsBaseName
Base names of the object.
std::shared_ptr< TCanvas > fpCanvas
Comparison canvas.
std::shared_ptr< ObjectDB > fpObjDB
Pointer to object database.
std::vector< TNamed * > fvpObjects
Vector of objects.
Specification of the handler for TProfile class.
void CreateCanvases(Option_t *opt="") override
Creates object comparison canvas.
Profile1DHandler(int iObject, int iFile, int iDataset)
Constructor.