Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
Femto1DCF.cxx
1/*
2 * Femto1DCF.cxx
3 *
4 * Created on: 12-03-2015
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "Femto1DCF.h"
11
12#include "CorrFit1DCF.h"
13#include "Cout.h"
14#include "Femto1DCFPainter.h"
15#include "FemtoPair.h"
16#include "FemtoSerializationInterface1D.h"
17#include "HtmlCore.h"
18#include "Std.h"
19
20#include <Rtypes.h>
21#include <RtypesCore.h>
22#include <TAttLine.h>
23#include <TAxis.h>
24#include <TBrowser.h>
25#include <TCanvas.h>
26#include <TH1.h>
27#include <TLegend.h>
28#include <TMathBase.h>
29#include <TObjArray.h>
30#include <TObjString.h>
31#include <TROOT.h>
32#include <TSystem.h>
33#include <iostream>
34
35
36namespace Hal {
37
38 void Femto1DCF::DrawScaled(Double_t scale, Option_t* opt) {
39 TH1* h = GetHist(kFALSE);
40 h->Scale(scale);
41 h->Draw(opt);
42 }
43
44 TString Femto1DCF::GetPic() const {
45 TObjArray* clones = new TObjArray();
46 clones->SetOwner(kTRUE);
47 clones->AddLast(new TObjString("divided_0.png"));
48 clones->AddLast(new TObjString("divided_1.png"));
49 clones->AddLast(new TObjString("divided_2.png"));
50 clones->AddLast(new TObjString("divided_3.png"));
51 TString pic_text = HtmlCore::ClickablePic("femto_pic", clones, 996, 1472);
52 delete clones;
53 return pic_text;
54 }
55
57 if (h) {
58 h->GetXaxis()->SetTitle(Femto::KinematicsToAxisLabel(fFrame, 0, 1));
59 h->GetYaxis()->SetTitle(Femto::KinematicsToAxisLabel(fFrame, 1, 1));
60 SetAxisName(Femto::KinematicsToAxisLabel(fFrame, 3, 1));
61 }
62 }
63
64 void Femto1DCF::FillNumObj(TObject* obj) {
65 FemtoPair* pair = (FemtoPair*) obj;
66 if (pair->IsAbs()) {
67 fNum->Fill(TMath::Abs(pair->GetT()), pair->GetWeight());
68 } else {
69 fNum->Fill(pair->GetT(), pair->GetWeight());
70 }
71 }
72
73 void Femto1DCF::FillDenObj(TObject* obj) {
74 FemtoPair* pair = (FemtoPair*) obj;
75 if (pair->IsAbs()) {
76 fDen->Fill(TMath::Abs(pair->GetT()), pair->GetWeight());
77 } else {
78 fDen->Fill(pair->GetT(), pair->GetWeight());
79 }
80 }
81
82 Femto1DCF::Femto1DCF() : DividedHisto1D(), fFrame(Femto::EKinematics::kLCMS) {}
83
84 Femto1DCF::Femto1DCF(const Femto1DCF& other) : DividedHisto1D(other), fFrame(other.fFrame) {}
85
86 Femto1DCF::~Femto1DCF() {
87 if (fPainter) delete fPainter;
88 // TODO Auto-generated destructor stub
89 }
90
91 Femto1DCF::Femto1DCF(TString name, Femto::EKinematics frame) : DividedHisto1D(name, 1), fFrame(frame) {
92 SetName(name);
93 AddLabel(Femto::KinematicsToLabel(fFrame));
94 }
95
96 Femto1DCF::Femto1DCF(TString name, Int_t bins, Double_t min, Double_t max, Femto::EKinematics frame) :
97 DividedHisto1D(name, bins, min, max), fFrame(frame) {
100 fNum->Sumw2();
101 fDen->Sumw2();
102 AddLabel(Femto::KinematicsToLabel(fFrame));
103 }
104
105 void Femto1DCF::AddNum(TH1* h, Option_t* opt) {
108 }
109
110 void Femto1DCF::AddDen(TH1* h, Option_t* opt) {
113 }
114
115 TString Femto1DCF::HTMLExtract(Int_t counter, TString dir) const {
116 TString path = Form("%s/divided_%i", dir.Data(), counter);
117 gSystem->MakeDirectory(path);
118 Bool_t batch = gROOT->IsBatch();
119 gROOT->SetBatch(kTRUE);
120 TCanvas* c1 = new TCanvas("canvas", "canvas", 0, 0, 1000, 1500);
121 c1->Divide(2, 1);
122 c1->cd(1);
123 TH1* cf = GetHist(kTRUE);
124 cf->Draw();
125 c1->cd(2);
126 TH1* num = (TH1*) GetNum()->Clone();
127 TH1* den = (TH1*) GetDen()->Clone();
128 den->SetLineColor(kRed);
129 num->SetLineColor(kGreen);
130 den->Draw();
131 num->Draw("SAME");
132 TLegend* leg = new TLegend(0.4, 0.8, 0.9, 0.9);
133 leg->AddEntry(num, "Numerator", "LPM");
134 leg->AddEntry(den, "Denominator", "LPM");
135 leg->Draw("SAME");
136 c1->SaveAs(Form("%s/divided.root", path.Data()));
137 delete c1;
138 delete num;
139 delete den;
140 delete cf;
141 delete leg;
142 gROOT->SetBatch(batch);
143 TString page = CommonExtract(counter, dir);
144 return page;
145 }
146
147 void Femto1DCF::Browse(TBrowser* b) {
148 gPad->Clear();
149 TVirtualPad* c1 = gPad;
150 Draw("all");
151 gPad = c1;
152 b->Add(fNum);
153 b->Add(fDen);
154 }
155
156 void Femto1DCF::Print(Option_t* opt) const {
157 DividedHisto1D::Print(opt);
158 TString text = Form("Frame : %s", Femto::KinematicsToLabel(fFrame).Data());
159 Cout::Text(text, "L", kWhite);
160 }
161
162 void Femto1DCF::Fit(CorrFit1DCF* fit) { fit->Fit(this); }
163
164 void Femto1DCF::FitDummy(CorrFit1DCF* fit) { fit->FitDummy(this); }
165
166 TObject* Femto1DCF::GetSpecial(TString opt) const {
167 if (opt == "serialization") return new FemtoSerializationInterface1D();
168 if (opt == "painter") return fPainter;
169 return nullptr;
170 }
171
172 void Femto1DCF::Draw(Option_t* option) {
173 TString options = option;
174 if (fPainter) {
175 fPainter->SetOption(option);
176 fPainter->Paint();
177 } else {
178 fPainter = new Hal::Femto1DCFPainter(this);
179 fPainter->SetOption(option);
180 fPainter->Paint();
181 }
182 }
183
184} // namespace Hal
virtual void Fit(TObject *histo)
virtual void FitDummy(TObject *histo)
static void Text(TString text, TString option="L", Color_t color=-1)
Definition Cout.cxx:92
void SetAxisName(TString name)
virtual TString CommonExtract(Int_t counter, TString dir) const
virtual void AddNum(TH1 *num, Option_t *opt="")
void AddLabel(TString label)
TH1 * GetHist(Bool_t normalized=kTRUE) const
virtual void AddDen(TH1 *den, Option_t *opt="")
void FillDenObj(TObject *obj)
Definition Femto1DCF.cxx:73
virtual void AddNum(TH1 *h, Option_t *opt="")
virtual void AddDen(TH1 *h, Option_t *opt="")
virtual void Browse(TBrowser *b)
virtual TString HTMLExtract(Int_t counter=0, TString dir=" ") const
virtual void DrawScaled(Double_t scale, Option_t *opt)
Definition Femto1DCF.cxx:38
void Fit(CorrFit1DCF *fit)
void FillNumObj(TObject *obj)
Definition Femto1DCF.cxx:64
void FitDummy(CorrFit1DCF *fit)
virtual TObject * GetSpecial(TString opt) const
virtual void SetAxisNames(TH1 *h)
Definition Femto1DCF.cxx:56
virtual TString GetPic() const
Definition Femto1DCF.cxx:44
virtual void Draw(Option_t *option="")
Double_t GetT() const
Definition FemtoPair.h:322
Double_t GetWeight() const
Definition FemtoPair.h:282
static TString ClickablePic(TString id_name, TObjArray *strings, Int_t width=796, Int_t height=572)
Definition HtmlCore.cxx:166
void Paint()
Definition Painter.cxx:42
virtual void SetOption(TString option)
Definition Painter.cxx:124