]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG0/AliRecoVsGeneCheck.cxx
commented logging message
[u/mrichter/AliRoot.git] / PWG0 / AliRecoVsGeneCheck.cxx
CommitLineData
fe646d50 1#include "AliRecoVsGeneCheck.h"
2
3#include <Riostream.h>
4
5//____________________________________________________________________
6ClassImp(AliRecoVsGeneCheck);
7
8//____________________________________________________________________
9AliRecoVsGeneCheck::AliRecoVsGeneCheck() {
10
11 fhVtzZRecoVsMC = new TH2F("vtx_z_reco_vs_mc","",100,-20,20,100,-20,20);
12
13 fhVtxZRes = new TH1F("vtx_z_res","",100,0,0.05);
14 fhVtxZResVsZ = new TH2F("vtx_z_res_vs_z","",100,-20,20,100,0,0.05);
15 fhVtxZResVsNPart = new TH2F("vtx_z_res_vs_npart","",100,-0.5,99.5,100,0,0.05);
16
17 fhVtxDzNorm = new TH1F("vtx_dz_norm","",100,-10,10);
18 fhVtxDzNormVsZ = new TH2F("vtx_dz_norm_vs_z","",100,-20,20,100,-10,10);
19 fhVtxDzNormVsNPart = new TH2F("vtx_dz_norm_vs_npart","",100,-0.5,99.5,100,-10,10);
20
21 fhVtxZMC = new TH1F("vtx_z_mc","",100,-20,20);
22 fhVtxZReco = new TH1F("vtx_z_reco","",100,-20,20);
23
24 fhNPart = new TH1F("n_part","",100,-0.5,99.5);
25
26 fhDPtVsPtVsEta = new TH3F("dpt_vs_pt_vs_eta","", 30,-1.5,1.5,30,0,3,120,-0.3,0.3);
27 fhDEtaVsPtVsEta = new TH3F("deta_vs_pt_vs_eta","",30,-1.5,1.5,30,0,3,120,-0.3,0.3);
28
29 fhVtxZMC ->Sumw2();
30 fhVtxZReco->Sumw2();
31 fhNPart ->Sumw2();
32
33 fhVtzZRecoVsMC ->SetXTitle("vertex z_{mc}");
34 fhVtzZRecoVsMC ->SetYTitle("vertex z_{reco}");
35
36 fhVtxZRes ->SetXTitle("vertex #sigma_{z_{reco}}");
37 fhVtxZResVsZ ->SetXTitle("vertex z_{mc}");
38 fhVtxZResVsZ ->SetYTitle("vertex #sigma_{z_{reco}}");
39 fhVtxZResVsNPart ->SetXTitle("n part");
40 fhVtxZResVsNPart ->SetYTitle("vertex #sigma_{z_{reco}}");
41
42 fhVtxDzNorm ->SetXTitle("vertex (z_{mc} - z_{reco}) / #sigma_{z_{reco}}");
43 fhVtxDzNormVsZ ->SetXTitle("vertex z_{mc}");
44 fhVtxDzNormVsZ ->SetYTitle("vertex (z_{mc} - z_{reco}) / #sigma_{z_{reco}}");
45 fhVtxDzNormVsNPart->SetXTitle("npart");
46 fhVtxDzNormVsNPart->SetYTitle("vertex (z_{mc} - z_{reco}) / #sigma_{z_{reco}}");
47
48 fhNPart ->SetXTitle("n part");
49
50 fhVtxZMC ->SetXTitle("vertex z_{mc}");
51 fhVtxZReco ->SetXTitle("vertex z_{reco}");
52
53 fhDPtVsPtVsEta ->SetXTitle("#eta");
54 fhDPtVsPtVsEta ->SetYTitle("p_{T} [GeV/c]");
55 fhDPtVsPtVsEta ->SetZTitle("#deltap_{T} [GeV/c]");
56
57 fhDEtaVsPtVsEta ->SetXTitle("#eta");
58 fhDEtaVsPtVsEta ->SetYTitle("p_{T} [GeV/c]");
59 fhDEtaVsPtVsEta ->SetZTitle("#delta#eta");
60
61}
62
63//____________________________________________________________________
64void
65AliRecoVsGeneCheck::Event(Double_t* vtx, Double_t* vtx_res, Double_t* mcvtx, Int_t n_part) {
66
67 fhVtzZRecoVsMC ->Fill(mcvtx[2],vtx[2]);
68
69 fhVtxZRes ->Fill(vtx_res[2]);
70 fhVtxZResVsZ ->Fill(mcvtx[2], vtx_res[2]);
71 fhVtxZResVsNPart ->Fill(n_part, vtx_res[2]);
72
73 if (vtx_res[2]!=0) {
74 Float_t dzNorm = (mcvtx[2] - vtx[2])/vtx_res[2];
75
76 fhVtxDzNorm ->Fill(dzNorm);
77 fhVtxDzNormVsZ ->Fill(mcvtx[2], dzNorm);
78 fhVtxDzNormVsNPart->Fill(n_part, dzNorm);
79
80 }
81 fhVtxZMC ->Fill(mcvtx[2]);
82 fhVtxZReco ->Fill(vtx[2]);
83 fhNPart ->Fill(n_part);
84}
85
86
87//____________________________________________________________________
88void
89AliRecoVsGeneCheck::Track(AliESDtrack* esdTrack, TParticle* mcParticle) {
90
91 Double_t p[3];
92 esdTrack->GetConstrainedPxPyPz(p); // ### TODO or GetInnerPxPyPy / GetOuterPxPyPy
93 TVector3 vector(p);
94
95 Float_t reco_p = esdTrack->GetP();
96 Float_t reco_theta = vector.Theta();
97 // Float_t reco_phi =
98 Float_t reco_eta = -TMath::Log(TMath::Tan(reco_theta/2.));
99 Float_t reco_pt = reco_p*TMath::Sin(reco_theta);
100
101 Float_t mc_eta = mcParticle->Eta();
102 Float_t mc_pt = mcParticle->Pt();
103
104 fhDPtVsPtVsEta->Fill(mc_eta, mc_pt, mc_pt - reco_pt);
105
106 fhDEtaVsPtVsEta->Fill(mc_eta, mc_pt, mc_eta - reco_eta);
107
108}
109
110//____________________________________________________________________
111void
112AliRecoVsGeneCheck::SaveHistograms(Char_t* dir) {
113
114 gDirectory->mkdir(dir);
115 gDirectory->cd(dir);
116
117 gDirectory->mkdir("event");
118 gDirectory->cd("event");
119
120 fhVtzZRecoVsMC ->Write();
121 fhVtxZRes ->Write();
122 fhVtxZResVsZ ->Write();
123 fhVtxZResVsNPart ->Write();
124 fhVtxDzNorm ->Write();
125 fhVtxDzNormVsZ ->Write();
126 fhVtxDzNormVsNPart ->Write();
127 fhVtxZMC ->Write();
128 fhVtxZReco ->Write();
129 fhNPart ->Write();
130
131 gDirectory->cd("../");
132
133
134 gDirectory->mkdir("track");
135 gDirectory->cd("track");
136
137 fhDPtVsPtVsEta ->Write();
138 fhDEtaVsPtVsEta ->Write();
139
140 gDirectory->cd("../");
141
142 gDirectory->cd("../");
143}