Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoCorrFunc1D.cxx
1/*
2 * FemtoCorrFunc1D.cxx
3 *
4 * Created on: 26 cze 2023
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "FemtoCorrFunc1D.h"
11
12#include "Array.h"
13#include "Cout.h"
14#include "DividedHisto.h"
15#include "FemtoPair.h"
16#include "HtmlCore.h"
17#include "HtmlFile.h"
18#include "HtmlTable.h"
19#include "ObjectMatrix.h"
20#include "Std.h"
21
22#include <initializer_list>
23
24#include <RtypesCore.h>
25#include <TCanvas.h>
26#include <TH1.h>
27#include <TNamed.h>
28#include <TROOT.h>
29#include <TString.h>
30#include <TSystem.h>
31
32
33namespace Hal {
34 TString FemtoCorrFunc1D::HTMLExtract(Int_t counter, TString dir) const {
35 TString path = Form("%s/corrfunc_%i", dir.Data(), counter);
36 gSystem->MakeDirectory(path);
37 TString filename = Form("%s/corrfunc.html", path.Data());
38 HtmlFile file(filename, kFALSE);
39 HtmlTable table("", "haltable", "");
40 HtmlRow row("", "dark_blue", "");
42 row.AddContent(HtmlCell("value"));
43 table.AddContent(row);
44 for (int i = 0; i < fArray->GetSize(); i++) {
45 DividedHisto1D* h = (DividedHisto1D*) fArray->At(i);
46 HtmlRow row1("", "light_blue", "");
47 if (!Angle()) {
48 row1.AddContent(HtmlCell(Form("%4.2f", fRange.Get(i))));
49 row1.AddContent(HtmlCell(Form("%4.2f", fRange.Get(i + 1))));
50 } else {
51 row1.AddContent(HtmlCell(Form("%4.2f&deg;", fRange.Get(i) * TMath::RadToDeg())));
52 row1.AddContent(HtmlCell(Form("%4.2f&deg;", fRange.Get(i + 1) * TMath::RadToDeg())));
53 }
54 row1.AddContent(HtmlCell(h->HTMLExtract(i, path)));
55 table.AddContent(row1);
56 }
57 if (fNumProp) {
58 Bool_t batch = gROOT->IsBatch();
59 gROOT->SetBatch(kTRUE);
60 TCanvas* c = new TCanvas("canvas", "canvas", 0, 0, 1000, 1500);
61 c->Divide(2, 1);
62 c->cd(1);
63 TH1D* hKt = (TH1D*) fNumProp->Clone();
64 hKt->Draw();
65 c->cd(2);
66 hKt = (TH1D*) fDenProp->Clone();
67 hKt->Draw();
68
69 c->SaveAs(Form("%s/corrfunc_%i/kt.root", dir.Data(), counter));
70 delete c;
71 gROOT->SetBatch(batch);
72 }
73 HtmlRow row2("", "light_blue", "");
74 row2.AddContent(HtmlCellCol("Comment:", 2));
76 table.AddContent(row2);
77 file.AddContent(table);
78 file.AddStringContent(HtmlCore::GetJsDiv("kt.root", "canvas;1"));
79 AddHTMLCode(file);
80 file.Save();
81 return HtmlCore::GetUrl(Form("corrfunc_%i/corrfunc.html", counter), this->ClassName());
82 }
83
85 if (!fNumProp) {
86 Cout::PrintInfo(Form("Missing histograms in %s, did you forget to set MakeHisto in c-tor?", ClassName()), EInfo::kError);
87 }
88 Int_t x_bins = fRange.GetSize() - 1;
89 for (int i = 0; i < x_bins; i++) {
90 DividedHisto1D* h = (DividedHisto1D*) fArray->At(i);
91 TString histname = h->GetName();
92 h->SetName(Form("%s[%i]", histname.Data(), i));
93 }
94 for (int i = 1; i < fRange.GetSize(); i++) {
95 if (fRange[i] < fRange[i - 1]) {
96 Cout::PrintInfo(Form("Wrong order in %s %4.2f<%4.2f", this->ClassName(), fRange[i], fRange[i - 1]),
97 EInfo::kCriticalError);
98 return kFALSE;
99 }
100 }
101 return kTRUE;
102 }
103
104 void FemtoCorrFunc1D::SetBins(const std::initializer_list<double>& init) {
105 if (fRange.GetSize() != (Int_t) init.size()) {
106 Cout::PrintInfo(Form("Wrong size of array during call SetKtBins %i %i", fRange.GetSize(), (int) init.size()),
107 EInfo::kError);
108 return;
109 }
110 std::initializer_list<double>::iterator it;
111 Int_t pos = 0;
112 for (it = init.begin(); it != init.end(); ++it) {
113 fRange[pos++] = *it;
114 }
115 }
116
117 FemtoCorrFunc1D::FemtoCorrFunc1D(const DividedHisto1D& h,
118 std::initializer_list<Double_t> cuts,
119 Int_t bins,
120 TString nameTitle,
121 TString xAxis,
122 TString yAxis) :
123 FemtoCorrFunc(h, cuts.size() - 1) {
124 auto vec = Hal::Std::GetVector(cuts);
125 fRange.MakeBigger(vec.size());
126 for (unsigned int i = 0; i < vec.size(); i++) {
127 fRange.Set(i, vec[i]);
128 }
129 fNumProp = new TH1D(nameTitle + " num", nameTitle + " num", bins, vec[0], vec[vec.size() - 1]);
130 fDenProp = new TH1D(nameTitle + " den", nameTitle + " den", bins, vec[0], vec[vec.size() - 1]);
131 fNumProp->SetXTitle(xAxis);
132 fDenProp->SetXTitle(xAxis);
133 fNumProp->SetYTitle(yAxis);
134 fDenProp->SetYTitle(yAxis);
135 }
136
137 FemtoCorrFunc1D::~FemtoCorrFunc1D() {
138 if (fNumProp) delete fNumProp;
139 if (fDenProp) delete fDenProp;
140 }
141
143 Double_t val = GetPairValNum(pair);
144 Int_t bin = ConvertVal(val);
145 if (bin < 0) return;
146 fNumProp->Fill(val);
147 ((DividedHisto1D*) fArray->At(bin))->FillNumObj(pair);
148 }
149
150 void FemtoCorrFunc1D::FillDenPerfect(FemtoPair* pair) {
151 Double_t val = GetPairValDenPerf(pair);
152 Int_t bin = ConvertVal(val);
153 if (bin < 0) return;
154 fDenProp->Fill(val);
155 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
156 }
157
158 void FemtoCorrFunc1D::FillDenRotated(FemtoPair* pair) {
159 Double_t val = GetPairValDenRot(pair);
160 Int_t bin = ConvertVal(val);
161 if (bin < 0) return;
162 fDenProp->Fill(val);
163 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
164 }
165
166 void FemtoCorrFunc1D::FillDenMixed(FemtoPair* pair) {
167 Double_t val = GetPairValDenMix(pair);
168 Int_t bin = ConvertVal(val);
169 if (bin < 0) return;
170 fDenProp->Fill(val);
171 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
172 }
173
174 void FemtoCorrFunc1D::FillDenHemisphere(FemtoPair* pair) {
175 Double_t val = GetPairValDenHemi(pair);
176 Int_t bin = ConvertVal(val);
177 if (bin < 0) return;
178 fDenProp->Fill(val);
179 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
180 }
181
182 void FemtoCorrFunc1D::FillDenCharged(FemtoPair* pair) {
183 Double_t val = GetPairValDenChar(pair);
184 Int_t bin = ConvertVal(val);
185 if (bin < 0) return;
186 fDenProp->Fill(val);
187 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
188 }
189
190 void FemtoCorrFunc1D::Add(const Object* pack) {
191 FemtoCorrFunc::Add(pack);
192 auto fcf = static_cast<const FemtoCorrFunc1D*>(pack);
193 if (fNumProp) {
194 fNumProp->Add(fcf->fNumProp);
195 fDenProp->Add(fcf->fDenProp);
196 }
197 }
198
199 Int_t FemtoCorrFunc1D::ConvertVal(Double_t val) const {
200 const Int_t size = fRange.GetSize();
201 if (val < fRange.Get(0)) return -1;
202 for (int i = 1; i < size; i++) {
203 if (val < fRange.Get(i)) { return i - 1; }
204 }
205 return -1;
206 }
207
208} /* namespace Hal */
void Set(Int_t index, T val)
Definition Array.h:103
T Get(Int_t i) const
Definition Array.h:97
void MakeBigger(Int_t new_dim)
Definition Array.cxx:5
Int_t GetSize() const
Definition Array.h:50
static void PrintInfo(TString text, Hal::EInfo status)
Definition Cout.cxx:370
virtual TString HTMLExtract(Int_t counter=0, TString dir=" ") const
virtual Bool_t Check()
virtual void Add(const Object *pack)
void FillNum(FemtoPair *pair)
Int_t ConvertVal(Double_t val) const
virtual TString HTMLExtract(Int_t no, TString dir="") const
TString GetComment() const
virtual void Add(const Object *pack)
static TString GetUrl(TString adress, TString text)
Definition HtmlCore.cxx:164
static TString GetJsDiv(TString root_file, TString object_name, TString draw_opt="colz", TString draw_div_name="drawing")
Definition HtmlCore.cxx:226
static TString CheckBr(TString text)
Definition HtmlCore.cxx:202
virtual void AddContent(const HtmlObject &obj)
Definition HtmlTable.cxx:29
virtual void AddContent(const HtmlObject &obj)
Definition HtmlTable.cxx:17
Int_t GetSize() const
TObject * At(Int_t i) const