Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoDebug2DCF.cxx
1/*
2 * FemtoDebug2DCF.cxx
3 *
4 * Created on: 25 lis 2023
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "FemtoDebug2DCF.h"
11#include "FemtoConst.h"
12#include "FemtoPair.h"
13
14#include "StdString.h"
15
16#include <RtypesCore.h>
17#include <TAxis.h>
18#include <TBrowser.h>
19#include <TCanvas.h>
20#include <TH1.h>
21#include <TH2.h>
22#include <TMath.h>
23#include <TROOT.h>
24#include <TString.h>
25#include <iostream>
26#include <stddef.h>
27
28namespace Hal {
29 std::function<void(TH1*, Hal::FemtoPair*)>* FemtoDebug2DCF::fMagicFunction = {nullptr};
31 SetAxisNames(GetNum());
32 SetAxisNames(GetDen());
33 }
34
35 void FemtoDebug2DCF::Browse(TBrowser* b) {
36 gPad->Clear();
37 TVirtualPad* c1 = gPad;
38 Draw("all");
39 gPad = c1;
40 b->Add(fNum);
41 b->Add(fDen);
42 }
43
44 TString FemtoDebug2DCF::HTMLExtract(Int_t counter, TString dir) const {
45 Bool_t batch = gROOT->IsBatch();
46 gROOT->SetBatch(kTRUE);
47 TCanvas* c1 = new TCanvas("canvas", "canvas", 0, 0, 1000, 1500);
48 c1->Divide(2, 2);
49 c1->cd(1);
50 TH2* cf = (TH2*) GetHist(kTRUE);
51 cf->Draw("SURF1");
52 c1->cd(2);
53 TH2* num = (TH2*) GetNum()->Clone("temp_n");
54 num->Draw("SURF1");
55 c1->cd(3);
56 TH2* den = (TH2*) GetDen()->Clone("temp_d");
57 den->Draw("SURF1");
58 c1->SaveAs(Form("%s/divided_%i/divided.root", dir.Data(), counter));
59 gROOT->SetBatch(batch);
60 delete num;
61 delete den;
62 delete cf;
63 TString page = CommonExtract(counter, dir);
64 return page;
65 }
66
67 void FemtoDebug2DCF::FillNumObj(TObject* obj) {
68 FemtoPair* pair = (FemtoPair*) obj;
69 (*fMagicFunction)(fNum, pair);
70 }
71
72 void FemtoDebug2DCF::FillDenObj(TObject* obj) {
73 FemtoPair* pair = (FemtoPair*) obj;
74 (*fMagicFunction)(fDen, pair);
75 }
76
77 FemtoDebug2DCF::FemtoDebug2DCF(TString name, Int_t phibins, Int_t eta_bins) :
78 DividedHisto2D(name, phibins + 1, -0.5, double(phibins) + 0.5, eta_bins + 1, -0.5, double(eta_bins) - 1, 'D') {
79 SetAxisNames(GetNum());
80 SetAxisNames(GetDen());
81 }
82
83 FemtoDebug2DCF::~FemtoDebug2DCF() {}
84
85 void FemtoDebug2DCF::Draw(Option_t* opt) {
86 TString option = opt;
87 if (Hal::Std::FindParam(option, "num", kTRUE)) {
88 GetNum()->Draw("SURF1");
89 } else if (Hal::Std::FindParam(option, "den", kTRUE)) {
90 GetDen()->Draw("SURF1");
91 } else if (Hal::Std::FindParam(option, "all", kTRUE)) {
92 TVirtualPad* c1 = gPad;
93 c1->Divide(2, 2);
94 c1->cd(1);
95 Bool_t rebin = kTRUE;
96 if (Hal::Std::FindParam(option, "hd", kTRUE)) { rebin = kFALSE; }
97 TH2* h = (TH2*) GetHist(kTRUE);
98 Int_t rebinX = 1;
99 Int_t rebinY = 1;
100 if (rebin) {
101 if (h->GetNbinsX() > 20) { rebinX = h->GetNbinsX() / 20; }
102 if (h->GetNbinsY() > 20) { rebinY = h->GetNbinsY() / 20; }
103 h->Rebin2D(rebinX, rebinY);
104 h->Scale(1.0 / ((Double_t) rebinX * rebinY));
105 }
106 h->Draw("COL");
107 c1->cd(2);
108 TH2* num = (TH2*) GetNum()->Clone("temp2dN");
109 if (rebin) num->Rebin2D(rebinX, rebinY);
110 num->Draw("COL");
111 c1->cd(3);
112 TH2* den = (TH2*) GetDen()->Clone("temp2dD");
113 if (rebin) den->Rebin2D(rebinX, rebinY);
114 den->Draw("COL");
115 gPad = c1;
116 } else {
117 GetHist(kFALSE)->Draw("COL");
118 }
119 }
120
121 void FemtoDebug2DCF::SetAxisNames(TH1* h) {
122 if (h == NULL) return;
123 h->GetXaxis()->SetTitle("Xval");
124 h->GetYaxis()->SetTitle("Yval");
125 }
126
127 TObject* FemtoDebug2DCF::Clone(const char* /*c*/) const { return new FemtoDebug2DCF(*this); }
128
129 FemtoDebug2DCF::FemtoDebug2DCF(const FemtoDebug2DCF& other) : Hal::DividedHisto2D(other) {}
130
131} /* namespace Hal */
virtual TString CommonExtract(Int_t counter, TString dir) const
TH1 * GetHist(Bool_t normalized=kTRUE) const
virtual void FillDenObj(TObject *obj)
virtual void FillNumObj(TObject *obj)
virtual TString HTMLExtract(Int_t counter=0, TString dir=" ") const
virtual void Draw(Option_t *opt)
virtual void Browse(TBrowser *b)