20 void tprofile2polyRealistic(
Int_t numEvents=100000)
24 TCanvas *
c2 =
new TCanvas(
"c2",
"Merge Individual moving charge plots", 800, 400);
45 for (
float i = minx; i < maxx; i += binsz) {
46 for (
float j = miny; j < maxy; j += binsz) {
47 tot_merge->AddBin(i, j, i + binsz, j + binsz);
48 for (
int l=0;
l<NUM_LS; ++
l) {
49 tot_avg_ls[
l].AddBin(i, j, i + binsz, j + binsz);
58 TString dir =
gROOT->GetTutorialDir();
59 dir.Append(
"/hist/data/tprofile2poly_tutorial.data");
60 infile.open(dir.Data());
64 std::cerr << dir.Data() << std::endl;
65 std::cerr <<
"Error code: " << strerror(errno) << std::endl;
68 std::cout <<
" WE ARE AFTER LOADING DATA " << std::endl;
70 vector<pair<Double_t, Double_t>> allCoords;
72 while (infile >> a >> b) {
73 pair<Double_t, Double_t> coord(a, b);
74 allCoords.push_back(coord);
77 if (allCoords.size() % 3 != 0) {
78 cout <<
"[ERROR] Bad file" << endl;
83 for (
Int_t i = 0; i < allCoords.size(); i += 3) {
84 x[0] = allCoords[i + 0].first;
85 y[0] = allCoords[i + 0].second;
86 x[1] = allCoords[i + 1].first;
87 y[1] = allCoords[i + 1].second;
88 x[2] = allCoords[i + 2].first;
89 y[2] = allCoords[i + 2].second;
91 det_avg_merge->AddBin(3, x, y);
92 det_err_merge->AddBin(3, x, y);
94 for (
int l=0;
l<NUM_LS; ++
l) {
95 det_avg_ls[
l].AddBin(3, x, y);
99 std::cout <<
" WE ARE AFTER ADDING BINS " << std::endl;
110 for (
int i = 0; i <= NUM_LS-1; ++i) {
111 std::cout <<
"[In Progress] LumiSection " << i << std::endl;
112 for (
int j = 0; j < numEvents; ++j) {
127 if (r2 > 3. - yoffset1 && r2 < 8. - yoffset1 &&
128 r1 > 1. + xoffset1 && r1 < 5. + xoffset1 ) {
132 if (r2 > -10 + yoffset2 && r2 < -8 + yoffset2 &&
133 r1 > -6 + xoffset2 && r1 < 8 + xoffset2 ) {
137 tot_avg_ls[i].Fill(r1, r2, val);
138 det_avg_ls[i].Fill(r1, r2, val);
145 title =
"Global View: Avg in LS " + to_string(i);
146 tot_avg_ls[i].SetTitle(title.c_str());
147 tot_avg_ls[i].SetStats(
false);
148 tot_avg_ls[i].Draw(
"COLZ");
151 c1->
cd((i+1)+NUM_LS);
152 title =
"Detector View: Avg in LS " + to_string(i);
153 det_avg_ls[i].SetTitle(title.c_str());
154 det_avg_ls[i].SetStats(
false);
155 det_avg_ls[i].Draw(
"COLZ");
159 c1->
cd((i+1)+(NUM_LS*2));
160 title =
"Detector View: Error in LS " + to_string(i);
161 det_avg_ls[i].SetTitle(title.c_str());
162 det_avg_ls[i].SetStats(
false);
163 det_avg_ls[i].SetContentToError();
164 det_avg_ls[i].Draw(
"COLZ");
168 std::vector<TProfile2Poly*> tot_avg_v;
169 std::vector<TProfile2Poly*> det_avg_v;
170 for (
Int_t t=0; t<NUM_LS; t++){
171 tot_avg_v.push_back(&tot_avg_ls[t]);
172 det_avg_v.push_back(&det_avg_ls[t]);
175 std::cout <<
"[In Progress] Merging" << std::endl;
179 tot_merge->Merge(tot_avg_v);
181 title =
"Total average merge";
183 tot_merge->Draw(
"COLZ");
185 det_avg_merge->Merge(det_avg_v);
187 title =
"Detector average merge";
188 det_avg_merge->
SetTitle(title.c_str());
189 det_avg_merge->SetContentToAverage();
190 det_avg_merge->Draw(
"COLZ");
192 det_err_merge->Merge(det_avg_v);
194 title =
"Detector error merge";
195 det_err_merge->
SetTitle(title.c_str());
196 det_err_merge->SetContentToError();
197 det_err_merge->Draw(
"COLZ");
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
This is the base class for the ROOT Random number generators.
2D Profile Histogram with Polygonal Bins.
virtual void SetTitle(const char *title="")=0
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
virtual void Update()
Update canvas pad buffers.
2D Histogram with Polygonal Bins