Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
tprofile2polyRealistic.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_hist
3/// \notebook
4/// Different charges depending on region
5///
6/// \macro_image
7/// \macro_code
8///
9/// \author Filip Ilic
10
11#include <iostream>
12#include <fstream>
13#include <TProfile2Poly.h>
14#include <TProfile2D.h>
15#include <TCanvas.h>
16#include <TRandom.h>
17#include <TROOT.h>
18
19void tprofile2polyRealistic(Int_t numEvents=100000)
20{
21 int NUM_LS = 8;
22 TCanvas *c1 = new TCanvas("c1", "moving charge", 900, 400);
23 TCanvas *c2 = new TCanvas("c2", "Merge Individual moving charge plots", 800, 400);
24
25 c1->Divide(NUM_LS, 3);
26 c2->Divide(3,1);
27
28 // -------------------- Construct Reference plot bins ------------------------
29 auto new_avg = new TProfile2Poly();
30
31 auto tot_avg_ls = new TProfile2Poly[NUM_LS];
32 auto det_avg_ls = new TProfile2Poly[NUM_LS];
33 auto det_err_ls = new TProfile2Poly[NUM_LS];
34 auto tot_merge = new TProfile2Poly();
35 auto det_avg_merge = new TProfile2Poly();
36 auto det_err_merge = new TProfile2Poly();
37
38 float minx = -15;
39 float maxx = 15;
40 float miny = -15;
41 float maxy = 15;
42 float binsz = 0.5;
43
44 for (float i = minx; i < maxx; i += binsz) {
45 for (float j = miny; j < maxy; j += binsz) {
46 tot_merge->AddBin(i, j, i + binsz, j + binsz);
47 for (int l = 0; l < NUM_LS; ++l) {
48 tot_avg_ls[l].AddBin(i, j, i + binsz, j + binsz);
49 }
50 }
51 }
52
53 // -------------------- Construct detector bins ------------------------
54 auto h2p = new TH2Poly();
55 auto tp2p = new TProfile2Poly();
56 std::ifstream infile;
57 TString dir = gROOT->GetTutorialDir();
58 dir.Append("/hist/data/tprofile2poly_tutorial.data");
59 infile.open(dir.Data());
60
61 if (!infile) // Verify that the file was open successfully
62 {
63 std::cerr << dir.Data() << std::endl; // Report error
64 std::cerr << "Error code: " << std::strerror(errno) << std::endl; // Get some info as to why
65 return;
66 }
67 std::cout << " WE ARE AFTER LOADING DATA " << std::endl;
68
69 std::vector<std::pair<Double_t, Double_t>> allCoords;
70 Double_t a, b;
71 while (infile >> a >> b)
72 allCoords.emplace_back(a, b);
73
74 if (allCoords.size() % 3 != 0) {
75 cout << "[ERROR] Bad file" << endl;
76 return;
77 }
78
79 Double_t x[3], y[3];
80 for (int i = 0; i < allCoords.size(); i += 3) {
81 x[0] = allCoords[i + 0].first;
82 y[0] = allCoords[i + 0].second;
83 x[1] = allCoords[i + 1].first;
84 y[1] = allCoords[i + 1].second;
85 x[2] = allCoords[i + 2].first;
86 y[2] = allCoords[i + 2].second;
87
88 det_avg_merge->AddBin(3, x, y);
89 det_err_merge->AddBin(3, x, y);
90
91 for (int l = 0; l < NUM_LS; ++l) {
92 det_avg_ls[l].AddBin(3, x, y);
93 det_err_ls[l].AddBin(3, x, y);
94 }
95 }
96
97 std::cout << " WE ARE AFTER ADDING BINS " << std::endl;
98
99 // -------------------- Simulate particles ------------------------
100 TRandom ran;
101
102 // moving error
103 Double_t xoffset1 = 0;
104 Double_t yoffset1 = 0;
105 Double_t xoffset2 = 0;
106 Double_t yoffset2 = 0;
107
108 for (int i = 0; i <= NUM_LS-1; ++i) { // LumiSection
109 std::cout << "[In Progress] LumiSection " << i << std::endl;
110 for (int j = 0; j < numEvents; ++j) { // Events
111 Double_t r1 = ran.Gaus(0, 10);
112 Double_t r2 = ran.Gaus(0, 8);
113 Double_t rok = ran.Gaus(10, 1);
114 Double_t rbad1 = ran.Gaus(8, 5);
115 Double_t rbad2 = ran.Gaus(-8, 5);
116
117 Double_t val = rok;
118
119 xoffset1 += 0.00002;
120 yoffset1 += 0.00002;
121
122 xoffset2 += 0.00003;
123 yoffset2 += 0.00004;
124
125 if (r2 > 3. - yoffset1 && r2 < 8. - yoffset1 &&
126 r1 > 1. + xoffset1 && r1 < 5. + xoffset1 ) {
127 val -= rbad1;
128 }
129
130 if (r2 > -10 + yoffset2 && r2 < -8 + yoffset2 &&
131 r1 > -6 + xoffset2 && r1 < 8 + xoffset2 ) {
132 val -= rbad2;
133 }
134
135 tot_avg_ls[i].Fill(r1, r2, val);
136 det_avg_ls[i].Fill(r1, r2, val);
137 det_err_ls[i].Fill(r1, r2, val);
138 }
139
140 std::string title;
141
142 c1->cd(i+1);
143 title = "Global View: Avg in LS " + std::to_string(i);
144 tot_avg_ls[i].SetTitle(title.c_str());
145 tot_avg_ls[i].SetStats(false);
146 tot_avg_ls[i].Draw("COLZ");
147 c1->Update();
148
149 c1->cd((i+1)+NUM_LS);
150 title = "Detector View: Avg in LS " + std::to_string(i);
151 det_avg_ls[i].SetTitle(title.c_str());
152 det_avg_ls[i].SetStats(false);
153 det_avg_ls[i].Draw("COLZ");
154 c1->Update();
155
156 c1->cd((i+1)+(NUM_LS*2));
157 title = "Detector View: Error in LS " + std::to_string(i);
158 det_err_ls[i].SetTitle(title.c_str());
159 det_err_ls[i].SetStats(false);
160 det_err_ls[i].SetContentToError();
161 det_err_ls[i].Draw("COLZ");
162 c1->Update();
163 }
164
165 std::vector<TProfile2Poly*> tot_avg_v;
166 std::vector<TProfile2Poly*> det_avg_v;
167 for (int t = 0; t < NUM_LS; t++){
168 tot_avg_v.push_back(&tot_avg_ls[t]);
169 det_avg_v.push_back(&det_avg_ls[t]);
170 }
171
172 std::cout << "[In Progress] Merging" << std::endl;
173
174 tot_merge->Merge(tot_avg_v);
175 c2->cd(1);
176 tot_merge->SetTitle("Total average merge");
177 tot_merge->Draw("COLZ");
178
179 det_avg_merge->Merge(det_avg_v);
180 c2->cd(2);
181 det_avg_merge->SetTitle("Detector average merge");
182 det_avg_merge->SetContentToAverage(); // implicit
183 det_avg_merge->Draw("COLZ");
184
185 det_err_merge->Merge(det_avg_v);
186 c2->cd(3);
187 det_err_merge->SetTitle("Detector error merge");
188 det_err_merge->SetContentToError();
189 det_err_merge->Draw("COLZ");
190}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
#define gROOT
Definition TROOT.h:407
The Canvas class.
Definition TCanvas.h:23
2D Histogram with Polygonal Bins
Definition TH2Poly.h:66
2D Profile Histogram with Polygonal Bins.
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
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...
Definition TRandom.cxx:274
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
TString & Append(const char *cs)
Definition TString.h:576
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
TLine l
Definition textangle.C:4