Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
FemtoCorrFunc2D.cxx
1/*
2 * FemtoCorrFunc2Db.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 "FemtoCorrFunc2D.h"
11
12#include "Cout.h"
13#include "FemtoPair.h"
14#include "HtmlCore.h"
15#include "HtmlFile.h"
16#include "HtmlTable.h"
17#include "ObjectMatrix.h"
18#include "Std.h"
19
20#include <RtypesCore.h>
21#include <TCanvas.h>
22#include <TNamed.h>
23#include <TROOT.h>
24#include <TStyle.h>
25#include <TSystem.h>
26
27namespace Hal {
28 FemtoCorrFunc2D::FemtoCorrFunc2D(const DividedHisto1D* h, TString labelX, TString labelY, Int_t sizeX, Int_t sizeY) :
29 FemtoCorrFunc(h, sizeX * sizeY), fLabelX(labelX), fLabelY(labelY) {
30 fRangeX.MakeBigger(sizeX + 1);
31 fRangeY.MakeBigger(sizeY + 1);
32 }
33
34 FemtoCorrFunc2D::FemtoCorrFunc2D(const DividedHisto1D& h, TString labelX, TString labelY, Int_t sizeX, Int_t sizeY) :
35 FemtoCorrFunc(h, sizeX * sizeY), fLabelX(labelX), fLabelY(labelY) {
36 fRangeX.MakeBigger(sizeX + 1);
37 fRangeY.MakeBigger(sizeY + 1);
38 }
39
40 TString FemtoCorrFunc2D::HTMLExtract(Int_t counter, TString dir) const {
41 TString path = Form("%s/corrfunc_%i", dir.Data(), counter);
42 gSystem->MakeDirectory(path);
43 TString filename = Form("%s/corrfunc.html", path.Data());
44 HtmlFile file(filename, kFALSE);
45 HtmlTable table("", "haltable", "");
46 table.SetDefaultClass();
47 HtmlRow row1;
48 row1.AddContent(HtmlCell(Form("%s/%s", fLabelX.Data(), fLabelY.Data())));
49 row1.SetClass(HtmlTableRowClass::DarkBlue());
50 auto GetStringContent = [&](Double_t a, Double_t b, Char_t flag) {
51 TString res = "";
52 if (Angle(flag)) {
53 res = ((Form("%4.2f - %4.2f&deg; ", a * TMath::RadToDeg(), b * TMath::RadToDeg())));
54 } else {
55 res = ((Form("%4.2f - %4.2f ", a, b)));
56 }
57 return res;
58 };
59 for (int i = 0; i < fRangeX.GetSize() - 1; i++) {
60 HtmlCell cell;
61 cell.SetStringContent(GetStringContent(fRangeX.Get(i), fRangeX.Get(i + 1), 'x'));
62 cell.SetClass(Hal::HtmlTableRowClass::DarkBlue());
63 row1.AddContent(cell);
64 }
65 table.AddContent(row1);
66 const Int_t y_size = fRangeY.GetSize() - 1;
67 const Int_t x_size = fRangeX.GetSize() - 1;
68 for (int i = 0; i < y_size; i++) {
69 HtmlRow row2;
70 HtmlCell celltmp;
71 celltmp.SetStringContent(GetStringContent(fRangeY.Get(i), fRangeY.Get(i + 1), 'y'));
72 celltmp.SetClass(HtmlTableRowClass::DarkBlue());
73 row2.AddContent(celltmp);
74 for (int j = 0; j < x_size; j++) {
75 Int_t pos = i * x_size + j;
76 DividedHisto1D* h = (DividedHisto1D*) fArray->At(pos);
77 HtmlCell cell(h->HTMLExtract(pos, path));
78 row2.AddContent(cell);
79 }
80 table.AddContent(row2);
81 }
82 HtmlRow row2("", "light_blue", "");
83 row2.AddContent(HtmlCellCol("Comment:", 2));
84 row2.AddContent(HtmlCellCol(HtmlCore::CheckBr(GetComment()), x_size - 1));
85 table.AddContent(row2);
86 file.AddContent(table);
87 HtmlTable table2("", "haltable", "");
88 table2.SetDefaultClass();
89 HtmlRow row3;
90 row3.SetClass(HtmlTableRowClass::DarkBlue());
91 row3.AddSimpleCells({"Numerator", "Denominator"});
92 table2.AddContent(row3);
93 HtmlRow row4;
94 row4.SetClass(HtmlTableRowClass::LightBlue());
95 HtmlCell cell1, cell2;
96 cell1.SetStringContent(HtmlCore::GetJsDiv("num.root", "Cnum;1", "colz", "dd1"));
97 cell2.SetStringContent(HtmlCore::GetJsDiv("den.root", "Cden;1", "colz", "dd2"));
98 Bool_t batch = gROOT->IsBatch();
99 gStyle->SetOptStat(0);
100 gROOT->SetBatch();
101 TCanvas* cnum = new TCanvas("Cnum");
102 fNumProp->DrawClone("colz");
103 cnum->SaveAs(Form("%s/num.root", path.Data()));
104 TCanvas* cden = new TCanvas("Cden");
105 fDenProp->DrawClone("colz");
106 cden->SaveAs(Form("%s/den.root", path.Data()));
107 delete cnum;
108 delete cden;
109 gROOT->SetBatch(batch);
110 row4.AddContent(cell1);
111 row4.AddContent(cell2);
112 table2.AddContent(row4);
113 file.AddContent(table2);
114 file.Save();
115 return HtmlCore::GetUrl(Form("corrfunc_%i/corrfunc.html", counter), this->ClassName());
116 }
117
118 Bool_t FemtoCorrFunc2D::Check() {
119 for (int i = 1; i < fRangeX.GetSize(); i++) {
120 if (fRangeX[i] < fRangeX[i - 1]) {
121 Cout::PrintInfo(Form("Wrong order in X-axis in %s", this->ClassName()), EInfo::kCriticalError);
122 return kFALSE;
123 }
124 }
125 for (int i = 1; i < fRangeY.GetSize(); i++) {
126 if (fRangeY[i] < fRangeY[i - 1]) {
127 Cout::PrintInfo(Form("Wrong order in Y-axis in %s", this->ClassName()), EInfo::kCriticalError);
128 return kFALSE;
129 }
130 }
131 Int_t x_bins = fRangeX.GetSize() - 1;
132 for (int i = 0; i < x_bins; i++) {
133 for (int j = 0; j < fRangeY.GetSize() - 1; j++) {
134 DividedHisto1D* h = (DividedHisto1D*) fArray->At(i + j * x_bins);
135 TString histname = h->GetName();
136 histname.ReplaceAll(Form("[%i]", i + j * x_bins), Form("[%i][%i]", i, j));
137 }
138 }
139 return InitPropMon();
140 }
141
142 FemtoCorrFunc2D::~FemtoCorrFunc2D() {
143 if (fNumProp) {
144 delete fNumProp;
145 delete fDenProp;
146 }
147 }
148
149 DividedHisto1D* FemtoCorrFunc2D::GetCF(Int_t i, Int_t j) const {
150 Int_t pos = i + j * (fRangeX.GetSize() - 1);
151 { return (DividedHisto1D*) fArray->At(pos); };
152 }
153
154 FemtoCorrFunc::FemtoCorrFunc() : Object(), fArray(nullptr) {}
155
157 auto par = GetPairValNum(pair);
158 Int_t bin = ConvertVal(par);
159 if (bin < 0) return;
160 fNumProp->Fill(par.first, par.second);
161 ((DividedHisto1D*) fArray->At(bin))->FillNumObj(pair);
162 }
163
164 void FemtoCorrFunc2D::FillDenMixed(FemtoPair* pair) {
165 auto par = GetPairValDenMix(pair);
166 Int_t bin = ConvertVal(par);
167 if (bin < 0) return;
168 fDenProp->Fill(par.first, par.second);
169 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
170 }
171
172 void FemtoCorrFunc2D::FillDenRotated(FemtoPair* pair) {
173 auto par = GetPairValDenMix(pair);
174 Int_t bin = ConvertVal(par);
175 if (bin < 0) return;
176 fDenProp->Fill(par.first, par.second);
177 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
178 }
179
180 void FemtoCorrFunc2D::FillDenHemisphere(FemtoPair* pair) {
181 auto par = GetPairValDenHemi(pair);
182 Int_t bin = ConvertVal(par);
183 if (bin < 0) return;
184 fDenProp->Fill(par.first, par.second);
185 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
186 }
187
188 void FemtoCorrFunc2D::FillDenCharged(Hal::FemtoPair* pair) {
189 auto par = GetPairValDenChar(pair);
190 Int_t bin = ConvertVal(par);
191 if (bin < 0) return;
192 fDenProp->Fill(par.first, par.second);
193 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
194 }
195
196 void FemtoCorrFunc2D::FillDenPerfect(Hal::FemtoPair* pair) {
197 auto par = GetPairValDenPerf(pair);
198 Int_t bin = ConvertVal(par);
199 if (bin < 0) return;
200 fDenProp->Fill(par.first, par.second);
201 ((DividedHisto1D*) fArray->At(bin))->FillDenObj(pair);
202 }
203
204 Int_t FemtoCorrFunc2D::ConvertVal(std::pair<Double_t, Double_t> vals) const {
205 const Int_t sizeX = fRangeX.GetSize();
206 int xbin = -1;
207 if (vals.first < fRangeX.Get(0)) return -1;
208 for (int i = 1; i < sizeX; i++) {
209 if (vals.first < fRangeX.Get(i)) {
210 xbin = i - 1;
211 break;
212 }
213 }
214 int ybin = -1;
215 const Int_t sizeY = fRangeY.GetSize();
216 if (vals.second < fRangeY.Get(0)) return -1;
217 for (int i = 1; i < sizeY; i++) {
218 if (vals.second < fRangeY.Get(i)) {
219 ybin = i - 1;
220 break;
221 }
222 }
223 return xbin + ybin * (fRangeX.GetSize() - 1);
224 }
225
226 FemtoCorrFunc2D::FemtoCorrFunc2D(const FemtoCorrFunc2D& other) :
227 Hal::FemtoCorrFunc(other),
228 fLabelX(other.fLabelX),
229 fLabelY(other.fLabelY),
230 fRangeX(other.fRangeX),
231 fRangeY(other.fRangeY),
232 fNumProp(nullptr),
233 fDenProp(nullptr) {
234 if (other.fNumProp) {
235 fNumProp = (TH2D*) other.fNumProp;
236 fDenProp = (TH2D*) other.fDenProp;
237 }
238 }
239
240 Bool_t FemtoCorrFunc2D::InitPropMon() {
241 Double_t mins[2] = {fRangeX[0], fRangeY[0]};
242 Double_t maxs[2] = {fRangeX[fRangeX.GetSize() - 1], fRangeY[fRangeY.GetSize() - 1]};
243 fNumProp = new TH2D("NumProp", "NumProp", fBinsX, mins[0], maxs[0], fBinsY, mins[1], maxs[1]);
244 fDenProp = new TH2D("DenProp", "DenProp", fBinsX, mins[0], maxs[0], fBinsY, mins[1], maxs[1]);
245 fNumProp->SetXTitle(fLabelX);
246 fDenProp->SetYTitle(fLabelY);
247 return kTRUE;
248 }
249
250} /* namespace Hal */
T Get(Int_t i) const
Definition Array.h:97
Int_t GetSize() const
Definition Array.h:50
virtual TString HTMLExtract(Int_t counter=0, TString dir=" ") const
void FillNum(FemtoPair *pair)
void SetClass(TString className)
Definition HtmlObject.h:40
virtual void AddContent(const HtmlObject &obj)
Definition HtmlTable.cxx:29
void SetDefaultClass()
Definition HtmlTable.h:40
virtual void AddContent(const HtmlObject &obj)
Definition HtmlTable.cxx:17
TObject * At(Int_t i) const