]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ANALYSIS/AliD0toKpiPlots.C
Coding conventions
[u/mrichter/AliRoot.git] / ANALYSIS / AliD0toKpiPlots.C
CommitLineData
3a9a3487 1void AliD0toKpiPlots(const Char_t *inName="AliD0toKpi.root",
2 const Char_t *outName="D0histograms.root") {
3 //--------------------------------------------------------------------------
4 // This macro histograms many variables of D0->Kpi candidates
5 //
6 // Andrea Dainese, andrea.dainese@pd.infn.it
7 //--------------------------------------------------------------------------
8
9 // set of cuts
10 Double_t D0Cuts[9] = {0.1, // mass [GeV]
11 1000000., // dca [micron]
12 1.1, // cosThetaStar
13 0., // pT K [GeV/c]
14 0., // pT Pi [GeV/c]
15 100000., // d0K upper [micron]
16 100000., // d0Pi upper [micron]
17 10000000000., // d0d0 [micron^2]
18 -1.1}; // cosThetaPointing
19
20 // number of events (for normalization)
21 Bool_t normalize = kFALSE;
22 Double_t events = 1.;
23
24
25
26 // define histograms
27 TH1F *hptK = new TH1F("hptK","\"K\" p_{t} distribution",50,0,10);
28 hptK->SetXTitle("p_{t} [GeV]");
29
30 TH1F *hptPi = new TH1F("hptPi","\"#pi\" p_{t} distribution",50,0,10);
31 hptPi->SetXTitle("p_{t} [GeV]");
32
33 TH1F *hDCA = new TH1F("hDCA","DCA",50,0,1000);
34 hDCA->SetXTitle("dca [#mu m]");
35
36 TH1F *hptD0 = new TH1F("hptD0","D^{0} p_{t} distribution",40,0,40);
37 hptD0->SetXTitle("p_{t} [GeV]");
38
39 TH1F *hyD0 = new TH1F("hyD0","D^{0} rapidity distribution",50,-2,2);
40 hyD0->SetXTitle("y");
41
42 TH1F *hCPtaD0 = new TH1F("hCPtaD0","cosine of pointing angle distribution",100,-1,1);
43 hCPtaD0->SetXTitle("cos #theta_{point}");
44
45 TH1F *hCPtaXY = new TH1F("hCPtaXY","cosine of pointing angle in (x,y) plane",100,-1,1);
46 hCPtaXY->SetXTitle("cos #theta_{point}");
47
48 TH1F *hCts = new TH1F("hCts","cosine of decay angle",50,-1.2,1.2);
49 hCts->SetXTitle("cos #theta^{*}");
50
51 TH2F *hCtsVsPtK = new TH2F("hCtsVsPtK","cosine of decay angle VS \"K\" p_{t}",50,0,5,50,-1,1);
52 hCtsVsPtK->SetYTitle("cos #theta^{*}");
53 hCtsVsPtK->SetXTitle("p_{t} [GeV]");
54
55 TH1F *hd0d0 = new TH1F("hd0d0","Product of the impact parameters",100,-100000,100000);
56 hd0d0->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
57
58 TH1F *hd0K = new TH1F("hd0K","Impact parameter of \"K\"",100,-5000,5000);
59 hd0K->SetXTitle("d_{0}^{K} [#mu m]");
60
61 TH1F *hd0Pi = new TH1F("hd0Pi","Impact parameter of \"#pi\"",100,-5000,5000);
62 hd0Pi->SetXTitle("d_{0}^{#pi} [#mu m]");
63
64 TH2F *hCPtaVsd0d0 = new TH2F("hCPtaVsd0d0","cos #theta_{point} vs d_{0}^{K} #times d_{0}^{#pi}",100,-100000,100000,100,-1,1);
65 hCPtaVsd0d0->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
66 hCPtaVsd0d0->SetYTitle("cos #theta_{point}");
67
68 TH2F *hCPtaVsd0d0zoom = new TH2F("hCPtaVsd0d0zoom","cos #theta_{point} vs d_{0}^{K} #times d_{0}^{#pi}",100,-100000,0,100,.9,1);
69 hCPtaVsd0d0zoom->SetXTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
70 hCPtaVsd0d0zoom->SetYTitle("cos #theta_{point}");
71
72 TH2F *hd0d0VSptD0 = new TH2F("hd0d0VSptD0","d_{0}^{K} #times d_{0}^{#pi} VS D^{0} p_{t}",50,0,25,100,-120000,120000);
73 hd0d0VSptD0->SetYTitle("d_{0}^{K} #times d_{0}^{#pi} [#mu m^{2}]");
74 hd0d0VSptD0->SetXTitle("D^{0} p_{t} [GeV]");
75
76 TH1F *hMass = new TH1F("hMass","Invariant mass distribution",50,1.765,1.965);
77 hMass->SetXTitle("M[K,#pi] [GeV]");
78
79 TH2F *hArm = new TH2F("hArm","Armenteros plot",50,-2,2,50,0,1);
80 hArm->SetXTitle("#alpha");
81 hArm->SetYTitle("q_{t}");
82
83 // open input file and get tree
84 TFile *inFile = TFile::Open(inName);
85
86 TTree *treeD0 = (TTree*)inFile->Get("TreeD0");
87 AliD0toKpi *D = 0;
88 treeD0->SetBranchAddress("D0toKpi",&D);
89 Int_t entries = (Int_t)treeD0->GetEntries();
90
91 printf("+++\n+++ Number of D0 in tree: %d\n+++\n",entries);
92
93 Double_t MD0,MD0bar,ctsD0,ctsD0bar,ctsPiD0,ctsPiD0bar;
94 Double_t WgtD0,WgtD0bar;
95 Double_t sampleABC=0.;
96 Int_t okD0=0,okD0bar=0;
97 Int_t nSel = 0;
98 Int_t ptbin;
99
100 // loop on D0
101 for(Int_t i=0; i<entries; i++) {
102 if(i%10000==0) printf(" candidate %d of %d\n",i,entries);
103
104 // get event from tree
105 treeD0->GetEvent(i);
106 //--- select the PID strategy & compute weights
107 // D->ApplyPID("TOFparam_PbPb");
108 // D->ComputeWgts();
109 // get weights for the three samples A+B+C
110 // D->GetWgts(WgtD0,WgtD0bar,"ABC");
111 WgtD0 = 1.; WgtD0bar = 1.;
112
113 // normalize to 1 event
114 if(normalize) { WgtD0 /= events; WgtD0bar /= events; }
115
116 // check if candidate passes selection (as D0 or D0bar)
117 D->Select(D0Cuts,okD0,okD0bar);
118
119 // set weights to 0 if the candidate doesn't pass selection
120 if(!okD0) WgtD0=0.;
121 if(!okD0bar) WgtD0bar=0.;
122 if(okD0 || okD0bar) nSel++;
123
124 // count selected candidates
125 sampleABC += WgtD0 + WgtD0bar;
126
127 // inv mass and cosThetaStar
128 D->InvMass(MD0,MD0bar);
129 D->CosThetaStar(ctsD0,ctsD0bar);
130
131 // fill histograms
132 hptK->Fill(D->PtChild(1),WgtD0);
133 hptK->Fill(D->PtChild(0),WgtD0bar);
134 hptPi->Fill(D->PtChild(0),WgtD0);
135 hptPi->Fill(D->PtChild(1),WgtD0bar);
136 hd0K->Fill(D->Getd0Child(1),WgtD0);
137 hd0K->Fill(D->Getd0Child(0),WgtD0bar);
138 hd0Pi->Fill(D->Getd0Child(0),WgtD0);
139 hd0Pi->Fill(D->Getd0Child(1),WgtD0bar);
140 hMass->Fill(MD0,WgtD0);
141 hMass->Fill(MD0bar,WgtD0bar);
142 hCts->Fill(ctsD0,WgtD0);
143 hCts->Fill(ctsD0bar,WgtD0bar);
144 hCtsVsPtK->Fill(D->PtChild(1),ctsD0,WgtD0);
145 hCtsVsPtK->Fill(D->PtChild(0),ctsD0bar,WgtD0bar);
146 hDCA->Fill(D->GetDCA(),WgtD0+WgtD0bar);
147 hptD0->Fill(D->Pt(),WgtD0+WgtD0bar);
148 hyD0->Fill(D->Rapidity(),WgtD0+WgtD0bar);
149 hd0d0->Fill(D->ProdImpParams(),WgtD0+WgtD0bar);
150 hCPtaD0->Fill(D->CosPointing(),WgtD0+WgtD0bar);
151 hCPtaXY->Fill(D->CosPointingXY(),WgtD0+WgtD0bar);
152 hCPtaVsd0d0->Fill(D->ProdImpParams(),D->CosPointing(),WgtD0+WgtD0bar);
153 hd0d0VSptD0->Fill(D->Pt(),D->ProdImpParams(),WgtD0+WgtD0bar);
154 hCPtaVsd0d0zoom->Fill(D->ProdImpParams(),D->CosPointing(),WgtD0+WgtD0bar);
155 hArm->Fill(D->Alpha(),D->Qt(),WgtD0+WgtD0bar);
156
157
158 } // end loop on D0 candidates
159
160 inFile->Close();
161
162 printf("\n\n --- Total number of candidates passing selection: %d\n\n --- Sum of weights sample A+B+C: %f\n\n",nSel,sampleABC);
163
164 // draw histograms
165 TCanvas *c1 = new TCanvas("c1","pt K & pi",0,0,700,700);
166 c1->SetLogy();
167 hptK->Draw();
168 hptPi->Draw("same");
169
170 TCanvas *c2 = new TCanvas("c2","pt D0",0,0,700,700);
171 c2->SetLogy();
172 hptD0->Draw();
173
174 TCanvas *c3 = new TCanvas("c3","rapidity D0",0,0,700,700);
175 hyD0->Draw();
176
177 TCanvas *c4 = new TCanvas("c4","pointing angle",0,0,700,700);
178 hCPtaD0->Draw();
179
180 TCanvas *c5 = new TCanvas("c5","d0 x d0",0,0,700,700);
181 c5->SetLogy();
182 hd0d0->Draw();
183
184 TCanvas *c6 = new TCanvas("c6","pointing angle VS d0d0",0,0,700,700);
185 c6->SetLogz();
186 hCPtaVsd0d0->Draw("box");
187
188 TCanvas *c7 = new TCanvas("c7","mass",0,0,700,700);
189 hMass->Draw();
190
191 TCanvas *c8 = new TCanvas("c8","armenteros",0,0,700,700);
192 hArm->Draw("box");
193
194 TCanvas *c9 = new TCanvas("c9","decay angle",0,0,700,700);
195 hCts->Draw();
196
197 TCanvas *c10 = new TCanvas("c10","dca",0,0,700,700);
198 c10->SetLogy();
199 hDCA->Draw();
200
201 TCanvas *c11 = new TCanvas("c11","d0 K & pi",0,0,700,700);
202 c11->SetLogy();
203 hd0K->Draw();
204 hd0Pi->Draw("same");
205
206 // write all histograms to file
207 TFile *outFile = new TFile(outName,"recreate");
208 hMass->Write();
209 hDCA->Write();
210 hCts->Write();
211 hCtsVsPtK->Write();
212 hArm->Write();
213 hCPtaVsd0d0->Write();
214 hd0d0VSptD0->Write();
215 hCPtaVsd0d0zoom->Write();
216 hd0d0->Write();
217 hCPtaD0->Write();
218 hCPtaXY->Write();
219 hptK->Write();
220 hptPi->Write();
221 hptD0->Write();
222 hyD0->Write();
223 hd0K->Write();
224 hd0Pi->Write();
225 outFile->Close();
226
227 return;
228}
229
230
231
232