CbmRoot
Loading...
Searching...
No Matches
CbmRichGeoTestOpt.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2021 UGiessen/JINR-LIT, Giessen/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
10#include "CbmRichGeoTestOpt.h"
11
12#include "CbmDrawHist.h"
13#include "TBox.h"
14#include "TFile.h"
15#include "TLine.h"
16#include "TStyle.h"
17#include "TText.h"
18#include "utils/CbmUtils.h"
19
20#include <iostream>
21#include <utility>
22
23using namespace std;
24
25CbmRichGeoTestOpt::CbmRichGeoTestOpt() : fReferenceInd(-1), fDrawReference(false), fHM(nullptr), fOutputDir("") {}
26
28
29void CbmRichGeoTestOpt::SetFilePathes(vector<string> geoTestBoxPathes, vector<string> geoTestOmega3Pathes,
30 vector<string> geoTestOmega8Pathes, vector<string> urqmdTestPathes,
31 vector<string> recoQaBoxPathes, vector<string> recoQaUrqmdPathes)
32{
33 fGeoTestBoxPathes = geoTestBoxPathes;
34 fGeoTestOmega3Pathes = geoTestOmega3Pathes;
35 fGeoTestOmega8Pathes = geoTestOmega8Pathes;
36 fUrqmdTestPathes = urqmdTestPathes;
37 fRecoQaBoxPathes = recoQaBoxPathes;
38 fRecoQaUrqmdPathes = recoQaUrqmdPathes;
39}
40
42
44{
45 string path = "";
46 if (fileEnum == kGeoTestBoxFile) path = fGeoTestBoxPathes[iFile];
47 if (fileEnum == kGeoTestOmega3File) path = fGeoTestOmega3Pathes[iFile];
48 if (fileEnum == kGeoTestOmega8File) path = fGeoTestOmega8Pathes[iFile];
49 if (fileEnum == kUrqmdTestFile) path = fUrqmdTestPathes[iFile];
50 if (fileEnum == kRecoQaBoxFile) path = fRecoQaBoxPathes[iFile];
51 if (fileEnum == kRecoQaUrqmdFile) path = fRecoQaUrqmdPathes[iFile];
52
53 return path;
54}
55
57{
58 string path = "";
59 if (fileEnum == kGeoTestBoxFile) return "kGeoTestBoxFile";
60 if (fileEnum == kGeoTestOmega3File) return "kGeoTestOmega3File";
61 if (fileEnum == kGeoTestOmega8File) return "kGeoTestOmega8File";
62 if (fileEnum == kUrqmdTestFile) return "kUrqmdTestFile";
63 if (fileEnum == kRecoQaBoxFile) return "kRecoQaBoxFile";
64 if (fileEnum == kRecoQaUrqmdFile) return "kRecoQaUrqmdFile";
65
66 return path;
67}
68
69pair<double, double> CbmRichGeoTestOpt::H1MeanRms(CbmRichGeoTestOptFileEnum fileEnum, int iFile, const string& histName)
70{
71 string path = GetFilePath(fileEnum, iFile);
72 if (path == "") return make_pair(0., 0.);
73
75 TFile* oldFile = gFile;
76 TDirectory* oldDir = gDirectory;
77
78 TFile* file = new TFile(path.c_str(), "READ");
79 if (file == nullptr) return make_pair(0., 0.);
80 TH1D* hist = file->Get<TH1D>(histName.c_str());
81 if (hist == nullptr) return make_pair(0., 0.);
82 double mean = hist->GetMean();
83 double rms = hist->GetRMS();
84
86 gFile = oldFile;
87 gDirectory = oldDir;
88
89 file->Close();
90 delete file;
91 return make_pair(mean, rms);
92}
93
94pair<double, double> CbmRichGeoTestOpt::H2ProjYMeanRms(CbmRichGeoTestOptFileEnum fileEnum, int iFile,
95 const string& histName)
96{
97 string path = GetFilePath(fileEnum, iFile);
98 if (path == "") return make_pair(0., 0.);
99
101 TFile* oldFile = gFile;
102 TDirectory* oldDir = gDirectory;
103
104 TFile* file = new TFile(path.c_str(), "READ");
105 if (file == nullptr) return make_pair(0., 0.);
106 TH2D* hist = file->Get<TH2D>(histName.c_str());
107 if (hist == nullptr) return make_pair(0., 0.);
108 TH1D* py = hist->ProjectionY((histName + to_string(iFile) + "_py").c_str());
109 double mean = py->GetMean();
110 double rms = py->GetRMS();
111
113 gFile = oldFile;
114 gDirectory = oldDir;
115
116 file->Close();
117 delete file;
118 return make_pair(mean, rms);
119}
120
121double CbmRichGeoTestOpt::HEntries(CbmRichGeoTestOptFileEnum fileEnum, int iFile, const string& histName)
122{
123 string path = GetFilePath(fileEnum, iFile);
124 if (path == "") return 0.;
125
127 TFile* oldFile = gFile;
128 TDirectory* oldDir = gDirectory;
129
130 TFile* file = new TFile(path.c_str(), "READ");
131 if (file == nullptr) return 0.;
132 TH1* hist = file->Get<TH1>(histName.c_str());
133 if (hist == nullptr) return 0.;
134 double entries = hist->GetEntries(); //hist->Integral();
135
137 gFile = oldFile;
138 gDirectory = oldDir;
139
140 file->Close();
141 delete file;
142 return entries;
143}
144
145
146void CbmRichGeoTestOpt::DrawMeanRms(CbmRichGeoTestOptFileEnum fileEnum, const string& histName,
147 CbmRichGeoTestOptHistEnum histEnum, const string& titleY, double minY, double maxY,
148 int nofParts, int nofFilesPart)
149{
150 string canvasName = "richgeoopt_sim_" + GetFileEnumName(fileEnum) + "_" + histName + to_string(nofParts);
151 fHM->CreateCanvas(canvasName, canvasName, 900, 600);
152 vector<TH1*> hists;
153 vector<string> legend;
154 double refValue = 0.;
155 for (int iP = 0; iP < nofParts; iP++) {
156 string histNameInd =
157 GetFileEnumName(fileEnum) + "_sim_" + histName + "_" + to_string(nofParts) + "_" + to_string(iP);
158 fHM->Create1<TH1D>(histNameInd, histNameInd + ";Geometry;" + titleY, nofFilesPart, .5, nofFilesPart + .5);
159 for (int iF = 0; iF < nofFilesPart; iF++) {
160 double value = 0.;
161 int ind = iP * nofFilesPart + iF;
162 if (histEnum == kH1MeanHist)
163 value = H1MeanRms(fileEnum, ind, histName).first;
164 else if (histEnum == kH1RmsHist)
165 value = H1MeanRms(fileEnum, ind, histName).second;
166 else if (histEnum == kH2MeanHist)
167 value = H2ProjYMeanRms(fileEnum, ind, histName).first;
168 else if (histEnum == kH2RmsHist)
169 value = H2ProjYMeanRms(fileEnum, ind, histName).second;
170 else if (histEnum == kHEntriesHist)
171 value = HEntries(fileEnum, ind, histName) / 1000.; // TODO: get nof events
172
173 if (ind + 1 == fReferenceInd) refValue = value;
174 fHM->H1(histNameInd)->SetBinContent(iF + 1, value);
175 }
176 hists.push_back(fHM->H1(histNameInd));
177 legend.push_back(to_string(iP));
178 }
179
180 DrawManyH1(hists, legend, minY, maxY);
181 DrawReferenceLineH1(refValue);
182}
183
184
185void CbmRichGeoTestOpt::DrawMeanEff(CbmRichGeoTestOptFileEnum fileEnum, const string& histName1,
186 const string& histName2, const string& titleY, double minY, double maxY,
187 int nofParts, int nofFilesPart, double effCoeff)
188{
189 string canvasName =
190 "richgeoopt_sim_" + GetFileEnumName(fileEnum) + "_" + histName1 + "_" + histName2 + to_string(nofParts);
191 fHM->CreateCanvas(canvasName, canvasName, 900, 600);
192 vector<TH1*> hists;
193 vector<string> legend;
194 double refValue = 0.;
195 for (int iP = 0; iP < nofParts; iP++) {
196 string histEffName = GetFileEnumName(fileEnum) + "_sim_" + histName1 + "_" + histName2 + "_" + to_string(nofParts)
197 + "_" + to_string(iP);
198 //int nofFilesPart = (GetNofFiles() / nofParts) + 1;
199 fHM->Create1<TH1D>(histEffName, histEffName + ";Geometry;" + titleY, nofFilesPart, .5, nofFilesPart + .5);
200 for (int iF = 0; iF < nofFilesPart; iF++) {
201 int ind = iP * nofFilesPart + iF;
202 double entries1 = HEntries(fileEnum, ind, histName1);
203 double entries2 = HEntries(fileEnum, ind, histName2);
204 double value = (entries2 != 0) ? effCoeff * entries1 / entries2 : 0.;
205 if (ind + 1 == fReferenceInd) refValue = value;
206 fHM->H1(histEffName)->SetBinContent(iF + 1, value);
207 }
208 hists.push_back(fHM->H1(histEffName));
209 legend.push_back(to_string(iP));
210 }
211 DrawManyH1(hists, legend, minY, maxY);
212 DrawReferenceLineH1(refValue);
213}
214
216 CbmRichGeoTestOptHistEnum histEnum, const string& titleZ, double minZ,
217 double maxZ, int precision)
218{
219 string canvasName = "richgeoopt_sim_" + GetFileEnumName(fileEnum) + "_" + histName + "_2d";
220 TCanvas* c = fHM->CreateCanvas(canvasName, canvasName, 2000, 1000);
221 c->Divide(4, 2);
222 vector<string> camTilt = {"0 deg", "3 deg", "6 deg", "9 deg", "12 deg", "15 deg", "18 deg", "21 deg"};
223 int nofParts = camTilt.size();
224 int nofFilesPart = GetNofFiles() / nofParts;
225 int nofX = 9;
226 int nofY = 9;
227 double shift = 0.5 * 200. / 10.;
228 double refBinCenterX = -1., refBinCenterY = -1.;
229 double refBinWidthX = -1., refBinWidthY = -1.;
230 for (int iP = 0; iP < nofParts; iP++) {
231 string histNameInd = GetFileEnumName(fileEnum) + "_sim_" + histName + "_" + to_string(iP) + "_2d";
232 fHM->Create2<TH2D>(histNameInd, histNameInd + ";Shift CamY [mm];Shift CamZ [mm];" + titleZ, nofX, -100. - shift,
233 100 + shift, nofY, -100. - shift, 100 + shift);
234 for (int iX = 0; iX < nofX; iX++) {
235 for (int iY = 0; iY < nofY; iY++) {
236 double value = 0.;
237 int ind = iP * nofFilesPart + iX * nofX + iY;
238 //cout << iP << " " << iX << " " << iY << " " << ind << endl;
239 if (histEnum == kH1MeanHist)
240 value = H1MeanRms(fileEnum, ind, histName).first;
241 else if (histEnum == kH1RmsHist)
242 value = H1MeanRms(fileEnum, ind, histName).second;
243 else if (histEnum == kH2MeanHist)
244 value = H2ProjYMeanRms(fileEnum, ind, histName).first;
245 else if (histEnum == kH2RmsHist)
246 value = H2ProjYMeanRms(fileEnum, ind, histName).second;
247 else if (histEnum == kHEntriesHist)
248 value = HEntries(fileEnum, ind, histName) / 1000.; // TODO: get nof events
249 if (ind + 1 == fReferenceInd) {
250 refBinCenterX = fHM->H2(histNameInd)->GetXaxis()->GetBinCenter(iX + 1);
251 refBinCenterY = fHM->H2(histNameInd)->GetYaxis()->GetBinCenter(iY + 1);
252 refBinWidthX = fHM->H2(histNameInd)->GetXaxis()->GetBinWidth(iX + 1);
253 refBinWidthY = fHM->H2(histNameInd)->GetYaxis()->GetBinWidth(iY + 1);
254 }
255 fHM->H2(histNameInd)->SetBinContent(iX + 1, iY + 1, value);
256 }
257 }
258 c->cd(iP + 1);
259 DrawH2(fHM->H2(histNameInd), kLinear, kLinear, kLinear, "COLZ");
260 DrawTextLabelsH2(fHM->H2(histNameInd), precision);
261 fHM->H2(histNameInd)->GetZaxis()->SetRangeUser(minZ, maxZ);
262 fHM->H2(histNameInd)->SetMarkerSize(1.2);
263 DrawTextOnPad(camTilt[iP], 0.1, 0.9, 0.9, 0.99);
264 if (refBinWidthX > 0. && refBinWidthY > 0.) {
265 DrawReferenceBoxH2(refBinCenterX, refBinCenterY, refBinWidthX, refBinWidthY);
266 refBinWidthX = -1.;
267 refBinWidthY = -1.;
268 }
269 }
270}
271
272void CbmRichGeoTestOpt::DrawMeanEff2D(CbmRichGeoTestOptFileEnum fileEnum, const string& histName1,
273 const string& histName2, const string& titleZ, double minZ, double maxZ,
274 double effCoeff, int precision)
275{
276 string canvasName = "richgeoopt_sim_" + GetFileEnumName(fileEnum) + "_" + histName1 + "_" + histName2 + "_2d";
277 TCanvas* c = fHM->CreateCanvas(canvasName, canvasName, 2000, 1000);
278 c->Divide(4, 2);
279 vector<string> camTilt = {"0 deg", "3 deg", "6 deg", "9 deg", "12 deg", "15 deg", "18 deg", "21 deg"};
280 int nofParts = camTilt.size();
281 int nofFilesPart = GetNofFiles() / nofParts;
282 int nofX = 9;
283 int nofY = 9;
284 double shift = 0.5 * 200. / 10.;
285 double refBinCenterX = -1., refBinCenterY = -1.;
286 double refBinWidthX = -1., refBinWidthY = -1.;
287 for (int iP = 0; iP < nofParts; iP++) {
288 string histNameInd =
289 GetFileEnumName(fileEnum) + "_sim_" + histName1 + "_" + histName2 + "_" + to_string(iP) + "_2d";
290 fHM->Create2<TH2D>(histNameInd, histNameInd + ";Shift CamY [mm];Shift CamZ [mm];" + titleZ, nofX, -100. - shift,
291 100 + shift, nofY, -100. - shift, 100 + shift);
292 for (int iX = 0; iX < nofX; iX++) {
293 for (int iY = 0; iY < nofY; iY++) {
294 int ind = iP * nofFilesPart + iX * nofX + iY;
295 //cout << iP << " " << iX << " " << iY << " " << ind << endl;
296 double entries1 = HEntries(fileEnum, ind, histName1);
297 double entries2 = HEntries(fileEnum, ind, histName2);
298 double value = (entries2 != 0) ? effCoeff * entries1 / entries2 : 0.;
299 if (ind + 1 == fReferenceInd) {
300 refBinCenterX = fHM->H2(histNameInd)->GetXaxis()->GetBinCenter(iX + 1);
301 refBinCenterY = fHM->H2(histNameInd)->GetYaxis()->GetBinCenter(iY + 1);
302 refBinWidthX = fHM->H2(histNameInd)->GetXaxis()->GetBinWidth(iX + 1);
303 refBinWidthY = fHM->H2(histNameInd)->GetYaxis()->GetBinWidth(iY + 1);
304 }
305 fHM->H2(histNameInd)->SetBinContent(iX + 1, iY + 1, value);
306 }
307 }
308 c->cd(iP + 1);
309 DrawH2(fHM->H2(histNameInd), kLinear, kLinear, kLinear, "COLZ");
310 DrawTextLabelsH2(fHM->H2(histNameInd), precision);
311 fHM->H2(histNameInd)->SetMarkerSize(1.2);
312 fHM->H2(histNameInd)->GetZaxis()->SetRangeUser(minZ, maxZ);
313 DrawTextOnPad(camTilt[iP], 0.1, 0.9, 0.9, 0.99);
314 if (refBinWidthX > 0. && refBinWidthY > 0.) {
315 DrawReferenceBoxH2(refBinCenterX, refBinCenterY, refBinWidthX, refBinWidthY);
316 refBinWidthX = -1.;
317 refBinWidthY = -1.;
318 }
319 }
320}
321
322void CbmRichGeoTestOpt::DrawTextLabelsH2(TH2* h, int precision)
323{
324 for (int i = 1; i <= h->GetNbinsX(); i++) {
325 for (int j = 1; j <= h->GetNbinsY(); j++) {
326 auto t = new TText(h->GetXaxis()->GetBinCenter(i), h->GetYaxis()->GetBinCenter(j),
327 Cbm::NumberToString(h->GetBinContent(i, j), precision).c_str());
328 t->SetTextAlign(22);
329 t->SetTextFont(43);
330 t->SetTextSize(12);
331 t->SetTextColor(kBlack);
332 t->Draw();
333 }
334 }
335}
336
337void CbmRichGeoTestOpt::DrawManyH1(const vector<TH1*>& hist, const vector<string>& legend, double minY, double maxY)
338{
339 vector<string> camTilt = {"0 deg", "3 deg", "6 deg", "9 deg", "12 deg", "15 deg", "18 deg", "21 deg"};
340 bool drawLegend = (hist.size() == 1) ? false : true;
341 DrawH1(hist, (hist.size() == camTilt.size()) ? camTilt : legend, kLinear, kLinear, drawLegend, 0.9, 0.75, 0.99, 0.99);
342 gPad->SetLeftMargin(0.1);
343 gPad->SetRightMargin(0.01);
344 hist[0]->GetYaxis()->SetTitleOffset(0.8);
345 hist[0]->GetYaxis()->SetRangeUser(minY, maxY);
346 DrawLines(false, true, minY, maxY);
347}
348
349void CbmRichGeoTestOpt::DrawLines(bool drawCamTilt, bool drawCamY, double minY, double maxY)
350{
351 double nMinY = minY + 0.8 * (maxY - minY);
352 double nMaxY = maxY;
353 // cam Y
354 if (drawCamY) {
355 for (int i = 1; i <= 72; i++) {
356 TLine* line = new TLine(i * 9 + .5, nMinY, i * 9 + .5, nMaxY);
357 line->SetLineColor(kGreen + 2);
358 line->SetLineWidth(2);
359 line->Draw();
360 }
361 }
362
363 // cam Tilt
364 if (drawCamTilt) {
365 for (int i = 1; i <= 8; i++) {
366 TLine* line = new TLine(i * 81 + .5, nMinY, i * 81 + .5, nMaxY);
367 line->SetLineColor(kBlue + 2);
368 line->SetLineWidth(2);
369 line->Draw();
370 }
371 }
372}
373
375{
376 if (!fDrawReference) return;
377
378 TLine* line = new TLine(0., value, 1e5, value);
379 line->SetLineColor(kGreen + 4);
380 line->SetLineWidth(1);
381 line->Draw();
382}
383
384void CbmRichGeoTestOpt::DrawReferenceBoxH2(double centerX, double centerY, double widthX, double widthY)
385{
386 if (!fDrawReference) return;
387
388 TBox* box = new TBox(centerX - 0.5 * widthX, centerY - 0.5 * widthY, centerX + 0.5 * widthX, centerY + 0.5 * widthY);
389 box->SetLineColor(kGreen + 4);
390 box->SetLineWidth(2);
391 box->SetFillStyle(0);
392 box->Draw();
393}
394
396{
397 fHM = new CbmHistManager();
398
400 // int nofFilesPartAll = GetNofFiles() / 1.;
401 int nofFilesPartCamTilt = GetNofFiles() / 8;
402
403 //1D
404 DrawMeanRms(kGeoTestBoxFile, "fhNofHits_hits", kH1MeanHist, "#hits in ring", 27, 31, 8, nofFilesPartCamTilt);
405 DrawMeanRms(kGeoTestBoxFile, "fhBoverAVsMom_points", kH2MeanHist, "B/A (hit fit)", 0.88, 0.95, 8,
406 nofFilesPartCamTilt);
407 DrawMeanRms(kGeoTestBoxFile, "fhDRVsNofHits", kH2RmsHist, "dR.RMS [cm]", 0.29, 0.4, 8, nofFilesPartCamTilt);
408 DrawMeanRms(kGeoTestBoxFile, "fhRadiusVsNofHits", kH2RmsHist, "Radius.RMS [cm]", 0.3, 0.6, 8, nofFilesPartCamTilt);
409 DrawMeanEff(kGeoTestBoxFile, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 88, 92, 8, nofFilesPartCamTilt,
410 100.);
411 DrawMeanRms(kUrqmdTestFile, "fh_nof_hits_per_event", kH1MeanHist, "Nof hits per event", 600, 850, 8,
412 nofFilesPartCamTilt);
413 DrawMeanRms(kUrqmdTestFile, "fh_gamma_nontarget_mom", kHEntriesHist, "e^{#pm}_{not target} from #gamma per event", 20,
414 24, 8, nofFilesPartCamTilt);
415 DrawMeanEff(kRecoQaBoxFile, "hte_Rich_Rich_Electron_Rec_p", "hte_Rich_Rich_Electron_Acc_p",
416 "Ring reco efficiency [%]", 84, 98, 8, nofFilesPartCamTilt, 100.);
417 DrawMeanEff(kRecoQaBoxFile, "hte_StsRich_StsRich_Electron_Rec_p", "hte_StsRich_StsRich_Electron_Acc_p",
418 "STS-RICH matching efficiency [%]", 80, 90, 8, nofFilesPartCamTilt, 100.);
419 DrawMeanEff(kRecoQaUrqmdFile, "hte_Rich_Rich_Electron_Rec_p", "hte_Rich_Rich_Electron_Acc_p",
420 "Ring reco efficiency [%]", 86, 100, 8, nofFilesPartCamTilt, 100.);
421 DrawMeanEff(kRecoQaUrqmdFile, "hte_StsRich_StsRich_Electron_Rec_p", "hte_StsRich_StsRich_Electron_Acc_p",
422 "STS-RICH matching efficiency [%]", 82, 94, 8, nofFilesPartCamTilt, 100.);
423 DrawMeanEff(kRecoQaUrqmdFile, "hps_Rich_All_RecPions_p", "hps_Rich_All_Rec_p", "Pion Suppression (RICH)", 100, 240, 8,
424 nofFilesPartCamTilt, 1.);
425 DrawMeanEff(kRecoQaUrqmdFile, "hpe_StsRich_StsRich_Electron_Rec_p", "hpe_StsRich_StsRich_Electron_Acc_p",
426 "El id (RICH) [%]", 80, 92, 8, nofFilesPartCamTilt, 100.);
427 DrawMeanEff(kGeoTestOmega3File, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 45, 52, 8,
428 nofFilesPartCamTilt, 100.);
429 DrawMeanEff(kGeoTestOmega8File, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 52, 58, 8,
430 nofFilesPartCamTilt, 100.);
431
432 //2D
433 DrawMeanRms2D(kGeoTestBoxFile, "fhNofHits_hits", kH1MeanHist, "#hits in ring", 27, 31, 1);
434 DrawMeanRms2D(kGeoTestBoxFile, "fhBoverAVsMom_points", kH2MeanHist, "B/A (point fit)", 0.88, 0.95, 1);
435 DrawMeanRms2D(kGeoTestBoxFile, "fhDRVsNofHits", kH2RmsHist, "dR.RMS [cm]", 0.29, 0.4, 1);
436 DrawMeanRms2D(kGeoTestBoxFile, "fhRadiusVsNofHits", kH2RmsHist, "Radius.RMS [cm]", 0.3, 0.6, 1);
437 DrawMeanEff2D(kGeoTestBoxFile, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 88, 92, 100., 1);
438 DrawMeanRms2D(kUrqmdTestFile, "fh_nof_hits_per_event", kH1MeanHist, "Nof hits per event", 600, 850, 0);
439 DrawMeanRms2D(kUrqmdTestFile, "fh_gamma_nontarget_mom", kHEntriesHist, "e^{#pm}_{not target} from #gamma per event",
440 20, 24, 1);
441 DrawMeanEff2D(kRecoQaBoxFile, "hte_Rich_Rich_Electron_Rec_p", "hte_Rich_Rich_Electron_Acc_p",
442 "Ring reco efficiency [%]", 84, 98, 100., 1);
443 DrawMeanEff2D(kRecoQaBoxFile, "hte_StsRich_StsRich_Electron_Rec_p", "hte_StsRich_StsRich_Electron_Acc_p",
444 "STS-RICH matching efficiency [%]", 80, 90, 100., 1);
445 DrawMeanEff2D(kRecoQaUrqmdFile, "hte_Rich_Rich_Electron_Rec_p", "hte_Rich_Rich_Electron_Acc_p",
446 "Ring reco efficiency [%]", 86, 100, 100., 1);
447 DrawMeanEff2D(kRecoQaUrqmdFile, "hte_StsRich_StsRich_Electron_Rec_p", "hte_StsRich_StsRich_Electron_Acc_p",
448 "STS-RICH matching efficiency [%]", 82, 94, 100., 1);
449 DrawMeanEff2D(kRecoQaUrqmdFile, "hps_Rich_All_RecPions_p", "hps_Rich_All_Rec_p", "Pion Suppression (RICH)", 100, 240,
450 1., 0);
451 DrawMeanEff2D(kRecoQaUrqmdFile, "hpe_StsRich_StsRich_Electron_Rec_p", "hpe_StsRich_StsRich_Electron_Acc_p",
452 "El id (RICH) [%]", 80, 92, 100., 1);
453 DrawMeanEff2D(kGeoTestOmega3File, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 45, 52, 100., 1);
454 DrawMeanEff2D(kGeoTestOmega8File, "fhAccMomEl", "fhMcMomEl", "Geometrical acceptance [%]", 52, 58, 100., 1);
455
457}
458
459
ClassImp(CbmConverterManager)
void DrawTextOnPad(const string &text, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
void SetDefaultDrawStyle()
void DrawH1(TH1 *hist, HistScale logx, HistScale logy, const string &drawOpt, Int_t color, Int_t lineWidth, Int_t lineStyle, Int_t markerSize, Int_t markerStyle)
void DrawH2(TH2 *hist, HistScale logx, HistScale logy, HistScale logz, const string &drawOpt)
Helper functions for drawing 1D and 2D histograms and graphs.
@ kLinear
Definition CbmDrawHist.h:69
Creates comparison plots for RICH geometry testing.
CbmRichGeoTestOptFileEnum
@ kGeoTestOmega3File
@ kRecoQaBoxFile
@ kUrqmdTestFile
@ kRecoQaUrqmdFile
@ kGeoTestBoxFile
@ kGeoTestOmega8File
CbmRichGeoTestOptHistEnum
@ kH2RmsHist
@ kHEntriesHist
@ kH1RmsHist
@ kH1MeanHist
@ kH2MeanHist
Generates beam ions for transport simulation.
Histogram manager.
TCanvas * CreateCanvas(const std::string &name, const std::string &title, Int_t width, Int_t height)
Create and draw TCanvas and store pointer to it.
void SaveCanvasToImage(const std::string &outputDir, const std::string &options="png,eps")
Save all stored canvases to images.
TH2 * H2(const std::string &name) const
Return pointer to TH2 histogram.
void Create2(const std::string &name, const std::string &title, Int_t nofBinsX, Double_t minBinX, Double_t maxBinX, Int_t nofBinsY, Double_t minBinY, Double_t maxBinY)
Helper function for creation of 2-dimensional histograms and profiles. Template argument is a real ob...
void Create1(const std::string &name, const std::string &title, Int_t nofBins, Double_t minBin, Double_t maxBin)
Helper function for creation of 1-dimensional histograms and profiles. Template argument is a real ob...
TH1 * H1(const std::string &name) const
Return pointer to TH1 histogram.
void DrawMeanRms(CbmRichGeoTestOptFileEnum fileEnum, const string &histName, CbmRichGeoTestOptHistEnum histEnum, const string &titleY, double minY, double maxY, int nofParts, int nofFilesPart)
void SetFilePathes(vector< string > geoTestPathes, vector< string > geoTestOmega3Pathes, vector< string > geoTestOmega8Pathes, vector< string > urqmdTestPathes, vector< string > recoQaBoxPathes, vector< string > recoQaUrqmdPathes)
CbmHistManager * fHM
vector< string > fGeoTestBoxPathes
void DrawReferenceLineH1(double value)
void Draw(Option_t *option="")
string GetFilePath(CbmRichGeoTestOptFileEnum fileType, int iFile)
vector< string > fUrqmdTestPathes
double HEntries(CbmRichGeoTestOptFileEnum fileEnum, int iFile, const string &histName)
void DrawTextLabelsH2(TH2 *h, int precision)
CbmRichGeoTestOpt()
Constructor.
pair< double, double > H2ProjYMeanRms(CbmRichGeoTestOptFileEnum fileType, int iFile, const string &histName)
string GetFileEnumName(CbmRichGeoTestOptFileEnum fileEnum)
vector< string > fRecoQaBoxPathes
void DrawLines(bool drawCamTilt, bool drawCamY, double minY, double maxY)
virtual ~CbmRichGeoTestOpt()
Destructor.
void DrawMeanEff2D(CbmRichGeoTestOptFileEnum fileEnum, const string &histName1, const string &histName2, const string &titleZ, double minZ, double maxZ, double effCoeff, int precision)
void DrawMeanRms2D(CbmRichGeoTestOptFileEnum fileEnum, const string &histName, CbmRichGeoTestOptHistEnum histEnum, const string &titleZ, double minZ, double maxZ, int precision)
vector< string > fRecoQaUrqmdPathes
vector< string > fGeoTestOmega3Pathes
void DrawMeanEff(CbmRichGeoTestOptFileEnum fileEnum, const string &histName1, const string &histName2, const string &titleY, double minY, double maxY, int nofParts, int nofFilesPart, double effCoeff)
void DrawManyH1(const vector< TH1 * > &hist, const vector< string > &legend, double minY, double maxY)
pair< double, double > H1MeanRms(CbmRichGeoTestOptFileEnum fileType, int iFile, const string &histName)
void DrawReferenceBoxH2(double centerX, double centerY, double widthX, double widthY)
vector< string > fGeoTestOmega8Pathes
std::string NumberToString(const T &value, int precision=1)
Definition CbmUtils.h:34
Hash for CbmL1LinkKey.