Newer
Older
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
////////////////////////////////////////
// Read reconstruction ROOT output file
// Plot variables
////////////////////////////////////////
#include "ROOT/RDataFrame.hxx"
#include <iostream>
#include "dd4pod/Geant4ParticleCollection.h"
#include "eicd/ClusterCollection.h"
#include "eicd/ClusterData.h"
#include "TCanvas.h"
#include "TStyle.h"
#include "TMath.h"
#include "TH1.h"
#include "TH1D.h"
// If root cannot find a dd4pod header make sure to add the include dir to ROOT_INCLUDE_PATH:
// export ROOT_INCLUDE_PATH=$HOME/include:$ROOT_INCLUDE_PATH
// export ROOT_INCLUDE_PATH=/home/jihee/stow/development/include:$ROOT_INCLUDE_PATH
using ROOT::RDataFrame;
using namespace ROOT::VecOps;
void emcal_electrons_analysis(const char* input_fname = "rec_electron_10kEvt.root")
{
// Setting for graphs
gROOT->SetStyle("Plain");
gStyle->SetOptFit(1);
gStyle->SetLineWidth(2);
gStyle->SetPadTickX(1);
gStyle->SetPadTickY(1);
gStyle->SetPadGridX(1);
gStyle->SetPadGridY(1);
gStyle->SetPadLeftMargin(0.14);
ROOT::EnableImplicitMT();
ROOT::RDataFrame d0("events", input_fname);
// Number of Clusters
auto ncluster = [] (const std::vector<eic::ClusterData>& evt) {return (int) evt.size(); };
// Cluster Energy [GeV]
auto clusterE = [] (const std::vector<eic::ClusterData>& evt) {
std::vector<float> result;
for (const auto& i: evt)
result.push_back(i.energy);
return result;
};
// Scattering Angle [degree]
auto theta = [] (const std::vector<eic::ClusterData>& evt) {
std::vector<float> result;
for(const auto& i: evt) {
auto r = TMath::Sqrt(i.position.x*i.position.x + i.position.y*i.position.y + i.position.z*i.position.z);
result.push_back(TMath::ACos(i.position.z/r)*TMath::RadToDeg());
}
return result;
};
// Pseudo-rapidity
auto eta = [] (const std::vector<eic::ClusterData>& evt) {
std::vector<float> result;
for(const auto& i: evt) {
auto r = TMath::Sqrt(i.position.x*i.position.x + i.position.y*i.position.y + i.position.z*i.position.z);
auto th = TMath::ACos(i.position.z/r)*TMath::RadToDeg();
result.push_back(-1.0 * TMath::Log(TMath::Tan((th*TMath::DegToRad())/2.0)));
}
return result;
};
// Mass [GeV]
auto mass2 = [](std::vector<dd4pod::Geant4ParticleData> const& input) {
std::vector<float> result;
result.push_back(input[2].mass*input[2].mass);
return result;
};
// Momentum [GeV]
auto p_rec = [](const std::vector<float>& energy_term, const std::vector<float>& mass_term) {
std::vector<float> result;
for (const auto& e : energy_term) {
for (const auto& m2 : mass_term) {
result.push_back(TMath::Sqrt((e*e) - m2));
}
}
return result;
};
// Thrown Momentum [GeV]
auto p_thr = [](std::vector<dd4pod::Geant4ParticleData> const& input) {
std::vector<float> result;
result.push_back(TMath::Sqrt(input[2].psx*input[2].psx + input[2].psy*input[2].psy + input[2].psz*input[2].psz));
return result;
};
// Thrown Energy [GeV]
auto E_thr = [](std::vector<dd4pod::Geant4ParticleData> const& input) {
std::vector<float> result;
result.push_back(TMath::Sqrt(input[2].psx*input[2].psx + input[2].psy*input[2].psy + input[2].psz*input[2].psz + input[2].mass*input[2].mass));
return result;
};
// Energy Resolution
auto E_res = [] (const std::vector<float>& Erec, const std::vector<float>& Ethr) {
std::vector<float> result;
for (const auto& E1 : Ethr) {
for (const auto& E2 : Erec) {
result.push_back((E2-E1)/E1);
}
}
return result;
};
// Momentum Ratio
auto p_ratio = [] (const std::vector<float>& Prec, const std::vector<float>& Pthr) {
std::vector<float> result;
for (const auto& P1 : Pthr) {
for (const auto& P2 : Prec) {
result.push_back(P2/P1);
}
}
return result;
};
// Define variables
auto d1 = d0.Define("ncluster", ncluster, {"EcalClusters"})
.Define("clusterE", clusterE, {"EcalClusters"})
.Define("theta", theta, {"EcalClusters"})
.Define("eta", eta, {"EcalClusters"})
.Define("mass2", mass2, {"mcparticles2"})
.Define("p_rec", p_rec, {"clusterE","mass2"})
.Define("p_thr", p_thr, {"mcparticles2"})
.Define("E_thr", E_thr, {"mcparticles2"})
.Define("E_res", E_res, {"clusterE","E_thr"})
.Define("p_ratio", p_ratio, {"p_rec","p_thr"})
;
// Define Histograms
auto hNCluster = d1.Histo1D({"hNCluster", "Number of Clusters; # of Clusters; Events", 5, -0.5, 4.5}, "ncluster");
auto hClusterE = d1.Histo1D({"hClusterE", "Cluster Energy; Cluster Energy [GeV]; Events", 100, -0.5, 30.5}, "clusterE");
auto hTheta = d1.Histo1D({"hTheta", "Scattering Angle; #theta [degree]; Events", 100, 130.0, 180.0}, "theta");
auto hEta = d1.Histo1D({"hEta", "Pseudo-rapidity; #eta; Events", 100, -5.0, 0.0}, "eta");
auto hPrec = d1.Histo1D({"hPrec", "Reconstructed Momentum; p_{rec}[GeV]; Events", 100, -0.5, 30.5}, "p_rec");
auto hPthr = d1.Histo1D({"hPthr", "Thrown Momentum; p_{thr}[GeV]; Events", 100, -0.5, 30.5}, "p_thr");
auto hEthr = d1.Histo1D({"hEthr", "Thrown Energy; E_{thr}[GeV]; Events", 100, -0.5, 30.5}, "E_thr");
// Select Events with One Cluster
auto d2 = d1.Filter("ncluster==1");
auto hClusterE1 = d2.Histo1D({"hClusterE1", "One Cluster Energy; Cluster Energy [GeV]; Events", 100, -0.5, 30.5}, "clusterE");
auto hEres = d2.Histo1D({"hEres", "Energy Resolution; #DeltaE/E; Events", 100, -1.0, 1.0}, "E_res");
auto hPratio = d2.Histo1D({"hPratio", "Momentum ratio; p_{rec}/p_{thr}; Events", 100, 0.0, 1.0}, "p_ratio");
auto hPthr_accepted = d2.Filter([=] (const std::vector<float>& Prec, const std::vector<float>& Pthr) {
for (const auto& P1 : Pthr) {
for (const auto& P2 : Prec) {
auto Pratio = P2/P1;
if (Pratio > 0.70) {
return true;
}
}
}
return false;}, {"p_rec","p_thr"})
.Histo1D({"hPthr_accepted", "Thrown momentum for reconstructed particle; p_{thr} [GeV]; Events", 100, -0.5, 30.5}, "p_thr");
// Cut on Radial Distance
auto d3 = d2.Filter([=] (const std::vector<eic::ClusterData>& evt) {
for(const auto& i: evt) {
auto pos_x = i.position.x;
auto pos_y = i.position.y;
auto radial_dist = TMath::Sqrt(pos_x*pos_x + pos_y*pos_y);
if (radial_dist > 18.0 && radial_dist < 30.0)
return true;
}
return false;}, {"EcalClusters"});
auto hEres_cut = d3.Histo1D({"hEres_cut", "Energy Resolution; #DeltaE/E; Events", 100, -1.0, 1.0}, "E_res");
auto hTheta_cut = d3.Histo1D({"hTheta_cut", "Scattering Angle; #theta [degree]; Events", 100, 130.0, 180.0}, "theta");
auto hEta_cut = d3.Histo1D({"hEta_cut", "Pseudo-rapidity; #eta; Events", 100, -5.0, 0.0}, "eta");
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
// Event Counts
auto nevents_thrown = d1.Count();
auto nevents_cluster1 = d2.Count();
std::cout << "Number of Events: " << (*nevents_cluster1) << " / " << (*nevents_thrown) << " = single cluster events/all \n";
// Draw Histograms
TCanvas *c1 = new TCanvas("c1", "c1", 500, 500);
c1->SetLogy(1);
hNCluster->GetYaxis()->SetTitleOffset(1.6);
hNCluster->SetLineWidth(2);
hNCluster->SetLineColor(kBlue);
hNCluster->DrawClone();
c1->SaveAs("results/emcal_electrons_ncluster.png");
c1->SaveAs("results/emcal_electrons_ncluster.pdf");
TCanvas *c2 = new TCanvas("c2", "c2", 500, 500);
hClusterE->GetYaxis()->SetTitleOffset(1.4);
hClusterE->SetLineWidth(2);
hClusterE->SetLineColor(kBlue);
hClusterE->DrawClone();
c2->SaveAs("results/emcal_electrons_clusterE.png");
c2->SaveAs("results/emcal_electrons_clusterE.pdf");
TCanvas *c3 = new TCanvas("c3", "c3", 500, 500);
hTheta->GetYaxis()->SetTitleOffset(1.4);
hTheta->SetLineWidth(2);
hTheta->SetLineColor(kBlue);
hTheta->DrawClone();
c3->SaveAs("results/emal_electrons_theta.png");
c3->SaveAs("results/emal_electrons_theta.pdf");
TCanvas *c4 = new TCanvas("c4", "c4", 500, 500);
hEta->GetYaxis()->SetTitleOffset(1.4);
hEta->SetLineWidth(2);
hEta->SetLineColor(kBlue);
hEta->DrawClone();
c4->SaveAs("results/emcal_electrons_eta.png");
c4->SaveAs("results/emcal_electrons_eta.pdf");
TCanvas *c5 = new TCanvas("c5", "c5", 500, 500);
hPrec->GetYaxis()->SetTitleOffset(1.4);
hPrec->SetLineWidth(2);
hPrec->SetLineColor(kBlue);
hPrec->DrawClone();
c5->SaveAs("results/emcal_electrons_Prec.png");
c5->SaveAs("results/emcal_electrons_Prec.pdf");
TCanvas *c6 = new TCanvas("c6", "c6", 500, 500);
hPthr->GetYaxis()->SetTitleOffset(1.4);
hPthr->SetLineWidth(2);
hPthr->SetLineColor(kBlue);
hPthr->DrawClone();
c6->SaveAs("results/emcal_electrons_Pthr.png");
c6->SaveAs("results/emcal_electrons_Pthr.pdf");
TCanvas *c7 = new TCanvas("c7", "c7", 500, 500);
hEthr->GetYaxis()->SetTitleOffset(1.4);
hEthr->SetLineWidth(2);
hEthr->SetLineColor(kBlue);
hEthr->DrawClone();
c7->SaveAs("results/emcal_electrons_Ethr.png");
c7->SaveAs("results/emcal_electrons_Ethr.pdf");
TCanvas *c8 = new TCanvas("c8", "c8", 500, 500);
hClusterE->GetYaxis()->SetTitleOffset(1.4);
hClusterE->SetLineWidth(2);
hClusterE->SetLineColor(kBlue);
hClusterE->DrawClone();
c8->SaveAs("results/emcal_electrons_clusterE1.png");
c8->SaveAs("results/emcal_electrons_clusterE1.pdf");
TCanvas *c9 = new TCanvas("c9", "c9", 500, 500);
hEres->GetYaxis()->SetTitleOffset(1.4);
hEres->SetLineWidth(2);
hEres->SetLineColor(kBlue);
hEres->Fit("gaus");
hEres->DrawClone();
c9->SaveAs("results/emcal_electrons_Eres.png");
c9->SaveAs("results/emcal_electrons_Eres.pdf");
TCanvas *c10 = new TCanvas("c10", "c10", 500, 500);
hPthr_accepted->GetYaxis()->SetTitleOffset(1.4);
hPthr_accepted->SetLineWidth(2);
hPthr_accepted->SetLineColor(kBlue);
hPthr_accepted->DrawClone();
c10->SaveAs("results/emcal_electrons_Pthr_accepted.png");
c10->SaveAs("results/emcal_electrons_Pthr_accepted.pdf");
TCanvas *c11 = new TCanvas("c11", "c11", 500, 500);
hPratio->GetYaxis()->SetTitleOffset(1.4);
hPratio->SetLineWidth(2);
hPratio->SetLineColor(kBlue);
hPratio->DrawClone();
c11->SaveAs("results/emcal_electrons_Pratio.png");
c11->SaveAs("results/emcal_electrons_Pratio.pdf");
TCanvas *c12 = new TCanvas("c12", "c12", 500, 500);
TH1D *hPacceptance = new TH1D("hPacceptance","Ratio p_{rec}/p_{thr}", 100, -0.5, 30.5);
hPacceptance = (TH1D*) hPthr_accepted->Clone();
hPacceptance->Divide(hPthr.GetPtr());
hPacceptance->SetTitle("Ratio p_{rec}/p_{thr}");
hPacceptance->SetLineColor(kBlue);
hPacceptance->SetLineWidth(2);
hPacceptance->GetXaxis()->SetTitle("p [GeV]");
hPacceptance->GetYaxis()->SetTitle("p_{rec}/p_{thr}");
hPacceptance->GetYaxis()->SetTitleOffset(1.4);
hPacceptance->DrawClone();
c12->SaveAs("results/emcal_electrons_Pacceptance.png");
c12->SaveAs("results/emcal_electrons_Pacceptance.pdf");
TCanvas *c13 = new TCanvas("c13", "c13", 500, 500);
hEres_cut->GetYaxis()->SetTitleOffset(1.4);
hEres_cut->SetLineWidth(2);
hEres_cut->SetLineColor(kBlue);
hEres_cut->Fit("gaus");
hEres_cut->DrawClone();
c13->SaveAs("results/emcal_electrons_Eres_cut.png");
c13->SaveAs("results/emcal_electrons_Eres_cut.pdf");
TCanvas *c14 = new TCanvas("c14", "c14", 500, 500);
hTheta_cut->GetYaxis()->SetTitleOffset(1.4);
hTheta_cut->SetLineWidth(2);
hTheta_cut->SetLineColor(kBlue);
hTheta_cut->DrawClone();
c14->SaveAs("results/emal_electrons_theta_cut.png");
c14->SaveAs("results/emal_electrons_theta_cut.pdf");
TCanvas *c15 = new TCanvas("c15", "c15", 500, 500);
hEta_cut->GetYaxis()->SetTitleOffset(1.4);
hEta_cut->SetLineWidth(2);
hEta_cut->SetLineColor(kBlue);
hEta_cut->DrawClone();
c15->SaveAs("results/emcal_electrons_eta_cut.png");
c15->SaveAs("results/emcal_electrons_eta_cut.pdf");