CbmRoot
Loading...
Searching...
No Matches
CbmQaHist.h
Go to the documentation of this file.
1/* Copyright (C) 2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergey Gorbunov [committer] */
4
9
10#ifndef CbmQaHist_H
11#define CbmQaHist_H
12
13#include "CbmQaCanvas.h"
14
15#include <Logger.h>
16
17#include <TFitResultPtr.h>
18#include <TPaveStats.h>
19#include <TStyle.h>
20#include <TVirtualPad.h>
21
22// The TH* headers are needed here for the ROOT linker
23
24#include <TF1.h>
25#include <TH1D.h>
26#include <TH1F.h>
27#include <TH1I.h>
28#include <TProfile.h>
29#include <TProfile2D.h>
30
36template<class HistTypeT>
37class CbmQaHist : public HistTypeT {
38 public:
40 CbmQaHist() : HistTypeT()
41 {
42 // we don't call SetOptStatFit() here since it will call Clone()
43 // which calls this default constructor in an endless recursion
44 if (gStyle) {
45 fOptStat = gStyle->GetOptStat();
46 fOptFit = gStyle->GetOptFit();
47 }
48 this->SetLineWidth(2);
49 }
50
52 CbmQaHist(const CbmQaHist& h) : HistTypeT(h) { SetOptStatFit(h.fOptStat, h.fOptFit); }
53
56 template<typename... Types>
57 CbmQaHist(Types... args) : HistTypeT(args...)
58 {
59 if (gStyle) {
60 SetOptStatFit(gStyle->GetOptStat(), gStyle->GetOptFit());
61 }
62 this->SetLineWidth(2);
63 }
64
67
71 template<typename... Types>
72 TFitResultPtr Fit(Types... args)
73 {
74 TVirtualPad* padsav = gPad;
76 this->Sumw2();
77 auto ret = HistTypeT::Fit(args...);
79
80 // make the output look nice
81
82 if (padsav) padsav->cd();
83 auto* f = this->GetFunction("gaus");
84 if (f) {
85 f->SetParNames("Peak", "#mu", "#sigma");
86 f->SetLineColor(kRed);
87 f->SetLineWidth(3);
88 TPaveStats* st = (TPaveStats*) this->FindObject("stats");
89 if (!st) {
90 LOG(fatal) << "CbmQaHist: can not access statistics of histogram with name \"" << this->GetName() << '\"';
91 }
92 else {
93 st->SetX1NDC(0.6);
94 st->SetX2NDC(0.940);
95 st->SetY1NDC(0.5);
96 st->SetY2NDC(0.930);
97 st->SetOptStat(111110);
98 st->SetOptFit(10001);
99 }
100 }
101 return ret;
102 }
103
105 void SetOptStat(Int_t stat = 1) { SetOptStatFit(stat, fOptFit); }
106
108 void SetOptFit(Int_t fit = 1) { SetOptStatFit(fOptStat, fit); }
109
111 void SetOptStatFit(int stat, int fit)
112 {
113 // the only way to create and auto-size the stat window is to draw the histogram
114
115 fOptStat = stat;
116 fOptFit = fit;
117 if (!gStyle) {
118 return;
119 } // should not happen
120
121 TVirtualPad* savePad = gPad;
122 int saveStat = gStyle->GetOptStat();
123 int saveFit = gStyle->GetOptFit();
124
126 gStyle->SetOptStat(fOptStat);
127 gStyle->SetOptFit(fOptFit);
128
129 this->SetStats(0); // remove the old stat window
130 this->SetStats(1); // set the flag to create the stat window during Draw()
131 this->Draw();
134
135 // restore the environment
136 gStyle->SetOptStat(saveStat);
137 gStyle->SetOptFit(saveFit);
138 if (savePad) savePad->cd();
139 }
140
141 private:
142 int fOptStat = 1;
143 int fOptFit = 0;
144
146};
147
148#endif
Definition of the CbmQaCanvas class.
static CbmQaCanvas & GetDummyCanvas()
a static canvas for temporary drawing
Definition CbmQaCanvas.h:51
void SetOptStatFit(int stat, int fit)
Set stat & fit drawing options and autoresize the stat window.
Definition CbmQaHist.h:111
~CbmQaHist()
Destructor.
Definition CbmQaHist.h:66
CbmQaHist()
Default constructor only for the streamer.
Definition CbmQaHist.h:40
CbmQaHist(Types... args)
Definition CbmQaHist.h:57
TFitResultPtr Fit(Types... args)
Definition CbmQaHist.h:72
int fOptStat
Definition CbmQaHist.h:142
void SetOptStat(Int_t stat=1)
Set stat drawing options and autoresize the stat window.
Definition CbmQaHist.h:105
int fOptFit
Definition CbmQaHist.h:143
ClassDefNV(CbmQaHist, 1)
void SetOptFit(Int_t fit=1)
Set fit drawing options and autoresize the stat window.
Definition CbmQaHist.h:108
CbmQaHist(const CbmQaHist &h)
Copy constructor.
Definition CbmQaHist.h:52
Data class with information on a STS local track.