CbmRoot
Loading...
Searching...
No Matches
CbmLitTofQaReport.cxx
Go to the documentation of this file.
1/* Copyright (C) 2013-2016 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer] */
4
10#include "CbmLitTofQaReport.h"
11
12#include "CbmDrawHist.h"
13#include "CbmHistManager.h"
14#include "CbmReportElement.h"
15#include "CbmUtils.h"
16#include "TCanvas.h"
17#include "TDirectory.h"
18#include "TF1.h"
19#include "TH1.h"
20#include "TLatex.h"
21#include "TMarker.h"
22#include "TPad.h"
23
24#include <boost/assign/list_of.hpp>
25
26#include <iostream>
27#include <map>
28
29using boost::assign::list_of;
31using Cbm::Split;
32using std::cout;
33using std::endl;
34using std::map;
35using std::string;
36using std::vector;
37
39
41
43{
44 Out().precision(3);
45 Out() << R()->DocumentBegin();
46 Out() << R()->Title(0, GetTitle());
47
49
50 Out() << R()->DocumentEnd();
51}
52
54{
55 HM()->NormalizeToIntegralByPattern("hmp_TofTrack_All_.+");
56
57 DrawH2ByPattern("hmp_Tof_Reco_.+_m2p", kLinear, kLinear, kLog, "colz");
58 DrawH2ByPattern("hmp_Tof_RecoMCID_.+_m2p", kLinear, kLinear, kLog, "colz");
59 DrawH2ByPattern("hmp_Tof_RecoAccTof_.+_m2p", kLinear, kLinear, kLog, "colz");
60 DrawH2ByPattern("hmp_Tof_RecoMCIDAccTof_.+_m2p", kLinear, kLinear, kLog, "colz");
61 DrawH1ByPattern("hmp_Tof_TimeZero_.+");
62 DrawH1ByPattern("hmp_Tof_dTime");
63 DrawH1ByPattern("hmp_Tof_Time_FirstTrack");
64
65 DrawH1ByPattern("hmp_TofTrack_All_.+");
66
67 DrawH1ByPattern("hmp_TofTrack_Z");
68 gPad->SetLogy(true);
69 DrawH1ByPattern("hmp_Tof_Z");
70 gPad->SetLogy(true);
71
73}
74
76{
77 vector<string> categories = list_of("Electron")("Muon")("Pion")("Proton")("Kaon");
78 Int_t nofCategories = categories.size();
79 for (Int_t iCat = 0; iCat < nofCategories; iCat++) {
80 string canvasName = "tof_qa_fit_histograms_" + categories[iCat];
81 TCanvas* canvas = CreateCanvas(canvasName.c_str(), canvasName.c_str(), 1200, 1000);
82 canvas->cd(1);
83 string histName = "hmp_Tof_RecoMCID_" + categories[iCat] + "_m2p";
84 TH2* hist = HM()->H2(histName);
85 DrawH2(hist);
86 hist->FitSlicesY();
87 TH1* meanHist = gDirectory->Get<TH1>(string(histName + "_1").c_str()); // mean
88 TH1* sigmaHist = gDirectory->Get<TH1>(string(histName + "_2").c_str()); // sigma
89 Int_t nofBins = meanHist->GetNbinsX();
90 TGraph* upGraph = new TGraph(nofBins);
91 upGraph->GetXaxis()->SetRangeUser(meanHist->GetXaxis()->GetXmin(), meanHist->GetXaxis()->GetXmax());
92 TGraph* downGraph = new TGraph(nofBins);
93 downGraph->GetXaxis()->SetRangeUser(meanHist->GetXaxis()->GetXmin(), meanHist->GetXaxis()->GetXmax());
94 TGraph* meanGraph = new TGraph(nofBins);
95 meanGraph->GetXaxis()->SetRangeUser(meanHist->GetXaxis()->GetXmin(), meanHist->GetXaxis()->GetXmax());
96 for (Int_t iBin = 1; iBin <= nofBins; iBin++) {
97 Double_t p = meanHist->GetBinCenter(iBin);
98 Double_t mean = meanHist->GetBinContent(iBin);
99 Double_t sigma = sigmaHist->GetBinContent(iBin);
100 upGraph->SetPoint(iBin - 1, p, mean + sigma);
101 downGraph->SetPoint(iBin - 1, p, mean - sigma);
102 meanGraph->SetPoint(iBin - 1, p, mean);
103 }
104 std::cout << "Upper function for " << categories[iCat] << std::endl;
105 FitFunction(upGraph);
106 std::cout << "Lower function for " << categories[iCat] << std::endl;
107 FitFunction(downGraph);
108 DrawGraph(meanGraph, kLinear, kLinear, "PSAME");
109 }
110}
111
113{
114 //TF1* f1 = new TF1("f1", "[0]*x*x*x*x+[1]*x*x*x+[2]*x*x+[3]*x+[4]", graph->GetXaxis()->GetXmin(), graph->GetXaxis()->GetXmax());
115 TF1* f1 = new TF1("f1", "[0]*x*x*x*x+[1]*x*x*x+[2]*x*x+[3]*x+[4]", 1., 8.);
116 f1->SetLineColor(kRed);
117 DrawGraph(graph, kLinear, kLinear, "PSAME", kRed, 2, 1, 1, 21);
118 graph->Fit(f1, "RQ");
119 Double_t p0 = f1->GetParameter(0);
120 Double_t p1 = f1->GetParameter(1);
121 Double_t p2 = f1->GetParameter(2);
122 Double_t p3 = f1->GetParameter(3);
123 Double_t p4 = f1->GetParameter(4);
124 std::cout << "Function: " << p0 << "*x^4" << string((p1 > 0) ? "+" : "") << p1 << "*x^3"
125 << string((p2 > 0) ? "+" : "") << p2 << "*x^2" << string((p3 > 0) ? "+" : "") << p3 << "*x"
126 << string((p4 > 0) ? "+" : "") << p4 << std::endl;
127}
128
ClassImp(CbmConverterManager)
void DrawH2(TH2 *hist, HistScale logx, HistScale logy, HistScale logz, const string &drawOpt)
void DrawGraph(TGraph *graph, HistScale logx, HistScale logy, const string &drawOpt, Int_t color, Int_t lineWidth, Int_t lineStyle, Int_t markerSize, Int_t markerStyle)
Helper functions for drawing 1D and 2D histograms and graphs.
@ kLinear
Definition CbmDrawHist.h:69
@ kLog
Definition CbmDrawHist.h:68
Histogram manager.
Create report for TOF QA.
Abstract class for basic report elements (headers, tables, images etc.).
TH2 * H2(const std::string &name) const
Return pointer to TH2 histogram.
void NormalizeToIntegralByPattern(const std::string &pattern)
Normalize histograms to integral which name matches specified pattern.
Create report for TOF QA.
void FitFunction(TGraph *graph)
CbmLitTofQaReport()
Constructor.
virtual void Draw()
Inherited from CbmSimulationReport.
virtual void Create()
Inherited from CbmSimulationReport.
virtual ~CbmLitTofQaReport()
Destructor.
virtual std::string Title(int size, const std::string &title) const =0
Return string with title.
virtual std::string DocumentBegin() const =0
Return string with open tags for document.
virtual std::string DocumentEnd() const =0
Return string with close tags of the document.
std::ostream & Out() const
All text output goes to this stream.
Definition CbmReport.h:66
void SetReportName(const std::string &name)
Definition CbmReport.h:69
void PrintCanvases() const
Print images created from canvases in the report.
const CbmReportElement * R() const
Accessor to CbmReportElement object. User has to write the report using available tags from CbmReport...
Definition CbmReport.h:61
TCanvas * CreateCanvas(const char *name, const char *title, Int_t ww, Int_t wh)
Create canvas and put it to vector of TCanvases. Canvases created with this function will be automati...
Definition CbmReport.cxx:94
Base class for simulation reports.
void DrawH1ByPattern(const std::string &histNamePattern)
Select by pattern TH1 histograms and draw each histogram on separate canvas.
CbmHistManager * HM() const
Return pointer to Histogram manager.
void DrawH2ByPattern(const std::string &histNamePattern, HistScale logx=kLinear, HistScale logy=kLinear, HistScale logz=kLinear, const std::string &drawOpt="")
Select by pattern TH2 histograms and draw each histogram on separate canvas.
vector< string > Split(const string &name, char delimiter)
Definition CbmUtils.cxx:67
std::string NumberToString(const T &value, int precision=1)
Definition CbmUtils.h:34