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