forked from Shyam-Uniba/Shyam-Uniba
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathResolution_Tracker.C
295 lines (242 loc) · 11.4 KB
/
Resolution_Tracker.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
// Momentum resolution for charged particles
// Shyam Kumar
// Code based on https://arxiv.org/abs/1805.12014 formula's by Werner Reigler
#include<TMath.h>
#include <TString.h>
void MultipleScatteringComponent(Double_t mp, Double_t momI, Double_t &Ptresol, Double_t effradlen, Int_t N);
Double_t mpi = 139.57018; // mass of charged pions
Double_t speed= 3.0e8; //speed of light
Int_t charge=1;
Double_t min_effradlen=50.e-4/9.37; // Silicon of thickness 50mum
Double_t N = 6; // 6 silicon layers from 2,3,4,5,6 and 7 cm
Double_t resolution_rphi = 10*1.0e-6/sqrt(12); // pixel resolution_rphi in meter
Double_t resolution_z = 10*1.0e-6/sqrt(12); // pixel resolution_z in meter
Double_t min_length = 0.05; // Radius of outer Barrel layer-first layer
Double_t magfield = 3.0; // Mag field Eic 3.0 T
Double_t eta = 0.01;
Double_t r0 = 0.02; // 2cm radius of first layer
Double_t N_factor_Pt_Arc = sqrt(720*pow(N,3)/((N-1)*(N+1)*(N+2)*(N+3)));
Double_t N_factor_Pt_MS = N/sqrt(pow(N,2)-1);
Double_t N_2 = pow(N,2); Double_t N_3 = pow(N,3);
Double_t r0_l0 = r0/min_length; Double_t r0_l0_2 = pow(r0_l0,2); Double_t r0_l0_3 = pow(r0_l0,3); Double_t r0_l0_4 = pow(r0_l0,4);
Double_t N_factor_Pt_TransPA_1 = 1/sqrt((N-1)*(N+1)*(N+2)*(N+3));
Double_t N_factor_Pt_TransPA_2 = sqrt((N_3-N/3.0-2.0/3)+(4*(2*N_3-N_2-N)*r0_l0)+(4*(7*N_3-N_2-N)*r0_l0_2)+40*N_3*r0_l0_3+20*N_3*r0_l0_4);
Double_t N_factor_Pt_TransPA_SR = N_factor_Pt_TransPA_1*N_factor_Pt_TransPA_2; // Spatial resolution_rphi
Double_t N_factor_Pt_TransPA_MS = sqrt((N-3./4)/(N-1)+N/(2.0*(N-1))*r0_l0+N_2/(4.0*(N-1))*r0_l0_2);
Double_t N_factor_Pt_LongPA_1 = 1./sqrt((N+1)*(N+2));
Double_t N_factor_Pt_LongPA_2 = sqrt((N+1./2)+3.0*N*r0_l0+3.0*N*r0_l0_2);
Double_t N_factor_Pt_LongPA_SR = N_factor_Pt_LongPA_1*N_factor_Pt_LongPA_2; // Spatial resolution_z
void Resolution_Tracker()
{
gStyle->SetTitleSize(0.04,"");
gStyle->SetTitleSize(0.04,"X");
gStyle->SetTitleSize(0.04,"Y");
gStyle->SetTitleOffset(1.05,"Y");
TCanvas *c1 = new TCanvas("c1", "c1",0,52,1400,1000);
c1->SetGridy();
c1->SetMargin(0.12, 0.01 ,0.10,0.07);
TCanvas *c2 = new TCanvas("c2", "c2",0,52,1400,1000);
c2->SetGridy();
c2->SetMargin(0.12, 0.01 ,0.10,0.07);
TCanvas *c3 = new TCanvas("c3", "c3",0,52,1400,1000);
c3->SetGridy();
c3->SetMargin(0.12, 0.01 ,0.10,0.07);
std::vector<Double_t> x,y,c,d,e,f;
std::vector<Double_t> pt_resol_arc,pt_resol_ms,Trans_PA_resol_arc,Trans_PA_resol_ms,Long_PA_resol_arc,Long_PA_resol_ms;
x.clear(); y.clear(); c.clear(); d.clear(); e.clear(); f.clear();
pt_resol_arc.clear(); pt_resol_ms.clear(); Trans_PA_resol_arc.clear(); Trans_PA_resol_ms.clear();
Long_PA_resol_arc.clear(); Long_PA_resol_ms.clear();
//--------Muon Bethe Bloch---------------------
for (Int_t pt=100; pt<=10000.;pt=pt+100)
{
x.push_back(pt*0.001);
Double_t p = pt*TMath::CosH(eta);
Double_t theta = 2.0*TMath::ATan(TMath::Exp(-1.0*eta));
Double_t sin_theta = fabs(TMath::Sin(theta));
Double_t length = fabs(min_length/sin_theta);
Double_t effradlen = fabs(min_effradlen/sin_theta);
//----Energy of incident particles----------
Double_t ms_el=0.,ms_mu=0., ms_pi=0., ms_k=0. , ms_p=0. ;
MultipleScatteringComponent(mpi,p,ms_pi,effradlen,N);
// pT resolution
Double_t pi_pT_resol_Arc = 100.*pt*0.001*resolution_rphi*N_factor_Pt_Arc/(0.3*magfield*length*length);
Double_t pi_pT_MS = 100.*N_factor_Pt_MS*ms_pi/(0.3*magfield*length);
// Transverse Pointing Angle Resolution
Double_t Trans_PA_resol_Arc = 3.0*resolution_rphi*N_factor_Pt_TransPA_SR;
Double_t Trans_PA_resol_MS = r0/(pt*0.001)*ms_pi*N_factor_Pt_TransPA_MS;
// Longitudinal Pointing Angle Resolution
Double_t Long_PA_resol_Arc = 2.0*resolution_z*N_factor_Pt_LongPA_SR;
Double_t Long_PA_resol_MS = r0/(pt*0.001*sin_theta)*ms_pi;
pt_resol_arc.push_back(pi_pT_resol_Arc);
pt_resol_ms.push_back(pi_pT_MS);
Trans_PA_resol_arc.push_back(Trans_PA_resol_Arc);
Trans_PA_resol_ms.push_back(Trans_PA_resol_MS);
Long_PA_resol_arc.push_back(Long_PA_resol_Arc);
Long_PA_resol_ms.push_back(Long_PA_resol_MS);
}
//----------Pion Pt Resolution----------------------------------
const Int_t n=x.size();
Double_t a[n], b[n], g[n],h[n];
for(Int_t i=0;i<n;i++)
{
a[i]=x[i];
b[i]= pt_resol_ms[i];
g[i] = pt_resol_arc[i];
h[i] = TMath::Sqrt(b[i]*b[i]+g[i]*g[i]);
}
TGraph *gr1 = new TGraph(n,a,b);
gr1->SetLineColor(2);
gr1->GetYaxis()->SetRangeUser(0.0, 20.0);
gr1->SetLineWidth(1);
gr1->SetMarkerColor(2);
gr1->SetMarkerStyle(6);
gr1->SetTitle(Form("Transverse Momentum Resolution for #eta = %1.2f",eta));
gr1->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr1->GetXaxis()->CenterTitle(true);
gr1->GetYaxis()->SetTitle("#frac{#sigma_{p_{T}}}{p_{T}}"); // #sigma_{x} (cm)
gr1->GetYaxis()->CenterTitle(true);
c1->cd();
gr1->Draw("ACP");
TGraph *gr2 = new TGraph(n,a,g);
gr2->SetLineWidth(1);
gr2->SetMarkerColor(kMagenta);
gr2->SetLineColor(kMagenta);
gr2->SetMarkerStyle(6);
gr2->Draw("same");
TGraph *gr1sumpi = new TGraph(n,a,h);
gr1sumpi->SetLineWidth(1);
gr1sumpi->SetMarkerColor(kBlack);
gr1sumpi->SetLineColor(kBlack);
gr1sumpi->SetMarkerStyle(6);
gr1sumpi-> SetLineStyle(9);
gr1sumpi->Draw("same");
//--------------Legend Draw----------------------
TLegend *leg_hist = new TLegend(0.65,0.7,0.99,0.93);
leg_hist->SetHeader("Particle in Silicon Tracker");
leg_hist->SetTextFont(42);
leg_hist->SetTextSize(0.03);
leg_hist->AddEntry(gr1,"p_{T} Resol. ( M.S.)","l");
leg_hist->AddEntry(gr2,"p_{T} Resol. (Sagitta)","l");
leg_hist->AddEntry(gr1sumpi,"Sum p_{T} Resol. (Pion)","l");
leg_hist->Draw();
TGraph *gr_Pt_Resol = (TGraph*)gr1sumpi->Clone();
gr_Pt_Resol->SetName("Pt_Resol");
gr_Pt_Resol->SetTitle(Form("Transverse Momentum Resolution for #eta = %1.2f",eta));
gr_Pt_Resol->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr_Pt_Resol->GetXaxis()->CenterTitle(true);
gr_Pt_Resol->GetYaxis()->SetTitle("#frac{#sigma_{p_{T}}}{p_{T}}"); // #sigma_{x} (cm)
gr_Pt_Resol->GetYaxis()->CenterTitle(true);
//----------Pion Transverse Pointing Resolution----------------------------------
for(Int_t i=0;i<n;i++)
{
a[i]=x[i];
b[i]= Trans_PA_resol_ms[i]*1.0e+6;
g[i] = Trans_PA_resol_arc[i]*1.0e+6;
h[i] = TMath::Sqrt(b[i]*b[i]+g[i]*g[i]);
}
gr1 = new TGraph(n,a,b);
gr1->SetLineColor(2);
gr1->GetYaxis()->SetRangeUser(0.0, 100.0);
gr1->SetLineWidth(1);
gr1->SetMarkerColor(2);
gr1->SetMarkerStyle(6);
gr1->SetTitle(Form("Transverse Pointing Resolution for #eta = %1.2f",eta));
gr1->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr1->GetXaxis()->CenterTitle(true);
gr1->GetYaxis()->SetTitle("Transverse Pointing Resolution (#mum)"); // #sigma_{x} (cm)
gr1->GetYaxis()->CenterTitle(true);
c2->cd();
gr1->Draw("ACP");
gr2 = new TGraph(n,a,g);
gr2->SetLineWidth(1);
gr2->SetMarkerColor(kMagenta);
gr2->SetLineColor(kMagenta);
gr2->SetMarkerStyle(6);
gr2->Draw("same");
gr1sumpi = new TGraph(n,a,h);
gr1sumpi->SetLineWidth(1);
gr1sumpi->SetMarkerColor(kBlack);
gr1sumpi->SetLineColor(kBlack);
gr1sumpi->SetMarkerStyle(6);
gr1sumpi-> SetLineStyle(9);
gr1sumpi->Draw("same");
TGraph *gr_TPA_Resol = (TGraph*)gr1sumpi->Clone();
gr_TPA_Resol->SetName("Trans_PA");
gr_TPA_Resol->SetTitle(Form("Transverse Pointing Resolution for #eta = %1.2f",eta));
gr_TPA_Resol->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr_TPA_Resol->GetXaxis()->CenterTitle(true);
gr_TPA_Resol->GetYaxis()->SetTitle("Transverse Pointing Resolution (#mum)"); // #sigma_{x} (cm)
gr_TPA_Resol->GetYaxis()->CenterTitle(true);
//--------------Legend Draw----------------------
leg_hist = new TLegend(0.60,0.7,0.99,0.93);
leg_hist->SetHeader("Particle in Silicon Tracker");
leg_hist->SetTextFont(42);
leg_hist->SetTextSize(0.03);
leg_hist->AddEntry(gr1,"Transverse PA Resol. ( M.S.)","l");
leg_hist->AddEntry(gr2,"Transverse PA Resol. (Sagitta)","l");
leg_hist->AddEntry(gr1sumpi,"Sum Transverse PA Resol. (Pion)","l");
leg_hist->Draw();
//----------Longitudinal Pointing Resolution----------------------------------
for(Int_t i=0;i<n;i++)
{
a[i]=x[i];
b[i]= Long_PA_resol_ms[i]*1.0e+6;
g[i] = Long_PA_resol_arc[i]*1.0e+6;
h[i] = TMath::Sqrt(b[i]*b[i]+g[i]*g[i]);
}
gr1 = new TGraph(n,a,b);
gr1->SetLineColor(2);
gr1->GetYaxis()->SetRangeUser(0.0, 100.0);
gr1->SetLineWidth(1);
gr1->SetMarkerColor(2);
gr1->SetMarkerStyle(6);
gr1->SetTitle(Form("Longitudinal Pointing Resolution for #eta = %1.2f",eta));
gr1->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr1->GetXaxis()->CenterTitle(true);
gr1->GetYaxis()->SetTitle("Longitudinal Pointing Resolution (#mum)"); // #sigma_{x} (cm)
gr1->GetYaxis()->CenterTitle(true);
c3->cd();
gr1->Draw("ACP");
gr2 = new TGraph(n,a,g);
gr2->SetLineWidth(1);
gr2->SetMarkerColor(kMagenta);
gr2->SetLineColor(kMagenta);
gr2->SetMarkerStyle(6);
gr2->Draw("same");
gr1sumpi = new TGraph(n,a,h);
gr1sumpi->SetLineWidth(1);
gr1sumpi->SetMarkerColor(kBlack);
gr1sumpi->SetLineColor(kBlack);
gr1sumpi->SetMarkerStyle(6);
gr1sumpi-> SetLineStyle(9);
gr1sumpi->Draw("same");
TGraph *gr_LPA_Resol = (TGraph*)gr1sumpi->Clone();
gr_LPA_Resol->SetName("Long_PA");
gr_LPA_Resol->SetTitle(Form("Longitudinal Pointing Resolution for #eta = %1.2f",eta));
gr_LPA_Resol->GetXaxis()->SetTitle("p_{T} (GeV/c)");
gr_LPA_Resol->GetXaxis()->CenterTitle(true);
gr_LPA_Resol->GetYaxis()->SetTitle("Longitudinal Pointing Resolution (#mum)"); // #sigma_{x} (cm)
gr_LPA_Resol->GetYaxis()->CenterTitle(true);
//--------------Legend Draw----------------------
leg_hist = new TLegend(0.60,0.7,0.99,0.93);
leg_hist->SetHeader("Particle in Silicon Tracker");
leg_hist->SetTextFont(42);
leg_hist->SetTextSize(0.03);
leg_hist->AddEntry(gr1,"Longitudinal PA Resol. ( M.S.)","l");
leg_hist->AddEntry(gr2,"Longitudinal PA Resol. (Sagitta)","l");
leg_hist->AddEntry(gr1sumpi,"Sum Longitudinal PA Resol. (Pion)","l");
leg_hist->Draw();
TFile *fout = new TFile(Form("TrackerResolution_eta_%1.2f.root",eta),"recreate");
fout->cd();
gr_Pt_Resol->Write();
gr_TPA_Resol->Write();
gr_LPA_Resol->Write();
fout->Close();
}
// This function is basically 1/beta*f(d/X0Sintheta)
void MultipleScatteringComponent(Double_t mp, Double_t momI, Double_t &Ptresol, Double_t effradlen, Int_t N){
Double_t En = TMath::Sqrt(mp*mp+momI*momI);
Double_t beta = momI/En;
Ptresol=(0.0136/beta)*(TMath::Sqrt(effradlen))*(1+0.038*TMath::Log(effradlen)); //*(1+0.038*TMath::Log(effradlen))
}