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